Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
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//
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
48#include "G4Track.hh"
49#include "G4Step.hh"
51
52class G4Navigator;
54class G4SafetyHelper;
56
58{
59 // Concrete class that does the geometrical transport
60
61 public: // with description
62
63 G4Transportation( G4int verbosityLevel= 1, const G4String& aName = "Transportation");
65
67 const G4Track& track,
68 G4double previousStepSize,
69 G4double currentMinimumStep,
70 G4double& currentSafety,
71 G4GPILSelection* selection
72 ); // override;
73
75 const G4Track& track,
76 const G4Step& stepData
77 ); // override;
78
80 const G4Track& track,
81 const G4Step& stepData
82 ); // override;
83 // Responsible for the relocation
84
86 const G4Track& ,
87 G4double previousStepSize,
88 G4ForceCondition* pForceCond
89 ); // override;
90 // Forces the PostStepDoIt action to be called,
91 // but does not limit the step
92
94
96 void SetPropagatorInField( G4PropagatorInField* pFieldPropagator);
97 // Access/set the assistant class that Propagate in a Field
98
101 inline G4int GetThresholdTrials() const;
102
103 inline void SetThresholdWarningEnergy( G4double newEnWarn );
104 inline void SetThresholdImportantEnergy( G4double newEnImp );
105 inline void SetThresholdTrials(G4int newMaxTrials );
106 // Get/Set parameters for killing loopers:
107 // Above 'important' energy a 'looping' particle in field will
108 // *NOT* be abandoned, except after fThresholdTrials attempts.
109 // Below Warning energy, no verbosity for looping particles is issued
110
111 void SetHighLooperThresholds(); // Shortcut method - old values (meant for HEP)
112 void SetLowLooperThresholds(); // Set low thresholds - for low-E applications
113 void PushThresholdsToLogger(); // Inform logger of current thresholds
114 void ReportLooperThresholds(); // Print values of looper thresholds
115
118 inline void ResetKilledStatistics( G4int report = 1);
119 // Statistics for tracks killed (currently due to looping in field)
120
121 inline void EnableShortStepOptimisation(G4bool optimise=true);
122 // Whether short steps < safety will avoid to call Navigator (if field=0)
123
124 static G4bool EnableMagneticMoment(G4bool useMoment=true);
125 // Whether to enable particles to be deflected with force due to magnetic moment
126
127 static G4bool EnableGravity(G4bool useGravity=true);
128 // Whether to enable particles to be deflected with force due to gravity
129
130 static void SetSilenceLooperWarnings( G4bool val);
131 // Do not warn (or throw exception) about 'looping' particles
133
134 public: // without description
135 static G4bool EnableUseMagneticMoment(G4bool useMoment=true)
136 { return EnableMagneticMoment(useMoment); } // Old name - will be deprecated
137
138 public: // without description
139
142 { return -1.0; } // No operation in AtRestGPIL
143
145 { return 0; } // No operation in AtRestDoIt
146
147 void StartTracking(G4Track* aTrack);
148 // Reset state for new (potentially resumed) track
149
150 virtual void ProcessDescription(std::ostream& outFile) const; // override;
151 void PrintStatistics( std::ostream& outStr) const;
152
153 protected:
154
155 void SetTouchableInformation(const G4TouchableHandle& touchable);
156
157 void ReportMissingLogger(const char * methodName);
158
159 protected:
160
162 // The navigator for the 'mass' geometry
163 // (the real one, that physics occurs in)
165 // The Propagators used to transport the particle
166
174 // The particle's state after this Step, Store for DoIt
175
177
179 G4bool fNewTrack= true; // Flag from StartTracking
181 G4bool fLastStepInVolume= false; // Last step - almost same as next flag
182 // (temporary redundancy for checking)
184 // Flag to determine whether a boundary was reached
185
186 G4bool fFieldExertedForce= false; // During current step
187
189
192 // Remember last safety origin & value.
193
195 // New ParticleChange
196
198
199 // Thresholds for looping particles:
200 //
201 G4double fThreshold_Warning_Energy = 1.0 * CLHEP::keV; // Warn above this energy
202 G4double fThreshold_Important_Energy = 1.0 * CLHEP::MeV; // Give a few trial above this E
203 G4int fThresholdTrials = 10; // Number of trials an important looper survives
204 // Above 'important' energy a 'looping' particle in field will
205 // *NOT* be abandoned, except after fThresholdTrials attempts.
206 G4int fAbandonUnstableTrials = 0; // Number of trials after which to abandon
207 // unstable loopers ( 0 = never )
208 // Counter for steps in which particle reports 'looping',
209 // ( Used if it is above 'Important' Energy. )
211
212 // Statistics for tracks abandoned due to looping - and 'saved' despite looping
213 //
218 unsigned long fNumLoopersKilled= 0;
227 // Whether to avoid calling G4Navigator for short step ( < safety)
228 // If using it, the safety estimate for endpoint will likely be smaller.
229 //
231
232 G4SafetyHelper* fpSafetyHelper; // To pass it the safety value obtained
233 G4TransportationLogger* fpLogger; // Reports issues / raises warnings
234
235 protected:
236
239 static G4bool fSilenceLooperWarnings; // Flag to *Supress* all 'looper' warnings
240
241};
242
243#include "G4Transportation.icc"
244
245#endif
G4ForceCondition
G4GPILSelection
CLHEP::Hep3Vector G4ThreeVector
G4ReferenceCountedHandle< G4VTouchable > G4TouchableHandle
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void SetThresholdImportantEnergy(G4double newEnImp)
G4ThreeVector fPreviousSftOrigin
G4PropagatorInField * GetPropagatorInField()
static void SetSilenceLooperWarnings(G4bool val)
G4double fThreshold_Important_Energy
G4Transportation(G4int verbosityLevel=1, const G4String &aName="Transportation")
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
G4double GetThresholdImportantEnergy() const
void SetTouchableInformation(const G4TouchableHandle &touchable)
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 *)
static G4bool fUseGravity
G4double GetThresholdWarningEnergy() const
G4Navigator * fLinearNavigator
G4TransportationLogger * fpLogger
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
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition G4VProcess.cc:46