Geant4  10.02.p03
G4Polyhedra Class Reference

#include <G4Polyhedra.hh>

Inheritance diagram for G4Polyhedra:
Collaboration diagram for G4Polyhedra:

Public Member Functions

 G4Polyhedra (const G4String &name, G4double phiStart, G4double phiTotal, G4int numSide, G4int numZPlanes, const G4double zPlane[], const G4double rInner[], const G4double rOuter[])
 
 G4Polyhedra (const G4String &name, G4double phiStart, G4double phiTotal, G4int numSide, G4int numRZ, const G4double r[], const G4double z[])
 
virtual ~G4Polyhedra ()
 
EInside Inside (const G4ThreeVector &p) const
 
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const
 
G4double DistanceToIn (const G4ThreeVector &p) const
 
void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
G4GeometryType GetEntityType () const
 
G4VSolidClone () const
 
G4ThreeVector GetPointOnSurface () const
 
std::ostream & StreamInfo (std::ostream &os) const
 
G4PolyhedronCreatePolyhedron () const
 
G4bool Reset ()
 
G4int GetNumSide () const
 
G4double GetStartPhi () const
 
G4double GetEndPhi () const
 
G4bool IsOpen () const
 
G4bool IsGeneric () const
 
G4int GetNumRZCorner () const
 
G4PolyhedraSideRZ GetCorner (const G4int index) const
 
G4PolyhedraHistoricalGetOriginalParameters () const
 
void SetOriginalParameters (G4PolyhedraHistorical *pars)
 
 G4Polyhedra (__void__ &)
 
 G4Polyhedra (const G4Polyhedra &source)
 
G4Polyhedraoperator= (const G4Polyhedra &source)
 
- Public Member Functions inherited from G4VCSGfaceted
 G4VCSGfaceted (const G4String &name)
 
virtual ~G4VCSGfaceted ()
 
 G4VCSGfaceted (const G4VCSGfaceted &source)
 
G4VCSGfacetedoperator= (const G4VCSGfaceted &source)
 
virtual G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
 
virtual G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const
 
virtual G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
 
virtual G4double DistanceToOut (const G4ThreeVector &p) const
 
virtual void DescribeYourselfTo (G4VGraphicsScene &scene) const
 
virtual G4VisExtent GetExtent () const
 
virtual G4PolyhedronGetPolyhedron () const
 
G4int GetCubVolStatistics () const
 
G4double GetCubVolEpsilon () const
 
void SetCubVolStatistics (G4int st)
 
void SetCubVolEpsilon (G4double ep)
 
G4int GetAreaStatistics () const
 
G4double GetAreaAccuracy () const
 
void SetAreaStatistics (G4int st)
 
void SetAreaAccuracy (G4double ep)
 
virtual G4double GetCubicVolume ()
 
virtual G4double GetSurfaceArea ()
 
 G4VCSGfaceted (__void__ &)
 
- Public Member Functions inherited from G4VSolid
 G4VSolid (const G4String &name)
 
virtual ~G4VSolid ()
 
G4bool operator== (const G4VSolid &s) const
 
G4String GetName () const
 
void SetName (const G4String &name)
 
G4double GetTolerance () const
 
void DumpInfo () const
 
virtual const G4VSolidGetConstituentSolid (G4int no) const
 
virtual G4VSolidGetConstituentSolid (G4int no)
 
virtual const G4DisplacedSolidGetDisplacedSolidPtr () const
 
virtual G4DisplacedSolidGetDisplacedSolidPtr ()
 
 G4VSolid (__void__ &)
 
 G4VSolid (const G4VSolid &rhs)
 
G4VSolidoperator= (const G4VSolid &rhs)
 
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
 
G4double EstimateSurfaceArea (G4int nStat, G4double ell) const
 

Protected Member Functions

void SetOriginalParameters (G4ReduciblePolygon *rz)
 
void Create (G4double phiStart, G4double phiTotal, G4int numSide, G4ReduciblePolygon *rz)
 
void CopyStuff (const G4Polyhedra &source)
 
void DeleteStuff ()
 
G4ThreeVector GetPointOnPlane (G4ThreeVector p0, G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3) const
 
G4ThreeVector GetPointOnTriangle (G4ThreeVector p0, G4ThreeVector p1, G4ThreeVector p2) const
 
G4ThreeVector GetPointOnSurfaceCorners () const
 
- Protected Member Functions inherited from G4VCSGfaceted
virtual G4double DistanceTo (const G4ThreeVector &p, const G4bool outgoing) const
 
G4ThreeVector GetPointOnSurfaceGeneric () const
 
void CopyStuff (const G4VCSGfaceted &source)
 
void DeleteStuff ()
 
- Protected Member Functions inherited from G4VSolid
void CalculateClippedPolygonExtent (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipCrossSection (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipBetweenSections (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipPolygon (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis) const
 

Protected Attributes

G4int numSide
 
G4double startPhi
 
G4double endPhi
 
G4bool phiIsOpen
 
G4bool genericPgon
 
G4int numCorner
 
G4PolyhedraSideRZcorners
 
G4PolyhedraHistoricaloriginal_parameters
 
G4EnclosingCylinderenclosingCylinder
 
- Protected Attributes inherited from G4VCSGfaceted
G4int numFace
 
G4VCSGface ** faces
 
G4double fCubicVolume
 
G4double fSurfaceArea
 
G4bool fRebuildPolyhedron
 
G4PolyhedronfpPolyhedron
 
- Protected Attributes inherited from G4VSolid
G4double kCarTolerance
 

Detailed Description

Definition at line 81 of file G4Polyhedra.hh.

Constructor & Destructor Documentation

◆ G4Polyhedra() [1/4]

G4Polyhedra::G4Polyhedra ( const G4String name,
G4double  phiStart,
G4double  phiTotal,
G4int  numSide,
G4int  numZPlanes,
const G4double  zPlane[],
const G4double  rInner[],
const G4double  rOuter[] 
)

Definition at line 79 of file G4Polyhedra.cc.

87  : G4VCSGfaceted( name ), genericPgon(false)
88 {
89  if (theNumSide <= 0)
90  {
91  std::ostringstream message;
92  message << "Solid must have at least one side - " << GetName() << G4endl
93  << " No sides specified !";
94  G4Exception("G4Polyhedra::G4Polyhedra()", "GeomSolids0002",
95  FatalErrorInArgument, message);
96  }
97 
98  //
99  // Calculate conversion factor from G3 radius to G4 radius
100  //
101  G4double phiTotal = thePhiTotal;
102  if ( (phiTotal <=0) || (phiTotal >= twopi*(1-DBL_EPSILON)) )
103  { phiTotal = twopi; }
104  G4double convertRad = std::cos(0.5*phiTotal/theNumSide);
105 
106  //
107  // Some historical stuff
108  //
110 
111  original_parameters->numSide = theNumSide;
112  original_parameters->Start_angle = phiStart;
114  original_parameters->Num_z_planes = numZPlanes;
115  original_parameters->Z_values = new G4double[numZPlanes];
116  original_parameters->Rmin = new G4double[numZPlanes];
117  original_parameters->Rmax = new G4double[numZPlanes];
118 
119  G4int i;
120  for (i=0; i<numZPlanes; i++)
121  {
122  if (( i < numZPlanes-1) && ( zPlane[i] == zPlane[i+1] ))
123  {
124  if( (rInner[i] > rOuter[i+1])
125  ||(rInner[i+1] > rOuter[i]) )
126  {
127  DumpInfo();
128  std::ostringstream message;
129  message << "Cannot create a Polyhedra with no contiguous segments."
130  << G4endl
131  << " Segments are not contiguous !" << G4endl
132  << " rMin[" << i << "] = " << rInner[i]
133  << " -- rMax[" << i+1 << "] = " << rOuter[i+1] << G4endl
134  << " rMin[" << i+1 << "] = " << rInner[i+1]
135  << " -- rMax[" << i << "] = " << rOuter[i];
136  G4Exception("G4Polyhedra::G4Polyhedra()", "GeomSolids0002",
137  FatalErrorInArgument, message);
138  }
139  }
140  original_parameters->Z_values[i] = zPlane[i];
141  original_parameters->Rmin[i] = rInner[i]/convertRad;
142  original_parameters->Rmax[i] = rOuter[i]/convertRad;
143  }
144 
145 
146  //
147  // Build RZ polygon using special PCON/PGON GEANT3 constructor
148  //
149  G4ReduciblePolygon *rz =
150  new G4ReduciblePolygon( rInner, rOuter, zPlane, numZPlanes );
151  rz->ScaleA( 1/convertRad );
152 
153  //
154  // Do the real work
155  //
156  Create( phiStart, phiTotal, theNumSide, rz );
157 
158  delete rz;
159 }
G4PolyhedraHistorical * original_parameters
Definition: G4Polyhedra.hh:187
G4bool genericPgon
Definition: G4Polyhedra.hh:184
G4VCSGfaceted(const G4String &name)
int G4int
Definition: G4Types.hh:78
void Create(G4double phiStart, G4double phiTotal, G4int numSide, G4ReduciblePolygon *rz)
Definition: G4Polyhedra.cc:201
G4String GetName() const
void DumpInfo() const
void ScaleA(G4double scale)
static const double twopi
Definition: G4SIunits.hh:75
#define DBL_EPSILON
Definition: templates.hh:87
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ G4Polyhedra() [2/4]

G4Polyhedra::G4Polyhedra ( const G4String name,
G4double  phiStart,
G4double  phiTotal,
G4int  numSide,
G4int  numRZ,
const G4double  r[],
const G4double  z[] 
)

Definition at line 165 of file G4Polyhedra.cc.

172  : G4VCSGfaceted( name ), genericPgon(true)
173 {
174  if (theNumSide <= 0)
175  {
176  std::ostringstream message;
177  message << "Solid must have at least one side - " << GetName() << G4endl
178  << " No sides specified !";
179  G4Exception("G4Polyhedra::G4Polyhedra()", "GeomSolids0002",
180  FatalErrorInArgument, message);
181  }
182 
183  G4ReduciblePolygon *rz = new G4ReduciblePolygon( r, z, numRZ );
184 
185  Create( phiStart, phiTotal, theNumSide, rz );
186 
187  // Set original_parameters struct for consistency
188  //
190 
191  delete rz;
192 }
G4bool genericPgon
Definition: G4Polyhedra.hh:184
G4VCSGfaceted(const G4String &name)
void Create(G4double phiStart, G4double phiTotal, G4int numSide, G4ReduciblePolygon *rz)
Definition: G4Polyhedra.cc:201
G4String GetName() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
void SetOriginalParameters(G4PolyhedraHistorical *pars)
Here is the call graph for this function:

◆ ~G4Polyhedra()

G4Polyhedra::~G4Polyhedra ( )
virtual

Definition at line 390 of file G4Polyhedra.cc.

391 {
392  delete [] corners;
394 
395  delete enclosingCylinder;
396 }
G4PolyhedraHistorical * original_parameters
Definition: G4Polyhedra.hh:187
G4PolyhedraSideRZ * corners
Definition: G4Polyhedra.hh:186
G4EnclosingCylinder * enclosingCylinder
Definition: G4Polyhedra.hh:189

◆ G4Polyhedra() [3/4]

G4Polyhedra::G4Polyhedra ( __void__ &  a)

Definition at line 379 of file G4Polyhedra.cc.

380  : G4VCSGfaceted(a), numSide(0), startPhi(0.), endPhi(0.),
381  phiIsOpen(false), genericPgon(false), numCorner(0), corners(0),
383 {
384 }
G4PolyhedraHistorical * original_parameters
Definition: G4Polyhedra.hh:187
G4bool genericPgon
Definition: G4Polyhedra.hh:184
G4PolyhedraSideRZ * corners
Definition: G4Polyhedra.hh:186
G4VCSGfaceted(const G4String &name)
G4double endPhi
Definition: G4Polyhedra.hh:182
G4EnclosingCylinder * enclosingCylinder
Definition: G4Polyhedra.hh:189
G4int numCorner
Definition: G4Polyhedra.hh:185
G4double startPhi
Definition: G4Polyhedra.hh:181
G4bool phiIsOpen
Definition: G4Polyhedra.hh:183

◆ G4Polyhedra() [4/4]

G4Polyhedra::G4Polyhedra ( const G4Polyhedra source)

Definition at line 402 of file G4Polyhedra.cc.

403  : G4VCSGfaceted( source )
404 {
405  CopyStuff( source );
406 }
void CopyStuff(const G4Polyhedra &source)
Definition: G4Polyhedra.cc:432
G4VCSGfaceted(const G4String &name)
Here is the call graph for this function:

Member Function Documentation

◆ Clone()

G4VSolid * G4Polyhedra::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 591 of file G4Polyhedra.cc.

592 {
593  return new G4Polyhedra(*this);
594 }
G4Polyhedra(const G4String &name, G4double phiStart, G4double phiTotal, G4int numSide, G4int numZPlanes, const G4double zPlane[], const G4double rInner[], const G4double rOuter[])
Definition: G4Polyhedra.cc:79
Here is the call graph for this function:

◆ ComputeDimensions()

void G4Polyhedra::ComputeDimensions ( G4VPVParameterisation p,
const G4int  n,
const G4VPhysicalVolume pRep 
)
virtual

Reimplemented from G4VSolid.

Definition at line 571 of file G4Polyhedra.cc.

574 {
575  p->ComputeDimensions(*this,n,pRep);
576 }
Char_t n[5]
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
Here is the call graph for this function:

◆ CopyStuff()

void G4Polyhedra::CopyStuff ( const G4Polyhedra source)
protected

Definition at line 432 of file G4Polyhedra.cc.

433 {
434  //
435  // Simple stuff
436  //
437  numSide = source.numSide;
438  startPhi = source.startPhi;
439  endPhi = source.endPhi;
440  phiIsOpen = source.phiIsOpen;
441  numCorner = source.numCorner;
442  genericPgon= source.genericPgon;
443 
444  //
445  // The corner array
446  //
448 
449  G4PolyhedraSideRZ *corn = corners,
450  *sourceCorn = source.corners;
451  do // Loop checking, 13.08.2015, G.Cosmo
452  {
453  *corn = *sourceCorn;
454  } while( ++sourceCorn, ++corn < corners+numCorner );
455 
456  //
457  // Original parameters
458  //
459  if (source.original_parameters)
460  {
463  }
464 
465  //
466  // Enclosing cylinder
467  //
469 
470  fRebuildPolyhedron = false;
471  fpPolyhedron = 0;
472 }
G4bool fRebuildPolyhedron
G4Polyhedron * fpPolyhedron
G4PolyhedraHistorical * original_parameters
Definition: G4Polyhedra.hh:187
G4bool genericPgon
Definition: G4Polyhedra.hh:184
G4PolyhedraSideRZ * corners
Definition: G4Polyhedra.hh:186
G4double endPhi
Definition: G4Polyhedra.hh:182
G4EnclosingCylinder * enclosingCylinder
Definition: G4Polyhedra.hh:189
G4int numCorner
Definition: G4Polyhedra.hh:185
G4double startPhi
Definition: G4Polyhedra.hh:181
G4bool phiIsOpen
Definition: G4Polyhedra.hh:183
Here is the caller graph for this function:

◆ Create()

void G4Polyhedra::Create ( G4double  phiStart,
G4double  phiTotal,
G4int  numSide,
G4ReduciblePolygon rz 
)
protected

Definition at line 201 of file G4Polyhedra.cc.

205 {
206  //
207  // Perform checks of rz values
208  //
209  if (rz->Amin() < 0.0)
210  {
211  std::ostringstream message;
212  message << "Illegal input parameters - " << GetName() << G4endl
213  << " All R values must be >= 0 !";
214  G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
215  FatalErrorInArgument, message);
216  }
217 
218  G4double rzArea = rz->Area();
219  if (rzArea < -kCarTolerance)
220  {
221  rz->ReverseOrder();
222  }
223  else if (rzArea < kCarTolerance)
224  {
225  std::ostringstream message;
226  message << "Illegal input parameters - " << GetName() << G4endl
227  << " R/Z cross section is zero or near zero: " << rzArea;
228  G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
229  FatalErrorInArgument, message);
230  }
231 
233  || (!rz->RemoveRedundantVertices( kCarTolerance )) )
234  {
235  std::ostringstream message;
236  message << "Illegal input parameters - " << GetName() << G4endl
237  << " Too few unique R/Z values !";
238  G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
239  FatalErrorInArgument, message);
240  }
241 
242  if (rz->CrossesItself( 1/kInfinity ))
243  {
244  std::ostringstream message;
245  message << "Illegal input parameters - " << GetName() << G4endl
246  << " R/Z segments cross !";
247  G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
248  FatalErrorInArgument, message);
249  }
250 
251  numCorner = rz->NumVertices();
252 
253 
254  startPhi = phiStart;
255  while( startPhi < 0 ) // Loop checking, 13.08.2015, G.Cosmo
256  startPhi += twopi;
257  //
258  // Phi opening? Account for some possible roundoff, and interpret
259  // nonsense value as representing no phi opening
260  //
261  if ( (phiTotal <= 0) || (phiTotal > twopi*(1-DBL_EPSILON)) )
262  {
263  phiIsOpen = false;
264  endPhi = phiStart+twopi;
265  }
266  else
267  {
268  phiIsOpen = true;
269 
270  //
271  // Convert phi into our convention
272  //
273  endPhi = phiStart+phiTotal;
274  while( endPhi < startPhi ) // Loop checking, 13.08.2015, G.Cosmo
275  endPhi += twopi;
276  }
277 
278  //
279  // Save number sides
280  //
281  numSide = theNumSide;
282 
283  //
284  // Allocate corner array.
285  //
287 
288  //
289  // Copy corners
290  //
291  G4ReduciblePolygonIterator iterRZ(rz);
292 
293  G4PolyhedraSideRZ *next = corners;
294  iterRZ.Begin();
295  do // Loop checking, 13.08.2015, G.Cosmo
296  {
297  next->r = iterRZ.GetA();
298  next->z = iterRZ.GetB();
299  } while( ++next, iterRZ.Next() );
300 
301  //
302  // Allocate face pointer array
303  //
305  faces = new G4VCSGface*[numFace];
306 
307  //
308  // Construct side faces
309  //
310  // To do so properly, we need to keep track of four successive RZ
311  // corners.
312  //
313  // But! Don't construct a face if both points are at zero radius!
314  //
315  G4PolyhedraSideRZ *corner = corners,
316  *prev = corners + numCorner-1,
317  *nextNext;
318  G4VCSGface **face = faces;
319  do // Loop checking, 13.08.2015, G.Cosmo
320  {
321  next = corner+1;
322  if (next >= corners+numCorner) next = corners;
323  nextNext = next+1;
324  if (nextNext >= corners+numCorner) nextNext = corners;
325 
326  if (corner->r < 1/kInfinity && next->r < 1/kInfinity) continue;
327 /*
328  // We must decide here if we can dare declare one of our faces
329  // as having a "valid" normal (i.e. allBehind = true). This
330  // is never possible if the face faces "inward" in r *unless*
331  // we have only one side
332  //
333  G4bool allBehind;
334  if ((corner->z > next->z) && (numSide > 1))
335  {
336  allBehind = false;
337  }
338  else
339  {
340  //
341  // Otherwise, it is only true if the line passing
342  // through the two points of the segment do not
343  // split the r/z cross section
344  //
345  allBehind = !rz->BisectedBy( corner->r, corner->z,
346  next->r, next->z, kCarTolerance );
347  }
348 */
349  *face++ = new G4PolyhedraSide( prev, corner, next, nextNext,
351  } while( prev=corner, corner=next, corner > corners );
352 
353  if (phiIsOpen)
354  {
355  //
356  // Construct phi open edges
357  //
358  *face++ = new G4PolyPhiFace( rz, startPhi, phiTotal/numSide, endPhi );
359  *face++ = new G4PolyPhiFace( rz, endPhi, phiTotal/numSide, startPhi );
360  }
361 
362  //
363  // We might have dropped a face or two: recalculate numFace
364  //
365  numFace = face-faces;
366 
367  //
368  // Make enclosingCylinder
369  //
371  new G4EnclosingCylinder( rz, phiIsOpen, phiStart, phiTotal );
372 }
G4bool CrossesItself(G4double tolerance)
static const G4double kInfinity
Definition: geomdefs.hh:42
G4PolyhedraSideRZ * corners
Definition: G4Polyhedra.hh:186
G4String GetName() const
G4double endPhi
Definition: G4Polyhedra.hh:182
G4bool RemoveDuplicateVertices(G4double tolerance)
G4EnclosingCylinder * enclosingCylinder
Definition: G4Polyhedra.hh:189
G4bool RemoveRedundantVertices(G4double tolerance)
G4int numCorner
Definition: G4Polyhedra.hh:185
G4double startPhi
Definition: G4Polyhedra.hh:181
static const double twopi
Definition: G4SIunits.hh:75
#define DBL_EPSILON
Definition: templates.hh:87
G4VCSGface ** faces
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double Amin() const
G4int NumVertices() const
#define G4endl
Definition: G4ios.hh:61
G4double kCarTolerance
Definition: G4VSolid.hh:304
double G4double
Definition: G4Types.hh:76
G4bool phiIsOpen
Definition: G4Polyhedra.hh:183
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CreatePolyhedron()

G4Polyhedron * G4Polyhedra::CreatePolyhedron ( ) const
virtual

Creates user defined polyhedron. This function allows to the user to define arbitrary polyhedron. The faces of the polyhedron should be either triangles or planar quadrilateral. Nodes of a face are defined by indexes pointing to the elements in the xyz array. Numeration of the elements in the array starts from 1 (like in fortran). The indexes can be positive or negative. Negative sign means that the corresponding edge is invisible. The normal of the face should be directed to exterior of the polyhedron.

Parameters
Nnodesnumber of nodes
Nfacesnumber of faces
xyznodes
faces_vecfaces (quadrilaterals or triangles)
Returns
status of the operation - is non-zero in case of problem

Implements G4VCSGfaceted.

Definition at line 920 of file G4Polyhedra.cc.

921 {
922  if (!genericPgon)
923  {
931  }
932  else
933  {
934  // The following code prepares for:
935  // HepPolyhedron::createPolyhedron(int Nnodes, int Nfaces,
936  // const double xyz[][3],
937  // const int faces_vec[][4])
938  // Here is an extract from the header file HepPolyhedron.h:
956  G4int nNodes;
957  G4int nFaces;
958  typedef G4double double3[3];
959  double3* xyz;
960  typedef G4int int4[4];
961  int4* faces_vec;
962  if (phiIsOpen)
963  {
964  // Triangulate open ends. Simple ear-chopping algorithm...
965  // I'm not sure how robust this algorithm is (J.Allison).
966  //
967  std::vector<G4bool> chopped(numCorner, false);
968  std::vector<G4int*> triQuads;
969  G4int remaining = numCorner;
970  G4int iStarter = 0;
971  while (remaining >= 3) // Loop checking, 13.08.2015, G.Cosmo
972  {
973  // Find unchopped corners...
974  //
975  G4int A = -1, B = -1, C = -1;
976  G4int iStepper = iStarter;
977  do // Loop checking, 13.08.2015, G.Cosmo
978  {
979  if (A < 0) { A = iStepper; }
980  else if (B < 0) { B = iStepper; }
981  else if (C < 0) { C = iStepper; }
982  do // Loop checking, 13.08.2015, G.Cosmo
983  {
984  if (++iStepper >= numCorner) iStepper = 0;
985  }
986  while (chopped[iStepper]);
987  }
988  while (C < 0 && iStepper != iStarter);
989 
990  // Check triangle at B is pointing outward (an "ear").
991  // Sign of z cross product determines...
992 
993  G4double BAr = corners[A].r - corners[B].r;
994  G4double BAz = corners[A].z - corners[B].z;
995  G4double BCr = corners[C].r - corners[B].r;
996  G4double BCz = corners[C].z - corners[B].z;
997  if (BAr * BCz - BAz * BCr < kCarTolerance)
998  {
999  G4int* tq = new G4int[3];
1000  tq[0] = A + 1;
1001  tq[1] = B + 1;
1002  tq[2] = C + 1;
1003  triQuads.push_back(tq);
1004  chopped[B] = true;
1005  --remaining;
1006  }
1007  else
1008  {
1009  do // Loop checking, 13.08.2015, G.Cosmo
1010  {
1011  if (++iStarter >= numCorner) { iStarter = 0; }
1012  }
1013  while (chopped[iStarter]);
1014  }
1015  }
1016 
1017  // Transfer to faces...
1018 
1019  nNodes = (numSide + 1) * numCorner;
1020  nFaces = numSide * numCorner + 2 * triQuads.size();
1021  faces_vec = new int4[nFaces];
1022  G4int iface = 0;
1023  G4int addition = numCorner * numSide;
1024  G4int d = numCorner - 1;
1025  for (G4int iEnd = 0; iEnd < 2; ++iEnd)
1026  {
1027  for (size_t i = 0; i < triQuads.size(); ++i)
1028  {
1029  // Negative for soft/auxiliary/normally invisible edges...
1030  //
1031  G4int a, b, c;
1032  if (iEnd == 0)
1033  {
1034  a = triQuads[i][0];
1035  b = triQuads[i][1];
1036  c = triQuads[i][2];
1037  }
1038  else
1039  {
1040  a = triQuads[i][0] + addition;
1041  b = triQuads[i][2] + addition;
1042  c = triQuads[i][1] + addition;
1043  }
1044  G4int ab = std::abs(b - a);
1045  G4int bc = std::abs(c - b);
1046  G4int ca = std::abs(a - c);
1047  faces_vec[iface][0] = (ab == 1 || ab == d)? a: -a;
1048  faces_vec[iface][1] = (bc == 1 || bc == d)? b: -b;
1049  faces_vec[iface][2] = (ca == 1 || ca == d)? c: -c;
1050  faces_vec[iface][3] = 0;
1051  ++iface;
1052  }
1053  }
1054 
1055  // Continue with sides...
1056 
1057  xyz = new double3[nNodes];
1058  const G4double dPhi = (endPhi - startPhi) / numSide;
1059  G4double phi = startPhi;
1060  G4int ixyz = 0;
1061  for (G4int iSide = 0; iSide < numSide; ++iSide)
1062  {
1063  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1064  {
1065  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1066  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1067  xyz[ixyz][2] = corners[iCorner].z;
1068  if (iCorner < numCorner - 1)
1069  {
1070  faces_vec[iface][0] = ixyz + 1;
1071  faces_vec[iface][1] = ixyz + numCorner + 1;
1072  faces_vec[iface][2] = ixyz + numCorner + 2;
1073  faces_vec[iface][3] = ixyz + 2;
1074  }
1075  else
1076  {
1077  faces_vec[iface][0] = ixyz + 1;
1078  faces_vec[iface][1] = ixyz + numCorner + 1;
1079  faces_vec[iface][2] = ixyz + 2;
1080  faces_vec[iface][3] = ixyz - numCorner + 2;
1081  }
1082  ++iface;
1083  ++ixyz;
1084  }
1085  phi += dPhi;
1086  }
1087 
1088  // Last corners...
1089 
1090  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1091  {
1092  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1093  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1094  xyz[ixyz][2] = corners[iCorner].z;
1095  ++ixyz;
1096  }
1097  }
1098  else // !phiIsOpen - i.e., a complete 360 degrees.
1099  {
1100  nNodes = numSide * numCorner;
1101  nFaces = numSide * numCorner;;
1102  xyz = new double3[nNodes];
1103  faces_vec = new int4[nFaces];
1104  // const G4double dPhi = (endPhi - startPhi) / numSide;
1105  const G4double dPhi = twopi / numSide; // !phiIsOpen endPhi-startPhi = 360 degrees.
1106  G4double phi = startPhi;
1107  G4int ixyz = 0, iface = 0;
1108  for (G4int iSide = 0; iSide < numSide; ++iSide)
1109  {
1110  for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1111  {
1112  xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1113  xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1114  xyz[ixyz][2] = corners[iCorner].z;
1115  if (iSide < numSide - 1)
1116  {
1117  if (iCorner < numCorner - 1)
1118  {
1119  faces_vec[iface][0] = ixyz + 1;
1120  faces_vec[iface][1] = ixyz + numCorner + 1;
1121  faces_vec[iface][2] = ixyz + numCorner + 2;
1122  faces_vec[iface][3] = ixyz + 2;
1123  }
1124  else
1125  {
1126  faces_vec[iface][0] = ixyz + 1;
1127  faces_vec[iface][1] = ixyz + numCorner + 1;
1128  faces_vec[iface][2] = ixyz + 2;
1129  faces_vec[iface][3] = ixyz - numCorner + 2;
1130  }
1131  }
1132  else // Last side joins ends...
1133  {
1134  if (iCorner < numCorner - 1)
1135  {
1136  faces_vec[iface][0] = ixyz + 1;
1137  faces_vec[iface][1] = ixyz + numCorner - nFaces + 1;
1138  faces_vec[iface][2] = ixyz + numCorner - nFaces + 2;
1139  faces_vec[iface][3] = ixyz + 2;
1140  }
1141  else
1142  {
1143  faces_vec[iface][0] = ixyz + 1;
1144  faces_vec[iface][1] = ixyz - nFaces + numCorner + 1;
1145  faces_vec[iface][2] = ixyz - nFaces + 2;
1146  faces_vec[iface][3] = ixyz - numCorner + 2;
1147  }
1148  }
1149  ++ixyz;
1150  ++iface;
1151  }
1152  phi += dPhi;
1153  }
1154  }
1155  G4Polyhedron* polyhedron = new G4Polyhedron;
1156  G4int problem = polyhedron->createPolyhedron(nNodes, nFaces, xyz, faces_vec);
1157  delete [] faces_vec;
1158  delete [] xyz;
1159  if (problem)
1160  {
1161  std::ostringstream message;
1162  message << "Problem creating G4Polyhedron for: " << GetName();
1163  G4Exception("G4Polyhedra::CreatePolyhedron()", "GeomSolids1002",
1164  JustWarning, message);
1165  delete polyhedron;
1166  return 0;
1167  }
1168  else
1169  {
1170  return polyhedron;
1171  }
1172  }
1173 }
Float_t d
G4int createPolyhedron(G4int Nnodes, G4int Nfaces, const G4double xyz[][3], const G4int faces[][4])
G4PolyhedraHistorical * original_parameters
Definition: G4Polyhedra.hh:187
G4bool genericPgon
Definition: G4Polyhedra.hh:184
G4PolyhedraSideRZ * corners
Definition: G4Polyhedra.hh:186
double C(double temp)
int G4int
Definition: G4Types.hh:78
G4String GetName() const
G4double endPhi
Definition: G4Polyhedra.hh:182
double A(double temperature)
G4int numCorner
Definition: G4Polyhedra.hh:185
G4double startPhi
Definition: G4Polyhedra.hh:181
static const double twopi
Definition: G4SIunits.hh:75
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const G4double ab
G4double kCarTolerance
Definition: G4VSolid.hh:304
double G4double
Definition: G4Types.hh:76
G4bool phiIsOpen
Definition: G4Polyhedra.hh:183
Here is the call graph for this function:

◆ DeleteStuff()

void G4Polyhedra::DeleteStuff ( )
protected

◆ DistanceToIn() [1/2]

G4double G4Polyhedra::DistanceToIn ( const G4ThreeVector p,
const G4ThreeVector v 
) const
virtual

Reimplemented from G4VCSGfaceted.

Definition at line 543 of file G4Polyhedra.cc.

545 {
546  //
547  // Quick test
548  //
549  if (enclosingCylinder->ShouldMiss(p,v))
550  return kInfinity;
551 
552  //
553  // Long answer
554  //
555  return G4VCSGfaceted::DistanceToIn( p, v );
556 }
static const G4double kInfinity
Definition: geomdefs.hh:42
G4bool ShouldMiss(const G4ThreeVector &p, const G4ThreeVector &v) const
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
G4EnclosingCylinder * enclosingCylinder
Definition: G4Polyhedra.hh:189
Here is the call graph for this function:

◆ DistanceToIn() [2/2]

G4double G4Polyhedra::DistanceToIn ( const G4ThreeVector p) const
virtual

Reimplemented from G4VCSGfaceted.

Definition at line 562 of file G4Polyhedra.cc.

563 {
564  return G4VCSGfaceted::DistanceToIn(p);
565 }
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Here is the call graph for this function:

◆ GetCorner()

G4PolyhedraSideRZ G4Polyhedra::GetCorner ( const G4int  index) const
inline
Here is the caller graph for this function:

◆ GetEndPhi()

G4double G4Polyhedra::GetEndPhi ( ) const
inline
Here is the caller graph for this function:

◆ GetEntityType()

G4GeometryType G4Polyhedra::GetEntityType ( ) const
virtual

Reimplemented from G4VCSGfaceted.

Definition at line 582 of file G4Polyhedra.cc.

583 {
584  return G4String("G4Polyhedra");
585 }

◆ GetNumRZCorner()

G4int G4Polyhedra::GetNumRZCorner ( ) const
inline
Here is the caller graph for this function:

◆ GetNumSide()

G4int G4Polyhedra::GetNumSide ( ) const
inline
Here is the caller graph for this function:

◆ GetOriginalParameters()

G4PolyhedraHistorical* G4Polyhedra::GetOriginalParameters ( ) const
inline
Here is the caller graph for this function:

◆ GetPointOnPlane()

G4ThreeVector G4Polyhedra::GetPointOnPlane ( G4ThreeVector  p0,
G4ThreeVector  p1,
G4ThreeVector  p2,
G4ThreeVector  p3 
) const
protected

Definition at line 654 of file G4Polyhedra.cc.

656 {
657  G4double lambda1, lambda2, chose,aOne,aTwo;
658  G4ThreeVector t, u, v, w, Area, normal;
659  aOne = 1.;
660  aTwo = 1.;
661 
662  t = p1 - p0;
663  u = p2 - p1;
664  v = p3 - p2;
665  w = p0 - p3;
666 
667  chose = G4RandFlat::shoot(0.,aOne+aTwo);
668  if( (chose>=0.) && (chose < aOne) )
669  {
670  lambda1 = G4RandFlat::shoot(0.,1.);
671  lambda2 = G4RandFlat::shoot(0.,lambda1);
672  return (p2+lambda1*v+lambda2*w);
673  }
674 
675  lambda1 = G4RandFlat::shoot(0.,1.);
676  lambda2 = G4RandFlat::shoot(0.,lambda1);
677  return (p0+lambda1*t+lambda2*u);
678 }
ThreeVector shoot(const G4int Ap, const G4int Af)
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetPointOnSurface()

G4ThreeVector G4Polyhedra::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 703 of file G4Polyhedra.cc.

704 {
705  if( !genericPgon ) // Polyhedra by faces
706  {
707  G4int j, numPlanes = original_parameters->Num_z_planes, Flag=0;
708  G4double chose, totArea=0., Achose1, Achose2,
709  rad1, rad2, sinphi1, sinphi2, cosphi1, cosphi2;
710  G4double a, b, l2, rang, totalPhi, ksi,
711  area, aTop=0., aBottom=0., zVal=0.;
712 
713  G4ThreeVector p0, p1, p2, p3;
714  std::vector<G4double> aVector1;
715  std::vector<G4double> aVector2;
716  std::vector<G4double> aVector3;
717 
718  totalPhi= (phiIsOpen) ? (endPhi-startPhi) : twopi;
719  ksi = totalPhi/numSide;
720  G4double cosksi = std::cos(ksi/2.);
721 
722  // Below we generate the areas relevant to our solid
723  //
724  for(j=0; j<numPlanes-1; j++)
725  {
726  a = original_parameters->Rmax[j+1];
727  b = original_parameters->Rmax[j];
729  -original_parameters->Z_values[j+1]) + sqr(b-a);
730  area = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
731  aVector1.push_back(area);
732  }
733 
734  for(j=0; j<numPlanes-1; j++)
735  {
736  a = original_parameters->Rmin[j+1];//*cosksi;
737  b = original_parameters->Rmin[j];//*cosksi;
739  -original_parameters->Z_values[j+1]) + sqr(b-a);
740  area = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
741  aVector2.push_back(area);
742  }
743 
744  for(j=0; j<numPlanes-1; j++)
745  {
746  if(phiIsOpen == true)
747  {
748  aVector3.push_back(0.5*(original_parameters->Rmax[j]
751  -original_parameters->Rmin[j+1])
752  *std::fabs(original_parameters->Z_values[j+1]
754  }
755  else { aVector3.push_back(0.); }
756  }
757 
758  for(j=0; j<numPlanes-1; j++)
759  {
760  totArea += numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j];
761  }
762 
763  // Must include top and bottom areas
764  //
765  if(original_parameters->Rmax[numPlanes-1] != 0.)
766  {
767  a = original_parameters->Rmax[numPlanes-1];
768  b = original_parameters->Rmin[numPlanes-1];
769  l2 = sqr(a-b);
770  aTop = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
771  }
772 
773  if(original_parameters->Rmax[0] != 0.)
774  {
775  a = original_parameters->Rmax[0];
776  b = original_parameters->Rmin[0];
777  l2 = sqr(a-b);
778  aBottom = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
779  }
780 
781  Achose1 = 0.;
782  Achose2 = numSide*(aVector1[0]+aVector2[0])+2.*aVector3[0];
783 
784  chose = G4RandFlat::shoot(0.,totArea+aTop+aBottom);
785  if( (chose >= 0.) && (chose < aTop + aBottom) )
786  {
787  chose = G4RandFlat::shoot(startPhi,startPhi+totalPhi);
788  rang = std::floor((chose-startPhi)/ksi-0.01);
789  if(rang<0) { rang=0; }
790  rang = std::fabs(rang);
791  sinphi1 = std::sin(startPhi+rang*ksi);
792  sinphi2 = std::sin(startPhi+(rang+1)*ksi);
793  cosphi1 = std::cos(startPhi+rang*ksi);
794  cosphi2 = std::cos(startPhi+(rang+1)*ksi);
795  chose = G4RandFlat::shoot(0., aTop + aBottom);
796  if(chose>=0. && chose<aTop)
797  {
798  rad1 = original_parameters->Rmin[numPlanes-1];
799  rad2 = original_parameters->Rmax[numPlanes-1];
800  zVal = original_parameters->Z_values[numPlanes-1];
801  }
802  else
803  {
804  rad1 = original_parameters->Rmin[0];
805  rad2 = original_parameters->Rmax[0];
806  zVal = original_parameters->Z_values[0];
807  }
808  p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal);
809  p1 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal);
810  p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal);
811  p3 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal);
812  return GetPointOnPlane(p0,p1,p2,p3);
813  }
814  else
815  {
816  for (j=0; j<numPlanes-1; j++)
817  {
818  if( ((chose >= Achose1) && (chose < Achose2)) || (j == numPlanes-2) )
819  {
820  Flag = j; break;
821  }
822  Achose1 += numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j];
823  Achose2 = Achose1 + numSide*(aVector1[j+1]+aVector2[j+1])
824  + 2.*aVector3[j+1];
825  }
826  }
827 
828  // At this point we have chosen a subsection
829  // between to adjacent plane cuts...
830 
831  j = Flag;
832 
833  totArea = numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j];
834  chose = G4RandFlat::shoot(0.,totArea);
835 
836  if( (chose>=0.) && (chose<numSide*aVector1[j]) )
837  {
838  chose = G4RandFlat::shoot(startPhi,startPhi+totalPhi);
839  rang = std::floor((chose-startPhi)/ksi-0.01);
840  if(rang<0) { rang=0; }
841  rang = std::fabs(rang);
842  rad1 = original_parameters->Rmax[j];
843  rad2 = original_parameters->Rmax[j+1];
844  sinphi1 = std::sin(startPhi+rang*ksi);
845  sinphi2 = std::sin(startPhi+(rang+1)*ksi);
846  cosphi1 = std::cos(startPhi+rang*ksi);
847  cosphi2 = std::cos(startPhi+(rang+1)*ksi);
848  zVal = original_parameters->Z_values[j];
849 
850  p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal);
851  p1 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal);
852 
853  zVal = original_parameters->Z_values[j+1];
854 
855  p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal);
856  p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal);
857  return GetPointOnPlane(p0,p1,p2,p3);
858  }
859  else if ( (chose >= numSide*aVector1[j])
860  && (chose <= numSide*(aVector1[j]+aVector2[j])) )
861  {
862  chose = G4RandFlat::shoot(startPhi,startPhi+totalPhi);
863  rang = std::floor((chose-startPhi)/ksi-0.01);
864  if(rang<0) { rang=0; }
865  rang = std::fabs(rang);
866  rad1 = original_parameters->Rmin[j];
867  rad2 = original_parameters->Rmin[j+1];
868  sinphi1 = std::sin(startPhi+rang*ksi);
869  sinphi2 = std::sin(startPhi+(rang+1)*ksi);
870  cosphi1 = std::cos(startPhi+rang*ksi);
871  cosphi2 = std::cos(startPhi+(rang+1)*ksi);
872  zVal = original_parameters->Z_values[j];
873 
874  p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal);
875  p1 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal);
876 
877  zVal = original_parameters->Z_values[j+1];
878 
879  p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal);
880  p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal);
881  return GetPointOnPlane(p0,p1,p2,p3);
882  }
883 
884  chose = G4RandFlat::shoot(0.,2.2);
885  if( (chose>=0.) && (chose < 1.) )
886  {
887  rang = startPhi;
888  }
889  else
890  {
891  rang = endPhi;
892  }
893 
894  cosphi1 = std::cos(rang); rad1 = original_parameters->Rmin[j];
895  sinphi1 = std::sin(rang); rad2 = original_parameters->Rmax[j];
896 
897  p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,
899  p1 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,
901 
902  rad1 = original_parameters->Rmax[j+1];
903  rad2 = original_parameters->Rmin[j+1];
904 
905  p2 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,
907  p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,
909  return GetPointOnPlane(p0,p1,p2,p3);
910  }
911  else // Generic polyhedra
912  {
913  return GetPointOnSurfaceGeneric();
914  }
915 }
ThreeVector shoot(const G4int Ap, const G4int Af)
CLHEP::Hep3Vector G4ThreeVector
G4PolyhedraHistorical * original_parameters
Definition: G4Polyhedra.hh:187
G4bool genericPgon
Definition: G4Polyhedra.hh:184
int G4int
Definition: G4Types.hh:78
G4ThreeVector GetPointOnPlane(G4ThreeVector p0, G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3) const
Definition: G4Polyhedra.cc:654
G4double endPhi
Definition: G4Polyhedra.hh:182
G4double startPhi
Definition: G4Polyhedra.hh:181
static const double twopi
Definition: G4SIunits.hh:75
G4ThreeVector GetPointOnSurfaceGeneric() const
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
G4bool phiIsOpen
Definition: G4Polyhedra.hh:183
Here is the call graph for this function:

◆ GetPointOnSurfaceCorners()

G4ThreeVector G4Polyhedra::GetPointOnSurfaceCorners ( ) const
protected

◆ GetPointOnTriangle()

G4ThreeVector G4Polyhedra::GetPointOnTriangle ( G4ThreeVector  p0,
G4ThreeVector  p1,
G4ThreeVector  p2 
) const
protected

Definition at line 686 of file G4Polyhedra.cc.

689 {
690  G4double lambda1,lambda2;
691  G4ThreeVector v=p3-p1, w=p1-p2;
692 
693  lambda1 = G4RandFlat::shoot(0.,1.);
694  lambda2 = G4RandFlat::shoot(0.,lambda1);
695 
696  return (p2 + lambda1*w + lambda2*v);
697 }
ThreeVector shoot(const G4int Ap, const G4int Af)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ GetStartPhi()

G4double G4Polyhedra::GetStartPhi ( ) const
inline
Here is the caller graph for this function:

◆ Inside()

EInside G4Polyhedra::Inside ( const G4ThreeVector p) const
virtual

Reimplemented from G4VCSGfaceted.

Definition at line 523 of file G4Polyhedra.cc.

524 {
525  //
526  // Quick test
527  //
528  if (enclosingCylinder->MustBeOutside(p)) return kOutside;
529 
530  //
531  // Long answer
532  //
533  return G4VCSGfaceted::Inside(p);
534 }
virtual EInside Inside(const G4ThreeVector &p) const
G4EnclosingCylinder * enclosingCylinder
Definition: G4Polyhedra.hh:189
G4bool MustBeOutside(const G4ThreeVector &p) const
Here is the call graph for this function:

◆ IsGeneric()

G4bool G4Polyhedra::IsGeneric ( ) const
inline
Here is the caller graph for this function:

◆ IsOpen()

G4bool G4Polyhedra::IsOpen ( ) const
inline
Here is the caller graph for this function:

◆ operator=()

G4Polyhedra & G4Polyhedra::operator= ( const G4Polyhedra source)

Definition at line 412 of file G4Polyhedra.cc.

413 {
414  if (this == &source) return *this;
415 
416  G4VCSGfaceted::operator=( source );
417 
418  delete [] corners;
420 
421  delete enclosingCylinder;
422 
423  CopyStuff( source );
424 
425  return *this;
426 }
void CopyStuff(const G4Polyhedra &source)
Definition: G4Polyhedra.cc:432
G4VCSGfaceted & operator=(const G4VCSGfaceted &source)
G4PolyhedraHistorical * original_parameters
Definition: G4Polyhedra.hh:187
G4PolyhedraSideRZ * corners
Definition: G4Polyhedra.hh:186
G4EnclosingCylinder * enclosingCylinder
Definition: G4Polyhedra.hh:189
Here is the call graph for this function:

◆ Reset()

G4bool G4Polyhedra::Reset ( )

Definition at line 481 of file G4Polyhedra.cc.

482 {
483  if (genericPgon)
484  {
485  std::ostringstream message;
486  message << "Solid " << GetName() << " built using generic construct."
487  << G4endl << "Not applicable to the generic construct !";
488  G4Exception("G4Polyhedra::Reset()", "GeomSolids1001",
489  JustWarning, message, "Parameters NOT resetted.");
490  return 1;
491  }
492 
493  //
494  // Clear old setup
495  //
497  delete [] corners;
498  delete enclosingCylinder;
499 
500  //
501  // Rebuild polyhedra
502  //
503  G4ReduciblePolygon *rz =
511  delete rz;
512 
513  return 0;
514 }
G4PolyhedraHistorical * original_parameters
Definition: G4Polyhedra.hh:187
G4bool genericPgon
Definition: G4Polyhedra.hh:184
G4PolyhedraSideRZ * corners
Definition: G4Polyhedra.hh:186
void Create(G4double phiStart, G4double phiTotal, G4int numSide, G4ReduciblePolygon *rz)
Definition: G4Polyhedra.cc:201
G4String GetName() const
G4EnclosingCylinder * enclosingCylinder
Definition: G4Polyhedra.hh:189
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetOriginalParameters() [1/2]

void G4Polyhedra::SetOriginalParameters ( G4PolyhedraHistorical pars)
inline
Here is the caller graph for this function:

◆ SetOriginalParameters() [2/2]

void G4Polyhedra::SetOriginalParameters ( G4ReduciblePolygon rz)
protected

Definition at line 1176 of file G4Polyhedra.cc.

1177 {
1178  G4int numPlanes = (G4int)numCorner;
1179  G4bool isConvertible=true;
1180  G4double Zmax=rz->Bmax();
1181  rz->StartWithZMin();
1182 
1183  // Prepare vectors for storage
1184  //
1185  std::vector<G4double> Z;
1186  std::vector<G4double> Rmin;
1187  std::vector<G4double> Rmax;
1188 
1189  G4int countPlanes=1;
1190  G4int icurr=0;
1191  G4int icurl=0;
1192 
1193  // first plane Z=Z[0]
1194  //
1195  Z.push_back(corners[0].z);
1196  G4double Zprev=Z[0];
1197  if (Zprev == corners[1].z)
1198  {
1199  Rmin.push_back(corners[0].r);
1200  Rmax.push_back (corners[1].r);icurr=1;
1201  }
1202  else if (Zprev == corners[numPlanes-1].z)
1203  {
1204  Rmin.push_back(corners[numPlanes-1].r);
1205  Rmax.push_back (corners[0].r);
1206  icurl=numPlanes-1;
1207  }
1208  else
1209  {
1210  Rmin.push_back(corners[0].r);
1211  Rmax.push_back (corners[0].r);
1212  }
1213 
1214  // next planes until last
1215  //
1216  G4int inextr=0, inextl=0;
1217  for (G4int i=0; i < numPlanes-2; i++)
1218  {
1219  inextr=1+icurr;
1220  inextl=(icurl <= 0)? numPlanes-1 : icurl-1;
1221 
1222  if((corners[inextr].z >= Zmax) & (corners[inextl].z >= Zmax)) { break; }
1223 
1224  G4double Zleft = corners[inextl].z;
1225  G4double Zright = corners[inextr].z;
1226  if(Zright>Zleft)
1227  {
1228  Z.push_back(Zleft);
1229  countPlanes++;
1230  G4double difZr=corners[inextr].z - corners[icurr].z;
1231  G4double difZl=corners[inextl].z - corners[icurl].z;
1232 
1233  if(std::fabs(difZl) < kCarTolerance)
1234  {
1235  if(std::fabs(difZr) < kCarTolerance)
1236  {
1237  Rmin.push_back(corners[inextl].r);
1238  Rmax.push_back(corners[icurr].r);
1239  }
1240  else
1241  {
1242  Rmin.push_back(corners[inextl].r);
1243  Rmax.push_back(corners[icurr].r + (Zleft-corners[icurr].z)/difZr
1244  *(corners[inextr].r - corners[icurr].r));
1245  }
1246  }
1247  else if (difZl >= kCarTolerance)
1248  {
1249  if(std::fabs(difZr) < kCarTolerance)
1250  {
1251  Rmin.push_back(corners[icurl].r);
1252  Rmax.push_back(corners[icurr].r);
1253  }
1254  else
1255  {
1256  Rmin.push_back(corners[icurl].r);
1257  Rmax.push_back(corners[icurr].r + (Zleft-corners[icurr].z)/difZr
1258  *(corners[inextr].r - corners[icurr].r));
1259  }
1260  }
1261  else
1262  {
1263  isConvertible=false; break;
1264  }
1265  icurl=(icurl == 0)? numPlanes-1 : icurl-1;
1266  }
1267  else if(std::fabs(Zright-Zleft)<kCarTolerance) // Zright=Zleft
1268  {
1269  Z.push_back(Zleft);
1270  countPlanes++;
1271  icurr++;
1272 
1273  icurl=(icurl == 0)? numPlanes-1 : icurl-1;
1274 
1275  Rmin.push_back(corners[inextl].r);
1276  Rmax.push_back (corners[inextr].r);
1277  }
1278  else // Zright<Zleft
1279  {
1280  Z.push_back(Zright);
1281  countPlanes++;
1282 
1283  G4double difZr=corners[inextr].z - corners[icurr].z;
1284  G4double difZl=corners[inextl].z - corners[icurl].z;
1285  if(std::fabs(difZr) < kCarTolerance)
1286  {
1287  if(std::fabs(difZl) < kCarTolerance)
1288  {
1289  Rmax.push_back(corners[inextr].r);
1290  Rmin.push_back(corners[icurr].r);
1291  }
1292  else
1293  {
1294  Rmin.push_back(corners[icurl].r + (Zright-corners[icurl].z)/difZl
1295  * (corners[inextl].r - corners[icurl].r));
1296  Rmax.push_back(corners[inextr].r);
1297  }
1298  icurr++;
1299  } // plate
1300  else if (difZr >= kCarTolerance)
1301  {
1302  if(std::fabs(difZl) < kCarTolerance)
1303  {
1304  Rmax.push_back(corners[inextr].r);
1305  Rmin.push_back (corners[icurr].r);
1306  }
1307  else
1308  {
1309  Rmax.push_back(corners[inextr].r);
1310  Rmin.push_back (corners[icurl].r+(Zright-corners[icurl].z)/difZl
1311  * (corners[inextl].r - corners[icurl].r));
1312  }
1313  icurr++;
1314  }
1315  else
1316  {
1317  isConvertible=false; break;
1318  }
1319  }
1320  } // end for loop
1321 
1322  // last plane Z=Zmax
1323  //
1324  Z.push_back(Zmax);
1325  countPlanes++;
1326  inextr=1+icurr;
1327  inextl=(icurl <= 0)? numPlanes-1 : icurl-1;
1328 
1329  Rmax.push_back(corners[inextr].r);
1330  Rmin.push_back(corners[inextl].r);
1331 
1332  // Set original parameters Rmin,Rmax,Z
1333  //
1334  if(isConvertible)
1335  {
1338  original_parameters->Z_values = new G4double[countPlanes];
1339  original_parameters->Rmin = new G4double[countPlanes];
1340  original_parameters->Rmax = new G4double[countPlanes];
1341 
1342  for(G4int j=0; j < countPlanes; j++)
1343  {
1344  original_parameters->Z_values[j] = Z[j];
1345  original_parameters->Rmax[j] = Rmax[j];
1346  original_parameters->Rmin[j] = Rmin[j];
1347  }
1350  original_parameters->Num_z_planes = countPlanes;
1351 
1352  }
1353  else // Set parameters(r,z) with Rmin==0 as convention
1354  {
1355 #ifdef G4SPECSDEBUG
1356  std::ostringstream message;
1357  message << "Polyhedra " << GetName() << G4endl
1358  << "cannot be converted to Polyhedra with (Rmin,Rmaz,Z) parameters!";
1359  G4Exception("G4Polyhedra::SetOriginalParameters()",
1360  "GeomSolids0002", JustWarning, message);
1361 #endif
1364  original_parameters->Z_values = new G4double[numPlanes];
1365  original_parameters->Rmin = new G4double[numPlanes];
1366  original_parameters->Rmax = new G4double[numPlanes];
1367 
1368  for(G4int j=0; j < numPlanes; j++)
1369  {
1371  original_parameters->Rmax[j] = corners[j].r;
1372  original_parameters->Rmin[j] = 0.0;
1373  }
1376  original_parameters->Num_z_planes = numPlanes;
1377  }
1378  //return isConvertible;
1379 }
G4PolyhedraHistorical * original_parameters
Definition: G4Polyhedra.hh:187
G4PolyhedraSideRZ * corners
Definition: G4Polyhedra.hh:186
int G4int
Definition: G4Types.hh:78
G4String GetName() const
G4double Bmax() const
G4double endPhi
Definition: G4Polyhedra.hh:182
G4int numCorner
Definition: G4Polyhedra.hh:185
Float_t Z
bool G4bool
Definition: G4Types.hh:79
G4double startPhi
Definition: G4Polyhedra.hh:181
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
int Zmax
#define G4endl
Definition: G4ios.hh:61
G4double kCarTolerance
Definition: G4VSolid.hh:304
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ StreamInfo()

std::ostream & G4Polyhedra::StreamInfo ( std::ostream &  os) const
virtual

Reimplemented from G4VCSGfaceted.

Definition at line 600 of file G4Polyhedra.cc.

601 {
602  G4int oldprc = os.precision(16);
603  os << "-----------------------------------------------------------\n"
604  << " *** Dump for solid - " << GetName() << " ***\n"
605  << " ===================================================\n"
606  << " Solid type: G4Polyhedra\n"
607  << " Parameters: \n"
608  << " starting phi angle : " << startPhi/degree << " degrees \n"
609  << " ending phi angle : " << endPhi/degree << " degrees \n"
610  << " number of sides : " << numSide << " \n";
611  G4int i=0;
612  if (!genericPgon)
613  {
615  os << " number of Z planes: " << numPlanes << "\n"
616  << " Z values: \n";
617  for (i=0; i<numPlanes; i++)
618  {
619  os << " Z plane " << i << ": "
620  << original_parameters->Z_values[i] << "\n";
621  }
622  os << " Tangent distances to inner surface (Rmin): \n";
623  for (i=0; i<numPlanes; i++)
624  {
625  os << " Z plane " << i << ": "
626  << original_parameters->Rmin[i] << "\n";
627  }
628  os << " Tangent distances to outer surface (Rmax): \n";
629  for (i=0; i<numPlanes; i++)
630  {
631  os << " Z plane " << i << ": "
632  << original_parameters->Rmax[i] << "\n";
633  }
634  }
635  os << " number of RZ points: " << numCorner << "\n"
636  << " RZ values (corners): \n";
637  for (i=0; i<numCorner; i++)
638  {
639  os << " "
640  << corners[i].r << ", " << corners[i].z << "\n";
641  }
642  os << "-----------------------------------------------------------\n";
643  os.precision(oldprc);
644 
645  return os;
646 }
G4PolyhedraHistorical * original_parameters
Definition: G4Polyhedra.hh:187
G4bool genericPgon
Definition: G4Polyhedra.hh:184
G4PolyhedraSideRZ * corners
Definition: G4Polyhedra.hh:186
int G4int
Definition: G4Types.hh:78
G4String GetName() const
G4double endPhi
Definition: G4Polyhedra.hh:182
G4int numCorner
Definition: G4Polyhedra.hh:185
G4double startPhi
Definition: G4Polyhedra.hh:181
static const double degree
Definition: G4SIunits.hh:143
Here is the call graph for this function:

Member Data Documentation

◆ corners

G4PolyhedraSideRZ* G4Polyhedra::corners
protected

Definition at line 186 of file G4Polyhedra.hh.

◆ enclosingCylinder

G4EnclosingCylinder* G4Polyhedra::enclosingCylinder
protected

Definition at line 189 of file G4Polyhedra.hh.

◆ endPhi

G4double G4Polyhedra::endPhi
protected

Definition at line 182 of file G4Polyhedra.hh.

◆ genericPgon

G4bool G4Polyhedra::genericPgon
protected

Definition at line 184 of file G4Polyhedra.hh.

◆ numCorner

G4int G4Polyhedra::numCorner
protected

Definition at line 185 of file G4Polyhedra.hh.

◆ numSide

G4int G4Polyhedra::numSide
protected

Definition at line 180 of file G4Polyhedra.hh.

◆ original_parameters

G4PolyhedraHistorical* G4Polyhedra::original_parameters
protected

Definition at line 187 of file G4Polyhedra.hh.

◆ phiIsOpen

G4bool G4Polyhedra::phiIsOpen
protected

Definition at line 183 of file G4Polyhedra.hh.

◆ startPhi

G4double G4Polyhedra::startPhi
protected

Definition at line 181 of file G4Polyhedra.hh.


The documentation for this class was generated from the following files: