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