Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4Transportation.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 //
26 //
27 // $Id: G4Transportation.hh 87829 2015-01-14 17:19:59Z gcosmo $
28 //
29 //
30 // ------------------------------------------------------------
31 // GEANT 4 include file implementation
32 // ------------------------------------------------------------
33 //
34 // Class description:
35 //
36 // G4Transportation is a process responsible for the transportation of
37 // a particle, i.e. the geometrical propagation encountering the
38 // geometrical sub-volumes of the detectors.
39 // It is also tasked with part of updating the "safety".
40 
41 // =======================================================================
42 // Created: 19 March 1997, J. Apostolakis
43 // =======================================================================
44 #ifndef G4Transportation_hh
45 #define G4Transportation_hh 1
46 
47 #include "G4VProcess.hh"
48 #include "G4FieldManager.hh"
49 
50 #include "G4Navigator.hh"
52 #include "G4PropagatorInField.hh"
53 #include "G4Track.hh"
54 #include "G4Step.hh"
56 class G4SafetyHelper;
58 
60 {
61  // Concrete class that does the geometrical transport
62 
63  public: // with description
64 
65  G4Transportation( G4int verbosityLevel= 1);
67 
69  const G4Track& track,
70  G4double previousStepSize,
71  G4double currentMinimumStep,
72  G4double& currentSafety,
73  G4GPILSelection* selection
74  );
75 
77  const G4Track& track,
78  const G4Step& stepData
79  );
80 
82  const G4Track& track,
83  const G4Step& stepData
84  );
85  // Responsible for the relocation.
86 
88  const G4Track& ,
89  G4double previousStepSize,
90  G4ForceCondition* pForceCond
91  );
92  // Forces the PostStepDoIt action to be called,
93  // but does not limit the step.
94 
95  G4bool FieldExertedForce() { return fFieldExertedForce; }
96 
98  void SetPropagatorInField( G4PropagatorInField* pFieldPropagator);
99  // Access/set the assistant class that Propagate in a Field.
100 
101  inline void SetVerboseLevel( G4int verboseLevel );
102  inline G4int GetVerboseLevel() const;
103  // Level of warnings regarding eg energy conservation
104  // in field integration.
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  static G4bool EnableUseMagneticMoment(G4bool useMoment=true);
128  // Whether to deflect particles with force due to magnetic moment
129 
130  public: // without description
131 
133  const G4Track& ,
135  ) { return -1.0; };
136  // No operation in AtRestDoIt.
137 
139  const G4Track& ,
140  const G4Step&
141  ) {return 0;};
142  // No operation in AtRestDoIt.
143 
144  void StartTracking(G4Track* aTrack);
145  // Reset state for new (potentially resumed) track
146 
147  protected:
148 
150  // Checks whether a field exists for the "global" field manager.
151 
152  private:
153 
154  G4Navigator* fLinearNavigator;
155  G4PropagatorInField* fFieldPropagator;
156  // The Propagators used to transport the particle
157 
158  // G4FieldManager* fGlobalFieldMgr; // Used MagneticField CC
159  // Field Manager for the whole Detector
160 
161  G4ThreeVector fTransportEndPosition;
162  G4ThreeVector fTransportEndMomentumDir;
163  G4double fTransportEndKineticEnergy;
164  G4ThreeVector fTransportEndSpin;
165  G4bool fMomentumChanged;
166  G4bool fEndGlobalTimeComputed;
167  G4double fCandidateEndGlobalTime;
168  // The particle's state after this Step, Store for DoIt
169 
170  G4bool fParticleIsLooping;
171  G4bool fNewTrack; // Flag from StartTracking
172  G4bool fFirstStepInVolume;
173  G4bool fLastStepInVolume; // Last step - almost same as next flag
174  // (temporary redundancy for checking)
175  G4bool fGeometryLimitedStep; // Flag to determine whether a boundary was reached.
176 
177  G4bool fFieldExertedForce; // During current step
178 
179  G4TouchableHandle fCurrentTouchableHandle;
180 
181  G4ThreeVector fPreviousSftOrigin;
182  G4double fPreviousSafety;
183  // Remember last safety origin & value.
184 
185  G4ParticleChangeForTransport fParticleChange;
186  // New ParticleChange
187 
188  G4double fEndPointDistance;
189 
190  // Thresholds for looping particles:
191  //
192  G4double fThreshold_Warning_Energy; // Warn above this energy
193  G4double fThreshold_Important_Energy; // Hesitate above this
194  G4int fThresholdTrials; // for this no of trials
195  // Above 'important' energy a 'looping' particle in field will
196  // *NOT* be abandoned, except after fThresholdTrials attempts.
197 
198  // Counter for steps in which particle reports 'looping',
199  // if it is above 'Important' Energy
200  G4int fNoLooperTrials;
201  // Statistics for tracks abandoned
202  G4double fSumEnergyKilled;
203  G4double fMaxEnergyKilled;
204 
205  // Whether to avoid calling G4Navigator for short step ( < safety)
206  // If using it, the safety estimate for endpoint will likely be smaller.
207  G4bool fShortStepOptimisation;
208 
209  G4SafetyHelper* fpSafetyHelper; // To pass it the safety value obtained
210 
211  // Verbosity
212  G4int fVerboseLevel;
213  // Verbosity level for warnings
214  // eg about energy non-conservation in magnetic field.
215 
216  // Whether to track state change from magnetic moment in a B-field
217 
218  private:
220  static G4bool fUseMagneticMoment;
221 
222 };
223 
224 #include "G4Transportation.icc"
225 
226 #endif
void SetPropagatorInField(G4PropagatorInField *pFieldPropagator)
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
G4double GetThresholdWarningEnergy() const
G4int verboseLevel
Definition: G4VProcess.hh:368
void StartTracking(G4Track *aTrack)
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
G4PropagatorInField * GetPropagatorInField()
int G4int
Definition: G4Types.hh:78
void SetThresholdImportantEnergy(G4double newEnImp)
void SetThresholdWarningEnergy(G4double newEnWarn)
void SetThresholdTrials(G4int newMaxTrials)
bool G4bool
Definition: G4Types.hh:79
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
G4bool DoesGlobalFieldExist()
Definition: G4Step.hh:76
G4bool FieldExertedForce()
G4Transportation(G4int verbosityLevel=1)
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
G4int GetVerboseLevel() const
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
G4int GetThresholdTrials() const
void ResetKilledStatistics(G4int report=1)
double G4double
Definition: G4Types.hh:76
G4ForceCondition
void EnableShortStepOptimisation(G4bool optimise=true)
void SetVerboseLevel(G4int verboseLevel)
G4GPILSelection
G4double GetSumEnergyKilled() const
G4double GetThresholdImportantEnergy() const
static G4bool EnableUseMagneticMoment(G4bool useMoment=true)
G4double GetMaxEnergyKilled() const