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

#include <G4ITSafetyHelper.hh>

+ Inheritance diagram for G4ITSafetyHelper:

Classes

class  State
 

Public Member Functions

 G4ITSafetyHelper ()
 
 ~G4ITSafetyHelper ()
 
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 ()
 
- Public Member Functions inherited from G4TrackStateDependent< G4ITSafetyHelper >
virtual ~G4TrackStateDependent ()
 
virtual void SetTrackState (G4shared_ptr< StateType > state)
 
virtual G4VTrackStateHandle PopTrackState ()
 
virtual G4VTrackStateHandle GetTrackState () const
 
virtual StateTypeHandle GetConcreteTrackState () const
 
virtual void LoadTrackState (G4TrackStateManager &manager)
 
virtual void SaveTrackState (G4TrackStateManager &manager)
 
virtual void NewTrackState ()
 
virtual StateTypeHandle CreateTrackState () const
 
virtual void ResetTrackState ()
 
- Public Member Functions inherited from G4VTrackStateDependent
 G4VTrackStateDependent ()
 
virtual ~G4VTrackStateDependent ()
 
virtual void NewTrackState ()=0
 
virtual void LoadTrackState (G4TrackStateManager &)=0
 
virtual void SaveTrackState (G4TrackStateManager &)=0
 
virtual G4VTrackStateHandle GetTrackState () const =0
 
virtual G4VTrackStateHandle PopTrackState ()=0
 
virtual void ResetTrackState ()=0
 

Additional Inherited Members

- Public Types inherited from G4TrackStateDependent< G4ITSafetyHelper >
typedef G4ITSafetyHelper ClassType
 
typedef G4TrackState< G4ITSafetyHelperStateType
 
typedef G4shared_ptr< StateTypeStateTypeHandle
 
- Protected Member Functions inherited from G4TrackStateDependent< G4ITSafetyHelper >
 G4TrackStateDependent ()
 
- Protected Attributes inherited from G4TrackStateDependent< G4ITSafetyHelper >
StateTypeHandle fpTrackState
 

Detailed Description

Definition at line 54 of file G4ITSafetyHelper.hh.

Constructor & Destructor Documentation

◆ G4ITSafetyHelper()

G4ITSafetyHelper::G4ITSafetyHelper ( )

Definition at line 41 of file G4ITSafetyHelper.cc.

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

◆ ~G4ITSafetyHelper()

G4ITSafetyHelper::~G4ITSafetyHelper ( )

Definition at line 89 of file G4ITSafetyHelper.cc.

90{
91}

Member Function Documentation

◆ CheckNextStep()

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

Definition at line 93 of file G4ITSafetyHelper.cc.

97{
98 // Distance in the Mass geometry
99 //
100 G4double linstep = fpMassNavigator->CheckNextStep(position, direction,
101 currentMaxStep, newSafety);
102 fpTrackState->fLastSafetyPosition = position;
103 fpTrackState->fLastSafety = newSafety;
104
105 // TODO: Can replace this with a call to PathFinder
106 // giving id of Mass Geometry --> this avoid doing the work twice
107
108 return linstep;
109}
double G4double
Definition: G4Types.hh:83

Referenced by G4DNABrownianTransportation::ComputeGeomLimit().

◆ ComputeSafety()

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

Definition at line 111 of file G4ITSafetyHelper.cc.

113{
114 G4double newSafety;
115
116 // Only recompute (calling Navigator/PathFinder) if 'position'
117 // is *not* the safety location and has moved 'significantly'
118 //
119 G4double moveLengthSq = (position - fpTrackState->fLastSafetyPosition).mag2();
120 if ((moveLengthSq > 0.0))
121 {
122 if (!fUseParallelGeometries)
123 {
124 // Safety for mass geometry
125 newSafety = fpMassNavigator->ComputeSafety(position, maxLength, true);
126 }
127 else
128 {
129 // Safety for all geometries
130 newSafety = fpPathFinder->ComputeSafety(position);
131 }
132
133 // We can only store a 'true' safety - one that was not restricted by maxLength
134 if (newSafety < maxLength)
135 {
136 fpTrackState->fLastSafety = newSafety;
137 fpTrackState->fLastSafetyPosition = position;
138 }
139 }
140 else
141 {
142 // return last value if position is not (significantly) changed
143 //
144 // G4double moveLength = 0;
145 // if( moveLengthSq > 0.0 ) { moveLength= std::sqrt(moveLengthSq); }
146 newSafety = fpTrackState->fLastSafety; // -moveLength;
147 }
148 return newSafety;
149}
G4double ComputeSafety(const G4ThreeVector &globalPoint)

◆ EnableParallelNavigation()

void G4ITSafetyHelper::EnableParallelNavigation ( G4bool  parallel)
inline

◆ GetWorldVolume()

G4VPhysicalVolume * G4ITSafetyHelper::GetWorldVolume ( )
inline

Definition at line 154 of file G4ITSafetyHelper.hh.

155{
156 return fpMassNavigator->GetWorldVolume();
157}

Referenced by G4DNABrownianTransportation::ComputeGeomLimit().

◆ InitialiseHelper()

void G4ITSafetyHelper::InitialiseHelper ( )

Definition at line 79 of file G4ITSafetyHelper.cc.

80{
82 if (fFirstCall)
83 {
85 }
86 fFirstCall = false;
87}

Referenced by G4ITTransportation::BuildPhysicsTable(), G4DNABrownianTransportation::BuildPhysicsTable(), and G4ITPathFinder::PrepareNewTrack().

◆ InitialiseNavigator()

void G4ITSafetyHelper::InitialiseNavigator ( )

Definition at line 55 of file G4ITSafetyHelper.cc.

56{
57 fpPathFinder = G4PathFinder::GetInstance();
58
59 G4ITTransportationManager* pTransportMgr =
61
62 fpMassNavigator = pTransportMgr->GetNavigatorForTracking();
63
64 if(fpMassNavigator == 0) abort();
65
66 // Check
67 //
68 G4VPhysicalVolume* worldPV = fpMassNavigator->GetWorldVolume();
69 if (worldPV == 0)
70 {
71 G4Exception("G4ITSafetyHelper::InitialiseNavigator",
72 "InvalidNavigatorWorld", FatalException,
73 "Found that existing tracking Navigator has NULL world");
74 }
75
76 // fMassNavigatorId = pTransportMgr->ActivateNavigator( fpMassNavigator );
77}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
static G4ITTransportationManager * GetTransportationManager()
G4ITNavigator * GetNavigatorForTracking() const
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:52

Referenced by InitialiseHelper().

◆ Locate()

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

Definition at line 183 of file G4ITSafetyHelper.cc.

185{
186 if (!fUseParallelGeometries)
187 {
188 fpMassNavigator->LocateGlobalPointAndSetup(newPosition, &newDirection, true,
189 false);
190 }
191 else
192 {
193 fpPathFinder->Locate(newPosition, newDirection);
194 }
195}
void Locate(const G4ThreeVector &position, const G4ThreeVector &direction, G4bool relativeSearch=true)

◆ ReLocateWithinVolume()

void G4ITSafetyHelper::ReLocateWithinVolume ( const G4ThreeVector pGlobalPoint)

Definition at line 151 of file G4ITSafetyHelper.cc.

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

◆ SetCurrentSafety()

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

Definition at line 160 of file G4ITSafetyHelper.hh.

161{
162 fpTrackState->fLastSafety = val;
163 fpTrackState->fLastSafetyPosition = pos;
164}

Referenced by G4ITTransportation::AlongStepGetPhysicalInteractionLength().

◆ SetVerboseLevel()

G4int G4ITSafetyHelper::SetVerboseLevel ( G4int  lev)
inline

Definition at line 97 of file G4ITSafetyHelper.hh.

98 { G4int oldlv= fVerbose; fVerbose= lev; return oldlv;}
int G4int
Definition: G4Types.hh:85

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