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