Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4GeomTools.cc
이 파일의 문서화 페이지로 가기
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: $
28 //
29 //
30 // class G4GeomTools implementation
31 //
32 // Author: evgueni.tcherniaev@cern.ch
33 //
34 // 10.10.2016 E.Tcherniaev: initial version.
35 // --------------------------------------------------------------------
36 
37 #include "G4GeomTools.hh"
38 
39 #include "geomdefs.hh"
40 #include "G4SystemOfUnits.hh"
41 #include "G4GeometryTolerance.hh"
42 
44 //
45 // Calculate area of a triangle in 2D
46 
48  G4double Bx, G4double By,
49  G4double Cx, G4double Cy)
50 {
51  return ((Bx-Ax)*(Cy-Ay) - (By-Ay)*(Cx-Ax))*0.5;
52 }
53 
55 //
56 // Calculate area of a triangle in 2D
57 
59  const G4TwoVector& B,
60  const G4TwoVector& C)
61 {
62  G4double Ax = A.x(), Ay = A.y();
63  return ((B.x()-Ax)*(C.y()-Ay) - (B.y()-Ay)*(C.x()-Ax))*0.5;
64 }
65 
67 //
68 // Calculate area of a quadrilateral in 2D
69 
71  const G4TwoVector& B,
72  const G4TwoVector& C,
73  const G4TwoVector& D)
74 {
75  return ((C.x()-A.x())*(D.y()-B.y()) - (C.y()-A.y())*(D.x()-B.x()))*0.5;
76 }
77 
79 //
80 // Calculate area of a polygon in 2D
81 
83 {
84  G4int n = p.size();
85  if (n < 3) return 0; // degerate polygon
86  G4double area = p[n-1].x()*p[0].y() - p[0].x()*p[n-1].y();
87  for(G4int i=1; i<n; ++i)
88  {
89  area += p[i-1].x()*p[i].y() - p[i].x()*p[i-1].y();
90  }
91  return area*0.5;
92 }
93 
95 //
96 // Point inside 2D triangle
97 
99  G4double Bx, G4double By,
100  G4double Cx, G4double Cy,
101  G4double Px, G4double Py)
102 
103 {
104  if ((Bx-Ax)*(Cy-Ay) - (By-Ay)*(Cx-Ax) > 0.)
105  {
106  if ((Ax-Cx)*(Py-Cy) - (Ay-Cy)*(Px-Cx) < 0.) return false;
107  if ((Bx-Ax)*(Py-Ay) - (By-Ay)*(Px-Ax) < 0.) return false;
108  if ((Cx-Bx)*(Py-By) - (Cy-By)*(Px-Bx) < 0.) return false;
109  }
110  else
111  {
112  if ((Ax-Cx)*(Py-Cy) - (Ay-Cy)*(Px-Cx) > 0.) return false;
113  if ((Bx-Ax)*(Py-Ay) - (By-Ay)*(Px-Ax) > 0.) return false;
114  if ((Cx-Bx)*(Py-By) - (Cy-By)*(Px-Bx) > 0.) return false;
115  }
116  return true;
117 }
118 
120 //
121 // Point inside 2D triangle
122 
124  const G4TwoVector& B,
125  const G4TwoVector& C,
126  const G4TwoVector& P)
127 {
128  G4double Ax = A.x(), Ay = A.y();
129  G4double Bx = B.x(), By = B.y();
130  G4double Cx = C.x(), Cy = C.y();
131  G4double Px = P.x(), Py = P.y();
132  if ((Bx-Ax)*(Cy-Ay) - (By-Ay)*(Cx-Ax) > 0.)
133  {
134  if ((Ax-Cx)*(Py-Cy) - (Ay-Cy)*(Px-Cx) < 0.) return false;
135  if ((Bx-Ax)*(Py-Ay) - (By-Ay)*(Px-Ax) < 0.) return false;
136  if ((Cx-Bx)*(Py-By) - (Cy-By)*(Px-Bx) < 0.) return false;
137  }
138  else
139  {
140  if ((Ax-Cx)*(Py-Cy) - (Ay-Cy)*(Px-Cx) > 0.) return false;
141  if ((Bx-Ax)*(Py-Ay) - (By-Ay)*(Px-Ax) > 0.) return false;
142  if ((Cx-Bx)*(Py-By) - (Cy-By)*(Px-Bx) > 0.) return false;
143  }
144  return true;
145 }
146 
148 //
149 // Point inside 2D polygon
150 
152  const G4TwoVectorList& v)
153 {
154  G4int Nv = v.size();
155  G4bool in = false;
156  for (G4int i = 0, k = Nv - 1; i < Nv; k = i++)
157  {
158  if ((v[i].y() > p.y()) != (v[k].y() > p.y()))
159  {
160  G4double ctg = (v[k].x()-v[i].x())/(v[k].y()-v[i].y());
161  in ^= (p.x() < (p.y()-v[i].y())*ctg + v[i].x());
162  }
163  }
164  return in;
165 }
166 
168 //
169 // Detemine whether 2D polygon is convex or not
170 
172 {
173  static const G4double kCarTolerance =
175 
176  G4bool gotNegative = false;
177  G4bool gotPositive = false;
178  G4int n = polygon.size();
179  if (n <= 0) return false;
180  for (G4int icur=0; icur<n; ++icur)
181  {
182  G4int iprev = (icur == 0) ? n-1 : icur-1;
183  G4int inext = (icur == n-1) ? 0 : icur+1;
184  G4TwoVector e1 = polygon[icur] - polygon[iprev];
185  G4TwoVector e2 = polygon[inext] - polygon[icur];
186  G4double cross = e1.x()*e2.y() - e1.y()*e2.x();
187  if (std::abs(cross) < kCarTolerance) return false;
188  if (cross < 0) gotNegative = true;
189  if (cross > 0) gotPositive = true;
190  if (gotNegative && gotPositive) return false;
191  }
192  return true;
193 }
194 
196 //
197 // Triangulate simple polygon
198 
201 {
202  result.resize(0);
203  std::vector<G4int> triangles;
204  G4bool reply = TriangulatePolygon(polygon,triangles);
205 
206  G4int n = triangles.size();
207  for (G4int i=0; i<n; ++i) result.push_back(polygon[triangles[i]]);
208  return reply;
209 }
210 
212 //
213 // Triangulation of a simple polygon by "ear clipping"
214 
216  std::vector<G4int>& result)
217 {
218  result.resize(0);
219 
220  // allocate and initialize list of Vertices in polygon
221  //
222  G4int n = polygon.size();
223  if (n < 3) return false;
224 
225  // we want a counter-clockwise polygon in V
226  //
227  G4double area = G4GeomTools::PolygonArea(polygon);
228  G4int* V = new G4int[n];
229  if (area > 0.)
230  for (G4int i=0; i<n; ++i) V[i] = i;
231  else
232  for (G4int i=0; i<n; ++i) V[i] = (n-1)-i;
233 
234  // Triangulation: remove nv-2 Vertices, creating 1 triangle every time
235  //
236  G4int nv = n;
237  G4int count = 2*nv; // error detection counter
238  for(G4int b=nv-1; nv>2; )
239  {
240  // ERROR: if we loop, it is probably a non-simple polygon
241  if ((count--) <= 0)
242  {
243  delete[] V;
244  if (area < 0.) std::reverse(result.begin(),result.end());
245  return false;
246  }
247 
248  // three consecutive vertices in current polygon, <a,b,c>
249  G4int a = (b < nv) ? b : 0; // previous
250  b = (a+1 < nv) ? a+1 : 0; // current
251  G4int c = (b+1 < nv) ? b+1 : 0; // next
252 
253  if (CheckSnip(polygon, a,b,c, nv,V))
254  {
255  // output Triangle
256  result.push_back(V[a]);
257  result.push_back(V[b]);
258  result.push_back(V[c]);
259 
260  // remove vertex b from remaining polygon
261  nv--;
262  for(G4int i=b; i<nv; ++i) V[i] = V[i+1];
263 
264  count = 2*nv; // resest error detection counter
265  }
266  }
267  delete[] V;
268  if (area < 0.) std::reverse(result.begin(),result.end());
269  return true;
270 }
271 
273 //
274 // Helper function for "ear clipping" polygon triangulation.
275 // Check for a valid snip
276 
278  G4int a, G4int b, G4int c,
279  G4int n, const G4int* V)
280 {
281  static const G4double kCarTolerance =
283 
284  // check orientation of Triangle
285  G4double Ax = contour[V[a]].x(), Ay = contour[V[a]].y();
286  G4double Bx = contour[V[b]].x(), By = contour[V[b]].y();
287  G4double Cx = contour[V[c]].x(), Cy = contour[V[c]].y();
288  if ((Bx-Ax)*(Cy-Ay) - (By-Ay)*(Cx-Ax) < kCarTolerance) return false;
289 
290  // check that there is no point inside Triangle
291  G4double xmin = std::min(std::min(Ax,Bx),Cx);
292  G4double xmax = std::max(std::max(Ax,Bx),Cx);
293  G4double ymin = std::min(std::min(Ay,By),Cy);
294  G4double ymax = std::max(std::max(Ay,By),Cy);
295  for (G4int i=0; i<n; ++i)
296  {
297  if((i == a) || (i == b) || (i == c)) continue;
298  G4double Px = contour[V[i]].x();
299  if (Px < xmin || Px > xmax) continue;
300  G4double Py = contour[V[i]].y();
301  if (Py < ymin || Py > ymax) continue;
302  if (PointInTriangle(Ax,Ay,Bx,By,Cx,Cy,Px,Py)) return false;
303  }
304  return true;
305 }
306 
308 //
309 // Remove collinear and coincident points from 2D polygon
310 
312  std::vector<G4int>& iout,
313  G4double tolerance)
314 {
315  iout.resize(0);
316  // set tolerance squared
317  G4double delta = tolerance*tolerance;
318  // set special value to mark vertices for removal
319  G4double removeIt = kInfinity;
320 
321  G4int nv = polygon.size();
322 
323  // Main loop: check every three consecutive points, if the points
324  // are collinear then mark middle point for removal
325  //
326  G4int icur = 0, iprev = 0, inext = 0, nout = 0;
327  for (G4int i=0; i<nv; ++i)
328  {
329  icur = i; // index of current point
330 
331  for (G4int k=1; k<nv+1; ++k) // set index of previous point
332  {
333  iprev = icur - k;
334  if (iprev < 0) iprev += nv;
335  if (polygon[iprev].x() != removeIt) break;
336  }
337 
338  for (G4int k=1; k<nv+1; ++k) // set index of next point
339  {
340  inext = icur + k;
341  if (inext >= nv) inext -= nv;
342  if (polygon[inext].x() != removeIt) break;
343  }
344 
345  if (iprev == inext) break; // degenerate polygon, stop
346 
347  // Calculate parameters of triangle (iprev->icur->inext),
348  // if triangle is too small or too narrow then mark current
349  // point for removal
350  G4TwoVector e1 = polygon[iprev] - polygon[icur];
351  G4TwoVector e2 = polygon[inext] - polygon[icur];
352 
353  // Check length of edges, then check height of the triangle
354  G4double leng1 = e1.mag2();
355  G4double leng2 = e2.mag2();
356  G4double leng3 = (e2-e1).mag2();
357  if (leng1 <= delta || leng2 <= delta || leng3 <= delta)
358  {
359  polygon[icur].setX(removeIt); nout++;
360  }
361  else
362  {
363  G4double lmax = std::max(std::max(leng1,leng2),leng3);
364  G4double area = std::abs(e1.x()*e2.y()-e1.y()*e2.x())*0.5;
365  if (area/std::sqrt(lmax) <= std::abs(tolerance))
366  {
367  polygon[icur].setX(removeIt); nout++;
368  }
369  }
370  }
371 
372  // Remove marked points
373  //
374  icur = 0;
375  if (nv - nout < 3) // degenerate polygon, remove all points
376  {
377  for (G4int i=0; i<nv; ++i) iout.push_back(i);
378  polygon.resize(0);
379  nv = 0;
380  }
381  for (G4int i=0; i<nv; ++i) // move points, if required
382  {
383  if (polygon[i].x() != removeIt)
384  polygon[icur++] = polygon[i];
385  else
386  iout.push_back(i);
387  }
388  if (icur < nv) polygon.resize(icur);
389  return;
390 }
391 
393 //
394 // Find bounding rectangle of a disk sector
395 
397  G4double startPhi, G4double delPhi,
398  G4TwoVector& pmin, G4TwoVector& pmax)
399 {
400  static const G4double kCarTolerance =
402 
403  // check parameters
404  //
405  pmin.set(0,0);
406  pmax.set(0,0);
407  if (rmin < 0) return false;
408  if (rmax <= rmin + kCarTolerance) return false;
409  if (delPhi <= 0 + kCarTolerance) return false;
410 
411  // calculate extent
412  //
413  pmin.set(-rmax,-rmax);
414  pmax.set( rmax, rmax);
415  if (delPhi >= CLHEP::twopi) return true;
416 
417  DiskExtent(rmin,rmax,
418  std::sin(startPhi),std::cos(startPhi),
419  std::sin(startPhi+delPhi),std::cos(startPhi+delPhi),
420  pmin,pmax);
421  return true;
422 }
423 
425 //
426 // Find bounding rectangle of a disk sector, fast version.
427 // No check of parameters !!!
428 
430  G4double sinStart, G4double cosStart,
431  G4double sinEnd, G4double cosEnd,
432  G4TwoVector& pmin, G4TwoVector& pmax)
433 {
434  static const G4double kCarTolerance =
436 
437  // check if 360 degrees
438  //
439  pmin.set(-rmax,-rmax);
440  pmax.set( rmax, rmax);
441 
442  if (std::abs(sinEnd-sinStart) < kCarTolerance &&
443  std::abs(cosEnd-cosStart) < kCarTolerance) return;
444 
445  // get start and end quadrants
446  //
447  // 1 | 0
448  // ---+---
449  // 3 | 2
450  //
451  G4int icase = (cosEnd < 0) ? 1 : 0;
452  if (sinEnd < 0) icase += 2;
453  if (cosStart < 0) icase += 4;
454  if (sinStart < 0) icase += 8;
455 
456  switch (icase)
457  {
458  // start quadrant 0
459  case 0: // start->end : 0->0
460  if (sinEnd < sinStart) break;
461  pmin.set(rmin*cosEnd,rmin*sinStart);
462  pmax.set(rmax*cosStart,rmax*sinEnd );
463  break;
464  case 1: // start->end : 0->1
465  pmin.set(rmax*cosEnd,std::min(rmin*sinStart,rmin*sinEnd));
466  pmax.set(rmax*cosStart,rmax );
467  break;
468  case 2: // start->end : 0->2
469  pmin.set(-rmax,-rmax);
470  pmax.set(std::max(rmax*cosStart,rmax*cosEnd),rmax);
471  break;
472  case 3: // start->end : 0->3
473  pmin.set(-rmax,rmax*sinEnd);
474  pmax.set(rmax*cosStart,rmax);
475  break;
476  // start quadrant 1
477  case 4: // start->end : 1->0
478  pmin.set(-rmax,-rmax);
479  pmax.set(rmax,std::max(rmax*sinStart,rmax*sinEnd));
480  break;
481  case 5: // start->end : 1->1
482  if (sinEnd > sinStart) break;
483  pmin.set(rmax*cosEnd,rmin*sinEnd );
484  pmax.set(rmin*cosStart,rmax*sinStart);
485  break;
486  case 6: // start->end : 1->2
487  pmin.set(-rmax,-rmax);
488  pmax.set(rmax*cosEnd,rmax*sinStart);
489  break;
490  case 7: // start->end : 1->3
491  pmin.set(-rmax,rmax*sinEnd);
492  pmax.set(std::max(rmin*cosStart,rmin*cosEnd),rmax*sinStart);
493  break;
494  // start quadrant 2
495  case 8: // start->end : 2->0
496  pmin.set(std::min(rmin*cosStart,rmin*cosEnd),rmax*sinStart);
497  pmax.set(rmax,rmax*sinEnd);
498  break;
499  case 9: // start->end : 2->1
500  pmin.set(rmax*cosEnd,rmax*sinStart);
501  pmax.set(rmax,rmax);
502  break;
503  case 10: // start->end : 2->2
504  if (sinEnd < sinStart) break;
505  pmin.set(rmin*cosStart,rmax*sinStart);
506  pmax.set(rmax*cosEnd,rmin*sinEnd );
507  break;
508  case 11: // start->end : 2->3
509  pmin.set(-rmax,std::min(rmax*sinStart,rmax*sinEnd));
510  pmax.set(rmax,rmax);
511  break;
512  // start quadrant 3
513  case 12: // start->end : 3->0
514  pmin.set(rmax*cosStart,-rmax);
515  pmax.set(rmax,rmax*sinEnd);
516  break;
517  case 13: // start->end : 3->1
518  pmin.set(std::min(rmax*cosStart,rmax*cosEnd),-rmax);
519  pmax.set(rmax,rmax);
520  break;
521  case 14: // start->end : 3->2
522  pmin.set(rmax*cosStart,-rmax);
523  pmax.set(rmax*cosEnd,std::max(rmin*sinStart,rmin*sinEnd));
524  break;
525  case 15: // start->end : 3->3
526  if (sinEnd > sinStart) break;
527  pmin.set(rmax*cosStart,rmax*sinEnd);
528  pmax.set(rmin*cosEnd,rmin*sinStart);
529  break;
530  }
531  return;
532 }
533 
535 //
536 // Compute the circumference (perimeter) of an ellipse
537 
539 {
540  G4double x = std::abs(pA);
541  G4double y = std::abs(pB);
542  G4double a = std::max(x,y);
543  G4double b = std::min(x,y);
544  G4double e = std::sqrt((1. - b/a)*(1. + b/a));
545  return 4. * a * comp_ellint_2(e);
546 }
547 
549 //
550 // Compute the lateral surface area of an elliptic cone
551 
553  G4double pB,
554  G4double pH)
555 {
556  G4double x = std::abs(pA);
557  G4double y = std::abs(pB);
558  G4double h = std::abs(pH);
559  G4double a = std::max(x,y);
560  G4double b = std::min(x,y);
561  G4double e = std::sqrt((1. - b/a)*(1. + b/a)) / std::hypot(1.,b/h);
562  return 2. * a * std::hypot(b,h) * comp_ellint_2(e);
563 }
564 
566 //
567 // Compute Elliptical Integral of the Second Kind
568 //
569 // The algorithm is based upon Carlson B.C., "Computation of real
570 // or complex elliptic integrals", Numerical Algorithms,
571 // Volume 10, Issue 1, 1995 (see equations 2.36 - 2.39)
572 //
573 // The code was adopted from C code at:
574 // http://paulbourke.net/geometry/ellipsecirc/
575 
577 {
578  const G4double eps = 1. / 134217728.; // 1 / 2^27
579 
580  G4double a = 1.;
581  G4double b = std::sqrt((1. - e)*(1. + e));
582  if (b == 1.) return CLHEP::halfpi;
583  if (b == 0.) return 1.;
584 
585  G4double x = 1.;
586  G4double y = b;
587  G4double S = 0.;
588  G4double M = 1.;
589  while (x - y > eps*y) {
590  G4double tmp = (x + y) * 0.5;
591  y = std::sqrt(x*y);
592  x = tmp;
593  M += M;
594  S += M * (x - y)*(x - y);
595  }
596  return 0.5 * CLHEP::halfpi * ((a + b)*(a + b) - S) / (x + y);
597 }
598 
600 //
601 // Calcuate area of a triangle in 3D
602 
604  const G4ThreeVector& B,
605  const G4ThreeVector& C)
606 {
607  return ((B-A).cross(C-A))*0.5;
608 }
609 
611 //
612 // Calcuate area of a quadrilateral in 3D
613 
615  const G4ThreeVector& B,
616  const G4ThreeVector& C,
617  const G4ThreeVector& D)
618 {
619  return ((C-A).cross(D-B))*0.5;
620 }
621 
623 //
624 // Calculate area of a polygon in 3D
625 
627 {
628  G4int n = p.size();
629  if (n < 3) return G4ThreeVector(0,0,0); // degerate polygon
630  G4ThreeVector normal = p[n-1].cross(p[0]);
631  for(G4int i=1; i<n; ++i)
632  {
633  normal += p[i-1].cross(p[i]);
634  }
635  return normal*0.5;
636 }
637 
639 //
640 // Calculate distance between point P and line segment AB in 3D
641 
643  const G4ThreeVector& A,
644  const G4ThreeVector& B)
645 {
646  G4ThreeVector AP = P - A;
647  G4ThreeVector AB = B - A;
648 
649  G4double u = AP.dot(AB);
650  if (u <= 0) return AP.mag(); // closest point is A
651 
652  G4double len2 = AB.mag2();
653  if (u >= len2) return (B-P).mag(); // closest point is B
654 
655  return ((u/len2)*AB - AP).mag(); // distance to line
656 }
657 
659 //
660 // Find closest point on line segment in 3D
661 
664  const G4ThreeVector& A,
665  const G4ThreeVector& B)
666 {
667  G4ThreeVector AP = P - A;
668  G4ThreeVector AB = B - A;
669 
670  G4double u = AP.dot(AB);
671  if (u <= 0) return A; // closest point is A
672 
673  G4double len2 = AB.mag2();
674  if (u >= len2) return B; // closest point is B
675 
676  G4double t = u/len2;
677  return A + t*AB; // closest point on segment
678 }
679 
681 //
682 // Find closest point on triangle in 3D.
683 //
684 // The implementation is based on the algorithm published in
685 // "Geometric Tools for Computer Graphics", Philip J Scheider and
686 // David H Eberly, Elsevier Science (USA), 2003.
687 //
688 // The algorithm is also available at:
689 // http://www.geometrictools.com/Documentation/DistancePoint3Triangle3.pdf
690 
693  const G4ThreeVector& A,
694  const G4ThreeVector& B,
695  const G4ThreeVector& C)
696 {
697  G4ThreeVector diff = A - P;
698  G4ThreeVector edge0 = B - A;
699  G4ThreeVector edge1 = C - A;
700 
701  G4double a = edge0.mag2();
702  G4double b = edge0.dot(edge1);
703  G4double c = edge1.mag2();
704  G4double d = diff.dot(edge0);
705  G4double e = diff.dot(edge1);
706 
707  G4double det = a*c - b*b;
708  G4double t0 = b*e - c*d;
709  G4double t1 = b*d - a*e;
710 
711  /*
712  ^ t1
713  \ 2 |
714  \ |
715  \ | regions
716  \|
717  C
718  |\
719  3 | \ 1
720  | \
721  | 0 \
722  | \
723  ---- A --- B ----> t0
724  | \
725  4 | 5 \ 6
726  | \
727  */
728 
729  G4int region = -1;
730  if (t0+t1 <= det)
731  region = (t0 < 0) ? ((t1 < 0) ? 4 : 3) : ((t1 < 0) ? 5 : 0);
732  else
733  region = (t0 < 0) ? 2 : ((t1 < 0) ? 6 : 1);
734 
735  switch (region)
736  {
737  case 0: // interior of triangle
738  {
739  G4double invDet = 1./det;
740  return A + (t0*invDet)*edge0 + (t1*invDet)*edge1;
741  }
742  case 1: // edge BC
743  {
744  G4double numer = c + e - b - d;
745  if (numer <= 0) return C;
746  G4double denom = a - 2*b + c;
747  return (numer >= denom) ? B : C + (numer/denom)*(edge0-edge1);
748  }
749  case 2: // edge AC or BC
750  {
751  G4double tmp0 = b + d;
752  G4double tmp1 = c + e;
753  if (tmp1 > tmp0)
754  {
755  G4double numer = tmp1 - tmp0;
756  G4double denom = a - 2*b + c;
757  return (numer >= denom) ? B : C + (numer/denom)*(edge0-edge1);
758  }
759  // same: (e >= 0) ? A : ((-e >= c) ? C : A + (-e/c)*edge1)
760  return (tmp1 <= 0) ? C : (( e >= 0) ? A : A + (-e/c)*edge1);
761  }
762  case 3: // edge AC
763  return (e >= 0) ? A : ((-e >= c) ? C : A + (-e/c)*edge1);
764 
765  case 4: // edge AB or AC
766  if (d < 0) return (-d >= a) ? B : A + (-d/a)*edge0;
767  return (e >= 0) ? A : ((-e >= c) ? C : A + (-e/c)*edge1);
768 
769  case 5: // edge AB
770  return (d >= 0) ? A : ((-d >= a) ? B : A + (-d/a)*edge0);
771 
772  case 6: // edge AB or BC
773  {
774  G4double tmp0 = b + e;
775  G4double tmp1 = a + d;
776  if (tmp1 > tmp0)
777  {
778  G4double numer = tmp1 - tmp0;
779  G4double denom = a - 2*b + c;
780  return (numer >= denom) ? C : B + (numer/denom)*(edge1-edge0);
781  }
782  // same: (d >= 0) ? A : ((-d >= a) ? B : A + (-d/a)*edge0)
783  return (tmp1 <= 0) ? B : (( d >= 0) ? A : A + (-d/a)*edge0);
784  }
785  default: // impossible case
787  }
788 }
789 
791 //
792 // Calculate bounding box of a spherical sector
793 
794 G4bool
796  G4double startTheta, G4double delTheta,
797  G4double startPhi, G4double delPhi,
798  G4ThreeVector& pmin, G4ThreeVector& pmax)
799 {
800  static const G4double kCarTolerance =
802 
803  // check parameters
804  //
805  pmin.set(0,0,0);
806  pmax.set(0,0,0);
807  if (rmin < 0) return false;
808  if (rmax <= rmin + kCarTolerance) return false;
809  if (delTheta <= 0 + kCarTolerance) return false;
810  if (delPhi <= 0 + kCarTolerance) return false;
811 
812  G4double stheta = startTheta;
813  G4double dtheta = delTheta;
814  if (stheta < 0 && stheta > CLHEP::pi) return false;
815  if (stheta + dtheta > CLHEP::pi) dtheta = CLHEP::pi - stheta;
816  if (dtheta <= 0 + kCarTolerance) return false;
817 
818  // calculate extent
819  //
820  pmin.set(-rmax,-rmax,-rmax);
821  pmax.set( rmax, rmax, rmax);
822  if (dtheta >= CLHEP::pi && delPhi >= CLHEP::twopi) return true;
823 
824  G4double etheta = stheta + dtheta;
825  G4double sinStart = std::sin(stheta);
826  G4double cosStart = std::cos(stheta);
827  G4double sinEnd = std::sin(etheta);
828  G4double cosEnd = std::cos(etheta);
829 
830  G4double rhomin = rmin*std::min(sinStart,sinEnd);
831  G4double rhomax = rmax;
832  if (stheta > CLHEP::halfpi) rhomax = rmax*sinStart;
833  if (etheta < CLHEP::halfpi) rhomax = rmax*sinEnd;
834 
835  G4TwoVector xymin,xymax;
836  DiskExtent(rhomin,rhomax,
837  std::sin(startPhi),std::cos(startPhi),
838  std::sin(startPhi+delPhi),std::cos(startPhi+delPhi),
839  xymin,xymax);
840 
841  G4double zmin = std::min(rmin*cosEnd,rmax*cosEnd);
842  G4double zmax = std::max(rmin*cosStart,rmax*cosStart);
843  pmin.set(xymin.x(),xymin.y(),zmin);
844  pmax.set(xymax.x(),xymax.y(),zmax);
845  return true;
846 }
Float_t x
Definition: compare.C:6
void set(double x, double y, double z)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static G4bool TriangulatePolygon(const G4TwoVectorList &polygon, G4TwoVectorList &result)
Definition: G4GeomTools.cc:199
CLHEP::Hep3Vector G4ThreeVector
static G4double comp_ellint_2(G4double e)
Definition: G4GeomTools.cc:576
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
TTree * t1
Definition: plottest35.C:26
static const G4double kInfinity
Definition: geomdefs.hh:42
static G4ThreeVector ClosestPointOnTriangle(const G4ThreeVector &P, const G4ThreeVector &A, const G4ThreeVector &B, const G4ThreeVector &C)
Definition: G4GeomTools.cc:692
static G4double DistancePointSegment(const G4ThreeVector &P, const G4ThreeVector &A, const G4ThreeVector &B)
Definition: G4GeomTools.cc:642
Float_t y
Definition: compare.C:6
double D(double temp)
const char * p
Definition: xmltok.h:285
double mag2() const
static G4ThreeVector PolygonAreaNormal(const G4ThreeVectorList &polygon)
Definition: G4GeomTools.cc:626
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
double dot(const Hep3Vector &) const
const G4double kCarTolerance
Float_t tmp
static G4bool DiskExtent(G4double rmin, G4double rmax, G4double startPhi, G4double delPhi, G4TwoVector &pmin, G4TwoVector &pmax)
Definition: G4GeomTools.cc:396
double S(double temp)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
std::vector< G4ThreeVector > G4ThreeVectorList
double x() const
static G4ThreeVector QuadAreaNormal(const G4ThreeVector &A, const G4ThreeVector &B, const G4ThreeVector &C, const G4ThreeVector &D)
Definition: G4GeomTools.cc:614
double A(double temperature)
Float_t d
double mag2() const
static G4bool IsConvex(const G4TwoVectorList &polygon)
Definition: G4GeomTools.cc:171
static G4double TriangleArea(G4double Ax, G4double Ay, G4double Bx, G4double By, G4double Cx, G4double Cy)
Definition: G4GeomTools.cc:47
static G4double EllipsePerimeter(G4double a, G4double b)
Definition: G4GeomTools.cc:538
Hep3Vector cross(const Hep3Vector &) const
G4double G4ParticleHPJENDLHEData::G4double result
static double P[]
static constexpr double twopi
Definition: SystemOfUnits.h:55
static G4double EllipticConeLateralArea(G4double a, G4double b, G4double h)
Definition: G4GeomTools.cc:552
static G4bool PointInTriangle(G4double Px, G4double Py, G4double Ax, G4double Ay, G4double Bx, G4double By, G4double Cx, G4double Cy)
Definition: G4GeomTools.cc:98
double mag() const
static G4ThreeVector ClosestPointOnSegment(const G4ThreeVector &P, const G4ThreeVector &A, const G4ThreeVector &B)
Definition: G4GeomTools.cc:663
static G4bool CheckSnip(const G4TwoVectorList &contour, G4int a, G4int b, G4int c, G4int n, const G4int *V)
Definition: G4GeomTools.cc:277
int G4int
Definition: G4Types.hh:78
std::vector< G4TwoVector > G4TwoVectorList
Definition: G4GeomTools.hh:50
ifstream in
Definition: comparison.C:7
void set(double x, double y)
double C(double temp)
static G4double QuadArea(const G4TwoVector &A, const G4TwoVector &B, const G4TwoVector &C, const G4TwoVector &D)
Definition: G4GeomTools.cc:70
double y() const
static void RemoveRedundantVertices(G4TwoVectorList &polygon, std::vector< G4int > &iout, G4double tolerance=0)
Definition: G4GeomTools.cc:311
static G4GeometryTolerance * GetInstance()
Char_t n[5]
static G4bool SphereExtent(G4double rmin, G4double rmax, G4double startTheta, G4double delTheta, G4double startPhi, G4double delPhi, G4ThreeVector &pmin, G4ThreeVector &pmax)
Definition: G4GeomTools.cc:795
static const G4double eps
static G4ThreeVector TriangleAreaNormal(const G4ThreeVector &A, const G4ThreeVector &B, const G4ThreeVector &C)
Definition: G4GeomTools.cc:603
static constexpr double halfpi
Definition: SystemOfUnits.h:56
double B(double temperature)
static G4bool PointInPolygon(const G4TwoVector &P, const G4TwoVectorList &Polygon)
Definition: G4GeomTools.cc:151
G4double GetSurfaceTolerance() const
static constexpr double pi
Definition: SystemOfUnits.h:54
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static G4double PolygonArea(const G4TwoVectorList &polygon)
Definition: G4GeomTools.cc:82