Geant4  10.00.p01
G4MonopoleTransportation.hh
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
28 //
29 // $Id: G4MonopoleTransportation.hh 69705 2013-05-13 09:09:52Z gcosmo $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 //
34 // ------------------------------------------------------------
35 // GEANT 4 include file implementation
36 // ------------------------------------------------------------
37 //
38 // Class description:
39 //
40 // G4MonopoleTransportation is a process responsible for the transportation of
41 // magnetic monopoles, i.e. the geometrical propagation encountering the
42 // geometrical sub-volumes of the detectors.
43 // It is also tasked with part of updating the "safety".
44 
45 // =======================================================================
46 // Created: 3 May 2010, J. Apostolakis, B. Bozsogi
47 // =======================================================================
48 
49 #ifndef G4MonopoleTransportation_hh
50 #define G4MonopoleTransportation_hh 1
51 
52 #include "G4VProcess.hh"
53 #include "G4FieldManager.hh"
54 
55 #include "G4Navigator.hh"
57 #include "G4PropagatorInField.hh"
58 #include "G4Track.hh"
59 #include "G4Step.hh"
61 #include "G4MonopoleFieldSetup.hh"
62 
63 class G4SafetyHelper;
64 class G4Monopole;
65 
67 {
68  // Concrete class that does the geometrical transport
69 
70  public: // with description
71 
72  G4MonopoleTransportation(const G4Monopole* p, G4int verbosityLevel= 1);
74 
76  const G4Track& track,
77  G4double previousStepSize,
78  G4double currentMinimumStep,
79  G4double& currentSafety,
80  G4GPILSelection* selection
81  );
82 
84  const G4Track& track,
85  const G4Step& stepData
86  );
87 
89  const G4Track& track,
90  const G4Step& stepData
91  );
92  // Responsible for the relocation.
93 
95  const G4Track& ,
96  G4double previousStepSize,
97  G4ForceCondition* pForceCond
98  );
99  // Forces the PostStepDoIt action to be called,
100  // but does not limit the step.
101 
103  void SetPropagatorInField( G4PropagatorInField* pFieldPropagator);
104  // Access/set the assistant class that Propagate in a Field.
105 
106  inline G4double GetThresholdWarningEnergy() const;
107  inline G4double GetThresholdImportantEnergy() const;
108  inline G4int GetThresholdTrials() const;
109 
110  inline void SetThresholdWarningEnergy( G4double newEnWarn );
111  inline void SetThresholdImportantEnergy( G4double newEnImp );
112  inline void SetThresholdTrials(G4int newMaxTrials );
113 
114  // Get/Set parameters for killing loopers:
115  // Above 'important' energy a 'looping' particle in field will
116  // *NOT* be abandoned, except after fThresholdTrials attempts.
117  // Below Warning energy, no verbosity for looping particles is issued
118 
119  inline G4double GetMaxEnergyKilled() const;
120  inline G4double GetSumEnergyKilled() const;
121  inline void ResetKilledStatistics( G4int report = 1);
122  // Statistics for tracks killed (currently due to looping in field)
123 
124  inline void EnableShortStepOptimisation(G4bool optimise=true);
125  // Whether short steps < safety will avoid to call Navigator (if field=0)
126 
127  public: // without description
128 
130  const G4Track& ,
132  ) { return -1.0; };
133  // No operation in AtRestDoIt.
134 
136  const G4Track& ,
137  const G4Step&
138  ) {return 0;};
139  // No operation in AtRestDoIt.
140 
141  virtual void StartTracking(G4Track* aTrack);
142  // Reset state for new (potentially resumed) track
143 
144  protected:
145 
147  // Checks whether a field exists for the "global" field manager.
148 
149  private:
150 
152 
154 
157  // The Propagators used to transport the particle
158 
164  // G4bool fEnergyChanged;
167  // The particle's state after this Step, Store for DoIt
168 
170 
172 
174  // Flag to determine whether a boundary was reached.
175 
178  // Remember last safety origin & value.
179 
181  // New ParticleChange
182 
184 
185  // Thresholds for looping particles:
186  //
187  G4double fThreshold_Warning_Energy; // Warn above this energy
188  G4double fThreshold_Important_Energy; // Hesitate above this
189  G4int fThresholdTrials; // for this no of trials
190  // Above 'important' energy a 'looping' particle in field will
191  // *NOT* be abandoned, except after fThresholdTrials attempts.
192  // G4double fUnimportant_Energy;
193  // Below this energy, no verbosity for looping particles is issued
194 
195  // Counter for steps in which particle reports 'looping',
196  // if it is above 'Important' Energy
198  // Statistics for tracks abandoned
201 
202  // Whether to avoid calling G4Navigator for short step ( < safety)
203  // If using it, the safety estimate for endpoint will likely be smaller.
205 
206  G4SafetyHelper* fpSafetyHelper; // To pass it the safety value obtained
207 
208 };
209 
211 
212 #endif
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
CLHEP::Hep3Vector G4ThreeVector
void SetThresholdImportantEnergy(G4double newEnImp)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
virtual void StartTracking(G4Track *aTrack)
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
G4PropagatorInField * fFieldPropagator
G4ParticleChangeForTransport fParticleChange
void SetPropagatorInField(G4PropagatorInField *pFieldPropagator)
int G4int
Definition: G4Types.hh:78
Definition of the G4MonopoleFieldSetup class.
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
G4double GetThresholdWarningEnergy() const
G4int GetThresholdTrials() const
void SetThresholdWarningEnergy(G4double newEnWarn)
bool G4bool
Definition: G4Types.hh:79
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
Definition: G4Step.hh:76
G4PropagatorInField * GetPropagatorInField()
void SetThresholdTrials(G4int newMaxTrials)
void ResetKilledStatistics(G4int report=1)
G4MonopoleTransportation(const G4Monopole *p, G4int verbosityLevel=1)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
G4double GetThresholdImportantEnergy() const
void EnableShortStepOptimisation(G4bool optimise=true)
double G4double
Definition: G4Types.hh:76
G4double GetSumEnergyKilled() const
G4ForceCondition
G4MonopoleFieldSetup * fMagSetup
G4GPILSelection
G4double GetMaxEnergyKilled() const