Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ITTransportation.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// $Id: G4ITTransportation.hh 64374 2012-10-31 16:37:23Z gcosmo $
27//
28/// \brief This class is a slightly modified version of G4Transportation
29/// initially written by John Apostolakis and colleagues (1997)
30/// But it should use the exact same algorithm
31//
32// Original Author : John Apostolakis
33//
34// Contact : Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
35//
36// WARNING : This class is released as a prototype.
37// It might strongly evolve or even disapear in the next releases.
38//
39// -------------------------------------------------------------------
40
41#ifndef G4ITTransportation_H
42#define G4ITTransportation_H
43
45
46#include "G4IT.hh"
47#include "G4VITProcess.hh"
48#include "G4Track.hh"
49#include "G4Step.hh"
51
52class G4ITNavigator;
53//class G4Navigator;
54class G4SafetyHelper;
56
58{
59 // Concrete class that does the geometrical transport
60public: // with description
61
62 G4ITTransportation(const G4String& aName = "ITTransportation",G4int verbosityLevel= 1);
63 virtual ~G4ITTransportation();
64
66
68
70
71 virtual void ComputeStep(const G4Track&,
72 const G4Step&,
73 const double timeStep,
74 double& spaceStep) ;
75
76 virtual void StartTracking(G4Track* aTrack);
77 // Give to the track a pointer to the transportation state
78
79//________________________________________________________
80public: // without description
81
83 const G4Track& ,
85 )
86 {
87 return -1.0;
88 }
89 // No operation in AtRestDoIt.
90
92 const G4Track& ,
93 const G4Step&
94 )
95 {
96 return 0;
97 }
98 // No operation in AtRestDoIt.
99
101 const G4Track& track,
102 G4double , // previousStepSize
103 G4double currentMinimumStep,
104 G4double& currentSafety,
105 G4GPILSelection* selection);
106
108 const G4Track& , // track
109 G4double , // previousStepSize
110 G4ForceCondition* pForceCond);
111
112 virtual G4VParticleChange* AlongStepDoIt( const G4Track& track,
113 const G4Step& stepData );
114
115 virtual G4VParticleChange* PostStepDoIt( const G4Track& track,
116 const G4Step& ) ;
117
118//________________________________________________________
119// inline virtual G4double GetTransportationTime() ;
120
123 // Access/set the assistant class that Propagate in a Field.
124
126 inline G4int GetVerboseLevel() const;
127 // Level of warnings regarding eg energy conservation
128 // in field integration.
129
132 inline G4int GetThresholdTrials() const;
133
134 inline void SetThresholdWarningEnergy( G4double newEnWarn );
135 inline void SetThresholdImportantEnergy( G4double newEnImp );
136 inline void SetThresholdTrials(G4int newMaxTrials );
137
138 // Get/Set parameters for killing loopers:
139 // Above 'important' energy a 'looping' particle in field will
140 // *NOT* be abandoned, except after fThresholdTrials attempts.
141 // Below Warning energy, no verbosity for looping particles is issued
142
145 inline void ResetKilledStatistics( G4int report = 1);
146 // Statistics for tracks killed (currently due to looping in field)
147
148 inline void EnableShortStepOptimisation(G4bool optimise=true);
149 // Whether short steps < safety will avoid to call Navigator (if field=0)
150
151protected :
152 //________________________________________________________________
153 // Protected methods
155 // Checks whether a field exists for the "global" field manager.
156
157 //________________________________________________________________
158 // Process information
160 {
161 public :
163 virtual ~G4ITTransportationState();
164
174 // The particle's state after this Step, Store for DoIt
175
178 // Flag to determine whether a boundary was reached.
179
182 // Remember last safety origin & value.
183
184 // Counter for steps in which particle reports 'looping',
185 // if it is above 'Important' Energy
187
188 // G4bool fFieldExists;
189 // Whether a magnetic field exists ...
190 // A data member for this is problematic: it is useful only if it
191 // can be initialised and updated -- and a scheme is not yet possible.
192
194 };
195
197
198 //________________________________________________________________
199 // Informations relative to the process only (meaning no information
200 // relative to the treated particle)
203 // The Propagators used to transport the particle
204
205// static const G4TouchableHandle nullTouchableHandle; // Points to (G4VTouchable*) 0
206
208 // New ParticleChange
209
210 // Thresholds for looping particles:
211 //
212 G4double fThreshold_Warning_Energy; // Warn above this energy
213 G4double fThreshold_Important_Energy; // Hesitate above this
214 G4int fThresholdTrials; // for this no of trials
215 // Above 'important' energy a 'looping' particle in field will
216 // *NOT* be abandoned, except after fThresholdTrials attempts.
218 // Below this energy, no verbosity for looping particles is issued
219
220
221 // Statistics for tracks abandoned
224
225 // Whether to avoid calling G4Navigator for short step ( < safety)
226 // If using it, the safety estimate for endpoint will likely be smaller.
228
229 G4SafetyHelper* fpSafetyHelper; // To pass it the safety value obtained
230
231 // Verbosity
233 // Verbosity level for warnings
234 // eg about energy non-conservation in magnetic field.
235
237 { fInstantiateProcessState = flag; }
238
239 G4bool InstantiateProcessState() { return fInstantiateProcessState; }
240
241private :
242 G4bool fInstantiateProcessState;
243 G4ITTransportation& operator=(const G4ITTransportation&);
244};
245
246#include "G4ITTransportation.icc"
247#endif // G4ITTransportation_H
#define G4IT_ADD_CLONE(parent_class, kid_class)
Definition: AddClone_def.hh:45
G4ForceCondition
G4GPILSelection
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4ParticleChangeForTransport fParticleChange
void SetVerboseLevel(G4int verboseLevel)
G4double fThreshold_Important_Energy
G4SafetyHelper * fpSafetyHelper
G4double GetMaxEnergyKilled() const
void SetThresholdWarningEnergy(G4double newEnWarn)
G4int GetVerboseLevel() const
virtual void ComputeStep(const G4Track &, const G4Step &, const double timeStep, double &spaceStep)
void EnableShortStepOptimisation(G4bool optimise=true)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double, G4ForceCondition *pForceCond)
void SetInstantiateProcessState(G4bool flag)
G4double GetThresholdImportantEnergy() const
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
G4PropagatorInField * fFieldPropagator
void SetThresholdTrials(G4int newMaxTrials)
void SetPropagatorInField(G4PropagatorInField *pFieldPropagator)
void SetThresholdImportantEnergy(G4double newEnImp)
G4ITTransportationState *const & fTransportationState
G4double GetThresholdWarningEnergy() const
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
G4double GetSumEnergyKilled() const
G4ITNavigator * fLinearNavigator
G4int GetThresholdTrials() const
G4PropagatorInField * GetPropagatorInField()
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &)
virtual void StartTracking(G4Track *aTrack)
void ResetKilledStatistics(G4int report=1)
Definition: G4Step.hh:78
G4int verboseLevel
Definition: G4VProcess.hh:368
#define const
Definition: zconf.h:118