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