Geant4  10.00.p03
G4VoxelNavigation.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: G4VoxelNavigation.cc 66356 2012-12-18 09:02:32Z gcosmo $
28 //
29 //
30 // class G4VoxelNavigation Implementation
31 //
32 // Author: P.Kent, 1996
33 //
34 // --------------------------------------------------------------------
35 
36 #include "G4VoxelNavigation.hh"
37 #include "G4GeometryTolerance.hh"
38 #include "G4VoxelSafety.hh"
39 
40 // ********************************************************************
41 // Constructor
42 // ********************************************************************
43 //
45  : fBList(), fVoxelDepth(-1),
46  fVoxelAxisStack(kNavigatorVoxelStackMax,kXAxis),
47  fVoxelNoSlicesStack(kNavigatorVoxelStackMax,0),
48  fVoxelSliceWidthStack(kNavigatorVoxelStackMax,0.),
49  fVoxelNodeNoStack(kNavigatorVoxelStackMax,0),
50  fVoxelHeaderStack(kNavigatorVoxelStackMax,(G4SmartVoxelHeader*)0),
51  fVoxelNode(0), fpVoxelSafety(0), fCheck(false), fBestSafety(false)
52 {
53  fLogger = new G4NavigationLogger("G4VoxelNavigation");
54  fpVoxelSafety = new G4VoxelSafety ();
55 }
56 
57 // ********************************************************************
58 // Destructor
59 // ********************************************************************
60 //
62 {
63  delete fpVoxelSafety;
64  delete fLogger;
65 }
66 
67 // ********************************************************************
68 // ComputeStep
69 // ********************************************************************
70 //
73  const G4ThreeVector& localDirection,
74  const G4double currentProposedStepLength,
75  G4double& newSafety,
76  G4NavigationHistory& history,
77  G4bool& validExitNormal,
78  G4ThreeVector& exitNormal,
79  G4bool& exiting,
80  G4bool& entering,
81  G4VPhysicalVolume *(*pBlockedPhysical),
82  G4int& blockedReplicaNo )
83 {
84  G4VPhysicalVolume *motherPhysical, *samplePhysical, *blockedExitedVol=0;
85  G4LogicalVolume *motherLogical;
86  G4VSolid *motherSolid;
87  G4ThreeVector sampleDirection;
88  G4double ourStep=currentProposedStepLength, motherSafety, ourSafety;
89  G4int localNoDaughters, sampleNo;
90 
91  G4bool initialNode, noStep;
92  G4SmartVoxelNode *curVoxelNode;
93  G4int curNoVolumes, contentNo;
94  G4double voxelSafety;
95 
96  motherPhysical = history.GetTopVolume();
97  motherLogical = motherPhysical->GetLogicalVolume();
98  motherSolid = motherLogical->GetSolid();
99 
100  //
101  // Compute mother safety
102  //
103 
104  motherSafety = motherSolid->DistanceToOut(localPoint);
105  ourSafety = motherSafety; // Working isotropic safety
106 
107 #ifdef G4VERBOSE
108  if ( fCheck )
109  {
110  fLogger->PreComputeStepLog (motherPhysical, motherSafety, localPoint);
111  }
112 #endif
113 
114  //
115  // Compute daughter safeties & intersections
116  //
117 
118  // Exiting normal optimisation
119  //
120  if ( exiting && validExitNormal )
121  {
122  if ( localDirection.dot(exitNormal)>=kMinExitingNormalCosine )
123  {
124  // Block exited daughter volume
125  //
126  blockedExitedVol = *pBlockedPhysical;
127  ourSafety = 0;
128  }
129  }
130  exiting = false;
131  entering = false;
132 
133  localNoDaughters = motherLogical->GetNoDaughters();
134 
135  fBList.Enlarge(localNoDaughters);
136  fBList.Reset();
137 
138  initialNode = true;
139  noStep = true;
140 
141  while (noStep)
142  {
143  curVoxelNode = fVoxelNode;
144  curNoVolumes = curVoxelNode->GetNoContained();
145  for (contentNo=curNoVolumes-1; contentNo>=0; contentNo--)
146  {
147  sampleNo = curVoxelNode->GetVolume(contentNo);
148  if ( !fBList.IsBlocked(sampleNo) )
149  {
150  fBList.BlockVolume(sampleNo);
151  samplePhysical = motherLogical->GetDaughter(sampleNo);
152  if ( samplePhysical!=blockedExitedVol )
153  {
154  G4AffineTransform sampleTf(samplePhysical->GetRotation(),
155  samplePhysical->GetTranslation());
156  sampleTf.Invert();
157  const G4ThreeVector samplePoint =
158  sampleTf.TransformPoint(localPoint);
159  const G4VSolid *sampleSolid =
160  samplePhysical->GetLogicalVolume()->GetSolid();
161  const G4double sampleSafety =
162  sampleSolid->DistanceToIn(samplePoint);
163 #ifdef G4VERBOSE
164  if( fCheck )
165  {
166  fLogger->PrintDaughterLog(sampleSolid,samplePoint,sampleSafety,0);
167  }
168 #endif
169  if ( sampleSafety<ourSafety )
170  {
171  ourSafety = sampleSafety;
172  }
173  if ( sampleSafety<=ourStep )
174  {
175  sampleDirection = sampleTf.TransformAxis(localDirection);
176  G4double sampleStep =
177  sampleSolid->DistanceToIn(samplePoint, sampleDirection);
178 #ifdef G4VERBOSE
179  if( fCheck )
180  {
181  fLogger->PrintDaughterLog(sampleSolid, samplePoint,
182  sampleSafety, sampleStep);
183  }
184 #endif
185  if ( sampleStep<=ourStep )
186  {
187  ourStep = sampleStep;
188  entering = true;
189  exiting = false;
190  *pBlockedPhysical = samplePhysical;
191  blockedReplicaNo = -1;
192 #ifdef G4VERBOSE
193  // Check to see that the resulting point is indeed in/on volume.
194  // This check could eventually be made only for successful
195  // candidate.
196 
197  if ( fCheck )
198  {
199  fLogger->AlongComputeStepLog (sampleSolid, samplePoint,
200  sampleDirection, localDirection, sampleSafety, sampleStep);
201  }
202 #endif
203  }
204  }
205  }
206  }
207  }
208  if (initialNode)
209  {
210  initialNode = false;
211  voxelSafety = ComputeVoxelSafety(localPoint);
212  if ( voxelSafety<ourSafety )
213  {
214  ourSafety = voxelSafety;
215  }
216  if ( currentProposedStepLength<ourSafety )
217  {
218  // Guaranteed physics limited
219  //
220  noStep = false;
221  entering = false;
222  exiting = false;
223  *pBlockedPhysical = 0;
224  ourStep = kInfinity;
225  }
226  else
227  {
228  //
229  // Compute mother intersection if required
230  //
231  if ( motherSafety<=ourStep )
232  {
233  G4double motherStep =
234  motherSolid->DistanceToOut(localPoint,
235  localDirection,
236  true, &validExitNormal, &exitNormal);
237 #ifdef G4VERBOSE
238  if ( fCheck )
239  {
240  fLogger->PostComputeStepLog(motherSolid, localPoint, localDirection,
241  motherStep, motherSafety);
242  }
243 #endif
244  if ( motherStep<=ourStep )
245  {
246  ourStep = motherStep;
247  exiting = true;
248  entering = false;
249  if ( validExitNormal )
250  {
251  const G4RotationMatrix *rot = motherPhysical->GetRotation();
252  if (rot)
253  {
254  exitNormal *= rot->inverse();
255  }
256  }
257  }
258  else
259  {
260  validExitNormal = false;
261  }
262  }
263  }
264  newSafety = ourSafety;
265  }
266  if (noStep)
267  {
268  noStep = LocateNextVoxel(localPoint, localDirection, ourStep);
269  }
270  } // end -while (noStep)- loop
271 
272  return ourStep;
273 }
274 
275 // ********************************************************************
276 // ComputeVoxelSafety
277 //
278 // Computes safety from specified point to voxel boundaries
279 // using already located point
280 // o collected boundaries for most derived level
281 // o adjacent boundaries for previous levels
282 // ********************************************************************
283 //
284 G4double
286 {
287  G4SmartVoxelHeader *curHeader;
288  G4double voxelSafety, curNodeWidth;
289  G4double curNodeOffset, minCurCommonDelta, maxCurCommonDelta;
290  G4int minCurNodeNoDelta, maxCurNodeNoDelta;
291  G4int localVoxelDepth, curNodeNo;
292  EAxis curHeaderAxis;
293 
294  localVoxelDepth = fVoxelDepth;
295 
296  curHeader = fVoxelHeaderStack[localVoxelDepth];
297  curHeaderAxis = fVoxelAxisStack[localVoxelDepth];
298  curNodeNo = fVoxelNodeNoStack[localVoxelDepth];
299  curNodeWidth = fVoxelSliceWidthStack[localVoxelDepth];
300 
301  // Compute linear intersection distance to boundaries of max/min
302  // to collected nodes at current level
303  //
304  curNodeOffset = curNodeNo*curNodeWidth;
305  maxCurNodeNoDelta = fVoxelNode->GetMaxEquivalentSliceNo()-curNodeNo;
306  minCurNodeNoDelta = curNodeNo-fVoxelNode->GetMinEquivalentSliceNo();
307  minCurCommonDelta = localPoint(curHeaderAxis)
308  - curHeader->GetMinExtent() - curNodeOffset;
309  maxCurCommonDelta = curNodeWidth-minCurCommonDelta;
310 
311  if ( minCurNodeNoDelta<maxCurNodeNoDelta )
312  {
313  voxelSafety = minCurNodeNoDelta*curNodeWidth;
314  voxelSafety += minCurCommonDelta;
315  }
316  else if (maxCurNodeNoDelta < minCurNodeNoDelta)
317  {
318  voxelSafety = maxCurNodeNoDelta*curNodeWidth;
319  voxelSafety += maxCurCommonDelta;
320  }
321  else // (maxCurNodeNoDelta == minCurNodeNoDelta)
322  {
323  voxelSafety = minCurNodeNoDelta*curNodeWidth;
324  voxelSafety += std::min(minCurCommonDelta,maxCurCommonDelta);
325  }
326 
327  // Compute isotropic safety to boundaries of previous levels
328  // [NOT to collected boundaries]
329  //
330  while ( (localVoxelDepth>0) && (voxelSafety>0) )
331  {
332  localVoxelDepth--;
333  curHeader = fVoxelHeaderStack[localVoxelDepth];
334  curHeaderAxis = fVoxelAxisStack[localVoxelDepth];
335  curNodeNo = fVoxelNodeNoStack[localVoxelDepth];
336  curNodeWidth = fVoxelSliceWidthStack[localVoxelDepth];
337  curNodeOffset = curNodeNo*curNodeWidth;
338  minCurCommonDelta = localPoint(curHeaderAxis)
339  - curHeader->GetMinExtent() - curNodeOffset;
340  maxCurCommonDelta = curNodeWidth-minCurCommonDelta;
341 
342  if ( minCurCommonDelta<voxelSafety )
343  {
344  voxelSafety = minCurCommonDelta;
345  }
346  if ( maxCurCommonDelta<voxelSafety )
347  {
348  voxelSafety = maxCurCommonDelta;
349  }
350  }
351  if ( voxelSafety<0 )
352  {
353  voxelSafety = 0;
354  }
355 
356  return voxelSafety;
357 }
358 
359 // ********************************************************************
360 // LocateNextVoxel
361 //
362 // Finds the next voxel from the current voxel and point
363 // in the specified direction
364 //
365 // Returns false if all voxels considered
366 // [current Step ends inside same voxel or leaves all voxels]
367 // true otherwise
368 // [the information on the next voxel is put into the set of
369 // fVoxel* variables & "stacks"]
370 // ********************************************************************
371 //
372 G4bool
374  const G4ThreeVector& localDirection,
375  const G4double currentStep)
376 {
377  G4SmartVoxelHeader *workHeader=0, *newHeader=0;
378  G4SmartVoxelProxy *newProxy=0;
379  G4SmartVoxelNode *newVoxelNode=0;
380  G4ThreeVector targetPoint, voxelPoint;
381  G4double workNodeWidth, workMinExtent, workCoord;
382  G4double minVal, maxVal, newDistance=0.;
383  G4double newHeaderMin, newHeaderNodeWidth;
384  G4int depth=0, newDepth=0, workNodeNo=0, newNodeNo=0, newHeaderNoSlices=0;
385  EAxis workHeaderAxis, newHeaderAxis;
386  G4bool isNewVoxel=false;
387 
388  G4double currentDistance = currentStep;
389  static const G4double sigma = 0.5*G4GeometryTolerance::GetInstance()
391 
392  // Determine if end of Step within current voxel
393  //
394  for (depth=0; depth<fVoxelDepth; depth++)
395  {
396  targetPoint = localPoint+localDirection*currentDistance;
397  newDistance = currentDistance;
398  workHeader = fVoxelHeaderStack[depth];
399  workHeaderAxis = fVoxelAxisStack[depth];
400  workNodeNo = fVoxelNodeNoStack[depth];
401  workNodeWidth = fVoxelSliceWidthStack[depth];
402  workMinExtent = workHeader->GetMinExtent();
403  workCoord = targetPoint(workHeaderAxis);
404  minVal = workMinExtent+workNodeNo*workNodeWidth;
405 
406  if ( minVal<=workCoord+sigma )
407  {
408  maxVal = minVal+workNodeWidth;
409  if ( maxVal<=workCoord-sigma )
410  {
411  // Must consider next voxel
412  //
413  newNodeNo = workNodeNo+1;
414  newHeader = workHeader;
415  newDistance = (maxVal-localPoint(workHeaderAxis))
416  / localDirection(workHeaderAxis);
417  isNewVoxel = true;
418  newDepth = depth;
419  }
420  }
421  else
422  {
423  newNodeNo = workNodeNo-1;
424  newHeader = workHeader;
425  newDistance = (minVal-localPoint(workHeaderAxis))
426  / localDirection(workHeaderAxis);
427  isNewVoxel = true;
428  newDepth = depth;
429  }
430  currentDistance = newDistance;
431  }
432  targetPoint = localPoint+localDirection*currentDistance;
433 
434  // Check if end of Step within collected boundaries of current voxel
435  //
436  depth = fVoxelDepth;
437  {
438  workHeader = fVoxelHeaderStack[depth];
439  workHeaderAxis = fVoxelAxisStack[depth];
440  workNodeNo = fVoxelNodeNoStack[depth];
441  workNodeWidth = fVoxelSliceWidthStack[depth];
442  workMinExtent = workHeader->GetMinExtent();
443  workCoord = targetPoint(workHeaderAxis);
444  minVal = workMinExtent+fVoxelNode->GetMinEquivalentSliceNo()*workNodeWidth;
445 
446  if ( minVal<=workCoord+sigma )
447  {
448  maxVal = workMinExtent+(fVoxelNode->GetMaxEquivalentSliceNo()+1)
449  *workNodeWidth;
450  if ( maxVal<=workCoord-sigma )
451  {
452  newNodeNo = fVoxelNode->GetMaxEquivalentSliceNo()+1;
453  newHeader = workHeader;
454  newDistance = (maxVal-localPoint(workHeaderAxis))
455  / localDirection(workHeaderAxis);
456  isNewVoxel = true;
457  newDepth = depth;
458  }
459  }
460  else
461  {
462  newNodeNo = fVoxelNode->GetMinEquivalentSliceNo()-1;
463  newHeader = workHeader;
464  newDistance = (minVal-localPoint(workHeaderAxis))
465  / localDirection(workHeaderAxis);
466  isNewVoxel = true;
467  newDepth = depth;
468  }
469  currentDistance = newDistance;
470  }
471  if (isNewVoxel)
472  {
473  // Compute new voxel & adjust voxel stack
474  //
475  // newNodeNo=Candidate node no at
476  // newDepth =refinement depth of crossed voxel boundary
477  // newHeader=Header for crossed voxel
478  // newDistance=distance to crossed voxel boundary (along the track)
479  //
480  if ( (newNodeNo<0) || (newNodeNo>=newHeader->GetNoSlices()))
481  {
482  // Leaving mother volume
483  //
484  isNewVoxel = false;
485  }
486  else
487  {
488  // Compute intersection point on the least refined
489  // voxel boundary that is hit
490  //
491  voxelPoint = localPoint+localDirection*newDistance;
492  fVoxelNodeNoStack[newDepth] = newNodeNo;
493  fVoxelDepth = newDepth;
494  newVoxelNode = 0;
495  while ( !newVoxelNode )
496  {
497  newProxy = newHeader->GetSlice(newNodeNo);
498  if (newProxy->IsNode())
499  {
500  newVoxelNode = newProxy->GetNode();
501  }
502  else
503  {
504  fVoxelDepth++;
505  newHeader = newProxy->GetHeader();
506  newHeaderAxis = newHeader->GetAxis();
507  newHeaderNoSlices = newHeader->GetNoSlices();
508  newHeaderMin = newHeader->GetMinExtent();
509  newHeaderNodeWidth = (newHeader->GetMaxExtent()-newHeaderMin)
510  / newHeaderNoSlices;
511  newNodeNo = G4int( (voxelPoint(newHeaderAxis)-newHeaderMin)
512  / newHeaderNodeWidth );
513  // Rounding protection
514  //
515  if ( newNodeNo<0 )
516  {
517  newNodeNo=0;
518  }
519  else if ( newNodeNo>=newHeaderNoSlices )
520  {
521  newNodeNo = newHeaderNoSlices-1;
522  }
523  // Stack info for stepping
524  //
525  fVoxelAxisStack[fVoxelDepth] = newHeaderAxis;
526  fVoxelNoSlicesStack[fVoxelDepth] = newHeaderNoSlices;
527  fVoxelSliceWidthStack[fVoxelDepth] = newHeaderNodeWidth;
528  fVoxelNodeNoStack[fVoxelDepth] = newNodeNo;
529  fVoxelHeaderStack[fVoxelDepth] = newHeader;
530  }
531  }
532  fVoxelNode = newVoxelNode;
533  }
534  }
535  return isNewVoxel;
536 }
537 
538 // ********************************************************************
539 // ComputeSafety
540 //
541 // Calculates the isotropic distance to the nearest boundary from the
542 // specified point in the local coordinate system.
543 // The localpoint utilised must be within the current volume.
544 // ********************************************************************
545 //
546 G4double
548  const G4NavigationHistory& history,
549  const G4double maxLength)
550 {
551  G4VPhysicalVolume *motherPhysical, *samplePhysical;
552  G4LogicalVolume *motherLogical;
553  G4VSolid *motherSolid;
554  G4double motherSafety, ourSafety;
555  G4int sampleNo;
556  G4SmartVoxelNode *curVoxelNode;
557  G4int curNoVolumes, contentNo;
558  G4double voxelSafety;
559 
560  motherPhysical = history.GetTopVolume();
561  motherLogical = motherPhysical->GetLogicalVolume();
562  motherSolid = motherLogical->GetSolid();
563 
564  if( fBestSafety )
565  {
566  return fpVoxelSafety->ComputeSafety( localPoint,*motherPhysical,maxLength );
567  }
568 
569  //
570  // Compute mother safety
571  //
572 
573  motherSafety = motherSolid->DistanceToOut(localPoint);
574  ourSafety = motherSafety; // Working isotropic safety
575 
576  if( motherSafety == 0.0 )
577  {
578 #ifdef G4DEBUG_NAVIGATION
579  // Check that point is inside mother volume
580  EInside insideMother= motherSolid->Inside(localPoint);
581 
582  if( insideMother == kOutside )
583  {
584  G4ExceptionDescription message;
585  message << "Safety method called for location outside current Volume." << G4endl
586  << "Location for safety is Outside this volume. " << G4endl
587  << "The approximate distance to the solid "
588  << "(safety from outside) is: "
589  << motherSolid->DistanceToIn( localPoint ) << G4endl;
590  message << " Problem occurred with physical volume: "
591  << " Name: " << motherPhysical->GetName()
592  << " Copy No: " << motherPhysical->GetCopyNo() << G4endl
593  << " Local Point = " << localPoint << G4endl;
594  message << " Description of solid: " << G4endl
595  << *motherSolid << G4endl;
596  G4Exception("G4VoxelNavigation::ComputeSafety()", "GeomNav0003",
597  JustWarning, // FatalException,
598  message);
599  }
600 
601  // Following check is NOT for an issue - it is only for information
602  // It is allowed that a solid gives approximate safety - even zero.
603  //
604  if( insideMother == kInside ) // && fVerbose )
605  {
606  G4ExceptionDescription messageIn;
607 
608  messageIn << " Point is Inside, but safety is Zero ." << G4endl;
609  messageIn << " Inexact safety for volume " << motherPhysical->GetName() << G4endl
610  << " Solid: Name= " << motherSolid->GetName()
611  << " Type= " << motherSolid->GetEntityType() << G4endl;
612  messageIn << " Local point= " << localPoint << G4endl;
613  messageIn << " Solid parameters: " << G4endl << *motherSolid << G4endl;
614  G4Exception("G4VoxelNavigation::ComputeSafety()", "GeomNav0003",
615  JustWarning, messageIn);
616  }
617 #endif
618  // if( insideMother != kInside )
619  return 0.0;
620  }
621 
622 #ifdef G4VERBOSE
623  if( fCheck )
624  {
625  fLogger->ComputeSafetyLog (motherSolid, localPoint, motherSafety, true);
626  }
627 #endif
628  //
629  // Compute daughter safeties
630  //
631  // Look only inside the current Voxel only (in the first version).
632  //
633  curVoxelNode = fVoxelNode;
634  curNoVolumes = curVoxelNode->GetNoContained();
635 
636  for ( contentNo=curNoVolumes-1; contentNo>=0; contentNo-- )
637  {
638  sampleNo = curVoxelNode->GetVolume(contentNo);
639  samplePhysical = motherLogical->GetDaughter(sampleNo);
640 
641  G4AffineTransform sampleTf(samplePhysical->GetRotation(),
642  samplePhysical->GetTranslation());
643  sampleTf.Invert();
644  const G4ThreeVector samplePoint =
645  sampleTf.TransformPoint(localPoint);
646  const G4VSolid *sampleSolid =
647  samplePhysical->GetLogicalVolume()->GetSolid();
648  G4double sampleSafety = sampleSolid->DistanceToIn(samplePoint);
649  if ( sampleSafety<ourSafety )
650  {
651  ourSafety = sampleSafety;
652  }
653 #ifdef G4VERBOSE
654  if( fCheck )
655  {
656  fLogger->ComputeSafetyLog (sampleSolid,samplePoint,sampleSafety,false);
657  }
658 #endif
659  }
660  voxelSafety = ComputeVoxelSafety(localPoint);
661  if ( voxelSafety<ourSafety )
662  {
663  ourSafety = voxelSafety;
664  }
665  return ourSafety;
666 }
667 
668 // ********************************************************************
669 // SetVerboseLevel
670 // ********************************************************************
671 //
673 {
674  if( fLogger ) fLogger->SetVerboseLevel(level);
676 }
void SetVerboseLevel(G4int level)
void AlongComputeStepLog(const G4VSolid *sampleSolid, const G4ThreeVector &samplePoint, const G4ThreeVector &sampleDirection, const G4ThreeVector &localDirection, G4double sampleSafety, G4double sampleStep) const
G4String GetName() const
G4SmartVoxelHeader * GetHeader() const
G4VPhysicalVolume * GetTopVolume() const
G4int GetMinEquivalentSliceNo() const
const G4ThreeVector & GetTranslation() const
static const G4double kInfinity
Definition: geomdefs.hh:42
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
CLHEP::Hep3Vector G4ThreeVector
CLHEP::HepRotation G4RotationMatrix
void PrintDaughterLog(const G4VSolid *sampleSolid, const G4ThreeVector &samplePoint, G4double sampleSafety, G4double sampleStep) const
G4int GetVolume(G4int pVolumeNo) const
G4double ComputeVoxelSafety(const G4ThreeVector &localPoint) const
G4double GetSurfaceTolerance() const
void ComputeSafetyLog(const G4VSolid *solid, const G4ThreeVector &point, G4double safety, G4bool banner) const
G4VPhysicalVolume * GetDaughter(const G4int i) const
G4double ComputeSafety(const G4ThreeVector &localPoint, const G4VPhysicalVolume &currentPhysical, G4double maxLength=DBL_MAX)
virtual G4GeometryType GetEntityType() const =0
std::vector< G4SmartVoxelHeader * > fVoxelHeaderStack
int G4int
Definition: G4Types.hh:78
std::vector< EAxis > fVoxelAxisStack
G4int GetNoContained() const
G4int GetMaxEquivalentSliceNo() const
G4bool IsNode() const
std::vector< G4double > fVoxelSliceWidthStack
G4AffineTransform & Invert()
EAxis GetAxis() const
G4BlockingList fBList
G4double GetMinExtent() const
const G4String & GetName() const
virtual EInside Inside(const G4ThreeVector &p) const =0
bool G4bool
Definition: G4Types.hh:79
G4bool LocateNextVoxel(const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4double currentStep)
std::vector< G4int > fVoxelNodeNoStack
void BlockVolume(const G4int v)
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4NavigationHistory &history, const G4double pMaxLength=DBL_MAX)
G4bool IsBlocked(const G4int v) const
std::vector< G4int > fVoxelNoSlicesStack
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
static const double kMinExitingNormalCosine
Definition: geomdefs.hh:46
G4NavigationLogger * fLogger
G4int GetNoDaughters() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void SetVerboseLevel(G4int level)
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
virtual ~G4VoxelNavigation()
G4LogicalVolume * GetLogicalVolume() const
G4SmartVoxelNode * GetNode() const
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
void SetVerboseLevel(G4int level)
void PreComputeStepLog(const G4VPhysicalVolume *motherPhysical, G4double motherSafety, const G4ThreeVector &localPoint) const
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
G4VoxelSafety * fpVoxelSafety
static const G4int kNavigatorVoxelStackMax
Definition: geomdefs.hh:80
G4SmartVoxelNode * fVoxelNode
#define G4endl
Definition: G4ios.hh:61
void PostComputeStepLog(const G4VSolid *motherSolid, const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, G4double motherStep, G4double motherSafety) 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
virtual G4double ComputeStep(const G4ThreeVector &globalPoint, const G4ThreeVector &globalDirection, const G4double currentProposedStepLength, G4double &newSafety, G4NavigationHistory &history, G4bool &validExitNormal, G4ThreeVector &exitNormal, G4bool &exiting, G4bool &entering, G4VPhysicalVolume *(*pBlockedPhysical), G4int &blockedReplicaNo)
static G4GeometryTolerance * GetInstance()
G4VSolid * GetSolid() const
void Enlarge(const G4int nv)