Geant4  10.00.p03
G4VSolid.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: G4VSolid.cc 72936 2013-08-14 13:17:11Z gcosmo $
28 //
29 // class G4VSolid
30 //
31 // Implementation for solid base class
32 //
33 // History:
34 //
35 // 06.12.02 V.Grichine, restored original conditions in ClipPolygon()
36 // 10.05.02 V.Grichine, ClipPolygon(): clip only other axis and limited voxels
37 // 15.04.02 V.Grichine, bug fixed in ClipPolygon(): clip only one axis
38 // 13.03.02 V.Grichine, cosmetics of voxel limit functions
39 // 15.11.00 D.Williams, V.Grichine, fix in CalculateClippedPolygonExtent()
40 // 10.07.95 P.Kent, Added == operator, solid Store entry
41 // 30.06.95 P.Kent, Created.
42 // --------------------------------------------------------------------
43 
44 #include "G4VSolid.hh"
45 #include "G4SolidStore.hh"
46 #include "globals.hh"
47 #include "Randomize.hh"
48 #include "G4GeometryTolerance.hh"
49 
50 #include "G4VoxelLimits.hh"
51 #include "G4AffineTransform.hh"
52 #include "G4VisExtent.hh"
53 
55 //
56 // Constructor
57 // - Copies name
58 // - Add ourselves to solid Store
59 
61  : fshapeName(name)
62 {
64 
65  // Register to store
66  //
68 }
69 
71 //
72 // Copy constructor
73 //
74 
76  : kCarTolerance(rhs.kCarTolerance), fshapeName(rhs.fshapeName)
77 {
78  // Register to store
79  //
81 }
82 
84 //
85 // Fake default constructor - sets only member data and allocates memory
86 // for usage restricted to object persistency.
87 //
88 G4VSolid::G4VSolid( __void__& )
89  : fshapeName("")
90 {
91  // Register to store
92  //
94 }
95 
97 //
98 // Destructor (virtual)
99 // - Remove ourselves from solid Store
100 
102 {
104 }
105 
107 //
108 // Assignment operator
109 
111 {
112  // Check assignment to self
113  //
114  if (this == &rhs) { return *this; }
115 
116  // Copy data
117  //
119  fshapeName = rhs.fshapeName;
120 
121  return *this;
122 }
123 
125 //
126 // Streaming operator dumping solid contents
127 
128 std::ostream& operator<< ( std::ostream& os, const G4VSolid& e )
129 {
130  return e.StreamInfo(os);
131 }
132 
134 //
135 // Throw exception if ComputeDimensions called for illegal derived class
136 
138  const G4int,
139  const G4VPhysicalVolume*)
140 {
141  std::ostringstream message;
142  message << "Illegal call to G4VSolid::ComputeDimensions()" << G4endl
143  << "Method not overloaded by derived class !";
144  G4Exception("G4VSolid::ComputeDimensions()", "GeomMgt0003",
145  FatalException, message);
146 }
147 
149 //
150 // Throw exception (warning) for solids not implementing the method
151 
153 {
154  std::ostringstream message;
155  message << "Not implemented for solid: "
156  << this->GetEntityType() << " !" << G4endl
157  << "Returning origin.";
158  G4Exception("G4VSolid::GetPointOnSurface()", "GeomMgt1001",
159  JustWarning, message);
160  return G4ThreeVector(0,0,0);
161 }
162 
164 //
165 // Dummy implementations ...
166 
168 { return 0; }
169 
171 { return 0; }
172 
174 { return 0; }
175 
177 { return 0; }
178 
180 //
181 // Returns an estimation of the solid volume in internal units.
182 // The number of statistics and error accuracy is fixed.
183 // This method may be overloaded by derived classes to compute the
184 // exact geometrical quantity for solids where this is possible.
185 // or anyway to cache the computed value.
186 // This implementation does NOT cache the computed value.
187 
189 {
190  G4int cubVolStatistics = 1000000;
191  G4double cubVolEpsilon = 0.001;
192  return EstimateCubicVolume(cubVolStatistics, cubVolEpsilon);
193 }
194 
196 //
197 // Calculate cubic volume based on Inside() method.
198 // Accuracy is limited by the second argument or the statistics
199 // expressed by the first argument.
200 // Implementation is courtesy of Vasiliki Despoina Mitsou,
201 // University of Athens.
202 
204 {
205  G4int iInside=0;
206  G4double px,py,pz,minX,maxX,minY,maxY,minZ,maxZ,volume,halfepsilon;
207  G4ThreeVector p;
208  EInside in;
209 
210  // values needed for CalculateExtent signature
211 
212  G4VoxelLimits limit; // Unlimited
213  G4AffineTransform origin;
214 
215  // min max extents of pSolid along X,Y,Z
216 
217  this->CalculateExtent(kXAxis,limit,origin,minX,maxX);
218  this->CalculateExtent(kYAxis,limit,origin,minY,maxY);
219  this->CalculateExtent(kZAxis,limit,origin,minZ,maxZ);
220 
221  // limits
222 
223  if(nStat < 100) nStat = 100;
224  if(epsilon > 0.01) epsilon = 0.01;
225  halfepsilon = 0.5*epsilon;
226 
227  for(G4int i = 0; i < nStat; i++ )
228  {
229  px = minX-halfepsilon+(maxX-minX+epsilon)*G4UniformRand();
230  py = minY-halfepsilon+(maxY-minY+epsilon)*G4UniformRand();
231  pz = minZ-halfepsilon+(maxZ-minZ+epsilon)*G4UniformRand();
232  p = G4ThreeVector(px,py,pz);
233  in = this->Inside(p);
234  if(in != kOutside) iInside++;
235  }
236  volume = (maxX-minX+epsilon)*(maxY-minY+epsilon)
237  * (maxZ-minZ+epsilon)*iInside/nStat;
238  return volume;
239 }
240 
242 //
243 // Returns an estimation of the solid surface area in internal units.
244 // The number of statistics and error accuracy is fixed.
245 // This method may be overloaded by derived classes to compute the
246 // exact geometrical quantity for solids where this is possible.
247 // or anyway to cache the computed value.
248 // This implementation does NOT cache the computed value.
249 
251 {
252  G4int stat = 1000000;
253  G4double ell = -1.;
254  return EstimateSurfaceArea(stat,ell);
255 }
256 
258 //
259 // Estimate surface area based on Inside(), DistanceToIn(), and
260 // DistanceToOut() methods. Accuracy is limited by the statistics
261 // defined by the first argument. Implemented by Mikhail Kosov.
262 
264 {
265  G4int inside=0;
266  G4double px,py,pz,minX,maxX,minY,maxY,minZ,maxZ,surf;
267  G4ThreeVector p;
268  EInside in;
269 
270  // values needed for CalculateExtent signature
271 
272  G4VoxelLimits limit; // Unlimited
273  G4AffineTransform origin;
274 
275  // min max extents of pSolid along X,Y,Z
276 
277  this->CalculateExtent(kXAxis,limit,origin,minX,maxX);
278  this->CalculateExtent(kYAxis,limit,origin,minY,maxY);
279  this->CalculateExtent(kZAxis,limit,origin,minZ,maxZ);
280 
281  // limits
282 
283  if(nStat < 100) { nStat = 100; }
284 
285  G4double dX=maxX-minX;
286  G4double dY=maxY-minY;
287  G4double dZ=maxZ-minZ;
288  if(ell<=0.) // Automatic definition of skin thickness
289  {
290  G4double minval=dX;
291  if(dY<dX) { minval=dY; }
292  if(dZ<minval) { minval=dZ; }
293  ell=.01*minval;
294  }
295 
296  G4double dd=2*ell;
297  minX-=ell; minY-=ell; minZ-=ell; dX+=dd; dY+=dd; dZ+=dd;
298 
299  for(G4int i = 0; i < nStat; i++ )
300  {
301  px = minX+dX*G4UniformRand();
302  py = minY+dY*G4UniformRand();
303  pz = minZ+dZ*G4UniformRand();
304  p = G4ThreeVector(px,py,pz);
305  in = this->Inside(p);
306  if(in != kOutside)
307  {
308  if (DistanceToOut(p)<ell) { inside++; }
309  }
310  else if(DistanceToIn(p)<ell) { inside++; }
311  }
312  // @@ The conformal correction can be upgraded
313  surf = dX*dY*dZ*inside/dd/nStat;
314  return surf;
315 }
316 
318 //
319 // Returns a pointer of a dynamically allocated copy of the solid.
320 // Returns NULL pointer with warning in case the concrete solid does not
321 // implement this method. The caller has responsibility for ownership.
322 //
323 
325 {
326  std::ostringstream message;
327  message << "Clone() method not implemented for type: "
328  << GetEntityType() << "!" << G4endl
329  << "Returning NULL pointer!";
330  G4Exception("G4VSolid::Clone()", "GeomMgt1001", JustWarning, message);
331  return 0;
332 }
333 
335 //
336 // Calculate the maximum and minimum extents of the polygon described
337 // by the vertices: pSectionIndex->pSectionIndex+1->
338 // pSectionIndex+2->pSectionIndex+3->pSectionIndex
339 // in the List pVertices
340 //
341 // If the minimum is <pMin pMin is set to the new minimum
342 // If the maximum is >pMax pMax is set to the new maximum
343 //
344 // No modifications are made to pVertices
345 //
346 
348  const G4int pSectionIndex,
349  const G4VoxelLimits& pVoxelLimit,
350  const EAxis pAxis,
351  G4double& pMin, G4double& pMax) const
352 {
353 
354  G4ThreeVectorList polygon;
355  polygon.reserve(4);
356  polygon.push_back((*pVertices)[pSectionIndex]);
357  polygon.push_back((*pVertices)[pSectionIndex+1]);
358  polygon.push_back((*pVertices)[pSectionIndex+2]);
359  polygon.push_back((*pVertices)[pSectionIndex+3]);
360  // G4cout<<"ClipCrossSection: 0-1-2-3"<<G4endl;
361  CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
362  return;
363 }
364 
366 //
367 // Calculate the maximum and minimum extents of the polygons
368 // joining the CrossSections at pSectionIndex->pSectionIndex+3 and
369 // pSectionIndex+4->pSectionIndex7
370 //
371 // in the List pVertices, within the boundaries of the voxel limits pVoxelLimit
372 //
373 // If the minimum is <pMin pMin is set to the new minimum
374 // If the maximum is >pMax pMax is set to the new maximum
375 //
376 // No modifications are made to pVertices
377 
379  const G4int pSectionIndex,
380  const G4VoxelLimits& pVoxelLimit,
381  const EAxis pAxis,
382  G4double& pMin, G4double& pMax) const
383 {
384  G4ThreeVectorList polygon;
385  polygon.reserve(4);
386  polygon.push_back((*pVertices)[pSectionIndex]);
387  polygon.push_back((*pVertices)[pSectionIndex+4]);
388  polygon.push_back((*pVertices)[pSectionIndex+5]);
389  polygon.push_back((*pVertices)[pSectionIndex+1]);
390  // G4cout<<"ClipBetweenSections: 0-4-5-1"<<G4endl;
391  CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
392  polygon.clear();
393 
394  polygon.push_back((*pVertices)[pSectionIndex+1]);
395  polygon.push_back((*pVertices)[pSectionIndex+5]);
396  polygon.push_back((*pVertices)[pSectionIndex+6]);
397  polygon.push_back((*pVertices)[pSectionIndex+2]);
398  // G4cout<<"ClipBetweenSections: 1-5-6-2"<<G4endl;
399  CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
400  polygon.clear();
401 
402  polygon.push_back((*pVertices)[pSectionIndex+2]);
403  polygon.push_back((*pVertices)[pSectionIndex+6]);
404  polygon.push_back((*pVertices)[pSectionIndex+7]);
405  polygon.push_back((*pVertices)[pSectionIndex+3]);
406  // G4cout<<"ClipBetweenSections: 2-6-7-3"<<G4endl;
407  CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
408  polygon.clear();
409 
410  polygon.push_back((*pVertices)[pSectionIndex+3]);
411  polygon.push_back((*pVertices)[pSectionIndex+7]);
412  polygon.push_back((*pVertices)[pSectionIndex+4]);
413  polygon.push_back((*pVertices)[pSectionIndex]);
414  // G4cout<<"ClipBetweenSections: 3-7-4-0"<<G4endl;
415  CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
416  return;
417 }
418 
419 
421 //
422 // Calculate the maximum and minimum extents of the convex polygon pPolygon
423 // along the axis pAxis, within the limits pVoxelLimit
424 //
425 
426 void
428  const G4VoxelLimits& pVoxelLimit,
429  const EAxis pAxis,
430  G4double& pMin,
431  G4double& pMax) const
432 {
433  G4int noLeft,i;
434  G4double component;
435  /*
436  G4cout<<G4endl;
437  for(i = 0 ; i < pPolygon.size() ; i++ )
438  {
439  G4cout << i << "\t"
440  << "p.x = " << pPolygon[i].operator()(pAxis) << "\t"
441  // << "p.y = " << pPolygon[i].y() << "\t"
442  // << "p.z = " << pPolygon[i].z() << "\t"
443  << G4endl;
444  }
445  G4cout<<G4endl;
446  */
447  ClipPolygon(pPolygon,pVoxelLimit,pAxis);
448  noLeft = pPolygon.size();
449 
450  if ( noLeft )
451  {
452  // G4cout<<G4endl;
453  for (i=0;i<noLeft;i++)
454  {
455  component = pPolygon[i].operator()(pAxis);
456  // G4cout <<i<<"\t"<<component<<G4endl;
457 
458  if (component < pMin)
459  {
460  // G4cout <<i<<"\t"<<"Pmin = "<<component<<G4endl;
461  pMin = component;
462  }
463  if (component > pMax)
464  {
465  // G4cout <<i<<"\t"<<"PMax = "<<component<<G4endl;
466  pMax = component;
467  }
468  }
469  // G4cout<<G4endl;
470  }
471  // G4cout<<"pMin = "<<pMin<<"\t"<<"pMax = "<<pMax<<G4endl;
472 }
473 
475 //
476 // Clip the convex polygon described by the vertices at
477 // pSectionIndex ->pSectionIndex+3 within pVertices to the limits pVoxelLimit
478 //
479 // Set pMin to the smallest
480 //
481 // Calculate the extent of the polygon along pAxis, when clipped to the
482 // limits pVoxelLimit. If the polygon exists after clippin, set pMin to
483 // the polygon's minimum extent along the axis if <pMin, and set pMax to
484 // the polygon's maximum extent along the axis if >pMax.
485 //
486 // The polygon is described by a set of vectors, where each vector represents
487 // a vertex, so that the polygon is described by the vertex sequence:
488 // 0th->1st 1st->2nd 2nd->... nth->0th
489 //
490 // Modifications to the polygon are made
491 //
492 // NOTE: Execessive copying during clipping
493 
495  const G4VoxelLimits& pVoxelLimit,
496  const EAxis ) const
497 {
498  G4ThreeVectorList outputPolygon;
499 
500  if ( pVoxelLimit.IsLimited() )
501  {
502  if (pVoxelLimit.IsXLimited() ) // && pAxis != kXAxis)
503  {
504  G4VoxelLimits simpleLimit1;
505  simpleLimit1.AddLimit(kXAxis,pVoxelLimit.GetMinXExtent(),kInfinity);
506  // G4cout<<"MinXExtent()"<<G4endl;
507  ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
508 
509  pPolygon.clear();
510 
511  if ( !outputPolygon.size() ) return;
512 
513  G4VoxelLimits simpleLimit2;
514  // G4cout<<"MaxXExtent()"<<G4endl;
515  simpleLimit2.AddLimit(kXAxis,-kInfinity,pVoxelLimit.GetMaxXExtent());
516  ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
517 
518  if ( !pPolygon.size() ) return;
519  else outputPolygon.clear();
520  }
521  if ( pVoxelLimit.IsYLimited() ) // && pAxis != kYAxis)
522  {
523  G4VoxelLimits simpleLimit1;
524  simpleLimit1.AddLimit(kYAxis,pVoxelLimit.GetMinYExtent(),kInfinity);
525  ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
526 
527  // Must always clear pPolygon - for clip to simpleLimit2 and in case of
528  // early exit
529 
530  pPolygon.clear();
531 
532  if ( !outputPolygon.size() ) return;
533 
534  G4VoxelLimits simpleLimit2;
535  simpleLimit2.AddLimit(kYAxis,-kInfinity,pVoxelLimit.GetMaxYExtent());
536  ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
537 
538  if ( !pPolygon.size() ) return;
539  else outputPolygon.clear();
540  }
541  if ( pVoxelLimit.IsZLimited() ) // && pAxis != kZAxis)
542  {
543  G4VoxelLimits simpleLimit1;
544  simpleLimit1.AddLimit(kZAxis,pVoxelLimit.GetMinZExtent(),kInfinity);
545  ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
546 
547  // Must always clear pPolygon - for clip to simpleLimit2 and in case of
548  // early exit
549 
550  pPolygon.clear();
551 
552  if ( !outputPolygon.size() ) return;
553 
554  G4VoxelLimits simpleLimit2;
555  simpleLimit2.AddLimit(kZAxis,-kInfinity,pVoxelLimit.GetMaxZExtent());
556  ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
557 
558  // Return after final clip - no cleanup
559  }
560  }
561 }
562 
564 //
565 // pVoxelLimits must be only limited along one axis, and either the maximum
566 // along the axis must be +kInfinity, or the minimum -kInfinity
567 
568 void
570  G4ThreeVectorList& outputPolygon,
571  const G4VoxelLimits& pVoxelLimit ) const
572 {
573  G4int i;
574  G4int noVertices=pPolygon.size();
575  G4ThreeVector vEnd,vStart;
576 
577  for (i = 0 ; i < noVertices ; i++ )
578  {
579  vStart = pPolygon[i];
580  // G4cout << "i = " << i << G4endl;
581  if ( i == noVertices-1 ) vEnd = pPolygon[0];
582  else vEnd = pPolygon[i+1];
583 
584  if ( pVoxelLimit.Inside(vStart) )
585  {
586  if (pVoxelLimit.Inside(vEnd))
587  {
588  // vStart and vEnd inside -> output end point
589  //
590  outputPolygon.push_back(vEnd);
591  }
592  else
593  {
594  // vStart inside, vEnd outside -> output crossing point
595  //
596  // G4cout << "vStart inside, vEnd outside" << G4endl;
597  pVoxelLimit.ClipToLimits(vStart,vEnd);
598  outputPolygon.push_back(vEnd);
599  }
600  }
601  else
602  {
603  if (pVoxelLimit.Inside(vEnd))
604  {
605  // vStart outside, vEnd inside -> output inside section
606  //
607  // G4cout << "vStart outside, vEnd inside" << G4endl;
608  pVoxelLimit.ClipToLimits(vStart,vEnd);
609  outputPolygon.push_back(vStart);
610  outputPolygon.push_back(vEnd);
611  }
612  else // Both point outside -> no output
613  {
614  // outputPolygon.push_back(vStart);
615  // outputPolygon.push_back(vEnd);
616  }
617  }
618  }
619 }
620 
622 {
623  G4VisExtent extent;
624  G4VoxelLimits voxelLimits; // Defaults to "infinite" limits.
625  G4AffineTransform affineTransform;
626  G4double vmin, vmax;
627  CalculateExtent(kXAxis,voxelLimits,affineTransform,vmin,vmax);
628  extent.SetXmin (vmin);
629  extent.SetXmax (vmax);
630  CalculateExtent(kYAxis,voxelLimits,affineTransform,vmin,vmax);
631  extent.SetYmin (vmin);
632  extent.SetYmax (vmax);
633  CalculateExtent(kZAxis,voxelLimits,affineTransform,vmin,vmax);
634  extent.SetZmin (vmin);
635  extent.SetZmax (vmax);
636  return extent;
637 }
638 
640 {
641  return 0;
642 }
643 
645 {
646  return 0;
647 }
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4VSolid.cc:644
void ClipCrossSection(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:347
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
static const G4double kInfinity
Definition: geomdefs.hh:42
G4double GetMinYExtent() const
static void Register(G4VSolid *pSolid)
CLHEP::Hep3Vector G4ThreeVector
static void DeRegister(G4VSolid *pSolid)
void ClipPolygonToSimpleLimits(G4ThreeVectorList &pPolygon, G4ThreeVectorList &outputPolygon, const G4VoxelLimits &pVoxelLimit) const
Definition: G4VSolid.cc:569
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:137
G4bool IsYLimited() const
G4String name
Definition: TRTMaterials.hh:40
G4double EstimateCubicVolume(G4int nStat, G4double epsilon) const
Definition: G4VSolid.cc:203
G4double GetSurfaceTolerance() const
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188
void SetXmax(G4double xmax)
Definition: G4VisExtent.hh:102
G4bool IsXLimited() const
virtual G4GeometryType GetEntityType() const =0
int G4int
Definition: G4Types.hh:78
void SetYmax(G4double ymax)
Definition: G4VisExtent.hh:106
G4double GetMaxXExtent() const
void AddLimit(const EAxis pAxis, const G4double pMin, const G4double pMax)
G4double GetMinZExtent() const
virtual std::ostream & StreamInfo(std::ostream &os) const =0
#define G4UniformRand()
Definition: Randomize.hh:87
G4bool IsLimited() const
void CalculateClippedPolygonExtent(G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:427
virtual EInside Inside(const G4ThreeVector &p) const =0
G4bool Inside(const G4ThreeVector &pVec) const
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:79
static G4SolidStore * GetInstance()
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
virtual G4Polyhedron * CreatePolyhedron() const
Definition: G4VSolid.cc:639
G4double EstimateSurfaceArea(G4int nStat, G4double ell) const
Definition: G4VSolid.cc:263
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual ~G4VSolid()
Definition: G4VSolid.cc:101
G4double GetMinXExtent() const
virtual G4VisExtent GetExtent() const
Definition: G4VSolid.cc:621
G4double GetMaxZExtent() const
void SetZmin(G4double zmin)
Definition: G4VisExtent.hh:108
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:152
virtual const G4VSolid * GetConstituentSolid(G4int no) const
Definition: G4VSolid.cc:167
#define G4endl
Definition: G4ios.hh:61
G4VSolid(const G4String &name)
Definition: G4VSolid.cc:60
G4double GetMaxYExtent() const
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:110
G4double kCarTolerance
Definition: G4VSolid.hh:305
virtual G4VSolid * Clone() const
Definition: G4VSolid.cc:324
void SetZmax(G4double zmax)
Definition: G4VisExtent.hh:110
void SetXmin(G4double xmin)
Definition: G4VisExtent.hh:100
double G4double
Definition: G4Types.hh:76
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
G4String fshapeName
Definition: G4VSolid.hh:318
std::ostream & operator<<(std::ostream &os, const G4VSolid &e)
Definition: G4VSolid.cc:128
void SetYmin(G4double ymin)
Definition: G4VisExtent.hh:104
G4bool IsZLimited() const
virtual G4double GetSurfaceArea()
Definition: G4VSolid.cc:250
void ClipBetweenSections(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:378
virtual const G4DisplacedSolid * GetDisplacedSolidPtr() const
Definition: G4VSolid.cc:173
static G4GeometryTolerance * GetInstance()
void ClipPolygon(G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis) const
Definition: G4VSolid.cc:494
G4bool ClipToLimits(G4ThreeVector &pStart, G4ThreeVector &pEnd) const