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