Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4Transportation.hh
이 파일의 문서화 페이지로 가기
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: G4Transportation.hh 110805 2018-06-15 06:52:15Z gcosmo $
28 //
29 //
30 // ------------------------------------------------------------
31 // GEANT 4 include file implementation
32 // ------------------------------------------------------------
33 //
34 // Class description:
35 //
36 // G4Transportation is a process responsible for the transportation of
37 // a particle, i.e. the geometrical propagation encountering the
38 // geometrical sub-volumes of the detectors.
39 // It is also tasked with part of updating the "safety".
40 
41 // =======================================================================
42 // Created: 19 March 1997, J. Apostolakis
43 // =======================================================================
44 #ifndef G4Transportation_hh
45 #define G4Transportation_hh 1
46 
47 #include "G4VProcess.hh"
48 #include "G4FieldManager.hh"
49 
50 #include "G4Navigator.hh"
52 #include "G4PropagatorInField.hh"
53 #include "G4Track.hh"
54 #include "G4Step.hh"
56 class G4SafetyHelper;
59 
61 {
62  // Concrete class that does the geometrical transport
63 
64  public: // with description
65 
66  G4Transportation( G4int verbosityLevel= 1);
68 
70  const G4Track& track,
71  G4double previousStepSize,
72  G4double currentMinimumStep,
73  G4double& currentSafety,
74  G4GPILSelection* selection
75  );
76 
78  const G4Track& track,
79  const G4Step& stepData
80  );
81 
83  const G4Track& track,
84  const G4Step& stepData
85  );
86  // Responsible for the relocation
87 
89  const G4Track& ,
90  G4double previousStepSize,
91  G4ForceCondition* pForceCond
92  );
93  // Forces the PostStepDoIt action to be called,
94  // but does not limit the step
95 
97 
99  void SetPropagatorInField( G4PropagatorInField* pFieldPropagator);
100  // Access/set the assistant class that Propagate in a Field
101 
102  inline G4double GetThresholdWarningEnergy() const;
103  inline G4double GetThresholdImportantEnergy() const;
104  inline G4int GetThresholdTrials() const;
105 
106  inline void SetThresholdWarningEnergy( G4double newEnWarn );
107  inline void SetThresholdImportantEnergy( G4double newEnImp );
108  inline void SetThresholdTrials(G4int newMaxTrials );
109 
110  // Get/Set parameters for killing loopers:
111  // Above 'important' energy a 'looping' particle in field will
112  // *NOT* be abandoned, except after fThresholdTrials attempts.
113  // Below Warning energy, no verbosity for looping particles is issued
114 
115  inline G4double GetMaxEnergyKilled() const;
116  inline G4double GetSumEnergyKilled() const;
117  inline void ResetKilledStatistics( G4int report = 1);
118  // Statistics for tracks killed (currently due to looping in field)
119 
120  inline void EnableShortStepOptimisation(G4bool optimise=true);
121  // Whether short steps < safety will avoid to call Navigator (if field=0)
122 
123  static G4bool EnableUseMagneticMoment(G4bool useMoment=true);
124  // Whether to deflect particles with force due to magnetic moment
125 
126  public: // without description
127 
130  { return -1.0; } // No operation in AtRestGPIL
131 
133  { return 0; } // No operation in AtRestDoIt
134 
135  void StartTracking(G4Track* aTrack);
136  // Reset state for new (potentially resumed) track
137 
138  protected:
139 
141  // Checks whether a field exists for the "global" field manager.
142 
143  // void ReportLoopingTrack( ... ) --> moved to Logger class
144  // Warn about dropping of tracks that repeatedly exceed number of
145  // iterations for propagaton in field
146 
147  private:
148 
151  // The Propagators used to transport the particle
152 
160  // The particle's state after this Step, Store for DoIt
161 
163  G4bool fNewTrack; // Flag from StartTracking
165  G4bool fLastStepInVolume; // Last step - almost same as next flag
166  // (temporary redundancy for checking)
167  G4bool fGeometryLimitedStep; // Flag to determine whether a boundary was reached.
168 
169  G4bool fFieldExertedForce; // During current step
170 
172 
175  // Remember last safety origin & value.
176 
178  // New ParticleChange
179 
181 
182  // Thresholds for looping particles:
183  //
184  G4double fThreshold_Warning_Energy; // Warn above this energy
185  G4double fThreshold_Important_Energy; // Hesitate above this
186  G4int fThresholdTrials; // for this no of trials
187  // Above 'important' energy a 'looping' particle in field will
188  // *NOT* be abandoned, except after fThresholdTrials attempts.
189 
190  // Counter for steps in which particle reports 'looping',
191  // if it is above 'Important' Energy
192  //
194 
195  // Statistics for tracks abandoned
196  //
199 
200  // Whether to avoid calling G4Navigator for short step ( < safety)
201  // If using it, the safety estimate for endpoint will likely be smaller.
202  //
204 
205  G4SafetyHelper* fpSafetyHelper; // To pass it the safety value obtained
206  G4TransportationLogger* fpLogger; // Reports issues / raises warnings
207 
208  private:
209 
212 };
213 
214 #include "G4Transportation.icc"
215 
216 #endif
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
static G4bool EnableUseMagneticMoment(G4bool useMoment=true)
G4ThreeVector fPreviousSftOrigin
G4ParticleChangeForTransport fParticleChange
G4bool DoesGlobalFieldExist()
G4double fTransportEndKineticEnergy
void SetThresholdTrials(G4int newMaxTrials)
void SetPropagatorInField(G4PropagatorInField *pFieldPropagator)
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
void ResetKilledStatistics(G4int report=1)
void SetThresholdWarningEnergy(G4double newEnWarn)
G4ThreeVector fTransportEndSpin
G4double GetThresholdWarningEnergy() const
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
G4ThreeVector fTransportEndPosition
G4double GetSumEnergyKilled() const
G4GPILSelection
static G4bool fUseMagneticMoment
G4ThreeVector fTransportEndMomentumDir
Definition: G4Step.hh:76
G4Transportation(G4int verbosityLevel=1)
G4TouchableHandle fCurrentTouchableHandle
G4double fCandidateEndGlobalTime
void SetThresholdImportantEnergy(G4double newEnImp)
int G4int
Definition: G4Types.hh:78
G4PropagatorInField * fFieldPropagator
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
void EnableShortStepOptimisation(G4bool optimise=true)
G4double GetMaxEnergyKilled() const
G4ForceCondition
G4PropagatorInField * GetPropagatorInField()
G4double GetThresholdImportantEnergy() const
G4SafetyHelper * fpSafetyHelper
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
G4bool FieldExertedForce()
G4Navigator * fLinearNavigator
void StartTracking(G4Track *aTrack)
G4double fThreshold_Warning_Energy
G4int GetThresholdTrials() const
G4TransportationLogger * fpLogger
G4double fThreshold_Important_Energy
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)