Geant4-11
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//
28//
29// ------------------------------------------------------------
30// GEANT 4 include file implementation
31// ------------------------------------------------------------
32//
33// Class description:
34//
35// G4Transportation is a process responsible for the transportation of
36// a particle, i.e. the geometrical propagation encountering the
37// geometrical sub-volumes of the detectors.
38// It is also tasked with part of updating the "safety".
39
40// =======================================================================
41// Created: 19 March 1997, J. Apostolakis
42// =======================================================================
43#ifndef G4Transportation_hh
44#define G4Transportation_hh 1
45
46#include "G4VProcess.hh"
47#include "G4FieldManager.hh"
48
49#include "G4Navigator.hh"
52#include "G4Track.hh"
53#include "G4Step.hh"
55
56class G4SafetyHelper;
59
61{
62 // Concrete class that does the geometrical transport
63
64 public: // with description
65
68
70 const G4Track& track,
71 G4double previousStepSize,
72 G4double currentMinimumStep,
73 G4double& currentSafety,
74 G4GPILSelection* selection
75 ); // override;
76
78 const G4Track& track,
79 const G4Step& stepData
80 ); // override;
81
83 const G4Track& track,
84 const G4Step& stepData
85 ); // override;
86 // Responsible for the relocation
87
89 const G4Track& ,
90 G4double previousStepSize,
91 G4ForceCondition* pForceCond
92 ); // override;
93 // Forces the PostStepDoIt action to be called,
94 // but does not limit the step
95
97
99 void SetPropagatorInField( G4PropagatorInField* pFieldPropagator);
100 // Access/set the assistant class that Propagate in a Field
101
104 inline G4int GetThresholdTrials() const;
105
106 inline void SetThresholdWarningEnergy( G4double newEnWarn );
107 inline void SetThresholdImportantEnergy( G4double newEnImp );
108 inline void SetThresholdTrials(G4int newMaxTrials );
109 // Get/Set parameters for killing loopers:
110 // Above 'important' energy a 'looping' particle in field will
111 // *NOT* be abandoned, except after fThresholdTrials attempts.
112 // Below Warning energy, no verbosity for looping particles is issued
113
114 void SetHighLooperThresholds(); // Shortcut method - old values (meant for HEP)
115 void SetLowLooperThresholds(); // Set low thresholds - for low-E applications
116 void PushThresholdsToLogger(); // Inform logger of current thresholds
117 void ReportLooperThresholds(); // Print values of looper thresholds
118
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 EnableMagneticMoment(G4bool useMoment=true);
128 // Whether to enable particles to be deflected with force due to magnetic moment
129
130 static G4bool EnableGravity(G4bool useGravity=true);
131 // Whether to enable particles to be deflected with force due to gravity
132
133 static void SetSilenceLooperWarnings( G4bool val);
134 // Do not warn (or throw exception) about 'looping' particles
136
137 public: // without description
138 static G4bool EnableUseMagneticMoment(G4bool useMoment=true)
139 { return EnableMagneticMoment(useMoment); } // Old name - will be deprecated
140
141 public: // without description
142
145 { return -1.0; } // No operation in AtRestGPIL
146
148 { return 0; } // No operation in AtRestDoIt
149
150 void StartTracking(G4Track* aTrack);
151 // Reset state for new (potentially resumed) track
152
153 virtual void ProcessDescription(std::ostream& outFile) const; // override;
154 void PrintStatistics( std::ostream& outStr) const;
155
156 protected:
157
159 // Check whether any field exists in the geometry
160 // - replaces method that checked only whether a field for the world volume
161
162 void ReportMissingLogger(const char * methodName);
163
164 private:
165
168 // The Propagators used to transport the particle
169
177 // The particle's state after this Step, Store for DoIt
178
180
182 G4bool fNewTrack= true; // Flag from StartTracking
184 G4bool fLastStepInVolume= false; // Last step - almost same as next flag
185 // (temporary redundancy for checking)
187 // Flag to determine whether a boundary was reached
188
189 G4bool fFieldExertedForce= false; // During current step
190
192
195 // Remember last safety origin & value.
196
198 // New ParticleChange
199
201
202 // Thresholds for looping particles:
203 //
204 G4double fThreshold_Warning_Energy = 1.0 * CLHEP::keV; // Warn above this energy
205 G4double fThreshold_Important_Energy = 1.0 * CLHEP::MeV; // Give a few trial above this E
206 G4int fThresholdTrials = 10; // Number of trials an important looper survives
207 // Above 'important' energy a 'looping' particle in field will
208 // *NOT* be abandoned, except after fThresholdTrials attempts.
209 G4int fAbandonUnstableTrials = 0; // Number of trials after which to abandon
210 // unstable loopers ( 0 = never )
211 // Counter for steps in which particle reports 'looping',
212 // ( Used if it is above 'Important' Energy. )
214
215 // Statistics for tracks abandoned due to looping - and 'saved' despite looping
216 //
221 unsigned long fNumLoopersKilled= 0;
230 // Whether to avoid calling G4Navigator for short step ( < safety)
231 // If using it, the safety estimate for endpoint will likely be smaller.
232 //
234
235 G4SafetyHelper* fpSafetyHelper; // To pass it the safety value obtained
236 G4TransportationLogger* fpLogger; // Reports issues / raises warnings
237
238 private:
239
243 static G4bool fSilenceLooperWarnings; // Flag to *Supress* all 'looper' warnings
244
245};
246
247#include "G4Transportation.icc"
248
249#endif
G4ForceCondition
G4GPILSelection
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
Definition: G4Step.hh:62
void SetThresholdImportantEnergy(G4double newEnImp)
G4ThreeVector fPreviousSftOrigin
G4PropagatorInField * GetPropagatorInField()
static void SetSilenceLooperWarnings(G4bool val)
G4double fThreshold_Important_Energy
G4double fMaxEnergyKilled_NonElectron
G4ThreeVector fTransportEndSpin
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
static G4bool EnableGravity(G4bool useGravity=true)
G4double fTransportEndKineticEnergy
void PushThresholdsToLogger()
void PrintStatistics(std::ostream &outStr) const
void EnableShortStepOptimisation(G4bool optimise=true)
G4int GetThresholdTrials() const
unsigned long fNumLoopersKilled_NonElectron
G4bool FieldExertedForce()
G4double GetThresholdImportantEnergy() const
static G4bool fSilenceLooperWarnings
G4PropagatorInField * fFieldPropagator
void SetPropagatorInField(G4PropagatorInField *pFieldPropagator)
G4ThreeVector fTransportEndPosition
G4double fThreshold_Warning_Energy
G4double fSumEnerSqKilled_NonElectron
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
void SetThresholdTrials(G4int newMaxTrials)
G4double GetSumEnergyKilled() const
G4ParticleChangeForTransport fParticleChange
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
static G4bool EnableUseMagneticMoment(G4bool useMoment=true)
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
G4double fSumEnergyUnstableSaved
static G4bool EnableMagneticMoment(G4bool useMoment=true)
void ResetKilledStatistics(G4int report=1)
G4double GetMaxEnergyKilled() const
G4TouchableHandle fCurrentTouchableHandle
void StartTracking(G4Track *aTrack)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
G4Transportation(G4int verbosityLevel=1)
static G4bool fUseGravity
G4double GetThresholdWarningEnergy() const
G4Navigator * fLinearNavigator
G4TransportationLogger * fpLogger
G4bool DoesAnyFieldExist()
static G4bool fUseMagneticMoment
void SetThresholdWarningEnergy(G4double newEnWarn)
G4double fSumEnergyKilled_NonElectron
static G4bool GetSilenceLooperWarnings()
void ReportMissingLogger(const char *methodName)
G4ThreeVector fTransportEndMomentumDir
unsigned long fNumLoopersKilled
G4SafetyHelper * fpSafetyHelper
G4double fCandidateEndGlobalTime
virtual void ProcessDescription(std::ostream &outFile) const
static constexpr double keV
static constexpr double MeV