CMS 3D CMS Logo

MonopoleTransportation.h
Go to the documentation of this file.
1 //
2 // =======================================================================
3 //
4 // Class MonopoleTransportation
5 //
6 // Created: 3 May 2010, J. Apostolakis, B. Bozsogi
7 // G4MonopoleTransportation class for
8 // Geant4 extended example "monopole"
9 //
10 // Adopted for CMSSW by V.Ivanchenko 30 April 2018
11 // from Geant4 global tag geant4-10-04-ref-03
12 //
13 // =======================================================================
14 //
15 // Class description:
16 //
17 // G4MonopoleTransportation is a process responsible for the transportation of
18 // magnetic monopoles, i.e. the geometrical propagation encountering the
19 // geometrical sub-volumes of the detectors.
20 // It is also tasked with part of updating the "safety".
21 // For monopoles, uses a different equation of motion and ignores energy
22 // conservation.
23 //
24 
25 #ifndef SimG4Core_PhysicsLists_MonopoleTransportation_h
26 #define SimG4Core_PhysicsLists_MonopoleTransportation_h 1
27 
28 #include "G4VProcess.hh"
29 #include "G4FieldManager.hh"
30 
31 #include "G4Navigator.hh"
32 #include "G4TransportationManager.hh"
33 #include "G4PropagatorInField.hh"
34 #include "G4Track.hh"
35 #include "G4Step.hh"
36 #include "G4ParticleChangeForTransport.hh"
37 #include "CLHEP/Units/SystemOfUnits.h"
38 
39 #include <memory>
40 
41 class G4SafetyHelper;
42 class Monopole;
43 class CMSFieldManager;
44 
45 class MonopoleTransportation : public G4VProcess
46 {
47 public:
48 
50  ~MonopoleTransportation() override;
51 
53  const G4Track& track,
54  G4double previousStepSize,
55  G4double currentMinimumStep,
56  G4double& currentSafety,
57  G4GPILSelection* selection
58  ) override;
59 
60  G4VParticleChange* AlongStepDoIt(
61  const G4Track& track,
62  const G4Step& stepData
63  ) override;
64 
65  G4VParticleChange* PostStepDoIt(
66  const G4Track& track,
67  const G4Step& stepData
68  ) override;
69  // Responsible for the relocation.
70 
72  const G4Track& ,
73  G4double previousStepSize,
74  G4ForceCondition* pForceCond
75  ) override;
76  // Forces the PostStepDoIt action to be called,
77  // but does not limit the step.
78 
79  G4PropagatorInField* GetPropagatorInField();
80  void SetPropagatorInField( G4PropagatorInField* pFieldPropagator);
81  // Access/set the assistant class that Propagate in a Field.
82 
83  inline G4double GetThresholdWarningEnergy() const;
84  inline G4double GetThresholdImportantEnergy() const;
85  inline G4int GetThresholdTrials() const;
86 
87  inline void SetThresholdWarningEnergy( G4double newEnWarn );
88  inline void SetThresholdImportantEnergy( G4double newEnImp );
89  inline void SetThresholdTrials(G4int newMaxTrials );
90 
91  // Get/Set parameters for killing loopers:
92  // Above 'important' energy a 'looping' particle in field will
93  // *NOT* be abandoned, except after fThresholdTrials attempts.
94  // Below Warning energy, no verbosity for looping particles is issued
95 
96  inline G4double GetMaxEnergyKilled() const;
97  inline G4double GetSumEnergyKilled() const;
98  void ResetKilledStatistics( G4int report = 1);
99  // Statistics for tracks killed (currently due to looping in field)
100 
101  inline void EnableShortStepOptimisation(G4bool optimise=true);
102  // Whether short steps < safety will avoid to call Navigator (if field=0)
103 
104  G4double AtRestGetPhysicalInteractionLength(const G4Track&,
105  G4ForceCondition*) override;
106 
107  G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&) override;
108  // No operation in AtRestDoIt.
109 
110  void StartTracking(G4Track* aTrack) override;
111  // Reset state for new (potentially resumed) track
112 
113 protected:
114 
115  G4bool DoesGlobalFieldExist();
116  // Checks whether a field exists for the "global" field manager.
117 
118 private:
119 
121 
123 
124  G4Navigator* fLinearNavigator;
125  G4PropagatorInField* fFieldPropagator;
126  // The Propagators used to transport the particle
127 
128  G4ThreeVector fTransportEndPosition;
131  G4ThreeVector fTransportEndSpin;
133 
136  // The particle's state after this Step, Store for DoIt
137 
139 
140  G4TouchableHandle fCurrentTouchableHandle;
141 
143  // Flag to determine whether a boundary was reached.
144 
145  G4ThreeVector fPreviousSftOrigin;
146  G4double fPreviousSafety;
147  // Remember last safety origin & value.
148 
149  G4ParticleChangeForTransport fParticleChange;
150  // New ParticleChange
151 
153 
154  // Thresholds for looping particles:
155  //
156  G4double fThreshold_Warning_Energy; // Warn above this energy
157  G4double fThreshold_Important_Energy; // Hesitate above this
158  G4int fThresholdTrials; // for this no of trials
159  // Above 'important' energy a 'looping' particle in field will
160  // *NOT* be abandoned, except after fThresholdTrials attempts.
161 
162  // Counter for steps in which particle reports 'looping',
163  // if it is above 'Important' Energy
165  // Statistics for tracks abandoned
168 
169  // Whether to avoid calling G4Navigator for short step ( < safety)
170  // If using it, the safety estimate for endpoint will likely be smaller.
172 
173  G4SafetyHelper* fpSafetyHelper; // To pass it the safety value obtained
174 
175 };
176 
177 inline void
178 MonopoleTransportation::SetPropagatorInField( G4PropagatorInField* pFieldPropagator)
179 {
180  fFieldPropagator= pFieldPropagator;
181 }
182 
183 inline G4PropagatorInField* MonopoleTransportation::GetPropagatorInField()
184 {
185  return fFieldPropagator;
186 }
187 
189 {
190  G4TransportationManager* transportMgr =
191  G4TransportationManager::GetTransportationManager();
192  return transportMgr->GetFieldManager()->DoesFieldExist();
193 }
194 
196 {
198 }
199 
201 {
203 }
204 
206 {
207  return fThresholdTrials;
208 }
209 
210 inline void MonopoleTransportation::SetThresholdWarningEnergy( G4double newEnWarn )
211 {
212  fThreshold_Warning_Energy= newEnWarn;
213 }
214 
216 {
217  fThreshold_Important_Energy = newEnImp;
218 }
219 
220 inline void MonopoleTransportation::SetThresholdTrials(G4int newMaxTrials )
221 {
222  fThresholdTrials = newMaxTrials;
223 }
224 
225 // Get parameters for killing loopers:
226 // Above 'important' energy a 'looping' particle in field will
227 // *NOT* be abandoned, except after fThresholdTrials attempts.
228 // Below Warning energy, no verbosity for looping particles is issued
229 
231 {
232  return fMaxEnergyKilled;
233 }
234 
236 {
237  return fSumEnergyKilled;
238 }
239 
240 inline void
242 {
243  fShortStepOptimisation=optimiseShortStep;
244 }
245 
246 #endif
G4PropagatorInField * GetPropagatorInField()
void SetThresholdTrials(G4int newMaxTrials)
void StartTracking(G4Track *aTrack) override
selection
main part
Definition: corrVsCorr.py:99
G4double GetThresholdWarningEnergy() const
G4double GetMaxEnergyKilled() const
void ResetKilledStatistics(G4int report=1)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *) override
G4double GetSumEnergyKilled() const
void SetThresholdWarningEnergy(G4double newEnWarn)
G4double GetThresholdImportantEnergy() const
G4ParticleChangeForTransport fParticleChange
void SetThresholdImportantEnergy(G4double newEnImp)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond) override
MonopoleTransportation(const Monopole *p, G4int verbosityLevel=1)
G4PropagatorInField * fFieldPropagator
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &) override
void SetPropagatorInField(G4PropagatorInField *pFieldPropagator)
G4TouchableHandle fCurrentTouchableHandle
void EnableShortStepOptimisation(G4bool optimise=true)
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData) override
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData) override