Geant4  10.00.p02
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 82326 2014-06-16 09:19:18Z gcosmo $
27 //
28 // Author: Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
29 //
30 // WARNING : This class is released as a prototype.
31 // It might strongly evolve or even disapear in the next releases.
32 //
33 // History:
34 // -----------
35 // 10 Oct 2011 M.Karamitros created
36 //
37 // -------------------------------------------------------------------
38 
39 #ifndef G4VITProcess_H
40 #define G4VITProcess_H
41 
42 #include <G4VProcess.hh>
43 #include "AddClone_def.hh"
44 #include "G4ReferenceCast.hh"
45 #include "G4shared_ptr.hh"
46 
47 class G4IT ;
49 
51  inline virtual ~G4ProcessState_Lock(){;}
52 };
53 
54 #define DowncastProcessState(destinationType) \
55  G4::dynamic_pointer_cast<destinationType>(G4VITProcess::fpState)
56 
57 #define UpcastProcessState(destinationType) \
58  G4::dynamic_pointer_cast<destinationType>(G4VITProcess::fpState)
59 
60 #define DowncastState(destinationType,source) \
61  G4::dynamic_pointer_cast<destinationType>(source)
62 
63 #define UpcastState(destinationType,source) \
64  G4::dynamic_pointer_cast<destinationType>(source)
65 
66 
74 class G4VITProcess : public G4VProcess
75 {
76 public:
77  //__________________________________
78  // Constructors & destructors
80 
81  virtual ~G4VITProcess();
82  G4VITProcess(const G4VITProcess& other);
83  G4VITProcess& operator=(const G4VITProcess& other);
84 
85  // equal opperators
86  G4int operator==(const G4VITProcess &right) const;
87  G4int operator!=(const G4VITProcess &right) const;
88 
90 
91  size_t GetProcessID() const
92  {
93  return fProcessID;
94  }
95 
96  G4::shared_ptr<G4ProcessState_Lock> GetProcessState()
97  {
99  }
100 
101  void SetProcessState(G4::shared_ptr<G4ProcessState_Lock> aProcInfo)
102  {
103  fpState = DowncastState(G4ProcessState,aProcInfo);
104  }
105 
107  {
108  fpState.reset();
109  }
110 
111 
112  //__________________________________
113  // Initialize and Save process info
114 
115  virtual void StartTracking(G4Track*);
116 
118 
120 
124  virtual void ResetNumberOfInteractionLengthLeft();
125 
126  inline G4bool ProposesTimeStep() const;
127 
128  inline static const size_t& GetMaxProcessIndex();
129 
130 protected: // with description
131 
132  void RetrieveProcessInfo();
133  void CreateInfo();
134 
135  //__________________________________
136  // Process info
137  // friend class G4TrackingInformation ;
138 
140  {
141  public:
142  G4ProcessState();
143  virtual ~G4ProcessState();
144 
146  // The flight length left for the current tracking particle
147  // in unit of "Interaction length".
148 
150  // Time left before the interaction : for at rest processes
151 
153  // The InteractionLength in the current material
154 
155  template<typename T> T* GetState()
156  {
157  return dynamic_cast<T*>(this);
158  }
159  };
160 
161  G4::shared_ptr<G4ProcessState> fpState;
162 
163  template<typename T> T* GetState()
164  {
165  return fpState->GetState<T>();
166  }
167 
168  inline virtual void ClearInteractionTimeLeft();
169 
170  inline virtual void ClearNumberOfInteractionLengthLeft();
171  // clear NumberOfInteractionLengthLeft
172  // !!! This method should be at the end of PostStepDoIt()
173  // !!! and AtRestDoIt
174  //_________________________________________________
175 
176 
178  { fInstantiateProcessState = flag; }
179 
181 
183 
184  private :
185 
186  size_t fProcessID;
187  // During all the simulation will identify a process, so if two identical
188  // processes are created using a copy constructor they will have the same
189  // fProcessID. NOTE: due to MT, this cannot be "const".
190 
191  static /*G4ThreadLocal*/ size_t *fNbProcess;
192 
194  //_________________________________________________
195  // Redefine needed members and method of G4VProcess
199 };
200 
202 {
203  fpState->theInteractionTimeLeft = -1.0;
204 }
205 
207 {
208  fpState->theNumberOfInteractionLengthLeft = -1.0;
209 }
210 
212 {
213  fpState->theNumberOfInteractionLengthLeft = -std::log( G4UniformRand() );
214 }
215 
217 {
218  if(fpState)
219  return fpState->theInteractionTimeLeft ;
220 
221  return -1 ;
222 }
223 
225 {
226  return fProposesTimeStep;
227 }
228 
229 inline const size_t& G4VITProcess::GetMaxProcessIndex()
230 {
231  if (!fNbProcess) fNbProcess = new size_t ( 0);
232  return *fNbProcess ;
233 }
234 #endif // G4VITProcess_H
G4IT is a interface which allows the inheriting object :
Definition: G4IT.hh:82
virtual ~G4VITProcess()
Definition: G4VITProcess.cc:58
#define DowncastState(destinationType, source)
Definition: G4VITProcess.hh:60
virtual void ClearNumberOfInteractionLengthLeft()
static const size_t & GetMaxProcessIndex()
void RetrieveProcessInfo()
#define G4IT_TO_BE_CLONED(parent_class)
Definition: AddClone_def.hh:42
virtual void ResetNumberOfInteractionLengthLeft()
WARNING : Redefine the method of G4VProcess reset (determine the value of)NumberOfInteractionLengthLe...
size_t fProcessID
G4String name
Definition: TRTMaterials.hh:40
G4VITProcess(const G4String &name, G4ProcessType type=fNotDefined)
Definition: G4VITProcess.cc:34
The class G4TrackingInformation (hold by G4IT) emcompasses processes informations computed at the PS/...
#define UpcastProcessState(destinationType)
Definition: G4VITProcess.hh:57
int G4int
Definition: G4Types.hh:78
void SetInstantiateProcessState(G4bool flag)
virtual void ClearInteractionTimeLeft()
G4int operator==(const G4VITProcess &right) const
#define G4UniformRand()
Definition: Randomize.hh:87
G4VITProcess & operator=(const G4VITProcess &other)
Definition: G4VITProcess.cc:74
virtual ~G4ProcessState_Lock()
Definition: G4VITProcess.hh:51
bool G4bool
Definition: G4Types.hh:79
virtual void StartTracking(G4Track *)
Definition: G4VITProcess.cc:81
G4double GetInteractionTimeLeft()
G4::shared_ptr< G4ProcessState > fpState
G4double * theNumberOfInteractionLengthLeft
void SetProcessState(G4::shared_ptr< G4ProcessState_Lock > aProcInfo)
G4bool fInstantiateProcessState
G4::shared_ptr< G4ProcessState_Lock > GetProcessState()
Definition: G4VITProcess.hh:96
G4bool fProposesTimeStep
G4bool InstantiateProcessState()
static size_t * fNbProcess
G4int operator!=(const G4VITProcess &right) const
double G4double
Definition: G4Types.hh:76
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4bool ProposesTimeStep() const
G4VITProcess inherits from G4VProcess.
Definition: G4VITProcess.hh:74
G4double * currentInteractionLength
void ResetProcessState()
size_t GetProcessID() const
Definition: G4VITProcess.hh:91
G4double * theInteractionTimeLeft
void CreateInfo()
G4ProcessType