Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Polyhedra.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 G4Polyhedra, a CSG polyhedra,
27// as an inherited class of G4VCSGfaceted.
28//
29// Utility classes:
30// * G4EnclosingCylinder: decided a quick check of geometry would be a
31// good idea (for CPU speed). If the quick check fails, the regular
32// full-blown G4VCSGfaceted version is invoked.
33// * G4ReduciblePolygon: Really meant as a check of input parameters,
34// this utility class also "converts" the GEANT3-like PGON/PCON
35// arguments into the newer ones.
36// Both these classes are implemented outside this file because they are
37// shared with G4Polycone.
38//
39// Author: David C. Williams ([email protected])
40// --------------------------------------------------------------------
41
42#include "G4Polyhedra.hh"
43
44#if !defined(G4GEOM_USE_UPOLYHEDRA)
45
46#include "G4PolyhedraSide.hh"
47#include "G4PolyPhiFace.hh"
48
49#include "G4GeomTools.hh"
50#include "G4VoxelLimits.hh"
51#include "G4AffineTransform.hh"
52#include "G4BoundingEnvelope.hh"
53
54#include "G4QuickRand.hh"
55
57#include "G4ReduciblePolygon.hh"
59
60namespace
61{
62 G4Mutex surface_elementsMutex = G4MUTEX_INITIALIZER;
63}
64
65using namespace CLHEP;
66
67// Constructor (GEANT3 style parameters)
68//
69// GEANT3 PGON radii are specified in the distance to the norm of each face.
70//
72 G4double phiStart,
73 G4double thePhiTotal,
74 G4int theNumSide,
75 G4int numZPlanes,
76 const G4double zPlane[],
77 const G4double rInner[],
78 const G4double rOuter[] )
79 : G4VCSGfaceted( name )
80{
81 if (theNumSide <= 0)
82 {
83 std::ostringstream message;
84 message << "Solid must have at least one side - " << GetName() << G4endl
85 << " No sides specified !";
86 G4Exception("G4Polyhedra::G4Polyhedra()", "GeomSolids0002",
87 FatalErrorInArgument, message);
88 }
89
90 //
91 // Calculate conversion factor from G3 radius to G4 radius
92 //
93 G4double phiTotal = thePhiTotal;
94 if ( (phiTotal <=0) || (phiTotal >= twopi*(1-DBL_EPSILON)) )
95 { phiTotal = twopi; }
96 G4double convertRad = std::cos(0.5*phiTotal/theNumSide);
97
98 //
99 // Some historical stuff
100 //
102
103 original_parameters->numSide = theNumSide;
106 original_parameters->Num_z_planes = numZPlanes;
107 original_parameters->Z_values = new G4double[numZPlanes];
108 original_parameters->Rmin = new G4double[numZPlanes];
109 original_parameters->Rmax = new G4double[numZPlanes];
110
111 for (G4int i=0; i<numZPlanes; ++i)
112 {
113 if (( i < numZPlanes-1) && ( zPlane[i] == zPlane[i+1] ))
114 {
115 if( (rInner[i] > rOuter[i+1])
116 ||(rInner[i+1] > rOuter[i]) )
117 {
118 DumpInfo();
119 std::ostringstream message;
120 message << "Cannot create a Polyhedra with no contiguous segments."
121 << G4endl
122 << " Segments are not contiguous !" << G4endl
123 << " rMin[" << i << "] = " << rInner[i]
124 << " -- rMax[" << i+1 << "] = " << rOuter[i+1] << G4endl
125 << " rMin[" << i+1 << "] = " << rInner[i+1]
126 << " -- rMax[" << i << "] = " << rOuter[i];
127 G4Exception("G4Polyhedra::G4Polyhedra()", "GeomSolids0002",
128 FatalErrorInArgument, message);
129 }
130 }
131 original_parameters->Z_values[i] = zPlane[i];
132 original_parameters->Rmin[i] = rInner[i]/convertRad;
133 original_parameters->Rmax[i] = rOuter[i]/convertRad;
134 }
135
136
137 //
138 // Build RZ polygon using special PCON/PGON GEANT3 constructor
139 //
141 new G4ReduciblePolygon( rInner, rOuter, zPlane, numZPlanes );
142 rz->ScaleA( 1/convertRad );
143
144 //
145 // Do the real work
146 //
147 Create( phiStart, phiTotal, theNumSide, rz );
148
149 delete rz;
150}
151
152// Constructor (generic parameters)
153//
155 G4double phiStart,
156 G4double phiTotal,
157 G4int theNumSide,
158 G4int numRZ,
159 const G4double r[],
160 const G4double z[] )
161 : G4VCSGfaceted( name ), genericPgon(true)
162{
163 if (theNumSide <= 0)
164 {
165 std::ostringstream message;
166 message << "Solid must have at least one side - " << GetName() << G4endl
167 << " No sides specified !";
168 G4Exception("G4Polyhedra::G4Polyhedra()", "GeomSolids0002",
169 FatalErrorInArgument, message);
170 }
171
172 G4ReduciblePolygon* rz = new G4ReduciblePolygon( r, z, numRZ );
173
174 Create( phiStart, phiTotal, theNumSide, rz );
175
176 // Set original_parameters struct for consistency
177 //
179
180 delete rz;
181}
182
183// Create
184//
185// Generic create routine, called by each constructor
186// after conversion of arguments
187//
189 G4double phiTotal,
190 G4int theNumSide,
192{
193 //
194 // Perform checks of rz values
195 //
196 if (rz->Amin() < 0.0)
197 {
198 std::ostringstream message;
199 message << "Illegal input parameters - " << GetName() << G4endl
200 << " All R values must be >= 0 !";
201 G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
202 FatalErrorInArgument, message);
203 }
204
205 G4double rzArea = rz->Area();
206 if (rzArea < -kCarTolerance)
207 {
208 rz->ReverseOrder();
209 }
210 else if (rzArea < kCarTolerance)
211 {
212 std::ostringstream message;
213 message << "Illegal input parameters - " << GetName() << G4endl
214 << " R/Z cross section is zero or near zero: " << rzArea;
215 G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
216 FatalErrorInArgument, message);
217 }
218
221 {
222 std::ostringstream message;
223 message << "Illegal input parameters - " << GetName() << G4endl
224 << " Too few unique R/Z values !";
225 G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
226 FatalErrorInArgument, message);
227 }
228
229 if (rz->CrossesItself( 1/kInfinity ))
230 {
231 std::ostringstream message;
232 message << "Illegal input parameters - " << GetName() << G4endl
233 << " R/Z segments cross !";
234 G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
235 FatalErrorInArgument, message);
236 }
237
238 numCorner = rz->NumVertices();
239
240
241 startPhi = phiStart;
242 while( startPhi < 0 ) // Loop checking, 13.08.2015, G.Cosmo
243 startPhi += twopi;
244 //
245 // Phi opening? Account for some possible roundoff, and interpret
246 // nonsense value as representing no phi opening
247 //
248 if ( (phiTotal <= 0) || (phiTotal > twopi*(1-DBL_EPSILON)) )
249 {
250 phiIsOpen = false;
251 endPhi = startPhi + twopi;
252 }
253 else
254 {
255 phiIsOpen = true;
256 endPhi = startPhi + phiTotal;
257 }
258
259 //
260 // Save number sides
261 //
262 numSide = theNumSide;
263
264 //
265 // Allocate corner array.
266 //
268
269 //
270 // Copy corners
271 //
273
275 iterRZ.Begin();
276 do // Loop checking, 13.08.2015, G.Cosmo
277 {
278 next->r = iterRZ.GetA();
279 next->z = iterRZ.GetB();
280 } while( ++next, iterRZ.Next() );
281
282 //
283 // Allocate face pointer array
284 //
286 faces = new G4VCSGface*[numFace];
287
288 //
289 // Construct side faces
290 //
291 // To do so properly, we need to keep track of four successive RZ
292 // corners.
293 //
294 // But! Don't construct a face if both points are at zero radius!
295 //
296 G4PolyhedraSideRZ* corner = corners,
297 * prev = corners + numCorner-1,
298 * nextNext;
299 G4VCSGface** face = faces;
300 do // Loop checking, 13.08.2015, G.Cosmo
301 {
302 next = corner+1;
303 if (next >= corners+numCorner) next = corners;
304 nextNext = next+1;
305 if (nextNext >= corners+numCorner) nextNext = corners;
306
307 if (corner->r < 1/kInfinity && next->r < 1/kInfinity) continue;
308/*
309 // We must decide here if we can dare declare one of our faces
310 // as having a "valid" normal (i.e. allBehind = true). This
311 // is never possible if the face faces "inward" in r *unless*
312 // we have only one side
313 //
314 G4bool allBehind;
315 if ((corner->z > next->z) && (numSide > 1))
316 {
317 allBehind = false;
318 }
319 else
320 {
321 //
322 // Otherwise, it is only true if the line passing
323 // through the two points of the segment do not
324 // split the r/z cross section
325 //
326 allBehind = !rz->BisectedBy( corner->r, corner->z,
327 next->r, next->z, kCarTolerance );
328 }
329*/
330 *face++ = new G4PolyhedraSide( prev, corner, next, nextNext,
332 } while( prev=corner, corner=next, corner > corners );
333
334 if (phiIsOpen)
335 {
336 //
337 // Construct phi open edges
338 //
339 *face++ = new G4PolyPhiFace( rz, startPhi, phiTotal/numSide, endPhi );
340 *face++ = new G4PolyPhiFace( rz, endPhi, phiTotal/numSide, startPhi );
341 }
342
343 //
344 // We might have dropped a face or two: recalculate numFace
345 //
346 numFace = (G4int)(face-faces);
347
348 //
349 // Make enclosingCylinder
350 //
352 new G4EnclosingCylinder( rz, phiIsOpen, phiStart, phiTotal );
353}
354
355// Fake default constructor - sets only member data and allocates memory
356// for usage restricted to object persistency.
357//
359 : G4VCSGfaceted(a), startPhi(0.), endPhi(0.)
360{
361}
362
363// Destructor
364//
366{
367 delete [] corners;
368 delete original_parameters;
369 delete enclosingCylinder;
370 delete fElements;
371 delete fpPolyhedron;
372 corners = nullptr;
373 original_parameters = nullptr;
374 enclosingCylinder = nullptr;
375 fElements = nullptr;
376 fpPolyhedron = nullptr;
377}
378
379// Copy constructor
380//
382 : G4VCSGfaceted( source )
383{
384 CopyStuff( source );
385}
386
387// Assignment operator
388//
390{
391 if (this == &source) return *this;
392
393 G4VCSGfaceted::operator=( source );
394
395 delete [] corners;
396 delete original_parameters;
397 delete enclosingCylinder;
398
399 CopyStuff( source );
400
401 return *this;
402}
403
404// CopyStuff
405//
407{
408 //
409 // Simple stuff
410 //
411 numSide = source.numSide;
412 startPhi = source.startPhi;
413 endPhi = source.endPhi;
414 phiIsOpen = source.phiIsOpen;
415 numCorner = source.numCorner;
416 genericPgon= source.genericPgon;
417
418 //
419 // The corner array
420 //
422
424 * sourceCorn = source.corners;
425 do // Loop checking, 13.08.2015, G.Cosmo
426 {
427 *corn = *sourceCorn;
428 } while( ++sourceCorn, ++corn < corners+numCorner );
429
430 //
431 // Original parameters
432 //
433 if (source.original_parameters)
434 {
437 }
438
439 //
440 // Enclosing cylinder
441 //
443
444 //
445 // Surface elements
446 //
447 delete fElements;
448 fElements = nullptr;
449
450 //
451 // Polyhedron
452 //
453 fRebuildPolyhedron = false;
454 delete fpPolyhedron;
455 fpPolyhedron = nullptr;
456}
457
458// Reset
459//
460// Recalculates and reshapes the solid, given pre-assigned scaled
461// original_parameters.
462//
464{
465 if (genericPgon)
466 {
467 std::ostringstream message;
468 message << "Solid " << GetName() << " built using generic construct."
469 << G4endl << "Not applicable to the generic construct !";
470 G4Exception("G4Polyhedra::Reset()", "GeomSolids1001",
471 JustWarning, message, "Parameters NOT resetted.");
472 return true;
473 }
474
475 //
476 // Clear old setup
477 //
479 delete [] corners;
480 delete enclosingCylinder;
481 delete fElements;
482 corners = nullptr;
483 fElements = nullptr;
484 enclosingCylinder = nullptr;
485
486 //
487 // Rebuild polyhedra
488 //
497 delete rz;
498
499 return false;
500}
501
502// Inside
503//
504// This is an override of G4VCSGfaceted::Inside, created in order
505// to speed things up by first checking with G4EnclosingCylinder.
506//
508{
509 //
510 // Quick test
511 //
513
514 //
515 // Long answer
516 //
517 return G4VCSGfaceted::Inside(p);
518}
519
520// DistanceToIn
521//
522// This is an override of G4VCSGfaceted::Inside, created in order
523// to speed things up by first checking with G4EnclosingCylinder.
524//
526 const G4ThreeVector& v ) const
527{
528 //
529 // Quick test
530 //
532 return kInfinity;
533
534 //
535 // Long answer
536 //
537 return G4VCSGfaceted::DistanceToIn( p, v );
538}
539
540// DistanceToIn
541//
543{
545}
546
547// Get bounding box
548//
550 G4ThreeVector& pMax) const
551{
552 G4double rmin = kInfinity, rmax = -kInfinity;
553 G4double zmin = kInfinity, zmax = -kInfinity;
554 for (G4int i=0; i<GetNumRZCorner(); ++i)
555 {
556 G4PolyhedraSideRZ corner = GetCorner(i);
557 if (corner.r < rmin) rmin = corner.r;
558 if (corner.r > rmax) rmax = corner.r;
559 if (corner.z < zmin) zmin = corner.z;
560 if (corner.z > zmax) zmax = corner.z;
561 }
562
563 G4double sphi = GetStartPhi();
564 G4double ephi = GetEndPhi();
565 G4double dphi = IsOpen() ? ephi-sphi : twopi;
566 G4int ksteps = GetNumSide();
567 G4double astep = dphi/ksteps;
568 G4double sinStep = std::sin(astep);
569 G4double cosStep = std::cos(astep);
570
571 G4double sinCur = GetSinStartPhi();
572 G4double cosCur = GetCosStartPhi();
573 if (!IsOpen()) rmin = 0.;
574 G4double xmin = rmin*cosCur, xmax = xmin;
575 G4double ymin = rmin*sinCur, ymax = ymin;
576 for (G4int k=0; k<ksteps+1; ++k)
577 {
578 G4double x = rmax*cosCur;
579 if (x < xmin) xmin = x;
580 if (x > xmax) xmax = x;
581 G4double y = rmax*sinCur;
582 if (y < ymin) ymin = y;
583 if (y > ymax) ymax = y;
584 if (rmin > 0)
585 {
586 G4double xx = rmin*cosCur;
587 if (xx < xmin) xmin = xx;
588 if (xx > xmax) xmax = xx;
589 G4double yy = rmin*sinCur;
590 if (yy < ymin) ymin = yy;
591 if (yy > ymax) ymax = yy;
592 }
593 G4double sinTmp = sinCur;
594 sinCur = sinCur*cosStep + cosCur*sinStep;
595 cosCur = cosCur*cosStep - sinTmp*sinStep;
596 }
597 pMin.set(xmin,ymin,zmin);
598 pMax.set(xmax,ymax,zmax);
599
600 // Check correctness of the bounding box
601 //
602 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
603 {
604 std::ostringstream message;
605 message << "Bad bounding box (min >= max) for solid: "
606 << GetName() << " !"
607 << "\npMin = " << pMin
608 << "\npMax = " << pMax;
609 G4Exception("G4Polyhedra::BoundingLimits()", "GeomMgt0001",
610 JustWarning, message);
611 DumpInfo();
612 }
613}
614
615// Calculate extent under transform and specified limit
616//
618 const G4VoxelLimits& pVoxelLimit,
619 const G4AffineTransform& pTransform,
620 G4double& pMin, G4double& pMax) const
621{
622 G4ThreeVector bmin, bmax;
623 G4bool exist;
624
625 // Check bounding box (bbox)
626 //
627 BoundingLimits(bmin,bmax);
628 G4BoundingEnvelope bbox(bmin,bmax);
629#ifdef G4BBOX_EXTENT
630 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
631#endif
632 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
633 {
634 return exist = (pMin < pMax) ? true : false;
635 }
636
637 // To find the extent, RZ contour of the polycone is subdivided
638 // in triangles. The extent is calculated as cumulative extent of
639 // all sub-polycones formed by rotation of triangles around Z
640 //
641 G4TwoVectorList contourRZ;
642 G4TwoVectorList triangles;
643 std::vector<G4int> iout;
644 G4double eminlim = pVoxelLimit.GetMinExtent(pAxis);
645 G4double emaxlim = pVoxelLimit.GetMaxExtent(pAxis);
646
647 // get RZ contour, ensure anticlockwise order of corners
648 for (G4int i=0; i<GetNumRZCorner(); ++i)
649 {
650 G4PolyhedraSideRZ corner = GetCorner(i);
651 contourRZ.push_back(G4TwoVector(corner.r,corner.z));
652 }
654 G4double area = G4GeomTools::PolygonArea(contourRZ);
655 if (area < 0.) std::reverse(contourRZ.begin(),contourRZ.end());
656
657 // triangulate RZ countour
658 if (!G4GeomTools::TriangulatePolygon(contourRZ,triangles))
659 {
660 std::ostringstream message;
661 message << "Triangulation of RZ contour has failed for solid: "
662 << GetName() << " !"
663 << "\nExtent has been calculated using boundary box";
664 G4Exception("G4Polyhedra::CalculateExtent()",
665 "GeomMgt1002",JustWarning,message);
666 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
667 }
668
669 // set trigonometric values
670 G4double sphi = GetStartPhi();
671 G4double ephi = GetEndPhi();
672 G4double dphi = IsOpen() ? ephi-sphi : twopi;
673 G4int ksteps = GetNumSide();
674 G4double astep = dphi/ksteps;
675 G4double sinStep = std::sin(astep);
676 G4double cosStep = std::cos(astep);
677 G4double sinStart = GetSinStartPhi();
678 G4double cosStart = GetCosStartPhi();
679
680 // allocate vector lists
681 std::vector<const G4ThreeVectorList *> polygons;
682 polygons.resize(ksteps+1);
683 for (G4int k=0; k<ksteps+1; ++k)
684 {
685 polygons[k] = new G4ThreeVectorList(3);
686 }
687
688 // main loop along triangles
689 pMin = kInfinity;
690 pMax = -kInfinity;
691 G4int ntria = (G4int)triangles.size()/3;
692 for (G4int i=0; i<ntria; ++i)
693 {
694 G4double sinCur = sinStart;
695 G4double cosCur = cosStart;
696 G4int i3 = i*3;
697 for (G4int k=0; k<ksteps+1; ++k) // rotate triangle
698 {
699 G4ThreeVectorList* ptr = const_cast<G4ThreeVectorList*>(polygons[k]);
700 G4ThreeVectorList::iterator iter = ptr->begin();
701 iter->set(triangles[i3+0].x()*cosCur,
702 triangles[i3+0].x()*sinCur,
703 triangles[i3+0].y());
704 iter++;
705 iter->set(triangles[i3+1].x()*cosCur,
706 triangles[i3+1].x()*sinCur,
707 triangles[i3+1].y());
708 iter++;
709 iter->set(triangles[i3+2].x()*cosCur,
710 triangles[i3+2].x()*sinCur,
711 triangles[i3+2].y());
712
713 G4double sinTmp = sinCur;
714 sinCur = sinCur*cosStep + cosCur*sinStep;
715 cosCur = cosCur*cosStep - sinTmp*sinStep;
716 }
717
718 // set sub-envelope and adjust extent
719 G4double emin,emax;
720 G4BoundingEnvelope benv(polygons);
721 if (!benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,emin,emax)) continue;
722 if (emin < pMin) pMin = emin;
723 if (emax > pMax) pMax = emax;
724 if (eminlim > pMin && emaxlim < pMax) break; // max possible extent
725 }
726 // free memory
727 for (G4int k=0; k<ksteps+1; ++k) { delete polygons[k]; polygons[k]=0;}
728 return (pMin < pMax);
729}
730
731// ComputeDimensions
732//
734 const G4int n,
735 const G4VPhysicalVolume* pRep )
736{
737 p->ComputeDimensions(*this,n,pRep);
738}
739
740// GetEntityType
741//
743{
744 return G4String("G4Polyhedra");
745}
746
747// Make a clone of the object
748//
750{
751 return new G4Polyhedra(*this);
752}
753
754// Stream object contents to an output stream
755//
756std::ostream& G4Polyhedra::StreamInfo( std::ostream& os ) const
757{
758 G4long oldprc = os.precision(16);
759 os << "-----------------------------------------------------------\n"
760 << " *** Dump for solid - " << GetName() << " ***\n"
761 << " ===================================================\n"
762 << " Solid type: G4Polyhedra\n"
763 << " Parameters: \n"
764 << " starting phi angle : " << startPhi/degree << " degrees \n"
765 << " ending phi angle : " << endPhi/degree << " degrees \n"
766 << " number of sides : " << numSide << " \n";
767 G4int i=0;
768 if (!genericPgon)
769 {
771 os << " number of Z planes: " << numPlanes << "\n"
772 << " Z values: \n";
773 for (i=0; i<numPlanes; ++i)
774 {
775 os << " Z plane " << i << ": "
776 << original_parameters->Z_values[i] << "\n";
777 }
778 os << " Tangent distances to inner surface (Rmin): \n";
779 for (i=0; i<numPlanes; ++i)
780 {
781 os << " Z plane " << i << ": "
782 << original_parameters->Rmin[i] << "\n";
783 }
784 os << " Tangent distances to outer surface (Rmax): \n";
785 for (i=0; i<numPlanes; ++i)
786 {
787 os << " Z plane " << i << ": "
788 << original_parameters->Rmax[i] << "\n";
789 }
790 }
791 os << " number of RZ points: " << numCorner << "\n"
792 << " RZ values (corners): \n";
793 for (i=0; i<numCorner; ++i)
794 {
795 os << " "
796 << corners[i].r << ", " << corners[i].z << "\n";
797 }
798 os << "-----------------------------------------------------------\n";
799 os.precision(oldprc);
800
801 return os;
802}
803
804//////////////////////////////////////////////////////////////////////////
805//
806// Return volume
807
809{
810 if (fCubicVolume == 0.)
811 {
812 G4double total = 0.;
813 G4int nrz = GetNumRZCorner();
814 G4PolyhedraSideRZ a = GetCorner(nrz - 1);
815 for (G4int i=0; i<nrz; ++i)
816 {
818 total += (b.r*b.r + b.r*a.r + a.r*a.r)*(b.z - a.z);
819 a = b;
820 }
821 fCubicVolume = std::abs(total)*
822 std::sin((GetEndPhi() - GetStartPhi())/GetNumSide())*GetNumSide()/6.;
823 }
824 return fCubicVolume;
825}
826
827//////////////////////////////////////////////////////////////////////////
828//
829// Return surface area
830
832{
833 if (fSurfaceArea == 0.)
834 {
835 G4double total = 0.;
836 G4int nrz = GetNumRZCorner();
837 if (IsOpen())
838 {
839 G4PolyhedraSideRZ a = GetCorner(nrz - 1);
840 for (G4int i=0; i<nrz; ++i)
841 {
843 total += a.r*b.z - a.z*b.r;
844 a = b;
845 }
846 total = std::abs(total);
847 }
849 G4double cosa = std::cos(alp);
850 G4double sina = std::sin(alp);
851 G4PolyhedraSideRZ a = GetCorner(nrz - 1);
852 for (G4int i=0; i<nrz; ++i)
853 {
855 G4ThreeVector p1(a.r, 0, a.z);
856 G4ThreeVector p2(a.r*cosa, a.r*sina, a.z);
857 G4ThreeVector p3(b.r*cosa, b.r*sina, b.z);
858 G4ThreeVector p4(b.r, 0, b.z);
859 total += GetNumSide()*(G4GeomTools::QuadAreaNormal(p1, p2, p3, p4)).mag();
860 a = b;
861 }
862 fSurfaceArea = total;
863 }
864 return fSurfaceArea;
865}
866
867//////////////////////////////////////////////////////////////////////////
868//
869// Set vector of surface elements, auxiliary method for sampling
870// random points on surface
871
873{
874 fElements = new std::vector<G4Polyhedra::surface_element>;
875 G4double total = 0.;
876 G4int nrz = GetNumRZCorner();
877
878 // set lateral surface elements
879 G4double dphi = (GetEndPhi() - GetStartPhi())/GetNumSide();
880 G4double cosa = std::cos(dphi);
881 G4double sina = std::sin(dphi);
882 G4int ia = nrz - 1;
883 for (G4int ib=0; ib<nrz; ++ib)
884 {
888 selem.i0 = ia;
889 selem.i1 = ib;
890 ia = ib;
891 if (a.r == 0. && b.r == 0.) continue;
892 G4ThreeVector p1(a.r, 0, a.z);
893 G4ThreeVector p2(a.r*cosa, a.r*sina, a.z);
894 G4ThreeVector p3(b.r*cosa, b.r*sina, b.z);
895 G4ThreeVector p4(b.r, 0, b.z);
896 if (a.r > 0.)
897 {
898 selem.i2 = -1;
899 total += GetNumSide()*(G4GeomTools::TriangleAreaNormal(p1, p2, p3)).mag();
900 selem.area = total;
901 fElements->push_back(selem);
902 }
903 if (b.r > 0.)
904 {
905 selem.i2 = -2;
906 total += GetNumSide()*(G4GeomTools::TriangleAreaNormal(p1, p3, p4)).mag();
907 selem.area = total;
908 fElements->push_back(selem);
909 }
910 }
911
912 // set elements for phi cuts
913 if (IsOpen())
914 {
915 G4TwoVectorList contourRZ;
916 std::vector<G4int> triangles;
917 for (G4int i=0; i<nrz; ++i)
918 {
919 G4PolyhedraSideRZ corner = GetCorner(i);
920 contourRZ.push_back(G4TwoVector(corner.r, corner.z));
921 }
922 G4GeomTools::TriangulatePolygon(contourRZ, triangles);
923 G4int ntria = (G4int)triangles.size();
924 for (G4int i=0; i<ntria; i+=3)
925 {
927 selem.i0 = triangles[i];
928 selem.i1 = triangles[i+1];
929 selem.i2 = triangles[i+2];
930 G4PolyhedraSideRZ a = GetCorner(selem.i0);
931 G4PolyhedraSideRZ b = GetCorner(selem.i1);
932 G4PolyhedraSideRZ c = GetCorner(selem.i2);
933 G4double stria =
934 std::abs(G4GeomTools::TriangleArea(a.r, a.z, b.r, b.z, c.r, c.z));
935 total += stria;
936 selem.area = total;
937 fElements->push_back(selem); // start phi
938 total += stria;
939 selem.area = total;
940 selem.i0 += nrz;
941 fElements->push_back(selem); // end phi
942 }
943 }
944}
945
946//////////////////////////////////////////////////////////////////////////
947//
948// Generate random point on surface
949
951{
952 // Set surface elements
953 if (!fElements)
954 {
955 G4AutoLock l(&surface_elementsMutex);
957 l.unlock();
958 }
959
960 // Select surface element
962 selem = fElements->back();
963 G4double select = selem.area*G4QuickRand();
964 auto it = std::lower_bound(fElements->begin(), fElements->end(), select,
965 [](const G4Polyhedra::surface_element& x, G4double val)
966 -> G4bool { return x.area < val; });
967
968 // Generate random point
969 G4double x = 0, y = 0, z = 0;
970 G4double u = G4QuickRand();
971 G4double v = G4QuickRand();
972 if (u + v > 1.) { u = 1. - u; v = 1. - v; }
973 G4int i0 = (*it).i0;
974 G4int i1 = (*it).i1;
975 G4int i2 = (*it).i2;
976 if (i2 < 0) // lateral surface
977 {
978 // sample point
979 G4int nside = GetNumSide();
980 G4double dphi = (GetEndPhi() - GetStartPhi())/nside;
981 G4double cosa = std::cos(dphi);
982 G4double sina = std::sin(dphi);
985 G4ThreeVector p0(a.r, 0, a.z);
986 G4ThreeVector p1(b.r, 0, b.z);
987 G4ThreeVector p2(b.r*cosa, b.r*sina, b.z);
988 if (i2 == -1) p1.set(a.r*cosa, a.r*sina, a.z);
989 p0 += (p1 - p0)*u + (p2 - p0)*v;
990 // find selected side and rotate point
991 G4double scurr = (*it).area;
992 G4double sprev = (it == fElements->begin()) ? 0. : (*(--it)).area;
993 G4int iside = nside*(select - sprev)/(scurr - sprev);
994 if (iside == 0 && GetStartPhi() == 0.)
995 {
996 x = p0.x();
997 y = p0.y();
998 z = p0.z();
999 }
1000 else
1001 {
1002 if (iside == nside) --iside; // iside must be less then nside
1003 G4double phi = iside*dphi + GetStartPhi();
1004 G4double cosphi = std::cos(phi);
1005 G4double sinphi = std::sin(phi);
1006 x = p0.x()*cosphi - p0.y()*sinphi;
1007 y = p0.x()*sinphi + p0.y()*cosphi;
1008 z = p0.z();
1009 }
1010 }
1011 else // phi cut
1012 {
1013 G4int nrz = GetNumRZCorner();
1014 G4double phi = (i0 < nrz) ? GetStartPhi() : GetEndPhi();
1015 if (i0 >= nrz) { i0 -= nrz; }
1019 G4double r = (p1.r - p0.r)*u + (p2.r - p0.r)*v + p0.r;
1020 x = r*std::cos(phi);
1021 y = r*std::sin(phi);
1022 z = (p1.z - p0.z)*u + (p2.z - p0.z)*v + p0.z;
1023 }
1024 return G4ThreeVector(x, y, z);
1025}
1026
1027//////////////////////////////////////////////////////////////////////////
1028//
1029// CreatePolyhedron
1030
1032{
1033 std::vector<G4TwoVector> rz(numCorner);
1034 for (G4int i = 0; i < numCorner; ++i)
1035 rz[i].set(corners[i].r, corners[i].z);
1036 return new G4PolyhedronPgon(startPhi, endPhi - startPhi, numSide, rz);
1037}
1038
1039// SetOriginalParameters
1040//
1042{
1043 G4int numPlanes = numCorner;
1044 G4bool isConvertible = true;
1045 G4double Zmax=rz->Bmax();
1046 rz->StartWithZMin();
1047
1048 // Prepare vectors for storage
1049 //
1050 std::vector<G4double> Z;
1051 std::vector<G4double> Rmin;
1052 std::vector<G4double> Rmax;
1053
1054 G4int countPlanes=1;
1055 G4int icurr=0;
1056 G4int icurl=0;
1057
1058 // first plane Z=Z[0]
1059 //
1060 Z.push_back(corners[0].z);
1061 G4double Zprev=Z[0];
1062 if (Zprev == corners[1].z)
1063 {
1064 Rmin.push_back(corners[0].r);
1065 Rmax.push_back (corners[1].r);icurr=1;
1066 }
1067 else if (Zprev == corners[numPlanes-1].z)
1068 {
1069 Rmin.push_back(corners[numPlanes-1].r);
1070 Rmax.push_back (corners[0].r);
1071 icurl=numPlanes-1;
1072 }
1073 else
1074 {
1075 Rmin.push_back(corners[0].r);
1076 Rmax.push_back (corners[0].r);
1077 }
1078
1079 // next planes until last
1080 //
1081 G4int inextr=0, inextl=0;
1082 for (G4int i=0; i < numPlanes-2; ++i)
1083 {
1084 inextr=1+icurr;
1085 inextl=(icurl <= 0)? numPlanes-1 : icurl-1;
1086
1087 if((corners[inextr].z >= Zmax) & (corners[inextl].z >= Zmax)) { break; }
1088
1089 G4double Zleft = corners[inextl].z;
1090 G4double Zright = corners[inextr].z;
1091 if(Zright>Zleft)
1092 {
1093 Z.push_back(Zleft);
1094 countPlanes++;
1095 G4double difZr=corners[inextr].z - corners[icurr].z;
1096 G4double difZl=corners[inextl].z - corners[icurl].z;
1097
1098 if(std::fabs(difZl) < kCarTolerance)
1099 {
1100 if(std::fabs(difZr) < kCarTolerance)
1101 {
1102 Rmin.push_back(corners[inextl].r);
1103 Rmax.push_back(corners[icurr].r);
1104 }
1105 else
1106 {
1107 Rmin.push_back(corners[inextl].r);
1108 Rmax.push_back(corners[icurr].r + (Zleft-corners[icurr].z)/difZr
1109 *(corners[inextr].r - corners[icurr].r));
1110 }
1111 }
1112 else if (difZl >= kCarTolerance)
1113 {
1114 if(std::fabs(difZr) < kCarTolerance)
1115 {
1116 Rmin.push_back(corners[icurl].r);
1117 Rmax.push_back(corners[icurr].r);
1118 }
1119 else
1120 {
1121 Rmin.push_back(corners[icurl].r);
1122 Rmax.push_back(corners[icurr].r + (Zleft-corners[icurr].z)/difZr
1123 *(corners[inextr].r - corners[icurr].r));
1124 }
1125 }
1126 else
1127 {
1128 isConvertible=false; break;
1129 }
1130 icurl=(icurl == 0)? numPlanes-1 : icurl-1;
1131 }
1132 else if(std::fabs(Zright-Zleft)<kCarTolerance) // Zright=Zleft
1133 {
1134 Z.push_back(Zleft);
1135 ++countPlanes;
1136 ++icurr;
1137
1138 icurl=(icurl == 0)? numPlanes-1 : icurl-1;
1139
1140 Rmin.push_back(corners[inextl].r);
1141 Rmax.push_back (corners[inextr].r);
1142 }
1143 else // Zright<Zleft
1144 {
1145 Z.push_back(Zright);
1146 ++countPlanes;
1147
1148 G4double difZr=corners[inextr].z - corners[icurr].z;
1149 G4double difZl=corners[inextl].z - corners[icurl].z;
1150 if(std::fabs(difZr) < kCarTolerance)
1151 {
1152 if(std::fabs(difZl) < kCarTolerance)
1153 {
1154 Rmax.push_back(corners[inextr].r);
1155 Rmin.push_back(corners[icurr].r);
1156 }
1157 else
1158 {
1159 Rmin.push_back(corners[icurl].r + (Zright-corners[icurl].z)/difZl
1160 * (corners[inextl].r - corners[icurl].r));
1161 Rmax.push_back(corners[inextr].r);
1162 }
1163 ++icurr;
1164 } // plate
1165 else if (difZr >= kCarTolerance)
1166 {
1167 if(std::fabs(difZl) < kCarTolerance)
1168 {
1169 Rmax.push_back(corners[inextr].r);
1170 Rmin.push_back (corners[icurr].r);
1171 }
1172 else
1173 {
1174 Rmax.push_back(corners[inextr].r);
1175 Rmin.push_back (corners[icurl].r+(Zright-corners[icurl].z)/difZl
1176 * (corners[inextl].r - corners[icurl].r));
1177 }
1178 ++icurr;
1179 }
1180 else
1181 {
1182 isConvertible=false; break;
1183 }
1184 }
1185 } // end for loop
1186
1187 // last plane Z=Zmax
1188 //
1189 Z.push_back(Zmax);
1190 ++countPlanes;
1191 inextr=1+icurr;
1192 inextl=(icurl <= 0)? numPlanes-1 : icurl-1;
1193
1194 Rmax.push_back(corners[inextr].r);
1195 Rmin.push_back(corners[inextl].r);
1196
1197 // Set original parameters Rmin,Rmax,Z
1198 //
1199 if(isConvertible)
1200 {
1203 original_parameters->Z_values = new G4double[countPlanes];
1204 original_parameters->Rmin = new G4double[countPlanes];
1205 original_parameters->Rmax = new G4double[countPlanes];
1206
1207 for(G4int j=0; j < countPlanes; ++j)
1208 {
1210 original_parameters->Rmax[j] = Rmax[j];
1211 original_parameters->Rmin[j] = Rmin[j];
1212 }
1215 original_parameters->Num_z_planes = countPlanes;
1216
1217 }
1218 else // Set parameters(r,z) with Rmin==0 as convention
1219 {
1220#ifdef G4SPECSDEBUG
1221 std::ostringstream message;
1222 message << "Polyhedra " << GetName() << G4endl
1223 << "cannot be converted to Polyhedra with (Rmin,Rmaz,Z) parameters!";
1224 G4Exception("G4Polyhedra::SetOriginalParameters()",
1225 "GeomSolids0002", JustWarning, message);
1226#endif
1229 original_parameters->Z_values = new G4double[numPlanes];
1230 original_parameters->Rmin = new G4double[numPlanes];
1231 original_parameters->Rmax = new G4double[numPlanes];
1232
1233 for(G4int j=0; j < numPlanes; ++j)
1234 {
1237 original_parameters->Rmin[j] = 0.0;
1238 }
1241 original_parameters->Num_z_planes = numPlanes;
1242 }
1243}
1244
1245#endif
std::vector< G4ThreeVector > G4ThreeVectorList
@ JustWarning
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::vector< G4TwoVector > G4TwoVectorList
Definition: G4GeomTools.hh:42
G4double G4QuickRand()
Definition: G4QuickRand.hh:34
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
std::mutex G4Mutex
Definition: G4Threading.hh:81
CLHEP::Hep3Vector G4ThreeVector
CLHEP::Hep2Vector G4TwoVector
Definition: G4TwoVector.hh:36
double G4double
Definition: G4Types.hh:83
long G4long
Definition: G4Types.hh:87
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
double z() const
double x() const
double y() const
void set(double x, double y, double z)
G4bool BoundingBoxVsVoxelLimits(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4bool ShouldMiss(const G4ThreeVector &p, const G4ThreeVector &v) const
G4bool MustBeOutside(const G4ThreeVector &p) const
static G4ThreeVector QuadAreaNormal(const G4ThreeVector &A, const G4ThreeVector &B, const G4ThreeVector &C, const G4ThreeVector &D)
Definition: G4GeomTools.cc:609
static G4bool TriangulatePolygon(const G4TwoVectorList &polygon, G4TwoVectorList &result)
Definition: G4GeomTools.cc:193
static void RemoveRedundantVertices(G4TwoVectorList &polygon, std::vector< G4int > &iout, G4double tolerance=0.0)
Definition: G4GeomTools.cc:305
static G4double TriangleArea(G4double Ax, G4double Ay, G4double Bx, G4double By, G4double Cx, G4double Cy)
Definition: G4GeomTools.cc:41
static G4double PolygonArea(const G4TwoVectorList &polygon)
Definition: G4GeomTools.cc:76
static G4ThreeVector TriangleAreaNormal(const G4ThreeVector &A, const G4ThreeVector &B, const G4ThreeVector &C)
Definition: G4GeomTools.cc:598
G4ThreeVector GetPointOnSurface() const
Definition: G4Polyhedra.cc:950
G4Polyhedra & operator=(const G4Polyhedra &source)
Definition: G4Polyhedra.cc:389
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Polyhedra.cc:756
G4PolyhedraHistorical * original_parameters
Definition: G4Polyhedra.hh:189
G4double GetCubicVolume()
Definition: G4Polyhedra.cc:808
G4PolyhedraSideRZ * corners
Definition: G4Polyhedra.hh:188
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Polyhedra.cc:617
G4Polyhedra(const G4String &name, G4double phiStart, G4double phiTotal, G4int numSide, G4int numZPlanes, const G4double zPlane[], const G4double rInner[], const G4double rOuter[])
Definition: G4Polyhedra.cc:71
G4int GetNumRZCorner() const
G4bool genericPgon
Definition: G4Polyhedra.hh:186
void CopyStuff(const G4Polyhedra &source)
Definition: G4Polyhedra.cc:406
G4double GetSinStartPhi() const
G4VSolid * Clone() const
Definition: G4Polyhedra.cc:749
void Create(G4double phiStart, G4double phiTotal, G4int numSide, G4ReduciblePolygon *rz)
Definition: G4Polyhedra.cc:188
G4bool IsOpen() const
G4GeometryType GetEntityType() const
Definition: G4Polyhedra.cc:742
G4double endPhi
Definition: G4Polyhedra.hh:184
G4double GetSurfaceArea()
Definition: G4Polyhedra.cc:831
G4double GetEndPhi() const
std::vector< surface_element > * fElements
Definition: G4Polyhedra.hh:194
void SetOriginalParameters(G4PolyhedraHistorical *pars)
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Polyhedra.cc:733
G4bool phiIsOpen
Definition: G4Polyhedra.hh:185
G4double startPhi
Definition: G4Polyhedra.hh:183
G4int GetNumSide() const
G4int numCorner
Definition: G4Polyhedra.hh:187
virtual ~G4Polyhedra()
Definition: G4Polyhedra.cc:365
G4PolyhedraSideRZ GetCorner(const G4int index) const
G4bool Reset()
Definition: G4Polyhedra.cc:463
G4double GetStartPhi() const
G4EnclosingCylinder * enclosingCylinder
Definition: G4Polyhedra.hh:191
G4double GetCosStartPhi() const
void SetSurfaceElements() const
Definition: G4Polyhedra.cc:872
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Polyhedra.cc:525
EInside Inside(const G4ThreeVector &p) const
Definition: G4Polyhedra.cc:507
G4Polyhedron * CreatePolyhedron() const
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4Polyhedra.cc:549
G4double Amin() const
G4bool RemoveDuplicateVertices(G4double tolerance)
G4int NumVertices() const
G4bool RemoveRedundantVertices(G4double tolerance)
G4bool CrossesItself(G4double tolerance)
G4double Bmax() const
void ScaleA(G4double scale)
G4VCSGfaceted & operator=(const G4VCSGfaceted &source)
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
G4double fCubicVolume
virtual EInside Inside(const G4ThreeVector &p) const
G4VCSGface ** faces
G4Polyhedron * fpPolyhedron
G4bool fRebuildPolyhedron
G4double fSurfaceArea
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
G4String GetName() const
void DumpInfo() const
G4double kCarTolerance
Definition: G4VSolid.hh:299
G4double GetMinExtent(const EAxis pAxis) const
G4double GetMaxExtent(const EAxis pAxis) const
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:67
@ kOutside
Definition: geomdefs.hh:68
Definition: DoubConv.h:17
#define DBL_EPSILON
Definition: templates.hh:66