Geant4  10.03
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 94983 2016-01-13 11:02:33Z 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;
72 
73 #include <map>
74 
76 class G4Track
78 {
79 
80 //--------
81 public: // With description
82 
83 // Constructor
84  G4Track();
85  G4Track(G4DynamicParticle* apValueDynamicParticle,
86  G4double aValueTime,
87  const G4ThreeVector& aValuePosition);
88  // aValueTime is a global time
89  G4Track(const G4Track&);
90  // Copy Constructor copys members other than tracking information
91 
92 private:
93  // Hide assignment operator as private
94  G4Track& operator=(const G4Track&);
95 
96 //--------
97 public: // With description
98 
99 // Destrcutor
100  ~G4Track();
101 
102 // Operators
103  inline void *operator new(size_t);
104  // Override "new" for "G4Allocator".
105  inline void operator delete(void *aTrack);
106  // Override "delete" for "G4Allocator".
107 
108  G4bool operator==( const G4Track& );
109 
110 //--------
111 public: // With description
112 // Copy information of the track (w/o tracking information)
113  void CopyTrackInfo(const G4Track&);
114 
115 // Get/Set functions
116  // track ID
117  G4int GetTrackID() const;
118  void SetTrackID(const G4int aValue);
119 
120  G4int GetParentID() const;
121  void SetParentID(const G4int aValue);
122 
123  // dynamic particle
124  const G4DynamicParticle* GetDynamicParticle() const;
125 
126  // particle definition
128  // following method of GetDefinition remains
129  // because of backward compatiblity. It will be removed in future
131 
132  // position, time
133  const G4ThreeVector& GetPosition() const;
134  void SetPosition(const G4ThreeVector& aValue);
135 
136  G4double GetGlobalTime() const;
137  void SetGlobalTime(const G4double aValue);
138  // Time since the event in which the track belongs is created.
139 
140  G4double GetLocalTime() const;
141  void SetLocalTime(const G4double aValue);
142  // Time since the current track is created.
143 
144  G4double GetProperTime() const;
145  void SetProperTime(const G4double aValue);
146  // Proper time of the current track
147 
148  // volume, material, touchable
149  G4VPhysicalVolume* GetVolume() const;
151 
152  G4Material* GetMaterial() const;
153  G4Material* GetNextMaterial() const;
154 
157 
158  const G4VTouchable* GetTouchable() const;
159  const G4TouchableHandle& GetTouchableHandle() const;
160  void SetTouchableHandle( const G4TouchableHandle& apValue);
161 
162  const G4VTouchable* GetNextTouchable() const;
164  void SetNextTouchableHandle( const G4TouchableHandle& apValue);
165 
166  const G4VTouchable* GetOriginTouchable() const;
168  void SetOriginTouchableHandle( const G4TouchableHandle& apValue);
169 
170  // energy
171  G4double GetKineticEnergy() const;
172  void SetKineticEnergy(const G4double aValue);
173 
174  G4double GetTotalEnergy() const;
175 
176 
177  // moemtnum
178  const G4ThreeVector& GetMomentumDirection() const;
179  void SetMomentumDirection(const G4ThreeVector& aValue);
180 
181  G4ThreeVector GetMomentum() const;
182 
183  // velocity
184  G4double GetVelocity() const;
185  void SetVelocity(G4double val);
186 
187  G4double CalculateVelocity() const;
189 
190  G4bool UseGivenVelocity() const;
191  void UseGivenVelocity(G4bool val);
192 
193  // polarization
194  const G4ThreeVector& GetPolarization() const;
195  void SetPolarization(const G4ThreeVector& aValue);
196 
197  // track status, flags for tracking
199  void SetTrackStatus(const G4TrackStatus aTrackStatus);
200 
201  G4bool IsBelowThreshold() const;
202  void SetBelowThresholdFlag(G4bool value = true);
203  // The flag of "BelowThreshold" is set to true
204  // if this track energy is below threshold energy
205  // in this material determined by the range cut value
206 
207  G4bool IsGoodForTracking() const;
208  void SetGoodForTrackingFlag(G4bool value = true);
209  // The flag of "GoodForTracking" is set by processes
210  // if this track should be tracked
211  // even if the energy is below threshold
212 
213  // track length
214  G4double GetTrackLength() const;
215  void AddTrackLength(const G4double aValue);
216  // Accumulated the track length
217 
218  // step information
219  const G4Step* GetStep() const;
220  void SetStep(const G4Step* aValue);
221 
222  G4int GetCurrentStepNumber() const;
224 
225  G4double GetStepLength() const;
226  void SetStepLength(G4double value);
227  // Before the end of the AlongStepDoIt loop,StepLength keeps
228  // the initial value which is determined by the shortest geometrical Step
229  // proposed by a physics process. After finishing the AlongStepDoIt,
230  // it will be set equal to 'StepLength' in G4Step.
231 
232  // vertex (,where this track was created) information
233  const G4ThreeVector& GetVertexPosition() const;
234  void SetVertexPosition(const G4ThreeVector& aValue);
235 
237  void SetVertexMomentumDirection(const G4ThreeVector& aValue);
238 
240  void SetVertexKineticEnergy(const G4double aValue);
241 
244 
245  const G4VProcess* GetCreatorProcess() const;
246  void SetCreatorProcess(const G4VProcess* aValue);
247 
248  inline void SetCreatorModelIndex(G4int idx);
249 
250  inline const G4String& GetCreatorModelName() const;
251 
252  inline G4int GetCreatorModelID() const;
253 
254  // track weight
255  // These are methods for manipulating a weight for this track.
256  G4double GetWeight() const;
257  void SetWeight(G4double aValue);
258 
259  // User information
261  void SetUserInformation(G4VUserTrackInformation* aValue) const;
262 
263  // Velocity table
264  static void SetVelocityTableProperties(G4double t_max, G4double t_min, G4int nbin);
267  static G4int GetNbinOfVelocityTable();
268 
269 //---------
270  private:
271 //---------
272  // Member data
273  G4int fCurrentStepNumber; // Total steps number up to now
274  G4ThreeVector fPosition; // Current positon
275  G4double fGlobalTime; // Time since the event is created
276  G4double fLocalTime; // Time since the track is created
277  G4double fTrackLength; // Accumulated track length
281 
285  // Touchable Handle
286 
289 
291  // This flag is set to true if this track energy is below
292  // threshold energy in this material determined by the range cut value
294  // This flag is set by processes if this track should be tracked
295  // even if the energy is below threshold
296 
298  // Before the end of the AlongStepDoIt loop, this keeps the initial
299  // Step length which is determined by the shortest geometrical Step
300  // proposed by a physics process. After finishing the AlongStepDoIt,
301  // this will be set equal to 'StepLength' in G4Step.
302 
304  // This is a weight for this track
305 
306  const G4Step* fpStep;
307 
308  G4ThreeVector fVtxPosition; // (x,y,z) of the vertex
309  G4ThreeVector fVtxMomentumDirection; // Momentum direction at the vertex
310  G4double fVtxKineticEnergy; // Kinetic energy at the vertex
311  const G4LogicalVolume* fpLVAtVertex; //Logical Volume at the vertex
312  const G4VProcess* fpCreatorProcess; // Process which created the track
313  G4int fCreatorModelIndex; // Index of the physics model which created the track
314 
316 
317  // cached values for CalculateVelocity
322 
325 
327  // do not calclulate velocity and just use current fVelocity
328  // if this flag is set
329 
330  mutable std::map<G4int,G4VAuxiliaryTrackInformation*>* fpAuxiliaryTrackInformationMap;
331 
332 //--------
333 public:
334 //--------
335 
338  std::map<G4int,G4VAuxiliaryTrackInformation*>* GetAuxiliaryTrackInformationMap() const
340 
343  // Note: G4VAuxiliaryTrackInformation object itself is *NOT* deleted
344 
345 //--------
346 private:
347 //--------
348 
350 };
351 
352 #include "G4Track.icc"
353 
354 #endif
static G4ThreadLocal G4VelocityTable * velTable
Definition: G4Track.hh:324
G4int GetCreatorModelID() const
G4int fCurrentStepNumber
Definition: G4Track.hh:273
const G4TouchableHandle & GetOriginTouchableHandle() const
void SetTrackStatus(const G4TrackStatus aTrackStatus)
std::map< G4int, G4VAuxiliaryTrackInformation * > * fpAuxiliaryTrackInformationMap
Definition: G4Track.hh:330
G4ParticleDefinition * GetDefinition() const
G4double fGlobalTime
Definition: G4Track.hh:275
G4int GetParentID() const
const G4ThreeVector & GetPolarization() const
void SetVelocity(G4double val)
void SetVertexMomentumDirection(const G4ThreeVector &aValue)
G4double GetLocalTime() const
G4int fTrackID
Definition: G4Track.hh:279
G4TouchableHandle fpOriginTouchable
Definition: G4Track.hh:284
G4double GetProperTime() const
void SetOriginTouchableHandle(const G4TouchableHandle &apValue)
CLHEP::Hep3Vector G4ThreeVector
G4double GetVelocity() const
const G4LogicalVolume * GetLogicalVolumeAtVertex() const
const G4DynamicParticle * GetDynamicParticle() const
void RemoveAuxiliaryTrackInformation(G4int idx)
Definition: G4Track.cc:351
const G4String & GetCreatorModelName() const
G4DynamicParticle * fpDynamicParticle
Definition: G4Track.hh:287
G4MaterialPropertyVector * groupvel
Definition: G4Track.hh:319
void SetAuxiliaryTrackInformation(G4int idx, G4VAuxiliaryTrackInformation *info) const
Definition: G4Track.cc:324
const G4VTouchable * GetOriginTouchable() const
const G4ThreeVector & GetPosition() const
G4TrackStatus GetTrackStatus() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void SetNextTouchableHandle(const G4TouchableHandle &apValue)
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4bool is_OpticalPhoton
Definition: G4Track.hh:323
void SetBelowThresholdFlag(G4bool value=true)
G4double fVtxKineticEnergy
Definition: G4Track.hh:310
std::map< G4int, G4VAuxiliaryTrackInformation * > * GetAuxiliaryTrackInformationMap() const
Definition: G4Track.hh:338
G4TouchableHandle fpNextTouchable
Definition: G4Track.hh:283
const char * name(G4int ptype)
#define G4ThreadLocal
Definition: tls.hh:89
const G4Step * GetStep() const
int G4int
Definition: G4Types.hh:78
G4double fLocalTime
Definition: G4Track.hh:276
void SetWeight(G4double aValue)
G4double prev_velocity
Definition: G4Track.hh:320
G4bool useGivenVelocity
Definition: G4Track.hh:326
void SetCreatorModelIndex(G4int idx)
void SetCreatorProcess(const G4VProcess *aValue)
G4VPhysicalVolume * GetNextVolume() const
G4VUserTrackInformation * GetUserInformation() const
const G4VProcess * fpCreatorProcess
Definition: G4Track.hh:312
G4VUserTrackInformation * fpUserInformation
Definition: G4Track.hh:315
G4double CalculateVelocityForOpticalPhoton() const
Definition: G4Track.cc:255
const G4VProcess * GetCreatorProcess() const
const G4MaterialCutsCouple * GetNextMaterialCutsCouple() const
G4bool fBelowThreshold
Definition: G4Track.hh:290
G4double GetKineticEnergy() const
G4double fWeight
Definition: G4Track.hh:303
void SetPosition(const G4ThreeVector &aValue)
G4ThreeVector fVtxMomentumDirection
Definition: G4Track.hh:309
G4int GetCurrentStepNumber() const
G4double fTrackLength
Definition: G4Track.hh:277
G4int fParentID
Definition: G4Track.hh:278
G4double GetVertexKineticEnergy() const
G4Material * GetNextMaterial() const
void SetGlobalTime(const G4double aValue)
bool G4bool
Definition: G4Types.hh:79
void SetStepLength(G4double value)
G4Material * prev_mat
Definition: G4Track.hh:318
G4ThreeVector fPosition
Definition: G4Track.hh:274
void SetVertexKineticEnergy(const G4double aValue)
const G4TouchableHandle & GetNextTouchableHandle() const
const G4ParticleDefinition * GetParticleDefinition() const
Definition: G4Step.hh:76
G4int GetTrackID() const
G4ThreeVector fVtxPosition
Definition: G4Track.hh:308
G4double GetGlobalTime() const
G4bool UseGivenVelocity() const
G4double CalculateVelocity() const
Definition: G4Track.cc:222
const G4TouchableHandle & GetTouchableHandle() const
const G4ThreeVector & GetVertexPosition() const
G4double GetTrackLength() const
G4Material * GetMaterial() const
void SetUserInformation(G4VUserTrackInformation *aValue) const
G4VAuxiliaryTrackInformation * GetAuxiliaryTrackInformation(G4int idx) const
Definition: G4Track.cc:340
static G4int GetNbinOfVelocityTable()
Definition: G4Track.cc:319
G4TouchableHandle fpTouchable
Definition: G4Track.hh:282
G4ThreeVector GetMomentum() const
G4double fStepLength
Definition: G4Track.hh:297
const G4ThreeVector & GetMomentumDirection() const
void IncrementCurrentStepNumber()
const G4VTouchable * GetNextTouchable() const
G4Track & operator=(const G4Track &)
Definition: G4Track.cc:158
G4double fVelocity
Definition: G4Track.hh:280
static G4double GetMinTOfVelocityTable()
Definition: G4Track.cc:314
G4TrackStatus fTrackStatus
Definition: G4Track.hh:288
void SetPolarization(const G4ThreeVector &aValue)
void SetVertexPosition(const G4ThreeVector &aValue)
const G4VTouchable * GetTouchable() const
void SetParentID(const G4int aValue)
G4bool IsBelowThreshold() const
G4double prev_momentum
Definition: G4Track.hh:321
G4int fCreatorModelIndex
Definition: G4Track.hh:313
static void SetVelocityTableProperties(G4double t_max, G4double t_min, G4int nbin)
Definition: G4Track.cc:301
G4VPhysicalVolume * GetVolume() const
const G4LogicalVolume * fpLVAtVertex
Definition: G4Track.hh:311
void ClearAuxiliaryTrackInformation()
Definition: G4Track.cc:370
G4double GetWeight() const
G4double GetTotalEnergy() const
G4bool IsGoodForTracking() const
static G4double GetMaxTOfVelocityTable()
Definition: G4Track.cc:309
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
G4bool operator==(const G4Track &)
G4TrackStatus
const G4ThreeVector & GetVertexMomentumDirection() const
G4Track()
Definition: G4Track.cc:98
void SetStep(const G4Step *aValue)
const G4Step * fpStep
Definition: G4Track.hh:306
void CopyTrackInfo(const G4Track &)
Definition: G4Track.cc:215
G4double GetStepLength() const
void SetMomentumDirection(const G4ThreeVector &aValue)
void SetTrackID(const G4int aValue)
void SetGoodForTrackingFlag(G4bool value=true)
G4bool fGoodForTracking
Definition: G4Track.hh:293
void SetLogicalVolumeAtVertex(const G4LogicalVolume *)
~G4Track()
Definition: G4Track.cc:149