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