Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4Track.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 // $Id: G4Track.hh 74276 2013-10-02 14:59:07Z gcosmo $
28 //
29 //
30 //---------------------------------------------------------------
31 //
32 // G4Track.hh
33 //
34 // Class Description:
35 // This class represents the partilce under tracking.
36 // It includes information related to tracking for examples:
37 // 1) current position/time of the particle,
38 // 2) static particle information,
39 // 3) the pointer to the physical volume where currently
40 // the particle exists
41 //
42 //---------------------------------------------------------------
43 // Modification for G4TouchableHandle 22 Oct. 2001 R.Chytracek
44 // Add MaterialCutCouple 08 Oct. 2002 H.Kurashige
45 // Add SetVelocityTableProperties 02 Apr. 2011 H.Kurashige
46 // Add fVelocity and Set/GetVelocity 29 Apr. 2011 H.Kurashige
47 // Use G4VelocityTable 17 AUg. 2011 H.Kurashige
48 
49 #ifndef G4Track_h
50 #define G4Track_h 1
51 
52 #include <cmath> // Include from 'system'
53 
54 #include "globals.hh" // Include from 'global'
55 #include "trkdefs.hh" // Include DLL defs...
56 #include "G4ThreeVector.hh" // Include from 'geometry'
57 #include "G4LogicalVolume.hh" // Include from 'geometry'
58 #include "G4VPhysicalVolume.hh" // Include from 'geometry'
59 #include "G4Allocator.hh" // Include from 'particle+matter'
60 #include "G4DynamicParticle.hh" // Include from 'particle+matter'
61 #include "G4TrackStatus.hh" // Include from 'tracking'
62 #include "G4TouchableHandle.hh" // Include from 'geometry'
64 #include "G4PhysicsModelCatalog.hh"
65 
66 #include "G4Material.hh"
67 
68 class G4Step; // Forward declaration
70 class G4VelocityTable;
71 
72 //////////////
73 class G4Track
74 //////////////
75 {
76 
77 //--------
78 public: // With description
79 
80 // Constructor
81  G4Track();
82  G4Track(G4DynamicParticle* apValueDynamicParticle,
83  G4double aValueTime,
84  const G4ThreeVector& aValuePosition);
85  // aValueTime is a global time
86  G4Track(const G4Track&);
87  // Copy Constructor copys members other than tracking information
88 
89 private:
90  // Hide assignment operator as private
91  G4Track& operator=(const G4Track&);
92 
93 //--------
94 public: // With description
95 
96 // Destrcutor
97  ~G4Track();
98 
99 // Operators
100  inline void *operator new(size_t);
101  // Override "new" for "G4Allocator".
102  inline void operator delete(void *aTrack);
103  // Override "delete" for "G4Allocator".
104 
105  G4bool operator==( const G4Track& );
106 
107 //--------
108 public: // With description
109 // Copy information of the track (w/o tracking information)
110  void CopyTrackInfo(const G4Track&);
111 
112 // Get/Set functions
113  // track ID
114  G4int GetTrackID() const;
115  void SetTrackID(const G4int aValue);
116 
117  G4int GetParentID() const;
118  void SetParentID(const G4int aValue);
119 
120  // dynamic particle
121  const G4DynamicParticle* GetDynamicParticle() const;
122 
123  // particle definition
125  // following method of GetDefinition remains
126  // because of backward compatiblity. It will be removed in future
128 
129  // position, time
130  const G4ThreeVector& GetPosition() const;
131  void SetPosition(const G4ThreeVector& aValue);
132 
133  G4double GetGlobalTime() const;
134  void SetGlobalTime(const G4double aValue);
135  // Time since the event in which the track belongs is created.
136 
137  G4double GetLocalTime() const;
138  void SetLocalTime(const G4double aValue);
139  // Time since the current track is created.
140 
141  G4double GetProperTime() const;
142  void SetProperTime(const G4double aValue);
143  // Proper time of the current track
144 
145  // volume, material, touchable
146  G4VPhysicalVolume* GetVolume() const;
148 
149  G4Material* GetMaterial() const;
150  G4Material* GetNextMaterial() const;
151 
154 
155  const G4VTouchable* GetTouchable() const;
156  const G4TouchableHandle& GetTouchableHandle() const;
157  void SetTouchableHandle( const G4TouchableHandle& apValue);
158 
159  const G4VTouchable* GetNextTouchable() const;
161  void SetNextTouchableHandle( const G4TouchableHandle& apValue);
162 
163  const G4VTouchable* GetOriginTouchable() const;
165  void SetOriginTouchableHandle( const G4TouchableHandle& apValue);
166 
167  // energy
168  G4double GetKineticEnergy() const;
169  void SetKineticEnergy(const G4double aValue);
170 
171  G4double GetTotalEnergy() const;
172 
173 
174  // moemtnum
175  const G4ThreeVector& GetMomentumDirection() const;
176  void SetMomentumDirection(const G4ThreeVector& aValue);
177 
178  G4ThreeVector GetMomentum() const;
179 
180  // velocity
181  G4double GetVelocity() const;
182  void SetVelocity(G4double val);
183 
184  G4double CalculateVelocity() const;
186 
187  G4bool UseGivenVelocity() const;
188  void UseGivenVelocity(G4bool val);
189 
190  // polarization
191  const G4ThreeVector& GetPolarization() const;
192  void SetPolarization(const G4ThreeVector& aValue);
193 
194  // track status, flags for tracking
196  void SetTrackStatus(const G4TrackStatus aTrackStatus);
197 
198  G4bool IsBelowThreshold() const;
199  void SetBelowThresholdFlag(G4bool value = true);
200  // The flag of "BelowThreshold" is set to true
201  // if this track energy is below threshold energy
202  // in this material determined by the range cut value
203 
204  G4bool IsGoodForTracking() const;
205  void SetGoodForTrackingFlag(G4bool value = true);
206  // The flag of "GoodForTracking" is set by processes
207  // if this track should be tracked
208  // even if the energy is below threshold
209 
210  // track length
211  G4double GetTrackLength() const;
212  void AddTrackLength(const G4double aValue);
213  // Accumulated the track length
214 
215  // step information
216  const G4Step* GetStep() const;
217  void SetStep(const G4Step* aValue);
218 
219  G4int GetCurrentStepNumber() const;
221 
222  G4double GetStepLength() const;
224  // Before the end of the AlongStepDoIt loop,StepLength keeps
225  // the initial value which is determined by the shortest geometrical Step
226  // proposed by a physics process. After finishing the AlongStepDoIt,
227  // it will be set equal to 'StepLength' in G4Step.
228 
229  // vertex (,where this track was created) information
230  const G4ThreeVector& GetVertexPosition() const;
231  void SetVertexPosition(const G4ThreeVector& aValue);
232 
234  void SetVertexMomentumDirection(const G4ThreeVector& aValue);
235 
237  void SetVertexKineticEnergy(const G4double aValue);
238 
241 
242  const G4VProcess* GetCreatorProcess() const;
243  void SetCreatorProcess(const G4VProcess* aValue);
244 
245  inline void SetCreatorModelIndex(G4int idx)
246  { fCreatorModelIndex = idx; }
248  { return G4PhysicsModelCatalog::GetModelName(fCreatorModelIndex); }
250  { return fCreatorModelIndex; }
251 
252  // track weight
253  // These are methods for manipulating a weight for this track.
254  G4double GetWeight() const;
255  void SetWeight(G4double aValue);
256 
257  // User information
260 
261  // Velocity table
262  static void SetVelocityTableProperties(G4double t_max, G4double t_min, G4int nbin);
265  static G4int GetNbinOfVelocityTable();
266 
267 //---------
268  private:
269 //---------
270  // Member data
271  G4int fCurrentStepNumber; // Total steps number up to now
272  G4ThreeVector fPosition; // Current positon
273  G4double fGlobalTime; // Time since the event is created
274  G4double fLocalTime; // Time since the track is created
275  G4double fTrackLength; // Accumulated track length
276  G4int fParentID;
277  G4int fTrackID;
278  G4double fVelocity;
279 
280  G4TouchableHandle fpTouchable;
281  G4TouchableHandle fpNextTouchable;
282  G4TouchableHandle fpOriginTouchable;
283  // Touchable Handle
284 
285  G4DynamicParticle* fpDynamicParticle;
286  G4TrackStatus fTrackStatus;
287 
288  G4bool fBelowThreshold;
289  // This flag is set to true if this track energy is below
290  // threshold energy in this material determined by the range cut value
291  G4bool fGoodForTracking;
292  // This flag is set by processes if this track should be tracked
293  // even if the energy is below threshold
294 
295  G4double fStepLength;
296  // Before the end of the AlongStepDoIt loop, this keeps the initial
297  // Step length which is determined by the shortest geometrical Step
298  // proposed by a physics process. After finishing the AlongStepDoIt,
299  // this will be set equal to 'StepLength' in G4Step.
300 
301  G4double fWeight;
302  // This is a weight for this track
303 
304  const G4Step* fpStep;
305 
306  G4ThreeVector fVtxPosition; // (x,y,z) of the vertex
307  G4ThreeVector fVtxMomentumDirection; // Momentum direction at the vertex
308  G4double fVtxKineticEnergy; // Kinetic energy at the vertex
309  const G4LogicalVolume* fpLVAtVertex; //Logical Volume at the vertex
310  const G4VProcess* fpCreatorProcess; // Process which created the track
311  G4int fCreatorModelIndex; // Index of the physics model which created the track
312 
313  G4VUserTrackInformation* fpUserInformation;
314 
315  // cached values for CalculateVelocity
316  mutable G4Material* prev_mat;
317  mutable G4MaterialPropertyVector* groupvel;
318  mutable G4double prev_velocity;
319  mutable G4double prev_momentum;
320 
321  G4bool is_OpticalPhoton;
322  static G4ThreadLocal G4VelocityTable* velTable;
323 
324  G4bool useGivenVelocity;
325  // do not calclulate velocity and just use current fVelocity
326  // if this flag is set
327 };
328 
329 #include "G4Track.icc"
330 
331 #endif
const G4TouchableHandle & GetOriginTouchableHandle() const
void SetTrackStatus(const G4TrackStatus aTrackStatus)
G4ParticleDefinition * GetDefinition() const
G4int GetParentID() const
const G4ThreeVector & GetPolarization() const
void SetVelocity(G4double val)
void SetVertexMomentumDirection(const G4ThreeVector &aValue)
G4double GetLocalTime() const
G4double GetProperTime() const
void SetOriginTouchableHandle(const G4TouchableHandle &apValue)
G4double GetVelocity() const
const G4LogicalVolume * GetLogicalVolumeAtVertex() const
const G4DynamicParticle * GetDynamicParticle() const
const G4VTouchable * GetOriginTouchable() const
const G4ThreeVector & GetPosition() const
G4TrackStatus GetTrackStatus() const
G4int GetCreatorModelID()
Definition: G4Track.hh:249
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void SetNextTouchableHandle(const G4TouchableHandle &apValue)
void SetTouchableHandle(const G4TouchableHandle &apValue)
void SetBelowThresholdFlag(G4bool value=true)
#define G4ThreadLocal
Definition: tls.hh:52
const G4Step * GetStep() const
int G4int
Definition: G4Types.hh:78
void SetWeight(G4double aValue)
void SetCreatorModelIndex(G4int idx)
Definition: G4Track.hh:245
void SetCreatorProcess(const G4VProcess *aValue)
G4VPhysicalVolume * GetNextVolume() const
G4VUserTrackInformation * GetUserInformation() const
G4double CalculateVelocityForOpticalPhoton() const
Definition: G4Track.cc:247
const G4VProcess * GetCreatorProcess() const
const G4MaterialCutsCouple * GetNextMaterialCutsCouple() const
G4double GetKineticEnergy() const
void SetPosition(const G4ThreeVector &aValue)
G4int GetCurrentStepNumber() const
G4double GetVertexKineticEnergy() const
G4Material * GetNextMaterial() const
void SetGlobalTime(const G4double aValue)
bool G4bool
Definition: G4Types.hh:79
void SetStepLength(G4double value)
static G4String & GetModelName(G4int)
void SetVertexKineticEnergy(const G4double aValue)
const G4TouchableHandle & GetNextTouchableHandle() const
const G4ParticleDefinition * GetParticleDefinition() const
Definition: G4Step.hh:76
G4int GetTrackID() const
G4double GetGlobalTime() const
G4String & GetCreatorModelName()
Definition: G4Track.hh:247
G4bool UseGivenVelocity() const
G4double CalculateVelocity() const
Definition: G4Track.cc:214
const G4TouchableHandle & GetTouchableHandle() const
const G4ThreeVector & GetVertexPosition() const
G4double GetTrackLength() const
G4Material * GetMaterial() const
static G4int GetNbinOfVelocityTable()
Definition: G4Track.cc:311
G4ThreeVector GetMomentum() const
const G4ThreeVector & GetMomentumDirection() const
void IncrementCurrentStepNumber()
const G4VTouchable * GetNextTouchable() const
static G4double GetMinTOfVelocityTable()
Definition: G4Track.cc:306
void SetPolarization(const G4ThreeVector &aValue)
void SetVertexPosition(const G4ThreeVector &aValue)
const G4VTouchable * GetTouchable() const
void SetParentID(const G4int aValue)
G4bool IsBelowThreshold() const
static void SetVelocityTableProperties(G4double t_max, G4double t_min, G4int nbin)
Definition: G4Track.cc:293
G4VPhysicalVolume * GetVolume() const
G4double GetWeight() const
G4double GetTotalEnergy() const
const XML_Char int const XML_Char * value
G4bool IsGoodForTracking() const
static G4double GetMaxTOfVelocityTable()
Definition: G4Track.cc:301
void SetLocalTime(const G4double aValue)
void AddTrackLength(const G4double aValue)
void SetProperTime(const G4double aValue)
void SetKineticEnergy(const G4double aValue)
double G4double
Definition: G4Types.hh:76
void SetUserInformation(G4VUserTrackInformation *aValue)
G4bool operator==(const G4Track &)
G4TrackStatus
const G4ThreeVector & GetVertexMomentumDirection() const
G4Track()
Definition: G4Track.cc:95
void SetStep(const G4Step *aValue)
void CopyTrackInfo(const G4Track &)
Definition: G4Track.cc:207
G4double GetStepLength() const
void SetMomentumDirection(const G4ThreeVector &aValue)
void SetTrackID(const G4int aValue)
void SetGoodForTrackingFlag(G4bool value=true)
void SetLogicalVolumeAtVertex(const G4LogicalVolume *)
~G4Track()
Definition: G4Track.cc:144