Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/dna/management/include/G4VITProcess.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 // Author: Mathieu Karamitros
 28 
 29 // The code is developed in the framework of the ESA AO7146
 30 //
 31 // We would be very happy hearing from you, send us your feedback! :)
 32 //
 33 // In order for Geant4-DNA to be maintained and still open-source,
 34 // article citations are crucial. 
 35 // If you use Geant4-DNA chemistry and you publish papers about your software, 
 36 // in addition to the general paper on Geant4-DNA:
 37 //
 38 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178
 39 //
 40 // we would be very happy if you could please also cite the following
 41 // reference papers on chemistry:
 42 //
 43 // J. Comput. Phys. 274 (2014) 841-882
 44 // Prog. Nucl. Sci. Tec. 2 (2011) 503-508 
 45 
 46 #ifndef G4VITProcess_H
 47 #define G4VITProcess_H
 48 
 49 #include <G4VProcess.hh>
 50 #include "AddClone_def.hh"
 51 #include "G4ReferenceCast.hh"
 52 #include "G4memory.hh"
 53 #include <typeinfo>
 54 
 55 class G4IT;
 56 class G4TrackingInformation;
 57 
 58 struct G4ProcessState_Lock
 59 {
 60   inline virtual ~G4ProcessState_Lock()
 61   {
 62     ;
 63   }
 64 };
 65 
 66 /*
 67  class G4ProcessStateHandle_Lock : public G4shared_ptr<G4ProcessState_Lock>
 68  {
 69  public:
 70  G4ProcessStateHandle_Lock(G4ProcessState_Lock* plock) : G4shared_ptr<G4ProcessState_Lock>(plock)
 71  {}
 72  virtual ~G4ProcessStateHandle_Lock(){}
 73  };
 74  */
 75 
 76 #define InitProcessState(destinationType,source) \
 77     reference_cast<destinationType>(source)
 78 
 79 #define DowncastProcessState(destinationType) \
 80   G4dynamic_pointer_cast<destinationType>(G4VITProcess::fpState)
 81 
 82 #define UpcastProcessState(destinationType) \
 83   G4dynamic_pointer_cast<destinationType>(G4VITProcess::fpState)
 84 
 85 #define DowncastState(destinationType,source) \
 86   G4dynamic_pointer_cast<destinationType>(source)
 87 
 88 #define UpcastState(destinationType,source) \
 89   G4dynamic_pointer_cast<destinationType>(source)
 90 
 91 /**
 92  * G4VITProcess inherits from G4VProcess.
 93  * A G4VITProcess is able to save its current state for a given track into G4IT.
 94  * This state may be retrieve latter on to be used by the G4VITProcess.
 95  * Each G4VITProcess is tagged.
 96  */
 97 
 98 class G4VITProcess : public G4VProcess
 99 {
100 public:
101   //__________________________________
102   // Constructors & destructors
103   G4VITProcess(const G4String& name, G4ProcessType type = fNotDefined);
104 
105   ~G4VITProcess() override;
106   G4VITProcess(const G4VITProcess& other);
107   G4VITProcess& operator=(const G4VITProcess& other);
108 
109   // equal opperators
110   G4bool operator==(const G4VITProcess &right) const;
111   G4bool operator!=(const G4VITProcess &right) const;
112 
113   G4IT_TO_BE_CLONED(G4VITProcess)
114 
115   size_t GetProcessID() const
116   {
117     return fProcessID;
118   }
119 
120   G4shared_ptr<G4ProcessState_Lock> GetProcessState()
121   {
122     return UpcastProcessState(G4ProcessState_Lock);
123   }
124 
125   void SetProcessState(G4shared_ptr<G4ProcessState_Lock> aProcInfo)
126   {
127     fpState = DowncastState(G4ProcessState, aProcInfo);
128   }
129 
130   void ResetProcessState()
131   {
132     fpState.reset();
133   }
134 
135   //__________________________________
136   // Initialize and Save process info
137 
138   void StartTracking(G4Track*) override;
139 
140   void BuildPhysicsTable(const G4ParticleDefinition&) override
141   {
142   }
143 
144   inline G4double GetInteractionTimeLeft();
145 
146   /** WARNING : Redefine the method of G4VProcess
147    * reset (determine the value of)NumberOfInteractionLengthLeft
148    */
149   void ResetNumberOfInteractionLengthLeft() override;
150 
151   inline G4bool ProposesTimeStep() const;
152 
153   inline static const size_t& GetMaxProcessIndex();
154 
155 protected:
156   // with description
157 
158   void RetrieveProcessInfo();
159   void CreateInfo();
160 
161   //__________________________________
162   // Process info
163   // friend class G4TrackingInformation ;
164 
165   struct G4ProcessState : public G4ProcessState_Lock
166   {
167   public:
168     G4ProcessState();
169     ~G4ProcessState() override;
170 
171     virtual G4String GetType()
172     {
173       return "G4ProcessState";
174     }
175 
176     G4double theNumberOfInteractionLengthLeft;
177     // The flight length left for the current tracking particle
178     // in unit of "Interaction length".
179 
180     G4double theInteractionTimeLeft;
181     // Time left before the interaction : for at rest processes
182 
183     G4double currentInteractionLength;
184     // The InteractionLength in the current material
185 
186     template<typename T>
187       T* GetState()
188       {
189         return dynamic_cast<T*>(this);
190       }
191   };
192 
193   template<typename T>
194     class G4ProcessStateBase : public G4ProcessState
195     {
196     public:
197       G4ProcessStateBase() :
198           G4ProcessState()
199       {
200       }
201       ~G4ProcessStateBase() override
202       = default;
203 
204       G4String GetType() override
205       {
206         return typeid(T).name();
207       }
208     };
209 
210   template<typename T>
211     T* GetState()
212     {
213       return fpState->GetState<T>();
214     }
215 
216   G4shared_ptr<G4ProcessState> fpState;
217 
218   void virtual SubtractNumberOfInteractionLengthLeft(G4double previousStepSize);
219 
220   inline virtual void ClearInteractionTimeLeft();
221 
222   inline virtual void ClearNumberOfInteractionLengthLeft();
223   // clear NumberOfInteractionLengthLeft
224   // !!! This method should be at the end of PostStepDoIt()
225   // !!! and AtRestDoIt
226   //_________________________________________________
227 
228   void SetInstantiateProcessState(G4bool flag)
229   {
230     fInstantiateProcessState = flag;
231   }
232 
233   G4bool InstantiateProcessState()
234   {
235     return fInstantiateProcessState;
236   }
237 
238   G4bool fProposesTimeStep;
239 
240 private:
241 
242   size_t fProcessID;
243   // During all the simulation will identify a process, so if two identical
244   // processes are created using a copy constructor they will have the same
245   // fProcessID. NOTE: due to MT, this cannot be "const".
246 
247   static/*G4ThreadLocal*/size_t *fNbProcess;
248 
249   G4bool fInstantiateProcessState;
250   //_________________________________________________
251   // Redefine needed members and method of G4VProcess
252   G4double* theNumberOfInteractionLengthLeft;
253   G4double* currentInteractionLength;
254   G4double* theInteractionTimeLeft;
255 };
256 
257 inline void G4VITProcess::ClearInteractionTimeLeft()
258 {
259   fpState->theInteractionTimeLeft = -1.0;
260 }
261 
262 inline void G4VITProcess::ClearNumberOfInteractionLengthLeft()
263 {
264   fpState->theNumberOfInteractionLengthLeft = -1.0;
265 }
266 
267 inline void G4VITProcess::ResetNumberOfInteractionLengthLeft()
268 {
269   fpState->theNumberOfInteractionLengthLeft = -std::log( G4UniformRand());
270 }
271 
272 inline G4double G4VITProcess::GetInteractionTimeLeft()
273 {
274   if (fpState) return fpState->theInteractionTimeLeft;
275 
276   return -1;
277 }
278 
279 inline G4bool G4VITProcess::ProposesTimeStep() const
280 {
281   return fProposesTimeStep;
282 }
283 
284 inline const size_t& G4VITProcess::GetMaxProcessIndex()
285 {
286   if (fNbProcess == nullptr) fNbProcess = new size_t(0);
287   return *fNbProcess;
288 }
289 
290 inline
291 void G4VITProcess::SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
292 {
293   if (fpState->currentInteractionLength > 0.0)
294   {
295     fpState->theNumberOfInteractionLengthLeft -= previousStepSize
296         / fpState->currentInteractionLength;
297     if (fpState->theNumberOfInteractionLengthLeft < 0.)
298     {
299       fpState->theNumberOfInteractionLengthLeft = CLHEP::perMillion;
300     }
301 
302   }
303   else
304   {
305 #ifdef G4VERBOSE
306     if (verboseLevel > 0)
307     {
308       G4cerr << "G4VITProcess::SubtractNumberOfInteractionLengthLeft()";
309       G4cerr << " [" << theProcessName << "]" << G4endl;
310       G4cerr << " currentInteractionLength = "
311              << fpState->currentInteractionLength << " [mm]";
312       G4cerr << " previousStepSize = " << previousStepSize << " [mm]";
313       G4cerr << G4endl;
314     }
315 #endif
316              G4String msg = "Negative currentInteractionLength for ";
317              msg += theProcessName;
318              G4Exception("G4VITProcess::SubtractNumberOfInteractionLengthLeft()",
319                 "ProcMan201",EventMustBeAborted,
320                 msg);
321   }
322 }
323 #endif // G4VITProcess_H
324