Geant4 9.6.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)
 
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 53 of file G4SafetyHelper.hh.

Constructor & Destructor Documentation

◆ G4SafetyHelper()

G4SafetyHelper::G4SafetyHelper ( )

Definition at line 42 of file G4SafetyHelper.cc.

43 : fUseParallelGeometries(false), // By default, one geometry only
44 fFirstCall(true),
45 fVerbose(0),
46 fLastSafetyPosition(0.0,0.0,0.0),
47 fLastSafety(0.0),
48 fRecomputeFactor(0.0)
49{
50 fpPathFinder= 0; // Cannot initialise this yet - a loop results
51
52 // Initialization of the Navigator pointer is postponed, and must
53 // be undertaken by another class calling InitialiseHelper()
54 //
55 fpMassNavigator= 0;
56 fMassNavigatorId= -1;
57}

◆ ~G4SafetyHelper()

G4SafetyHelper::~G4SafetyHelper ( )

Definition at line 89 of file G4SafetyHelper.cc.

90{
91}

Member Function Documentation

◆ CheckNextStep()

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

Definition at line 94 of file G4SafetyHelper.cc.

98{
99 // Distance in the Mass geometry
100 //
101 G4double linstep = fpMassNavigator->CheckNextStep( position,
102 direction,
103 currentMaxStep,
104 newSafety);
105 fLastSafetyPosition = position;
106 fLastSafety = newSafety;
107
108 // TO-DO: Can replace this with a call to PathFinder
109 // giving id of Mass Geometry --> this avoid doing the work twice
110
111 return linstep;
112}
double G4double
Definition: G4Types.hh:64
G4double CheckNextStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
#define position
Definition: xmlparse.cc:605

Referenced by G4VMscModel::ComputeGeomLimit().

◆ ComputeSafety()

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

Definition at line 114 of file G4SafetyHelper.cc.

115{
116 G4double newSafety;
117
118 // Only recompute (calling Navigator/PathFinder) if 'position'
119 // is *not* the safety location and has moved 'significantly'
120 //
121 G4double moveLengthSq = (position-fLastSafetyPosition).mag2();
122 G4double safeDistance = fRecomputeFactor*fLastSafety;
123 if( (moveLengthSq > 0.0 )
124 && (moveLengthSq >= safeDistance*safeDistance))
125 {
126 fLastSafetyPosition = position;
127
128 if( !fUseParallelGeometries )
129 {
130 // Safety for mass geometry
131 fLastSafety = fpMassNavigator->ComputeSafety(position, maxLength, true);
132 }
133 else
134 {
135 // Safety for all geometries
136 fLastSafety = fpPathFinder->ComputeSafety(position);
137 }
138 newSafety = fLastSafety;
139 }
140 else
141 {
142 // return last value if position is not significantly changed
143 //
144 G4double moveLength = 0;
145 if( moveLengthSq > 0.0 )
146 {
147 moveLength= std::sqrt(moveLengthSq);
148 }
149 newSafety = fLastSafety-moveLength;
150 }
151 return newSafety;
152}
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
G4double ComputeSafety(const G4ThreeVector &globalPoint)

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

◆ EnableParallelNavigation()

void G4SafetyHelper::EnableParallelNavigation ( G4bool  parallel)
inline

Definition at line 130 of file G4SafetyHelper.hh.

131{
132 fUseParallelGeometries = parallel;
133}

Referenced by G4PathFinder::EnableParallelNavigation().

◆ GetWorldVolume()

G4VPhysicalVolume * G4SafetyHelper::GetWorldVolume ( )
inline

Definition at line 136 of file G4SafetyHelper.hh.

137{
138 return fpMassNavigator->GetWorldVolume();
139}
G4VPhysicalVolume * GetWorldVolume() const

Referenced by G4VMscModel::ComputeGeomLimit().

◆ InitialiseHelper()

void G4SafetyHelper::InitialiseHelper ( )

Definition at line 81 of file G4SafetyHelper.cc.

82{
83 fLastSafetyPosition = G4ThreeVector(0.0,0.0,0.0);
84 fLastSafety = 0.0;
85 if (fFirstCall) { InitialiseNavigator(); }
86 fFirstCall = false;
87}
CLHEP::Hep3Vector G4ThreeVector
void InitialiseNavigator()

Referenced by G4VEnergyLossProcess::BuildPhysicsTable(), G4VMscModel::GetParticleChangeForMSC(), G4PathFinder::PrepareNewTrack(), and G4VMultipleScattering::PreparePhysicsTable().

◆ InitialiseNavigator()

void G4SafetyHelper::InitialiseNavigator ( )

Definition at line 59 of file G4SafetyHelper.cc.

60{
61 fpPathFinder= G4PathFinder::GetInstance();
62
63 G4TransportationManager* pTransportMgr=
65
66 fpMassNavigator = pTransportMgr->GetNavigatorForTracking();
67
68 // Check
69 //
70 G4VPhysicalVolume* worldPV = fpMassNavigator->GetWorldVolume();
71 if( worldPV == 0 )
72 {
73 G4Exception("G4SafetyHelper::InitialiseNavigator",
74 "InvalidNavigatorWorld", FatalException,
75 "Found that existing tracking Navigator has NULL world");
76 }
77
78 fMassNavigatorId = pTransportMgr->ActivateNavigator( fpMassNavigator );
79}
@ FatalException
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:57
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
G4int ActivateNavigator(G4Navigator *aNavigator)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Referenced by InitialiseHelper().

◆ Locate()

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

Definition at line 184 of file G4SafetyHelper.cc.

186{
187 if( !fUseParallelGeometries)
188 {
189 fpMassNavigator->LocateGlobalPointAndSetup(newPosition, &newDirection,
190 true, false);
191 }
192 else
193 {
194 fpPathFinder->Locate( newPosition, newDirection );
195 }
196}
virtual G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=0, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
Definition: G4Navigator.cc:116
void Locate(const G4ThreeVector &position, const G4ThreeVector &direction, G4bool relativeSearch=true)

◆ ReLocateWithinVolume()

void G4SafetyHelper::ReLocateWithinVolume ( const G4ThreeVector pGlobalPoint)

Definition at line 154 of file G4SafetyHelper.cc.

155{
156#ifdef G4VERBOSE
157 if( fVerbose > 0 ) {
158 // There is an opportunity - and need - to check whether the proposed move is safe
159 G4ThreeVector moveVec= newPosition - fLastSafetyPosition;
160 if( moveVec.mag2() > sqr(fLastSafety) )
161 {
162 // A problem exists - we are proposing to move outside 'Safety Sphere'
164 ed << " Safety Sphere: Radius = " << fLastSafety;
165 ed << " Center = " << fLastSafetyPosition << G4endl;
166 ed << " New Location : Move = " << moveVec.mag2();
167 ed << " Position = " << newPosition << G4endl;
168 G4Exception("G4SafetyHelper::ReLocateWithinVolume", "GeomNav999", JustWarning,
169 "Unsafe Move> Asked to relocate beyond 'Safety sphere'.");
170 }
171 }
172#endif
173
174 if( !fUseParallelGeometries )
175 {
176 fpMassNavigator->LocateGlobalPointWithinVolume( newPosition );
177 }
178 else
179 {
180 fpPathFinder->ReLocate( newPosition );
181 }
182}
@ JustWarning
#define G4endl
Definition: G4ios.hh:52
double mag2() const
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
Definition: G4Navigator.cc:548
void ReLocate(const G4ThreeVector &position)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
T sqr(const T &x)
Definition: templates.hh:145

Referenced by G4VMultipleScattering::PostStepDoIt().

◆ SetCurrentSafety()

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

◆ SetVerboseLevel()

G4int G4SafetyHelper::SetVerboseLevel ( G4int  lev)
inline

Definition at line 96 of file G4SafetyHelper.hh.

96{ G4int oldlv= fVerbose; fVerbose= lev; return oldlv; }
int G4int
Definition: G4Types.hh:66

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