Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4SafetyCalculator.hh
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// G4SafetyCalculator
27//
28// Class description:
29//
30// A class that provides an estimate of the isotropic safety - the
31// minimum distance from a global point to the nearest boundary
32// of the current volume or the nearest daughter volumes.
33// This estimate can be an underestimate, either because a solid
34// provides an underestimate (for speed) or in order to avoid
35// substantial additional computations.
36//
37// Obtains from the navigator the current transformation history.
38
39// Author: John Apostolakis, CERN - February 2023
40// --------------------------------------------------------------------
41#ifndef G4SafetyCalculator_HH
42#define G4SafetyCalculator_HH 1
43
44#include "geomdefs.hh"
45
46#include "G4ThreeVector.hh"
47#include "G4AffineTransform.hh"
48#include "G4RotationMatrix.hh"
49
50#include "G4LogicalVolume.hh" // Used in inline methods
52
54#include "G4NormalNavigation.hh"
55#include "G4VoxelNavigation.hh"
60
61#include "G4VoxelSafety.hh"
62
63#include <iostream>
64
66
68{
69 public:
70
71 G4SafetyCalculator( const G4Navigator& navigator,
72 const G4NavigationHistory& navHistory );
73 // Constructor - initialisers and setup.
74
77 // Copy constructor & assignment operator not allowed.
78
80 // Destructor. No actions.
81
83 G4VPhysicalVolume* physicalVolume,
84 const G4double pProposedMaxLength = DBL_MAX,
85 G4bool verbose = false );
86 // Calculate the isotropic distance to the nearest boundary from the
87 // specified point in the global coordinate system.
88 // The globalpoint utilised *must* be located exactly within the
89 // current volume (it also must *not* be in a daughter volume).
90 // The value returned can be an underestimate (and typically will be
91 // if complex volumes are involved).
92 // The calculation will not look beyond the proposed maximum length
93 // to avoid extra volume safety calculations. The geometry must be closed.
94
97 // Accessor & modifier for custom external navigation.
98
99 void CompareSafetyValues( G4double oldSafety,
100 G4double newValue,
101 G4VPhysicalVolume* motherPhysical,
102 const G4ThreeVector &globalPoint,
103 G4bool keepState,
104 G4double maxLength,
105 G4bool enteredVolume,
106 G4bool exitedVolume );
107 // Compare estimates of the safety, and report if difference(s) found.
108
109 protected:
110
111 void QuickLocateWithinVolume(const G4ThreeVector& pointLocal,
112 G4VPhysicalVolume* motherPhysical);
113 // Prepare state of sub-navigators by informing them of current point.
114
115 inline G4ThreeVector ComputeLocalPoint(const G4ThreeVector& rGlobPoint) const;
116 // Return position vector in local coordinate system, given a position
117 // vector in world coordinate system.
118
119 inline G4ThreeVector ComputeLocalAxis(const G4ThreeVector& pVec) const;
120 // Return the local direction of the specified vector in the reference
121 // system of the volume that was found by LocalGlobalPointAndSetup().
122 // The Local Coordinates of point in world coordinate system.
123
124 inline EVolume CharacteriseDaughters(const G4LogicalVolume* pLog) const;
125 // Characterise daughter of logical volume.
126
128 // Get regular structure ID of first daughter.
129
130 private:
131
132 // BEGIN -- Tracking Invariants part 1
133 //
134 const G4Navigator& fNavigator;
135 // Associated navigator. Needed for details of current state,
136 // for optimisation
137
138 const G4NavigationHistory& fNavHistory;
139 // Associated navigator's navigation history. Transformation and history
140 // of the current path through the geometrical hierarchy.
141 //
142 // END -- Tracking Invariants part 1
143
144 G4double fkCarTolerance;
145 // Cached tolerance.
146
147 // BEGIN State information
148 //
149 G4ThreeVector fPreviousSftOrigin;
150 G4double fPreviousSafety = 0.0;
151 // Memory of last safety origin & value. Used in ComputeStep to ensure
152 // that origin of current Step is in the same volume as the point of the
153 // last relocation.
154
155 // Helpers/Utility classes - their state can change
156 //
157 G4NormalNavigation fnormalNav;
158 G4VoxelNavigation fvoxelNav;
160 G4ReplicaNavigation freplicaNav;
161 G4RegularNavigation fregularNav;
162 G4VExternalNavigation* fpExternalNav = nullptr;
163 G4VoxelSafety fVoxelSafety;
164};
165
166// Auxiliary inline methods -- copied from G4Navigator
167
168// Return local coordinates given point in the world coord system.
169//
170inline G4ThreeVector
172{
173 return fNavHistory.GetTopTransform().TransformPoint(pGlobalPoint);
174}
175
176// Returns local direction given vector direction in world coord system.
177//
178inline G4ThreeVector
180{
181 return fNavHistory.GetTopTransform().TransformAxis(pVec);
182}
183
184inline EVolume
189
190inline G4int
192{
193 G4int regId = 0;
194 G4VPhysicalVolume *pVol;
195
196 if ( pLog->GetNoDaughters() == 1 )
197 {
198 pVol = pLog->GetDaughter(0);
199 regId = pVol->GetRegularStructureId();
200 }
201 return regId;
202}
203
204#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
EVolume CharacteriseDaughters() const
std::size_t GetNoDaughters() const
G4VPhysicalVolume * GetDaughter(const std::size_t i) const
const G4AffineTransform & GetTopTransform() const
void SetExternalNavigation(G4VExternalNavigation *externalNav)
G4VExternalNavigation * GetExternalNavigation() const
G4SafetyCalculator & operator=(const G4SafetyCalculator &)=delete
G4SafetyCalculator(const G4Navigator &navigator, const G4NavigationHistory &navHistory)
G4int GetDaughtersRegularStructureId(const G4LogicalVolume *pLv) const
G4ThreeVector ComputeLocalAxis(const G4ThreeVector &pVec) const
void QuickLocateWithinVolume(const G4ThreeVector &pointLocal, G4VPhysicalVolume *motherPhysical)
void CompareSafetyValues(G4double oldSafety, G4double newValue, G4VPhysicalVolume *motherPhysical, const G4ThreeVector &globalPoint, G4bool keepState, G4double maxLength, G4bool enteredVolume, G4bool exitedVolume)
G4ThreeVector ComputeLocalPoint(const G4ThreeVector &rGlobPoint) const
~G4SafetyCalculator()=default
G4double SafetyInCurrentVolume(const G4ThreeVector &globalpoint, G4VPhysicalVolume *physicalVolume, const G4double pProposedMaxLength=DBL_MAX, G4bool verbose=false)
G4SafetyCalculator(const G4SafetyCalculator &)=delete
EVolume CharacteriseDaughters(const G4LogicalVolume *pLog) const
virtual G4int GetRegularStructureId() const =0
EVolume
Definition geomdefs.hh:83
#define DBL_MAX
Definition templates.hh:62