Belle II Software development
G4LongLivedNeutralTransportation.h
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9// modified from GEANT4 G4Transportation class
10
11
12#ifndef G4LongLivedNeutralTransportation_hh
13#define G4LongLivedNeutralTransportation_hh 1
14
15#include "G4VProcess.hh"
16#include "G4FieldManager.hh"
17
18#include "G4Navigator.hh"
19#include "G4TransportationManager.hh"
20#include "G4PropagatorInField.hh"
21#include "G4Track.hh"
22#include "G4Step.hh"
23#include "G4ParticleChangeForTransport.hh"
24
25class G4SafetyHelper;
26class G4CoupledTransportation;
27namespace Belle2 {
36 class G4LongLivedNeutralTransportation : public G4VProcess {
37
38 public:
39
44 G4LongLivedNeutralTransportation(G4int verbosityLevel = 1);
45
50
62 const G4Track& track,
63 G4double previousStepSize,
64 G4double currentMinimumStep,
65 G4double& currentSafety,
66 G4GPILSelection* selection
67 );
68
75 G4VParticleChange* AlongStepDoIt(
76 const G4Track& track,
77 const G4Step& stepData
78 );
79
86 G4VParticleChange* PostStepDoIt(
87 const G4Track& track,
88 const G4Step& stepData
89 );
90
102 const G4Track& track,
103 G4double previousStepSize,
104 G4ForceCondition* pForceCond
105 );
106
107
108
109 inline G4bool FieldExertedForce() { return fFieldExertedForce; }
112 G4PropagatorInField* GetPropagatorInField();
114 void SetPropagatorInField(G4PropagatorInField* pFieldPropagator);
117 inline G4double GetThresholdWarningEnergy() const;
118 inline G4double GetThresholdImportantEnergy() const;
119 inline G4int GetThresholdTrials() const;
123 inline void SetThresholdWarningEnergy(G4double newEnWarn)
124 {
125 fThreshold_Warning_Energy = newEnWarn;
126 }
127 inline void SetThresholdImportantEnergy(G4double newEnImp)
128 {
130 }
131 inline void SetThresholdTrials(G4int newMaxTrials)
132 {
133 fThresholdTrials = newMaxTrials;
134 }
150 inline G4double GetMaxEnergyKilled() const;
152 inline G4double GetSumEnergyKilled() const;
154 inline void ResetKilledStatistics(G4int report = 1);
157 inline void EnableShortStepOptimisation(G4bool optimise = true);
160 static void SetSilenceLooperWarnings(G4bool val);
162 static G4bool GetSilenceLooperWarnings();
165 public: // without description
166
170 G4double AtRestGetPhysicalInteractionLength(const G4Track&,
171 G4ForceCondition*)
172 { return -1.0; }
173
177 G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&)
178 { return 0; }
179
180 void StartTracking(G4Track* aTrack);
188 virtual void ProcessDescription(std::ostream& outStream) const;
189
190 void PrintStatistics(std::ostream& outStr) const;
193 protected:
194
198 inline G4bool DoesGlobalFieldExist()
199 {
200 G4TransportationManager* transportMgr;
201 transportMgr = G4TransportationManager::GetTransportationManager();
202
203 // fFieldExists= transportMgr->GetFieldManager()->DoesFieldExist();
204 // return fFieldExists;
205 return transportMgr->GetFieldManager()->DoesFieldExist();
206 }
207
208 private:
209
210 G4Navigator* fLinearNavigator;
211 G4PropagatorInField* fFieldPropagator;
213 G4ThreeVector fTransportEndPosition = G4ThreeVector(0.0, 0.0, 0.0);
214 G4ThreeVector fTransportEndMomentumDir = G4ThreeVector(0.0, 0.0, 0.0);
216 G4ThreeVector fTransportEndSpin = G4ThreeVector(0.0, 0.0, 0.0);
217 G4bool fMomentumChanged = true;
218 G4bool fEndGlobalTimeComputed = false;
219 G4double fCandidateEndGlobalTime = 0.0;
221 G4bool fAnyFieldExists = false;
223 G4bool fNewTrack = true;
224 G4bool fFirstStepInVolume = true;
225 G4bool fLastStepInVolume = false;
226 G4bool fGeometryLimitedStep = true;
228 G4bool fFieldExertedForce = false;
230 G4TouchableHandle fCurrentTouchableHandle;
232 G4ThreeVector fPreviousSftOrigin;
235 G4ParticleChangeForTransport fParticleChange;
239 // Thresholds for looping particles:
240 //
241 G4double fThreshold_Warning_Energy = 1.0 * CLHEP::keV;
242 G4double fThreshold_Important_Energy = 1.0 * CLHEP::MeV;
243 G4int fThresholdTrials = 10;
244 // Above 'important' energy a 'looping' particle in field will
245 // *NOT* be abandoned, except after fThresholdTrials attempts.
247 // unstable loopers ( 0 = never )
248 // Counter for steps in which particle reports 'looping',
249 // ( Used if it is above 'Important' Energy. )
250 G4int fNoLooperTrials = 0;
252 // Statistics for tracks abandoned due to looping, only used as flag, no looping neutral particles
253 //
254 G4double fSumEnergyKilled = 0.0;
260
261 G4SafetyHelper* fpSafetyHelper;
263 private:
264
265 friend class G4CoupledTransportation;
266 static G4bool fUseMagneticMoment;
267 static G4bool fUseGravity;
270 };
272}
273
274#endif
Concrete class that does the geometrical transport.
G4TouchableHandle fCurrentTouchableHandle
Current touchable handle.
G4int GetThresholdTrials() const
Access fThresholdTrials.
G4double GetThresholdImportantEnergy() const
Access fThreshold_Important_Energy.
void PrintStatistics(std::ostream &outStr) const
returns current logging info of the algorithm
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
G4VProcess::PostStepDoIt() implementation,.
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
No operation in AtRestGPIL.
G4double fCandidateEndGlobalTime
The particle's state after this Step, Store for DoIt.
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
G4VProcess::AlongStepGetPhysicalInteractionLength() implementation,.
G4int fThresholdTrials
Number of trials an important looper survives.
G4PropagatorInField * GetPropagatorInField()
Access fFieldPropagator, the assistant class that Propagate in a Field.
G4double GetMaxEnergyKilled() const
Print values of looper thresholds.
void SetThresholdWarningEnergy(G4double newEnWarn)
Set fThreshold_Warning_Energy.
void SetThresholdImportantEnergy(G4double newEnImp)
Set fThreshold_Important_Energy.
virtual void ProcessDescription(std::ostream &outStream) const
G4LongLivedNeutralTransportation::ProcessDescription()
G4bool fGeometryLimitedStep
Flag to determine whether a boundary was reached.
G4ThreeVector fTransportEndMomentumDir
The particle's state after this Step, Store for DoIt.
void SetThresholdTrials(G4int newMaxTrials)
Set fThresholdTrials.
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
G4VProcess::AlongStepDoIt() implementation,.
G4int fAbandonUnstableTrials
Number of trials after which to abandon.
G4SafetyHelper * fpSafetyHelper
To pass it the safety value obtained.
void SetPropagatorInField(G4PropagatorInField *pFieldPropagator)
Access/set the assistant class that Propagate in a Field.
void SetHighLooperThresholds()
Get/Set parameters for killing loopers: Above 'important' energy a 'looping' particle in field will N...
G4double fSumEnergyKilled
Sum of abandoned looping tracks energies.
G4ThreeVector fPreviousSftOrigin
Remember last safety origin.
static G4bool fUseGravity
Flag take into account gravity.
void ReportLooperThresholds()
Set low thresholds - for low-E applications.
static void SetSilenceLooperWarnings(G4bool val)
Do not warn about 'looping' particles.
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *pForceCond)
G4VProcess::PostStepGetPhysicalInteractionLength() implementation.
G4ThreeVector fTransportEndSpin
The particle's state after this Step, Store for DoIt.
G4Navigator * fLinearNavigator
Propagator used to transport the particle.
static G4bool fSilenceLooperWarnings
Flag to Supress all 'looper' warnings.
G4double fThreshold_Important_Energy
Give a few trial above this E for looping particle.
G4ThreeVector fTransportEndPosition
The particle's state after this Step, Store for DoIt.
G4bool fEndGlobalTimeComputed
The particle's state after this Step, Store for DoIt.
G4bool fShortStepOptimisation
Whether to avoid calling G4Navigator for short step ( < safety) If using it, the safety estimate for ...
G4int fNoLooperTrials
Counter for steps in which particle reports 'looping', if it is above 'Important' Energy.
G4bool fMomentumChanged
The particle's state after this Step, Store for DoIt.
void StartTracking(G4Track *aTrack)
Reset state for new (potentially resumed) track.
static G4bool fUseMagneticMoment
Flag take into account magnetic moment.
G4ParticleChangeForTransport fParticleChange
New ParticleChange.
void SetLowLooperThresholds()
Shortcut method - old values (meant for HEP)
G4double fThreshold_Warning_Energy
Warn above this energy about looping particle.
void EnableShortStepOptimisation(G4bool optimise=true)
Whether short steps < safety will avoid to call Navigator (if field=0)
G4bool FieldExertedForce()
References fFieldExertedForce.
G4bool DoesGlobalFieldExist()
Checks whether a field exists for the "global" field manager.
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
No operation in AtRestDoIt.
static G4bool GetSilenceLooperWarnings()
Do not throw exception about 'looping' particles.
G4double GetSumEnergyKilled() const
Access fSumEnergyKilled.
void ResetKilledStatistics(G4int report=1)
Statistics for tracks killed (currently due to looping in field)
G4double GetThresholdWarningEnergy() const
Access fThreshold_Warning_Energy.
G4PropagatorInField * fFieldPropagator
Propagator used to transport the particle.
G4double fTransportEndKineticEnergy
The particle's state after this Step, Store for DoIt.
Abstract base class for different kinds of events.