Geant4  10.01.p01
G4PolyconeSide.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: G4PolyconeSide.cc 88948 2015-03-16 16:26:50Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class source file
32 //
33 //
34 // G4PolyconeSide.cc
35 //
36 // Implementation of the face representing one conical side of a polycone
37 //
38 // --------------------------------------------------------------------
39 
40 #include "G4PolyconeSide.hh"
41 #include "meshdefs.hh"
42 #include "G4PhysicalConstants.hh"
43 #include "G4IntersectingCone.hh"
44 #include "G4ClippablePolygon.hh"
45 #include "G4AffineTransform.hh"
46 #include "G4SolidExtentList.hh"
47 #include "G4GeometryTolerance.hh"
48 
49 #include "Randomize.hh"
50 
51 // This static member is thread local. For each thread, it points to the
52 // array of G4PlSideData instances.
53 //
54 template <class G4PlSideData> G4ThreadLocal
56 
57 // This new field helps to use the class G4PlSideManager.
58 //
60 
61 
62 // Returns the private data instance manager.
63 //
65 {
66  return subInstanceManager;
67 }
68 
69 
70 //
71 // Constructor
72 //
73 // Values for r1,z1 and r2,z2 should be specified in clockwise
74 // order in (r,z).
75 //
77  const G4PolyconeSideRZ *tail,
78  const G4PolyconeSideRZ *head,
79  const G4PolyconeSideRZ *nextRZ,
80  G4double thePhiStart,
81  G4double theDeltaPhi,
82  G4bool thePhiIsOpen,
83  G4bool isAllBehind )
84  : ncorners(0), corners(0)
85 {
86 
88 
90  fSurfaceArea = 0.0;
91  G4MT_pcphi.first = G4ThreeVector(0,0,0);
92  G4MT_pcphi.second= 0.0;
93 
94  //
95  // Record values
96  //
97  r[0] = tail->r; z[0] = tail->z;
98  r[1] = head->r; z[1] = head->z;
99 
100  phiIsOpen = thePhiIsOpen;
101  if (phiIsOpen)
102  {
103  deltaPhi = theDeltaPhi;
104  startPhi = thePhiStart;
105 
106  //
107  // Set phi values to our conventions
108  //
109  while (deltaPhi < 0.0) deltaPhi += twopi;
110  while (startPhi < 0.0) startPhi += twopi;
111 
112  //
113  // Calculate corner coordinates
114  //
115  ncorners = 4;
117 
118  corners[0] = G4ThreeVector( tail->r*std::cos(startPhi),
119  tail->r*std::sin(startPhi), tail->z );
120  corners[1] = G4ThreeVector( head->r*std::cos(startPhi),
121  head->r*std::sin(startPhi), head->z );
122  corners[2] = G4ThreeVector( tail->r*std::cos(startPhi+deltaPhi),
123  tail->r*std::sin(startPhi+deltaPhi), tail->z );
124  corners[3] = G4ThreeVector( head->r*std::cos(startPhi+deltaPhi),
125  head->r*std::sin(startPhi+deltaPhi), head->z );
126  }
127  else
128  {
129  deltaPhi = twopi;
130  startPhi = 0.0;
131  }
132 
133  allBehind = isAllBehind;
134 
135  //
136  // Make our intersecting cone
137  //
138  cone = new G4IntersectingCone( r, z );
139 
140  //
141  // Calculate vectors in r,z space
142  //
143  rS = r[1]-r[0]; zS = z[1]-z[0];
144  length = std::sqrt( rS*rS + zS*zS);
145  rS /= length; zS /= length;
146 
147  rNorm = +zS;
148  zNorm = -rS;
149 
150  G4double lAdj;
151 
152  prevRS = r[0]-prevRZ->r;
153  prevZS = z[0]-prevRZ->z;
154  lAdj = std::sqrt( prevRS*prevRS + prevZS*prevZS );
155  prevRS /= lAdj;
156  prevZS /= lAdj;
157 
158  rNormEdge[0] = rNorm + prevZS;
159  zNormEdge[0] = zNorm - prevRS;
160  lAdj = std::sqrt( rNormEdge[0]*rNormEdge[0] + zNormEdge[0]*zNormEdge[0] );
161  rNormEdge[0] /= lAdj;
162  zNormEdge[0] /= lAdj;
163 
164  nextRS = nextRZ->r-r[1];
165  nextZS = nextRZ->z-z[1];
166  lAdj = std::sqrt( nextRS*nextRS + nextZS*nextZS );
167  nextRS /= lAdj;
168  nextZS /= lAdj;
169 
170  rNormEdge[1] = rNorm + nextZS;
171  zNormEdge[1] = zNorm - nextRS;
172  lAdj = std::sqrt( rNormEdge[1]*rNormEdge[1] + zNormEdge[1]*zNormEdge[1] );
173  rNormEdge[1] /= lAdj;
174  zNormEdge[1] /= lAdj;
175 }
176 
177 
178 //
179 // Fake default constructor - sets only member data and allocates memory
180 // for usage restricted to object persistency.
181 //
183  : startPhi(0.), deltaPhi(0.), phiIsOpen(false), allBehind(false),
184  cone(0), rNorm(0.), zNorm(0.), rS(0.), zS(0.), length(0.),
185  prevRS(0.), prevZS(0.), nextRS(0.), nextZS(0.), ncorners(0), corners(0),
186  kCarTolerance(0.), fSurfaceArea(0.), instanceID(0)
187 {
188  r[0] = r[1] = 0.;
189  z[0] = z[1] = 0.;
190  rNormEdge[0]= rNormEdge[1] = 0.;
191  zNormEdge[0]= zNormEdge[1] = 0.;
192 }
193 
194 
195 //
196 // Destructor
197 //
199 {
200  delete cone;
201  if (phiIsOpen) { delete [] corners; }
202 }
203 
204 
205 //
206 // Copy constructor
207 //
209  : G4VCSGface(), ncorners(0), corners(0)
210 {
212 
213  CopyStuff( source );
214 }
215 
216 
217 //
218 // Assignment operator
219 //
221 {
222  if (this == &source) { return *this; }
223 
224  delete cone;
225  if (phiIsOpen) { delete [] corners; }
226 
227  CopyStuff( source );
228 
229  return *this;
230 }
231 
232 
233 //
234 // CopyStuff
235 //
237 {
238  r[0] = source.r[0];
239  r[1] = source.r[1];
240  z[0] = source.z[0];
241  z[1] = source.z[1];
242 
243  startPhi = source.startPhi;
244  deltaPhi = source.deltaPhi;
245  phiIsOpen = source.phiIsOpen;
246  allBehind = source.allBehind;
247 
248  kCarTolerance = source.kCarTolerance;
249  fSurfaceArea = source.fSurfaceArea;
250 
251  cone = new G4IntersectingCone( *source.cone );
252 
253  rNorm = source.rNorm;
254  zNorm = source.zNorm;
255  rS = source.rS;
256  zS = source.zS;
257  length = source.length;
258  prevRS = source.prevRS;
259  prevZS = source.prevZS;
260  nextRS = source.nextRS;
261  nextZS = source.nextZS;
262 
263  rNormEdge[0] = source.rNormEdge[0];
264  rNormEdge[1] = source.rNormEdge[1];
265  zNormEdge[0] = source.zNormEdge[0];
266  zNormEdge[1] = source.zNormEdge[1];
267 
268  if (phiIsOpen)
269  {
270  ncorners = 4;
272 
273  corners[0] = source.corners[0];
274  corners[1] = source.corners[1];
275  corners[2] = source.corners[2];
276  corners[3] = source.corners[3];
277  }
278 }
279 
280 
281 //
282 // Intersect
283 //
285  const G4ThreeVector &v,
286  G4bool outgoing,
287  G4double surfTolerance,
288  G4double &distance,
289  G4double &distFromSurface,
291  G4bool &isAllBehind )
292 {
293  G4double s1, s2;
294  G4double normSign = outgoing ? +1 : -1;
295 
296  isAllBehind = allBehind;
297 
298  //
299  // Check for two possible intersections
300  //
301  G4int nside = cone->LineHitsCone( p, v, &s1, &s2 );
302  if (nside == 0) return false;
303 
304  //
305  // Check the first side first, since it is (supposed to be) closest
306  //
307  G4ThreeVector hit = p + s1*v;
308 
309  if (PointOnCone( hit, normSign, p, v, normal ))
310  {
311  //
312  // Good intersection! What about the normal?
313  //
314  if (normSign*v.dot(normal) > 0)
315  {
316  //
317  // We have a valid intersection, but it could very easily
318  // be behind the point. To decide if we tolerate this,
319  // we have to see if the point p is on the surface near
320  // the intersecting point.
321  //
322  // What does it mean exactly for the point p to be "near"
323  // the intersection? It means that if we draw a line from
324  // p to the hit, the line remains entirely within the
325  // tolerance bounds of the cone. To test this, we can
326  // ask if the normal is correct near p.
327  //
328  G4double pr = p.perp();
329  if (pr < DBL_MIN) pr = DBL_MIN;
330  G4ThreeVector pNormal( rNorm*p.x()/pr, rNorm*p.y()/pr, zNorm );
331  if (normSign*v.dot(pNormal) > 0)
332  {
333  //
334  // p and intersection in same hemisphere
335  //
336  G4double distOutside2;
337  distFromSurface = -normSign*DistanceAway( p, false, distOutside2 );
338  if (distOutside2 < surfTolerance*surfTolerance)
339  {
340  if (distFromSurface > -surfTolerance)
341  {
342  //
343  // We are just inside or away from the
344  // surface. Accept *any* value of distance.
345  //
346  distance = s1;
347  return true;
348  }
349  }
350  }
351  else
352  distFromSurface = s1;
353 
354  //
355  // Accept positive distances
356  //
357  if (s1 > 0)
358  {
359  distance = s1;
360  return true;
361  }
362  }
363  }
364 
365  if (nside==1) return false;
366 
367  //
368  // Well, try the second hit
369  //
370  hit = p + s2*v;
371 
372  if (PointOnCone( hit, normSign, p, v, normal ))
373  {
374  //
375  // Good intersection! What about the normal?
376  //
377  if (normSign*v.dot(normal) > 0)
378  {
379  G4double pr = p.perp();
380  if (pr < DBL_MIN) pr = DBL_MIN;
381  G4ThreeVector pNormal( rNorm*p.x()/pr, rNorm*p.y()/pr, zNorm );
382  if (normSign*v.dot(pNormal) > 0)
383  {
384  G4double distOutside2;
385  distFromSurface = -normSign*DistanceAway( p, false, distOutside2 );
386  if (distOutside2 < surfTolerance*surfTolerance)
387  {
388  if (distFromSurface > -surfTolerance)
389  {
390  distance = s2;
391  return true;
392  }
393  }
394  }
395  else
396  distFromSurface = s2;
397 
398  if (s2 > 0)
399  {
400  distance = s2;
401  return true;
402  }
403  }
404  }
405 
406  //
407  // Better luck next time
408  //
409  return false;
410 }
411 
412 
413 //
414 // Distance
415 //
417 {
418  G4double normSign = outgoing ? -1 : +1;
419  G4double distFrom, distOut2;
420 
421  //
422  // We have two tries for each hemisphere. Try the closest first.
423  //
424  distFrom = normSign*DistanceAway( p, false, distOut2 );
425  if (distFrom > -0.5*kCarTolerance )
426  {
427  //
428  // Good answer
429  //
430  if (distOut2 > 0)
431  return std::sqrt( distFrom*distFrom + distOut2 );
432  else
433  return std::fabs(distFrom);
434  }
435 
436  //
437  // Try second side.
438  //
439  distFrom = normSign*DistanceAway( p, true, distOut2 );
440  if (distFrom > -0.5*kCarTolerance)
441  {
442 
443  if (distOut2 > 0)
444  return std::sqrt( distFrom*distFrom + distOut2 );
445  else
446  return std::fabs(distFrom);
447  }
448 
449  return kInfinity;
450 }
451 
452 
453 //
454 // Inside
455 //
458  G4double *bestDistance )
459 {
460  G4double distFrom, distOut2, dist2;
461  G4double edgeRZnorm;
462 
463  distFrom = DistanceAway( p, distOut2, &edgeRZnorm );
464  dist2 = distFrom*distFrom + distOut2;
465 
466  *bestDistance = std::sqrt( dist2);
467 
468  // Okay then, inside or out?
469  //
470  if ( (std::fabs(edgeRZnorm) < tolerance)
471  && (distOut2< tolerance*tolerance) )
472  return kSurface;
473  else if (edgeRZnorm < 0)
474  return kInside;
475  else
476  return kOutside;
477 }
478 
479 
480 //
481 // Normal
482 //
484  G4double *bestDistance )
485 {
486  if (p == G4ThreeVector(0.,0.,0.)) { return p; }
487 
488  G4double dFrom, dOut2;
489 
490  dFrom = DistanceAway( p, false, dOut2 );
491 
492  *bestDistance = std::sqrt( dFrom*dFrom + dOut2 );
493 
494  G4double rds = p.perp();
495  if (rds!=0.) { return G4ThreeVector(rNorm*p.x()/rds,rNorm*p.y()/rds,zNorm); }
496  return G4ThreeVector( 0.,0., zNorm ).unit();
497 }
498 
499 
500 //
501 // Extent
502 //
504 {
505  if (axis.perp2() < DBL_MIN)
506  {
507  //
508  // Special case
509  //
510  return axis.z() < 0 ? -cone->ZLo() : cone->ZHi();
511  }
512 
513  //
514  // Is the axis pointing inside our phi gap?
515  //
516  if (phiIsOpen)
517  {
518  G4double phi = GetPhi(axis);
519  while( phi < startPhi ) phi += twopi;
520 
521  if (phi > deltaPhi+startPhi)
522  {
523  //
524  // Yeah, looks so. Make four three vectors defining the phi
525  // opening
526  //
527  G4double cosP = std::cos(startPhi), sinP = std::sin(startPhi);
528  G4ThreeVector a( r[0]*cosP, r[0]*sinP, z[0] );
529  G4ThreeVector b( r[1]*cosP, r[1]*sinP, z[1] );
530  cosP = std::cos(startPhi+deltaPhi); sinP = std::sin(startPhi+deltaPhi);
531  G4ThreeVector c( r[0]*cosP, r[0]*sinP, z[0] );
532  G4ThreeVector d( r[1]*cosP, r[1]*sinP, z[1] );
533 
534  G4double ad = axis.dot(a),
535  bd = axis.dot(b),
536  cd = axis.dot(c),
537  dd = axis.dot(d);
538 
539  if (bd > ad) ad = bd;
540  if (cd > ad) ad = cd;
541  if (dd > ad) ad = dd;
542 
543  return ad;
544  }
545  }
546 
547  //
548  // Check either end
549  //
550  G4double aPerp = axis.perp();
551 
552  G4double a = aPerp*r[0] + axis.z()*z[0];
553  G4double b = aPerp*r[1] + axis.z()*z[1];
554 
555  if (b > a) a = b;
556 
557  return a;
558 }
559 
560 
561 //
562 // CalculateExtent
563 //
564 // See notes in G4VCSGface
565 //
567  const G4VoxelLimits &voxelLimit,
568  const G4AffineTransform &transform,
569  G4SolidExtentList &extentList )
570 {
571  G4ClippablePolygon polygon;
572 
573  //
574  // Here we will approximate (ala G4Cons) and divide our conical section
575  // into segments, like G4Polyhedra. When doing so, the radius
576  // is extented far enough such that the segments always lie
577  // just outside the surface of the conical section we are
578  // approximating.
579  //
580 
581  //
582  // Choose phi size of our segment(s) based on constants as
583  // defined in meshdefs.hh
584  //
585  G4int numPhi = (G4int)(deltaPhi/kMeshAngleDefault) + 1;
586  if (numPhi < kMinMeshSections)
587  numPhi = kMinMeshSections;
588  else if (numPhi > kMaxMeshSections)
589  numPhi = kMaxMeshSections;
590 
591  G4double sigPhi = deltaPhi/numPhi;
592 
593  //
594  // Determine radius factor to keep segments outside
595  //
596  G4double rFudge = 1.0/std::cos(0.5*sigPhi);
597 
598  //
599  // Decide which radius to use on each end of the side,
600  // and whether a transition mesh is required
601  //
602  // {r0,z0} - Beginning of this side
603  // {r1,z1} - Ending of this side
604  // {r2,z0} - Beginning of transition piece connecting previous
605  // side (and ends at beginning of this side)
606  //
607  // So, order is 2 --> 0 --> 1.
608  // -------
609  //
610  // r2 < 0 indicates that no transition piece is required
611  //
612  G4double r0, r1, r2, z0, z1;
613 
614  r2 = -1; // By default: no transition piece
615 
616  if (rNorm < -DBL_MIN)
617  {
618  //
619  // This side faces *inward*, and so our mesh has
620  // the same radius
621  //
622  r1 = r[1];
623  z1 = z[1];
624  z0 = z[0];
625  r0 = r[0];
626 
627  r2 = -1;
628 
629  if (prevZS > DBL_MIN)
630  {
631  //
632  // The previous side is facing outwards
633  //
634  if ( prevRS*zS - prevZS*rS > 0 )
635  {
636  //
637  // Transition was convex: build transition piece
638  //
639  if (r[0] > DBL_MIN) r2 = r[0]*rFudge;
640  }
641  else
642  {
643  //
644  // Transition was concave: short this side
645  //
646  FindLineIntersect( z0, r0, zS, rS,
647  z0, r0*rFudge, prevZS, prevRS*rFudge, z0, r0 );
648  }
649  }
650 
651  if ( nextZS > DBL_MIN && (rS*nextZS - zS*nextRS < 0) )
652  {
653  //
654  // The next side is facing outwards, forming a
655  // concave transition: short this side
656  //
657  FindLineIntersect( z1, r1, zS, rS,
658  z1, r1*rFudge, nextZS, nextRS*rFudge, z1, r1 );
659  }
660  }
661  else if (rNorm > DBL_MIN)
662  {
663  //
664  // This side faces *outward* and is given a boost to
665  // it radius
666  //
667  r0 = r[0]*rFudge;
668  z0 = z[0];
669  r1 = r[1]*rFudge;
670  z1 = z[1];
671 
672  if (prevZS < -DBL_MIN)
673  {
674  //
675  // The previous side is facing inwards
676  //
677  if ( prevRS*zS - prevZS*rS > 0 )
678  {
679  //
680  // Transition was convex: build transition piece
681  //
682  if (r[0] > DBL_MIN) r2 = r[0];
683  }
684  else
685  {
686  //
687  // Transition was concave: short this side
688  //
689  FindLineIntersect( z0, r0, zS, rS*rFudge,
690  z0, r[0], prevZS, prevRS, z0, r0 );
691  }
692  }
693 
694  if ( nextZS < -DBL_MIN && (rS*nextZS - zS*nextRS < 0) )
695  {
696  //
697  // The next side is facing inwards, forming a
698  // concave transition: short this side
699  //
700  FindLineIntersect( z1, r1, zS, rS*rFudge,
701  z1, r[1], nextZS, nextRS, z1, r1 );
702  }
703  }
704  else
705  {
706  //
707  // This side is perpendicular to the z axis (is a disk)
708  //
709  // Whether or not r0 needs a rFudge factor depends
710  // on the normal of the previous edge. Similar with r1
711  // and the next edge. No transition piece is required.
712  //
713  r0 = r[0];
714  r1 = r[1];
715  z0 = z[0];
716  z1 = z[1];
717 
718  if (prevZS > DBL_MIN) r0 *= rFudge;
719  if (nextZS > DBL_MIN) r1 *= rFudge;
720  }
721 
722  //
723  // Loop
724  //
725  G4double phi = startPhi,
726  cosPhi = std::cos(phi),
727  sinPhi = std::sin(phi);
728 
729  G4ThreeVector v0( r0*cosPhi, r0*sinPhi, z0 ),
730  v1( r1*cosPhi, r1*sinPhi, z1 ),
731  v2, w0, w1, w2;
732  transform.ApplyPointTransform( v0 );
733  transform.ApplyPointTransform( v1 );
734 
735  if (r2 >= 0)
736  {
737  v2 = G4ThreeVector( r2*cosPhi, r2*sinPhi, z0 );
738  transform.ApplyPointTransform( v2 );
739  }
740 
741  do
742  {
743  phi += sigPhi;
744  if (numPhi == 1) phi = startPhi+deltaPhi; // Try to avoid roundoff
745  cosPhi = std::cos(phi),
746  sinPhi = std::sin(phi);
747 
748  w0 = G4ThreeVector( r0*cosPhi, r0*sinPhi, z0 );
749  w1 = G4ThreeVector( r1*cosPhi, r1*sinPhi, z1 );
750  transform.ApplyPointTransform( w0 );
751  transform.ApplyPointTransform( w1 );
752 
753  G4ThreeVector deltaV = r0 > r1 ? w0-v0 : w1-v1;
754 
755  //
756  // Build polygon, taking special care to keep the vertices
757  // in order
758  //
759  polygon.ClearAllVertices();
760 
761  polygon.AddVertexInOrder( v0 );
762  polygon.AddVertexInOrder( v1 );
763  polygon.AddVertexInOrder( w1 );
764  polygon.AddVertexInOrder( w0 );
765 
766  //
767  // Get extent
768  //
769  if (polygon.PartialClip( voxelLimit, axis ))
770  {
771  //
772  // Get dot product of normal with target axis
773  //
774  polygon.SetNormal( deltaV.cross(v1-v0).unit() );
775 
776  extentList.AddSurface( polygon );
777  }
778 
779  if (r2 >= 0)
780  {
781  //
782  // Repeat, for transition piece
783  //
784  w2 = G4ThreeVector( r2*cosPhi, r2*sinPhi, z0 );
785  transform.ApplyPointTransform( w2 );
786 
787  polygon.ClearAllVertices();
788 
789  polygon.AddVertexInOrder( v2 );
790  polygon.AddVertexInOrder( v0 );
791  polygon.AddVertexInOrder( w0 );
792  polygon.AddVertexInOrder( w2 );
793 
794  if (polygon.PartialClip( voxelLimit, axis ))
795  {
796  polygon.SetNormal( deltaV.cross(v0-v2).unit() );
797 
798  extentList.AddSurface( polygon );
799  }
800 
801  v2 = w2;
802  }
803 
804  //
805  // Next vertex
806  //
807  v0 = w0;
808  v1 = w1;
809  } while( --numPhi > 0 );
810 
811  //
812  // We are almost done. But, it is important that we leave no
813  // gaps in the surface of our solid. By using rFudge, however,
814  // we've done exactly that, if we have a phi segment.
815  // Add two additional faces if necessary
816  //
817  if (phiIsOpen && rNorm > DBL_MIN)
818  {
819  cosPhi = std::cos(startPhi);
820  sinPhi = std::sin(startPhi);
821 
822  G4ThreeVector a0( r[0]*cosPhi, r[0]*sinPhi, z[0] ),
823  a1( r[1]*cosPhi, r[1]*sinPhi, z[1] ),
824  b0( r0*cosPhi, r0*sinPhi, z[0] ),
825  b1( r1*cosPhi, r1*sinPhi, z[1] );
826 
827  transform.ApplyPointTransform( a0 );
828  transform.ApplyPointTransform( a1 );
829  transform.ApplyPointTransform( b0 );
830  transform.ApplyPointTransform( b1 );
831 
832  polygon.ClearAllVertices();
833 
834  polygon.AddVertexInOrder( a0 );
835  polygon.AddVertexInOrder( a1 );
836  polygon.AddVertexInOrder( b0 );
837  polygon.AddVertexInOrder( b1 );
838 
839  if (polygon.PartialClip( voxelLimit , axis))
840  {
841  G4ThreeVector normal( sinPhi, -cosPhi, 0 );
842  polygon.SetNormal( transform.TransformAxis( normal ) );
843 
844  extentList.AddSurface( polygon );
845  }
846 
847  cosPhi = std::cos(startPhi+deltaPhi);
848  sinPhi = std::sin(startPhi+deltaPhi);
849 
850  a0 = G4ThreeVector( r[0]*cosPhi, r[0]*sinPhi, z[0] ),
851  a1 = G4ThreeVector( r[1]*cosPhi, r[1]*sinPhi, z[1] ),
852  b0 = G4ThreeVector( r0*cosPhi, r0*sinPhi, z[0] ),
853  b1 = G4ThreeVector( r1*cosPhi, r1*sinPhi, z[1] );
854  transform.ApplyPointTransform( a0 );
855  transform.ApplyPointTransform( a1 );
856  transform.ApplyPointTransform( b0 );
857  transform.ApplyPointTransform( b1 );
858 
859  polygon.ClearAllVertices();
860 
861  polygon.AddVertexInOrder( a0 );
862  polygon.AddVertexInOrder( a1 );
863  polygon.AddVertexInOrder( b0 );
864  polygon.AddVertexInOrder( b1 );
865 
866  if (polygon.PartialClip( voxelLimit, axis ))
867  {
868  G4ThreeVector normal( -sinPhi, cosPhi, 0 );
869  polygon.SetNormal( transform.TransformAxis( normal ) );
870 
871  extentList.AddSurface( polygon );
872  }
873  }
874 
875  return;
876 }
877 
878 
879 //
880 // GetPhi
881 //
882 // Calculate Phi for a given 3-vector (point), if not already cached for the
883 // same point, in the attempt to avoid consecutive computation of the same
884 // quantity
885 //
887 {
888  G4double val=0.;
889 
890  if (G4MT_pcphi.first != p)
891  {
892  val = p.phi();
893  G4MT_pcphi.first = p;
894  G4MT_pcphi.second = val;
895  }
896  else
897  {
898  val = G4MT_pcphi.second;
899  }
900  return val;
901 }
902 
903 
904 //
905 // DistanceAway
906 //
907 // Calculate distance of a point from our conical surface, including the effect
908 // of any phi segmentation
909 //
910 // Arguments:
911 // p - (in) Point to check
912 // opposite - (in) If true, check opposite hemisphere (see below)
913 // distOutside - (out) Additional distance outside the edges of the surface
914 // edgeRZnorm - (out) if negative, point is inside
915 //
916 // return value = distance from the conical plane, if extrapolated beyond edges,
917 // signed by whether the point is in inside or outside the shape
918 //
919 // Notes:
920 // * There are two answers, depending on which hemisphere is considered.
921 //
923  G4bool opposite,
924  G4double &distOutside2,
925  G4double *edgeRZnorm )
926 {
927  //
928  // Convert our point to r and z
929  //
930  G4double rx = p.perp(), zx = p.z();
931 
932  //
933  // Change sign of r if opposite says we should
934  //
935  if (opposite) rx = -rx;
936 
937  //
938  // Calculate return value
939  //
940  G4double deltaR = rx - r[0], deltaZ = zx - z[0];
941  G4double answer = deltaR*rNorm + deltaZ*zNorm;
942 
943  //
944  // Are we off the surface in r,z space?
945  //
946  G4double q = deltaR*rS + deltaZ*zS;
947  if (q < 0)
948  {
949  distOutside2 = q*q;
950  if (edgeRZnorm) *edgeRZnorm = deltaR*rNormEdge[0] + deltaZ*zNormEdge[0];
951  }
952  else if (q > length)
953  {
954  distOutside2 = sqr( q-length );
955  if (edgeRZnorm)
956  {
957  deltaR = rx - r[1];
958  deltaZ = zx - z[1];
959  *edgeRZnorm = deltaR*rNormEdge[1] + deltaZ*zNormEdge[1];
960  }
961  }
962  else
963  {
964  distOutside2 = 0;
965  if (edgeRZnorm) *edgeRZnorm = answer;
966  }
967 
968  if (phiIsOpen)
969  {
970  //
971  // Finally, check phi
972  //
973  G4double phi = GetPhi(p);
974  while( phi < startPhi ) phi += twopi;
975 
976  if (phi > startPhi+deltaPhi)
977  {
978  //
979  // Oops. Are we closer to the start phi or end phi?
980  //
982  while( phi > startPhi ) phi -= twopi;
983  G4double d2 = startPhi-phi;
984 
985  if (d2 < d1) d1 = d2;
986 
987  //
988  // Add result to our distance
989  //
990  G4double dist = d1*rx;
991 
992  distOutside2 += dist*dist;
993  if (edgeRZnorm)
994  {
995  *edgeRZnorm = std::max(std::fabs(*edgeRZnorm),std::fabs(dist));
996  }
997  }
998  }
999 
1000  return answer;
1001 }
1002 
1003 
1004 //
1005 // DistanceAway
1006 //
1007 // Special version of DistanceAway for Inside.
1008 // Opposite parameter is not used, instead use sign of rx for choosing the side
1009 //
1011  G4double &distOutside2,
1012  G4double *edgeRZnorm )
1013 {
1014  //
1015  // Convert our point to r and z
1016  //
1017  G4double rx = p.perp(), zx = p.z();
1018 
1019  //
1020  // Change sign of r if we should
1021  //
1022  G4int part = 1;
1023  if (rx < 0) part = -1;
1024 
1025  //
1026  // Calculate return value
1027  //
1028  G4double deltaR = rx - r[0]*part, deltaZ = zx - z[0];
1029  G4double answer = deltaR*rNorm*part + deltaZ*zNorm;
1030 
1031  //
1032  // Are we off the surface in r,z space?
1033  //
1034  G4double q = deltaR*rS*part + deltaZ*zS;
1035  if (q < 0)
1036  {
1037  distOutside2 = q*q;
1038  if (edgeRZnorm)
1039  {
1040  *edgeRZnorm = deltaR*rNormEdge[0]*part + deltaZ*zNormEdge[0];
1041  }
1042  }
1043  else if (q > length)
1044  {
1045  distOutside2 = sqr( q-length );
1046  if (edgeRZnorm)
1047  {
1048  deltaR = rx - r[1]*part;
1049  deltaZ = zx - z[1];
1050  *edgeRZnorm = deltaR*rNormEdge[1]*part + deltaZ*zNormEdge[1];
1051  }
1052  }
1053  else
1054  {
1055  distOutside2 = 0;
1056  if (edgeRZnorm) *edgeRZnorm = answer;
1057  }
1058 
1059  if (phiIsOpen)
1060  {
1061  //
1062  // Finally, check phi
1063  //
1064  G4double phi = GetPhi(p);
1065  while( phi < startPhi ) phi += twopi;
1066 
1067  if (phi > startPhi+deltaPhi)
1068  {
1069  //
1070  // Oops. Are we closer to the start phi or end phi?
1071  //
1072  G4double d1 = phi-startPhi-deltaPhi;
1073  while( phi > startPhi ) phi -= twopi;
1074  G4double d2 = startPhi-phi;
1075 
1076  if (d2 < d1) d1 = d2;
1077 
1078  //
1079  // Add result to our distance
1080  //
1081  G4double dist = d1*rx*part;
1082 
1083  distOutside2 += dist*dist;
1084  if (edgeRZnorm)
1085  {
1086  *edgeRZnorm = std::max(std::fabs(*edgeRZnorm),std::fabs(dist));
1087  }
1088  }
1089  }
1090 
1091  return answer;
1092 }
1093 
1094 
1095 //
1096 // PointOnCone
1097 //
1098 // Decide if a point is on a cone and return normal if it is
1099 //
1101  G4double normSign,
1102  const G4ThreeVector &p,
1103  const G4ThreeVector &v,
1105 {
1106  G4double rx = hit.perp();
1107  //
1108  // Check radial/z extent, as appropriate
1109  //
1110  if (!cone->HitOn( rx, hit.z() )) return false;
1111 
1112  if (phiIsOpen)
1113  {
1114  G4double phiTolerant = 2.0*kCarTolerance/(rx+kCarTolerance);
1115  //
1116  // Check phi segment. Here we have to be careful
1117  // to use the standard method consistent with
1118  // PolyPhiFace. See PolyPhiFace::InsideEdgesExact
1119  //
1120  G4double phi = GetPhi(hit);
1121  while( phi < startPhi-phiTolerant ) phi += twopi;
1122 
1123  if (phi > startPhi+deltaPhi+phiTolerant) return false;
1124 
1125  if (phi > startPhi+deltaPhi-phiTolerant)
1126  {
1127  //
1128  // Exact treatment
1129  //
1130  G4ThreeVector qx = p + v;
1131  G4ThreeVector qa = qx - corners[2],
1132  qb = qx - corners[3];
1133  G4ThreeVector qacb = qa.cross(qb);
1134 
1135  if (normSign*qacb.dot(v) < 0) return false;
1136  }
1137  else if (phi < phiTolerant)
1138  {
1139  G4ThreeVector qx = p + v;
1140  G4ThreeVector qa = qx - corners[1],
1141  qb = qx - corners[0];
1142  G4ThreeVector qacb = qa.cross(qb);
1143 
1144  if (normSign*qacb.dot(v) < 0) return false;
1145  }
1146  }
1147 
1148  //
1149  // We have a good hit! Calculate normal
1150  //
1151  if (rx < DBL_MIN)
1152  normal = G4ThreeVector( 0, 0, zNorm < 0 ? -1 : 1 );
1153  else
1154  normal = G4ThreeVector( rNorm*hit.x()/rx, rNorm*hit.y()/rx, zNorm );
1155  return true;
1156 }
1157 
1158 
1159 //
1160 // FindLineIntersect
1161 //
1162 // Decide the point at which two 2-dimensional lines intersect
1163 //
1164 // Equation of line: x = x1 + s*tx1
1165 // y = y1 + s*ty1
1166 //
1167 // It is assumed that the lines are *not* parallel
1168 //
1170  G4double tx1, G4double ty1,
1171  G4double x2, G4double y2,
1172  G4double tx2, G4double ty2,
1173  G4double &x, G4double &y )
1174 {
1175  //
1176  // The solution is a simple linear equation
1177  //
1178  G4double deter = tx1*ty2 - tx2*ty1;
1179 
1180  G4double s1 = ((x2-x1)*ty2 - tx2*(y2-y1))/deter;
1181  G4double s2 = ((x2-x1)*ty1 - tx1*(y2-y1))/deter;
1182 
1183  //
1184  // We want the answer to not depend on which order the
1185  // lines were specified. Take average.
1186  //
1187  x = 0.5*( x1+s1*tx1 + x2+s2*tx2 );
1188  y = 0.5*( y1+s1*ty1 + y2+s2*ty2 );
1189 }
1190 
1191 
1192 //
1193 // Calculate surface area for GetPointOnSurface()
1194 //
1196 {
1197  if(fSurfaceArea==0)
1198  {
1199  fSurfaceArea = (r[0]+r[1])* std::sqrt(sqr(r[0]-r[1])+sqr(z[0]-z[1]));
1200  fSurfaceArea *= 0.5*(deltaPhi);
1201  }
1202  return fSurfaceArea;
1203 }
1204 
1205 
1206 //
1207 // GetPointOnFace
1208 //
1210 {
1211  G4double x,y,zz;
1212  G4double rr,phi,dz,dr;
1213  dr=r[1]-r[0];dz=z[1]-z[0];
1215  rr=r[0]+dr*G4UniformRand();
1216 
1217  x=rr*std::cos(phi);
1218  y=rr*std::sin(phi);
1219 
1220  // PolyconeSide has a Ring Form
1221  //
1222  if (dz==0.)
1223  {
1224  zz=z[0];
1225  }
1226  else
1227  {
1228  if(dr==0.) // PolyconeSide has a Tube Form
1229  {
1230  zz = z[0]+dz*G4UniformRand();
1231  }
1232  else
1233  {
1234  zz = z[0]+(rr-r[0])*dz/dr;
1235  }
1236  }
1237 
1238  return G4ThreeVector(x,y,zz);
1239 }
const G4double a0
G4ThreeVector GetPointOnFace()
static G4GEOM_DLL G4PlSideManager subInstanceManager
static const G4double d1
static const G4double kInfinity
Definition: geomdefs.hh:42
G4double SurfaceArea()
#define G4MT_pcphi
CLHEP::Hep3Vector G4ThreeVector
G4double zNormEdge[2]
static const G4double a1
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, G4bool outgoing, G4double surfTolerance, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal, G4bool &isAllBehind)
G4double GetSurfaceTolerance() const
void SetNormal(const G4ThreeVector &newNormal)
virtual G4bool PartialClip(const G4VoxelLimits &voxelLimit, const EAxis IgnoreMe)
virtual void AddVertexInOrder(const G4ThreeVector vertex)
static const G4double tolerance
static const G4double cd
EInside Inside(const G4ThreeVector &p, G4double tolerance, G4double *bestDistance)
G4double a
Definition: TRTMaterials.hh:39
virtual void ClearAllVertices()
G4bool HitOn(const G4double r, const G4double z)
G4double Extent(const G4ThreeVector axis)
G4ThreeVector Normal(const G4ThreeVector &p, G4double *bestDistance)
G4int LineHitsCone(const G4ThreeVector &p, const G4ThreeVector &v, G4double *s1, G4double *s2)
#define G4ThreadLocal
Definition: tls.hh:89
const G4int kMinMeshSections
Definition: meshdefs.hh:45
G4int CreateSubInstance()
void CopyStuff(const G4PolyconeSide &source)
int G4int
Definition: G4Types.hh:78
static void FindLineIntersect(G4double x1, G4double y1, G4double tx1, G4double ty1, G4double x2, G4double y2, G4double tx2, G4double ty2, G4double &x, G4double &y)
G4double ZLo() const
G4ThreeVector * corners
virtual ~G4PolyconeSide()
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
G4double rNormEdge[2]
#define G4UniformRand()
Definition: Randomize.hh:95
void CalculateExtent(const EAxis axis, const G4VoxelLimits &voxelLimit, const G4AffineTransform &tranform, G4SolidExtentList &extentList)
G4PolyconeSide(const G4PolyconeSideRZ *prevRZ, const G4PolyconeSideRZ *tail, const G4PolyconeSideRZ *head, const G4PolyconeSideRZ *nextRZ, G4double phiStart, G4double deltaPhi, G4bool phiIsOpen, G4bool isAllBehind=false)
G4double ZHi() const
bool G4bool
Definition: G4Types.hh:79
void AddSurface(const G4ClippablePolygon &surface)
static const G4PlSideManager & GetSubInstanceManager()
G4double fSurfaceArea
G4IntersectingCone * cone
T max(const T t1, const T t2)
brief Return the largest of the two arguments
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
G4double GetPhi(const G4ThreeVector &p)
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
#define DBL_MIN
Definition: templates.hh:75
static const G4double b1
G4double DistanceAway(const G4ThreeVector &p, G4bool opposite, G4double &distOutside2, G4double *rzNorm=0)
G4bool PointOnCone(const G4ThreeVector &hit, G4double normSign, const G4ThreeVector &p, const G4ThreeVector &v, G4ThreeVector &normal)
G4double Distance(const G4ThreeVector &p, G4bool outgoing)
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
static const G4double d2
static const G4double b0
G4double kCarTolerance
G4PolyconeSide & operator=(const G4PolyconeSide &source)
const G4double kMeshAngleDefault
Definition: meshdefs.hh:42
static G4GeometryTolerance * GetInstance()
const G4double r0
void ApplyPointTransform(G4ThreeVector &vec) const
const G4int kMaxMeshSections
Definition: meshdefs.hh:46