Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
28//
29// ------------------------------------------------------------
30// GEANT 4 include file implementation
31// ------------------------------------------------------------
32//
33// Class description:
34//
35// G4CoupledTransportation is an optional process to transport
36// a particle, in case of coupled navigation in parallel geometries
37// i.e. the geometrical propagation will be done
38// encountering the geometrical volumes of the detectors and
39// those of parallel geometries (eg for biasing, scoring, fast simulation)
40// It is tasked with updating the "safety" to reflect the geometrical
41// distance to the nearest volume, and the time of flight of the particle.
42
43// =======================================================================
44// Created: 17 May 2006, J. Apostolakis
45// =======================================================================
46#ifndef G4CoupledTransportation_hh
47#define G4CoupledTransportation_hh 1
48
49#include "G4VProcess.hh"
50
51#include "G4FieldManager.hh"
52
53#include "G4Navigator.hh"
56#include "G4PathFinder.hh"
57
58#include "G4Track.hh"
59#include "G4Step.hh"
61
62class G4SafetyHelper;
64
66{
67 // Concrete class that does the geometrical transport
68
69 public: // with description
70
71 G4CoupledTransportation( G4int verbosityLevel= 0);
73
75 const G4Track& track,
76 G4double previousStepSize,
77 G4double currentMinimumStep,
78 G4double& currentSafety,
79 G4GPILSelection* selection
80 );
81
83 const G4Track& track,
84 const G4Step& stepData
85 );
86
88 const G4Track& track,
89 const G4Step& stepData
90 );
91 // Responsible for the relocation
92
94 const G4Track& ,
95 G4double previousStepSize,
96 G4ForceCondition* pForceCond
97 );
98 // Forces the PostStepDoIt action to be called,
99 // but does not limit the step
100
103 // Access/set the assistant class that Propagate in a Field
104
107 inline G4int GetThresholdTrials() const;
108
109 inline void SetThresholdWarningEnergy( G4double newEnWarn );
110 inline void SetThresholdImportantEnergy( G4double newEnImp );
111 inline void SetThresholdTrials(G4int newMaxTrials );
112
113 void SetHighLooperThresholds(); // Shortcut method - old values (meant for HEP)
114 void SetLowLooperThresholds(); // Set low thresholds - for low-E applications
115 void PushThresholdsToLogger(); // Inform logger of current thresholds
116 void ReportLooperThresholds(); // Print values of looper thresholds
117
118 // Get/Set parameters for killing loopers:
119 // Above 'important' energy a 'looping' particle in field will
120 // *NOT* be abandoned, except after fThresholdTrials attempts.
121 // Below Warning energy, no verbosity for looping particles is issued
122
125 inline void ResetKilledStatistics( G4int report = 1);
126 // Statistics for tracks killed (currently due to looping in field)
127
128 static G4bool EnableMagneticMoment(G4bool useMoment=true);
129 // Whether to deflect particles with force due to magnetic moment
130
131 static G4bool EnableGravity(G4bool useGravity);
132 // Turn on the capability to deflect particles with a gravity field
133
134 static void SetSilenceLooperWarnings( G4bool val);
135 // Do not warn (or throw exception) about 'looping' particles
137
138 static void SetSignifyStepsInAnyVolume( G4bool anyVol )
139 { fSignifyStepInAnyVolume = anyVol; }
141 { return fSignifyStepInAnyVolume; }
142 // Flag in step corresponds to first/last step in a volume 'any'
143 // geometry (if this is true) or refers to first/last step in mass
144 // geometry only (if false)
145
146 // The following methods give access to first/last step in particular
147 // geometry *independent* of the choice of the 'Signify' flag
148 //
149 G4bool IsFirstStepInAnyVolume() const { return fFirstStepInAnyVolume; }
150 G4bool IsLastStepInAnyVolume() const { return fAnyGeometryLimitedStep; }
151 G4bool IsFirstStepInMassVolume() const { return fFirstStepInMassVolume; }
152 G4bool IsLastStepInMassVolume() const { return fMassGeometryLimitedStep; }
153
154 public: // without description
155
156 void StartTracking(G4Track* aTrack);
157 void EndTracking();
158
159 static G4bool EnableUseMagneticMoment(G4bool useMoment=true)
160 { return EnableMagneticMoment(useMoment); }
161 // Old name ... obsolete
162
165 { return -1.0; } // No operation in AtRestGPIL
166
168 { return 0; } // No operation in AtRestDoIt
169
170 void PrintStatistics( std::ostream& outStr) const;
171
172 protected:
173
175 // Check whether any field exists in the geometry
176 // - replaces method that checked only whether a field for the world volume
177
178 void ReportInexactEnergy(G4double startEnergy, G4double endEnergy);
179 // Issue warning
180
181 void ReportMove( G4ThreeVector OldVector, G4ThreeVector NewVector,
182 const G4String& Quantity );
183
184 void ReportMissingLogger(const char * methodName);
185
186 private:
187
188 G4Navigator* fMassNavigator;
189 // The navigator for the 'mass' geometry
190 // (the real one, that physics occurs in)
191 G4PathFinder* fPathFinder;
192 G4int fNavigatorId;
193 // The PathFinder used to transport the particle
194
195 G4PropagatorInField* fFieldPropagator;
196 // Still required in order to find/set the fieldmanager
197
198 G4bool fAnyFieldExists;
199 // G4bool fStartedNewTrack; // True for first step or restarted tracking
200 // until first step's AlongStepGPIL
201
202 G4ThreeVector fTransportEndPosition;
203 G4ThreeVector fTransportEndMomentumDir;
204 G4double fTransportEndKineticEnergy;
205 G4ThreeVector fTransportEndSpin;
206 G4bool fMomentumChanged;
207 // The particle's state after this Step, Store for DoIt
208
209 G4bool fEndGlobalTimeComputed;
210 G4double fCandidateEndGlobalTime;
211
212 G4bool fParticleIsLooping;
213 G4bool fNewTrack;
214 G4ThreeVector fPreviousSftOrigin;
215 G4double fPreviousMassSafety;
216 G4double fPreviousFullSafety;
217
218 G4TouchableHandle fCurrentTouchableHandle;
219
220 // G4bool fFieldExists;
221 // Whether a magnetic field exists ...
222 // A data member for this is problematic: it is useful only if it
223 // can be initialised and updated -- and a scheme is not yet possible.
224
225 G4bool fMassGeometryLimitedStep;
226 // Flag to determine whether a 'mass' boundary was reached.
227 G4bool fAnyGeometryLimitedStep;
228 // Did any geometry limit the step ?
229
230 G4ParticleChangeForTransport fParticleChange;
231 // New ParticleChange
232
233 G4double fEndpointDistance;
234
235
236 // Thresholds for looping particles:
237 //
238 G4double fThreshold_Warning_Energy = 1.0 * CLHEP::keV; // Warn above this energy
239 G4double fThreshold_Important_Energy = 1.0 * CLHEP::MeV; // Give a few trial above this E
240 G4int fThresholdTrials = 10; // Number of trials an important looper survives
241 // Above 'important' energy a 'looping' particle in field will
242 // *NOT* be abandoned, except after fThresholdTrials attempts.
243
244 // Counter for steps in which particle reports 'looping',
245 // if it is above 'Important' Energy
246 //
247 G4int fNoLooperTrials=0;
248
249 // Statistics for tracks abandoned
250 //
251 G4double fSumEnergyKilled= 0.0;
252 G4double fSumEnerSqKilled= 0.0;
253 G4double fMaxEnergyKilled= -1.0;
254 G4int fMaxEnergyKilledPDG= 0;
255 unsigned long fNumLoopersKilled= 0;
256 G4double fSumEnergyKilled_NonElectron= 0.0;
257 G4double fSumEnerSqKilled_NonElectron= 0.0;
258 G4double fMaxEnergyKilled_NonElectron= -1.0;
259 G4int fMaxEnergyKilled_NonElecPDG= 0;
260 unsigned long fNumLoopersKilled_NonElectron= 0;
261 G4double fSumEnergySaved= 0.0;
262 G4double fMaxEnergySaved= -1.0;
263 G4double fSumEnergyUnstableSaved = 0.0;
264
265 G4SafetyHelper* fpSafetyHelper; // To pass it the safety value obtained
266 G4TransportationLogger* fpLogger; // Reports issues / raises warnings
267
268 private:
269
270 friend class G4Transportation;
271 static G4bool fUseMagneticMoment;
272 static G4bool fUseGravity;
273 static G4bool fSilenceLooperWarnings; // Flag to *Supress* all 'looper' warnings
274
275 private:
276
277 G4bool fFirstStepInMassVolume;
278 G4bool fFirstStepInAnyVolume;
279 // G4bool fLastStepInMassVolume; => use fMassGeometryLimitedStep
280 // G4bool fLastStepInAnyVolume; => use fAnyGeometryLimitedStep
281
282 static G4bool fSignifyStepInAnyVolume;
283 // True: First/Last step in any one of the geometries
284 // False: First/Last step in volume of 'mass' geometry
285};
286
287#include "G4CoupledTransportation.icc"
288
289#endif
G4ForceCondition
G4GPILSelection
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4PropagatorInField * GetPropagatorInField()
static void SetSignifyStepsInAnyVolume(G4bool anyVol)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
static G4bool GetSignifyStepsInAnyVolume()
static G4bool EnableUseMagneticMoment(G4bool useMoment=true)
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)
static G4bool EnableGravity(G4bool useGravity)
G4int GetThresholdTrials() const
void ReportMove(G4ThreeVector OldVector, G4ThreeVector NewVector, const G4String &Quantity)
static void SetSilenceLooperWarnings(G4bool val)
void ReportInexactEnergy(G4double startEnergy, G4double endEnergy)
static G4bool EnableMagneticMoment(G4bool useMoment=true)
G4double GetSumEnergyKilled() const
void PrintStatistics(std::ostream &outStr) 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)
void ReportMissingLogger(const char *methodName)
Definition: G4Step.hh:62