Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PolyPhiFace.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 G4PolyPhiFace, the face that bounds a polycone or
27// polyhedra at its phi opening.
28//
29// Author: David C. Williams (UCSC), 1998
30// --------------------------------------------------------------------
31
32#include "G4PolyPhiFace.hh"
33#include "G4ClippablePolygon.hh"
34#include "G4ReduciblePolygon.hh"
35#include "G4AffineTransform.hh"
36#include "G4SolidExtentList.hh"
38
39#include "G4QuickRand.hh"
40#include "G4TwoVector.hh"
41
42// Constructor
43//
44// Points r,z should be supplied in clockwise order in r,z. For example:
45//
46// [1]---------[2] ^ R
47// | | |
48// | | +--> z
49// [0]---------[3]
50//
52 G4double phi,
53 G4double deltaPhi,
54 G4double phiOther )
55{
57
58 numEdges = rz->NumVertices();
59
60 rMin = rz->Amin();
61 rMax = rz->Amax();
62 zMin = rz->Bmin();
63 zMax = rz->Bmax();
64
65 //
66 // Is this the "starting" phi edge of the two?
67 //
68 G4bool start = (phiOther > phi);
69
70 //
71 // Build radial vector
72 //
73 radial = G4ThreeVector( std::cos(phi), std::sin(phi), 0.0 );
74
75 //
76 // Build normal
77 //
78 G4double zSign = start ? 1 : -1;
79 normal = G4ThreeVector( zSign*radial.y(), -zSign*radial.x(), 0 );
80
81 //
82 // Is allBehind?
83 //
84 allBehind = (zSign*(std::cos(phiOther)*radial.y() - std::sin(phiOther)*radial.x()) < 0);
85
86 //
87 // Adjacent edges
88 //
89 G4double midPhi = phi + (start ? +0.5 : -0.5)*deltaPhi;
90 G4double cosMid = std::cos(midPhi),
91 sinMid = std::sin(midPhi);
92 //
93 // Allocate corners
94 //
95 const std::size_t maxEdges = numEdges>0 ? numEdges : 1;
96 corners = new G4PolyPhiFaceVertex[maxEdges];
97 //
98 // Fill them
99 //
101
102 G4PolyPhiFaceVertex* corn = corners;
103 G4PolyPhiFaceVertex* helper = corners;
104
105 iterRZ.Begin();
106 do // Loop checking, 13.08.2015, G.Cosmo
107 {
108 corn->r = iterRZ.GetA();
109 corn->z = iterRZ.GetB();
110 corn->x = corn->r*radial.x();
111 corn->y = corn->r*radial.y();
112
113 // Add pointer on prev corner
114 //
115 if( corn == corners )
116 { corn->prev = corners+maxEdges-1;}
117 else
118 { corn->prev = helper; }
119
120 // Add pointer on next corner
121 //
122 if( corn < corners+maxEdges-1 )
123 { corn->next = corn+1;}
124 else
125 { corn->next = corners; }
126
127 helper = corn;
128 } while( ++corn, iterRZ.Next() );
129
130 //
131 // Allocate edges
132 //
133 edges = new G4PolyPhiFaceEdge[maxEdges];
134
135 //
136 // Fill them
137 //
138 G4double rFact = std::cos(0.5*deltaPhi);
139 G4double rFactNormalize = 1.0/std::sqrt(1.0+rFact*rFact);
140
141 G4PolyPhiFaceVertex* prev = corners+maxEdges-1,
142 * here = corners;
143 G4PolyPhiFaceEdge* edge = edges;
144 do // Loop checking, 13.08.2015, G.Cosmo
145 {
146 G4ThreeVector sideNorm;
147
148 edge->v0 = prev;
149 edge->v1 = here;
150
151 G4double dr = here->r - prev->r,
152 dz = here->z - prev->z;
153
154 edge->length = std::sqrt( dr*dr + dz*dz );
155
156 edge->tr = dr/edge->length;
157 edge->tz = dz/edge->length;
158
159 if ((here->r < DBL_MIN) && (prev->r < DBL_MIN))
160 {
161 //
162 // Sigh! Always exceptions!
163 // This edge runs at r==0, so its adjoing surface is not a
164 // PolyconeSide or PolyhedraSide, but the opposite PolyPhiFace.
165 //
166 G4double zSignOther = start ? -1 : 1;
167 sideNorm = G4ThreeVector( zSignOther*std::sin(phiOther),
168 -zSignOther*std::cos(phiOther), 0 );
169 }
170 else
171 {
172 sideNorm = G4ThreeVector( edge->tz*cosMid,
173 edge->tz*sinMid,
174 -edge->tr*rFact );
175 sideNorm *= rFactNormalize;
176 }
177 sideNorm += normal;
178
179 edge->norm3D = sideNorm.unit();
180 } while( edge++, prev=here, ++here < corners+maxEdges );
181
182 //
183 // Go back and fill in corner "normals"
184 //
185 G4PolyPhiFaceEdge* prevEdge = edges+maxEdges-1;
186 edge = edges;
187 do // Loop checking, 13.08.2015, G.Cosmo
188 {
189 //
190 // Calculate vertex 2D normals (on the phi surface)
191 //
192 G4double rPart = prevEdge->tr + edge->tr;
193 G4double zPart = prevEdge->tz + edge->tz;
194 G4double norm = std::sqrt( rPart*rPart + zPart*zPart );
195 G4double rNorm = +zPart/norm;
196 G4double zNorm = -rPart/norm;
197
198 edge->v0->rNorm = rNorm;
199 edge->v0->zNorm = zNorm;
200
201 //
202 // Calculate the 3D normals.
203 //
204 // Find the vector perpendicular to the z axis
205 // that defines the plane that contains the vertex normal
206 //
207 G4ThreeVector xyVector;
208
209 if (edge->v0->r < DBL_MIN)
210 {
211 //
212 // This is a vertex at r==0, which is a special
213 // case. The normal we will construct lays in the
214 // plane at the center of the phi opening.
215 //
216 // We also know that rNorm < 0
217 //
218 G4double zSignOther = start ? -1 : 1;
219 G4ThreeVector normalOther( zSignOther*std::sin(phiOther),
220 -zSignOther*std::cos(phiOther), 0 );
221
222 xyVector = - normal - normalOther;
223 }
224 else
225 {
226 //
227 // This is a vertex at r > 0. The plane
228 // is the average of the normal and the
229 // normal of the adjacent phi face
230 //
231 xyVector = G4ThreeVector( cosMid, sinMid, 0 );
232 if (rNorm < 0)
233 {
234 xyVector -= normal;
235 }
236 else
237 {
238 xyVector += normal;
239 }
240 }
241
242 //
243 // Combine it with the r/z direction from the face
244 //
245 edge->v0->norm3D = rNorm*xyVector.unit() + G4ThreeVector( 0, 0, zNorm );
246 } while( prevEdge=edge, ++edge < edges+maxEdges );
247
248 //
249 // Build point on surface
250 //
251 G4double rAve = 0.5*(rMax-rMin),
252 zAve = 0.5*(zMax-zMin);
253 surface = G4ThreeVector( rAve*radial.x(), rAve*radial.y(), zAve );
254}
255
256// Diagnose
257//
258// Throw an exception if something is found inconsistent with
259// the solid.
260//
261// For debugging purposes only
262//
264{
265 G4PolyPhiFaceVertex* corner = corners;
266 do // Loop checking, 13.08.2015, G.Cosmo
267 {
268 G4ThreeVector test(corner->x, corner->y, corner->z);
269 test -= 1E-6*corner->norm3D;
270
271 if (owner->Inside(test) != kInside)
272 {
273 G4Exception( "G4PolyPhiFace::Diagnose()", "GeomSolids0002",
274 FatalException, "Bad vertex normal found." );
275 }
276 } while( ++corner < corners+numEdges );
277}
278
279// Fake default constructor - sets only member data and allocates memory
280// for usage restricted to object persistency.
281//
283 : rMin(0.), rMax(0.), zMin(0.), zMax(0.), kCarTolerance(0.)
284{
285}
286
287
288//
289// Destructor
290//
292{
293 delete [] edges;
294 delete [] corners;
295}
296
297
298//
299// Copy constructor
300//
302{
303 CopyStuff( source );
304}
305
306
307//
308// Assignment operator
309//
311{
312 if (this == &source) { return *this; }
313
314 delete [] edges;
315 delete [] corners;
316
317 CopyStuff( source );
318
319 return *this;
320}
321
322
323//
324// CopyStuff (protected)
325//
326void G4PolyPhiFace::CopyStuff( const G4PolyPhiFace& source )
327{
328 //
329 // The simple stuff
330 //
331 numEdges = source.numEdges;
332 normal = source.normal;
333 radial = source.radial;
334 surface = source.surface;
335 rMin = source.rMin;
336 rMax = source.rMax;
337 zMin = source.zMin;
338 zMax = source.zMax;
339 allBehind = source.allBehind;
340 triangles = nullptr;
341
342 kCarTolerance = source.kCarTolerance;
343 fSurfaceArea = source.fSurfaceArea;
344
345 const std::size_t maxEdges = (numEdges > 0) ? numEdges : 1;
346
347 //
348 // Corner dynamic array
349 //
350 corners = new G4PolyPhiFaceVertex[maxEdges];
351 G4PolyPhiFaceVertex *corn = corners,
352 *sourceCorn = source.corners;
353 do // Loop checking, 13.08.2015, G.Cosmo
354 {
355 *corn = *sourceCorn;
356 } while( ++sourceCorn, ++corn < corners+maxEdges );
357
358 //
359 // Edge dynamic array
360 //
361 edges = new G4PolyPhiFaceEdge[maxEdges];
362
363 G4PolyPhiFaceVertex* prev = corners+maxEdges-1,
364 * here = corners;
365 G4PolyPhiFaceEdge* edge = edges,
366 * sourceEdge = source.edges;
367 do // Loop checking, 13.08.2015, G.Cosmo
368 {
369 *edge = *sourceEdge;
370 edge->v0 = prev;
371 edge->v1 = here;
372 } while( ++sourceEdge, ++edge, prev=here, ++here < corners+maxEdges );
373}
374
375// Intersect
376//
378 const G4ThreeVector& v,
379 G4bool outgoing,
380 G4double surfTolerance,
381 G4double& distance,
382 G4double& distFromSurface,
383 G4ThreeVector& aNormal,
384 G4bool& isAllBehind )
385{
386 G4double normSign = outgoing ? +1 : -1;
387
388 //
389 // These don't change
390 //
391 isAllBehind = allBehind;
392 aNormal = normal;
393
394 //
395 // Correct normal? Here we have straight sides, and can safely ignore
396 // intersections where the dot product with the normal is zero.
397 //
398 G4double dotProd = normSign*normal.dot(v);
399
400 if (dotProd <= 0) { return false; }
401
402 //
403 // Calculate distance to surface. If the side is too far
404 // behind the point, we must reject it.
405 //
406 G4ThreeVector ps = p - surface;
407 distFromSurface = -normSign*ps.dot(normal);
408
409 if (distFromSurface < -surfTolerance) { return false; }
410
411 //
412 // Calculate precise distance to intersection with the side
413 // (along the trajectory, not normal to the surface)
414 //
415 distance = distFromSurface/dotProd;
416
417 //
418 // Calculate intersection point in r,z
419 //
420 G4ThreeVector ip = p + distance*v;
421
422 G4double r = radial.dot(ip);
423
424 //
425 // And is it inside the r/z extent?
426 //
427 return InsideEdgesExact( r, ip.z(), normSign, p, v );
428}
429
430// Distance
431//
433{
434 G4double normSign = outgoing ? +1 : -1;
435 //
436 // Correct normal?
437 //
438 G4ThreeVector ps = p - surface;
439 G4double distPhi = -normSign*normal.dot(ps);
440
441 if (distPhi < -0.5*kCarTolerance)
442 {
443 return kInfinity;
444 }
445 if (distPhi < 0)
446 {
447 distPhi = 0.0;
448 }
449
450 //
451 // Calculate projected point in r,z
452 //
453 G4double r = radial.dot(p);
454
455 //
456 // Are we inside the face?
457 //
458 G4double distRZ2;
459
460 if (InsideEdges( r, p.z(), &distRZ2, nullptr ))
461 {
462 //
463 // Yup, answer is just distPhi
464 //
465 return distPhi;
466 }
467
468 //
469 // Nope. Penalize by distance out
470 //
471 return std::sqrt( distPhi*distPhi + distRZ2 );
472}
473
474// Inside
475//
477 G4double tolerance,
478 G4double* bestDistance )
479{
480 //
481 // Get distance along phi, which if negative means the point
482 // is nominally inside the shape.
483 //
484 G4ThreeVector ps = p - surface;
485 G4double distPhi = normal.dot(ps);
486
487 //
488 // Calculate projected point in r,z
489 //
490 G4double r = radial.dot(p);
491
492 //
493 // Are we inside the face?
494 //
495 G4double distRZ2;
496 G4PolyPhiFaceVertex* base3Dnorm = nullptr;
497 G4ThreeVector* head3Dnorm = nullptr;
498
499 if (InsideEdges( r, p.z(), &distRZ2, &base3Dnorm, &head3Dnorm ))
500 {
501 //
502 // Looks like we're inside. Distance is distance in phi.
503 //
504 *bestDistance = std::fabs(distPhi);
505
506 //
507 // Use distPhi to decide fate
508 //
509 if (distPhi < -tolerance) { return kInside; }
510 if (distPhi < tolerance) { return kSurface; }
511 return kOutside;
512 }
513
514 //
515 // We're outside the extent of the face,
516 // so the distance is penalized by distance from edges in RZ
517 //
518 *bestDistance = std::sqrt( distPhi*distPhi + distRZ2 );
519
520 //
521 // Use edge normal to decide fate
522 //
523 G4ThreeVector cc( base3Dnorm->r*radial.x(),
524 base3Dnorm->r*radial.y(),
525 base3Dnorm->z );
526 cc = p - cc;
527 G4double normDist = head3Dnorm->dot(cc);
528 if ( distRZ2 > tolerance*tolerance )
529 {
530 //
531 // We're far enough away that kSurface is not possible
532 //
533 return normDist < 0 ? kInside : kOutside;
534 }
535
536 if (normDist < -tolerance) { return kInside; }
537 if (normDist < tolerance) { return kSurface; }
538 return kOutside;
539}
540
541// Normal
542//
543// This virtual member is simple for our planer shape,
544// which has only one normal
545//
547 G4double* bestDistance )
548{
549 //
550 // Get distance along phi, which if negative means the point
551 // is nominally inside the shape.
552 //
553 G4double distPhi = normal.dot(p);
554
555 //
556 // Calculate projected point in r,z
557 //
558 G4double r = radial.dot(p);
559
560 //
561 // Are we inside the face?
562 //
563 G4double distRZ2;
564
565 if (InsideEdges( r, p.z(), &distRZ2, nullptr ))
566 {
567 //
568 // Yup, answer is just distPhi
569 //
570 *bestDistance = std::fabs(distPhi);
571 }
572 else
573 {
574 //
575 // Nope. Penalize by distance out
576 //
577 *bestDistance = std::sqrt( distPhi*distPhi + distRZ2 );
578 }
579
580 return normal;
581}
582
583// Extent
584//
585// This actually isn't needed by polycone or polyhedra...
586//
588{
589 G4double max = -kInfinity;
590
591 G4PolyPhiFaceVertex* corner = corners;
592 do // Loop checking, 13.08.2015, G.Cosmo
593 {
594 G4double here = axis.x()*corner->r*radial.x()
595 + axis.y()*corner->r*radial.y()
596 + axis.z()*corner->z;
597 if (here > max) { max = here; }
598 } while( ++corner < corners + numEdges );
599
600 return max;
601}
602
603// CalculateExtent
604//
605// See notes in G4VCSGface
606//
608 const G4VoxelLimits& voxelLimit,
609 const G4AffineTransform& transform,
610 G4SolidExtentList& extentList )
611{
612 //
613 // Construct a (sometimes big) clippable polygon,
614 //
615 // Perform the necessary transformations while doing so
616 //
617 G4ClippablePolygon polygon;
618
619 G4PolyPhiFaceVertex* corner = corners;
620 do // Loop checking, 13.08.2015, G.Cosmo
621 {
622 G4ThreeVector point( 0, 0, corner->z );
623 point += radial*corner->r;
624
625 polygon.AddVertexInOrder( transform.TransformPoint( point ) );
626 } while( ++corner < corners + numEdges );
627
628 //
629 // Clip away
630 //
631 if (polygon.PartialClip( voxelLimit, axis ))
632 {
633 //
634 // Add it to the list
635 //
636 polygon.SetNormal( transform.TransformAxis(normal) );
637 extentList.AddSurface( polygon );
638 }
639}
640
641// InsideEdgesExact
642//
643// Decide if the point in r,z is inside the edges of our face,
644// **but** do so consistently with other faces.
645//
646// This routine has functionality similar to InsideEdges, but uses
647// an algorithm to decide if a trajectory falls inside or outside the
648// face that uses only the trajectory p,v values and the three dimensional
649// points representing the edges of the polygon. The objective is to plug up
650// any leaks between touching G4PolyPhiFaces (at r==0) and any other face
651// that uses the same convention.
652//
653// See: "Computational Geometry in C (Second Edition)"
654// http://cs.smith.edu/~orourke/
655//
656G4bool G4PolyPhiFace::InsideEdgesExact( G4double r, G4double z,
657 G4double normSign,
658 const G4ThreeVector& p,
659 const G4ThreeVector& v )
660{
661 //
662 // Quick check of extent
663 //
664 if ( (r < rMin-kCarTolerance)
665 || (r > rMax+kCarTolerance) ) { return false; }
666
667 if ( (z < zMin-kCarTolerance)
668 || (z > zMax+kCarTolerance) ) { return false; }
669
670 //
671 // Exact check: loop over all vertices
672 //
673 G4double qx = p.x() + v.x(),
674 qy = p.y() + v.y(),
675 qz = p.z() + v.z();
676
677 G4int answer = 0;
678 G4PolyPhiFaceVertex *corn = corners,
679 *prev = corners+numEdges-1;
680
681 G4double cornZ, prevZ;
682
683 prevZ = ExactZOrder( z, qx, qy, qz, v, normSign, prev );
684 do // Loop checking, 13.08.2015, G.Cosmo
685 {
686 //
687 // Get z order of this vertex, and compare to previous vertex
688 //
689 cornZ = ExactZOrder( z, qx, qy, qz, v, normSign, corn );
690
691 if (cornZ < 0)
692 {
693 if (prevZ < 0) { continue; }
694 }
695 else if (cornZ > 0)
696 {
697 if (prevZ > 0) { continue; }
698 }
699 else
700 {
701 //
702 // By chance, we overlap exactly (within precision) with
703 // the current vertex. Continue if the same happened previously
704 // (e.g. the previous vertex had the same z value)
705 //
706 if (prevZ == 0) { continue; }
707
708 //
709 // Otherwise, to decide what to do, we need to know what is
710 // coming up next. Specifically, we need to find the next vertex
711 // with a non-zero z order.
712 //
713 // One might worry about infinite loops, but the above conditional
714 // should prevent it
715 //
716 G4PolyPhiFaceVertex *next = corn;
717 G4double nextZ;
718 do // Loop checking, 13.08.2015, G.Cosmo
719 {
720 ++next;
721 if (next == corners+numEdges) { next = corners; }
722
723 nextZ = ExactZOrder( z, qx, qy, qz, v, normSign, next );
724 } while( nextZ == 0 );
725
726 //
727 // If we won't be changing direction, go to the next vertex
728 //
729 if (nextZ*prevZ < 0) { continue; }
730 }
731
732
733 //
734 // We overlap in z with the side of the face that stretches from
735 // vertex "prev" to "corn". On which side (left or right) do
736 // we lay with respect to this segment?
737 //
738 G4ThreeVector qa( qx - prev->x, qy - prev->y, qz - prev->z ),
739 qb( qx - corn->x, qy - corn->y, qz - corn->z );
740
741 G4double aboveOrBelow = normSign*qa.cross(qb).dot(v);
742
743 if (aboveOrBelow > 0)
744 {
745 ++answer;
746 }
747 else if (aboveOrBelow < 0)
748 {
749 --answer;
750 }
751 else
752 {
753 //
754 // A precisely zero answer here means we exactly
755 // intersect (within roundoff) the edge of the face.
756 // Return true in this case.
757 //
758 return true;
759 }
760 } while( prevZ = cornZ, prev=corn, ++corn < corners+numEdges );
761
762 return answer!=0;
763}
764
765// InsideEdges (don't care about distance)
766//
767// Decide if the point in r,z is inside the edges of our face
768//
769// This routine can be made a zillion times quicker by implementing
770// better code, for example:
771//
772// int pnpoly(int npol, float *xp, float *yp, float x, float y)
773// {
774// int i, j, c = 0;
775// for (i = 0, j = npol-1; i < npol; j = i++) {
776// if ((((yp[i]<=y) && (y<yp[j])) ||
777// ((yp[j]<=y) && (y<yp[i]))) &&
778// (x < (xp[j] - xp[i]) * (y - yp[i]) / (yp[j] - yp[i]) + xp[i]))
779//
780// c = !c;
781// }
782// return c;
783// }
784//
785// See "Point in Polyon Strategies", Eric Haines [Graphic Gems IV] pp. 24-46
786//
787// The algorithm below is rather unique, but is based on code needed to
788// calculate the distance to the shape. Left here for testing...
789//
790G4bool G4PolyPhiFace::InsideEdges( G4double r, G4double z )
791{
792 //
793 // Quick check of extent
794 //
795 if ( r < rMin || r > rMax ) { return false; }
796 if ( z < zMin || z > zMax ) { return false; }
797
798 //
799 // More thorough check
800 //
801 G4double notUsed;
802
803 return InsideEdges( r, z, &notUsed, nullptr );
804}
805
806// InsideEdges (care about distance)
807//
808// Decide if the point in r,z is inside the edges of our face
809//
810G4bool G4PolyPhiFace::InsideEdges( G4double r, G4double z,
811 G4double* bestDist2,
812 G4PolyPhiFaceVertex** base3Dnorm,
813 G4ThreeVector** head3Dnorm )
814{
815 G4double bestDistance2 = kInfinity;
816 G4bool answer = false;
817
818 G4PolyPhiFaceEdge* edge = edges;
819 do // Loop checking, 13.08.2015, G.Cosmo
820 {
821 G4PolyPhiFaceVertex* testMe = nullptr;
822 G4PolyPhiFaceVertex* v0_vertex = edge->v0;
823 G4PolyPhiFaceVertex* v1_vertex = edge->v1;
824 //
825 // Get distance perpendicular to the edge
826 //
827 G4double dr = (r-v0_vertex->r), dz = (z-v0_vertex->z);
828
829 G4double distOut = dr*edge->tz - dz*edge->tr;
830 G4double distance2 = distOut*distOut;
831 if (distance2 > bestDistance2) { continue; } // No hope!
832
833 //
834 // Check to see if normal intersects edge within the edge's boundary
835 //
836 G4double q = dr*edge->tr + dz*edge->tz;
837
838 //
839 // If it doesn't, penalize distance2 appropriately
840 //
841 if (q < 0)
842 {
843 distance2 += q*q;
844 testMe = v0_vertex;
845 }
846 else if (q > edge->length)
847 {
848 G4double s2 = q-edge->length;
849 distance2 += s2*s2;
850 testMe = v1_vertex;
851 }
852
853 //
854 // Closest edge so far?
855 //
856 if (distance2 < bestDistance2)
857 {
858 bestDistance2 = distance2;
859 if (testMe != nullptr)
860 {
861 G4double distNorm = dr*testMe->rNorm + dz*testMe->zNorm;
862 answer = (distNorm <= 0);
863 if (base3Dnorm != nullptr)
864 {
865 *base3Dnorm = testMe;
866 *head3Dnorm = &testMe->norm3D;
867 }
868 }
869 else
870 {
871 answer = (distOut <= 0);
872 if (base3Dnorm != nullptr)
873 {
874 *base3Dnorm = v0_vertex;
875 *head3Dnorm = &edge->norm3D;
876 }
877 }
878 }
879 } while( ++edge < edges + numEdges );
880
881 *bestDist2 = bestDistance2;
882 return answer;
883}
884
885// Calculation of Surface Area of a Triangle
886// In the same time Random Point in Triangle is given
887//
888G4double G4PolyPhiFace::SurfaceTriangle( const G4ThreeVector& p1,
889 const G4ThreeVector& p2,
890 const G4ThreeVector& p3,
891 G4ThreeVector* p4 )
892{
893 G4ThreeVector v, w;
894
895 v = p3 - p1;
896 w = p1 - p2;
897 G4double lambda1 = G4QuickRand();
898 G4double lambda2 = lambda1*G4QuickRand();
899
900 *p4=p2 + lambda1*w + lambda2*v;
901 return 0.5*(v.cross(w)).mag();
902}
903
904// Compute surface area
905//
907{
908 if ( fSurfaceArea==0. ) { Triangulate(); }
909 return fSurfaceArea;
910}
911
912// Return random point on face
913//
915{
916 Triangulate();
917 return surface_point;
918}
919
920//
921// Auxiliary Functions used for Finding the PointOnFace using Triangulation
922//
923
924// Calculation of 2*Area of Triangle with Sign
925//
926G4double G4PolyPhiFace::Area2( const G4TwoVector& a,
927 const G4TwoVector& b,
928 const G4TwoVector& c )
929{
930 return ((b.x()-a.x())*(c.y()-a.y())-
931 (c.x()-a.x())*(b.y()-a.y()));
932}
933
934// Boolean function for sign of Surface
935//
936G4bool G4PolyPhiFace::Left( const G4TwoVector& a,
937 const G4TwoVector& b,
938 const G4TwoVector& c )
939{
940 return Area2(a,b,c)>0;
941}
942
943// Boolean function for sign of Surface
944//
945G4bool G4PolyPhiFace::LeftOn( const G4TwoVector& a,
946 const G4TwoVector& b,
947 const G4TwoVector& c )
948{
949 return Area2(a,b,c)>=0;
950}
951
952// Boolean function for sign of Surface
953//
954G4bool G4PolyPhiFace::Collinear( const G4TwoVector& a,
955 const G4TwoVector& b,
956 const G4TwoVector& c )
957{
958 return Area2(a,b,c)==0;
959}
960
961// Boolean function for finding "Proper" Intersection
962// That means Intersection of two lines segments (a,b) and (c,d)
963//
964G4bool G4PolyPhiFace::IntersectProp( const G4TwoVector& a,
965 const G4TwoVector& b,
966 const G4TwoVector& c, const G4TwoVector& d )
967{
968 if( Collinear(a,b,c) || Collinear(a,b,d)||
969 Collinear(c,d,a) || Collinear(c,d,b) ) { return false; }
970
971 G4bool Positive;
972 Positive = !(Left(a,b,c))^!(Left(a,b,d));
973 return Positive && ((!Left(c,d,a)^!Left(c,d,b)) != 0);
974}
975
976// Boolean function for determining if Point c is between a and b
977// For the tree points(a,b,c) on the same line
978//
979G4bool G4PolyPhiFace::Between( const G4TwoVector& a, const G4TwoVector& b, const G4TwoVector& c )
980{
981 if( !Collinear(a,b,c) ) { return false; }
982
983 if(a.x()!=b.x())
984 {
985 return ((a.x()<=c.x())&&(c.x()<=b.x()))||
986 ((a.x()>=c.x())&&(c.x()>=b.x()));
987 }
988
989 return ((a.y()<=c.y())&&(c.y()<=b.y()))||
990 ((a.y()>=c.y())&&(c.y()>=b.y()));
991}
992
993// Boolean function for finding Intersection "Proper" or not
994// Between two line segments (a,b) and (c,d)
995//
997 const G4TwoVector& b,
998 const G4TwoVector& c, const G4TwoVector& d )
999{
1000 if( IntersectProp(a,b,c,d) ) { return true; }
1001
1002 if( Between(a,b,c) || Between(a,b,d)
1003 || Between(c,d,a) || Between(c,d,b) ) { return true; }
1004
1005 return false;
1006}
1007
1008// Boolean Diagonalie help to determine
1009// if diagonal s of segment (a,b) is convex or reflex
1010//
1011G4bool G4PolyPhiFace::Diagonalie( G4PolyPhiFaceVertex* a,
1013{
1014 G4PolyPhiFaceVertex *corner = triangles;
1015 G4PolyPhiFaceVertex *corner_next=triangles;
1016
1017 // For each Edge (corner,corner_next)
1018 do // Loop checking, 13.08.2015, G.Cosmo
1019 {
1020 corner_next=corner->next;
1021
1022 // Skip edges incident to a of b
1023 //
1024 if( (corner!=a)&&(corner_next!=a)
1025 &&(corner!=b)&&(corner_next!=b) )
1026 {
1027 G4TwoVector rz1,rz2,rz3,rz4;
1028 rz1 = G4TwoVector(a->r,a->z);
1029 rz2 = G4TwoVector(b->r,b->z);
1030 rz3 = G4TwoVector(corner->r,corner->z);
1031 rz4 = G4TwoVector(corner_next->r,corner_next->z);
1032 if( Intersect(rz1,rz2,rz3,rz4) ) { return false; }
1033 }
1034 corner=corner->next;
1035
1036 } while( corner != triangles );
1037
1038 return true;
1039}
1040
1041// Boolean function that determine if b is Inside Cone (a0,a,a1)
1042// being a the center of the Cone
1043//
1044G4bool G4PolyPhiFace::InCone( G4PolyPhiFaceVertex* a, G4PolyPhiFaceVertex* b )
1045{
1046 // a0,a and a1 are consecutive vertices
1047 //
1048 G4PolyPhiFaceVertex *a0,*a1;
1049 a1=a->next;
1050 a0=a->prev;
1051
1052 G4TwoVector arz,arz0,arz1,brz;
1053 arz=G4TwoVector(a->r,a->z);arz0=G4TwoVector(a0->r,a0->z);
1054 arz1=G4TwoVector(a1->r,a1->z);brz=G4TwoVector(b->r,b->z);
1055
1056
1057 if(LeftOn(arz,arz1,arz0)) // If a is convex vertex
1058 {
1059 return Left(arz,brz,arz0)&&Left(brz,arz,arz1);
1060 }
1061
1062 // Else a is reflex
1063 return !( LeftOn(arz,brz,arz1)&&LeftOn(brz,arz,arz0));
1064}
1065
1066// Boolean function finding if Diagonal is possible
1067// inside Polycone or PolyHedra
1068//
1069G4bool G4PolyPhiFace::Diagonal( G4PolyPhiFaceVertex* a, G4PolyPhiFaceVertex* b )
1070{
1071 return InCone(a,b) && InCone(b,a) && Diagonalie(a,b);
1072}
1073
1074// Initialisation for Triangulisation by ear tips
1075// For details see "Computational Geometry in C" by Joseph O'Rourke
1076//
1077void G4PolyPhiFace::EarInit()
1078{
1079 G4PolyPhiFaceVertex* corner = triangles;
1080 G4PolyPhiFaceVertex* c_prev, *c_next;
1081
1082 do // Loop checking, 13.08.2015, G.Cosmo
1083 {
1084 // We need to determine three consecutive vertices
1085 //
1086 c_next=corner->next;
1087 c_prev=corner->prev;
1088
1089 // Calculation of ears
1090 //
1091 corner->ear=Diagonal(c_prev,c_next);
1092 corner=corner->next;
1093
1094 } while( corner!=triangles );
1095}
1096
1097// Triangulisation by ear tips for Polycone or Polyhedra
1098// For details see "Computational Geometry in C" by Joseph O'Rourke
1099//
1100void G4PolyPhiFace::Triangulate()
1101{
1102 // The copy of Polycone is made and this copy is reordered in order to
1103 // have a list of triangles. This list is used for GetPointOnFace().
1104
1105 const std::size_t maxEdges = (numEdges > 0) ? numEdges : 1;
1106 auto tri_help = new G4PolyPhiFaceVertex[maxEdges];
1107 triangles = tri_help;
1108 G4PolyPhiFaceVertex* triang = triangles;
1109
1110 std::vector<G4double> areas;
1111 std::vector<G4ThreeVector> points;
1112 G4double area=0.;
1113 G4PolyPhiFaceVertex *v0,*v1,*v2,*v3,*v4;
1114 v2=triangles;
1115
1116 // Make copy for prev/next for triang=corners
1117 //
1118 G4PolyPhiFaceVertex* helper = corners;
1119 G4PolyPhiFaceVertex* helper2 = corners;
1120 do // Loop checking, 13.08.2015, G.Cosmo
1121 {
1122 triang->r = helper->r;
1123 triang->z = helper->z;
1124 triang->x = helper->x;
1125 triang->y = helper->y;
1126
1127 // add pointer on prev corner
1128 //
1129 if( helper==corners )
1130 { triang->prev=triangles+maxEdges-1; }
1131 else
1132 { triang->prev=helper2; }
1133
1134 // add pointer on next corner
1135 //
1136 if( helper<corners+maxEdges-1 )
1137 { triang->next=triang+1; }
1138 else
1139 { triang->next=triangles; }
1140 helper2=triang;
1141 helper=helper->next;
1142 triang=triang->next;
1143
1144 } while( helper!=corners );
1145
1146 EarInit();
1147
1148 G4int n=numEdges;
1149 G4int i=0;
1150 G4ThreeVector p1,p2,p3,p4;
1151 const G4int max_n_loops=numEdges*10000; // protection against infinite loop
1152
1153 // Each step of outer loop removes one ear
1154 //
1155 while(n>3) // Loop checking, 13.08.2015, G.Cosmo
1156 { // Inner loop searches for one ear
1157 v2=triangles;
1158 do // Loop checking, 13.08.2015, G.Cosmo
1159 {
1160 if(v2->ear) // Ear found. Fill variables
1161 {
1162 // (v1,v3) is diagonal
1163 //
1164 v3=v2->next; v4=v3->next;
1165 v1=v2->prev; v0=v1->prev;
1166
1167 // Calculate areas and points
1168
1169 p1=G4ThreeVector((v2)->x,(v2)->y,(v2)->z);
1170 p2=G4ThreeVector((v1)->x,(v1)->y,(v1)->z);
1171 p3=G4ThreeVector((v3)->x,(v3)->y,(v3)->z);
1172
1173 G4double result1 = SurfaceTriangle(p1,p2,p3,&p4 );
1174 points.push_back(p4);
1175 areas.push_back(result1);
1176 area=area+result1;
1177
1178 // Update earity of diagonal endpoints
1179 //
1180 v1->ear=Diagonal(v0,v3);
1181 v3->ear=Diagonal(v1,v4);
1182
1183 // Cut off the ear v2
1184 // Has to be done for a copy and not for real PolyPhiFace
1185 //
1186 v1->next=v3;
1187 v3->prev=v1;
1188 triangles=v3; // In case the head was v2
1189 --n;
1190
1191 break; // out of inner loop
1192 } // end if ear found
1193
1194 v2=v2->next;
1195
1196 } while( v2!=triangles );
1197
1198 ++i;
1199 if(i>=max_n_loops)
1200 {
1201 G4Exception( "G4PolyPhiFace::Triangulation()",
1202 "GeomSolids0003", FatalException,
1203 "Maximum number of steps is reached for triangulation!" );
1204 }
1205 } // end outer while loop
1206
1207 if(v2->next != nullptr)
1208 {
1209 // add last triangle
1210 //
1211 v2=v2->next;
1212 p1=G4ThreeVector((v2)->x,(v2)->y,(v2)->z);
1213 p2=G4ThreeVector((v2->next)->x,(v2->next)->y,(v2->next)->z);
1214 p3=G4ThreeVector((v2->prev)->x,(v2->prev)->y,(v2->prev)->z);
1215 G4double result1 = SurfaceTriangle(p1,p2,p3,&p4 );
1216 points.push_back(p4);
1217 areas.push_back(result1);
1218 area=area+result1;
1219 }
1220
1221 // Surface Area is stored
1222 //
1223 fSurfaceArea = area;
1224
1225 // Second Step: choose randomly one surface
1226 //
1227 G4double chose = area*G4QuickRand();
1228
1229 // Third Step: Get a point on choosen surface
1230 //
1231 G4double Achose1, Achose2;
1232 Achose1=0; Achose2=0.;
1233 i=0;
1234 do // Loop checking, 13.08.2015, G.Cosmo
1235 {
1236 Achose2+=areas[i];
1237 if(chose>=Achose1 && chose<Achose2)
1238 {
1239 G4ThreeVector point;
1240 point=points[i];
1241 surface_point=point;
1242 break;
1243 }
1244 ++i; Achose1=Achose2;
1245 } while( i<numEdges-2 );
1246
1247 delete [] tri_help;
1248 tri_help = nullptr;
1249}
const G4double kCarTolerance
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
const G4double a0
G4double G4QuickRand(uint32_t seed=0)
CLHEP::Hep3Vector G4ThreeVector
CLHEP::Hep2Vector G4TwoVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
double x() const
double y() const
double z() const
Hep3Vector unit() const
double x() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
double dot(const Hep3Vector &) const
double mag() const
G4AffineTransform is a class for geometric affine transformations. It supports efficient arbitrary ro...
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) 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)
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4PolyPhiFace is a face that bounds a polycone or polyhedra when it has a phi opening....
G4double Distance(const G4ThreeVector &p, G4bool outgoing) override
EInside Inside(const G4ThreeVector &p, G4double tolerance, G4double *bestDistance) override
void Diagnose(G4VSolid *solid)
void CalculateExtent(const EAxis axis, const G4VoxelLimits &voxelLimit, const G4AffineTransform &tranform, G4SolidExtentList &extentList) override
G4PolyPhiFace & operator=(const G4PolyPhiFace &source)
G4double Extent(const G4ThreeVector axis) override
G4ThreeVector Normal(const G4ThreeVector &p, G4double *bestDistance) override
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, G4bool outgoing, G4double surfTolerance, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal, G4bool &allBehind) override
~G4PolyPhiFace() override
G4double SurfaceArea() override
G4PolyPhiFace(const G4ReduciblePolygon *rz, G4double phi, G4double deltaPhi, G4double phiOther)
G4ReduciblePolygonIterator is companion class for iterating over the vertices of a polygon.
G4ReduciblePolygon is a utility class used to specify, test, reduce, and/or otherwise manipulate a 2D...
G4double Amin() const
G4double Bmin() const
G4double Amax() const
G4SolidExtentList is utility class designed for calculating the extent of a CSG-like solid for a voxe...
void AddSurface(const G4ClippablePolygon &surface)
virtual G4ThreeVector GetPointOnFace()=0
G4VSolid is an abstract base class for solids, physical shapes that can be tracked through....
Definition G4VSolid.hh:80
virtual EInside Inside(const G4ThreeVector &p) const =0
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 axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668
G4PolyPhiFaceVertex * v1
G4ThreeVector norm3D
G4PolyPhiFaceVertex * v0
G4ThreeVector norm3D
G4PolyPhiFaceVertex * next
G4PolyPhiFaceVertex * prev
#define DBL_MIN
Definition templates.hh:54