Geant4  10.00.p03
G4VoxelSafety.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 // $Id: G4VoxelSafety.cc,v $
27 //
28 // Author: John Apostolakis
29 // First version: 31 May 2010
30 //
31 // --------------------------------------------------------------------
32 #include "G4VoxelSafety.hh"
33 
34 #include "G4GeometryTolerance.hh"
35 
36 #include "G4SmartVoxelProxy.hh"
37 #include "G4SmartVoxelNode.hh"
38 #include "G4SmartVoxelHeader.hh"
39 
40 // ********************************************************************
41 // Constructor
42 // - copied from G4VoxelNavigation (1st version)
43 // ********************************************************************
44 //
46  : fBlockList(),
47  fpMotherLogical(0),
48  fVoxelDepth(-1),
49  fVoxelAxisStack(kNavigatorVoxelStackMax,kXAxis),
50  fVoxelNoSlicesStack(kNavigatorVoxelStackMax,0),
51  fVoxelSliceWidthStack(kNavigatorVoxelStackMax,0.),
52  fVoxelNodeNoStack(kNavigatorVoxelStackMax,0),
53  fVoxelHeaderStack(kNavigatorVoxelStackMax,(G4SmartVoxelHeader*)0),
54  fVoxelNode(0),
55  fCheck(false),
56  fVerbose(0)
57 {
59 }
60 
61 // ********************************************************************
62 // Destructor
63 // ********************************************************************
64 //
66 {
67 }
68 
69 // ********************************************************************
70 // ComputeSafety
71 //
72 // Calculates the isotropic distance to the nearest boundary from the
73 // specified point in the local coordinate system.
74 // The localpoint utilised must be within the current volume.
75 // ********************************************************************
76 //
79  const G4VPhysicalVolume& currentPhysical,
80  G4double maxLength)
81 {
82  G4LogicalVolume *motherLogical;
83  G4VSolid *motherSolid;
84  G4SmartVoxelHeader *motherVoxelHeader;
85  G4double motherSafety, ourSafety;
86  G4int localNoDaughters;
87  G4double daughterSafety;
88 
89  motherLogical = currentPhysical.GetLogicalVolume();
90  fpMotherLogical= motherLogical; // For use by the other methods
91  motherSolid = motherLogical->GetSolid();
92  motherVoxelHeader= motherLogical->GetVoxelHeader();
93 
94 #ifdef G4VERBOSE
95  if( fVerbose > 0 )
96  {
97  G4cout << "*** G4VoxelSafety::ComputeSafety(): ***" << G4endl;
98  }
99 #endif
100 
101  // Check that point is inside mother volume
102  //
103  EInside insideMother= motherSolid->Inside(localPoint);
104  if( insideMother != kInside )
105  {
106 #ifdef G4DEBUG_NAVIGATION
107  if( insideMother == kOutside )
108  {
109  std::ostringstream message;
110  message << "Safety method called for location outside current Volume."
111  << G4endl
112  << "Location for safety is Outside this volume. " << G4endl
113  << "The approximate distance to the solid "
114  << "(safety from outside) is: "
115  << motherSolid->DistanceToIn( localPoint ) << G4endl;
116  message << " Problem occurred with physical volume: "
117  << " Name: " << currentPhysical.GetName()
118  << " Copy No: " << currentPhysical.GetCopyNo() << G4endl
119  << " Local Point = " << localPoint << G4endl;
120  message << " Description of solid: " << G4endl
121  << *motherSolid << G4endl;
122  G4Exception("G4VoxelSafety::ComputeSafety()", "GeomNav0003",
123  FatalException, message);
124  }
125 #endif
126  return 0.0;
127  }
128 
129  // First limit: mother safety - distance to outer boundaries
130  //
131  motherSafety = motherSolid->DistanceToOut(localPoint);
132  ourSafety = motherSafety; // Working isotropic safety
133 
134 #ifdef G4VERBOSE
135  if(( fCheck ) ) // && ( fVerbose == 1 ))
136  {
137  G4cout << " Invoked DistanceToOut(p) for mother solid: "
138  << motherSolid->GetName()
139  << ". Solid replied: " << motherSafety << G4endl
140  << " For local point p: " << localPoint
141  << ", to be considered as 'mother safety'." << G4endl;
142  }
143 #endif
144  localNoDaughters = motherLogical->GetNoDaughters();
145 
146  fBlockList.Enlarge(localNoDaughters);
147  fBlockList.Reset();
148 
149  fVoxelDepth = -1; // Resets the depth -- must be done for next method
150  daughterSafety= SafetyForVoxelHeader(motherVoxelHeader, localPoint, maxLength,
151  currentPhysical, 0, ourSafety);
152  ourSafety= std::min( motherSafety, daughterSafety );
153 
154  return ourSafety;
155 }
156 
157 // ********************************************************************
158 // SafetyForVoxelNode
159 //
160 // Calculate the safety for volumes included in current Voxel Node
161 // ********************************************************************
162 //
163 G4double
165  const G4ThreeVector& localPoint )
166 {
167  G4double ourSafety= DBL_MAX;
168 
169  G4int curNoVolumes, contentNo, sampleNo;
170  G4VPhysicalVolume *samplePhysical;
171 
172  G4double sampleSafety=0.0;
173  G4ThreeVector samplePoint;
174  G4VSolid* ptrSolid=0;
175 
176  curNoVolumes = curVoxelNode->GetNoContained();
177 
178  for ( contentNo=curNoVolumes-1; contentNo>=0; contentNo-- )
179  {
180  sampleNo = curVoxelNode->GetVolume(contentNo);
181  if ( !fBlockList.IsBlocked(sampleNo) )
182  {
183  fBlockList.BlockVolume(sampleNo);
184 
185  samplePhysical = fpMotherLogical->GetDaughter(sampleNo);
186  G4AffineTransform sampleTf(samplePhysical->GetRotation(),
187  samplePhysical->GetTranslation());
188  sampleTf.Invert();
189  samplePoint = sampleTf.TransformPoint(localPoint);
190  ptrSolid = samplePhysical->GetLogicalVolume()->GetSolid();
191 
192  sampleSafety = ptrSolid->DistanceToIn(samplePoint);
193  ourSafety = std::min( sampleSafety, ourSafety );
194 #ifdef G4VERBOSE
195  if(( fCheck ) && ( fVerbose == 1 ))
196  {
197  // ReportSolidSafetyToIn( MethodName, solid, value, point );
198  G4cout << "*** G4VoxelSafety::SafetyForVoxelNode(): ***" << G4endl
199  << " Invoked DistanceToIn(p) for daughter solid: "
200  << ptrSolid->GetName()
201  << ". Solid replied: " << sampleSafety << G4endl
202  << " For local point p: " << samplePoint
203  << ", to be considered as 'daughter safety'." << G4endl;
204  }
205 #endif
206  }
207  } // end for contentNo
208 
209  return ourSafety;
210 }
211 
212 // ********************************************************************
213 // SafetyForVoxelHeader
214 //
215 // Cycles through levels of headers to process each node level
216 // Obtained by modifying VoxelLocate (to cycle through Node Headers)
217 // *********************************************************************
218 //
219 G4double
221  const G4ThreeVector& localPoint,
222  G4double maxLength,
223  const G4VPhysicalVolume& currentPhysical, //Debug
224  G4double distUpperDepth_Sq,
225  G4double previousMinSafety
226  )
227 {
228  const G4SmartVoxelHeader * const targetVoxelHeader=pHeader;
229  G4SmartVoxelNode *targetVoxelNode=0;
230 
231  const G4SmartVoxelProxy *sampleProxy;
232  EAxis targetHeaderAxis;
233  G4double targetHeaderMin, targetHeaderMax, targetHeaderNodeWidth;
234  G4int targetHeaderNoSlices;
235  G4int targetNodeNo;
236 
237  G4double minSafety= previousMinSafety;
238  G4double ourSafety= DBL_MAX;
239  unsigned int checkedNum= 0;
240 
241  fVoxelDepth++;
242  // fVoxelDepth set by ComputeSafety or previous level call
243 
244  targetHeaderAxis = targetVoxelHeader->GetAxis();
245  targetHeaderNoSlices = targetVoxelHeader->GetNoSlices();
246  targetHeaderMin = targetVoxelHeader->GetMinExtent();
247  targetHeaderMax = targetVoxelHeader->GetMaxExtent();
248 
249  targetHeaderNodeWidth = (targetHeaderMax-targetHeaderMin)
250  / targetHeaderNoSlices;
251 
252  G4double localCrd= localPoint(targetHeaderAxis);
253 
254  const G4int candNodeNo= G4int( (localCrd-targetHeaderMin)
255  / targetHeaderNodeWidth );
256  // Ensure that it is between 0 and targetHeader->GetMaxExtent() - 1
257 
258 #ifdef G4DEBUG_VOXELISATION
259  if( candNodeNo < 0 || candNodeNo > targetHeaderNoSlices-1 )
260  {
262  ed << " Potential ERROR."
263  << " Point is outside range of Voxel in current coordinate" << G4endl;
264  ed << " Node number of point " << localPoint
265  << "is outside the range. " << G4endl;
266  ed << " Voxel node Num= " << candNodeNo << " versus minimum= " << 0
267  << " and maximum= " << targetHeaderNoSlices-1 << G4endl;
268  ed << " Axis = " << targetHeaderAxis
269  << " No of slices = " << targetHeaderNoSlices << G4endl;
270  ed << " Local coord = " << localCrd
271  << " Voxel Min = " << targetHeaderMin
272  << " Max = " << targetHeaderMax << G4endl;
273  G4LogicalVolume *pLogical= currentPhysical.GetLogicalVolume();
274  ed << " Current volume (physical) = " << currentPhysical.GetName()
275  << " (logical) = " << pLogical->GetName() << G4endl;
276  G4VSolid* pSolid= pLogical->GetSolid();
277  ed << " Solid type = " << pSolid->GetEntityType() << G4endl;
278  ed << *pSolid << G4endl;
279  G4Exception("G4VoxelSafety::SafetyForVoxelHeader()", "GeomNav1003",
280  JustWarning, ed,
281  "Point is outside range of Voxel in current coordinate");
282  }
283 #endif
284 
285  const G4int pointNodeNo =
286  std::max( 0, std::min( candNodeNo, targetHeaderNoSlices-1 ) );
287 
288 #ifdef G4VERBOSE
289  if( fVerbose > 2 )
290  {
291  G4cout << G4endl;
292  G4cout << "**** G4VoxelSafety::SafetyForVoxelHeader " << G4endl;
293  G4cout << " Called at Depth = " << fVoxelDepth;
294  G4cout << " Calculated pointNodeNo= " << pointNodeNo
295  << " from position= " << localPoint(targetHeaderAxis)
296  << " min= " << targetHeaderMin
297  << " max= " << targetHeaderMax
298  << " width= " << targetHeaderNodeWidth
299  << " no-slices= " << targetHeaderNoSlices
300  << " axis= " << targetHeaderAxis << G4endl;
301  }
302  else if (fVerbose == 1)
303  {
304  G4cout << " VoxelSafety: Depth = " << fVoxelDepth
305  << " Number of Slices = " << targetHeaderNoSlices
306  << " Header (address) = " << targetVoxelHeader << G4endl;
307  }
308 #endif
309 
310  // Stack info for stepping
311  //
312  fVoxelAxisStack[fVoxelDepth] = targetHeaderAxis;
313  fVoxelNoSlicesStack[fVoxelDepth] = targetHeaderNoSlices;
314  fVoxelSliceWidthStack[fVoxelDepth] = targetHeaderNodeWidth;
315 
316  fVoxelHeaderStack[fVoxelDepth] = pHeader;
317 
318  G4int trialUp= -1, trialDown= -1;
319  G4double distUp= DBL_MAX, distDown= DBL_MAX;
320 
321  // Using Equivalent voxels - this is pre-initialisation only
322  //
323  G4int nextUp= pointNodeNo+1;
324  G4int nextDown= pointNodeNo-1;
325 
326  G4int nextNodeNo= pointNodeNo;
327  G4double distAxis; // Distance in current Axis
328  distAxis= 0.0; // Starting in node containing local Coordinate
329 
330  G4bool nextIsInside= false;
331 
332  G4double distMaxInterest= std::min( previousMinSafety, maxLength);
333  // We will not look beyond this distance.
334  // This distance will be updated to reflect the
335  // max ( minSafety, maxLength ) at each step
336 
337  targetNodeNo= pointNodeNo;
338  do
339  {
340  G4double nodeSafety= DBL_MAX, headerSafety= DBL_MAX;
341  fVoxelNodeNoStack[fVoxelDepth] = targetNodeNo;
342 
343  checkedNum++;
344 
345  sampleProxy = targetVoxelHeader->GetSlice(targetNodeNo);
346 
347 #ifdef G4DEBUG_NAVIGATION
348  assert( sampleProxy != 0);
349  if( fVerbose > 2 )
350  {
351  G4cout << " -Checking node " << targetNodeNo
352  << " is proxy with address " << sampleProxy << G4endl;
353  }
354 #endif
355 
356  if ( sampleProxy == 0 )
357  {
359  ed << " Problem for node number= " << targetNodeNo
360  << " Number of slides = " << targetHeaderNoSlices
361  << G4endl;
362  G4Exception( "G4VoxelSafety::SafetyForVoxelHeader()", "GeomNav0003",
363  FatalException, ed,
364  "Problem sampleProxy is Zero. Failure in loop.");
365  }
366  else if ( sampleProxy->IsNode() )
367  {
368  targetVoxelNode = sampleProxy->GetNode();
369 
370  // Deal with the node here [ i.e. the last level ]
371  //
372  nodeSafety= SafetyForVoxelNode( targetVoxelNode, localPoint);
373 #ifdef G4DEBUG_NAVIGATION
374  if( fVerbose > 2 )
375  {
376  G4cout << " -- It is a Node ";
377  G4cout << " its safety= " << nodeSafety
378  << " our level Saf = " << ourSafety
379  << " MinSaf= " << minSafety << G4endl;
380  }
381 #endif
382  ourSafety= std::min( ourSafety, nodeSafety );
383 
384  trialUp = targetVoxelNode->GetMaxEquivalentSliceNo()+1;
385  trialDown = targetVoxelNode->GetMinEquivalentSliceNo()-1;
386  }
387  else
388  {
389  const G4SmartVoxelHeader *pNewVoxelHeader = sampleProxy->GetHeader();
390 
391  G4double distCombined_Sq;
392  distCombined_Sq = distUpperDepth_Sq + distAxis*distAxis;
393 
394 #ifdef G4DEBUG_NAVIGATION
395  if( fVerbose > 2 )
396  {
397  G4double distCombined= std::sqrt( distCombined_Sq );
398  G4double distUpperDepth= std::sqrt ( distUpperDepth_Sq );
399  G4cout << " -- It is a Header " << G4endl;
400  G4cout << " Recurse to deal with next level, fVoxelDepth= "
401  << fVoxelDepth+1 << G4endl;
402  G4cout << " Distances: Upper= " << distUpperDepth
403  << " Axis= " << distAxis
404  << " Combined= " << distCombined << G4endl;
405  }
406 #endif
407 
408  // Recurse to deal with lower levels
409  //
410  headerSafety= SafetyForVoxelHeader( pNewVoxelHeader, localPoint,
411  maxLength, currentPhysical,
412  distCombined_Sq, minSafety);
413  ourSafety= std::min( ourSafety, headerSafety );
414 
415 #ifdef G4DEBUG_NAVIGATION
416  if( fVerbose > 2 )
417  {
418  G4cout << " >> Header safety = " << headerSafety
419  << " our level Saf = " << ourSafety << G4endl;
420  }
421 #endif
422  trialUp = pNewVoxelHeader->GetMaxEquivalentSliceNo()+1;
423  trialDown = pNewVoxelHeader->GetMinEquivalentSliceNo()-1;
424  }
425  minSafety= std::min( minSafety, ourSafety );
426 
427  // Find next closest Voxel
428  // - first try: by simple subtraction
429  // - later: using distance (TODO - tbc)
430  //
431  if( targetNodeNo >= pointNodeNo )
432  {
433  nextUp = trialUp;
434  // distUp = std::max( targetHeaderMax-localCrd, 0.0 );
435  G4double lowerEdgeOfNext = targetHeaderMin
436  + nextUp * targetHeaderNodeWidth;
437  distUp = lowerEdgeOfNext-localCrd ;
438  if( distUp < 0.0 )
439  {
440  distUp = DBL_MAX; // On the wrong side - must not be considered
441  }
442 #ifdef G4DEBUG_NAVIGATION
443  if( fVerbose > 2 )
444  {
445  G4cout << " > Updated nextUp= " << nextUp << G4endl;
446  }
447 #endif
448  }
449 
450  if( targetNodeNo <= pointNodeNo )
451  {
452  nextDown = trialDown;
453  // distDown = std::max( localCrd-targetHeaderMin, 0.0);
454  G4double upperEdgeOfNext = targetHeaderMin
455  + (nextDown+1) * targetHeaderNodeWidth;
456  distDown = localCrd-upperEdgeOfNext;
457  if( distDown < 0.0 )
458  {
459  distDown= DBL_MAX; // On the wrong side - must not be considered
460  }
461 #ifdef G4DEBUG_NAVIGATION
462  if( fVerbose > 2 )
463  {
464  G4cout << " > Updated nextDown= " << nextDown << G4endl;
465  }
466 #endif
467  }
468 
469 #ifdef G4DEBUG_NAVIGATION
470  if( fVerbose > 2 )
471  {
472  G4cout << " Node= " << pointNodeNo
473  << " Up: next= " << nextUp << " d# "
474  << nextUp - pointNodeNo
475  << " trialUp= " << trialUp << " d# "
476  << trialUp - pointNodeNo
477  << " Down: next= " << nextDown << " d# "
478  << targetNodeNo - nextDown
479  << " trialDwn= " << trialDown << " d# "
480  << targetNodeNo - trialDown
481  << " condition (next is Inside)= " << nextIsInside
482  << G4endl;
483  }
484 #endif
485 
486  G4bool UpIsClosest;
487  UpIsClosest= distUp < distDown;
488 
489  if( (nextUp < targetHeaderNoSlices)
490  && (UpIsClosest || (nextDown < 0)) )
491  {
492  nextNodeNo=nextUp;
493  distAxis = distUp;
494  ++nextUp; // Default
495 #ifdef G4VERBOSE
496  if( fVerbose > 2 )
497  {
498  G4cout << " > Chose Up. Depth= " << fVoxelDepth
499  << " Nodes: next= " << nextNodeNo
500  << " new nextUp= " << nextUp
501  << " Dist = " << distAxis << G4endl;
502  }
503 #endif
504  }
505  else
506  {
507  nextNodeNo=nextDown;
508  distAxis = distDown;
509  --nextDown; // A default value
510 #ifdef G4VERBOSE
511  if( fVerbose > 2 )
512  {
513  G4cout << " > Chose Down. Depth= " << fVoxelDepth
514  << " Nodes: next= " << nextNodeNo
515  << " new nextDown= " << nextDown
516  << " Dist = " << distAxis << G4endl;
517  }
518 #endif
519  }
520 
521  nextIsInside = (nextNodeNo >= 0) && (nextNodeNo < targetHeaderNoSlices);
522  if( nextIsInside )
523  {
524  targetNodeNo= nextNodeNo;
525 
526 #ifdef G4DEBUG_NAVIGATION
527  assert( targetVoxelHeader->GetSlice(nextNodeNo) != 0 );
528  G4bool bContinue= (distAxis<minSafety);
529  if( !bContinue )
530  {
531  if( fVerbose > 2 )
532  {
533  G4cout << " Can skip remaining at depth " << targetHeaderAxis
534  << " >> distAxis= " << distAxis
535  << " minSafety= " << minSafety << G4endl;
536  }
537  }
538 #endif
539  }
540  else
541  {
542 #ifdef G4DEBUG_NAVIGATION
543  if( fVerbose > 2)
544  {
545  G4cout << " VoxSaf> depth= " << fVoxelDepth << G4endl;
546  G4cout << " VoxSaf> No more candidates: nodeDown= " << nextDown
547  << " nodeUp= " << nextUp
548  << " lastSlice= " << targetHeaderNoSlices << G4endl;
549  }
550 #endif
551  }
552 
553  // This calculation can be 'hauled'-up to where minSafety is calculated
554  //
555  distMaxInterest = std::min( minSafety, maxLength );
556 
557  } while ( nextIsInside && ( distAxis*distAxis + distUpperDepth_Sq
558  < distMaxInterest*distMaxInterest ) );
559 
560 #ifdef G4VERBOSE
561  if( fVerbose > 0 )
562  {
563  G4cout << " Ended for targetNodeNo -- checked " << checkedNum << " out of "
564  << targetHeaderNoSlices << " slices." << G4endl;
565  G4cout << " ===== Returning from SafetyForVoxelHeader "
566  << " Depth= " << fVoxelDepth << G4endl
567  << G4endl;
568  }
569 #endif
570 
571  // Go back one level
572  fVoxelDepth--;
573 
574  return ourSafety;
575 }
G4String GetName() const
G4LogicalVolume * fpMotherLogical
G4SmartVoxelHeader * GetVoxelHeader() const
G4SmartVoxelHeader * GetHeader() const
G4int GetMinEquivalentSliceNo() const
G4int GetMaxEquivalentSliceNo() const
G4double SafetyForVoxelNode(const G4SmartVoxelNode *curVoxelNode, const G4ThreeVector &localPoint)
G4double SafetyForVoxelHeader(const G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint, G4double maxLength, const G4VPhysicalVolume &currentPhysical, G4double distUpperDepth=0.0, G4double previousMinSafety=DBL_MAX)
const G4ThreeVector & GetTranslation() const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4String GetName() const
CLHEP::Hep3Vector G4ThreeVector
G4int GetVolume(G4int pVolumeNo) const
#define assert(x)
Definition: mymalloc.cc:1309
std::vector< G4int > fVoxelNodeNoStack
G4double GetSurfaceTolerance() const
G4VPhysicalVolume * GetDaughter(const G4int i) const
G4double ComputeSafety(const G4ThreeVector &localPoint, const G4VPhysicalVolume &currentPhysical, G4double maxLength=DBL_MAX)
virtual G4GeometryType GetEntityType() const =0
int G4int
Definition: G4Types.hh:78
G4int GetNoContained() const
G4int GetMaxEquivalentSliceNo() const
G4bool IsNode() const
G4AffineTransform & Invert()
EAxis GetAxis() const
G4double GetMinExtent() const
G4GLOB_DLL std::ostream G4cout
std::vector< G4double > fVoxelSliceWidthStack
const G4String & GetName() const
virtual EInside Inside(const G4ThreeVector &p) const =0
bool G4bool
Definition: G4Types.hh:79
void BlockVolume(const G4int v)
G4bool IsBlocked(const G4int v) const
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
G4int GetNoSlices() const
G4double GetMaxExtent() const
std::vector< const G4SmartVoxelHeader * > fVoxelHeaderStack
G4int GetNoDaughters() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::vector< G4int > fVoxelNoSlicesStack
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4LogicalVolume * GetLogicalVolume() const
G4SmartVoxelNode * GetNode() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
EInside
Definition: geomdefs.hh:58
G4BlockingList fBlockList
G4double kCarTolerance
EAxis
Definition: geomdefs.hh:54
virtual G4int GetCopyNo() const =0
const G4RotationMatrix * GetRotation() const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static const G4int kNavigatorVoxelStackMax
Definition: geomdefs.hh:80
#define G4endl
Definition: G4ios.hh:61
G4SmartVoxelProxy * GetSlice(G4int n) const
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
G4int GetMinEquivalentSliceNo() const
#define DBL_MAX
Definition: templates.hh:83
static G4GeometryTolerance * GetInstance()
G4VSolid * GetSolid() const
void Enlarge(const G4int nv)
std::vector< EAxis > fVoxelAxisStack