Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VoxelSafety Class Reference

#include <G4VoxelSafety.hh>

Public Member Functions

 G4VoxelSafety ()
 
 ~G4VoxelSafety ()
 
G4SmartVoxelNodeVoxelLocate (G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint)
 
G4double ComputeSafety (const G4ThreeVector &localPoint, const G4VPhysicalVolume &currentPhysical, G4double maxLength=DBL_MAX)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int level)
 

Protected Member Functions

G4double SafetyForVoxelHeader (const G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint, G4double maxLength, const G4VPhysicalVolume &currentPhysical, G4double distUpperDepth=0.0, G4double previousMinSafety=DBL_MAX)
 
G4double SafetyForVoxelNode (const G4SmartVoxelNode *curVoxelNode, const G4ThreeVector &localPoint)
 
G4SmartVoxelNodeVoxelLocateLight (G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint) const
 

Detailed Description

Definition at line 55 of file G4VoxelSafety.hh.

Constructor & Destructor Documentation

◆ G4VoxelSafety()

G4VoxelSafety::G4VoxelSafety ( )

Definition at line 44 of file G4VoxelSafety.cc.

45 : fBlockList(),
46 fVoxelAxisStack(kNavigatorVoxelStackMax,kXAxis),
47 fVoxelNoSlicesStack(kNavigatorVoxelStackMax,0),
48 fVoxelSliceWidthStack(kNavigatorVoxelStackMax,0.),
49 fVoxelNodeNoStack(kNavigatorVoxelStackMax,0),
50 fVoxelHeaderStack(kNavigatorVoxelStackMax,(G4SmartVoxelHeader*)nullptr)
51{
53}
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
@ kXAxis
Definition: geomdefs.hh:55

◆ ~G4VoxelSafety()

G4VoxelSafety::~G4VoxelSafety ( )

Definition at line 59 of file G4VoxelSafety.cc.

60{
61}

Member Function Documentation

◆ ComputeSafety()

G4double G4VoxelSafety::ComputeSafety ( const G4ThreeVector localPoint,
const G4VPhysicalVolume currentPhysical,
G4double  maxLength = DBL_MAX 
)

Definition at line 72 of file G4VoxelSafety.cc.

75{
76 G4LogicalVolume *motherLogical;
77 G4VSolid *motherSolid;
78 G4SmartVoxelHeader *motherVoxelHeader;
79 G4double motherSafety, ourSafety;
80 G4int localNoDaughters;
81 G4double daughterSafety;
82
83 motherLogical = currentPhysical.GetLogicalVolume();
84 fpMotherLogical= motherLogical; // For use by the other methods
85 motherSolid = motherLogical->GetSolid();
86 motherVoxelHeader = motherLogical->GetVoxelHeader();
87
88#ifdef G4VERBOSE
89 if( fVerbose > 0 )
90 {
91 G4cout << "*** G4VoxelSafety::ComputeSafety(): ***" << G4endl;
92 }
93#endif
94
95 // Check that point is inside mother volume
96 //
97 EInside insideMother = motherSolid->Inside(localPoint);
98 if( insideMother != kInside )
99 {
100#ifdef G4DEBUG_NAVIGATION
101 if( insideMother == kOutside )
102 {
103 std::ostringstream message;
104 message << "Safety method called for location outside current Volume."
105 << G4endl
106 << "Location for safety is Outside this volume. " << G4endl
107 << "The approximate distance to the solid "
108 << "(safety from outside) is: "
109 << motherSolid->DistanceToIn( localPoint ) << G4endl;
110 message << " Problem occurred with physical volume: "
111 << " Name: " << currentPhysical.GetName()
112 << " Copy No: " << currentPhysical.GetCopyNo() << G4endl
113 << " Local Point = " << localPoint << G4endl;
114 message << " Description of solid: " << G4endl
115 << *motherSolid << G4endl;
116 G4Exception("G4VoxelSafety::ComputeSafety()", "GeomNav0003",
117 FatalException, message);
118 }
119#endif
120 return 0.0;
121 }
122
123 // First limit: mother safety - distance to outer boundaries
124 //
125 motherSafety = motherSolid->DistanceToOut(localPoint);
126 ourSafety = motherSafety; // Working isotropic safety
127
128#ifdef G4VERBOSE
129 if(( fCheck ) ) // && ( fVerbose == 1 ))
130 {
131 G4cout << " Invoked DistanceToOut(p) for mother solid: "
132 << motherSolid->GetName()
133 << ". Solid replied: " << motherSafety << G4endl
134 << " For local point p: " << localPoint
135 << ", to be considered as 'mother safety'." << G4endl;
136 }
137#endif
138 localNoDaughters = motherLogical->GetNoDaughters();
139
140 fBlockList.Enlarge(localNoDaughters);
141 fBlockList.Reset();
142
143 fVoxelDepth = -1; // Resets the depth -- must be done for next method
144 daughterSafety= SafetyForVoxelHeader(motherVoxelHeader, localPoint, maxLength,
145 currentPhysical, 0.0, ourSafety);
146 ourSafety= std::min( motherSafety, daughterSafety );
147
148 return ourSafety;
149}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void Enlarge(const G4int nv)
G4VSolid * GetSolid() const
size_t GetNoDaughters() const
G4SmartVoxelHeader * GetVoxelHeader() const
G4LogicalVolume * GetLogicalVolume() const
virtual G4int GetCopyNo() const =0
const G4String & GetName() const
G4String GetName() const
virtual EInside Inside(const G4ThreeVector &p) const =0
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
G4double SafetyForVoxelHeader(const G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint, G4double maxLength, const G4VPhysicalVolume &currentPhysical, G4double distUpperDepth=0.0, G4double previousMinSafety=DBL_MAX)
EInside
Definition: geomdefs.hh:67
@ kInside
Definition: geomdefs.hh:70
@ kOutside
Definition: geomdefs.hh:68

Referenced by G4Navigator::ComputeSafety(), G4ITNavigator1::ComputeSafety(), G4ITNavigator2::ComputeSafety(), and G4VoxelNavigation::ComputeSafety().

◆ GetVerboseLevel()

G4int G4VoxelSafety::GetVerboseLevel ( ) const
inline

Definition at line 69 of file G4VoxelSafety.hh.

69{ return fVerbose; }

◆ SafetyForVoxelHeader()

G4double G4VoxelSafety::SafetyForVoxelHeader ( const G4SmartVoxelHeader pHead,
const G4ThreeVector localPoint,
G4double  maxLength,
const G4VPhysicalVolume currentPhysical,
G4double  distUpperDepth = 0.0,
G4double  previousMinSafety = DBL_MAX 
)
protected

Definition at line 213 of file G4VoxelSafety.cc.

220{
221 const G4SmartVoxelHeader* const targetVoxelHeader = pHeader;
222 G4SmartVoxelNode* targetVoxelNode = nullptr;
223
224 const G4SmartVoxelProxy* sampleProxy;
225 EAxis targetHeaderAxis;
226 G4double targetHeaderMin, targetHeaderMax, targetHeaderNodeWidth;
227 G4int targetHeaderNoSlices;
228 G4int targetNodeNo;
229
230 G4double minSafety = previousMinSafety;
231 G4double ourSafety = DBL_MAX;
232 unsigned int checkedNum= 0;
233
234 ++fVoxelDepth;
235 // fVoxelDepth set by ComputeSafety or previous level call
236
237 targetHeaderAxis = targetVoxelHeader->GetAxis();
238 targetHeaderNoSlices = targetVoxelHeader->GetNoSlices();
239 targetHeaderMin = targetVoxelHeader->GetMinExtent();
240 targetHeaderMax = targetVoxelHeader->GetMaxExtent();
241
242 targetHeaderNodeWidth = (targetHeaderMax-targetHeaderMin)
243 / targetHeaderNoSlices;
244
245 G4double localCrd = localPoint(targetHeaderAxis);
246
247 const G4int candNodeNo = G4int( (localCrd-targetHeaderMin)
248 / targetHeaderNodeWidth );
249 // Ensure that it is between 0 and targetHeader->GetMaxExtent() - 1
250
251#ifdef G4DEBUG_VOXELISATION
252 if( candNodeNo < 0 || candNodeNo > targetHeaderNoSlices-1 )
253 {
255 ed << " Potential ERROR."
256 << " Point is outside range of Voxel in current coordinate" << G4endl;
257 ed << " Node number of point " << localPoint
258 << "is outside the range. " << G4endl;
259 ed << " Voxel node Num= " << candNodeNo << " versus minimum= " << 0
260 << " and maximum= " << targetHeaderNoSlices-1 << G4endl;
261 ed << " Axis = " << targetHeaderAxis
262 << " No of slices = " << targetHeaderNoSlices << G4endl;
263 ed << " Local coord = " << localCrd
264 << " Voxel Min = " << targetHeaderMin
265 << " Max = " << targetHeaderMax << G4endl;
266 G4LogicalVolume *pLogical= currentPhysical.GetLogicalVolume();
267 ed << " Current volume (physical) = " << currentPhysical.GetName()
268 << " (logical) = " << pLogical->GetName() << G4endl;
269 G4VSolid* pSolid= pLogical->GetSolid();
270 ed << " Solid type = " << pSolid->GetEntityType() << G4endl;
271 ed << *pSolid << G4endl;
272 G4Exception("G4VoxelSafety::SafetyForVoxelHeader()", "GeomNav1003",
273 JustWarning, ed,
274 "Point is outside range of Voxel in current coordinate");
275 }
276#endif
277
278 const G4int pointNodeNo =
279 std::max( 0, std::min( candNodeNo, targetHeaderNoSlices-1 ) );
280
281#ifdef G4VERBOSE
282 if( fVerbose > 2 )
283 {
284 G4cout << G4endl;
285 G4cout << "**** G4VoxelSafety::SafetyForVoxelHeader " << G4endl;
286 G4cout << " Called at Depth = " << fVoxelDepth;
287 G4cout << " Calculated pointNodeNo= " << pointNodeNo
288 << " from position= " << localPoint(targetHeaderAxis)
289 << " min= " << targetHeaderMin
290 << " max= " << targetHeaderMax
291 << " width= " << targetHeaderNodeWidth
292 << " no-slices= " << targetHeaderNoSlices
293 << " axis= " << targetHeaderAxis << G4endl;
294 }
295 else if (fVerbose == 1)
296 {
297 G4cout << " VoxelSafety: Depth = " << fVoxelDepth
298 << " Number of Slices = " << targetHeaderNoSlices
299 << " Header (address) = " << targetVoxelHeader << G4endl;
300 }
301#endif
302
303 // Stack info for stepping
304 //
305 fVoxelAxisStack[fVoxelDepth] = targetHeaderAxis;
306 fVoxelNoSlicesStack[fVoxelDepth] = targetHeaderNoSlices;
307 fVoxelSliceWidthStack[fVoxelDepth] = targetHeaderNodeWidth;
308
309 fVoxelHeaderStack[fVoxelDepth] = pHeader;
310
311 G4int trialUp = -1, trialDown = -1;
312 G4double distUp = DBL_MAX, distDown = DBL_MAX;
313
314 // Using Equivalent voxels - this is pre-initialisation only
315 //
316 G4int nextUp = pointNodeNo+1;
317 G4int nextDown = pointNodeNo-1;
318
319 G4int nextNodeNo = pointNodeNo;
320 G4double distAxis; // Distance in current Axis
321 distAxis = 0.0; // Starting in node containing local Coordinate
322
323 G4bool nextIsInside = false;
324
325 G4double distMaxInterest= std::min( previousMinSafety, maxLength);
326 // We will not look beyond this distance.
327 // This distance will be updated to reflect the
328 // max ( minSafety, maxLength ) at each step
329
330 targetNodeNo = pointNodeNo;
331 do
332 {
333 G4double nodeSafety = DBL_MAX, headerSafety = DBL_MAX;
334 fVoxelNodeNoStack[fVoxelDepth] = targetNodeNo;
335
336 ++checkedNum;
337
338 sampleProxy = targetVoxelHeader->GetSlice(targetNodeNo);
339
340#ifdef G4DEBUG_NAVIGATION
341 assert( sampleProxy != 0);
342 if( fVerbose > 2 )
343 {
344 G4cout << " -Checking node " << targetNodeNo
345 << " is proxy with address " << sampleProxy << G4endl;
346 }
347#endif
348
349 if ( sampleProxy == 0 )
350 {
352 ed << " Problem for node number= " << targetNodeNo
353 << " Number of slides = " << targetHeaderNoSlices
354 << G4endl;
355 G4Exception( "G4VoxelSafety::SafetyForVoxelHeader()", "GeomNav0003",
356 FatalException, ed,
357 "Problem sampleProxy is Zero. Failure in loop.");
358 }
359 else if ( sampleProxy->IsNode() )
360 {
361 targetVoxelNode = sampleProxy->GetNode();
362
363 // Deal with the node here [ i.e. the last level ]
364 //
365 nodeSafety= SafetyForVoxelNode( targetVoxelNode, localPoint);
366#ifdef G4DEBUG_NAVIGATION
367 if( fVerbose > 2 )
368 {
369 G4cout << " -- It is a Node ";
370 G4cout << " its safety= " << nodeSafety
371 << " our level Saf = " << ourSafety
372 << " MinSaf= " << minSafety << G4endl;
373 }
374#endif
375 ourSafety= std::min( ourSafety, nodeSafety );
376
377 trialUp = targetVoxelNode->GetMaxEquivalentSliceNo()+1;
378 trialDown = targetVoxelNode->GetMinEquivalentSliceNo()-1;
379 }
380 else
381 {
382 const G4SmartVoxelHeader* pNewVoxelHeader = sampleProxy->GetHeader();
383
384 G4double distCombined_Sq;
385 distCombined_Sq = distUpperDepth_Sq + distAxis*distAxis;
386
387#ifdef G4DEBUG_NAVIGATION
388 if( fVerbose > 2 )
389 {
390 G4double distCombined= std::sqrt( distCombined_Sq );
391 G4double distUpperDepth= std::sqrt ( distUpperDepth_Sq );
392 G4cout << " -- It is a Header " << G4endl;
393 G4cout << " Recurse to deal with next level, fVoxelDepth= "
394 << fVoxelDepth+1 << G4endl;
395 G4cout << " Distances: Upper= " << distUpperDepth
396 << " Axis= " << distAxis
397 << " Combined= " << distCombined << G4endl;
398 }
399#endif
400
401 // Recurse to deal with lower levels
402 //
403 headerSafety= SafetyForVoxelHeader( pNewVoxelHeader, localPoint,
404 maxLength, currentPhysical,
405 distCombined_Sq, minSafety);
406 ourSafety = std::min( ourSafety, headerSafety );
407
408#ifdef G4DEBUG_NAVIGATION
409 if( fVerbose > 2 )
410 {
411 G4cout << " >> Header safety = " << headerSafety
412 << " our level Saf = " << ourSafety << G4endl;
413 }
414#endif
415 trialUp = pNewVoxelHeader->GetMaxEquivalentSliceNo()+1;
416 trialDown = pNewVoxelHeader->GetMinEquivalentSliceNo()-1;
417 }
418 minSafety = std::min( minSafety, ourSafety );
419
420 // Find next closest Voxel
421 // - first try: by simple subtraction
422 // - later: using distance (TODO - tbc)
423 //
424 if( targetNodeNo >= pointNodeNo )
425 {
426 nextUp = trialUp;
427 // distUp = std::max( targetHeaderMax-localCrd, 0.0 );
428 G4double lowerEdgeOfNext = targetHeaderMin
429 + nextUp * targetHeaderNodeWidth;
430 distUp = lowerEdgeOfNext-localCrd ;
431 if( distUp < 0.0 )
432 {
433 distUp = DBL_MAX; // On the wrong side - must not be considered
434 }
435#ifdef G4DEBUG_NAVIGATION
436 if( fVerbose > 2 )
437 {
438 G4cout << " > Updated nextUp= " << nextUp << G4endl;
439 }
440#endif
441 }
442
443 if( targetNodeNo <= pointNodeNo )
444 {
445 nextDown = trialDown;
446 // distDown = std::max( localCrd-targetHeaderMin, 0.0);
447 G4double upperEdgeOfNext = targetHeaderMin
448 + (nextDown+1) * targetHeaderNodeWidth;
449 distDown = localCrd-upperEdgeOfNext;
450 if( distDown < 0.0 )
451 {
452 distDown= DBL_MAX; // On the wrong side - must not be considered
453 }
454#ifdef G4DEBUG_NAVIGATION
455 if( fVerbose > 2 )
456 {
457 G4cout << " > Updated nextDown= " << nextDown << G4endl;
458 }
459#endif
460 }
461
462#ifdef G4DEBUG_NAVIGATION
463 if( fVerbose > 2 )
464 {
465 G4cout << " Node= " << pointNodeNo
466 << " Up: next= " << nextUp << " d# "
467 << nextUp - pointNodeNo
468 << " trialUp= " << trialUp << " d# "
469 << trialUp - pointNodeNo
470 << " Down: next= " << nextDown << " d# "
471 << targetNodeNo - nextDown
472 << " trialDwn= " << trialDown << " d# "
473 << targetNodeNo - trialDown
474 << " condition (next is Inside)= " << nextIsInside
475 << G4endl;
476 }
477#endif
478
479 G4bool UpIsClosest;
480 UpIsClosest = distUp < distDown;
481
482 if( (nextUp < targetHeaderNoSlices)
483 && (UpIsClosest || (nextDown < 0)) )
484 {
485 nextNodeNo = nextUp;
486 distAxis = distUp;
487 ++nextUp; // Default
488#ifdef G4VERBOSE
489 if( fVerbose > 2 )
490 {
491 G4cout << " > Chose Up. Depth= " << fVoxelDepth
492 << " Nodes: next= " << nextNodeNo
493 << " new nextUp= " << nextUp
494 << " Dist = " << distAxis << G4endl;
495 }
496#endif
497 }
498 else
499 {
500 nextNodeNo = nextDown;
501 distAxis = distDown;
502 --nextDown; // A default value
503#ifdef G4VERBOSE
504 if( fVerbose > 2 )
505 {
506 G4cout << " > Chose Down. Depth= " << fVoxelDepth
507 << " Nodes: next= " << nextNodeNo
508 << " new nextDown= " << nextDown
509 << " Dist = " << distAxis << G4endl;
510 }
511#endif
512 }
513
514 nextIsInside = (nextNodeNo >= 0) && (nextNodeNo < targetHeaderNoSlices);
515 if( nextIsInside )
516 {
517 targetNodeNo= nextNodeNo;
518
519#ifdef G4DEBUG_NAVIGATION
520 assert( targetVoxelHeader->GetSlice(nextNodeNo) != 0 );
521 G4bool bContinue = (distAxis<minSafety);
522 if( !bContinue )
523 {
524 if( fVerbose > 2 )
525 {
526 G4cout << " Can skip remaining at depth " << targetHeaderAxis
527 << " >> distAxis= " << distAxis
528 << " minSafety= " << minSafety << G4endl;
529 }
530 }
531#endif
532 }
533 else
534 {
535#ifdef G4DEBUG_NAVIGATION
536 if( fVerbose > 2)
537 {
538 G4cout << " VoxSaf> depth= " << fVoxelDepth << G4endl;
539 G4cout << " VoxSaf> No more candidates: nodeDown= " << nextDown
540 << " nodeUp= " << nextUp
541 << " lastSlice= " << targetHeaderNoSlices << G4endl;
542 }
543#endif
544 }
545
546 // This calculation can be 'hauled'-up to where minSafety is calculated
547 //
548 distMaxInterest = std::min( minSafety, maxLength );
549
550 } while ( nextIsInside && ( distAxis*distAxis + distUpperDepth_Sq
551 < distMaxInterest*distMaxInterest ) );
552
553#ifdef G4VERBOSE
554 if( fVerbose > 0 )
555 {
556 G4cout << " Ended for targetNodeNo -- checked " << checkedNum << " out of "
557 << targetHeaderNoSlices << " slices." << G4endl;
558 G4cout << " ===== Returning from SafetyForVoxelHeader "
559 << " Depth= " << fVoxelDepth << G4endl
560 << G4endl;
561 }
562#endif
563
564 // Go back one level
565 fVoxelDepth--;
566
567 return ourSafety;
568}
@ JustWarning
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
bool G4bool
Definition: G4Types.hh:86
const G4String & GetName() const
G4int GetMaxEquivalentSliceNo() const
G4double GetMaxExtent() const
G4double GetMinExtent() const
EAxis GetAxis() const
G4SmartVoxelProxy * GetSlice(G4int n) const
size_t GetNoSlices() const
G4int GetMinEquivalentSliceNo() const
G4int GetMaxEquivalentSliceNo() const
G4int GetMinEquivalentSliceNo() const
G4bool IsNode() const
G4SmartVoxelNode * GetNode() const
G4SmartVoxelHeader * GetHeader() const
virtual G4GeometryType GetEntityType() const =0
G4double SafetyForVoxelNode(const G4SmartVoxelNode *curVoxelNode, const G4ThreeVector &localPoint)
EAxis
Definition: geomdefs.hh:54
#define DBL_MAX
Definition: templates.hh:62

Referenced by ComputeSafety(), and SafetyForVoxelHeader().

◆ SafetyForVoxelNode()

G4double G4VoxelSafety::SafetyForVoxelNode ( const G4SmartVoxelNode curVoxelNode,
const G4ThreeVector localPoint 
)
protected

Definition at line 158 of file G4VoxelSafety.cc.

160{
161 G4double ourSafety = DBL_MAX;
162
163 G4int curNoVolumes, contentNo, sampleNo;
164 G4VPhysicalVolume* samplePhysical;
165
166 G4double sampleSafety = 0.0;
167 G4ThreeVector samplePoint;
168 G4VSolid* ptrSolid = nullptr;
169
170 curNoVolumes = curVoxelNode->GetNoContained();
171
172 for ( contentNo=curNoVolumes-1; contentNo>=0; contentNo-- )
173 {
174 sampleNo = curVoxelNode->GetVolume(contentNo);
175 if ( !fBlockList.IsBlocked(sampleNo) )
176 {
177 fBlockList.BlockVolume(sampleNo);
178
179 samplePhysical = fpMotherLogical->GetDaughter(sampleNo);
180 G4AffineTransform sampleTf(samplePhysical->GetRotation(),
181 samplePhysical->GetTranslation());
182 sampleTf.Invert();
183 samplePoint = sampleTf.TransformPoint(localPoint);
184 ptrSolid = samplePhysical->GetLogicalVolume()->GetSolid();
185
186 sampleSafety = ptrSolid->DistanceToIn(samplePoint);
187 ourSafety = std::min( sampleSafety, ourSafety );
188#ifdef G4VERBOSE
189 if(( fCheck ) && ( fVerbose == 1 ))
190 {
191 G4cout << "*** G4VoxelSafety::SafetyForVoxelNode(): ***" << G4endl
192 << " Invoked DistanceToIn(p) for daughter solid: "
193 << ptrSolid->GetName()
194 << ". Solid replied: " << sampleSafety << G4endl
195 << " For local point p: " << samplePoint
196 << ", to be considered as 'daughter safety'." << G4endl;
197 }
198#endif
199 }
200 } // end for contentNo
201
202 return ourSafety;
203}
void BlockVolume(const G4int v)
G4bool IsBlocked(const G4int v) const
G4VPhysicalVolume * GetDaughter(const G4int i) const
G4int GetVolume(G4int pVolumeNo) const
size_t GetNoContained() const
const G4RotationMatrix * GetRotation() const
const G4ThreeVector GetTranslation() const

Referenced by SafetyForVoxelHeader().

◆ SetVerboseLevel()

void G4VoxelSafety::SetVerboseLevel ( G4int  level)
inline

Definition at line 70 of file G4VoxelSafety.hh.

70{ fVerbose = level; }

Referenced by G4VoxelNavigation::SetVerboseLevel().

◆ VoxelLocate()

G4SmartVoxelNode * G4VoxelSafety::VoxelLocate ( G4SmartVoxelHeader pHead,
const G4ThreeVector localPoint 
)

◆ VoxelLocateLight()

G4SmartVoxelNode * G4VoxelSafety::VoxelLocateLight ( G4SmartVoxelHeader pHead,
const G4ThreeVector localPoint 
) const
protected

The documentation for this class was generated from the following files: