Geant4  10.02.p01
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: G4CoupledTransportation.hh 86964 2014-11-21 11:47:44Z gcosmo $
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"
56 #include "G4PropagatorInField.hh"
57 #include "G4PathFinder.hh"
58 
59 #include "G4Track.hh"
60 #include "G4Step.hh"
62 class 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 
101  void SetPropagatorInField( G4PropagatorInField* pFieldPropagator);
102  // Access/set the assistant class that Propagate in a Field.
103 
104  inline void SetVerboseLevel( G4int verboseLevel );
105  inline G4int GetVerboseLevel() const;
106  // Level of warnings regarding eg energy conservation
107  // in field integration.
108 
109  inline G4double GetThresholdWarningEnergy() const;
110  inline G4double GetThresholdImportantEnergy() const;
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 
122  inline G4double GetMaxEnergyKilled() const;
123  inline G4double GetSumEnergyKilled() const;
124  inline void ResetKilledStatistics( G4int report = 1);
125  // Statistics for tracks killed (currently due to looping in field)
126 
127  static 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  void ReportMove( G4ThreeVector OldVector, G4ThreeVector NewVector,
156  const G4String& Quantity );
157 
158  private:
159 
161  // The navigator for the 'mass' geometry (the real one, that physics occurs in)
164  // The PathFinder used to transport the particle
165 
167  // Still required in order to find/set the fieldmanager
168 
170  // G4bool fStartedNewTrack; // True for first step or restarted tracking
171  // until first step's AlongStepGPIL
172 
178  // The particle's state after this Step, Store for DoIt
179 
182 
184 
188 
190 
191  // G4bool fFieldExists;
192  // Whether a magnetic field exists ...
193  // A data member for this is problematic: it is useful only if it
194  // can be initialised and updated -- and a scheme is not yet possible.
195 
197  // Flag to determine whether a 'mass' boundary was reached.
199  // Did any geometry limit the step ?
200 
202  // New ParticleChange
203 
205 
206 
207  // Thresholds for looping particles:
208  //
209  G4double fThreshold_Warning_Energy; // Warn above this energy
210  G4double fThreshold_Important_Energy; // Hesitate above this
211  G4int fThresholdTrials; // for this no of trials
212  // Above 'important' energy a 'looping' particle in field will
213  // *NOT* be abandoned, except after fThresholdTrials attempts.
214 
215  // Counter for steps in which particle reports 'looping',
216  // if it is above 'Important' Energy
218  // Statistics for tracks abandoned
221 
222  G4SafetyHelper* fpSafetyHelper; // To pass it the safety value obtained
223 
224  // Verbosity
226  // Verbosity level for warnings
227  // eg about energy non-conservation in magnetic field.
228 
229  // Whether to track state change from magnetic moment in a B-field
230  private:
231  friend class G4Transportation;
233 
234 };
235 
237 
238 #endif
G4double GetThresholdImportantEnergy() const
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
void SetPropagatorInField(G4PropagatorInField *pFieldPropagator)
G4int verboseLevel
Definition: G4VProcess.hh:368
void ReportMove(G4ThreeVector OldVector, G4ThreeVector NewVector, const G4String &Quantity)
void StartTracking(G4Track *aTrack)
CLHEP::Hep3Vector G4ThreeVector
G4TouchableHandle fCurrentTouchableHandle
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
G4int GetVerboseLevel() const
G4double GetMaxEnergyKilled() const
int G4int
Definition: G4Types.hh:78
G4double GetSumEnergyKilled() const
static G4bool EnableUseMagneticMoment(G4bool useMoment=true)
void ReportInexactEnergy(G4double startEnergy, G4double endEnergy)
G4int GetThresholdTrials() const
void SetThresholdImportantEnergy(G4double newEnImp)
bool G4bool
Definition: G4Types.hh:79
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
void SetThresholdWarningEnergy(G4double newEnWarn)
Definition: G4Step.hh:76
void SetVerboseLevel(G4int verboseLevel)
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
void ResetKilledStatistics(G4int report=1)
G4double GetThresholdWarningEnergy() const
void SetThresholdTrials(G4int newMaxTrials)
G4PropagatorInField * GetPropagatorInField()
double G4double
Definition: G4Types.hh:76
G4ForceCondition
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
G4PropagatorInField * fFieldPropagator
G4GPILSelection
G4ParticleChangeForTransport fParticleChange
G4CoupledTransportation(G4int verbosityLevel=0)