Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/dna/management/include/G4ITTransportation.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 
 61 #include <CLHEP/Units/SystemOfUnits.h>
 62 
 63 #include "G4VITProcess.hh"
 64 #include "G4Track.hh"
 65 #include "G4Step.hh"
 66 #include "G4ParticleChangeForTransport.hh"
 67 
 68 class G4ITNavigator;
 69 //class G4Navigator;
 70 class G4ITSafetyHelper;
 71 class G4PropagatorInField;
 72 
 73 class G4ITTransportation : public G4VITProcess
 74 {
 75   // Concrete class that does the geometrical transport
 76 public:
 77   // with description
 78 
 79   G4ITTransportation(const G4String& aName = "ITTransportation",
 80                      G4int verbosityLevel = 0);
 81   ~G4ITTransportation() override;
 82 
 83   G4ITTransportation(const G4ITTransportation&);
 84 
 85   G4IT_ADD_CLONE(G4VITProcess, G4ITTransportation)
 86 
 87   void BuildPhysicsTable(const G4ParticleDefinition&) override;
 88 
 89   virtual void ComputeStep(const G4Track&,
 90                            const G4Step&,
 91                            const double timeStep,
 92                            double& spaceStep);
 93 
 94   void StartTracking(G4Track* aTrack) override;
 95   // Give to the track a pointer to the transportation state
 96 
 97   G4bool IsStepLimitedByGeometry()
 98   {
 99     return GetState<G4ITTransportationState>()->fGeometryLimitedStep;
100   }
101 
102   //________________________________________________________
103 public:
104   // without description
105 
106   G4double AtRestGetPhysicalInteractionLength(const G4Track&,
107                                                       G4ForceCondition*) override
108   {
109     return -1.0;
110   }
111   // No operation in  AtRestDoIt.
112 
113   G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&) override
114   {
115     return nullptr;
116   }
117   // No operation in  AtRestDoIt.
118 
119   G4double AlongStepGetPhysicalInteractionLength(const G4Track& track,
120                                                          G4double, //   previousStepSize
121                                                          G4double currentMinimumStep,
122                                                          G4double& currentSafety,
123                                                          G4GPILSelection* selection) override;
124 
125   G4double PostStepGetPhysicalInteractionLength(const G4Track&, // track
126                                                         G4double, // previousStepSize
127                                                         G4ForceCondition* pForceCond) override;
128 
129   G4VParticleChange* AlongStepDoIt(const G4Track& track,
130                                            const G4Step& stepData) override;
131 
132   G4VParticleChange* PostStepDoIt(const G4Track& track, const G4Step&) override;
133 
134   //________________________________________________________
135   //    inline virtual G4double GetTransportationTime() ;
136 
137   G4PropagatorInField* GetPropagatorInField();
138   void SetPropagatorInField(G4PropagatorInField* pFieldPropagator);
139   // Access/set the assistant class that Propagate in a Field.
140 
141   inline void SetVerboseLevel(G4int verboseLevel);
142   inline G4int GetVerboseLevel() const;
143   // Level of warnings regarding eg energy conservation
144   // in field integration.
145 
146   inline G4double GetThresholdWarningEnergy() const;
147   inline G4double GetThresholdImportantEnergy() const;
148   inline G4int GetThresholdTrials() const;
149 
150   inline void SetThresholdWarningEnergy(G4double newEnWarn);
151   inline void SetThresholdImportantEnergy(G4double newEnImp);
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 
159   inline G4double GetMaxEnergyKilled() const;
160   inline G4double GetSumEnergyKilled() const;
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 
167 protected:
168   //________________________________________________________________
169   // Protected methods
170   G4bool DoesGlobalFieldExist();
171   // Checks whether a field exists for the "global" field manager.
172 
173   //________________________________________________________________
174   // Process information
175   struct G4ITTransportationState : public G4ProcessState
176   {
177   public:
178     G4ITTransportationState();
179     ~G4ITTransportationState() override;
180     G4String GetType() override
181     {
182       return "G4ITTransportationState";
183     }
184 
185     G4ThreeVector fTransportEndPosition;
186     G4ThreeVector fTransportEndMomentumDir;
187     G4double fTransportEndKineticEnergy;
188     G4ThreeVector fTransportEndSpin;
189     G4bool fMomentumChanged;
190     G4bool fEnergyChanged;
191     G4bool fEndGlobalTimeComputed;
192     G4double fCandidateEndGlobalTime;
193     G4bool fParticleIsLooping;
194     // The particle's state after this Step, Store for DoIt
195 
196     G4TouchableHandle fCurrentTouchableHandle;
197     G4bool fGeometryLimitedStep;
198     // Flag to determine whether a boundary was reached.
199 
200     G4ThreeVector fPreviousSftOrigin;
201     G4double fPreviousSafety;
202     // Remember last safety origin & value.
203 
204     // Counter for steps in which particle reports 'looping',
205     //   if it is above 'Important' Energy
206     G4int fNoLooperTrials;
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 
213     G4double fEndPointDistance;
214   };
215 
216   //________________________________________________________________
217   // Informations relative to the process only (meaning no information
218   // relative to the treated particle)
219   G4ITNavigator* fLinearNavigator;
220   G4PropagatorInField* fFieldPropagator;
221   // The Propagators used to transport the particle
222 
223   // static const G4TouchableHandle nullTouchableHandle;
224   // Points to (G4VTouchable*) 0
225 
226   G4ParticleChangeForTransport fParticleChange;
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{10}; //    for this no of trials
234   // Above 'important' energy a 'looping' particle in field will
235   //   *NOT* be abandoned, except after fThresholdTrials attempts.
236   G4double fUnimportant_Energy;
237   //  Below this energy, no verbosity for looping particles is issued
238 
239   // Statistics for tracks abandoned
240   G4double fSumEnergyKilled{0.0};
241   G4double fMaxEnergyKilled{0.0};
242 
243   // Whether to avoid calling G4Navigator for short step ( < safety)
244   //   If using it, the safety estimate for endpoint will likely be smaller.
245   G4bool fShortStepOptimisation{false};
246 
247   G4ITSafetyHelper* fpSafetyHelper; // To pass it the safety value obtained
248 
249   // Verbosity
250   G4int fVerboseLevel;
251   // Verbosity level for warnings
252   // eg about energy non-conservation in magnetic field.
253 
254   void SetInstantiateProcessState(G4bool flag)
255   {
256     fInstantiateProcessState = flag;
257   }
258 
259   G4bool InstantiateProcessState()
260   {
261     return fInstantiateProcessState;
262   }
263 
264 private:
265   G4bool fInstantiateProcessState;
266   G4ITTransportation& operator=(const G4ITTransportation&);
267 };
268 
269 #include "G4ITTransportation.icc"
270 #endif // G4ITTransportation_H
271