Geant4 11.4.0
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 (UCSC), 1998
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 "G4QuickRand.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 {
98 deltaPhi += twopi;
99 }
100 while (startPhi < 0.0) // Loop checking, 13.08.2015, G.Cosmo
101 {
102 startPhi += twopi;
103 }
104
105 //
106 // Calculate corner coordinates
107 //
108 ncorners = 4;
109 corners = new G4ThreeVector[ncorners];
110
111 corners[0] = G4ThreeVector( tail->r*std::cos(startPhi),
112 tail->r*std::sin(startPhi), tail->z );
113 corners[1] = G4ThreeVector( head->r*std::cos(startPhi),
114 head->r*std::sin(startPhi), head->z );
115 corners[2] = G4ThreeVector( tail->r*std::cos(startPhi+deltaPhi),
116 tail->r*std::sin(startPhi+deltaPhi), tail->z );
117 corners[3] = G4ThreeVector( head->r*std::cos(startPhi+deltaPhi),
118 head->r*std::sin(startPhi+deltaPhi), head->z );
119 }
120 else
121 {
122 deltaPhi = twopi;
123 startPhi = 0.0;
124 }
125
126 allBehind = isAllBehind;
127
128 //
129 // Make our intersecting cone
130 //
131 cone = new G4IntersectingCone( r, z );
132
133 //
134 // Calculate vectors in r,z space
135 //
136 rS = r[1]-r[0]; zS = z[1]-z[0];
137 length = std::sqrt( rS*rS + zS*zS);
138 rS /= length; zS /= length;
139
140 rNorm = +zS;
141 zNorm = -rS;
142
143 G4double lAdj;
144
145 prevRS = r[0]-prevRZ->r;
146 prevZS = z[0]-prevRZ->z;
147 lAdj = std::sqrt( prevRS*prevRS + prevZS*prevZS );
148 prevRS /= lAdj;
149 prevZS /= lAdj;
150
151 rNormEdge[0] = rNorm + prevZS;
152 zNormEdge[0] = zNorm - prevRS;
153 lAdj = std::sqrt( rNormEdge[0]*rNormEdge[0] + zNormEdge[0]*zNormEdge[0] );
154 rNormEdge[0] /= lAdj;
155 zNormEdge[0] /= lAdj;
156
157 nextRS = nextRZ->r-r[1];
158 nextZS = nextRZ->z-z[1];
159 lAdj = std::sqrt( nextRS*nextRS + nextZS*nextZS );
160 nextRS /= lAdj;
161 nextZS /= lAdj;
162
163 rNormEdge[1] = rNorm + nextZS;
164 zNormEdge[1] = zNorm - nextRS;
165 lAdj = std::sqrt( rNormEdge[1]*rNormEdge[1] + zNormEdge[1]*zNormEdge[1] );
166 rNormEdge[1] /= lAdj;
167 zNormEdge[1] /= lAdj;
168}
169
170// Fake default constructor - sets only member data and allocates memory
171// for usage restricted to object persistency.
172//
174 : startPhi(0.), deltaPhi(0.),
175 rNorm(0.), zNorm(0.), rS(0.), zS(0.), length(0.),
176 prevRS(0.), prevZS(0.), nextRS(0.), nextZS(0.),
177 kCarTolerance(0.), instanceID(0)
178{
179 r[0] = r[1] = 0.;
180 z[0] = z[1] = 0.;
181 rNormEdge[0]= rNormEdge[1] = 0.;
182 zNormEdge[0]= zNormEdge[1] = 0.;
183}
184
185// Destructor
186//
188{
189 delete cone;
190 if (phiIsOpen) { delete [] corners; }
191}
192
193// Copy constructor
194//
196{
197 instanceID = subInstanceManager.CreateSubInstance();
198
199 CopyStuff( source );
200}
201
202// Assignment operator
203//
205{
206 if (this == &source) { return *this; }
207
208 delete cone;
209 if (phiIsOpen) { delete [] corners; }
210
211 CopyStuff( source );
212
213 return *this;
214}
215
216// CopyStuff
217//
218void G4PolyconeSide::CopyStuff( const G4PolyconeSide& source )
219{
220 r[0] = source.r[0];
221 r[1] = source.r[1];
222 z[0] = source.z[0];
223 z[1] = source.z[1];
224
225 startPhi = source.startPhi;
226 deltaPhi = source.deltaPhi;
227 phiIsOpen = source.phiIsOpen;
228 allBehind = source.allBehind;
229
230 kCarTolerance = source.kCarTolerance;
231 fSurfaceArea = source.fSurfaceArea;
232
233 cone = new G4IntersectingCone( *source.cone );
234
235 rNorm = source.rNorm;
236 zNorm = source.zNorm;
237 rS = source.rS;
238 zS = source.zS;
239 length = source.length;
240 prevRS = source.prevRS;
241 prevZS = source.prevZS;
242 nextRS = source.nextRS;
243 nextZS = source.nextZS;
244
245 rNormEdge[0] = source.rNormEdge[0];
246 rNormEdge[1] = source.rNormEdge[1];
247 zNormEdge[0] = source.zNormEdge[0];
248 zNormEdge[1] = source.zNormEdge[1];
249
250 if (phiIsOpen)
251 {
252 ncorners = 4;
253 corners = new G4ThreeVector[ncorners];
254
255 corners[0] = source.corners[0];
256 corners[1] = source.corners[1];
257 corners[2] = source.corners[2];
258 corners[3] = source.corners[3];
259 }
260}
261
262// Intersect
263//
265 const G4ThreeVector& v,
266 G4bool outgoing,
267 G4double surfTolerance,
268 G4double& distance,
269 G4double& distFromSurface,
270 G4ThreeVector& normal,
271 G4bool& isAllBehind )
272{
273 G4double s1=0., s2=0.;
274 G4double normSign = outgoing ? +1 : -1;
275
276 isAllBehind = allBehind;
277
278 //
279 // Check for two possible intersections
280 //
281 G4int nside = cone->LineHitsCone( p, v, &s1, &s2 );
282 if (nside == 0) { return false; }
283
284 //
285 // Check the first side first, since it is (supposed to be) closest
286 //
287 G4ThreeVector hit = p + s1*v;
288
289 if (PointOnCone( hit, normSign, p, v, normal ))
290 {
291 //
292 // Good intersection! What about the normal?
293 //
294 if (normSign*v.dot(normal) > 0)
295 {
296 //
297 // We have a valid intersection, but it could very easily
298 // be behind the point. To decide if we tolerate this,
299 // we have to see if the point p is on the surface near
300 // the intersecting point.
301 //
302 // What does it mean exactly for the point p to be "near"
303 // the intersection? It means that if we draw a line from
304 // p to the hit, the line remains entirely within the
305 // tolerance bounds of the cone. To test this, we can
306 // ask if the normal is correct near p.
307 //
308 G4double pr = p.perp();
309 if (pr < DBL_MIN) { pr = DBL_MIN; }
310 G4ThreeVector pNormal( rNorm*p.x()/pr, rNorm*p.y()/pr, zNorm );
311 if (normSign*v.dot(pNormal) > 0)
312 {
313 //
314 // p and intersection in same hemisphere
315 //
316 G4double distOutside2;
317 distFromSurface = -normSign*DistanceAway( p, false, distOutside2 );
318 if (distOutside2 < surfTolerance*surfTolerance)
319 {
320 if (distFromSurface > -surfTolerance)
321 {
322 //
323 // We are just inside or away from the
324 // surface. Accept *any* value of distance.
325 //
326 distance = s1;
327 return true;
328 }
329 }
330 }
331 else
332 {
333 distFromSurface = s1;
334 }
335
336 //
337 // Accept positive distances
338 //
339 if (s1 > 0)
340 {
341 distance = s1;
342 return true;
343 }
344 }
345 }
346
347 if (nside==1) { return false;
348}
349
350 //
351 // Well, try the second hit
352 //
353 hit = p + s2*v;
354
355 if (PointOnCone( hit, normSign, p, v, normal ))
356 {
357 //
358 // Good intersection! What about the normal?
359 //
360 if (normSign*v.dot(normal) > 0)
361 {
362 G4double pr = p.perp();
363 if (pr < DBL_MIN) { pr = DBL_MIN; }
364 G4ThreeVector pNormal( rNorm*p.x()/pr, rNorm*p.y()/pr, zNorm );
365 if (normSign*v.dot(pNormal) > 0)
366 {
367 G4double distOutside2;
368 distFromSurface = -normSign*DistanceAway( p, false, distOutside2 );
369 if (distOutside2 < surfTolerance*surfTolerance)
370 {
371 if (distFromSurface > -surfTolerance)
372 {
373 distance = s2;
374 return true;
375 }
376 }
377 }
378 else
379 {
380 distFromSurface = s2;
381 }
382
383 if (s2 > 0)
384 {
385 distance = s2;
386 return true;
387 }
388 }
389 }
390
391 //
392 // Better luck next time
393 //
394 return false;
395}
396
397// Distance
398//
400{
401 G4double normSign = outgoing ? -1 : +1;
402 G4double distFrom, distOut2;
403
404 //
405 // We have two tries for each hemisphere. Try the closest first.
406 //
407 distFrom = normSign*DistanceAway( p, false, distOut2 );
408 if (distFrom > -0.5*kCarTolerance )
409 {
410 //
411 // Good answer
412 //
413 if (distOut2 > 0)
414 {
415 return std::sqrt( distFrom*distFrom + distOut2 );
416 }
417 return std::fabs(distFrom);
418 }
419
420 //
421 // Try second side.
422 //
423 distFrom = normSign*DistanceAway( p, true, distOut2 );
424 if (distFrom > -0.5*kCarTolerance)
425 {
426 if (distOut2 > 0)
427 {
428 return std::sqrt( distFrom*distFrom + distOut2 );
429 }
430 return std::fabs(distFrom);
431 }
432
433 return kInfinity;
434}
435
436// Inside
437//
439 G4double tolerance,
440 G4double* bestDistance )
441{
442 G4double distFrom, distOut2, dist2;
443 G4double edgeRZnorm;
444
445 distFrom = DistanceAway( p, distOut2, &edgeRZnorm );
446 dist2 = distFrom*distFrom + distOut2;
447
448 *bestDistance = std::sqrt( dist2);
449
450 // Okay then, inside or out?
451 //
452 if ( (std::fabs(edgeRZnorm) < tolerance)
453 && (distOut2< tolerance*tolerance) )
454 {
455 return kSurface;
456 }
457 if (edgeRZnorm < 0)
458 {
459 return kInside;
460 }
461
462 return kOutside;
463}
464
465// Normal
466//
468 G4double* bestDistance )
469{
470 if (p == G4ThreeVector(0.,0.,0.)) { return p; }
471
472 G4double dFrom, dOut2;
473
474 dFrom = DistanceAway( p, false, dOut2 );
475
476 *bestDistance = std::sqrt( dFrom*dFrom + dOut2 );
477
478 G4double rds = p.perp();
479 if (rds!=0.) { return {rNorm*p.x()/rds,rNorm*p.y()/rds,zNorm}; }
480 return G4ThreeVector( 0.,0., zNorm ).unit();
481}
482
483// Extent
484//
486{
487 if (axis.perp2() < DBL_MIN)
488 {
489 //
490 // Special case
491 //
492 return axis.z() < 0 ? -cone->ZLo() : cone->ZHi();
493 }
494
495 //
496 // Is the axis pointing inside our phi gap?
497 //
498 if (phiIsOpen)
499 {
500 G4double phi = GetPhi(axis);
501 while( phi < startPhi ) // Loop checking, 13.08.2015, G.Cosmo
502 {
503 phi += twopi;
504 }
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// CalculateExtent
546//
547// See notes in G4VCSGface
548//
550 const G4VoxelLimits& voxelLimit,
551 const G4AffineTransform& transform,
552 G4SolidExtentList& extentList )
553{
554 G4ClippablePolygon polygon;
555
556 //
557 // Here we will approximate (ala G4Cons) and divide our conical section
558 // into segments, like G4Polyhedra. When doing so, the radius
559 // is extented far enough such that the segments always lie
560 // just outside the surface of the conical section we are
561 // approximating.
562 //
563
564 //
565 // Choose phi size of our segment(s) based on constants as
566 // defined in meshdefs.hh
567 //
568 G4int numPhi = (G4int)(deltaPhi/kMeshAngleDefault) + 1;
569 if (numPhi < kMinMeshSections)
570 {
571 numPhi = kMinMeshSections;
572 }
573 else if (numPhi > kMaxMeshSections)
574 {
575 numPhi = kMaxMeshSections;
576 }
577
578 G4double sigPhi = deltaPhi/numPhi;
579
580 //
581 // Determine radius factor to keep segments outside
582 //
583 G4double rFudge = 1.0/std::cos(0.5*sigPhi);
584
585 //
586 // Decide which radius to use on each end of the side,
587 // and whether a transition mesh is required
588 //
589 // {r0,z0} - Beginning of this side
590 // {r1,z1} - Ending of this side
591 // {r2,z0} - Beginning of transition piece connecting previous
592 // side (and ends at beginning of this side)
593 //
594 // So, order is 2 --> 0 --> 1.
595 // -------
596 //
597 // r2 < 0 indicates that no transition piece is required
598 //
599 G4double r0, r1, r2, z0, z1;
600
601 r2 = -1; // By default: no transition piece
602
603 if (rNorm < -DBL_MIN)
604 {
605 //
606 // This side faces *inward*, and so our mesh has
607 // the same radius
608 //
609 r1 = r[1];
610 z1 = z[1];
611 z0 = z[0];
612 r0 = r[0];
613
614 r2 = -1;
615
616 if (prevZS > DBL_MIN)
617 {
618 //
619 // The previous side is facing outwards
620 //
621 if ( prevRS*zS - prevZS*rS > 0 )
622 {
623 //
624 // Transition was convex: build transition piece
625 //
626 if (r[0] > DBL_MIN) { r2 = r[0]*rFudge; }
627 }
628 else
629 {
630 //
631 // Transition was concave: short this side
632 //
633 FindLineIntersect( z0, r0, zS, rS,
634 z0, r0*rFudge, prevZS, prevRS*rFudge, z0, r0 );
635 }
636 }
637
638 if ( nextZS > DBL_MIN && (rS*nextZS - zS*nextRS < 0) )
639 {
640 //
641 // The next side is facing outwards, forming a
642 // concave transition: short this side
643 //
644 FindLineIntersect( z1, r1, zS, rS,
645 z1, r1*rFudge, nextZS, nextRS*rFudge, z1, r1 );
646 }
647 }
648 else if (rNorm > DBL_MIN)
649 {
650 //
651 // This side faces *outward* and is given a boost to
652 // it radius
653 //
654 r0 = r[0]*rFudge;
655 z0 = z[0];
656 r1 = r[1]*rFudge;
657 z1 = z[1];
658
659 if (prevZS < -DBL_MIN)
660 {
661 //
662 // The previous side is facing inwards
663 //
664 if ( prevRS*zS - prevZS*rS > 0 )
665 {
666 //
667 // Transition was convex: build transition piece
668 //
669 if (r[0] > DBL_MIN) { r2 = r[0]; }
670 }
671 else
672 {
673 //
674 // Transition was concave: short this side
675 //
676 FindLineIntersect( z0, r0, zS, rS*rFudge,
677 z0, r[0], prevZS, prevRS, z0, r0 );
678 }
679 }
680
681 if ( nextZS < -DBL_MIN && (rS*nextZS - zS*nextRS < 0) )
682 {
683 //
684 // The next side is facing inwards, forming a
685 // concave transition: short this side
686 //
687 FindLineIntersect( z1, r1, zS, rS*rFudge,
688 z1, r[1], nextZS, nextRS, z1, r1 );
689 }
690 }
691 else
692 {
693 //
694 // This side is perpendicular to the z axis (is a disk)
695 //
696 // Whether or not r0 needs a rFudge factor depends
697 // on the normal of the previous edge. Similar with r1
698 // and the next edge. No transition piece is required.
699 //
700 r0 = r[0];
701 r1 = r[1];
702 z0 = z[0];
703 z1 = z[1];
704
705 if (prevZS > DBL_MIN) { r0 *= rFudge; }
706 if (nextZS > DBL_MIN) { r1 *= rFudge; }
707 }
708
709 //
710 // Loop
711 //
712 G4double phi = startPhi,
713 cosPhi = std::cos(phi),
714 sinPhi = std::sin(phi);
715
716 G4ThreeVector v0( r0*cosPhi, r0*sinPhi, z0 ),
717 v1( r1*cosPhi, r1*sinPhi, z1 ),
718 v2, w0, w1, w2;
719 transform.ApplyPointTransform( v0 );
720 transform.ApplyPointTransform( v1 );
721
722 if (r2 >= 0)
723 {
724 v2 = G4ThreeVector( r2*cosPhi, r2*sinPhi, z0 );
725 transform.ApplyPointTransform( v2 );
726 }
727
728 do // Loop checking, 13.08.2015, G.Cosmo
729 {
730 phi += sigPhi;
731 if (numPhi == 1) { phi = startPhi+deltaPhi; } // Try to avoid roundoff
732 cosPhi = std::cos(phi),
733 sinPhi = std::sin(phi);
734
735 w0 = G4ThreeVector( r0*cosPhi, r0*sinPhi, z0 );
736 w1 = G4ThreeVector( r1*cosPhi, r1*sinPhi, z1 );
737 transform.ApplyPointTransform( w0 );
738 transform.ApplyPointTransform( w1 );
739
740 G4ThreeVector deltaV = r0 > r1 ? w0-v0 : w1-v1;
741
742 //
743 // Build polygon, taking special care to keep the vertices
744 // in order
745 //
746 polygon.ClearAllVertices();
747
748 polygon.AddVertexInOrder( v0 );
749 polygon.AddVertexInOrder( v1 );
750 polygon.AddVertexInOrder( w1 );
751 polygon.AddVertexInOrder( w0 );
752
753 //
754 // Get extent
755 //
756 if (polygon.PartialClip( voxelLimit, axis ))
757 {
758 //
759 // Get dot product of normal with target axis
760 //
761 polygon.SetNormal( deltaV.cross(v1-v0).unit() );
762
763 extentList.AddSurface( polygon );
764 }
765
766 if (r2 >= 0)
767 {
768 //
769 // Repeat, for transition piece
770 //
771 w2 = G4ThreeVector( r2*cosPhi, r2*sinPhi, z0 );
772 transform.ApplyPointTransform( w2 );
773
774 polygon.ClearAllVertices();
775
776 polygon.AddVertexInOrder( v2 );
777 polygon.AddVertexInOrder( v0 );
778 polygon.AddVertexInOrder( w0 );
779 polygon.AddVertexInOrder( w2 );
780
781 if (polygon.PartialClip( voxelLimit, axis ))
782 {
783 polygon.SetNormal( deltaV.cross(v0-v2).unit() );
784
785 extentList.AddSurface( polygon );
786 }
787
788 v2 = w2;
789 }
790
791 //
792 // Next vertex
793 //
794 v0 = w0;
795 v1 = w1;
796 } while( --numPhi > 0 );
797
798 //
799 // We are almost done. But, it is important that we leave no
800 // gaps in the surface of our solid. By using rFudge, however,
801 // we've done exactly that, if we have a phi segment.
802 // Add two additional faces if necessary
803 //
804 if (phiIsOpen && rNorm > DBL_MIN)
805 {
806 cosPhi = std::cos(startPhi);
807 sinPhi = std::sin(startPhi);
808
809 G4ThreeVector a0( r[0]*cosPhi, r[0]*sinPhi, z[0] ),
810 a1( r[1]*cosPhi, r[1]*sinPhi, z[1] ),
811 b0( r0*cosPhi, r0*sinPhi, z[0] ),
812 b1( r1*cosPhi, r1*sinPhi, z[1] );
813
814 transform.ApplyPointTransform( a0 );
815 transform.ApplyPointTransform( a1 );
816 transform.ApplyPointTransform( b0 );
817 transform.ApplyPointTransform( b1 );
818
819 polygon.ClearAllVertices();
820
821 polygon.AddVertexInOrder( a0 );
822 polygon.AddVertexInOrder( a1 );
823 polygon.AddVertexInOrder( b0 );
824 polygon.AddVertexInOrder( b1 );
825
826 if (polygon.PartialClip( voxelLimit , axis))
827 {
828 G4ThreeVector normal( sinPhi, -cosPhi, 0 );
829 polygon.SetNormal( transform.TransformAxis( normal ) );
830
831 extentList.AddSurface( polygon );
832 }
833
834 cosPhi = std::cos(startPhi+deltaPhi);
835 sinPhi = std::sin(startPhi+deltaPhi);
836
837 a0 = G4ThreeVector( r[0]*cosPhi, r[0]*sinPhi, z[0] ),
838 a1 = G4ThreeVector( r[1]*cosPhi, r[1]*sinPhi, z[1] ),
839 b0 = G4ThreeVector( r0*cosPhi, r0*sinPhi, z[0] ),
840 b1 = G4ThreeVector( r1*cosPhi, r1*sinPhi, z[1] );
841 transform.ApplyPointTransform( a0 );
842 transform.ApplyPointTransform( a1 );
843 transform.ApplyPointTransform( b0 );
844 transform.ApplyPointTransform( b1 );
845
846 polygon.ClearAllVertices();
847
848 polygon.AddVertexInOrder( a0 );
849 polygon.AddVertexInOrder( a1 );
850 polygon.AddVertexInOrder( b0 );
851 polygon.AddVertexInOrder( b1 );
852
853 if (polygon.PartialClip( voxelLimit, axis ))
854 {
855 G4ThreeVector normal( -sinPhi, cosPhi, 0 );
856 polygon.SetNormal( transform.TransformAxis( normal ) );
857
858 extentList.AddSurface( polygon );
859 }
860 }
861
862 return;
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//
871G4double G4PolyconeSide::GetPhi( const G4ThreeVector& p )
872{
873 G4double val=0.;
875
876 if (vphi != p)
877 {
878 val = p.phi();
879 G4MT_pcphix = p.x(); G4MT_pcphiy = p.y(); G4MT_pcphiz = p.z();
880 G4MT_pcphik = val;
881 }
882 else
883 {
884 val = G4MT_pcphik;
885 }
886 return val;
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//
906G4double G4PolyconeSide::DistanceAway( const G4ThreeVector& p,
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 != nullptr)
935 {
936 *edgeRZnorm = deltaR*rNormEdge[0] + deltaZ*zNormEdge[0];
937 }
938 }
939 else if (q > length)
940 {
941 distOutside2 = sqr( q-length );
942 if (edgeRZnorm != nullptr)
943 {
944 deltaR = rx - r[1];
945 deltaZ = zx - z[1];
946 *edgeRZnorm = deltaR*rNormEdge[1] + deltaZ*zNormEdge[1];
947 }
948 }
949 else
950 {
951 distOutside2 = 0.;
952 if (edgeRZnorm != nullptr) { *edgeRZnorm = answer; }
953 }
954
955 if (phiIsOpen)
956 {
957 //
958 // Finally, check phi
959 //
960 G4double phi = GetPhi(p);
961 while( phi < startPhi ) // Loop checking, 13.08.2015, G.Cosmo
962 {
963 phi += twopi;
964 }
965
966 if (phi > startPhi+deltaPhi)
967 {
968 //
969 // Oops. Are we closer to the start phi or end phi?
970 //
971 G4double d1 = phi-startPhi-deltaPhi;
972 while( phi > startPhi ) // Loop checking, 13.08.2015, G.Cosmo
973 {
974 phi -= twopi;
975 }
976 G4double d2 = startPhi-phi;
977
978 if (d2 < d1) { d1 = d2; }
979
980 //
981 // Add result to our distance
982 //
983 G4double dist = d1*rx;
984
985 distOutside2 += dist*dist;
986 if (edgeRZnorm != nullptr)
987 {
988 *edgeRZnorm = std::max(std::fabs(*edgeRZnorm),std::fabs(dist));
989 }
990 }
991 }
992
993 return answer;
994}
995
996// DistanceAway
997//
998// Special version of DistanceAway for Inside.
999// Opposite parameter is not used, instead use sign of rx for choosing the side
1000//
1001G4double G4PolyconeSide::DistanceAway( const G4ThreeVector& p,
1002 G4double& distOutside2,
1003 G4double* edgeRZnorm )
1004{
1005 //
1006 // Convert our point to r and z
1007 //
1008 G4double rx = p.perp(), zx = p.z();
1009
1010 //
1011 // Change sign of r if we should
1012 //
1013 G4int part = 1;
1014 if (rx < 0) { part = -1; }
1015
1016 //
1017 // Calculate return value
1018 //
1019 G4double deltaR = rx - r[0]*part, deltaZ = zx - z[0];
1020 G4double answer = deltaR*rNorm*part + deltaZ*zNorm;
1021
1022 //
1023 // Are we off the surface in r,z space?
1024 //
1025 G4double q = deltaR*rS*part + deltaZ*zS;
1026 if (q < 0)
1027 {
1028 distOutside2 = q*q;
1029 if (edgeRZnorm != nullptr)
1030 {
1031 *edgeRZnorm = deltaR*rNormEdge[0]*part + deltaZ*zNormEdge[0];
1032 }
1033 }
1034 else if (q > length)
1035 {
1036 distOutside2 = sqr( q-length );
1037 if (edgeRZnorm != nullptr)
1038 {
1039 deltaR = rx - r[1]*part;
1040 deltaZ = zx - z[1];
1041 *edgeRZnorm = deltaR*rNormEdge[1]*part + deltaZ*zNormEdge[1];
1042 }
1043 }
1044 else
1045 {
1046 distOutside2 = 0.;
1047 if (edgeRZnorm != nullptr) { *edgeRZnorm = answer; }
1048 }
1049
1050 if (phiIsOpen)
1051 {
1052 //
1053 // Finally, check phi
1054 //
1055 G4double phi = GetPhi(p);
1056 while( phi < startPhi ) // Loop checking, 13.08.2015, G.Cosmo
1057 {
1058 phi += twopi;
1059 }
1060
1061 if (phi > startPhi+deltaPhi)
1062 {
1063 //
1064 // Oops. Are we closer to the start phi or end phi?
1065 //
1066 G4double d1 = phi-startPhi-deltaPhi;
1067 while( phi > startPhi ) // Loop checking, 13.08.2015, G.Cosmo
1068 {
1069 phi -= twopi;
1070 }
1071 G4double d2 = startPhi-phi;
1072
1073 if (d2 < d1) { d1 = d2; }
1074
1075 //
1076 // Add result to our distance
1077 //
1078 G4double dist = d1*rx*part;
1079
1080 distOutside2 += dist*dist;
1081 if (edgeRZnorm != nullptr)
1082 {
1083 *edgeRZnorm = std::max(std::fabs(*edgeRZnorm),std::fabs(dist));
1084 }
1085 }
1086 }
1087
1088 return answer;
1089}
1090
1091// PointOnCone
1092//
1093// Decide if a point is on a cone and return normal if it is
1094//
1095G4bool G4PolyconeSide::PointOnCone( const G4ThreeVector& hit,
1096 G4double normSign,
1097 const G4ThreeVector& p,
1098 const G4ThreeVector& v,
1099 G4ThreeVector& normal )
1100{
1101 G4double rx = hit.perp();
1102 //
1103 // Check radial/z extent, as appropriate
1104 //
1105 if (!cone->HitOn( rx, hit.z() )) { return false; }
1106
1107 if (phiIsOpen)
1108 {
1109 G4double phiTolerant = 2.0*kCarTolerance/(rx+kCarTolerance);
1110 //
1111 // Check phi segment. Here we have to be careful
1112 // to use the standard method consistent with
1113 // PolyPhiFace. See PolyPhiFace::InsideEdgesExact
1114 //
1115 G4double phi = GetPhi(hit);
1116 while( phi < startPhi-phiTolerant ) // Loop checking, 13.08.2015, G.Cosmo
1117 {
1118 phi += twopi;
1119 }
1120
1121 if (phi > startPhi+deltaPhi+phiTolerant) { return false; }
1122
1123 if (phi > startPhi+deltaPhi-phiTolerant)
1124 {
1125 //
1126 // Exact treatment
1127 //
1128 G4ThreeVector qx = p + v;
1129 G4ThreeVector qa = qx - corners[2],
1130 qb = qx - corners[3];
1131 G4ThreeVector qacb = qa.cross(qb);
1132
1133 if (normSign*qacb.dot(v) < 0) { return false; }
1134 }
1135 else if (phi < phiTolerant)
1136 {
1137 G4ThreeVector qx = p + v;
1138 G4ThreeVector qa = qx - corners[1],
1139 qb = qx - corners[0];
1140 G4ThreeVector qacb = qa.cross(qb);
1141
1142 if (normSign*qacb.dot(v) < 0) { return false; }
1143 }
1144 }
1145
1146 //
1147 // We have a good hit! Calculate normal
1148 //
1149 if (rx < DBL_MIN)
1150 {
1151 normal = G4ThreeVector( 0, 0, zNorm < 0 ? -1 : 1 );
1152 }
1153 else
1154 {
1155 normal = G4ThreeVector( rNorm*hit.x()/rx, rNorm*hit.y()/rx, zNorm );
1156 }
1157 return true;
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//
1169void G4PolyconeSide::FindLineIntersect( G4double x1, G4double y1,
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// Calculate surface area for GetPointOnSurface()
1192//
1194{
1195 if(fSurfaceArea==0.)
1196 {
1197 fSurfaceArea = (r[0]+r[1])* std::sqrt(sqr(r[0]-r[1])+sqr(z[0]-z[1]));
1198 fSurfaceArea *= 0.5*(deltaPhi);
1199 }
1200 return fSurfaceArea;
1201}
1202
1203// GetPointOnFace
1204//
1206{
1207 G4double x,y,zz;
1208 G4double rr,phi,dz,dr;
1209 dr=r[1]-r[0];dz=z[1]-z[0];
1210 phi=startPhi+deltaPhi*G4QuickRand();
1211 rr=r[0]+dr*G4QuickRand();
1212
1213 x=rr*std::cos(phi);
1214 y=rr*std::sin(phi);
1215
1216 // PolyconeSide has a Ring Form
1217 //
1218 if (dz==0.)
1219 {
1220 zz=z[0];
1221 }
1222 else
1223 {
1224 if(dr==0.) // PolyconeSide has a Tube Form
1225 {
1226 zz = z[0]+dz*G4QuickRand();
1227 }
1228 else
1229 {
1230 zz = z[0]+(rr-r[0])*dz/dr;
1231 }
1232 }
1233
1234 return {x,y,zz};
1235}
const G4double kCarTolerance
const G4double a0
#define G4MT_pcphiz
#define G4MT_pcphik
#define G4MT_pcphix
#define G4MT_pcphiy
G4GeomSplitter< G4PlSideData > G4PlSideManager
G4double G4QuickRand(uint32_t seed=0)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
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 perp() const
G4AffineTransform is a class for geometric affine transformations. It supports efficient arbitrary ro...
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
void ApplyPointTransform(G4ThreeVector &vec) const
G4ClippablePolygon in a utility class defining a polygon that can be clipped by a voxel.
G4bool PartialClip(const G4VoxelLimits &voxelLimit, const EAxis IgnoreMe)
void AddVertexInOrder(const G4ThreeVector &vertex)
void SetNormal(const G4ThreeVector &newNormal)
G4int CreateSubInstance()
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4IntersectingCone is a utility class used to calculate the intersection of an arbitrary line with a ...
G4PolyconeSide is a utility class implementing a face that represents one conical side of a polycone.
void CalculateExtent(const EAxis axis, const G4VoxelLimits &voxelLimit, const G4AffineTransform &tranform, G4SolidExtentList &extentList) override
~G4PolyconeSide() override
G4ThreeVector GetPointOnFace() override
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
EInside Inside(const G4ThreeVector &p, G4double tolerance, G4double *bestDistance) override
G4double Distance(const G4ThreeVector &p, G4bool outgoing) override
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)
G4SolidExtentList is utility class designed for calculating the extent of a CSG-like solid for a voxe...
void AddSurface(const G4ClippablePolygon &surface)
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
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 axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668
T sqr(const T &x)
Definition templates.hh:128
#define DBL_MIN
Definition templates.hh:54