Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PolyhedraSide.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// $Id$
28//
29//
30// --------------------------------------------------------------------
31// GEANT 4 class source file
32//
33//
34// G4PolyhedraSide.cc
35//
36// Implementation of the face representing one segmented side of a Polyhedra
37//
38// --------------------------------------------------------------------
39
40#include "G4PolyhedraSide.hh"
42#include "G4IntersectingCone.hh"
43#include "G4ClippablePolygon.hh"
44#include "G4AffineTransform.hh"
45#include "G4SolidExtentList.hh"
47
48#include "Randomize.hh"
49
50//
51// Constructor
52//
53// Values for r1,z1 and r2,z2 should be specified in clockwise
54// order in (r,z).
55//
57 const G4PolyhedraSideRZ *tail,
58 const G4PolyhedraSideRZ *head,
59 const G4PolyhedraSideRZ *nextRZ,
60 G4int theNumSide,
61 G4double thePhiStart,
62 G4double thePhiTotal,
63 G4bool thePhiIsOpen,
64 G4bool isAllBehind )
65{
67 fSurfaceArea=0.;
68 fPhi.first = G4ThreeVector(0,0,0);
69 fPhi.second= 0.0;
70
71 //
72 // Record values
73 //
74 r[0] = tail->r; z[0] = tail->z;
75 r[1] = head->r; z[1] = head->z;
76
77 G4double phiTotal;
78
79 //
80 // Set phi to our convention
81 //
82 startPhi = thePhiStart;
83 while (startPhi < 0.0) startPhi += twopi;
84
85 phiIsOpen = thePhiIsOpen;
86 phiTotal = (phiIsOpen) ? thePhiTotal : twopi;
87
88 allBehind = isAllBehind;
89
90 //
91 // Make our intersecting cone
92 //
93 cone = new G4IntersectingCone( r, z );
94
95 //
96 // Construct side plane vector set
97 //
98 numSide = theNumSide;
99 deltaPhi = phiTotal/theNumSide;
100 endPhi = startPhi+phiTotal;
101
103
105
106 //
107 // ...this is where we start
108 //
109 G4double phi = startPhi;
110 G4ThreeVector a1( r[0]*std::cos(phi), r[0]*std::sin(phi), z[0] ),
111 b1( r[1]*std::cos(phi), r[1]*std::sin(phi), z[1] ),
112 c1( prevRZ->r*std::cos(phi), prevRZ->r*std::sin(phi), prevRZ->z ),
113 d1( nextRZ->r*std::cos(phi), nextRZ->r*std::sin(phi), nextRZ->z ),
114 a2, b2, c2, d2;
116
118 do
119 {
120 //
121 // ...this is where we are going
122 //
123 phi += deltaPhi;
124 a2 = G4ThreeVector( r[0]*std::cos(phi), r[0]*std::sin(phi), z[0] );
125 b2 = G4ThreeVector( r[1]*std::cos(phi), r[1]*std::sin(phi), z[1] );
126 c2 = G4ThreeVector( prevRZ->r*std::cos(phi), prevRZ->r*std::sin(phi), prevRZ->z );
127 d2 = G4ThreeVector( nextRZ->r*std::cos(phi), nextRZ->r*std::sin(phi), nextRZ->z );
128
129 G4ThreeVector tt;
130
131 //
132 // ...build some relevant vectors.
133 // the point is to sacrifice a little memory with precalcs
134 // to gain speed
135 //
136 vec->center = 0.25*( a1 + a2 + b1 + b2 );
137
138 tt = b2 + b1 - a2 - a1;
139 vec->surfRZ = tt.unit();
140 if (vec==vecs) lenRZ = 0.25*tt.mag();
141
142 tt = b2 - b1 + a2 - a1;
143 vec->surfPhi = tt.unit();
144 if (vec==vecs)
145 {
146 lenPhi[0] = 0.25*tt.mag();
147 tt = b2 - b1;
148 lenPhi[1] = (0.5*tt.mag()-lenPhi[0])/lenRZ;
149 }
150
151 tt = vec->surfPhi.cross(vec->surfRZ);
152 vec->normal = tt.unit();
153
154 //
155 // ...edge normals are the average of the normals of
156 // the two faces they connect.
157 //
158 // ...edge normals are necessary if we are to accurately
159 // decide if a point is "inside" a face. For non-convex
160 // shapes, it is absolutely necessary to know information
161 // on adjacent faces to accurate determine this.
162 //
163 // ...we don't need them for the phi edges, since that
164 // information is taken care of internally. The r/z edges,
165 // however, depend on the adjacent G4PolyhedraSide.
166 //
167 G4ThreeVector a12, adj;
168
169 a12 = a2-a1;
170
171 adj = 0.5*(c1+c2-a1-a2);
172 adj = adj.cross(a12);
173 adj = adj.unit() + vec->normal;
174 vec->edgeNorm[0] = adj.unit();
175
176 a12 = b1-b2;
177 adj = 0.5*(d1+d2-b1-b2);
178 adj = adj.cross(a12);
179 adj = adj.unit() + vec->normal;
180 vec->edgeNorm[1] = adj.unit();
181
182 //
183 // ...the corners are crucial. It is important that
184 // they are calculated consistently for adjacent
185 // G4PolyhedraSides, to avoid gaps caused by roundoff.
186 //
187 vec->edges[0] = edge;
188 edge->corner[0] = a1;
189 edge->corner[1] = b1;
190 edge++;
191 vec->edges[1] = edge;
192
193 a1 = a2;
194 b1 = b2;
195 c1 = c2;
196 d1 = d2;
197 } while( ++vec < vecs+numSide );
198
199 //
200 // Clean up hanging edge
201 //
202 if (phiIsOpen)
203 {
204 edge->corner[0] = a2;
205 edge->corner[1] = b2;
206 }
207 else
208 {
209 vecs[numSide-1].edges[1] = edges;
210 }
211
212 //
213 // Go back and fill in remaining fields in edges
214 //
215 vec = vecs;
217 do
218 {
219 edge = vec->edges[0]; // The edge between prev and vec
220
221 //
222 // Okay: edge normal is average of normals of adjacent faces
223 //
224 G4ThreeVector eNorm = vec->normal + prev->normal;
225 edge->normal = eNorm.unit();
226
227 //
228 // Vertex normal is average of norms of adjacent surfaces (all four)
229 // However, vec->edgeNorm is unit vector in some direction
230 // as the sum of normals of adjacent PolyhedraSide with vec.
231 // The normalization used for this vector should be the same
232 // for vec and prev.
233 //
234 eNorm = vec->edgeNorm[0] + prev->edgeNorm[0];
235 edge->cornNorm[0] = eNorm.unit();
236
237 eNorm = vec->edgeNorm[1] + prev->edgeNorm[1];
238 edge->cornNorm[1] = eNorm.unit();
239 } while( prev=vec, ++vec < vecs + numSide );
240
241 if (phiIsOpen)
242 {
243 // G4double rFact = std::cos(0.5*deltaPhi);
244 //
245 // If phi is open, we need to patch up normals of the
246 // first and last edges and their corresponding
247 // vertices.
248 //
249 // We use vectors that are in the plane of the
250 // face. This should be safe.
251 //
252 vec = vecs;
253
254 G4ThreeVector normvec = vec->edges[0]->corner[0]
255 - vec->edges[0]->corner[1];
256 normvec = normvec.cross(vec->normal);
257 if (normvec.dot(vec->surfPhi) > 0) normvec = -normvec;
258
259 vec->edges[0]->normal = normvec.unit();
260
261 vec->edges[0]->cornNorm[0] = (vec->edges[0]->corner[0]
262 - vec->center).unit();
263 vec->edges[0]->cornNorm[1] = (vec->edges[0]->corner[1]
264 - vec->center).unit();
265
266 //
267 // Repeat for ending phi
268 //
269 vec = vecs + numSide - 1;
270
271 normvec = vec->edges[1]->corner[0] - vec->edges[1]->corner[1];
272 normvec = normvec.cross(vec->normal);
273 if (normvec.dot(vec->surfPhi) < 0) normvec = -normvec;
274
275 vec->edges[1]->normal = normvec.unit();
276
277 vec->edges[1]->cornNorm[0] = (vec->edges[1]->corner[0]
278 - vec->center).unit();
279 vec->edges[1]->cornNorm[1] = (vec->edges[1]->corner[1]
280 - vec->center).unit();
281 }
282
283 //
284 // edgeNorm is the factor one multiplies the distance along vector phi
285 // on the surface of one of our sides in order to calculate the distance
286 // from the edge. (see routine DistanceAway)
287 //
288 edgeNorm = 1.0/std::sqrt( 1.0 + lenPhi[1]*lenPhi[1] );
289}
290
291
292//
293// Fake default constructor - sets only member data and allocates memory
294// for usage restricted to object persistency.
295//
297 : numSide(0), startPhi(0.), deltaPhi(0.), endPhi(0.),
298 phiIsOpen(false), allBehind(false), cone(0), vecs(0), edges(0),
299 lenRZ(0.), edgeNorm(0.), kCarTolerance(0.), fSurfaceArea(0.)
300{
301 r[0] = r[1] = 0.;
302 z[0] = z[1] = 0.;
303 lenPhi[0] = lenPhi[1] = 0.;
304}
305
306
307//
308// Destructor
309//
311{
312 delete cone;
313 delete [] vecs;
314 delete [] edges;
315}
316
317
318//
319// Copy constructor
320//
322 : G4VCSGface()
323{
324 CopyStuff( source );
325}
326
327
328//
329// Assignment operator
330//
332{
333 if (this == &source) return *this;
334
335 delete cone;
336 delete [] vecs;
337 delete [] edges;
338
339 CopyStuff( source );
340
341 return *this;
342}
343
344
345//
346// CopyStuff
347//
349{
350 //
351 // The simple stuff
352 //
353 numSide = source.numSide;
354 r[0] = source.r[0];
355 r[1] = source.r[1];
356 z[0] = source.z[0];
357 z[1] = source.z[1];
358 startPhi = source.startPhi;
359 deltaPhi = source.deltaPhi;
360 endPhi = source.endPhi;
361 phiIsOpen = source.phiIsOpen;
362 allBehind = source.allBehind;
363
364 lenRZ = source.lenRZ;
365 lenPhi[0] = source.lenPhi[0];
366 lenPhi[1] = source.lenPhi[1];
367 edgeNorm = source.edgeNorm;
368
369 kCarTolerance = source.kCarTolerance;
370 fSurfaceArea = source.fSurfaceArea;
371
372 cone = new G4IntersectingCone( *source.cone );
373
374 //
375 // Duplicate edges
376 //
377 G4int numEdges = phiIsOpen ? numSide+1 : numSide;
378 edges = new G4PolyhedraSideEdge[numEdges];
379
381 *sourceEdge = source.edges;
382 do
383 {
384 *edge = *sourceEdge;
385 } while( ++sourceEdge, ++edge < edges + numEdges);
386
387 //
388 // Duplicate vecs
389 //
391
393 *sourceVec = source.vecs;
394 do
395 {
396 *vec = *sourceVec;
397 vec->edges[0] = edges + (sourceVec->edges[0] - source.edges);
398 vec->edges[1] = edges + (sourceVec->edges[1] - source.edges);
399 } while( ++sourceVec, ++vec < vecs + numSide );
400}
401
402
403//
404// Intersect
405//
406// Decide if a line intersects the face.
407//
408// Arguments:
409// p = (in) starting point of line segment
410// v = (in) direction of line segment (assumed a unit vector)
411// A, B = (in) 2d transform variables (see note top of file)
412// normSign = (in) desired sign for dot product with normal (see below)
413// surfTolerance = (in) minimum distance from the surface
414// vecs = (in) Vector set array
415// distance = (out) distance to surface furfilling all requirements
416// distFromSurface = (out) distance from the surface
417// thisNormal = (out) normal vector of the intersecting surface
418//
419// Return value:
420// true if an intersection is found. Otherwise, output parameters are
421// undefined.
422//
423// Notes:
424// * normSign: if we are "inside" the shape and only want to find out how far
425// to leave the shape, we only want to consider intersections with surfaces in
426// which the trajectory is leaving the shape. Since the normal vectors to the
427// surface always point outwards from the inside, this means we want the dot
428// product of the trajectory direction v and the normal of the side normals[i]
429// to be positive. Thus, we should specify normSign as +1.0. Otherwise, if
430// we are outside and want to go in, normSign should be set to -1.0.
431// Don't set normSign to zero, or you will get no intersections!
432//
433// * surfTolerance: see notes on argument "surfTolerance" in routine
434// "IntersectSidePlane".
435// ----HOWEVER---- We should *not* apply this surface tolerance if the
436// starting point is not within phi or z of the surface. Specifically,
437// if the starting point p angle in x/y places it on a separate side from the
438// intersection or if the starting point p is outside the z bounds of the
439// segment, surfTolerance must be ignored or we should *always* accept the
440// intersection!
441// This is simply because the sides do not have infinite extent.
442//
443//
445 const G4ThreeVector &v,
446 G4bool outgoing,
447 G4double surfTolerance,
448 G4double &distance,
449 G4double &distFromSurface,
450 G4ThreeVector &normal,
451 G4bool &isAllBehind )
452{
453 G4double normSign = outgoing ? +1 : -1;
454
455 //
456 // ------------------TO BE IMPLEMENTED---------------------
457 // Testing the intersection of individual phi faces is
458 // pretty straight forward. The simple thing therefore is to
459 // form a loop and check them all in sequence.
460 //
461 // But, I worry about one day someone making
462 // a polygon with a thousands sides. A linear search
463 // would not be ideal in such a case.
464 //
465 // So, it would be nice to be able to quickly decide
466 // which face would be intersected. One can make a very
467 // good guess by using the intersection with a cone.
468 // However, this is only reliable in 99% of the cases.
469 //
470 // My solution: make a decent guess as to the one or
471 // two potential faces might get intersected, and then
472 // test them. If we have the wrong face, use the test
473 // to make a better guess.
474 //
475 // Since we might have two guesses, form a queue of
476 // potential intersecting faces. Keep an array of
477 // already tested faces to avoid doing one more than
478 // once.
479 //
480 // Result: at worst, an iterative search. On average,
481 // a little more than two tests would be required.
482 //
483 G4ThreeVector q = p + v;
484
485 G4int face = 0;
487 do
488 {
489 //
490 // Correct normal?
491 //
492 G4double dotProd = normSign*v.dot(vec->normal);
493 if (dotProd <= 0) continue;
494
495 //
496 // Is this face in front of the point along the trajectory?
497 //
498 G4ThreeVector delta = p - vec->center;
499 distFromSurface = -normSign*delta.dot(vec->normal);
500
501 if (distFromSurface < -surfTolerance) continue;
502
503 //
504 // phi
505 // c -------- d ^
506 // | | |
507 // a -------- b +---> r/z
508 //
509 //
510 // Do we remain on this particular segment?
511 //
512 G4ThreeVector qc = q - vec->edges[1]->corner[0];
513 G4ThreeVector qd = q - vec->edges[1]->corner[1];
514
515 if (normSign*qc.cross(qd).dot(v) < 0) continue;
516
517 G4ThreeVector qa = q - vec->edges[0]->corner[0];
518 G4ThreeVector qb = q - vec->edges[0]->corner[1];
519
520 if (normSign*qa.cross(qb).dot(v) > 0) continue;
521
522 //
523 // We found the one and only segment we might be intersecting.
524 // Do we remain within r/z bounds?
525 //
526
527 if (r[0] > 1/kInfinity && normSign*qa.cross(qc).dot(v) < 0) return false;
528 if (r[1] > 1/kInfinity && normSign*qb.cross(qd).dot(v) > 0) return false;
529
530 //
531 // We allow the face to be slightly behind the trajectory
532 // (surface tolerance) only if the point p is within
533 // the vicinity of the face
534 //
535 if (distFromSurface < 0)
536 {
537 G4ThreeVector ps = p - vec->center;
538
539 G4double rz = ps.dot(vec->surfRZ);
540 if (std::fabs(rz) > lenRZ+surfTolerance) return false;
541
542 G4double pp = ps.dot(vec->surfPhi);
543 if (std::fabs(pp) > lenPhi[0] + lenPhi[1]*rz + surfTolerance) return false;
544 }
545
546
547 //
548 // Intersection found. Return answer.
549 //
550 distance = distFromSurface/dotProd;
551 normal = vec->normal;
552 isAllBehind = allBehind;
553 return true;
554 } while( ++vec, ++face < numSide );
555
556 //
557 // Oh well. Better luck next time.
558 //
559 return false;
560}
561
562
564{
565 G4double normSign = outgoing ? -1 : +1;
566
567 //
568 // Try the closest phi segment first
569 //
570 G4int iPhi = ClosestPhiSegment( GetPhi(p) );
571
572 G4ThreeVector pdotc = p - vecs[iPhi].center;
573 G4double normDist = pdotc.dot(vecs[iPhi].normal);
574
575 if (normSign*normDist > -0.5*kCarTolerance)
576 {
577 return DistanceAway( p, vecs[iPhi], &normDist );
578 }
579
580 //
581 // Now we have an interesting problem... do we try to find the
582 // closest facing side??
583 //
584 // Considered carefully, the answer is no. We know that if we
585 // are asking for the distance out, we are supposed to be inside,
586 // and vice versa.
587 //
588
589 return kInfinity;
590}
591
592
593//
594// Inside
595//
597 G4double tolerance,
598 G4double *bestDistance )
599{
600 //
601 // Which phi segment is closest to this point?
602 //
603 G4int iPhi = ClosestPhiSegment( GetPhi(p) );
604
605 G4double norm;
606
607 //
608 // Get distance to this segment
609 //
610 *bestDistance = DistanceToOneSide( p, vecs[iPhi], &norm );
611
612 //
613 // Use distance along normal to decide return value
614 //
615 if ( (std::fabs(norm) < tolerance) && (*bestDistance < 2.0*tolerance) )
616 return kSurface;
617 else if (norm < 0)
618 return kInside;
619 else
620 return kOutside;
621}
622
623
624//
625// Normal
626//
628 G4double *bestDistance )
629{
630 //
631 // Which phi segment is closest to this point?
632 //
633 G4int iPhi = ClosestPhiSegment( GetPhi(p) );
634
635 //
636 // Get distance to this segment
637 //
638 G4double norm;
639 *bestDistance = DistanceToOneSide( p, vecs[iPhi], &norm );
640
641 return vecs[iPhi].normal;
642}
643
644
645//
646// Extent
647//
649{
650 if (axis.perp2() < DBL_MIN)
651 {
652 //
653 // Special case
654 //
655 return axis.z() < 0 ? -cone->ZLo() : cone->ZHi();
656 }
657
658 G4int iPhi, i1, i2;
659 G4double best;
660 G4ThreeVector *list[4];
661
662 //
663 // Which phi segment, if any, does the axis belong to
664 //
665 iPhi = PhiSegment( GetPhi(axis) );
666
667 if (iPhi < 0)
668 {
669 //
670 // No phi segment? Check front edge of first side and
671 // last edge of second side
672 //
673 i1 = 0; i2 = numSide-1;
674 }
675 else
676 {
677 //
678 // Check all corners of matching phi side
679 //
680 i1 = iPhi; i2 = iPhi;
681 }
682
683 list[0] = vecs[i1].edges[0]->corner;
684 list[1] = vecs[i1].edges[0]->corner+1;
685 list[2] = vecs[i2].edges[1]->corner;
686 list[3] = vecs[i2].edges[1]->corner+1;
687
688 //
689 // Who's biggest?
690 //
691 best = -kInfinity;
692 G4ThreeVector **vec = list;
693 do
694 {
695 G4double answer = (*vec)->dot(axis);
696 if (answer > best) best = answer;
697 } while( ++vec < list+4 );
698
699 return best;
700}
701
702
703//
704// CalculateExtent
705//
706// See notes in G4VCSGface
707//
709 const G4VoxelLimits &voxelLimit,
710 const G4AffineTransform &transform,
711 G4SolidExtentList &extentList )
712{
713 //
714 // Loop over all sides
715 //
717 do
718 {
719 //
720 // Fill our polygon with the four corners of
721 // this side, after the specified transformation
722 //
723 G4ClippablePolygon polygon;
724
725 polygon.AddVertexInOrder(transform.
726 TransformPoint(vec->edges[0]->corner[0]));
727 polygon.AddVertexInOrder(transform.
728 TransformPoint(vec->edges[0]->corner[1]));
729 polygon.AddVertexInOrder(transform.
730 TransformPoint(vec->edges[1]->corner[1]));
731 polygon.AddVertexInOrder(transform.
732 TransformPoint(vec->edges[1]->corner[0]));
733
734 //
735 // Get extent
736 //
737 if (polygon.PartialClip( voxelLimit, axis ))
738 {
739 //
740 // Get dot product of normal along target axis
741 //
742 polygon.SetNormal( transform.TransformAxis(vec->normal) );
743
744 extentList.AddSurface( polygon );
745 }
746 } while( ++vec < vecs+numSide );
747
748 return;
749}
750
751
752//
753// IntersectSidePlane
754//
755// Decide if a line correctly intersects one side plane of our segment.
756// It is assumed that the correct side has been chosen, and thus only
757// the z bounds (of the entire segment) are checked.
758//
759// normSign - To be multiplied against normal:
760// = +1.0 normal is unchanged
761// = -1.0 normal is reversed (now points inward)
762//
763// Arguments:
764// p - (in) Point
765// v - (in) Direction
766// vec - (in) Description record of the side plane
767// normSign - (in) Sign (+/- 1) to apply to normal
768// surfTolerance - (in) Surface tolerance (generally > 0, see below)
769// distance - (out) Distance along v to intersection
770// distFromSurface - (out) Distance from surface normal
771//
772// Notes:
773// surfTolerance - Used to decide if a point is behind the surface,
774// a point is allow to be -surfTolerance behind the
775// surface (as measured along the normal), but *only*
776// if the point is within the r/z bounds + surfTolerance
777// of the segment.
778//
780 const G4ThreeVector &v,
781 const G4PolyhedraSideVec& vec,
782 G4double normSign,
783 G4double surfTolerance,
784 G4double &distance,
785 G4double &distFromSurface )
786{
787 //
788 // Correct normal? Here we have straight sides, and can safely ignore
789 // intersections where the dot product with the normal is zero.
790 //
791 G4double dotProd = normSign*v.dot(vec.normal);
792
793 if (dotProd <= 0) return false;
794
795 //
796 // Calculate distance to surface. If the side is too far
797 // behind the point, we must reject it.
798 //
799 G4ThreeVector delta = p - vec.center;
800 distFromSurface = -normSign*delta.dot(vec.normal);
801
802 if (distFromSurface < -surfTolerance) return false;
803
804 //
805 // Calculate precise distance to intersection with the side
806 // (along the trajectory, not normal to the surface)
807 //
808 distance = distFromSurface/dotProd;
809
810 //
811 // Do we fall off the r/z extent of the segment?
812 //
813 // Calculate this very, very carefully! Why?
814 // 1. If a RZ end is at R=0, you can't miss!
815 // 2. If you just fall off in RZ, the answer must
816 // be consistent with adjacent G4PolyhedraSide faces.
817 // (2) implies that only variables used by other G4PolyhedraSide
818 // faces may be used, which includes only: p, v, and the edge corners.
819 // It also means that one side is a ">" or "<", which the other
820 // must be ">=" or "<=". Fortunately, this isn't a new problem.
821 // The solution below I borrowed from Joseph O'Rourke,
822 // "Computational Geometry in C (Second Edition)"
823 // See: http://cs.smith.edu/~orourke/
824 //
825 G4ThreeVector ic = p + distance*v - vec.center;
826 G4double atRZ = vec.surfRZ.dot(ic);
827
828 if (atRZ < 0)
829 {
830 if (r[0]==0) return true; // Can't miss!
831
832 if (atRZ < -lenRZ*1.2) return false; // Forget it! Missed by a mile.
833
834 G4ThreeVector q = p + v;
835 G4ThreeVector qa = q - vec.edges[0]->corner[0],
836 qb = q - vec.edges[1]->corner[0];
837 G4ThreeVector qacb = qa.cross(qb);
838 if (normSign*qacb.dot(v) < 0) return false;
839
840 if (distFromSurface < 0)
841 {
842 if (atRZ < -lenRZ-surfTolerance) return false;
843 }
844 }
845 else if (atRZ > 0)
846 {
847 if (r[1]==0) return true; // Can't miss!
848
849 if (atRZ > lenRZ*1.2) return false; // Missed by a mile
850
851 G4ThreeVector q = p + v;
852 G4ThreeVector qa = q - vec.edges[0]->corner[1],
853 qb = q - vec.edges[1]->corner[1];
854 G4ThreeVector qacb = qa.cross(qb);
855 if (normSign*qacb.dot(v) >= 0) return false;
856
857 if (distFromSurface < 0)
858 {
859 if (atRZ > lenRZ+surfTolerance) return false;
860 }
861 }
862
863 return true;
864}
865
866
867//
868// LineHitsSegments
869//
870// Calculate which phi segments a line intersects in three dimensions.
871// No check is made as to whether the intersections are within the z bounds of
872// the segment.
873//
875 const G4ThreeVector &v,
876 G4int *i1, G4int *i2 )
877{
878 G4double s1, s2;
879 //
880 // First, decide if and where the line intersects the cone
881 //
882 G4int n = cone->LineHitsCone( p, v, &s1, &s2 );
883
884 if (n==0) return 0;
885
886 //
887 // Try first intersection.
888 //
889 *i1 = PhiSegment( std::atan2( p.y() + s1*v.y(), p.x() + s1*v.x() ) );
890 if (n==1)
891 {
892 return (*i1 < 0) ? 0 : 1;
893 }
894
895 //
896 // Try second intersection
897 //
898 *i2 = PhiSegment( std::atan2( p.y() + s2*v.y(), p.x() + s2*v.x() ) );
899 if (*i1 == *i2) return 0;
900
901 if (*i1 < 0)
902 {
903 if (*i2 < 0) return 0;
904 *i1 = *i2;
905 return 1;
906 }
907
908 if (*i2 < 0) return 1;
909
910 return 2;
911}
912
913
914//
915// ClosestPhiSegment
916//
917// Decide which phi segment is closest in phi to the point.
918// The result is the same as PhiSegment if there is no phi opening.
919//
921{
922 G4int iPhi = PhiSegment( phi0 );
923 if (iPhi >= 0) return iPhi;
924
925 //
926 // Boogers! The points falls inside the phi segment.
927 // Look for the closest point: the start, or end
928 //
929 G4double phi = phi0;
930
931 while( phi < startPhi ) phi += twopi;
932 G4double d1 = phi-endPhi;
933
934 while( phi > startPhi ) phi -= twopi;
935 G4double d2 = startPhi-phi;
936
937 return (d2 < d1) ? 0 : numSide-1;
938}
939
940
941//
942// PhiSegment
943//
944// Decide which phi segment an angle belongs to, counting from zero.
945// A value of -1 indicates that the phi value is outside the shape
946// (only possible if phiTotal < 360 degrees).
947//
949{
950 //
951 // How far are we from phiStart? Come up with a positive answer
952 // that is less than 2*PI
953 //
954 G4double phi = phi0 - startPhi;
955 while( phi < 0 ) phi += twopi;
956 while( phi > twopi ) phi -= twopi;
957
958 //
959 // Divide
960 //
961 G4int answer = (G4int)(phi/deltaPhi);
962
963 if (answer >= numSide)
964 {
965 if (phiIsOpen)
966 {
967 return -1; // Looks like we missed
968 }
969 else
970 {
971 answer = numSide-1; // Probably just roundoff
972 }
973 }
974
975 return answer;
976}
977
978
979//
980// GetPhi
981//
982// Calculate Phi for a given 3-vector (point), if not already cached for the
983// same point, in the attempt to avoid consecutive computation of the same
984// quantity
985//
987{
988 G4double val=0.;
989
990 if (fPhi.first != p)
991 {
992 val = p.phi();
993 fPhi.first = p;
994 fPhi.second = val;
995 }
996 else
997 {
998 val = fPhi.second;
999 }
1000 return val;
1001}
1002
1003
1004//
1005// DistanceToOneSide
1006//
1007// Arguments:
1008// p - (in) Point to check
1009// vec - (in) vector set of this side
1010// normDist - (out) distance normal to the side or edge, as appropriate, signed
1011// Return value = total distance from the side
1012//
1014 const G4PolyhedraSideVec &vec,
1015 G4double *normDist )
1016{
1017 G4ThreeVector pct = p - vec.center;
1018
1019 //
1020 // Get normal distance
1021 //
1022 *normDist = vec.normal.dot(pct);
1023
1024 //
1025 // Add edge penalty
1026 //
1027 return DistanceAway( p, vec, normDist );
1028}
1029
1030
1031//
1032// DistanceAway
1033//
1034// Add distance from side edges, if necesssary, to total distance,
1035// and updates normDist appropriate depending on edge normals.
1036//
1038 const G4PolyhedraSideVec &vec,
1039 G4double *normDist )
1040{
1041 G4double distOut2;
1042 G4ThreeVector pct = p - vec.center;
1043 G4double distFaceNorm = *normDist;
1044
1045 //
1046 // Okay, are we inside bounds?
1047 //
1048 G4double pcDotRZ = pct.dot(vec.surfRZ);
1049 G4double pcDotPhi = pct.dot(vec.surfPhi);
1050
1051 //
1052 // Go through all permutations.
1053 // Phi
1054 // | | ^
1055 // B | H | E |
1056 // ------[1]------------[3]----- |
1057 // |XXXXXXXXXXXXXX| +----> RZ
1058 // C |XXXXXXXXXXXXXX| F
1059 // |XXXXXXXXXXXXXX|
1060 // ------[0]------------[2]----
1061 // A | G | D
1062 // | |
1063 //
1064 // It's real messy, but at least it's quick
1065 //
1066
1067 if (pcDotRZ < -lenRZ)
1068 {
1069 G4double lenPhiZ = lenPhi[0] - lenRZ*lenPhi[1];
1070 G4double distOutZ = pcDotRZ+lenRZ;
1071 //
1072 // Below in RZ
1073 //
1074 if (pcDotPhi < -lenPhiZ)
1075 {
1076 //
1077 // ...and below in phi. Find distance to point (A)
1078 //
1079 G4double distOutPhi = pcDotPhi+lenPhiZ;
1080 distOut2 = distOutPhi*distOutPhi + distOutZ*distOutZ;
1081 G4ThreeVector pa = p - vec.edges[0]->corner[0];
1082 *normDist = pa.dot(vec.edges[0]->cornNorm[0]);
1083 }
1084 else if (pcDotPhi > lenPhiZ)
1085 {
1086 //
1087 // ...and above in phi. Find distance to point (B)
1088 //
1089 G4double distOutPhi = pcDotPhi-lenPhiZ;
1090 distOut2 = distOutPhi*distOutPhi + distOutZ*distOutZ;
1091 G4ThreeVector pb = p - vec.edges[1]->corner[0];
1092 *normDist = pb.dot(vec.edges[1]->cornNorm[0]);
1093 }
1094 else
1095 {
1096 //
1097 // ...and inside in phi. Find distance to line (C)
1098 //
1099 G4ThreeVector pa = p - vec.edges[0]->corner[0];
1100 distOut2 = distOutZ*distOutZ;
1101 *normDist = pa.dot(vec.edgeNorm[0]);
1102 }
1103 }
1104 else if (pcDotRZ > lenRZ)
1105 {
1106 G4double lenPhiZ = lenPhi[0] + lenRZ*lenPhi[1];
1107 G4double distOutZ = pcDotRZ-lenRZ;
1108 //
1109 // Above in RZ
1110 //
1111 if (pcDotPhi < -lenPhiZ)
1112 {
1113 //
1114 // ...and below in phi. Find distance to point (D)
1115 //
1116 G4double distOutPhi = pcDotPhi+lenPhiZ;
1117 distOut2 = distOutPhi*distOutPhi + distOutZ*distOutZ;
1118 G4ThreeVector pd = p - vec.edges[0]->corner[1];
1119 *normDist = pd.dot(vec.edges[0]->cornNorm[1]);
1120 }
1121 else if (pcDotPhi > lenPhiZ)
1122 {
1123 //
1124 // ...and above in phi. Find distance to point (E)
1125 //
1126 G4double distOutPhi = pcDotPhi-lenPhiZ;
1127 distOut2 = distOutPhi*distOutPhi + distOutZ*distOutZ;
1128 G4ThreeVector pe = p - vec.edges[1]->corner[1];
1129 *normDist = pe.dot(vec.edges[1]->cornNorm[1]);
1130 }
1131 else
1132 {
1133 //
1134 // ...and inside in phi. Find distance to line (F)
1135 //
1136 distOut2 = distOutZ*distOutZ;
1137 G4ThreeVector pd = p - vec.edges[0]->corner[1];
1138 *normDist = pd.dot(vec.edgeNorm[1]);
1139 }
1140 }
1141 else
1142 {
1143 G4double lenPhiZ = lenPhi[0] + pcDotRZ*lenPhi[1];
1144 //
1145 // We are inside RZ bounds
1146 //
1147 if (pcDotPhi < -lenPhiZ)
1148 {
1149 //
1150 // ...and below in phi. Find distance to line (G)
1151 //
1152 G4double distOut = edgeNorm*(pcDotPhi+lenPhiZ);
1153 distOut2 = distOut*distOut;
1154 G4ThreeVector pd = p - vec.edges[0]->corner[1];
1155 *normDist = pd.dot(vec.edges[0]->normal);
1156 }
1157 else if (pcDotPhi > lenPhiZ)
1158 {
1159 //
1160 // ...and above in phi. Find distance to line (H)
1161 //
1162 G4double distOut = edgeNorm*(pcDotPhi-lenPhiZ);
1163 distOut2 = distOut*distOut;
1164 G4ThreeVector pe = p - vec.edges[1]->corner[1];
1165 *normDist = pe.dot(vec.edges[1]->normal);
1166 }
1167 else
1168 {
1169 //
1170 // Inside bounds! No penalty.
1171 //
1172 return std::fabs(distFaceNorm);
1173 }
1174 }
1175 return std::sqrt( distFaceNorm*distFaceNorm + distOut2 );
1176}
1177
1178
1179//
1180// Calculation of surface area of a triangle.
1181// At the same time a random point in the triangle is given
1182//
1184 G4ThreeVector p2,
1185 G4ThreeVector p3,
1186 G4ThreeVector *p4 )
1187{
1188 G4ThreeVector v, w;
1189
1190 v = p3 - p1;
1191 w = p1 - p2;
1192 G4double lambda1 = G4UniformRand();
1193 G4double lambda2 = lambda1*G4UniformRand();
1194
1195 *p4=p2 + lambda1*w + lambda2*v;
1196 return 0.5*(v.cross(w)).mag();
1197}
1198
1199
1200//
1201// GetPointOnPlane
1202//
1203// Auxiliary method for GetPointOnSurface()
1204//
1208 G4double *Area )
1209{
1210 G4double chose,aOne,aTwo;
1211 G4ThreeVector point1,point2;
1212 aOne = SurfaceTriangle(p0,p1,p2,&point1);
1213 aTwo = SurfaceTriangle(p2,p3,p0,&point2);
1214 *Area= aOne+aTwo;
1215
1216 chose = G4UniformRand()*(aOne+aTwo);
1217 if( (chose>=0.) && (chose < aOne) )
1218 {
1219 return (point1);
1220 }
1221 return (point2);
1222}
1223
1224
1225//
1226// SurfaceArea()
1227//
1229{
1230 if( fSurfaceArea==0. )
1231 {
1232 // Define the variables
1233 //
1234 G4double area,areas;
1235 G4ThreeVector point1;
1236 G4ThreeVector v1,v2,v3,v4;
1237 G4PolyhedraSideVec *vec = vecs;
1238 areas=0.;
1239
1240 // Do a loop on all SideEdge
1241 //
1242 do
1243 {
1244 // Define 4points for a Plane or Triangle
1245 //
1246 v1=vec->edges[0]->corner[0];
1247 v2=vec->edges[0]->corner[1];
1248 v3=vec->edges[1]->corner[1];
1249 v4=vec->edges[1]->corner[0];
1250 point1=GetPointOnPlane(v1,v2,v3,v4,&area);
1251 areas+=area;
1252 } while( ++vec < vecs + numSide);
1253
1254 fSurfaceArea=areas;
1255 }
1256 return fSurfaceArea;
1257}
1258
1259
1260//
1261// GetPointOnFace()
1262//
1264{
1265 // Define the variables
1266 //
1267 std::vector<G4double>areas;
1268 std::vector<G4ThreeVector>points;
1269 G4double area=0;
1270 G4double result1;
1271 G4ThreeVector point1;
1272 G4ThreeVector v1,v2,v3,v4;
1273 G4PolyhedraSideVec *vec = vecs;
1274
1275 // Do a loop on all SideEdge
1276 //
1277 do
1278 {
1279 // Define 4points for a Plane or Triangle
1280 //
1281 v1=vec->edges[0]->corner[0];
1282 v2=vec->edges[0]->corner[1];
1283 v3=vec->edges[1]->corner[1];
1284 v4=vec->edges[1]->corner[0];
1285 point1=GetPointOnPlane(v1,v2,v3,v4,&result1);
1286 points.push_back(point1);
1287 areas.push_back(result1);
1288 area+=result1;
1289 } while( ++vec < vecs+numSide );
1290
1291 // Choose randomly one of the surfaces and point on it
1292 //
1293 G4double chose = area*G4UniformRand();
1294 G4double Achose1,Achose2;
1295 Achose1=0;Achose2=0.;
1296 G4int i=0;
1297 do
1298 {
1299 Achose2+=areas[i];
1300 if(chose>=Achose1 && chose<Achose2)
1301 {
1302 point1=points[i] ; break;
1303 }
1304 i++; Achose1=Achose2;
1305 } while( i<numSide );
1306
1307 return point1;
1308}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4UniformRand()
Definition: Randomize.hh:53
double z() const
Hep3Vector unit() const
double phi() const
double x() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
double dot(const Hep3Vector &) const
double perp2() const
double mag() 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 ZHi() const
G4double ZLo() const
G4int LineHitsCone(const G4ThreeVector &p, const G4ThreeVector &v, G4double *s1, G4double *s2)
EInside Inside(const G4ThreeVector &p, G4double tolerance, G4double *bestDistance)
G4PolyhedraSideVec * vecs
G4double Distance(const G4ThreeVector &p, G4bool outgoing)
G4PolyhedraSideEdge * edges
G4ThreeVector GetPointOnPlane(G4ThreeVector p0, G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3, G4double *Area)
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, G4bool outgoing, G4double surfTolerance, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal, G4bool &allBehind)
G4PolyhedraSide & operator=(const G4PolyhedraSide &source)
G4double SurfaceTriangle(G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector *p4)
G4double GetPhi(const G4ThreeVector &p)
G4bool IntersectSidePlane(const G4ThreeVector &p, const G4ThreeVector &v, const G4PolyhedraSideVec &vec, G4double normSign, G4double surfTolerance, G4double &distance, G4double &distFromSurface)
G4PolyhedraSide(const G4PolyhedraSideRZ *prevRZ, const G4PolyhedraSideRZ *tail, const G4PolyhedraSideRZ *head, const G4PolyhedraSideRZ *nextRZ, G4int numSide, G4double phiStart, G4double phiTotal, G4bool phiIsOpen, G4bool isAllBehind=false)
G4ThreeVector GetPointOnFace()
G4double Extent(const G4ThreeVector axis)
G4int ClosestPhiSegment(G4double phi)
G4double lenPhi[2]
G4int PhiSegment(G4double phi)
G4int LineHitsSegments(const G4ThreeVector &p, const G4ThreeVector &v, G4int *i1, G4int *i2)
G4ThreeVector Normal(const G4ThreeVector &p, G4double *bestDistance)
G4double DistanceAway(const G4ThreeVector &p, const G4PolyhedraSideVec &vec, G4double *normDist)
G4double DistanceToOneSide(const G4ThreeVector &p, const G4PolyhedraSideVec &vec, G4double *normDist)
G4double SurfaceArea()
virtual ~G4PolyhedraSide()
G4IntersectingCone * cone
void CalculateExtent(const EAxis axis, const G4VoxelLimits &voxelLimit, const G4AffineTransform &tranform, G4SolidExtentList &extentList)
void CopyStuff(const G4PolyhedraSide &source)
void AddSurface(const G4ClippablePolygon &surface)
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:58
@ kInside
Definition: geomdefs.hh:58
@ kOutside
Definition: geomdefs.hh:58
@ kSurface
Definition: geomdefs.hh:58
#define DBL_MIN
Definition: templates.hh:75