Geant4  10.00.p02
G4USolid.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 //
27 // $Id:$
28 // GEANT4 tag $Name:$
29 //
30 //
31 // G4USolid implementation
32 //
33 // --------------------------------------------------------------------
34 
35 #include "G4USolid.hh"
36 #include "G4AffineTransform.hh"
37 #include "G4VoxelLimits.hh"
38 #include "G4VGraphicsScene.hh"
39 #include "G4Polyhedron.hh"
40 #include "G4PolyhedronArbitrary.hh"
41 #include "G4VisExtent.hh"
42 #include "G4PhysicalConstants.hh"
43 
45  G4VSolid(name), fShape(s), fPolyhedron(0)
46 {
47 }
48 
50  : G4VSolid(a), fShape(0), fPolyhedron(0)
51 {
52 }
53 
55 {
56  delete fPolyhedron;
57 }
58 
60 {
61  return (this == &s) ? true : false;
62 }
63 
65 {
66  UVector3 pt;
67  VUSolid::EnumInside in_temp;
68  EInside in = kOutside;
69  pt.x = p.x();
70  pt.y = p.y();
71  pt.z = p.z(); // better assign at construction
72 
73  in_temp = fShape->Inside(pt);
74 
75  if (in_temp == VUSolid::eSurface)return kSurface;
76  if (in_temp == VUSolid::eInside)return kInside;
77 
78  return in;
79 }
80 
82 {
83  UVector3 p;
84  p.x = pt.x();
85  p.y = pt.y();
86  p.z = pt.z();
87  UVector3 n;
88  fShape->Normal(p, n);
89  return G4ThreeVector(n.x, n.y, n.z);
90 }
91 
93  const G4ThreeVector& d)const
94 {
95  UVector3 p;
96  p.x = pt.x();
97  p.y = pt.y();
98  p.z = pt.z(); // better assign at construction
99  UVector3 v;
100  v.x = d.x();
101  v.y = d.y();
102  v.z = d.z(); // better assign at construction
103  G4double dist = fShape->DistanceToIn(p, v);
104  if (dist > kInfinity) dist = kInfinity;
105  return dist;
106 }
107 
109 {
110  UVector3 p;
111  p.x = pt.x();
112  p.y = pt.y();
113  p.z = pt.z(); // better assign at construction
114  G4double dist = fShape->SafetyFromOutside(p); // true?
115  if (dist > kInfinity) dist = kInfinity;
116  return dist;
117 }
118 
120  const G4ThreeVector& d,
121  const G4bool calcNorm,
122  G4bool* validNorm,
123  G4ThreeVector* norm) const
124 {
125  UVector3 p;
126  p.x = pt.x();
127  p.y = pt.y();
128  p.z = pt.z(); // better assign at construction
129  UVector3 v;
130  v.x = d.x();
131  v.y = d.y();
132  v.z = d.z(); // better assign at construction
133  UVector3 n;
134  bool valid;
135  G4double dist = fShape->DistanceToOut(p, v, n,valid); // should use local variable
136  if(calcNorm)
137  {
138  if(valid){ *validNorm = true;}
139  else {* validNorm =false;}
140  if(*validNorm)
141  { norm->setX(n.x);
142  norm->setY(n.y);
143  norm->setZ(n.z);
144  } // *norm = n, but only after calcNorm check
145  }
146  if (dist > kInfinity) dist = kInfinity;
147  return dist;
148 }
149 
151 {
152  UVector3 p;
153  p.x = pt.x();
154  p.y = pt.y();
155  p.z = pt.z(); // better assign at construction
156  return fShape->SafetyFromInside(p); // true?
157 }
158 
160 {
161  return fShape->Capacity();
162 }
163 
165 {
166  return fShape->SurfaceArea();
167 }
168 
170 {
171  UVector3 p;
172  p = fShape->GetPointOnSurface();
173  return G4ThreeVector(p.x, p.y, p.z);
174 }
175 
177  const G4VoxelLimits& pVoxelLimit,
178  const G4AffineTransform& pTransform,
179  G4double& pMin, G4double& pMax) const
180 {
181  if (!pTransform.IsRotated())
182  {
184  G4double offset = pTransform.NetTranslation().x();
185  if (pAxis == kYAxis)
186  {
187  eAxis = VUSolid::eYaxis;
188  offset = pTransform.NetTranslation().y();
189  }
190  if (pAxis == kZAxis)
191  {
192  eAxis = VUSolid::eZaxis;
193  offset = pTransform.NetTranslation().z();
194  }
195  fShape->ExtentAxis(eAxis, pMin, pMax);
196 
197  pMin += offset;
198  pMax += offset;
199 
200  if (pVoxelLimit.IsLimited())
201  {
202  switch (pAxis)
203  {
204  case kXAxis:
205  if ((pMin > pVoxelLimit.GetMaxXExtent() + kCarTolerance) ||
206  (pMax < pVoxelLimit.GetMinXExtent() - kCarTolerance))
207  {
208  return false;
209  }
210  else
211  {
212  pMin = std::max(pMin, pVoxelLimit.GetMinXExtent());
213  pMax = std::min(pMax, pVoxelLimit.GetMaxXExtent());
214  }
215  break;
216  case kYAxis:
217  if ((pMin > pVoxelLimit.GetMaxYExtent() + kCarTolerance) ||
218  (pMax < pVoxelLimit.GetMinYExtent() - kCarTolerance))
219  {
220  return false;
221  }
222  else
223  {
224  pMin = std::max(pMin, pVoxelLimit.GetMinYExtent());
225  pMax = std::min(pMax, pVoxelLimit.GetMaxYExtent());
226  }
227  break;
228  case kZAxis:
229  if ((pMin > pVoxelLimit.GetMaxZExtent() + kCarTolerance) ||
230  (pMax < pVoxelLimit.GetMinZExtent() - kCarTolerance))
231  {
232  return false;
233  }
234  else
235  {
236  pMin = std::max(pMin, pVoxelLimit.GetMinZExtent());
237  pMax = std::min(pMax, pVoxelLimit.GetMaxZExtent());
238  }
239  break;
240  default:
241  break;
242  }
243  pMin -= kCarTolerance ;
244  pMax += kCarTolerance ;
245  }
246  return true;
247  }
248  else // General rotated case - create and clip mesh to boundaries
249  {
250  // Rotate BoundingBox and Calculate Extent as for BREPS
251 
252  G4bool existsAfterClip = false ;
253  G4ThreeVectorList* vertices ;
254 
255  pMin = +kInfinity ;
256  pMax = -kInfinity ;
257 
258  // Calculate rotated vertex coordinates
259 
260  vertices = CreateRotatedVertices(pTransform) ;
261  ClipCrossSection(vertices, 0, pVoxelLimit, pAxis, pMin, pMax) ;
262  ClipCrossSection(vertices, 4, pVoxelLimit, pAxis, pMin, pMax) ;
263  ClipBetweenSections(vertices, 0, pVoxelLimit, pAxis, pMin, pMax) ;
264 
265  if (pVoxelLimit.IsLimited(pAxis) == false)
266  {
267  if ((pMin != kInfinity) || (pMax != -kInfinity))
268  {
269  existsAfterClip = true ;
270 
271  // Add 2*tolerance to avoid precision troubles
272 
273  pMin -= kCarTolerance;
274  pMax += kCarTolerance;
275  }
276  }
277  else
278  {
279  G4ThreeVector clipCentre(
280  (pVoxelLimit.GetMinXExtent() + pVoxelLimit.GetMaxXExtent()) * 0.5,
281  (pVoxelLimit.GetMinYExtent() + pVoxelLimit.GetMaxYExtent()) * 0.5,
282  (pVoxelLimit.GetMinZExtent() + pVoxelLimit.GetMaxZExtent()) * 0.5);
283 
284  if ((pMin != kInfinity) || (pMax != -kInfinity))
285  {
286  existsAfterClip = true ;
287 
288 
289  // Check to see if endpoints are in the solid
290 
291  clipCentre(pAxis) = pVoxelLimit.GetMinExtent(pAxis);
292 
293  if (Inside(pTransform.Inverse().TransformPoint(clipCentre)) != kOutside)
294  {
295  pMin = pVoxelLimit.GetMinExtent(pAxis);
296  }
297  else
298  {
299  pMin -= kCarTolerance;
300  }
301  clipCentre(pAxis) = pVoxelLimit.GetMaxExtent(pAxis);
302 
303  if (Inside(pTransform.Inverse().TransformPoint(clipCentre)) != kOutside)
304  {
305  pMax = pVoxelLimit.GetMaxExtent(pAxis);
306  }
307  else
308  {
309  pMax += kCarTolerance;
310  }
311  }
312 
313  // Check for case where completely enveloping clipping volume
314  // If point inside then we are confident that the solid completely
315  // envelopes the clipping volume. Hence set min/max extents according
316  // to clipping volume extents along the specified axis.
317 
318  else if (Inside(pTransform.Inverse().TransformPoint(clipCentre))
319  != kOutside)
320  {
321  existsAfterClip = true ;
322  pMin = pVoxelLimit.GetMinExtent(pAxis) ;
323  pMax = pVoxelLimit.GetMaxExtent(pAxis) ;
324  }
325  }
326  delete vertices;
327  return existsAfterClip;
328  }
329 }
330 
332  const G4int,
333  const G4VPhysicalVolume*)
334 {
335  std::ostringstream message;
336  message << "Illegal call to G4USolid::ComputeDimensions()" << G4endl
337  << "Method not overloaded by derived class !";
338  G4Exception("G4USolid::ComputeDimensions()", "GeomSolids0003",
339  FatalException, message);
340 }
341 
343 {
344  scene.AddSolid(*this);
345 }
346 
348 {
349 
350  G4String string = fShape->GetEntityType();
351  return "G4" + string;
352 }
353 
354 std::ostream& G4USolid::StreamInfo(std::ostream& os) const
355 {
356  os << "-----------------------------------------------------------\n"
357  << " *** Dump for solid - " << fShape->GetName() << " ***\n";
358  // << " ===================================================\n"
359  // << " Solid type: " << fShape->GetEntityType() << "\n"
360  // << "-----------------------------------------------------------\n";
361 
362  return os;//fShape->StreamInfo(os);
363 }
364 
366  : G4VSolid(rhs), fPolyhedron(0)
367 {
368  fShape = rhs.fShape->Clone();
370 }
371 
373 {
374  // Check assignment to self
375  //
376  if (this == &rhs)
377  {
378  return *this;
379  }
380 
381  // Copy base class data
382  //
383  G4VSolid::operator=(rhs);
384 
385  // Copy data
386  //
387  fShape = rhs.fShape->Clone();
389 
390  return *this;
391 }
392 
394 {
395  std::ostringstream message;
396  message << "Clone() method not implemented for type: "
397  << GetEntityType() << "!" << G4endl
398  << "Returning NULL pointer!";
399  G4Exception("G4USolid::Clone()", "GeomSolids1001", JustWarning, message);
400  return 0;
401 }
402 
405 {
406  G4double xMin, xMax, yMin, yMax, zMin, zMax;
407 
408  fShape->ExtentAxis(VUSolid::eXaxis, xMin, xMax);
409  fShape->ExtentAxis(VUSolid::eYaxis, yMin, yMax);
410  fShape->ExtentAxis(VUSolid::eZaxis, zMin, zMax);
411 
412  G4ThreeVectorList* vertices;
413  vertices = new G4ThreeVectorList();
414 
415  if (vertices)
416  {
417  vertices->reserve(8);
418  G4ThreeVector vertex0(xMin, yMin, zMin);
419  G4ThreeVector vertex1(xMax, yMin, zMin);
420  G4ThreeVector vertex2(xMax, yMax, zMin);
421  G4ThreeVector vertex3(xMin, yMax, zMin);
422  G4ThreeVector vertex4(xMin, yMin, zMax);
423  G4ThreeVector vertex5(xMax, yMin, zMax);
424  G4ThreeVector vertex6(xMax, yMax, zMax);
425  G4ThreeVector vertex7(xMin, yMax, zMax);
426 
427  vertices->push_back(pTransform.TransformPoint(vertex0));
428  vertices->push_back(pTransform.TransformPoint(vertex1));
429  vertices->push_back(pTransform.TransformPoint(vertex2));
430  vertices->push_back(pTransform.TransformPoint(vertex3));
431  vertices->push_back(pTransform.TransformPoint(vertex4));
432  vertices->push_back(pTransform.TransformPoint(vertex5));
433  vertices->push_back(pTransform.TransformPoint(vertex6));
434  vertices->push_back(pTransform.TransformPoint(vertex7));
435  }
436  else
437  {
438  G4Exception("G4VUSolid::CreateRotatedVertices()", "FatalError",
439  FatalException, "Out of memory - Cannot allocate vertices!");
440  }
441  return vertices;
442 }
443 
445 {
446  G4int index = 0;
447  if (fShape->GetEntityType() == "Box")
448  {
449  double array[3];
450  fShape->GetParametersList(index, array);
451  return new G4PolyhedronBox(array[0], array[1], array[2]);
452  }
453  if (fShape->GetEntityType() == "Tubs")
454  {
455  double array[5];
456  fShape->GetParametersList(index, array);
457  return new G4PolyhedronTubs(array[0], array[1], array[2], array[3], array[4]);
458  }
459  if (fShape->GetEntityType() == "Cons")
460  {
461  double array[7];
462  fShape->GetParametersList(index, array);
463  return new G4PolyhedronCons(array[0], array[1], array[2], array[3], array[4], array[5], array[6]);
464  }
465  if (fShape->GetEntityType() == "Orb")
466  {
467  double array[1];
468  fShape->GetParametersList(index, array);
469  return new G4PolyhedronSphere(0., array[0], 0., 2 * pi, 0., pi);
470  }
471  if (fShape->GetEntityType() == "Sphere")
472  {
473  double array[6];
474  fShape->GetParametersList(index, array);
475  return new G4PolyhedronSphere(array[0], array[1], array[2], array[3], array[4], array[5]);
476  }
477  if (fShape->GetEntityType() == "Tet")
478  {
479  double array[12];
480  fShape->GetParametersList(index, array);
481  G4Polyhedron* ph = new G4Polyhedron;
482  double xyz[4][3];
483  static int faces[4][4] = {{1, 3, 2, 0}, {1, 4, 3, 0}, {1, 2, 4, 0}, {2, 3, 4, 0}};
484  xyz[0][0] = array[0];
485  xyz[0][1] = array[1];
486  xyz[0][2] = array[2];
487  xyz[1][0] = array[3];
488  xyz[1][1] = array[4];
489  xyz[1][2] = array[5];
490  xyz[2][0] = array[6];
491  xyz[2][1] = array[7];
492  xyz[2][2] = array[8];
493  xyz[3][0] = array[9];
494  xyz[3][1] = array[10];
495  xyz[3][2] = array[11];
496 
497  ph->createPolyhedron(4, 4, xyz, faces);
498  return ph;
499  }
500  if (fShape->GetEntityType() == "Trd")
501  {
502  double array[5];
503  fShape->GetParametersList(index, array);
504  return new G4PolyhedronTrd2(array[0], array[1], array[2], array[3], array[4]);
505  }
506  if (fShape->GetEntityType() == "Trap")
507  {
508  double array[12];
509  fShape->GetParametersList(index, array);
510  double phi = (array[11] != 1.0) ? (std::atan(array[10] / array[9])) : (0.0);
511  double alpha1 = std::atan(array[4]);
512  double alpha2 = std::atan(array[8]);
513  double theta = std::acos(array[11]);
514 
515  return new G4PolyhedronTrap(array[0], theta, phi,
516  array[1], array[2], array[3], alpha1,
517  array[5], array[6], array[7], alpha2);
518  }
519 
520  /*
521  if(fShape->GetEntityType()=="TessellatedSolid"){
522 
523  G4Polyhedron *uPolyhedron=fShape->GetPolyhedron();
524  std::size_t nVertices = (*uPolyhedron).vertices.size();
525  std::size_t nFacets = (*uPolyhedron).facets.size();
526 
527  G4PolyhedronArbitrary *polyhedron =
528  new G4PolyhedronArbitrary (nVertices, nFacets);
529 
530  for (std::vector<UVector3>::const_iterator v = (*uPolyhedron).vertices.begin();
531  v!=(*uPolyhedron).vertices.end(); v++)
532  {
533  UVector3 p=(*v);
534  G4ThreeVector pt(p.x,p.y,p.z);
535 
536  polyhedron->AddVertex(pt);
537  }
538  for (std::vector<UFacet>::const_iterator f=(*uPolyhedron).facets.begin();
539  f != (*uPolyhedron).facets.end(); f++)
540  {
541  polyhedron->AddFacet((*f).f1,(*f).f2,(*f).f3,(*f).f4);
542  }
543 
544  return (G4Polyhedron*) polyhedron;
545  }
546  */
547 
548  return 0;
549 }
550 
552 {
553  if (!fPolyhedron ||
555  fPolyhedron->GetNumberOfRotationSteps())
556  {
557  delete fPolyhedron;
559  }
560  return fPolyhedron;
561 }
562 
564 {
565  G4VisExtent extent;
566  G4VoxelLimits voxelLimits; // Defaults to "infinite" limits.
567  G4AffineTransform affineTransform;
568  G4double vmin, vmax;
569  CalculateExtent(kXAxis, voxelLimits, affineTransform, vmin, vmax);
570  extent.SetXmin(vmin);
571  extent.SetXmax(vmax);
572  CalculateExtent(kYAxis, voxelLimits, affineTransform, vmin, vmax);
573  extent.SetYmin(vmin);
574  extent.SetYmax(vmax);
575  CalculateExtent(kZAxis, voxelLimits, affineTransform, vmin, vmax);
576  extent.SetZmin(vmin);
577  extent.SetZmax(vmax);
578  return extent;
579 }
virtual UVector3 GetPointOnSurface() const =0
void ClipCrossSection(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:347
G4ThreeVectorList * CreateRotatedVertices(const G4AffineTransform &pT) const
Definition: G4USolid.cc:404
virtual UGeometryType GetEntityType() const =0
virtual double SurfaceArea()=0
Definition: VUSolid.cc:255
virtual void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4USolid.cc:342
static const G4double kInfinity
Definition: geomdefs.hh:42
G4double GetMinYExtent() const
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4USolid.cc:119
CLHEP::Hep3Vector G4ThreeVector
G4AffineTransform Inverse() const
G4Polyhedron * CreatePolyhedron() const
Definition: G4USolid.cc:444
const std::string & GetName() const
Definition: VUSolid.hh:103
G4String name
Definition: TRTMaterials.hh:40
const G4double pi
G4bool IsRotated() const
virtual std::ostream & StreamInfo(std::ostream &os) const
Definition: G4USolid.cc:354
void SetXmax(G4double xmax)
Definition: G4VisExtent.hh:102
G4ThreeVector NetTranslation() const
virtual G4double GetCubicVolume()
Definition: G4USolid.cc:159
virtual G4double GetSurfaceArea()
Definition: G4USolid.cc:164
G4double a
Definition: TRTMaterials.hh:39
virtual bool Normal(const UVector3 &aPoint, UVector3 &aNormal) const =0
virtual EnumInside Inside(const UVector3 &aPoint) const =0
virtual EInside Inside(const G4ThreeVector &p) const
Definition: G4USolid.cc:64
virtual void AddSolid(const G4Box &)=0
int G4int
Definition: G4Types.hh:78
double x
Definition: UVector3.hh:136
G4bool operator==(const G4USolid &s) const
Definition: G4USolid.cc:59
void SetYmax(G4double ymax)
Definition: G4VisExtent.hh:106
static const double s
Definition: G4SIunits.hh:150
G4double GetMaxXExtent() const
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4USolid.cc:551
G4double GetMinZExtent() const
virtual double DistanceToIn(const UVector3 &aPoint, const UVector3 &aDirection, double aPstep=UUtils::kInfinity) const =0
G4bool IsLimited() const
VUSolid * fShape
Definition: G4USolid.hh:182
bool G4bool
Definition: G4Types.hh:79
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:79
EAxisType
Definition: VUSolid.hh:27
EnumInside
Definition: VUSolid.hh:23
const G4int n
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4USolid.cc:92
virtual double SafetyFromInside(const UVector3 &aPoint, bool aAccurate=false) const =0
G4double GetMinXExtent() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4double GetMaxZExtent() const
void SetZmin(G4double zmin)
Definition: G4VisExtent.hh:108
G4USolid(const G4String &pName, VUSolid *shape)
Definition: G4USolid.cc:44
G4Polyhedron * fPolyhedron
Definition: G4USolid.hh:183
T max(const T t1, const T t2)
brief Return the largest of the two arguments
EInside
Definition: geomdefs.hh:58
virtual void ExtentAxis(EAxisType aAxis, double &aMin, double &aMax) const
Definition: VUSolid.cc:318
EAxis
Definition: geomdefs.hh:54
virtual double SafetyFromOutside(const UVector3 &aPoint, bool aAccurate=false) const =0
virtual G4VisExtent GetExtent() const
Definition: G4USolid.cc:563
G4USolid & operator=(const G4USolid &rhs)
Definition: G4USolid.cc:372
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4USolid.cc:81
#define G4endl
Definition: G4ios.hh:61
G4double GetMaxYExtent() const
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:110
G4double kCarTolerance
Definition: G4VSolid.hh:305
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
virtual double DistanceToOut(const UVector3 &aPoint, const UVector3 &aDirection, UVector3 &aNormalVector, bool &aConvex, double aPstep=UUtils::kInfinity) const =0
void SetZmax(G4double zmax)
Definition: G4VisExtent.hh:110
double z
Definition: UVector3.hh:138
virtual VUSolid * Clone() const =0
virtual ~G4USolid()
Definition: G4USolid.cc:54
void SetXmin(G4double xmin)
Definition: G4VisExtent.hh:100
double G4double
Definition: G4Types.hh:76
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
Definition: G4USolid.cc:176
G4double GetMaxExtent(const EAxis pAxis) const
void SetYmin(G4double ymin)
Definition: G4VisExtent.hh:104
virtual G4VSolid * Clone() const
Definition: G4USolid.cc:393
double y
Definition: UVector3.hh:137
void ClipBetweenSections(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:378
G4double GetMinExtent(const EAxis pAxis) const
virtual double Capacity()=0
Definition: VUSolid.cc:199
virtual G4GeometryType GetEntityType() const
Definition: G4USolid.cc:347
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4USolid.cc:331
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4USolid.cc:169
virtual void GetParametersList(int aNumber, double *aArray) const =0