Belle II Software  release-08-01-10
G4MonopoleTransportation.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 exoticphysics/monopole/*
10 
11 #pragma once
12 
13 // inline include at the end of the file
14 
15 #include <G4VProcess.hh>
16 #include <G4FieldManager.hh>
17 #include <G4Navigator.hh>
18 #include <G4TransportationManager.hh>
19 #include <G4PropagatorInField.hh>
20 #include <G4Track.hh>
21 #include <G4Step.hh>
22 #include <G4ParticleChangeForTransport.hh>
23 
24 #include <simulation/monopoles/G4MonopoleFieldSetup.h>
25 
26 class G4SafetyHelper;
27 
28 namespace Belle2 {
34  namespace Monopoles {
35 
36  class G4Monopole;
37 
44  class G4MonopoleTransportation : public G4VProcess {
45 
46  public:
47 
53  explicit G4MonopoleTransportation(const G4Monopole* mpl, G4int verb = 1);
58 
75  const G4Track& track,
76  G4double previousStepSize,
77  G4double currentMinimumStep,
78  G4double& currentSafety,
79  G4GPILSelection* selection
80  );
81 
90  virtual G4VParticleChange* AlongStepDoIt(
91  const G4Track& track,
92  const G4Step& stepData
93  );
94 
104  virtual G4VParticleChange* PostStepDoIt(
105  const G4Track& track,
106  const G4Step& stepData
107  );
108 
120  virtual G4double PostStepGetPhysicalInteractionLength(
121  const G4Track& track,
122  G4double previousStepSize,
123  G4ForceCondition* pForceCond
124  );
125 
126  inline G4PropagatorInField* GetPropagatorInField();
127  inline void SetPropagatorInField(G4PropagatorInField*
128  pFieldPropagator);
130  inline G4double GetThresholdWarningEnergy() const;
131  inline G4double GetThresholdImportantEnergy() const;
132  inline G4int GetThresholdTrials() const;
134  inline void SetThresholdWarningEnergy(G4double newEnWarn);
135  inline void SetThresholdImportantEnergy(G4double newEnImp);
136  inline void SetThresholdTrials(G4int newMaxTrials);
138  // Get/Set parameters for killing loopers:
139  // Above 'important' energy a 'looping' particle in field will
140  // *NOT* be abandoned, except after fThresholdTrials attempts.
141  // Below Warning energy, no verbosity for looping particles is issued
142 
143  inline G4double GetMaxEnergyKilled() const;
144  inline G4double GetSumEnergyKilled() const;
145  inline void ResetKilledStatistics(G4int report = 1);
147  inline void EnableShortStepOptimisation(G4bool optimise =
148  true);
150  public:
151 
156  const G4Track&,
157  G4ForceCondition*
158  ) { return -1.0; };
159 
163  virtual G4VParticleChange* AtRestDoIt(
164  const G4Track&,
165  const G4Step&
166  ) {return 0;};
167 
168  virtual void StartTracking(G4Track* aTrack);
170  protected:
171 
172  G4bool DoesGlobalFieldExist();
174  private:
175 
180  G4Navigator* fLinearNavigator;
181  G4PropagatorInField* fFieldPropagator;
183  G4ThreeVector fTransportEndPosition;
184  G4ThreeVector fTransportEndMomentumDir;
186  G4ThreeVector fTransportEndSpin;
188  // G4bool fEnergyChanged; /**< The particle's state after this Step, Store for DoIt*/
194  G4TouchableHandle fCurrentTouchableHandle;
198  G4ThreeVector fPreviousSftOrigin;
199  G4double fPreviousSafety;
201  G4ParticleChangeForTransport fParticleChange;
203  G4double endpointDistance;
205  G4double fThreshold_Trap_Energy; //FIXME should be dependent on e.g. density
206 
207  // Thresholds for looping particles:
211  // Above 'important' energy a 'looping' particle in field will
212  // *NOT* be abandoned, except after fThresholdTrials attempts.
213  // G4double fUnimportant_Energy;
214  // Below this energy, no verbosity for looping particles is issued
215 
217  G4double fSumEnergyKilled;
218  G4double fMaxEnergyKilled;
225 
226  G4SafetyHelper* fpSafetyHelper;
228  };
229 
230  } //end Monopoles namespace
231 
233 } //end Belle2 namespace
234 
235 #include <simulation/monopoles/G4MonopoleTransportationInline.h>
Monopole field setup singleton class, that takes care of switching between conventional particle tran...
Concrete class that does the geometrical transport.
G4TouchableHandle fCurrentTouchableHandle
Current touchable handle.
G4int GetThresholdTrials() const
Access fThresholdTrials.
G4double GetThresholdImportantEnergy() const
Access fThreshold_Important_Energy.
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
G4VProcess::PostStepDoIt() implementation.
G4double fCandidateEndGlobalTime
The particle's state after this Step, Store for DoIt.
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
G4VProcess::AlongStepGetPhysicalInteractionLength() implementation.
G4int fThresholdTrials
Nubmer of trials for looping particles.
G4PropagatorInField * GetPropagatorInField()
Access fFieldPropagator, the assistant class that Propagate in a Field.
G4double GetMaxEnergyKilled() const
Access fMaxEnergyKilled.
void SetThresholdWarningEnergy(G4double newEnWarn)
Set fThreshold_Warning_Energy.
void SetThresholdImportantEnergy(G4double newEnImp)
Set fThreshold_Important_Energy.
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
No operation in AtRestDoIt.
G4double fMaxEnergyKilled
Max of abandoned looping tracks energies.
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.
G4double fPreviousSafety
Remember last safety value.
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
G4VProcess::AlongStepDoIt() implementation, Proposes changes during step to fParticleChange of this c...
G4SafetyHelper * fpSafetyHelper
To pass it the safety value obtained.
void SetPropagatorInField(G4PropagatorInField *pFieldPropagator)
Set fFieldPropagator, the assistant class that Propagate in a Field.
G4double fSumEnergyKilled
Sum of abandoned looping tracks energies.
G4ThreeVector fPreviousSftOrigin
Remember last safety origin.
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
No operation in AtRestDoIt.
virtual 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.
G4double fThreshold_Important_Energy
Hesitate above this about looping particle for a certain no of trials.
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.
virtual void StartTracking(G4Track *aTrack)
Reset state for new (potentially resumed) track.
G4ParticleChangeForTransport fParticleChange
New ParticleChange.
G4MonopoleTransportation(const G4Monopole *mpl, G4int verb=1)
Constructor.
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 DoesGlobalFieldExist()
Checks whether a field exists for the "global" field manager.
const G4Monopole * fParticleDef
Monopole definition for charge and mass reference.
G4bool fParticleIsLooping
Is the monopole stuck in looping.
G4double GetSumEnergyKilled() const
Access fSumEnergyKilled.
void ResetKilledStatistics(G4int report=1)
Statistics for tracks killed (currently due to looping in field)
G4double fThreshold_Trap_Energy
Assume monopoles below this can bound to material.
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.
G4MonopoleFieldSetup * fMagSetup
Monpole field setup.
A class to hold monopole description as a particle.
Definition: G4Monopole.h:33
Abstract base class for different kinds of events.