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