Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4UPolyhedra.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 G4UPolyhedra wrapper class
27//
28// 31.10.13 G.Cosmo, CERN
29// --------------------------------------------------------------------
30
31#include "G4Polyhedra.hh"
32#include "G4UPolyhedra.hh"
33
34#if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
35
36#include "G4GeomTools.hh"
38#include "G4AffineTransform.hh"
40#include "G4BoundingEnvelope.hh"
41
42using namespace CLHEP;
43
44////////////////////////////////////////////////////////////////////////
45//
46// Constructor (GEANT3 style parameters)
47//
48// GEANT3 PGON radii are specified in the distance to the norm of each face.
49//
50G4UPolyhedra::G4UPolyhedra(const G4String& name,
51 G4double phiStart,
52 G4double phiTotal,
53 G4int numSide,
54 G4int numZPlanes,
55 const G4double zPlane[],
56 const G4double rInner[],
57 const G4double rOuter[] )
58 : Base_t(name, phiStart, phiTotal, numSide,
59 numZPlanes, zPlane, rInner, rOuter)
60{
61 fGenericPgon = false;
62 SetOriginalParameters();
63 wrStart = phiStart;
64 while (wrStart < 0)
65 {
66 wrStart += twopi;
67 }
68 wrDelta = phiTotal;
69 if (wrDelta <= 0. || wrDelta >= twopi*(1-DBL_EPSILON))
70 {
71 wrDelta = twopi;
72 }
73 wrNumSide = numSide;
74 G4double convertRad = 1./std::cos(0.5*wrDelta/wrNumSide);
75 rzcorners.resize(0);
76 for (G4int i=0; i<numZPlanes; ++i)
77 {
78 G4double z = zPlane[i];
79 G4double r = rOuter[i]*convertRad;
80 rzcorners.emplace_back(r,z);
81 }
82 for (G4int i=numZPlanes-1; i>=0; --i)
83 {
84 G4double z = zPlane[i];
85 G4double r = rInner[i]*convertRad;
86 rzcorners.emplace_back(r,z);
87 }
88 std::vector<G4int> iout;
90}
91
92
93////////////////////////////////////////////////////////////////////////
94//
95// Constructor (generic parameters)
96//
97G4UPolyhedra::G4UPolyhedra(const G4String& name,
98 G4double phiStart,
99 G4double phiTotal,
100 G4int numSide,
101 G4int numRZ,
102 const G4double r[],
103 const G4double z[] )
104 : Base_t(name, phiStart, phiTotal, numSide, numRZ, r, z)
105{
106 fGenericPgon = true;
107 SetOriginalParameters();
108 wrStart = phiStart;
109 while (wrStart < 0.)
110 {
111 wrStart += twopi;
112 }
113 wrDelta = phiTotal;
114 if (wrDelta <= 0. || wrDelta >= twopi*(1-DBL_EPSILON))
115 {
116 wrDelta = twopi;
117 }
118 wrNumSide = numSide;
119 rzcorners.resize(0);
120 for (G4int i=0; i<numRZ; ++i)
121 {
122 rzcorners.emplace_back(r[i],z[i]);
123 }
124 std::vector<G4int> iout;
126}
127
128
129////////////////////////////////////////////////////////////////////////
130//
131// Copy constructor
132//
133G4UPolyhedra::G4UPolyhedra( const G4UPolyhedra& source )
134 : Base_t( source )
135{
136 fGenericPgon = source.fGenericPgon;
137 fOriginalParameters = source.fOriginalParameters;
138 wrStart = source.wrStart;
139 wrDelta = source.wrDelta;
140 wrNumSide = source.wrNumSide;
141 rzcorners = source.rzcorners;
142}
143
144
145////////////////////////////////////////////////////////////////////////
146//
147// Assignment operator
148//
149G4UPolyhedra& G4UPolyhedra::operator=( const G4UPolyhedra& source )
150{
151 if (this == &source) return *this;
152
153 Base_t::operator=( source );
154 fGenericPgon = source.fGenericPgon;
155 fOriginalParameters = source.fOriginalParameters;
156 wrStart = source.wrStart;
157 wrDelta = source.wrDelta;
158 wrNumSide = source.wrNumSide;
159 rzcorners = source.rzcorners;
160
161 return *this;
162}
163
164
165////////////////////////////////////////////////////////////////////////
166//
167// Accessors & modifiers
168//
169G4int G4UPolyhedra::GetNumSide() const
170{
171 return wrNumSide;
172}
173G4double G4UPolyhedra::GetStartPhi() const
174{
175 return wrStart;
176}
177G4double G4UPolyhedra::GetEndPhi() const
178{
179 return (wrStart + wrDelta);
180}
181G4double G4UPolyhedra::GetSinStartPhi() const
182{
183 G4double phi = GetStartPhi();
184 return std::sin(phi);
185}
186G4double G4UPolyhedra::GetCosStartPhi() const
187{
188 G4double phi = GetStartPhi();
189 return std::cos(phi);
190}
191G4double G4UPolyhedra::GetSinEndPhi() const
192{
193 G4double phi = GetEndPhi();
194 return std::sin(phi);
195}
196G4double G4UPolyhedra::GetCosEndPhi() const
197{
198 G4double phi = GetEndPhi();
199 return std::cos(phi);
200}
201G4bool G4UPolyhedra::IsOpen() const
202{
203 return (wrDelta < twopi);
204}
205G4bool G4UPolyhedra::IsGeneric() const
206{
207 return fGenericPgon;
208}
209G4int G4UPolyhedra::GetNumRZCorner() const
210{
211 return rzcorners.size();
212}
213G4PolyhedraSideRZ G4UPolyhedra::GetCorner(G4int index) const
214{
215 G4TwoVector rz = rzcorners.at(index);
216 G4PolyhedraSideRZ psiderz = { rz.x(), rz.y() };
217
218 return psiderz;
219}
220G4PolyhedraHistorical* G4UPolyhedra::GetOriginalParameters() const
221{
222 return new G4PolyhedraHistorical(fOriginalParameters);
223}
224void G4UPolyhedra::SetOriginalParameters()
225{
226 G4double startPhi = GetPhiStart();
227 G4double deltaPhi = GetPhiDelta();
228 G4int numPlanes = GetZSegmentCount() + 1;
229 G4int numSides = GetSideCount();
230
231 fOriginalParameters.Start_angle = startPhi;
232 fOriginalParameters.Opening_angle = deltaPhi;
233 fOriginalParameters.Num_z_planes = numPlanes;
234 fOriginalParameters.numSide = numSides;
235
236 delete [] fOriginalParameters.Z_values;
237 delete [] fOriginalParameters.Rmin;
238 delete [] fOriginalParameters.Rmax;
239 fOriginalParameters.Z_values = new G4double[numPlanes];
240 fOriginalParameters.Rmin = new G4double[numPlanes];
241 fOriginalParameters.Rmax = new G4double[numPlanes];
242
243 G4double convertRad = fGenericPgon
244 ? 1.0 : std::cos(0.5*deltaPhi/numSides);
245 for (G4int i=0; i<numPlanes; ++i)
246 {
247 fOriginalParameters.Z_values[i] = GetZPlanes()[i];
248 fOriginalParameters.Rmax[i] = GetRMax()[i]/convertRad;
249 fOriginalParameters.Rmin[i] = GetRMin()[i]/convertRad;
250 }
251}
252void G4UPolyhedra::SetOriginalParameters(G4PolyhedraHistorical* pars)
253{
254 fOriginalParameters = *pars;
255 fRebuildPolyhedron = true;
256 Reset();
257}
258
259G4bool G4UPolyhedra::Reset()
260{
261 if (fGenericPgon)
262 {
263 std::ostringstream message;
264 message << "Solid " << GetName() << " built using generic construct."
265 << G4endl << "Not applicable to the generic construct !";
266 G4Exception("G4UPolyhedra::Reset()", "GeomSolids1001",
267 JustWarning, message, "Parameters NOT reset.");
268 return true; // error code set
269 }
270
271 //
272 // Rebuild polyhedra based on original parameters
273 //
274 wrStart = fOriginalParameters.Start_angle;
275 while (wrStart < 0.)
276 {
277 wrStart += twopi;
278 }
279 wrDelta = fOriginalParameters.Opening_angle;
280 if (wrDelta <= 0. || wrDelta >= twopi*(1-DBL_EPSILON))
281 {
282 wrDelta = twopi;
283 }
284 wrNumSide = fOriginalParameters.numSide;
285 rzcorners.resize(0);
286 for (G4int i=0; i<fOriginalParameters.Num_z_planes; ++i)
287 {
288 G4double z = fOriginalParameters.Z_values[i];
289 G4double r = fOriginalParameters.Rmax[i];
290 rzcorners.emplace_back(r,z);
291 }
292 for (G4int i=fOriginalParameters.Num_z_planes-1; i>=0; --i)
293 {
294 G4double z = fOriginalParameters.Z_values[i];
295 G4double r = fOriginalParameters.Rmin[i];
296 rzcorners.emplace_back(r,z);
297 }
298 std::vector<G4int> iout;
300
301 return false; // error code unset
302}
303
304
305////////////////////////////////////////////////////////////////////////
306//
307// Dispatch to parameterisation for replication mechanism dimension
308// computation & modification.
309//
310void G4UPolyhedra::ComputeDimensions(G4VPVParameterisation* p,
311 const G4int n,
312 const G4VPhysicalVolume* pRep)
313{
314 p->ComputeDimensions(*(G4Polyhedra*)this,n,pRep);
315}
316
317
318//////////////////////////////////////////////////////////////////////////
319//
320// Make a clone of the object
321
322G4VSolid* G4UPolyhedra::Clone() const
323{
324 return new G4UPolyhedra(*this);
325}
326
327
328//////////////////////////////////////////////////////////////////////////
329//
330// Get bounding box
331
332void G4UPolyhedra::BoundingLimits(G4ThreeVector& pMin,
333 G4ThreeVector& pMax) const
334{
335 static G4bool checkBBox = true;
336 static G4bool checkPhi = true;
337
338 G4double rmin = kInfinity, rmax = -kInfinity;
339 G4double zmin = kInfinity, zmax = -kInfinity;
340 for (G4int i=0; i<GetNumRZCorner(); ++i)
341 {
342 G4PolyhedraSideRZ corner = GetCorner(i);
343 if (corner.r < rmin) rmin = corner.r;
344 if (corner.r > rmax) rmax = corner.r;
345 if (corner.z < zmin) zmin = corner.z;
346 if (corner.z > zmax) zmax = corner.z;
347 }
348
349 G4double sphi = GetStartPhi();
350 G4double ephi = GetEndPhi();
351 G4double dphi = IsOpen() ? ephi-sphi : twopi;
352 G4int ksteps = GetNumSide();
353 G4double astep = dphi/ksteps;
354 G4double sinStep = std::sin(astep);
355 G4double cosStep = std::cos(astep);
356
357 G4double sinCur = GetSinStartPhi();
358 G4double cosCur = GetCosStartPhi();
359 if (!IsOpen()) rmin = 0.;
360 G4double xmin = rmin*cosCur, xmax = xmin;
361 G4double ymin = rmin*sinCur, ymax = ymin;
362 for (G4int k=0; k<ksteps+1; ++k)
363 {
364 G4double x = rmax*cosCur;
365 if (x < xmin) xmin = x;
366 if (x > xmax) xmax = x;
367 G4double y = rmax*sinCur;
368 if (y < ymin) ymin = y;
369 if (y > ymax) ymax = y;
370 if (rmin > 0.)
371 {
372 G4double xx = rmin*cosCur;
373 if (xx < xmin) xmin = xx;
374 if (xx > xmax) xmax = xx;
375 G4double yy = rmin*sinCur;
376 if (yy < ymin) ymin = yy;
377 if (yy > ymax) ymax = yy;
378 }
379 G4double sinTmp = sinCur;
380 sinCur = sinCur*cosStep + cosCur*sinStep;
381 cosCur = cosCur*cosStep - sinTmp*sinStep;
382 }
383 pMin.set(xmin,ymin,zmin);
384 pMax.set(xmax,ymax,zmax);
385
386 // Check correctness of the bounding box
387 //
388 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
389 {
390 std::ostringstream message;
391 message << "Bad bounding box (min >= max) for solid: "
392 << GetName() << " !"
393 << "\npMin = " << pMin
394 << "\npMax = " << pMax;
395 G4Exception("G4UPolyhedra::BoundingLimits()", "GeomMgt0001",
396 JustWarning, message);
397 StreamInfo(G4cout);
398 }
399
400 // Check consistency of bounding boxes
401 //
402 if (checkBBox)
403 {
404 U3Vector vmin, vmax;
405 Extent(vmin,vmax);
406 if (std::abs(pMin.x()-vmin.x()) > kCarTolerance ||
407 std::abs(pMin.y()-vmin.y()) > kCarTolerance ||
408 std::abs(pMin.z()-vmin.z()) > kCarTolerance ||
409 std::abs(pMax.x()-vmax.x()) > kCarTolerance ||
410 std::abs(pMax.y()-vmax.y()) > kCarTolerance ||
411 std::abs(pMax.z()-vmax.z()) > kCarTolerance)
412 {
413 std::ostringstream message;
414 message << "Inconsistency in bounding boxes for solid: "
415 << GetName() << " !"
416 << "\nBBox min: wrapper = " << pMin << " solid = " << vmin
417 << "\nBBox max: wrapper = " << pMax << " solid = " << vmax;
418 G4Exception("G4UPolyhedra::BoundingLimits()", "GeomMgt0001",
419 JustWarning, message);
420 checkBBox = false;
421 }
422 }
423
424 // Check consistency of angles
425 //
426 if (checkPhi)
427 {
428 if (GetStartPhi() != GetPhiStart() ||
429 GetEndPhi() != GetPhiEnd() ||
430 GetNumSide() != GetSideCount() ||
431 IsOpen() != (Base_t::GetPhiDelta() < twopi))
432 {
433 std::ostringstream message;
434 message << "Inconsistency in Phi angles or # of sides for solid: "
435 << GetName() << " !"
436 << "\nPhi start : wrapper = " << GetStartPhi()
437 << " solid = " << GetPhiStart()
438 << "\nPhi end : wrapper = " << GetEndPhi()
439 << " solid = " << GetPhiEnd()
440 << "\nPhi # sides: wrapper = " << GetNumSide()
441 << " solid = " << GetSideCount()
442 << "\nPhi is open: wrapper = " << (IsOpen() ? "true" : "false")
443 << " solid = "
444 << ((Base_t::GetPhiDelta() < twopi) ? "true" : "false");
445 G4Exception("G4UPolyhedra::BoundingLimits()", "GeomMgt0001",
446 JustWarning, message);
447 checkPhi = false;
448 }
449 }
450}
451
452//////////////////////////////////////////////////////////////////////////
453//
454// Calculate extent under transform and specified limit
455
456G4bool
457G4UPolyhedra::CalculateExtent(const EAxis pAxis,
458 const G4VoxelLimits& pVoxelLimit,
459 const G4AffineTransform& pTransform,
460 G4double& pMin, G4double& pMax) const
461{
462 G4ThreeVector bmin, bmax;
463 G4bool exist;
464
465 // Check bounding box (bbox)
466 //
467 BoundingLimits(bmin,bmax);
468 G4BoundingEnvelope bbox(bmin,bmax);
469#ifdef G4BBOX_EXTENT
470 if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
471#endif
472 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
473 {
474 return exist = pMin < pMax;
475 }
476
477 // To find the extent, RZ contour of the polycone is subdivided
478 // in triangles. The extent is calculated as cumulative extent of
479 // all sub-polycones formed by rotation of triangles around Z
480 //
481 G4TwoVectorList contourRZ;
482 G4TwoVectorList triangles;
483 std::vector<G4int> iout;
484 G4double eminlim = pVoxelLimit.GetMinExtent(pAxis);
485 G4double emaxlim = pVoxelLimit.GetMaxExtent(pAxis);
486
487 // get RZ contour, ensure anticlockwise order of corners
488 for (G4int i=0; i<GetNumRZCorner(); ++i)
489 {
490 G4PolyhedraSideRZ corner = GetCorner(i);
491 contourRZ.emplace_back(corner.r,corner.z);
492 }
494 G4double area = G4GeomTools::PolygonArea(contourRZ);
495 if (area < 0.) std::reverse(contourRZ.begin(),contourRZ.end());
496
497 // triangulate RZ countour
498 if (!G4GeomTools::TriangulatePolygon(contourRZ,triangles))
499 {
500 std::ostringstream message;
501 message << "Triangulation of RZ contour has failed for solid: "
502 << GetName() << " !"
503 << "\nExtent has been calculated using boundary box";
504 G4Exception("G4UPolyhedra::CalculateExtent()",
505 "GeomMgt1002",JustWarning,message);
506 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
507 }
508
509 // set trigonometric values
510 G4double sphi = GetStartPhi();
511 G4double ephi = GetEndPhi();
512 G4double dphi = IsOpen() ? ephi-sphi : twopi;
513 G4int ksteps = GetNumSide();
514 G4double astep = dphi/ksteps;
515 G4double sinStep = std::sin(astep);
516 G4double cosStep = std::cos(astep);
517 G4double sinStart = GetSinStartPhi();
518 G4double cosStart = GetCosStartPhi();
519
520 // allocate vector lists
521 std::vector<const G4ThreeVectorList *> polygons;
522 polygons.resize(ksteps+1);
523 for (G4int k=0; k<ksteps+1; ++k)
524 {
525 polygons[k] = new G4ThreeVectorList(3);
526 }
527
528 // main loop along triangles
529 pMin = kInfinity;
530 pMax = -kInfinity;
531 G4int ntria = triangles.size()/3;
532 for (G4int i=0; i<ntria; ++i)
533 {
534 G4double sinCur = sinStart;
535 G4double cosCur = cosStart;
536 G4int i3 = i*3;
537 for (G4int k=0; k<ksteps+1; ++k) // rotate triangle
538 {
539 auto ptr = const_cast<G4ThreeVectorList*>(polygons[k]);
540 auto iter = ptr->begin();
541 iter->set(triangles[i3+0].x()*cosCur,
542 triangles[i3+0].x()*sinCur,
543 triangles[i3+0].y());
544 iter++;
545 iter->set(triangles[i3+1].x()*cosCur,
546 triangles[i3+1].x()*sinCur,
547 triangles[i3+1].y());
548 iter++;
549 iter->set(triangles[i3+2].x()*cosCur,
550 triangles[i3+2].x()*sinCur,
551 triangles[i3+2].y());
552
553 G4double sinTmp = sinCur;
554 sinCur = sinCur*cosStep + cosCur*sinStep;
555 cosCur = cosCur*cosStep - sinTmp*sinStep;
556 }
557
558 // set sub-envelope and adjust extent
559 G4double emin,emax;
560 G4BoundingEnvelope benv(polygons);
561 if (!benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,emin,emax)) continue;
562 if (emin < pMin) pMin = emin;
563 if (emax > pMax) pMax = emax;
564 if (eminlim > pMin && emaxlim < pMax) break; // max possible extent
565 }
566 // free memory
567 for (G4int k=0; k<ksteps+1; ++k) { delete polygons[k]; polygons[k]=nullptr;}
568 return (pMin < pMax);
569}
570
571
572////////////////////////////////////////////////////////////////////////
573//
574// CreatePolyhedron
575//
576G4Polyhedron* G4UPolyhedra::CreatePolyhedron() const
577{
578 // Check the validity of the delta phi
579 G4double deltaPhi = wrDelta;
580 if (deltaPhi <= 0. || deltaPhi >= twopi*(1-DBL_EPSILON))
581 {
582 deltaPhi = twopi;
583 }
584 return new G4PolyhedronPgon(wrStart, deltaPhi, wrNumSide, rzcorners);
585}
586
587#endif // G4GEOM_USE_USOLIDS
const G4double kCarTolerance
std::vector< G4ThreeVector > G4ThreeVectorList
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::vector< G4TwoVector > G4TwoVectorList
CLHEP::Hep3Vector G4ThreeVector
CLHEP::Hep2Vector G4TwoVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double x() const
double y() const
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...
static void RemoveRedundantVertices(G4TwoVectorList &polygon, std::vector< G4int > &iout, G4double tolerance=0.0)
static G4double PolygonArea(const G4TwoVectorList &polygon)
static G4bool TriangulatePolygon(const G4TwoVectorList &polygon, std::vector< G4int > &result)
G4PolyhedraHistorical is a data structure for use in G4Polyhedra.
G4Polyhedra represents a composed closed polyhedra (PGON) made of planar sizes along the Z axis,...
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....
G4VSolid is an abstract base class for solids, physical shapes that can be tracked through....
Definition G4VSolid.hh:80
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
const char * name(G4int ptype)
#define DBL_EPSILON
Definition templates.hh:66