Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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 64057 2012-10-30 15:04:49Z 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
46class G4IT ;
48
50 inline virtual ~G4ProcessState_Lock(){;}
51};
52
53#define InitProcessState(destination,source) \
54 destination(reference_cast(destination,source))
55
56/**
57 * G4VITProcess inherits from G4VProcess.
58 * A G4VITProcess is able to save its current state for a given track into G4IT.
59 * This state may be retrieve latter on to be used by the G4VITProcess.
60 * Each G4VITProcess is tagged.
61 */
62
64{
65public:
66 //__________________________________
67 // Constructors & destructors
69
70 virtual ~G4VITProcess();
71 G4VITProcess(const G4VITProcess& other);
72 G4VITProcess& operator=(const G4VITProcess& other);
73
74 // equal opperators
75 G4int operator==(const G4VITProcess &right) const;
76 G4int operator!=(const G4VITProcess &right) const;
77
79
80 size_t GetProcessID() const
81 {
82 return fProcessID;
83 }
84
86 {
87 return fpState;
88 }
89
91 {
92 fpState = (G4ProcessState*) aProcInfo;
93 }
94
95 //__________________________________
96 // Initialize and Save process info
97
98 virtual void StartTracking(G4Track*);
99
101
103
104 /** WARNING : Redefine the method of G4VProcess
105 * reset (determine the value of)NumberOfInteractionLengthLeft
106 */
108
109 inline G4bool ProposesTimeStep() const;
110
111 inline static const size_t& GetMaxProcessIndex();
112
113protected: // with description
114
117
118 //__________________________________
119 // Process info
120 // friend class G4TrackingInformation ;
121
123 {
124 public:
126 virtual ~G4ProcessState();
127
129 // The flight length left for the current tracking particle
130 // in unit of "Interaction length".
131
133 // Time left before the interaction : for at rest processes
134
136 // The InteractionLength in the current material
137 };
138
140
141 inline virtual void ClearInteractionTimeLeft();
142
143 //_________________________________________________
144 // Redefine needed members and method of G4VProcess
146 G4double previousStepSize
147 );
148 // subtract NumberOfInteractionLengthLeft by the value corresponding to
149 // previousStepSize
150
151 inline virtual void ClearNumberOfInteractionLengthLeft();
152 // clear NumberOfInteractionLengthLeft
153 // !!! This method should be at the end of PostStepDoIt()
154 // !!! and AtRestDoIt
155 //_________________________________________________
156
157
159 { fInstantiateProcessState = flag; }
160
161 G4bool InstantiateProcessState() { return fInstantiateProcessState; }
162
164
165private :
166 const size_t fProcessID; // During all the simulation will identify a
167 // process, so if two identical process are created using a copy constructor
168 // they will have the same fProcessID
169 static size_t fNbProcess ;
170
171 G4bool fInstantiateProcessState;
172 //_________________________________________________
173 // Redefine needed members and method of G4VProcess
174 G4double* theNumberOfInteractionLengthLeft;
175 G4double* currentInteractionLength;
176 G4double* theInteractionTimeLeft;
177};
178
180{
182}
183
185{
187}
188
190{
192}
193
195{
196 if(fpState)
198
199 return -1 ;
200}
201
203{
204 return fProposesTimeStep;
205}
206
208{
209 return fNbProcess ;
210}
211#endif // G4VITProcess_H
#define G4IT_TO_BE_CLONED(parent_class)
Definition: AddClone_def.hh:42
G4ProcessType
@ fNotDefined
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4UniformRand()
Definition: Randomize.hh:53
Definition: G4IT.hh:83
G4bool ProposesTimeStep() const
virtual void StartTracking(G4Track *)
Definition: G4VITProcess.cc:81
void SetInstantiateProcessState(G4bool flag)
size_t GetProcessID() const
Definition: G4VITProcess.hh:80
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
Definition: G4VITProcess.cc:96
G4int operator==(const G4VITProcess &right) const
void CreateInfo()
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4int operator!=(const G4VITProcess &right) const
G4ProcessState * fpState
void SetProcessState(G4ProcessState_Lock *aProcInfo)
Definition: G4VITProcess.hh:90
virtual void ResetNumberOfInteractionLengthLeft()
virtual void ClearInteractionTimeLeft()
G4VITProcess & operator=(const G4VITProcess &other)
Definition: G4VITProcess.cc:74
void RetrieveProcessInfo()
G4bool InstantiateProcessState()
G4ProcessState_Lock * GetProcessState()
Definition: G4VITProcess.hh:85
G4double GetInteractionTimeLeft()
static const size_t & GetMaxProcessIndex()
virtual void ClearNumberOfInteractionLengthLeft()
G4bool fProposesTimeStep
virtual ~G4VITProcess()
Definition: G4VITProcess.cc:57
virtual ~G4ProcessState_Lock()
Definition: G4VITProcess.hh:50