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

#include <G4SafetyHelper.hh>

Public Member Functions

 G4SafetyHelper ()
 
 ~G4SafetyHelper ()
 
G4double CheckNextStep (const G4ThreeVector &position, const G4ThreeVector &direction, const G4double currentMaxStep, G4double &newSafety)
 
G4double ComputeSafety (const G4ThreeVector &pGlobalPoint, G4double maxRadius=DBL_MAX)
 
void Locate (const G4ThreeVector &pGlobalPoint, const G4ThreeVector &direction)
 
void ReLocateWithinVolume (const G4ThreeVector &pGlobalPoint)
 
G4bool RecheckDistanceToCurrentBoundary (const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double *prDistance, G4double *prNewSafety=nullptr) const
 
void EnableParallelNavigation (G4bool parallel)
 
void InitialiseNavigator ()
 
G4int SetVerboseLevel (G4int lev)
 
G4VPhysicalVolumeGetWorldVolume ()
 
void SetCurrentSafety (G4double val, const G4ThreeVector &pos)
 
void InitialiseHelper ()
 

Detailed Description

Definition at line 46 of file G4SafetyHelper.hh.

Constructor & Destructor Documentation

◆ G4SafetyHelper()

G4SafetyHelper::G4SafetyHelper ( )

Definition at line 38 of file G4SafetyHelper.cc.

39 : fLastSafetyPosition(0.0,0.0,0.0)
40{
41}

◆ ~G4SafetyHelper()

G4SafetyHelper::~G4SafetyHelper ( )

Definition at line 73 of file G4SafetyHelper.cc.

74{
75}

Member Function Documentation

◆ CheckNextStep()

G4double G4SafetyHelper::CheckNextStep ( const G4ThreeVector position,
const G4ThreeVector direction,
const G4double  currentMaxStep,
G4double newSafety 
)

Definition at line 78 of file G4SafetyHelper.cc.

82{
83 // Distance in the Mass geometry
84 //
85 G4double linstep = fpMassNavigator->CheckNextStep( position,
86 direction,
87 currentMaxStep,
88 newSafety);
89 fLastSafetyPosition = position;
90 fLastSafety = newSafety;
91
92 // TO-DO: Can replace this with a call to PathFinder
93 // giving id of Mass Geometry --> this avoid doing the work twice
94
95 return linstep;
96}
double G4double
Definition: G4Types.hh:83
G4double CheckNextStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
#define position
Definition: xmlparse.cc:622

Referenced by G4VMscModel::ComputeGeomLimit().

◆ ComputeSafety()

G4double G4SafetyHelper::ComputeSafety ( const G4ThreeVector pGlobalPoint,
G4double  maxRadius = DBL_MAX 
)

Definition at line 98 of file G4SafetyHelper.cc.

100{
101 G4double newSafety;
102
103 // Only recompute (calling Navigator/PathFinder) if 'position'
104 // is *not* the safety location and has moved 'significantly'
105 //
106 G4double moveLengthSq = (position-fLastSafetyPosition).mag2();
107 if( (moveLengthSq > 0.0 ) )
108 {
109 if( !fUseParallelGeometries )
110 {
111 // Safety for mass geometry
112 newSafety = fpMassNavigator->ComputeSafety(position, maxLength, true);
113
114 // Only store a 'true' safety - one that was not restricted by maxLength
115 if( newSafety < maxLength )
116 {
117 fLastSafety= newSafety;
118 fLastSafetyPosition = position;
119 }
120 }
121 else
122 {
123 // Safety for all geometries
124 newSafety = fpPathFinder->ComputeSafety(position);
125
126 fLastSafety= newSafety;
127 fLastSafetyPosition = position;
128 }
129
130 }
131 else
132 {
133 // return last value if position is not (significantly) changed
134 //
135 // G4double moveLength = 0;
136 // if( moveLengthSq > 0.0 ) { moveLength= std::sqrt(moveLengthSq); }
137 newSafety = fLastSafety; // -moveLength;
138 }
139
140 return newSafety;
141}
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
G4double ComputeSafety(const G4ThreeVector &globalPoint)

Referenced by G4VEnergyLossProcess::AlongStepDoIt(), G4VMultipleScattering::AlongStepDoIt(), and G4VMscModel::ComputeSafety().

◆ EnableParallelNavigation()

void G4SafetyHelper::EnableParallelNavigation ( G4bool  parallel)
inline

Definition at line 143 of file G4SafetyHelper.hh.

144{
145 fUseParallelGeometries = parallel;
146}

Referenced by G4PathFinder::EnableParallelNavigation().

◆ GetWorldVolume()

G4VPhysicalVolume * G4SafetyHelper::GetWorldVolume ( )
inline

Definition at line 149 of file G4SafetyHelper.hh.

150{
151 return fpMassNavigator->GetWorldVolume();
152}
G4VPhysicalVolume * GetWorldVolume() const

◆ InitialiseHelper()

void G4SafetyHelper::InitialiseHelper ( )

◆ InitialiseNavigator()

void G4SafetyHelper::InitialiseNavigator ( )

Definition at line 43 of file G4SafetyHelper.cc.

44{
45 fpPathFinder = G4PathFinder::GetInstance();
46
47 G4TransportationManager* pTransportMgr=
49
50 fpMassNavigator = pTransportMgr->GetNavigatorForTracking();
51
52 // Check
53 //
54 G4VPhysicalVolume* worldPV = fpMassNavigator->GetWorldVolume();
55 if( worldPV == nullptr )
56 {
57 G4Exception("G4SafetyHelper::InitialiseNavigator",
58 "GeomNav0003", FatalException,
59 "Found that existing tracking Navigator has NULL world");
60 }
61
62 fMassNavigatorId = pTransportMgr->ActivateNavigator( fpMassNavigator );
63}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:52
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
G4int ActivateNavigator(G4Navigator *aNavigator)

Referenced by InitialiseHelper().

◆ Locate()

void G4SafetyHelper::Locate ( const G4ThreeVector pGlobalPoint,
const G4ThreeVector direction 
)

Definition at line 177 of file G4SafetyHelper.cc.

179{
180 if( !fUseParallelGeometries )
181 {
182 fpMassNavigator->LocateGlobalPointAndSetup(newPosition, &newDirection,
183 true, false);
184 }
185 else
186 {
187 fpPathFinder->Locate( newPosition, newDirection );
188 }
189}
virtual G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=nullptr, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
Definition: G4Navigator.cc:126
void Locate(const G4ThreeVector &position, const G4ThreeVector &direction, G4bool relativeSearch=true)

◆ RecheckDistanceToCurrentBoundary()

G4bool G4SafetyHelper::RecheckDistanceToCurrentBoundary ( const G4ThreeVector pGlobalPoint,
const G4ThreeVector pDirection,
const G4double  pCurrentProposedStepLength,
G4double prDistance,
G4double prNewSafety = nullptr 
) const

Definition at line 191 of file G4SafetyHelper.cc.

197{
198 G4bool retval;
199 if( !fUseParallelGeometries )
200 {
201 retval = fpMassNavigator->RecheckDistanceToCurrentBoundary(pGlobalPoint,
202 pDirection,
203 aProposedMove,
204 prDistance,
205 prNewSafety);
206 }
207 else
208 {
209 retval = fpPathFinder->RecheckDistanceToCurrentBoundary(pGlobalPoint,
210 pDirection,
211 aProposedMove,
212 prDistance,
213 prNewSafety);
214 }
215 return retval;
216}
bool G4bool
Definition: G4Types.hh:86
virtual G4bool RecheckDistanceToCurrentBoundary(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double CurrentProposedStepLength, G4double *prDistance, G4double *prNewSafety=nullptr) const
G4bool RecheckDistanceToCurrentBoundary(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double *prDistance, G4double *prNewSafety=nullptr) const

◆ ReLocateWithinVolume()

void G4SafetyHelper::ReLocateWithinVolume ( const G4ThreeVector pGlobalPoint)

Definition at line 143 of file G4SafetyHelper.cc.

144{
145#ifdef G4VERBOSE
146 if( fVerbose > 0 )
147 {
148 // There is an opportunity - and need - to check whether
149 // the proposed move is safe
150 G4ThreeVector moveVec = newPosition - fLastSafetyPosition;
151 if( moveVec.mag2() > sqr(fLastSafety) )
152 {
153 // A problem exists - we are proposing to move outside 'Safety Sphere'
155 ed << "Unsafe Move> Asked to relocate beyond 'Safety sphere'. Details: "
156 << G4endl;
157 ed << " Safety Sphere: Radius = " << fLastSafety;
158 ed << " Center = " << fLastSafetyPosition << G4endl;
159 ed << " New Location : Move = " << moveVec.mag();
160 ed << " Position = " << newPosition << G4endl;
161 G4Exception("G4SafetyHelper::ReLocateWithinVolume",
162 "GeomNav1001", JustWarning, ed);
163 }
164 }
165#endif
166
167 if( !fUseParallelGeometries )
168 {
169 fpMassNavigator->LocateGlobalPointWithinVolume( newPosition );
170 }
171 else
172 {
173 fpPathFinder->ReLocate( newPosition );
174 }
175}
@ JustWarning
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4endl
Definition: G4ios.hh:57
double mag2() const
double mag() const
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
Definition: G4Navigator.cc:608
void ReLocate(const G4ThreeVector &position)
T sqr(const T &x)
Definition: templates.hh:128

Referenced by G4VMultipleScattering::AlongStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), and G4MuNeutrinoNucleusProcess::PostStepDoIt().

◆ SetCurrentSafety()

void G4SafetyHelper::SetCurrentSafety ( G4double  val,
const G4ThreeVector pos 
)
inline

Definition at line 155 of file G4SafetyHelper.hh.

156{
157 fLastSafety = val;
158 fLastSafetyPosition = pos;
159}

Referenced by G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), and G4Transportation::AlongStepGetPhysicalInteractionLength().

◆ SetVerboseLevel()

G4int G4SafetyHelper::SetVerboseLevel ( G4int  lev)
inline

Definition at line 135 of file G4SafetyHelper.hh.

136{
137 G4int oldlv = fVerbose;
138 fVerbose = lev;
139 return oldlv;
140}
int G4int
Definition: G4Types.hh:85

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