Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PolyhedraSide.cc
Go to the documentation of this file.
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: G4PolyhedraSide.cc 100427 2016-10-21 12:58:28Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class source file
32 //
33 //
34 // G4PolyhedraSide.cc
35 //
36 // Implementation of the face representing one segmented side of a Polyhedra
37 //
38 // --------------------------------------------------------------------
39 
40 #include "G4PolyhedraSide.hh"
41 #include "G4PhysicalConstants.hh"
42 #include "G4IntersectingCone.hh"
43 #include "G4ClippablePolygon.hh"
44 #include "G4AffineTransform.hh"
45 #include "G4SolidExtentList.hh"
46 #include "G4GeometryTolerance.hh"
47 
48 #include "Randomize.hh"
49 
50 // This new field helps to use the class G4PhSideManager.
51 //
52 G4PhSideManager G4PolyhedraSide::subInstanceManager;
53 
54 // This macro changes the references to fields that are now encapsulated
55 // in the class G4PhSideData.
56 //
57 #define G4MT_phphi ((subInstanceManager.offset[instanceID]).fPhi)
58 
59 // Returns the private data instance manager.
60 //
62 {
63  return subInstanceManager;
64 }
65 
66 //
67 // Constructor
68 //
69 // Values for r1,z1 and r2,z2 should be specified in clockwise
70 // order in (r,z).
71 //
73  const G4PolyhedraSideRZ *tail,
74  const G4PolyhedraSideRZ *head,
75  const G4PolyhedraSideRZ *nextRZ,
76  G4int theNumSide,
77  G4double thePhiStart,
78  G4double thePhiTotal,
79  G4bool thePhiIsOpen,
80  G4bool isAllBehind )
81 {
82 
83  instanceID = subInstanceManager.CreateSubInstance();
84 
86  fSurfaceArea=0.;
87  G4MT_phphi.first = G4ThreeVector(0,0,0);
88  G4MT_phphi.second= 0.0;
89 
90  //
91  // Record values
92  //
93  r[0] = tail->r; z[0] = tail->z;
94  r[1] = head->r; z[1] = head->z;
95 
96  G4double phiTotal;
97 
98  //
99  // Set phi to our convention
100  //
101  startPhi = thePhiStart;
102  while (startPhi < 0.0) // Loop checking, 13.08.2015, G.Cosmo
103  startPhi += twopi;
104 
105  phiIsOpen = thePhiIsOpen;
106  phiTotal = (phiIsOpen) ? thePhiTotal : twopi;
107 
108  allBehind = isAllBehind;
109 
110  //
111  // Make our intersecting cone
112  //
113  cone = new G4IntersectingCone( r, z );
114 
115  //
116  // Construct side plane vector set
117  //
118  numSide = theNumSide;
119  deltaPhi = phiTotal/theNumSide;
120  endPhi = startPhi+phiTotal;
121 
123 
125 
126  //
127  // ...this is where we start
128  //
129  G4double phi = startPhi;
130  G4ThreeVector a1( r[0]*std::cos(phi), r[0]*std::sin(phi), z[0] ),
131  b1( r[1]*std::cos(phi), r[1]*std::sin(phi), z[1] ),
132  c1( prevRZ->r*std::cos(phi), prevRZ->r*std::sin(phi), prevRZ->z ),
133  d1( nextRZ->r*std::cos(phi), nextRZ->r*std::sin(phi), nextRZ->z ),
134  a2, b2, c2, d2;
135  G4PolyhedraSideEdge *edge = edges;
136 
137  G4PolyhedraSideVec *vec = vecs;
138  do // Loop checking, 13.08.2015, G.Cosmo
139  {
140  //
141  // ...this is where we are going
142  //
143  phi += deltaPhi;
144  a2 = G4ThreeVector( r[0]*std::cos(phi), r[0]*std::sin(phi), z[0] );
145  b2 = G4ThreeVector( r[1]*std::cos(phi), r[1]*std::sin(phi), z[1] );
146  c2 = G4ThreeVector( prevRZ->r*std::cos(phi), prevRZ->r*std::sin(phi), prevRZ->z );
147  d2 = G4ThreeVector( nextRZ->r*std::cos(phi), nextRZ->r*std::sin(phi), nextRZ->z );
148 
149  G4ThreeVector tt;
150 
151  //
152  // ...build some relevant vectors.
153  // the point is to sacrifice a little memory with precalcs
154  // to gain speed
155  //
156  vec->center = 0.25*( a1 + a2 + b1 + b2 );
157 
158  tt = b2 + b1 - a2 - a1;
159  vec->surfRZ = tt.unit();
160  if (vec==vecs) lenRZ = 0.25*tt.mag();
161 
162  tt = b2 - b1 + a2 - a1;
163  vec->surfPhi = tt.unit();
164  if (vec==vecs)
165  {
166  lenPhi[0] = 0.25*tt.mag();
167  tt = b2 - b1;
168  lenPhi[1] = (0.5*tt.mag()-lenPhi[0])/lenRZ;
169  }
170 
171  tt = vec->surfPhi.cross(vec->surfRZ);
172  vec->normal = tt.unit();
173 
174  //
175  // ...edge normals are the average of the normals of
176  // the two faces they connect.
177  //
178  // ...edge normals are necessary if we are to accurately
179  // decide if a point is "inside" a face. For non-convex
180  // shapes, it is absolutely necessary to know information
181  // on adjacent faces to accurate determine this.
182  //
183  // ...we don't need them for the phi edges, since that
184  // information is taken care of internally. The r/z edges,
185  // however, depend on the adjacent G4PolyhedraSide.
186  //
187  G4ThreeVector a12, adj;
188 
189  a12 = a2-a1;
190 
191  adj = 0.5*(c1+c2-a1-a2);
192  adj = adj.cross(a12);
193  adj = adj.unit() + vec->normal;
194  vec->edgeNorm[0] = adj.unit();
195 
196  a12 = b1-b2;
197  adj = 0.5*(d1+d2-b1-b2);
198  adj = adj.cross(a12);
199  adj = adj.unit() + vec->normal;
200  vec->edgeNorm[1] = adj.unit();
201 
202  //
203  // ...the corners are crucial. It is important that
204  // they are calculated consistently for adjacent
205  // G4PolyhedraSides, to avoid gaps caused by roundoff.
206  //
207  vec->edges[0] = edge;
208  edge->corner[0] = a1;
209  edge->corner[1] = b1;
210  edge++;
211  vec->edges[1] = edge;
212 
213  a1 = a2;
214  b1 = b2;
215  c1 = c2;
216  d1 = d2;
217  } while( ++vec < vecs+numSide );
218 
219  //
220  // Clean up hanging edge
221  //
222  if (phiIsOpen)
223  {
224  edge->corner[0] = a2;
225  edge->corner[1] = b2;
226  }
227  else
228  {
229  vecs[numSide-1].edges[1] = edges;
230  }
231 
232  //
233  // Go back and fill in remaining fields in edges
234  //
235  vec = vecs;
236  G4PolyhedraSideVec *prev = vecs+numSide-1;
237  do // Loop checking, 13.08.2015, G.Cosmo
238  {
239  edge = vec->edges[0]; // The edge between prev and vec
240 
241  //
242  // Okay: edge normal is average of normals of adjacent faces
243  //
244  G4ThreeVector eNorm = vec->normal + prev->normal;
245  edge->normal = eNorm.unit();
246 
247  //
248  // Vertex normal is average of norms of adjacent surfaces (all four)
249  // However, vec->edgeNorm is unit vector in some direction
250  // as the sum of normals of adjacent PolyhedraSide with vec.
251  // The normalization used for this vector should be the same
252  // for vec and prev.
253  //
254  eNorm = vec->edgeNorm[0] + prev->edgeNorm[0];
255  edge->cornNorm[0] = eNorm.unit();
256 
257  eNorm = vec->edgeNorm[1] + prev->edgeNorm[1];
258  edge->cornNorm[1] = eNorm.unit();
259  } while( prev=vec, ++vec < vecs + numSide );
260 
261  if (phiIsOpen)
262  {
263  // G4double rFact = std::cos(0.5*deltaPhi);
264  //
265  // If phi is open, we need to patch up normals of the
266  // first and last edges and their corresponding
267  // vertices.
268  //
269  // We use vectors that are in the plane of the
270  // face. This should be safe.
271  //
272  vec = vecs;
273 
274  G4ThreeVector normvec = vec->edges[0]->corner[0]
275  - vec->edges[0]->corner[1];
276  normvec = normvec.cross(vec->normal);
277  if (normvec.dot(vec->surfPhi) > 0) normvec = -normvec;
278 
279  vec->edges[0]->normal = normvec.unit();
280 
281  vec->edges[0]->cornNorm[0] = (vec->edges[0]->corner[0]
282  - vec->center).unit();
283  vec->edges[0]->cornNorm[1] = (vec->edges[0]->corner[1]
284  - vec->center).unit();
285 
286  //
287  // Repeat for ending phi
288  //
289  vec = vecs + numSide - 1;
290 
291  normvec = vec->edges[1]->corner[0] - vec->edges[1]->corner[1];
292  normvec = normvec.cross(vec->normal);
293  if (normvec.dot(vec->surfPhi) < 0) normvec = -normvec;
294 
295  vec->edges[1]->normal = normvec.unit();
296 
297  vec->edges[1]->cornNorm[0] = (vec->edges[1]->corner[0]
298  - vec->center).unit();
299  vec->edges[1]->cornNorm[1] = (vec->edges[1]->corner[1]
300  - vec->center).unit();
301  }
302 
303  //
304  // edgeNorm is the factor one multiplies the distance along vector phi
305  // on the surface of one of our sides in order to calculate the distance
306  // from the edge. (see routine DistanceAway)
307  //
308  edgeNorm = 1.0/std::sqrt( 1.0 + lenPhi[1]*lenPhi[1] );
309 }
310 
311 
312 //
313 // Fake default constructor - sets only member data and allocates memory
314 // for usage restricted to object persistency.
315 //
317  : numSide(0), startPhi(0.), deltaPhi(0.), endPhi(0.),
318  phiIsOpen(false), allBehind(false), cone(0), vecs(0), edges(0),
319  lenRZ(0.), edgeNorm(0.), kCarTolerance(0.), fSurfaceArea(0.), instanceID(0)
320 {
321  r[0] = r[1] = 0.;
322  z[0] = z[1] = 0.;
323  lenPhi[0] = lenPhi[1] = 0.;
324 }
325 
326 
327 //
328 // Destructor
329 //
331 {
332  delete cone;
333  delete [] vecs;
334  delete [] edges;
335 }
336 
337 
338 //
339 // Copy constructor
340 //
342  : G4VCSGface()
343 {
344  instanceID = subInstanceManager.CreateSubInstance();
345 
346  CopyStuff( source );
347 }
348 
349 
350 //
351 // Assignment operator
352 //
354 {
355  if (this == &source) return *this;
356 
357  delete cone;
358  delete [] vecs;
359  delete [] edges;
360 
361  CopyStuff( source );
362 
363  return *this;
364 }
365 
366 
367 //
368 // CopyStuff
369 //
371 {
372  //
373  // The simple stuff
374  //
375  numSide = source.numSide;
376  r[0] = source.r[0];
377  r[1] = source.r[1];
378  z[0] = source.z[0];
379  z[1] = source.z[1];
380  startPhi = source.startPhi;
381  deltaPhi = source.deltaPhi;
382  endPhi = source.endPhi;
383  phiIsOpen = source.phiIsOpen;
384  allBehind = source.allBehind;
385 
386  lenRZ = source.lenRZ;
387  lenPhi[0] = source.lenPhi[0];
388  lenPhi[1] = source.lenPhi[1];
389  edgeNorm = source.edgeNorm;
390 
391  kCarTolerance = source.kCarTolerance;
392  fSurfaceArea = source.fSurfaceArea;
393 
394  cone = new G4IntersectingCone( *source.cone );
395 
396  //
397  // Duplicate edges
398  //
399  G4int numEdges = phiIsOpen ? numSide+1 : numSide;
400  edges = new G4PolyhedraSideEdge[numEdges];
401 
402  G4PolyhedraSideEdge *edge = edges,
403  *sourceEdge = source.edges;
404  do // Loop checking, 13.08.2015, G.Cosmo
405  {
406  *edge = *sourceEdge;
407  } while( ++sourceEdge, ++edge < edges + numEdges);
408 
409  //
410  // Duplicate vecs
411  //
413 
414  G4PolyhedraSideVec *vec = vecs,
415  *sourceVec = source.vecs;
416  do // Loop checking, 13.08.2015, G.Cosmo
417  {
418  *vec = *sourceVec;
419  vec->edges[0] = edges + (sourceVec->edges[0] - source.edges);
420  vec->edges[1] = edges + (sourceVec->edges[1] - source.edges);
421  } while( ++sourceVec, ++vec < vecs + numSide );
422 }
423 
424 
425 //
426 // Intersect
427 //
428 // Decide if a line intersects the face.
429 //
430 // Arguments:
431 // p = (in) starting point of line segment
432 // v = (in) direction of line segment (assumed a unit vector)
433 // A, B = (in) 2d transform variables (see note top of file)
434 // normSign = (in) desired sign for dot product with normal (see below)
435 // surfTolerance = (in) minimum distance from the surface
436 // vecs = (in) Vector set array
437 // distance = (out) distance to surface furfilling all requirements
438 // distFromSurface = (out) distance from the surface
439 // thisNormal = (out) normal vector of the intersecting surface
440 //
441 // Return value:
442 // true if an intersection is found. Otherwise, output parameters are
443 // undefined.
444 //
445 // Notes:
446 // * normSign: if we are "inside" the shape and only want to find out how far
447 // to leave the shape, we only want to consider intersections with surfaces in
448 // which the trajectory is leaving the shape. Since the normal vectors to the
449 // surface always point outwards from the inside, this means we want the dot
450 // product of the trajectory direction v and the normal of the side normals[i]
451 // to be positive. Thus, we should specify normSign as +1.0. Otherwise, if
452 // we are outside and want to go in, normSign should be set to -1.0.
453 // Don't set normSign to zero, or you will get no intersections!
454 //
455 // * surfTolerance: see notes on argument "surfTolerance" in routine
456 // "IntersectSidePlane".
457 // ----HOWEVER---- We should *not* apply this surface tolerance if the
458 // starting point is not within phi or z of the surface. Specifically,
459 // if the starting point p angle in x/y places it on a separate side from the
460 // intersection or if the starting point p is outside the z bounds of the
461 // segment, surfTolerance must be ignored or we should *always* accept the
462 // intersection!
463 // This is simply because the sides do not have infinite extent.
464 //
465 //
467  const G4ThreeVector &v,
468  G4bool outgoing,
469  G4double surfTolerance,
470  G4double &distance,
471  G4double &distFromSurface,
473  G4bool &isAllBehind )
474 {
475  G4double normSign = outgoing ? +1 : -1;
476 
477  //
478  // ------------------TO BE IMPLEMENTED---------------------
479  // Testing the intersection of individual phi faces is
480  // pretty straight forward. The simple thing therefore is to
481  // form a loop and check them all in sequence.
482  //
483  // But, I worry about one day someone making
484  // a polygon with a thousands sides. A linear search
485  // would not be ideal in such a case.
486  //
487  // So, it would be nice to be able to quickly decide
488  // which face would be intersected. One can make a very
489  // good guess by using the intersection with a cone.
490  // However, this is only reliable in 99% of the cases.
491  //
492  // My solution: make a decent guess as to the one or
493  // two potential faces might get intersected, and then
494  // test them. If we have the wrong face, use the test
495  // to make a better guess.
496  //
497  // Since we might have two guesses, form a queue of
498  // potential intersecting faces. Keep an array of
499  // already tested faces to avoid doing one more than
500  // once.
501  //
502  // Result: at worst, an iterative search. On average,
503  // a little more than two tests would be required.
504  //
505  G4ThreeVector q = p + v;
506 
507  G4int face = 0;
508  G4PolyhedraSideVec *vec = vecs;
509  do // Loop checking, 13.08.2015, G.Cosmo
510  {
511  //
512  // Correct normal?
513  //
514  G4double dotProd = normSign*v.dot(vec->normal);
515  if (dotProd <= 0) continue;
516 
517  //
518  // Is this face in front of the point along the trajectory?
519  //
520  G4ThreeVector delta = p - vec->center;
521  distFromSurface = -normSign*delta.dot(vec->normal);
522 
523  if (distFromSurface < -surfTolerance) continue;
524 
525  //
526  // phi
527  // c -------- d ^
528  // | | |
529  // a -------- b +---> r/z
530  //
531  //
532  // Do we remain on this particular segment?
533  //
534  G4ThreeVector qc = q - vec->edges[1]->corner[0];
535  G4ThreeVector qd = q - vec->edges[1]->corner[1];
536 
537  if (normSign*qc.cross(qd).dot(v) < 0) continue;
538 
539  G4ThreeVector qa = q - vec->edges[0]->corner[0];
540  G4ThreeVector qb = q - vec->edges[0]->corner[1];
541 
542  if (normSign*qa.cross(qb).dot(v) > 0) continue;
543 
544  //
545  // We found the one and only segment we might be intersecting.
546  // Do we remain within r/z bounds?
547  //
548 
549  if (r[0] > 1/kInfinity && normSign*qa.cross(qc).dot(v) < 0) return false;
550  if (r[1] > 1/kInfinity && normSign*qb.cross(qd).dot(v) > 0) return false;
551 
552  //
553  // We allow the face to be slightly behind the trajectory
554  // (surface tolerance) only if the point p is within
555  // the vicinity of the face
556  //
557  if (distFromSurface < 0)
558  {
559  G4ThreeVector ps = p - vec->center;
560 
561  G4double rz = ps.dot(vec->surfRZ);
562  if (std::fabs(rz) > lenRZ+surfTolerance) return false;
563 
564  G4double pp = ps.dot(vec->surfPhi);
565  if (std::fabs(pp) > lenPhi[0] + lenPhi[1]*rz + surfTolerance) return false;
566  }
567 
568 
569  //
570  // Intersection found. Return answer.
571  //
572  distance = distFromSurface/dotProd;
573  normal = vec->normal;
574  isAllBehind = allBehind;
575  return true;
576  } while( ++vec, ++face < numSide );
577 
578  //
579  // Oh well. Better luck next time.
580  //
581  return false;
582 }
583 
584 
586 {
587  G4double normSign = outgoing ? -1 : +1;
588 
589  //
590  // Try the closest phi segment first
591  //
592  G4int iPhi = ClosestPhiSegment( GetPhi(p) );
593 
594  G4ThreeVector pdotc = p - vecs[iPhi].center;
595  G4double normDist = pdotc.dot(vecs[iPhi].normal);
596 
597  if (normSign*normDist > -0.5*kCarTolerance)
598  {
599  return DistanceAway( p, vecs[iPhi], &normDist );
600  }
601 
602  //
603  // Now we have an interesting problem... do we try to find the
604  // closest facing side??
605  //
606  // Considered carefully, the answer is no. We know that if we
607  // are asking for the distance out, we are supposed to be inside,
608  // and vice versa.
609  //
610 
611  return kInfinity;
612 }
613 
614 
615 //
616 // Inside
617 //
619  G4double tolerance,
620  G4double *bestDistance )
621 {
622  //
623  // Which phi segment is closest to this point?
624  //
625  G4int iPhi = ClosestPhiSegment( GetPhi(p) );
626 
627  G4double norm;
628 
629  //
630  // Get distance to this segment
631  //
632  *bestDistance = DistanceToOneSide( p, vecs[iPhi], &norm );
633 
634  //
635  // Use distance along normal to decide return value
636  //
637  if ( (std::fabs(norm) < tolerance) && (*bestDistance < 2.0*tolerance) )
638  return kSurface;
639  else if (norm < 0)
640  return kInside;
641  else
642  return kOutside;
643 }
644 
645 
646 //
647 // Normal
648 //
650  G4double *bestDistance )
651 {
652  //
653  // Which phi segment is closest to this point?
654  //
655  G4int iPhi = ClosestPhiSegment( GetPhi(p) );
656 
657  //
658  // Get distance to this segment
659  //
660  G4double norm;
661  *bestDistance = DistanceToOneSide( p, vecs[iPhi], &norm );
662 
663  return vecs[iPhi].normal;
664 }
665 
666 
667 //
668 // Extent
669 //
671 {
672  if (axis.perp2() < DBL_MIN)
673  {
674  //
675  // Special case
676  //
677  return axis.z() < 0 ? -cone->ZLo() : cone->ZHi();
678  }
679 
680  G4int iPhi, i1, i2;
681  G4double best;
682  G4ThreeVector *list[4];
683 
684  //
685  // Which phi segment, if any, does the axis belong to
686  //
687  iPhi = PhiSegment( GetPhi(axis) );
688 
689  if (iPhi < 0)
690  {
691  //
692  // No phi segment? Check front edge of first side and
693  // last edge of second side
694  //
695  i1 = 0; i2 = numSide-1;
696  }
697  else
698  {
699  //
700  // Check all corners of matching phi side
701  //
702  i1 = iPhi; i2 = iPhi;
703  }
704 
705  list[0] = vecs[i1].edges[0]->corner;
706  list[1] = vecs[i1].edges[0]->corner+1;
707  list[2] = vecs[i2].edges[1]->corner;
708  list[3] = vecs[i2].edges[1]->corner+1;
709 
710  //
711  // Who's biggest?
712  //
713  best = -kInfinity;
714  G4ThreeVector **vec = list;
715  do // Loop checking, 13.08.2015, G.Cosmo
716  {
717  G4double answer = (*vec)->dot(axis);
718  if (answer > best) best = answer;
719  } while( ++vec < list+4 );
720 
721  return best;
722 }
723 
724 
725 //
726 // CalculateExtent
727 //
728 // See notes in G4VCSGface
729 //
731  const G4VoxelLimits &voxelLimit,
732  const G4AffineTransform &transform,
733  G4SolidExtentList &extentList )
734 {
735  //
736  // Loop over all sides
737  //
738  G4PolyhedraSideVec *vec = vecs;
739  do // Loop checking, 13.08.2015, G.Cosmo
740  {
741  //
742  // Fill our polygon with the four corners of
743  // this side, after the specified transformation
744  //
745  G4ClippablePolygon polygon;
746 
747  polygon.AddVertexInOrder(transform.
748  TransformPoint(vec->edges[0]->corner[0]));
749  polygon.AddVertexInOrder(transform.
750  TransformPoint(vec->edges[0]->corner[1]));
751  polygon.AddVertexInOrder(transform.
752  TransformPoint(vec->edges[1]->corner[1]));
753  polygon.AddVertexInOrder(transform.
754  TransformPoint(vec->edges[1]->corner[0]));
755 
756  //
757  // Get extent
758  //
759  if (polygon.PartialClip( voxelLimit, axis ))
760  {
761  //
762  // Get dot product of normal along target axis
763  //
764  polygon.SetNormal( transform.TransformAxis(vec->normal) );
765 
766  extentList.AddSurface( polygon );
767  }
768  } while( ++vec < vecs+numSide );
769 
770  return;
771 }
772 
773 
774 //
775 // IntersectSidePlane
776 //
777 // Decide if a line correctly intersects one side plane of our segment.
778 // It is assumed that the correct side has been chosen, and thus only
779 // the z bounds (of the entire segment) are checked.
780 //
781 // normSign - To be multiplied against normal:
782 // = +1.0 normal is unchanged
783 // = -1.0 normal is reversed (now points inward)
784 //
785 // Arguments:
786 // p - (in) Point
787 // v - (in) Direction
788 // vec - (in) Description record of the side plane
789 // normSign - (in) Sign (+/- 1) to apply to normal
790 // surfTolerance - (in) Surface tolerance (generally > 0, see below)
791 // distance - (out) Distance along v to intersection
792 // distFromSurface - (out) Distance from surface normal
793 //
794 // Notes:
795 // surfTolerance - Used to decide if a point is behind the surface,
796 // a point is allow to be -surfTolerance behind the
797 // surface (as measured along the normal), but *only*
798 // if the point is within the r/z bounds + surfTolerance
799 // of the segment.
800 //
802  const G4ThreeVector &v,
803  const G4PolyhedraSideVec& vec,
804  G4double normSign,
805  G4double surfTolerance,
806  G4double &distance,
807  G4double &distFromSurface )
808 {
809  //
810  // Correct normal? Here we have straight sides, and can safely ignore
811  // intersections where the dot product with the normal is zero.
812  //
813  G4double dotProd = normSign*v.dot(vec.normal);
814 
815  if (dotProd <= 0) return false;
816 
817  //
818  // Calculate distance to surface. If the side is too far
819  // behind the point, we must reject it.
820  //
821  G4ThreeVector delta = p - vec.center;
822  distFromSurface = -normSign*delta.dot(vec.normal);
823 
824  if (distFromSurface < -surfTolerance) return false;
825 
826  //
827  // Calculate precise distance to intersection with the side
828  // (along the trajectory, not normal to the surface)
829  //
830  distance = distFromSurface/dotProd;
831 
832  //
833  // Do we fall off the r/z extent of the segment?
834  //
835  // Calculate this very, very carefully! Why?
836  // 1. If a RZ end is at R=0, you can't miss!
837  // 2. If you just fall off in RZ, the answer must
838  // be consistent with adjacent G4PolyhedraSide faces.
839  // (2) implies that only variables used by other G4PolyhedraSide
840  // faces may be used, which includes only: p, v, and the edge corners.
841  // It also means that one side is a ">" or "<", which the other
842  // must be ">=" or "<=". Fortunately, this isn't a new problem.
843  // The solution below I borrowed from Joseph O'Rourke,
844  // "Computational Geometry in C (Second Edition)"
845  // See: http://cs.smith.edu/~orourke/
846  //
847  G4ThreeVector ic = p + distance*v - vec.center;
848  G4double atRZ = vec.surfRZ.dot(ic);
849 
850  if (atRZ < 0)
851  {
852  if (r[0]==0) return true; // Can't miss!
853 
854  if (atRZ < -lenRZ*1.2) return false; // Forget it! Missed by a mile.
855 
856  G4ThreeVector q = p + v;
857  G4ThreeVector qa = q - vec.edges[0]->corner[0],
858  qb = q - vec.edges[1]->corner[0];
859  G4ThreeVector qacb = qa.cross(qb);
860  if (normSign*qacb.dot(v) < 0) return false;
861 
862  if (distFromSurface < 0)
863  {
864  if (atRZ < -lenRZ-surfTolerance) return false;
865  }
866  }
867  else if (atRZ > 0)
868  {
869  if (r[1]==0) return true; // Can't miss!
870 
871  if (atRZ > lenRZ*1.2) return false; // Missed by a mile
872 
873  G4ThreeVector q = p + v;
874  G4ThreeVector qa = q - vec.edges[0]->corner[1],
875  qb = q - vec.edges[1]->corner[1];
876  G4ThreeVector qacb = qa.cross(qb);
877  if (normSign*qacb.dot(v) >= 0) return false;
878 
879  if (distFromSurface < 0)
880  {
881  if (atRZ > lenRZ+surfTolerance) return false;
882  }
883  }
884 
885  return true;
886 }
887 
888 
889 //
890 // LineHitsSegments
891 //
892 // Calculate which phi segments a line intersects in three dimensions.
893 // No check is made as to whether the intersections are within the z bounds of
894 // the segment.
895 //
897  const G4ThreeVector &v,
898  G4int *i1, G4int *i2 )
899 {
900  G4double s1, s2;
901  //
902  // First, decide if and where the line intersects the cone
903  //
904  G4int n = cone->LineHitsCone( p, v, &s1, &s2 );
905 
906  if (n==0) return 0;
907 
908  //
909  // Try first intersection.
910  //
911  *i1 = PhiSegment( std::atan2( p.y() + s1*v.y(), p.x() + s1*v.x() ) );
912  if (n==1)
913  {
914  return (*i1 < 0) ? 0 : 1;
915  }
916 
917  //
918  // Try second intersection
919  //
920  *i2 = PhiSegment( std::atan2( p.y() + s2*v.y(), p.x() + s2*v.x() ) );
921  if (*i1 == *i2) return 0;
922 
923  if (*i1 < 0)
924  {
925  if (*i2 < 0) return 0;
926  *i1 = *i2;
927  return 1;
928  }
929 
930  if (*i2 < 0) return 1;
931 
932  return 2;
933 }
934 
935 
936 //
937 // ClosestPhiSegment
938 //
939 // Decide which phi segment is closest in phi to the point.
940 // The result is the same as PhiSegment if there is no phi opening.
941 //
943 {
944  G4int iPhi = PhiSegment( phi0 );
945  if (iPhi >= 0) return iPhi;
946 
947  //
948  // Boogers! The points falls inside the phi segment.
949  // Look for the closest point: the start, or end
950  //
951  G4double phi = phi0;
952 
953  while( phi < startPhi ) // Loop checking, 13.08.2015, G.Cosmo
954  phi += twopi;
955  G4double d1 = phi-endPhi;
956 
957  while( phi > startPhi ) // Loop checking, 13.08.2015, G.Cosmo
958  phi -= twopi;
959  G4double d2 = startPhi-phi;
960 
961  return (d2 < d1) ? 0 : numSide-1;
962 }
963 
964 
965 //
966 // PhiSegment
967 //
968 // Decide which phi segment an angle belongs to, counting from zero.
969 // A value of -1 indicates that the phi value is outside the shape
970 // (only possible if phiTotal < 360 degrees).
971 //
973 {
974  //
975  // How far are we from phiStart? Come up with a positive answer
976  // that is less than 2*PI
977  //
978  G4double phi = phi0 - startPhi;
979  while( phi < 0 ) // Loop checking, 13.08.2015, G.Cosmo
980  phi += twopi;
981  while( phi > twopi ) // Loop checking, 13.08.2015, G.Cosmo
982  phi -= twopi;
983 
984  //
985  // Divide
986  //
987  G4int answer = (G4int)(phi/deltaPhi);
988 
989  if (answer >= numSide)
990  {
991  if (phiIsOpen)
992  {
993  return -1; // Looks like we missed
994  }
995  else
996  {
997  answer = numSide-1; // Probably just roundoff
998  }
999  }
1000 
1001  return answer;
1002 }
1003 
1004 
1005 //
1006 // GetPhi
1007 //
1008 // Calculate Phi for a given 3-vector (point), if not already cached for the
1009 // same point, in the attempt to avoid consecutive computation of the same
1010 // quantity
1011 //
1013 {
1014  G4double val=0.;
1015 
1016  if (G4MT_phphi.first != p)
1017  {
1018  val = p.phi();
1019  G4MT_phphi.first = p;
1020  G4MT_phphi.second = val;
1021  }
1022  else
1023  {
1024  val = G4MT_phphi.second;
1025  }
1026  return val;
1027 }
1028 
1029 
1030 //
1031 // DistanceToOneSide
1032 //
1033 // Arguments:
1034 // p - (in) Point to check
1035 // vec - (in) vector set of this side
1036 // normDist - (out) distance normal to the side or edge, as appropriate, signed
1037 // Return value = total distance from the side
1038 //
1040  const G4PolyhedraSideVec &vec,
1041  G4double *normDist )
1042 {
1043  G4ThreeVector pct = p - vec.center;
1044 
1045  //
1046  // Get normal distance
1047  //
1048  *normDist = vec.normal.dot(pct);
1049 
1050  //
1051  // Add edge penalty
1052  //
1053  return DistanceAway( p, vec, normDist );
1054 }
1055 
1056 
1057 //
1058 // DistanceAway
1059 //
1060 // Add distance from side edges, if necesssary, to total distance,
1061 // and updates normDist appropriate depending on edge normals.
1062 //
1064  const G4PolyhedraSideVec &vec,
1065  G4double *normDist )
1066 {
1067  G4double distOut2;
1068  G4ThreeVector pct = p - vec.center;
1069  G4double distFaceNorm = *normDist;
1070 
1071  //
1072  // Okay, are we inside bounds?
1073  //
1074  G4double pcDotRZ = pct.dot(vec.surfRZ);
1075  G4double pcDotPhi = pct.dot(vec.surfPhi);
1076 
1077  //
1078  // Go through all permutations.
1079  // Phi
1080  // | | ^
1081  // B | H | E |
1082  // ------[1]------------[3]----- |
1083  // |XXXXXXXXXXXXXX| +----> RZ
1084  // C |XXXXXXXXXXXXXX| F
1085  // |XXXXXXXXXXXXXX|
1086  // ------[0]------------[2]----
1087  // A | G | D
1088  // | |
1089  //
1090  // It's real messy, but at least it's quick
1091  //
1092 
1093  if (pcDotRZ < -lenRZ)
1094  {
1095  G4double lenPhiZ = lenPhi[0] - lenRZ*lenPhi[1];
1096  G4double distOutZ = pcDotRZ+lenRZ;
1097  //
1098  // Below in RZ
1099  //
1100  if (pcDotPhi < -lenPhiZ)
1101  {
1102  //
1103  // ...and below in phi. Find distance to point (A)
1104  //
1105  G4double distOutPhi = pcDotPhi+lenPhiZ;
1106  distOut2 = distOutPhi*distOutPhi + distOutZ*distOutZ;
1107  G4ThreeVector pa = p - vec.edges[0]->corner[0];
1108  *normDist = pa.dot(vec.edges[0]->cornNorm[0]);
1109  }
1110  else if (pcDotPhi > lenPhiZ)
1111  {
1112  //
1113  // ...and above in phi. Find distance to point (B)
1114  //
1115  G4double distOutPhi = pcDotPhi-lenPhiZ;
1116  distOut2 = distOutPhi*distOutPhi + distOutZ*distOutZ;
1117  G4ThreeVector pb = p - vec.edges[1]->corner[0];
1118  *normDist = pb.dot(vec.edges[1]->cornNorm[0]);
1119  }
1120  else
1121  {
1122  //
1123  // ...and inside in phi. Find distance to line (C)
1124  //
1125  G4ThreeVector pa = p - vec.edges[0]->corner[0];
1126  distOut2 = distOutZ*distOutZ;
1127  *normDist = pa.dot(vec.edgeNorm[0]);
1128  }
1129  }
1130  else if (pcDotRZ > lenRZ)
1131  {
1132  G4double lenPhiZ = lenPhi[0] + lenRZ*lenPhi[1];
1133  G4double distOutZ = pcDotRZ-lenRZ;
1134  //
1135  // Above in RZ
1136  //
1137  if (pcDotPhi < -lenPhiZ)
1138  {
1139  //
1140  // ...and below in phi. Find distance to point (D)
1141  //
1142  G4double distOutPhi = pcDotPhi+lenPhiZ;
1143  distOut2 = distOutPhi*distOutPhi + distOutZ*distOutZ;
1144  G4ThreeVector pd = p - vec.edges[0]->corner[1];
1145  *normDist = pd.dot(vec.edges[0]->cornNorm[1]);
1146  }
1147  else if (pcDotPhi > lenPhiZ)
1148  {
1149  //
1150  // ...and above in phi. Find distance to point (E)
1151  //
1152  G4double distOutPhi = pcDotPhi-lenPhiZ;
1153  distOut2 = distOutPhi*distOutPhi + distOutZ*distOutZ;
1154  G4ThreeVector pe = p - vec.edges[1]->corner[1];
1155  *normDist = pe.dot(vec.edges[1]->cornNorm[1]);
1156  }
1157  else
1158  {
1159  //
1160  // ...and inside in phi. Find distance to line (F)
1161  //
1162  distOut2 = distOutZ*distOutZ;
1163  G4ThreeVector pd = p - vec.edges[0]->corner[1];
1164  *normDist = pd.dot(vec.edgeNorm[1]);
1165  }
1166  }
1167  else
1168  {
1169  G4double lenPhiZ = lenPhi[0] + pcDotRZ*lenPhi[1];
1170  //
1171  // We are inside RZ bounds
1172  //
1173  if (pcDotPhi < -lenPhiZ)
1174  {
1175  //
1176  // ...and below in phi. Find distance to line (G)
1177  //
1178  G4double distOut = edgeNorm*(pcDotPhi+lenPhiZ);
1179  distOut2 = distOut*distOut;
1180  G4ThreeVector pd = p - vec.edges[0]->corner[1];
1181  *normDist = pd.dot(vec.edges[0]->normal);
1182  }
1183  else if (pcDotPhi > lenPhiZ)
1184  {
1185  //
1186  // ...and above in phi. Find distance to line (H)
1187  //
1188  G4double distOut = edgeNorm*(pcDotPhi-lenPhiZ);
1189  distOut2 = distOut*distOut;
1190  G4ThreeVector pe = p - vec.edges[1]->corner[1];
1191  *normDist = pe.dot(vec.edges[1]->normal);
1192  }
1193  else
1194  {
1195  //
1196  // Inside bounds! No penalty.
1197  //
1198  return std::fabs(distFaceNorm);
1199  }
1200  }
1201  return std::sqrt( distFaceNorm*distFaceNorm + distOut2 );
1202 }
1203 
1204 
1205 //
1206 // Calculation of surface area of a triangle.
1207 // At the same time a random point in the triangle is given
1208 //
1210  G4ThreeVector p2,
1211  G4ThreeVector p3,
1212  G4ThreeVector *p4 )
1213 {
1214  G4ThreeVector v, w;
1215 
1216  v = p3 - p1;
1217  w = p1 - p2;
1218  G4double lambda1 = G4UniformRand();
1219  G4double lambda2 = lambda1*G4UniformRand();
1220 
1221  *p4=p2 + lambda1*w + lambda2*v;
1222  return 0.5*(v.cross(w)).mag();
1223 }
1224 
1225 
1226 //
1227 // GetPointOnPlane
1228 //
1229 // Auxiliary method for GetPointOnSurface()
1230 //
1234  G4double *Area )
1235 {
1236  G4double chose,aOne,aTwo;
1237  G4ThreeVector point1,point2;
1238  aOne = SurfaceTriangle(p0,p1,p2,&point1);
1239  aTwo = SurfaceTriangle(p2,p3,p0,&point2);
1240  *Area= aOne+aTwo;
1241 
1242  chose = G4UniformRand()*(aOne+aTwo);
1243  if( (chose>=0.) && (chose < aOne) )
1244  {
1245  return (point1);
1246  }
1247  return (point2);
1248 }
1249 
1250 
1251 //
1252 // SurfaceArea()
1253 //
1255 {
1256  if( fSurfaceArea==0. )
1257  {
1258  // Define the variables
1259  //
1260  G4double area,areas;
1261  G4ThreeVector point1;
1262  G4ThreeVector v1,v2,v3,v4;
1263  G4PolyhedraSideVec *vec = vecs;
1264  areas=0.;
1265 
1266  // Do a loop on all SideEdge
1267  //
1268  do // Loop checking, 13.08.2015, G.Cosmo
1269  {
1270  // Define 4points for a Plane or Triangle
1271  //
1272  v1=vec->edges[0]->corner[0];
1273  v2=vec->edges[0]->corner[1];
1274  v3=vec->edges[1]->corner[1];
1275  v4=vec->edges[1]->corner[0];
1276  point1=GetPointOnPlane(v1,v2,v3,v4,&area);
1277  areas+=area;
1278  } while( ++vec < vecs + numSide);
1279 
1280  fSurfaceArea=areas;
1281  }
1282  return fSurfaceArea;
1283 }
1284 
1285 
1286 //
1287 // GetPointOnFace()
1288 //
1290 {
1291  // Define the variables
1292  //
1293  std::vector<G4double>areas;
1294  std::vector<G4ThreeVector>points;
1295  G4double area=0;
1296  G4double result1;
1297  G4ThreeVector point1;
1298  G4ThreeVector v1,v2,v3,v4;
1299  G4PolyhedraSideVec *vec = vecs;
1300 
1301  // Do a loop on all SideEdge
1302  //
1303  do // Loop checking, 13.08.2015, G.Cosmo
1304  {
1305  // Define 4points for a Plane or Triangle
1306  //
1307  v1=vec->edges[0]->corner[0];
1308  v2=vec->edges[0]->corner[1];
1309  v3=vec->edges[1]->corner[1];
1310  v4=vec->edges[1]->corner[0];
1311  point1=GetPointOnPlane(v1,v2,v3,v4,&result1);
1312  points.push_back(point1);
1313  areas.push_back(result1);
1314  area+=result1;
1315  } while( ++vec < vecs+numSide );
1316 
1317  // Choose randomly one of the surfaces and point on it
1318  //
1319  G4double chose = area*G4UniformRand();
1320  G4double Achose1,Achose2;
1321  Achose1=0;Achose2=0.;
1322  G4int i=0;
1323  do // Loop checking, 13.08.2015, G.Cosmo
1324  {
1325  Achose2+=areas[i];
1326  if(chose>=Achose1 && chose<Achose2)
1327  {
1328  point1=points[i] ; break;
1329  }
1330  i++; Achose1=Achose2;
1331  } while( i<numSide );
1332 
1333  return point1;
1334 }
G4ThreeVector GetPointOnPlane(G4ThreeVector p0, G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3, G4double *Area)
static c2_factory< G4double > c2
G4double GetPhi(const G4ThreeVector &p)
G4double lenPhi[2]
G4double Extent(const G4ThreeVector axis)
G4PolyhedraSide(const G4PolyhedraSideRZ *prevRZ, const G4PolyhedraSideRZ *tail, const G4PolyhedraSideRZ *head, const G4PolyhedraSideRZ *nextRZ, G4int numSide, G4double phiStart, G4double phiTotal, G4bool phiIsOpen, G4bool isAllBehind=false)
double perp2() const
static const G4double kInfinity
Definition: geomdefs.hh:42
CLHEP::Hep3Vector G4ThreeVector
double x() const
EInside Inside(const G4ThreeVector &p, G4double tolerance, G4double *bestDistance)
double dot(const Hep3Vector &) const
virtual ~G4PolyhedraSide()
G4int ClosestPhiSegment(G4double phi)
G4double DistanceToOneSide(const G4ThreeVector &p, const G4PolyhedraSideVec &vec, G4double *normDist)
const char * p
Definition: xmltok.h:285
static const G4PhSideManager & GetSubInstanceManager()
static const G4double d2
G4double GetSurfaceTolerance() const
void SetNormal(const G4ThreeVector &newNormal)
G4ThreeVector Normal(const G4ThreeVector &p, G4double *bestDistance)
virtual G4bool PartialClip(const G4VoxelLimits &voxelLimit, const EAxis IgnoreMe)
virtual void AddVertexInOrder(const G4ThreeVector vertex)
G4int LineHitsSegments(const G4ThreeVector &p, const G4ThreeVector &v, G4int *i1, G4int *i2)
G4int LineHitsCone(const G4ThreeVector &p, const G4ThreeVector &v, G4double *s1, G4double *s2)
G4int CreateSubInstance()
G4PolyhedraSideVec * vecs
int G4int
Definition: G4Types.hh:78
G4double ZLo() const
double z() const
static constexpr double twopi
Definition: G4SIunits.hh:76
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
G4double DistanceAway(const G4ThreeVector &p, const G4PolyhedraSideVec &vec, G4double *normDist)
G4double Distance(const G4ThreeVector &p, G4bool outgoing)
G4double SurfaceTriangle(G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector *p4)
G4bool IntersectSidePlane(const G4ThreeVector &p, const G4ThreeVector &v, const G4PolyhedraSideVec &vec, G4double normSign, G4double surfTolerance, G4double &distance, G4double &distFromSurface)
#define G4UniformRand()
Definition: Randomize.hh:97
G4double ZHi() const
bool G4bool
Definition: G4Types.hh:79
void AddSurface(const G4ClippablePolygon &surface)
const G4double kCarTolerance
double phi() const
G4int PhiSegment(G4double phi)
const G4int n
static const G4double d1
tuple v
Definition: test.py:18
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, G4bool outgoing, G4double surfTolerance, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal, G4bool &allBehind)
G4ThreeVector GetPointOnFace()
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
Hep3Vector unit() const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
double y() const
G4double SurfaceArea()
#define DBL_MIN
Definition: templates.hh:75
Hep3Vector cross(const Hep3Vector &) const
void CopyStuff(const G4PolyhedraSide &source)
G4IntersectingCone * cone
G4PolyhedraSideEdge * edges
double G4double
Definition: G4Types.hh:76
static constexpr double ps
Definition: G4SIunits.hh:172
void CalculateExtent(const EAxis axis, const G4VoxelLimits &voxelLimit, const G4AffineTransform &tranform, G4SolidExtentList &extentList)
G4PolyhedraSide & operator=(const G4PolyhedraSide &source)
double mag() const
#define G4MT_phphi
static G4GeometryTolerance * GetInstance()
tuple c1
Definition: plottest35.py:14