Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
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
55class G4IT;
57
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
99{
100public:
101 //__________________________________
102 // Constructors & destructors
103 G4VITProcess(const G4String& name, G4ProcessType type = fNotDefined);
104
105 virtual ~G4VITProcess();
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
114
115 size_t GetProcessID() const
116 {
117 return fProcessID;
118 }
119
120 G4shared_ptr<G4ProcessState_Lock> GetProcessState()
121 {
123 }
124
125 void SetProcessState(G4shared_ptr<G4ProcessState_Lock> aProcInfo)
126 {
128 }
129
131 {
132 fpState.reset();
133 }
134
135 //__________________________________
136 // Initialize and Save process info
137
138 virtual void StartTracking(G4Track*);
139
141 {
142 }
143
145
146 /** WARNING : Redefine the method of G4VProcess
147 * reset (determine the value of)NumberOfInteractionLengthLeft
148 */
150
151 inline G4bool ProposesTimeStep() const;
152
153 inline static const size_t& GetMaxProcessIndex();
154
155protected:
156 // with description
157
160
161 //__________________________________
162 // Process info
163 // friend class G4TrackingInformation ;
164
166 {
167 public:
169 virtual ~G4ProcessState();
170
172 {
173 return "G4ProcessState";
174 }
175
177 // The flight length left for the current tracking particle
178 // in unit of "Interaction length".
179
181 // Time left before the interaction : for at rest processes
182
184 // The InteractionLength in the current material
185
186 template<typename T>
188 {
189 return dynamic_cast<T*>(this);
190 }
191 };
192
193 template<typename T>
195 {
196 public:
199 {
200 }
202 {
203 }
204
206 {
207 return typeid(T).name();
208 }
209 };
210
211 template<typename T>
213 {
214 return fpState->GetState<T>();
215 }
216
217 G4shared_ptr<G4ProcessState> fpState;
218
219 void virtual SubtractNumberOfInteractionLengthLeft(G4double previousStepSize);
220
221 inline virtual void ClearInteractionTimeLeft();
222
223 inline virtual void ClearNumberOfInteractionLengthLeft();
224 // clear NumberOfInteractionLengthLeft
225 // !!! This method should be at the end of PostStepDoIt()
226 // !!! and AtRestDoIt
227 //_________________________________________________
228
230 {
231 fInstantiateProcessState = flag;
232 }
233
235 {
236 return fInstantiateProcessState;
237 }
238
240
241private:
242
243 size_t fProcessID;
244 // During all the simulation will identify a process, so if two identical
245 // processes are created using a copy constructor they will have the same
246 // fProcessID. NOTE: due to MT, this cannot be "const".
247
248 static/*G4ThreadLocal*/size_t *fNbProcess;
249
250 G4bool fInstantiateProcessState;
251 //_________________________________________________
252 // Redefine needed members and method of G4VProcess
253 G4double* theNumberOfInteractionLengthLeft;
254 G4double* currentInteractionLength;
255 G4double* theInteractionTimeLeft;
256};
257
259{
260 fpState->theInteractionTimeLeft = -1.0;
261}
262
264{
265 fpState->theNumberOfInteractionLengthLeft = -1.0;
266}
267
269{
270 fpState->theNumberOfInteractionLengthLeft = -std::log( G4UniformRand());
271}
272
274{
275 if (fpState) return fpState->theInteractionTimeLeft;
276
277 return -1;
278}
279
281{
282 return fProposesTimeStep;
283}
284
286{
287 if (!fNbProcess) fNbProcess = new size_t(0);
288 return *fNbProcess;
289}
290
291inline
293{
294 if (fpState->currentInteractionLength > 0.0)
295 {
296 fpState->theNumberOfInteractionLengthLeft -= previousStepSize
297 / fpState->currentInteractionLength;
298 if (fpState->theNumberOfInteractionLengthLeft < 0.)
299 {
300 fpState->theNumberOfInteractionLengthLeft = CLHEP::perMillion;
301 }
302
303 }
304 else
305 {
306#ifdef G4VERBOSE
307 if (verboseLevel > 0)
308 {
309 G4cerr << "G4VITProcess::SubtractNumberOfInteractionLengthLeft()";
310 G4cerr << " [" << theProcessName << "]" << G4endl;
311 G4cerr << " currentInteractionLength = "
312 << fpState->currentInteractionLength << " [mm]";
313 G4cerr << " previousStepSize = " << previousStepSize << " [mm]";
314 G4cerr << G4endl;
315 }
316#endif
317 G4String msg = "Negative currentInteractionLength for ";
318 msg += theProcessName;
319 G4Exception("G4VITProcess::SubtractNumberOfInteractionLengthLeft()",
320 "ProcMan201",EventMustBeAborted,
321 msg);
322 }
323}
324#endif // G4VITProcess_H
#define G4IT_TO_BE_CLONED(parent_class)
Definition: AddClone_def.hh:49
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
G4ProcessType
@ fNotDefined
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#define DowncastState(destinationType, source)
Definition: G4VITProcess.hh:85
#define UpcastProcessState(destinationType)
Definition: G4VITProcess.hh:82
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
#define G4UniformRand()
Definition: Randomize.hh:52
Definition: G4IT.hh:88
G4bool ProposesTimeStep() const
virtual void StartTracking(G4Track *)
Definition: G4VITProcess.cc:85
G4bool operator!=(const G4VITProcess &right) const
void SetInstantiateProcessState(G4bool flag)
size_t GetProcessID() const
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
void CreateInfo()
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4shared_ptr< G4ProcessState > fpState
void ResetProcessState()
G4bool operator==(const G4VITProcess &right) const
G4shared_ptr< G4ProcessState_Lock > GetProcessState()
virtual void ResetNumberOfInteractionLengthLeft()
virtual void ClearInteractionTimeLeft()
void SetProcessState(G4shared_ptr< G4ProcessState_Lock > aProcInfo)
G4VITProcess & operator=(const G4VITProcess &other)
Definition: G4VITProcess.cc:78
void RetrieveProcessInfo()
G4bool InstantiateProcessState()
G4double GetInteractionTimeLeft()
static const size_t & GetMaxProcessIndex()
virtual void ClearNumberOfInteractionLengthLeft()
G4bool fProposesTimeStep
virtual ~G4VITProcess()
Definition: G4VITProcess.cc:60
G4int verboseLevel
Definition: G4VProcess.hh:360
G4String theProcessName
Definition: G4VProcess.hh:345
virtual ~G4ProcessState_Lock()
Definition: G4VITProcess.hh:60
virtual G4String GetType()