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 class G4SafetyHelper;
40 class Monopole;
41 class CMSFieldManager;
42 
43 namespace sim {
44  class ChordFinderSetter;
45 }
46 
47 #include <memory>
48 
49 class MonopoleTransportation : public G4VProcess
50 {
51 public:
52 
54  G4int verbosityLevel= 1);
55  ~MonopoleTransportation() override;
56 
57  G4double AlongStepGetPhysicalInteractionLength(
58  const G4Track& track,
59  G4double previousStepSize,
60  G4double currentMinimumStep,
61  G4double& currentSafety,
62  G4GPILSelection* selection
63  ) override;
64 
65  G4VParticleChange* AlongStepDoIt(
66  const G4Track& track,
67  const G4Step& stepData
68  ) override;
69 
70  G4VParticleChange* PostStepDoIt(
71  const G4Track& track,
72  const G4Step& stepData
73  ) override;
74  // Responsible for the relocation.
75 
76  G4double PostStepGetPhysicalInteractionLength(
77  const G4Track& ,
78  G4double previousStepSize,
79  G4ForceCondition* pForceCond
80  ) override;
81  // Forces the PostStepDoIt action to be called,
82  // but does not limit the step.
83 
84  G4PropagatorInField* GetPropagatorInField();
85  void SetPropagatorInField( G4PropagatorInField* pFieldPropagator);
86  // Access/set the assistant class that Propagate in a Field.
87 
88  inline G4double GetThresholdWarningEnergy() const;
89  inline G4double GetThresholdImportantEnergy() const;
90  inline G4int GetThresholdTrials() const;
91 
92  inline void SetThresholdWarningEnergy( G4double newEnWarn );
93  inline void SetThresholdImportantEnergy( G4double newEnImp );
94  inline void SetThresholdTrials(G4int newMaxTrials );
95 
96  // Get/Set parameters for killing loopers:
97  // Above 'important' energy a 'looping' particle in field will
98  // *NOT* be abandoned, except after fThresholdTrials attempts.
99  // Below Warning energy, no verbosity for looping particles is issued
100 
101  inline G4double GetMaxEnergyKilled() const;
102  inline G4double GetSumEnergyKilled() const;
103  inline void ResetKilledStatistics( G4int report = 1);
104  // Statistics for tracks killed (currently due to looping in field)
105 
106  inline void EnableShortStepOptimisation(G4bool optimise=true);
107  // Whether short steps < safety will avoid to call Navigator (if field=0)
108 
110  const G4Track& ,
111  G4ForceCondition*
112  ) override { return -1.0; };
113  // No operation in AtRestDoIt.
114 
115  G4VParticleChange* AtRestDoIt(
116  const G4Track& ,
117  const G4Step&
118  ) override {return nullptr;};
119  // No operation in AtRestDoIt.
120 
121  void StartTracking(G4Track* aTrack) override;
122  // Reset state for new (potentially resumed) track
123 
124 protected:
125 
126  G4bool DoesGlobalFieldExist();
127  // Checks whether a field exists for the "global" field manager.
128 
129 private:
130 
132 
133  sim::ChordFinderSetter *fChordFinderSetter;;
135 
136  G4Navigator* fLinearNavigator;
137  G4PropagatorInField* fFieldPropagator;
138  // The Propagators used to transport the particle
139 
140  G4ThreeVector fTransportEndPosition;
143  G4ThreeVector fTransportEndSpin;
145 
148  // The particle's state after this Step, Store for DoIt
149 
151 
152  G4TouchableHandle fCurrentTouchableHandle;
153 
155  // Flag to determine whether a boundary was reached.
156 
157  G4ThreeVector fPreviousSftOrigin;
158  G4double fPreviousSafety;
159  // Remember last safety origin & value.
160 
161  G4ParticleChangeForTransport fParticleChange;
162  // New ParticleChange
163 
165 
166  // Thresholds for looping particles:
167  //
168  G4double fThreshold_Warning_Energy; // Warn above this energy
169  G4double fThreshold_Important_Energy; // Hesitate above this
170  G4int fThresholdTrials; // for this no of trials
171  // Above 'important' energy a 'looping' particle in field will
172  // *NOT* be abandoned, except after fThresholdTrials attempts.
173 
174  // Counter for steps in which particle reports 'looping',
175  // if it is above 'Important' Energy
177  // Statistics for tracks abandoned
180 
181  // Whether to avoid calling G4Navigator for short step ( < safety)
182  // If using it, the safety estimate for endpoint will likely be smaller.
184 
185  G4SafetyHelper* fpSafetyHelper; // To pass it the safety value obtained
186 
187 };
188 
189 inline void
190 MonopoleTransportation::SetPropagatorInField( G4PropagatorInField* pFieldPropagator)
191 {
192  fFieldPropagator= pFieldPropagator;
193 }
194 
195 inline G4PropagatorInField* MonopoleTransportation::GetPropagatorInField()
196 {
197  return fFieldPropagator;
198 }
199 
201 {
202  G4TransportationManager* transportMgr;
203  transportMgr= G4TransportationManager::GetTransportationManager();
204 
205  // fFieldExists= transportMgr->GetFieldManager()->DoesFieldExist();
206  // return fFieldExists;
207  return transportMgr->GetFieldManager()->DoesFieldExist();
208 }
209 
211 {
212  return fThreshold_Warning_Energy;
213 }
214 
216 {
217  return fThreshold_Important_Energy;
218 }
219 
221 {
222  return fThresholdTrials;
223 }
224 
225 inline void MonopoleTransportation::SetThresholdWarningEnergy( G4double newEnWarn )
226 {
227  fThreshold_Warning_Energy= newEnWarn;
228 }
229 
231 {
232  fThreshold_Important_Energy = newEnImp;
233 }
234 
235 inline void MonopoleTransportation::SetThresholdTrials(G4int newMaxTrials )
236 {
237  fThresholdTrials = newMaxTrials;
238 }
239 
240 // Get/Set parameters for killing loopers:
241 // Above 'important' energy a 'looping' particle in field will
242 // *NOT* be abandoned, except after fThresholdTrials attempts.
243 // Below Warning energy, no verbosity for looping particles is issued
244 
246 {
247  return fMaxEnergyKilled;
248 }
249 
251 {
252  return fSumEnergyKilled;
253 }
254 
256 {
257  if( report ) {
258  G4cout << " MonopoleTransportation: Statistics for looping particles " << G4endl;
259  G4cout << " Sum of energy of loopers killed: " << fSumEnergyKilled << G4endl;
260  G4cout << " Max energy of loopers killed: " << fMaxEnergyKilled << G4endl;
261  }
262 
263  fSumEnergyKilled= 0;
264  fMaxEnergyKilled= -1.0*CLHEP::GeV;
265 }
266 // Statistics for tracks killed (currently due to looping in field)
267 
268 inline void
270 {
271  fShortStepOptimisation=optimiseShortStep;
272 }
273 
274 #endif
G4PropagatorInField * GetPropagatorInField()
const double GeV
Definition: MathUtil.h:16
void SetThresholdTrials(G4int newMaxTrials)
selection
main part
Definition: corrVsCorr.py:98
G4double GetThresholdWarningEnergy() const
G4double GetMaxEnergyKilled() const
void ResetKilledStatistics(G4int report=1)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *) override
G4double GetSumEnergyKilled() const
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &) override
void SetThresholdWarningEnergy(G4double newEnWarn)
G4double GetThresholdImportantEnergy() const
G4ParticleChangeForTransport fParticleChange
void SetThresholdImportantEnergy(G4double newEnImp)
Definition: RunManager.h:28
G4PropagatorInField * fFieldPropagator
void SetPropagatorInField(G4PropagatorInField *pFieldPropagator)
G4TouchableHandle fCurrentTouchableHandle
void EnableShortStepOptimisation(G4bool optimise=true)