Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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//
27/// \brief This class is a slightly modified version of G4Transportation
28/// initially written by John Apostolakis and colleagues (1997)
29/// But it should use the exact same algorithm
30//
31// Original Author : John Apostolakis
32//
33// Contact : Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
34//
35// WARNING : This class is released as a prototype.
36// It might strongly evolve or even disapear in the next releases.
37//
38// -------------------------------------------------------------------
39// Author: Mathieu Karamitros
40
41// The code is developed in the framework of the ESA AO7146
42//
43// We would be very happy hearing from you, send us your feedback! :)
44//
45// In order for Geant4-DNA to be maintained and still open-source,
46// article citations are crucial.
47// If you use Geant4-DNA chemistry and you publish papers about your software,
48// in addition to the general paper on Geant4-DNA:
49//
50// Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178
51//
52// we would be very happy if you could please also cite the following
53// reference papers on chemistry:
54//
55// J. Comput. Phys. 274 (2014) 841-882
56// Prog. Nucl. Sci. Tec. 2 (2011) 503-508
57
58#ifndef G4ITTransportation_H
59#define G4ITTransportation_H
60
62
63#include "G4VITProcess.hh"
64#include "G4Track.hh"
65#include "G4Step.hh"
67
68class G4ITNavigator;
69//class G4Navigator;
72
74{
75 // Concrete class that does the geometrical transport
76public:
77 // with description
78
79 G4ITTransportation(const G4String& aName = "ITTransportation",
80 G4int verbosityLevel = 0);
81 virtual ~G4ITTransportation();
82
84
86
87 virtual void BuildPhysicsTable(const G4ParticleDefinition&);
88
89 virtual void ComputeStep(const G4Track&,
90 const G4Step&,
91 const double timeStep,
92 double& spaceStep);
93
94 virtual void StartTracking(G4Track* aTrack);
95 // Give to the track a pointer to the transportation state
96
98 {
99 return GetState<G4ITTransportationState>()->fGeometryLimitedStep;
100 }
101
102 //________________________________________________________
103public:
104 // without description
105
108 {
109 return -1.0;
110 }
111 // No operation in AtRestDoIt.
112
114 {
115 return 0;
116 }
117 // No operation in AtRestDoIt.
118
120 G4double, // previousStepSize
121 G4double currentMinimumStep,
122 G4double& currentSafety,
123 G4GPILSelection* selection);
124
126 G4double, // previousStepSize
127 G4ForceCondition* pForceCond);
128
129 virtual G4VParticleChange* AlongStepDoIt(const G4Track& track,
130 const G4Step& stepData);
131
132 virtual G4VParticleChange* PostStepDoIt(const G4Track& track, const G4Step&);
133
134 //________________________________________________________
135 // inline virtual G4double GetTransportationTime() ;
136
139 // Access/set the assistant class that Propagate in a Field.
140
142 inline G4int GetVerboseLevel() const;
143 // Level of warnings regarding eg energy conservation
144 // in field integration.
145
148 inline G4int GetThresholdTrials() const;
149
150 inline void SetThresholdWarningEnergy(G4double newEnWarn);
152 inline void SetThresholdTrials(G4int newMaxTrials);
153
154 // Get/Set parameters for killing loopers:
155 // Above 'important' energy a 'looping' particle in field will
156 // *NOT* be abandoned, except after fThresholdTrials attempts.
157 // Below Warning energy, no verbosity for looping particles is issued
158
161 inline void ResetKilledStatistics(G4int report = 1);
162 // Statistics for tracks killed (currently due to looping in field)
163
164 inline void EnableShortStepOptimisation(G4bool optimise = true);
165 // Whether short steps < safety will avoid to call Navigator (if field=0)
166
167protected:
168 //________________________________________________________________
169 // Protected methods
171 // Checks whether a field exists for the "global" field manager.
172
173 //________________________________________________________________
174 // Process information
176 {
177 public:
179 virtual ~G4ITTransportationState();
181 {
182 return "G4ITTransportationState";
183 }
184
194 // The particle's state after this Step, Store for DoIt
195
198 // Flag to determine whether a boundary was reached.
199
202 // Remember last safety origin & value.
203
204 // Counter for steps in which particle reports 'looping',
205 // if it is above 'Important' Energy
207
208 // G4bool fFieldExists;
209 // Whether a magnetic field exists ...
210 // A data member for this is problematic: it is useful only if it
211 // can be initialised and updated -- and a scheme is not yet possible.
212
214 };
215
216 //________________________________________________________________
217 // Informations relative to the process only (meaning no information
218 // relative to the treated particle)
219 G4ITNavigator* fLinearNavigator;
221 // The Propagators used to transport the particle
222
223 // static const G4TouchableHandle nullTouchableHandle;
224 // Points to (G4VTouchable*) 0
225
227 // New ParticleChange
228
229 // Thresholds for looping particles:
230 //
231 G4double fThreshold_Warning_Energy; // Warn above this energy
232 G4double fThreshold_Important_Energy; // Hesitate above this
233 G4int fThresholdTrials; // for this no of trials
234 // Above 'important' energy a 'looping' particle in field will
235 // *NOT* be abandoned, except after fThresholdTrials attempts.
237 // Below this energy, no verbosity for looping particles is issued
238
239 // Statistics for tracks abandoned
242
243 // Whether to avoid calling G4Navigator for short step ( < safety)
244 // If using it, the safety estimate for endpoint will likely be smaller.
246
247 G4ITSafetyHelper* fpSafetyHelper; // To pass it the safety value obtained
248
249 // Verbosity
251 // Verbosity level for warnings
252 // eg about energy non-conservation in magnetic field.
253
255 {
256 fInstantiateProcessState = flag;
257 }
258
260 {
261 return fInstantiateProcessState;
262 }
263
264private:
265 G4bool fInstantiateProcessState;
266 G4ITTransportation& operator=(const G4ITTransportation&);
267};
268
269#include "G4ITTransportation.icc"
270#endif // G4ITTransportation_H
#define G4IT_ADD_CLONE(parent_class, kid_class)
Definition: AddClone_def.hh:52
G4ForceCondition
G4GPILSelection
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4ParticleChangeForTransport fParticleChange
void SetVerboseLevel(G4int verboseLevel)
G4double fThreshold_Important_Energy
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 void BuildPhysicsTable(const G4ParticleDefinition &)
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)
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)
G4ITSafetyHelper * fpSafetyHelper
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:62
G4int verboseLevel
Definition: G4VProcess.hh:360