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

G4SafetyCalculator is a class that provides an estimate of the isotropic safety (the minimum distance from a global point to the nearest boundary of the current volume or the nearest daughter volumes). More...

#include <G4SafetyCalculator.hh>

Public Member Functions

 G4SafetyCalculator (const G4Navigator &navigator, const G4NavigationHistory &navHistory)
 G4SafetyCalculator (const G4SafetyCalculator &)=delete
G4SafetyCalculatoroperator= (const G4SafetyCalculator &)=delete
 ~G4SafetyCalculator ()=default
G4double SafetyInCurrentVolume (const G4ThreeVector &globalpoint, G4VPhysicalVolume *physicalVolume, const G4double pProposedMaxLength=DBL_MAX, G4bool verbose=false)
G4VExternalNavigationGetExternalNavigation () const
void SetExternalNavigation (G4VExternalNavigation *externalNav)
void CompareSafetyValues (G4double oldSafety, G4double newValue, G4VPhysicalVolume *motherPhysical, const G4ThreeVector &globalPoint, G4bool keepState, G4double maxLength, G4bool enteredVolume, G4bool exitedVolume)

Protected Member Functions

void QuickLocateWithinVolume (const G4ThreeVector &pointLocal, G4VPhysicalVolume *motherPhysical)
G4ThreeVector ComputeLocalPoint (const G4ThreeVector &rGlobPoint) const
G4ThreeVector ComputeLocalAxis (const G4ThreeVector &pVec) const
EVolume CharacteriseDaughters (const G4LogicalVolume *pLog) const
G4int GetDaughtersRegularStructureId (const G4LogicalVolume *pLv) const

Detailed Description

G4SafetyCalculator is a class that provides an estimate of the isotropic safety (the minimum distance from a global point to the nearest boundary of the current volume or the nearest daughter volumes).

Definition at line 72 of file G4SafetyCalculator.hh.

Constructor & Destructor Documentation

◆ G4SafetyCalculator() [1/2]

G4SafetyCalculator::G4SafetyCalculator ( const G4Navigator & navigator,
const G4NavigationHistory & navHistory )

Constructor, initialisers and setup.

Definition at line 36 of file G4SafetyCalculator.cc.

39 : fNavigator(navigator), fNavHistory(navHistory)
40{
42}
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()

Referenced by G4SafetyCalculator(), and operator=().

◆ G4SafetyCalculator() [2/2]

G4SafetyCalculator::G4SafetyCalculator ( const G4SafetyCalculator & )
delete

Copy constructor & assignment operator not allowed.

◆ ~G4SafetyCalculator()

G4SafetyCalculator::~G4SafetyCalculator ( )
default

Destructor. No actions.

Member Function Documentation

◆ CharacteriseDaughters()

EVolume G4SafetyCalculator::CharacteriseDaughters ( const G4LogicalVolume * pLog) const
inlineprotected

Characterises the daughter of logical volume.

Referenced by QuickLocateWithinVolume(), and SafetyInCurrentVolume().

◆ CompareSafetyValues()

void G4SafetyCalculator::CompareSafetyValues ( G4double oldSafety,
G4double newValue,
G4VPhysicalVolume * motherPhysical,
const G4ThreeVector & globalPoint,
G4bool keepState,
G4double maxLength,
G4bool enteredVolume,
G4bool exitedVolume )

Compares estimates of the safety, and reports if found difference(s).

Definition at line 210 of file G4SafetyCalculator.cc.

219{
220 constexpr G4double reportThreshold= 3.0e-14;
221 // At least warn if rel-error exceeds it
222 constexpr G4double errorThreshold= 1.0e-08;
223 // Fatal if relative error is larger
224 constexpr G4double epsilonLen= 1.0e-20;
225 // Baseline minimal value for divisor
226
227 const G4double oldSafetyPlus = std::fabs(oldSafety)+epsilonLen;
228 if( std::fabs( newValue - oldSafety) > reportThreshold * oldSafetyPlus )
229 {
231 std::ostringstream msg;
232 G4double diff= (newValue-oldSafety);
233 G4double relativeDiff= diff / oldSafetyPlus;
234
235 msg << " New (G4SafetyCalculator) value *disagrees* by relative diff " << relativeDiff
236 << " in physical volume '" << motherPhysical->GetName() << "' "
237 << "copy-no = " << motherPhysical->GetCopyNo();
238 if( enteredDaughterVol ) { msg << " ( Just Entered new daughter volume. ) "; }
239 if( exitedMotherVol ) { msg << " ( Just Exited previous volume. ) "; }
240 msg << G4endl;
241 msg << " Safeties: old= " << std::setprecision(12) << oldSafety
242 << " trial " << newValue
243 << " new-old= " << std::setprecision(7) << diff << G4endl;
244
245 if( std::fabs(diff) < errorThreshold * ( std::fabs(oldSafety)+1.0e-20 ) )
246 {
247 msg << " (tiny difference) ";
248 severity= JustWarning;
249 }
250 else
251 {
252 msg << " (real difference) ";
253 severity= FatalException;
254
255 // Extra information -- for big errors
256 msg << " NOTE: keepState = " << keepState << G4endl;
257 msg << " Location - Global coordinates: " << globalPoint
258 << " volume= '" << motherPhysical->GetName() << "'"
259 << " copy-no= " << motherPhysical->GetCopyNo() << G4endl;
260 msg << " Argument maxLength= " << maxLength << G4endl;
261
262 std::size_t depth= fNavHistory.GetDepth();
263 msg << " Navigation History: depth = " << depth << G4endl;
264 for( G4int i=1; i<(G4int)depth; ++i )
265 {
266 msg << " d= " << i << " " << std::setw(32)
267 << fNavHistory.GetVolume(i)->GetName()
268 << " copyNo= " << fNavHistory.GetReplicaNo(i);
269 msg << G4endl;
270 }
271 }
272
273#ifdef G4DEBUG_NAVIGATION
274 G4double redo= SafetyInCurrentVolume(globalPoint, motherPhysical,
275 maxLength, true);
276 msg << " Redoing estimator: value = " << std::setprecision(16) << redo
277 << " diff/last= " << std::setprecision(7) << redo - newValue
278 << " diff/old= " << redo - oldSafety << G4endl;
279#endif
280
281 G4Exception("G4SafetyCalculator::CompareSafetyValues()", "GeomNav1007",
282 severity, msg);
283 }
284}
G4ExceptionSeverity
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4double SafetyInCurrentVolume(const G4ThreeVector &globalpoint, G4VPhysicalVolume *physicalVolume, const G4double pProposedMaxLength=DBL_MAX, G4bool verbose=false)
virtual G4int GetCopyNo() const =0
const G4String & GetName() const

◆ ComputeLocalAxis()

G4ThreeVector G4SafetyCalculator::ComputeLocalAxis ( const G4ThreeVector & pVec) const
inlineprotected

Computes the local direction of the specified vector in the reference system of the volume that was found by LocateGlobalPointAndSetup().

Parameters
[in]pVecVector in global coordinates.
Returns
The local direction of the specified vector in global coordinates.

◆ ComputeLocalPoint()

G4ThreeVector G4SafetyCalculator::ComputeLocalPoint ( const G4ThreeVector & rGlobPoint) const
inlineprotected

Computes point in local coordinates system, given a position vector in world coordinate system.

Parameters
[in]rGlobPointPoint in global coordinates.
Returns
The point in local coordinates system.

Referenced by SafetyInCurrentVolume().

◆ GetDaughtersRegularStructureId()

G4int G4SafetyCalculator::GetDaughtersRegularStructureId ( const G4LogicalVolume * pLv) const
inlineprotected

Gets regular structure ID of first daughter.

Referenced by QuickLocateWithinVolume(), and SafetyInCurrentVolume().

◆ GetExternalNavigation()

G4VExternalNavigation * G4SafetyCalculator::GetExternalNavigation ( ) const

Accessor & modifier for custom external navigation.

Definition at line 194 of file G4SafetyCalculator.cc.

195{
196 return fpExternalNav;
197}

◆ operator=()

G4SafetyCalculator & G4SafetyCalculator::operator= ( const G4SafetyCalculator & )
delete

◆ QuickLocateWithinVolume()

void G4SafetyCalculator::QuickLocateWithinVolume ( const G4ThreeVector & pointLocal,
G4VPhysicalVolume * motherPhysical )
protected

Prepare state of sub-navigators by informing them of current point.

Parameters
[in]pointLocalPoint in local coordinates.
[in]motherPhysicalPointer to current volume where to relocate in case an external custom navigator is used.

Definition at line 152 of file G4SafetyCalculator.cc.

155{
156 // For the case of Voxel (or Parameterised) volume the respective
157 // sub-navigator must be messaged to update its voxel information etc
158
159 // Update the state of the Sub Navigators
160 // - in particular any voxel information they store/cache
161 //
162 G4LogicalVolume* motherLogical = motherPhysical->GetLogicalVolume();
163 G4SmartVoxelHeader* pVoxelHeader = motherLogical->GetVoxelHeader();
164
165 switch( CharacteriseDaughters(motherLogical) )
166 {
167 case kNormal:
168 if ( pVoxelHeader != nullptr )
169 {
170 fvoxelNav.VoxelLocate( pVoxelHeader, pointLocal );
171 }
172 break;
173 case kParameterised:
174 if( GetDaughtersRegularStructureId(motherLogical) != 1 )
175 {
176 // Resets state & returns voxel node
177 //
178 fparamNav.ParamVoxelLocate( pVoxelHeader, pointLocal );
179 }
180 break;
181 case kReplica:
182 // Nothing to do
183 break;
184 case kExternal:
185 fpExternalNav->RelocateWithinVolume( motherPhysical,
186 pointLocal );
187 break;
188 }
189}
G4SmartVoxelHeader * GetVoxelHeader() const
G4int GetDaughtersRegularStructureId(const G4LogicalVolume *pLv) const
EVolume CharacteriseDaughters(const G4LogicalVolume *pLog) const
G4LogicalVolume * GetLogicalVolume() const
@ kNormal
Definition geomdefs.hh:84
@ kParameterised
Definition geomdefs.hh:86
@ kExternal
Definition geomdefs.hh:87
@ kReplica
Definition geomdefs.hh:85

Referenced by SafetyInCurrentVolume().

◆ SafetyInCurrentVolume()

G4double G4SafetyCalculator::SafetyInCurrentVolume ( const G4ThreeVector & globalpoint,
G4VPhysicalVolume * physicalVolume,
const G4double pProposedMaxLength = DBL_MAX,
G4bool verbose = false )

Calculates the isotropic distance to the nearest boundary from the specified point in the global coordinate system.

Parameters
[in]globalPointThe point in global coordinates; it must be located exactly within the current volume (it also must not be in a daughter volume.
[in]physicalVolumeCurrent volume.
[in]pProposedMaxLengthThe calculation will not look beyond the proposed maximum length to avoid extra volume safety calculations.
[in]verboseFlag to enable verbosity (default is false).
Returns
An underestimate of the safety distance (and typically will be if complex volumes are involved.

Definition at line 53 of file G4SafetyCalculator.cc.

58{
59 G4double safety = 0.0;
60 G4ThreeVector stepEndPoint = fNavigator.GetLastStepEndPoint();
61
62 G4ThreeVector localPoint = ComputeLocalPoint(pGlobalpoint);
63
64 G4double distEndpointSq = (pGlobalpoint-stepEndPoint).mag2();
65 G4bool stayedOnEndpoint = distEndpointSq < sqr(fkCarTolerance);
66 G4bool endpointOnSurface = fNavigator.EnteredDaughterVolume()
67 || fNavigator.ExitedMotherVolume();
68
69 G4VPhysicalVolume* motherPhysical = fNavHistory.GetTopVolume();
70 if( motherPhysical != physicalVolume )
71 {
72 std::ostringstream msg;
73 msg << " Current (navigation) phys-volume: " << motherPhysical
74 << " name= " << motherPhysical->GetName() << G4endl
75 << " Request made for phys-volume: " << physicalVolume
76 << " name= " << physicalVolume->GetName() << G4endl;
77 G4Exception("G4SafetyCalculator::SafetyInCurrentVolume", "GeomNav0001",
78 FatalException, msg,
79 "This method must be called only in the Current volume.");
80 }
81
82 if( !(endpointOnSurface && stayedOnEndpoint) )
83 {
84 G4LogicalVolume* motherLogical = motherPhysical->GetLogicalVolume();
85 G4SmartVoxelHeader* pVoxelHeader = motherLogical->GetVoxelHeader();
86
87 // Pseudo-relocate to this point (updates voxel information only)
88 //
89 QuickLocateWithinVolume( localPoint, motherPhysical );
90 //*********************
91
92 // switch(CharacteriseDaughters(motherLogical))
93 auto dtype= CharacteriseDaughters(motherLogical);
94 switch(dtype)
95 {
96 case kNormal:
97 if ( pVoxelHeader != nullptr )
98 {
99 // New way: best safety
100 safety = fVoxelSafety.ComputeSafety(localPoint,
101 *motherPhysical, pMaxLength);
102 }
103 else
104 {
105 safety=fnormalNav.ComputeSafety(localPoint,fNavHistory,pMaxLength);
106 }
107 break;
108 case kParameterised:
109 if( GetDaughtersRegularStructureId(motherLogical) != 1 )
110 {
111 safety=fparamNav.ComputeSafety(localPoint,fNavHistory,pMaxLength);
112 }
113 else // Regular structure
114 {
115 safety=fregularNav.ComputeSafety(localPoint,fNavHistory,pMaxLength);
116 }
117 break;
118 case kReplica:
119 safety = freplicaNav.ComputeSafety(pGlobalpoint, localPoint,
120 fNavHistory, pMaxLength);
121 break;
122 case kExternal:
123 safety = fpExternalNav->ComputeSafety(localPoint, fNavHistory,
124 pMaxLength);
125 break;
126 }
127
128 // Remember last safety origin & value
129 //
130 fPreviousSftOrigin = pGlobalpoint;
131 fPreviousSafety = safety;
132 }
133
134 return safety;
135}
CLHEP::Hep3Vector G4ThreeVector
bool G4bool
Definition G4Types.hh:86
void QuickLocateWithinVolume(const G4ThreeVector &pointLocal, G4VPhysicalVolume *motherPhysical)
G4ThreeVector ComputeLocalPoint(const G4ThreeVector &rGlobPoint) const
T sqr(const T &x)
Definition templates.hh:128

Referenced by CompareSafetyValues().

◆ SetExternalNavigation()

void G4SafetyCalculator::SetExternalNavigation ( G4VExternalNavigation * externalNav)

Definition at line 202 of file G4SafetyCalculator.cc.

203{
204 fpExternalNav = eNav;
205}

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