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