Belle II Software  release-05-01-25
G4LongLivedNeutralTransportation.h
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2018 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Sascha Dreyer *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 // modified from GEANT4 G4Transportation class
12 
13 
14 #ifndef G4LongLivedNeutralTransportation_hh
15 #define G4LongLivedNeutralTransportation_hh 1
16 
17 #include "G4VProcess.hh"
18 #include "G4FieldManager.hh"
19 
20 #include "G4Navigator.hh"
21 #include "G4TransportationManager.hh"
22 #include "G4PropagatorInField.hh"
23 #include "G4Track.hh"
24 #include "G4Step.hh"
25 #include "G4ParticleChangeForTransport.hh"
26 
27 class G4SafetyHelper;
28 class G4CoupledTransportation;
29 namespace Belle2 {
38  class G4LongLivedNeutralTransportation : public G4VProcess {
39 
40  public:
41 
46  G4LongLivedNeutralTransportation(G4int verbosityLevel = 1);
47 
52 
64  const G4Track& track,
65  G4double previousStepSize,
66  G4double currentMinimumStep,
67  G4double& currentSafety,
68  G4GPILSelection* selection
69  );
70 
77  G4VParticleChange* AlongStepDoIt(
78  const G4Track& track,
79  const G4Step& stepData
80  );
81 
88  G4VParticleChange* PostStepDoIt(
89  const G4Track& track,
90  const G4Step& stepData
91  );
92 
104  const G4Track&,
105  G4double previousStepSize,
106  G4ForceCondition* pForceCond
107  );
108 
109 
110 
111  inline G4bool FieldExertedForce() { return fFieldExertedForce; }
114  G4PropagatorInField* GetPropagatorInField();
116  void SetPropagatorInField(G4PropagatorInField* pFieldPropagator);
119  inline G4double GetThresholdWarningEnergy() const;
120  inline G4double GetThresholdImportantEnergy() const;
121  inline G4int GetThresholdTrials() const;
125  inline void SetThresholdWarningEnergy(G4double newEnWarn)
126  {
127  fThreshold_Warning_Energy = newEnWarn;
128  }
129  inline void SetThresholdImportantEnergy(G4double newEnImp)
130  {
131  fThreshold_Important_Energy = newEnImp;
132  }
133  inline void SetThresholdTrials(G4int newMaxTrials)
134  {
135  fThresholdTrials = newMaxTrials;
136  }
147  void SetLowLooperThresholds();
149  void ReportLooperThresholds();
152  inline G4double GetMaxEnergyKilled() const;
154  inline G4double GetSumEnergyKilled() const;
156  inline void ResetKilledStatistics(G4int report = 1);
159  inline void EnableShortStepOptimisation(G4bool optimise = true);
162  static void SetSilenceLooperWarnings(G4bool val);
164  static G4bool GetSilenceLooperWarnings();
167  public: // without description
168 
172  G4double AtRestGetPhysicalInteractionLength(const G4Track&,
173  G4ForceCondition*)
174  { return -1.0; }
175 
179  G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&)
180  { return 0; }
181 
182  void StartTracking(G4Track* aTrack);
190  virtual void ProcessDescription(std::ostream& outFile) const;
191 
192  void PrintStatistics(std::ostream& outStr) const;
195  protected:
196 
200  inline G4bool DoesGlobalFieldExist()
201  {
202  G4TransportationManager* transportMgr;
203  transportMgr = G4TransportationManager::GetTransportationManager();
204 
205  // fFieldExists= transportMgr->GetFieldManager()->DoesFieldExist();
206  // return fFieldExists;
207  return transportMgr->GetFieldManager()->DoesFieldExist();
208  }
209 
210  private:
211 
212  G4Navigator* fLinearNavigator;
213  G4PropagatorInField* fFieldPropagator;
215  G4ThreeVector fTransportEndPosition = G4ThreeVector(0.0, 0.0, 0.0);
216  G4ThreeVector fTransportEndMomentumDir = G4ThreeVector(0.0, 0.0, 0.0);
217  G4double fTransportEndKineticEnergy = 0.0;
218  G4ThreeVector fTransportEndSpin = G4ThreeVector(0.0, 0.0, 0.0);
219  G4bool fMomentumChanged = true;
220  G4bool fEndGlobalTimeComputed = false;
221  G4double fCandidateEndGlobalTime = 0.0;
223  G4bool fAnyFieldExists = false;
225  G4bool fNewTrack = true;
226  G4bool fFirstStepInVolume = true;
227  G4bool fLastStepInVolume = false;
228  G4bool fGeometryLimitedStep = true;
230  G4bool fFieldExertedForce = false;
232  G4TouchableHandle fCurrentTouchableHandle;
234  G4ThreeVector fPreviousSftOrigin;
235  G4double fPreviousSafety;
237  G4ParticleChangeForTransport fParticleChange;
239  G4double fEndPointDistance;
241  // Thresholds for looping particles:
242  //
243  G4double fThreshold_Warning_Energy = 1.0 * CLHEP::keV;
244  G4double fThreshold_Important_Energy = 1.0 * CLHEP::MeV;
245  G4int fThresholdTrials = 10;
246  // Above 'important' energy a 'looping' particle in field will
247  // *NOT* be abandoned, except after fThresholdTrials attempts.
249  // unstable loopers ( 0 = never )
250  // Counter for steps in which particle reports 'looping',
251  // ( Used if it is above 'Important' Energy. )
252  G4int fNoLooperTrials = 0;
254  // Statistics for tracks abandoned due to looping, only used as flag, no looping neutral particles
255  //
256  G4double fSumEnergyKilled = 0.0;
262 
263  G4SafetyHelper* fpSafetyHelper;
265  private:
266 
267  friend class G4CoupledTransportation;
268  static G4bool fUseMagneticMoment;
269  static G4bool fUseGravity;
270  static G4bool fSilenceLooperWarnings;
272  };
274 }
275 
276 #endif
Belle2::G4LongLivedNeutralTransportation::StartTracking
void StartTracking(G4Track *aTrack)
Reset state for new (potentially resumed) track.
Definition: G4LongLivedNeutralTransportation.cc:457
Belle2::G4LongLivedNeutralTransportation::SetLowLooperThresholds
void SetLowLooperThresholds()
Shortcut method - old values (meant for HEP)
Definition: G4LongLivedNeutralTransportation.cc:543
Belle2::G4LongLivedNeutralTransportation::fAbandonUnstableTrials
G4int fAbandonUnstableTrials
Number of trials after which to abandon.
Definition: G4LongLivedNeutralTransportation.h:248
Belle2::G4LongLivedNeutralTransportation::AtRestGetPhysicalInteractionLength
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
No operation in AtRestGPIL.
Definition: G4LongLivedNeutralTransportation.h:172
Belle2::G4LongLivedNeutralTransportation::fShortStepOptimisation
G4bool fShortStepOptimisation
Whether to avoid calling G4Navigator for short step ( < safety) If using it, the safety estimate for ...
Definition: G4LongLivedNeutralTransportation.h:261
Belle2::G4LongLivedNeutralTransportation::fLinearNavigator
G4Navigator * fLinearNavigator
Propagator used to transport the particle.
Definition: G4LongLivedNeutralTransportation.h:212
Belle2::G4LongLivedNeutralTransportation::fUseMagneticMoment
static G4bool fUseMagneticMoment
Flag take into account magnetic moment.
Definition: G4LongLivedNeutralTransportation.h:268
Belle2::G4LongLivedNeutralTransportation::G4LongLivedNeutralTransportation
G4LongLivedNeutralTransportation(G4int verbosityLevel=1)
Constructor.
Definition: G4LongLivedNeutralTransportation.cc:35
Belle2::G4LongLivedNeutralTransportation::fNoLooperTrials
G4int fNoLooperTrials
Counter for steps in which particle reports 'looping', if it is above 'Important' Energy.
Definition: G4LongLivedNeutralTransportation.h:252
Belle2::G4LongLivedNeutralTransportation::fPreviousSafety
G4double fPreviousSafety
Remember last safety value.
Definition: G4LongLivedNeutralTransportation.h:235
Belle2::G4LongLivedNeutralTransportation::AlongStepGetPhysicalInteractionLength
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
G4VProcess::AlongStepGetPhysicalInteractionLength() implementation,.
Definition: G4LongLivedNeutralTransportation.cc:113
Belle2::G4LongLivedNeutralTransportation::GetSumEnergyKilled
G4double GetSumEnergyKilled() const
Access fSumEnergyKilled.
Belle2::G4LongLivedNeutralTransportation::fGeometryLimitedStep
G4bool fGeometryLimitedStep
Flag to determine whether a boundary was reached.
Definition: G4LongLivedNeutralTransportation.h:228
Belle2::G4LongLivedNeutralTransportation::fLastStepInVolume
G4bool fLastStepInVolume
Flag last step in a geom.
Definition: G4LongLivedNeutralTransportation.h:227
Belle2::G4LongLivedNeutralTransportation::fThreshold_Important_Energy
G4double fThreshold_Important_Energy
Give a few trial above this E for looping particle.
Definition: G4LongLivedNeutralTransportation.h:244
Belle2::G4LongLivedNeutralTransportation::fTransportEndKineticEnergy
G4double fTransportEndKineticEnergy
The particle's state after this Step, Store for DoIt.
Definition: G4LongLivedNeutralTransportation.h:217
Belle2::G4LongLivedNeutralTransportation::fTransportEndSpin
G4ThreeVector fTransportEndSpin
The particle's state after this Step, Store for DoIt.
Definition: G4LongLivedNeutralTransportation.h:218
Belle2::G4LongLivedNeutralTransportation::GetThresholdWarningEnergy
G4double GetThresholdWarningEnergy() const
Access fThreshold_Warning_Energy.
Belle2::G4LongLivedNeutralTransportation::fpSafetyHelper
G4SafetyHelper * fpSafetyHelper
To pass it the safety value obtained.
Definition: G4LongLivedNeutralTransportation.h:263
Belle2::G4LongLivedNeutralTransportation::fThresholdTrials
G4int fThresholdTrials
Number of trials an important looper survives.
Definition: G4LongLivedNeutralTransportation.h:245
Belle2::G4LongLivedNeutralTransportation::fCandidateEndGlobalTime
G4double fCandidateEndGlobalTime
The particle's state after this Step, Store for DoIt.
Definition: G4LongLivedNeutralTransportation.h:221
Belle2::G4LongLivedNeutralTransportation::PrintStatistics
void PrintStatistics(std::ostream &outStr) const
returns current logging info of the algorithm
Definition: G4LongLivedNeutralTransportation.cc:98
Belle2::G4LongLivedNeutralTransportation::fPreviousSftOrigin
G4ThreeVector fPreviousSftOrigin
Remember last safety origin.
Definition: G4LongLivedNeutralTransportation.h:234
Belle2::G4LongLivedNeutralTransportation::ReportLooperThresholds
void ReportLooperThresholds()
Set low thresholds - for low-E applications.
Belle2::G4LongLivedNeutralTransportation::fEndPointDistance
G4double fEndPointDistance
Endpoint distance.
Definition: G4LongLivedNeutralTransportation.h:239
Belle2::G4LongLivedNeutralTransportation::~G4LongLivedNeutralTransportation
~G4LongLivedNeutralTransportation()
Destructor.
Definition: G4LongLivedNeutralTransportation.cc:87
Belle2::G4LongLivedNeutralTransportation::fParticleChange
G4ParticleChangeForTransport fParticleChange
New ParticleChange.
Definition: G4LongLivedNeutralTransportation.h:237
Belle2::G4LongLivedNeutralTransportation::fUseGravity
static G4bool fUseGravity
Flag take into account gravity.
Definition: G4LongLivedNeutralTransportation.h:269
Belle2::G4LongLivedNeutralTransportation::PostStepDoIt
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
G4VProcess::PostStepDoIt() implementation,.
Definition: G4LongLivedNeutralTransportation.cc:350
Belle2::G4LongLivedNeutralTransportation::fFieldExertedForce
G4bool fFieldExertedForce
During current step.
Definition: G4LongLivedNeutralTransportation.h:230
Belle2::G4LongLivedNeutralTransportation::GetPropagatorInField
G4PropagatorInField * GetPropagatorInField()
Access fFieldPropagator, the assistant class that Propagate in a Field.
Belle2::G4LongLivedNeutralTransportation::PostStepGetPhysicalInteractionLength
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
G4VProcess::PostStepGetPhysicalInteractionLength() implementation.
Definition: G4LongLivedNeutralTransportation.cc:338
Belle2::G4LongLivedNeutralTransportation::fFirstStepInVolume
G4bool fFirstStepInVolume
Flag first step in a geom.
Definition: G4LongLivedNeutralTransportation.h:226
Belle2::G4LongLivedNeutralTransportation::AtRestDoIt
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
No operation in AtRestDoIt.
Definition: G4LongLivedNeutralTransportation.h:179
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::G4LongLivedNeutralTransportation::DoesGlobalFieldExist
G4bool DoesGlobalFieldExist()
Checks whether a field exists for the "global" field manager.
Definition: G4LongLivedNeutralTransportation.h:200
Belle2::G4LongLivedNeutralTransportation::ProcessDescription
virtual void ProcessDescription(std::ostream &outFile) const
G4LongLivedNeutralTransportation::ProcessDescription()
Definition: G4LongLivedNeutralTransportation.cc:562
Belle2::G4LongLivedNeutralTransportation::fEndGlobalTimeComputed
G4bool fEndGlobalTimeComputed
The particle's state after this Step, Store for DoIt.
Definition: G4LongLivedNeutralTransportation.h:220
Belle2::G4LongLivedNeutralTransportation::fAnyFieldExists
G4bool fAnyFieldExists
Flag for existing fields.
Definition: G4LongLivedNeutralTransportation.h:223
Belle2::G4LongLivedNeutralTransportation
Concrete class that does the geometrical transport.
Definition: G4LongLivedNeutralTransportation.h:38
Belle2::G4LongLivedNeutralTransportation::SetSilenceLooperWarnings
static void SetSilenceLooperWarnings(G4bool val)
Do not warn about 'looping' particles.
Definition: G4LongLivedNeutralTransportation.cc:511
Belle2::G4LongLivedNeutralTransportation::EnableShortStepOptimisation
void EnableShortStepOptimisation(G4bool optimise=true)
Whether short steps < safety will avoid to call Navigator (if field=0)
Belle2::G4LongLivedNeutralTransportation::ResetKilledStatistics
void ResetKilledStatistics(G4int report=1)
Statistics for tracks killed (currently due to looping in field)
Belle2::G4LongLivedNeutralTransportation::GetSilenceLooperWarnings
static G4bool GetSilenceLooperWarnings()
Do not throw exception about 'looping' particles.
Definition: G4LongLivedNeutralTransportation.cc:519
Belle2::G4LongLivedNeutralTransportation::fFieldPropagator
G4PropagatorInField * fFieldPropagator
Propagator used to transport the particle.
Definition: G4LongLivedNeutralTransportation.h:213
Belle2::G4LongLivedNeutralTransportation::AlongStepDoIt
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
G4VProcess::AlongStepDoIt() implementation,.
Definition: G4LongLivedNeutralTransportation.cc:263
Belle2::G4LongLivedNeutralTransportation::fSumEnergyKilled
G4double fSumEnergyKilled
Sum of abandoned looping tracks energies.
Definition: G4LongLivedNeutralTransportation.h:256
Belle2::G4LongLivedNeutralTransportation::SetPropagatorInField
void SetPropagatorInField(G4PropagatorInField *pFieldPropagator)
Access/set the assistant class that Propagate in a Field.
Belle2::G4LongLivedNeutralTransportation::SetHighLooperThresholds
void SetHighLooperThresholds()
Get/Set parameters for killing loopers: Above 'important' energy a 'looping' particle in field will N...
Definition: G4LongLivedNeutralTransportation.cc:527
Belle2::G4LongLivedNeutralTransportation::fNewTrack
G4bool fNewTrack
Flag from StartTracking.
Definition: G4LongLivedNeutralTransportation.h:225
Belle2::G4LongLivedNeutralTransportation::GetThresholdTrials
G4int GetThresholdTrials() const
Access fThresholdTrials.
Belle2::G4LongLivedNeutralTransportation::fThreshold_Warning_Energy
G4double fThreshold_Warning_Energy
Warn above this energy about looping particle.
Definition: G4LongLivedNeutralTransportation.h:243
Belle2::G4LongLivedNeutralTransportation::fTransportEndPosition
G4ThreeVector fTransportEndPosition
The particle's state after this Step, Store for DoIt.
Definition: G4LongLivedNeutralTransportation.h:215
Belle2::G4LongLivedNeutralTransportation::fSilenceLooperWarnings
static G4bool fSilenceLooperWarnings
Flag to Supress all 'looper' warnings.
Definition: G4LongLivedNeutralTransportation.h:270
Belle2::G4LongLivedNeutralTransportation::fMomentumChanged
G4bool fMomentumChanged
The particle's state after this Step, Store for DoIt.
Definition: G4LongLivedNeutralTransportation.h:219
Belle2::G4LongLivedNeutralTransportation::fTransportEndMomentumDir
G4ThreeVector fTransportEndMomentumDir
The particle's state after this Step, Store for DoIt.
Definition: G4LongLivedNeutralTransportation.h:216
Belle2::G4LongLivedNeutralTransportation::GetThresholdImportantEnergy
G4double GetThresholdImportantEnergy() const
Access fThreshold_Important_Energy.
Belle2::G4LongLivedNeutralTransportation::GetMaxEnergyKilled
G4double GetMaxEnergyKilled() const
Print values of looper thresholds.
Belle2::G4LongLivedNeutralTransportation::fCurrentTouchableHandle
G4TouchableHandle fCurrentTouchableHandle
Current touchable handle.
Definition: G4LongLivedNeutralTransportation.h:232