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