Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4CoupledTransportation.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$
28//
29//
30// ------------------------------------------------------------
31// GEANT 4 include file implementation
32// ------------------------------------------------------------
33//
34// Class description:
35//
36// G4CoupledTransportation is an optional process to transport
37// a particle, in case of coupled navigation in parallel geometries
38// i.e. the geometrical propagation will be done
39// encountering the geometrical volumes of the detectors and
40// those of parallel geometries (eg for biasing, scoring, fast simulation)
41// It is tasked with updating the "safety" to reflect the geometrical
42// distance to the nearest volume, and the time of flight of the particle.
43
44// =======================================================================
45// Created: 17 May 2006, J. Apostolakis
46// =======================================================================
47#ifndef G4CoupledTransportation_hh
48#define G4CoupledTransportation_hh 1
49
50#include "G4VProcess.hh"
51
52#include "G4FieldManager.hh"
53
54#include "G4Navigator.hh"
57#include "G4PathFinder.hh"
58
59#include "G4Track.hh"
60#include "G4Step.hh"
62class G4SafetyHelper;
63
65{
66 // Concrete class that does the geometrical transport
67
68 public: // with description
69
70 G4CoupledTransportation( G4int verbosityLevel= 0);
72
74 const G4Track& track,
75 G4double previousStepSize,
76 G4double currentMinimumStep,
77 G4double& currentSafety,
78 G4GPILSelection* selection
79 );
80
82 const G4Track& track,
83 const G4Step& stepData
84 );
85
87 const G4Track& track,
88 const G4Step& stepData
89 );
90 // Responsible for the relocation.
91
93 const G4Track& ,
94 G4double previousStepSize,
95 G4ForceCondition* pForceCond
96 );
97 // Forces the PostStepDoIt action to be called,
98 // but does not limit the step.
99
102 // Access/set the assistant class that Propagate in a Field.
103
105 inline G4int GetVerboseLevel() const;
106 // Level of warnings regarding eg energy conservation
107 // in field integration.
108
111 inline G4int GetThresholdTrials() const;
112
113 inline void SetThresholdWarningEnergy( G4double newEnWarn );
114 inline void SetThresholdImportantEnergy( G4double newEnImp );
115 inline void SetThresholdTrials(G4int newMaxTrials );
116
117 // Get/Set parameters for killing loopers:
118 // Above 'important' energy a 'looping' particle in field will
119 // *NOT* be abandoned, except after fThresholdTrials attempts.
120 // Below Warning energy, no verbosity for looping particles is issued
121
124 inline void ResetKilledStatistics( G4int report = 1);
125 // Statistics for tracks killed (currently due to looping in field)
126
127 inline G4bool EnableUseMagneticMoment(G4bool useMoment=true);
128 // Whether to deflect particles with force due to magnetic moment
129
130 public: // without description
131
132 void StartTracking(G4Track* aTrack);
133 void EndTracking();
134
136 const G4Track& ,
138 ) { return -1.0; };
139 // No operation in AtRestDoIt.
140
142 const G4Track& ,
143 const G4Step&
144 ) {return 0;};
145 // No operation in AtRestDoIt.
146
147 protected:
148
150 // Checks whether a field exists for the "global" field manager.
151
152 void ReportInexactEnergy(G4double startEnergy, G4double endEnergy);
153 // Issue warning
154
155 private:
156
157 G4Navigator* fMassNavigator;
158 // The navigator for the 'mass' geometry (the real one, that physics occurs in)
159 G4PathFinder* fPathFinder;
160 G4int fNavigatorId;
161 // The PathFinder used to transport the particle
162
163 G4PropagatorInField* fFieldPropagator;
164 // Still required in order to find/set the fieldmanager
165
166 G4bool fGlobalFieldExists;
167 // G4bool fStartedNewTrack; // True for first step or restarted tracking
168 // until first step's AlongStepGPIL
169
170 G4ThreeVector fTransportEndPosition;
171 G4ThreeVector fTransportEndMomentumDir;
172 G4double fTransportEndKineticEnergy;
173 G4ThreeVector fTransportEndSpin;
174 G4bool fMomentumChanged;
175 G4bool fEnergyChanged;
176 // The particle's state after this Step, Store for DoIt
177
178 G4bool fEndGlobalTimeComputed;
179 G4double fCandidateEndGlobalTime;
180
181 G4bool fParticleIsLooping;
182
183 G4ThreeVector fPreviousSftOrigin;
184 G4double fPreviousMassSafety;
185 G4double fPreviousFullSafety;
186
187 G4TouchableHandle fCurrentTouchableHandle;
188
189 // G4bool fFieldExists;
190 // Whether a magnetic field exists ...
191 // A data member for this is problematic: it is useful only if it
192 // can be initialised and updated -- and a scheme is not yet possible.
193
194 G4bool fMassGeometryLimitedStep;
195 // Flag to determine whether a 'mass' boundary was reached.
196 G4bool fAnyGeometryLimitedStep;
197 // Did any geometry limit the step ?
198
199 G4ParticleChangeForTransport fParticleChange;
200 // New ParticleChange
201
202 G4double endpointDistance;
203
204
205 // Thresholds for looping particles:
206 //
207 G4double fThreshold_Warning_Energy; // Warn above this energy
208 G4double fThreshold_Important_Energy; // Hesitate above this
209 G4int fThresholdTrials; // for this no of trials
210 // Above 'important' energy a 'looping' particle in field will
211 // *NOT* be abandoned, except after fThresholdTrials attempts.
212 G4double fUnimportant_Energy;
213 // Below this energy, no verbosity for looping particles is issued
214
215 // Counter for steps in which particle reports 'looping',
216 // if it is above 'Important' Energy
217 G4int fNoLooperTrials;
218 // Statistics for tracks abandoned
219 G4double fSumEnergyKilled;
220 G4double fMaxEnergyKilled;
221
222 G4SafetyHelper* fpSafetyHelper; // To pass it the safety value obtained
223
224 // Whether to track state change from magnetic moment in a B-field
225 G4bool fUseMagneticMoment;
226
227 // Verbosity
228 G4int fVerboseLevel;
229 // Verbosity level for warnings
230 // eg about energy non-conservation in magnetic field.
231};
232
233#include "G4CoupledTransportation.icc"
234
235#endif
G4ForceCondition
G4GPILSelection
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4PropagatorInField * GetPropagatorInField()
G4bool EnableUseMagneticMoment(G4bool useMoment=true)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
G4int GetVerboseLevel() const
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
G4double GetMaxEnergyKilled() const
G4double GetThresholdImportantEnergy() const
void ResetKilledStatistics(G4int report=1)
void SetThresholdTrials(G4int newMaxTrials)
G4int GetThresholdTrials() const
void ReportInexactEnergy(G4double startEnergy, G4double endEnergy)
void SetVerboseLevel(G4int verboseLevel)
G4double GetSumEnergyKilled() const
G4double GetThresholdWarningEnergy() const
void SetPropagatorInField(G4PropagatorInField *pFieldPropagator)
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
void SetThresholdImportantEnergy(G4double newEnImp)
void SetThresholdWarningEnergy(G4double newEnWarn)
void StartTracking(G4Track *aTrack)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
Definition: G4Step.hh:78
G4int verboseLevel
Definition: G4VProcess.hh:368