Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// G4Track
27//
28// Class description:
29//
30// This class describes the particle under tracking.
31// It includes information related to tracking, i.e.:
32// 1) current position/time of the particle,
33// 2) static particle information,
34// 3) the pointer to the physical volume where currently
35// the particle exists
36
37// Author: Katsuya Amako, KEK - 1995
38// Revisions: Hisaya Kurashige, 1998-2011
39// --------------------------------------------------------------------
40#ifndef G4Track_hh
41#define G4Track_hh 1
42
43#include <cmath> // Include from 'system'
44#include <map>
46
47#include "globals.hh" // Include from 'global'
48#include "trkdefs.hh" // Include DLL defs...
49#include "G4ThreeVector.hh" // Include from 'geometry'
50#include "G4LogicalVolume.hh" // Include from 'geometry'
51#include "G4VPhysicalVolume.hh" // Include from 'geometry'
52#include "G4Allocator.hh" // Include from 'particle+matter'
53#include "G4DynamicParticle.hh" // Include from 'particle+matter'
54#include "G4TrackStatus.hh" // Include from 'tracking'
55#include "G4TouchableHandle.hh" // Include from 'geometry'
58#include "G4Material.hh"
59#include "G4Profiler.hh"
60
61class G4Step; // Forward declaration
64class G4VProcess;
65
67{
68 public:
69 // the profiler aliases are only used when compiled with GEANT4_USE_TIMEMORY
71
72 G4Track();
73 // Default constructor
74 G4Track(G4DynamicParticle* apValueDynamicParticle,
75 G4double aValueTime,
76 const G4ThreeVector& aValuePosition);
77 // Constructor - aValueTime is a global time
78
79 G4Track(const G4Track&);
80 // Copy Constructor - copies members other than tracking information
81
82 ~G4Track();
83 // Destructor
84
85 inline void* operator new(std::size_t);
86 // Override "new" for "G4Allocator".
87 inline void operator delete(void* aTrack);
88 // Override "delete" for "G4Allocator".
89
90 G4Track& operator=(const G4Track&);
91 // Assignment operator
92
93 inline G4bool operator==(const G4Track&);
94 inline G4bool operator!=(const G4Track&);
95 // Equality operators
96
97 void CopyTrackInfo(const G4Track&);
98 // Copy information of the track (w/o tracking information)
99
101 void SetTrackID(const G4int aValue);
102 // Get/Set functions track ID
103
105 void SetParentID(const G4int aValue);
106
108 // Dynamic particle
109
111 // Particle definition
113 // Obsolete, for backwards compatibility...
114
116 void SetPosition(const G4ThreeVector& aValue);
117 // Position, time
118
120 void SetGlobalTime(const G4double aValue);
121 // Time since the event in which the track belongs is created
122
124 void SetLocalTime(const G4double aValue);
125 // Time since the current track is created
126
128 void SetProperTime(const G4double aValue);
129 // Proper time of the current track
130
133 // Volume, material, touchable
134
137
140
144
148
152
154 void SetKineticEnergy(const G4double aValue);
155 // Energy
156
158
161 // Momentum
162
164
167 // Velocity
168
171
174
176 void SetPolarization(const G4ThreeVector& aValue);
177 // Polarization
178
180 void SetTrackStatus(const G4TrackStatus aTrackStatus);
181 // Track status, flags for tracking
182
184 void SetBelowThresholdFlag(G4bool value = true);
185 // The flag of "BelowThreshold" is set to true
186 // If this track energy is below threshold energy
187 // in this material is determined by the range cut value
188
190 void SetGoodForTrackingFlag(G4bool value = true);
191 // The flag of "GoodForTracking" is set by processes
192 // if this track should be tracked
193 // even if the energy is below threshold
194
196 void AddTrackLength(const G4double aValue);
197 // Accumulated track length
198
199 const G4Step* GetStep() const;
200 void SetStep(const G4Step* aValue);
201 // Step information
202
205
208 // Before the end of the AlongStepDoIt() loop, StepLength keeps
209 // the initial value which is determined by the shortest geometrical Step
210 // proposed by a physics process. After finishing the AlongStepDoIt(),
211 // it will be set equal to 'StepLength' in G4Step
212
214 void SetVertexPosition(const G4ThreeVector& aValue);
215 // Vertex (where this track was created) information
216
219
222
225
227 void SetCreatorProcess(const G4VProcess* aValue);
228
229 inline void SetCreatorModelID(const G4int id);
230 inline G4int GetCreatorModelID() const;
232 inline const G4String GetCreatorModelName() const;
233 // Identification of the physics model that created the track:
234 // each of the three information (ID, index, name) is unique
235 // (the model ID and its name are supposed to be used in Geant4
236 // code, whereas the index is meant for plotting in user code)
237
239 inline void SetParentResonanceDef(const G4ParticleDefinition* parent);
241 inline void SetParentResonanceID(const G4int parentID );
246 // Because short-lived resonances (e.g. omega, phi, rho, delta, etc.)
247 // do not have corresponding track objects, if the track is produced
248 // by a resonance parent, these methods allow to get/set information
249 // regarding this short-lived parent.
250 // The ID is a unique (integer) identifier for each resonance (which
251 // corresponds to the rounded integer of the mass of the resonance
252 // in keV), which allows to know if two (or more) tracks originated
253 // from the same parent resonance: this should not be confused with
254 // the parent-track-ID (fParentID) which corresponds to its closest
255 // ancestor which is not a short-lived resonance (and therefore has
256 // a corresponding track object).
257 // In the case of a track non originating from a resonance parent,
258 // the above "Get" methods return, respectively: nullptr, 0, false,
259 // 0, "", 0.
260
262 void SetWeight(G4double aValue);
263 // Track weight; methods for manipulating a weight for this track
264
267 // User information
268
270 G4VAuxiliaryTrackInformation* info) const;
272 inline std::map<G4int, G4VAuxiliaryTrackInformation*>*
274
277 // Note: G4VAuxiliaryTrackInformation object itself is *NOT* deleted
278
279 private:
280
281 void ClearAuxiliaryTrackInformation();
282
283 // Member data
284
285 G4ThreeVector fPosition;
286 // Current positon
287 G4double fGlobalTime = 0.0;
288 // Time since the event is created
289 G4double fLocalTime = 0.0;
290 // Time since the track is created
291 G4double fTrackLength = 0.0;
292 // Accumulated track length
293
294 G4double fVelocity = 0.0;
295
296 G4TouchableHandle fpTouchable;
297 G4TouchableHandle fpNextTouchable;
298 G4TouchableHandle fpOriginTouchable;
299 // Touchable Handle
300
301 G4DynamicParticle* fpDynamicParticle = nullptr;
302 mutable G4TrackStatus fTrackStatus = fAlive;
303
304 G4double fStepLength = 0.0;
305 // Before the end of the AlongStepDoIt loop, this keeps the initial
306 // Step length which is determined by the shortest geometrical Step
307 // proposed by a physics process. After finishing the AlongStepDoIt,
308 // this will be set equal to 'StepLength' in G4Step.
309
310 G4double fWeight = 1.0;
311 // This is a weight for this track
312
313 const G4Step* fpStep = nullptr;
314
315 G4ThreeVector fVtxPosition;
316 // (x,y,z) of the vertex
317 G4ThreeVector fVtxMomentumDirection;
318 // Momentum direction at the vertex
319 G4double fVtxKineticEnergy = 0.0;
320 // Kinetic energy at the vertex
321 const G4LogicalVolume* fpLVAtVertex = nullptr;
322 // Logical Volume at the vertex
323 const G4VProcess* fpCreatorProcess = nullptr;
324 // Process which created the track
325
326 mutable G4VUserTrackInformation* fpUserInformation = nullptr;
327
328 mutable G4Material* prev_mat = nullptr;
329 mutable G4MaterialPropertyVector* groupvel = nullptr;
330 mutable G4double prev_velocity = 0.0;
331 mutable G4double prev_momentum = 0.0;
332 // cached values for CalculateVelocity
333
334 mutable std::map<G4int, G4VAuxiliaryTrackInformation*>*
335 fpAuxiliaryTrackInformationMap = nullptr;
336
337 G4int fCurrentStepNumber = 0;
338 // Total steps number up to now
339
340 G4int fCreatorModelID = -1;
341 // ID of the physics model which created the track
342
343 const G4ParticleDefinition* fParentResonanceDef = nullptr;
344 // Pointer to the particle definition of a short-lived resonance,
345 // in the case that the track is produced by a resonance parent
346 // (which does not have a corresponding track object)
347 G4int fParentResonanceID = 0;
348 // Unique ID for the parent resonance, in the case that the track
349 // is produced by a resonance parent, else 0
350
351 G4int fParentID = 0;
352 G4int fTrackID = 0;
353
354 G4bool fBelowThreshold = false;
355 // This flag is set to true if this track energy is below
356 // threshold energy in this material determined by the range cut value
357 G4bool fGoodForTracking = false;
358 // This flag is set by processes if this track should be tracked
359 // even if the energy is below threshold
360
361 G4bool is_OpticalPhoton = false;
362
363 G4bool useGivenVelocity = false;
364 // do not calculate velocity and just use current fVelocity
365 // if this flag is set
366};
367
368#include "G4Track.icc"
369
370#endif
G4TrackStatus
@ fAlive
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void SetVelocity(G4double val)
G4Material * GetNextMaterial() const
G4TrackStatus GetTrackStatus() const
std::map< G4int, G4VAuxiliaryTrackInformation * > * GetAuxiliaryTrackInformationMap() const
const G4ParticleDefinition * GetParentResonanceDef() const
void SetTrackStatus(const G4TrackStatus aTrackStatus)
void SetStepLength(G4double value)
G4int GetParentResonancePDGEncoding() const
void SetPosition(const G4ThreeVector &aValue)
G4double GetVelocity() const
void SetPolarization(const G4ThreeVector &aValue)
G4double CalculateVelocityForOpticalPhoton() const
Definition G4Track.cc:160
const G4MaterialCutsCouple * GetNextMaterialCutsCouple() const
G4double GetVertexKineticEnergy() const
G4int GetTrackID() const
void SetStep(const G4Step *aValue)
const G4ParticleDefinition * GetParticleDefinition() const
G4VPhysicalVolume * GetVolume() const
void SetAuxiliaryTrackInformation(G4int id, G4VAuxiliaryTrackInformation *info) const
Definition G4Track.cc:209
void SetVertexPosition(const G4ThreeVector &aValue)
const G4TouchableHandle & GetNextTouchableHandle() const
const G4VProcess * GetCreatorProcess() const
void SetVertexMomentumDirection(const G4ThreeVector &aValue)
const G4LogicalVolume * GetLogicalVolumeAtVertex() const
G4VPhysicalVolume * GetNextVolume() const
void SetNextTouchableHandle(const G4TouchableHandle &apValue)
G4bool UseGivenVelocity() const
void CopyTrackInfo(const G4Track &)
Definition G4Track.cc:154
G4double GetWeight() const
G4int GetCreatorModelID() const
void SetWeight(G4double aValue)
void RemoveAuxiliaryTrackInformation(G4int id)
Definition G4Track.cc:241
const G4String GetCreatorModelName() const
void SetParentResonanceID(const G4int parentID)
const G4ThreeVector & GetPosition() const
G4double GetTrackLength() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
const G4ThreeVector & GetVertexMomentumDirection() const
G4double GetGlobalTime() const
G4double GetProperTime() const
G4int GetCurrentStepNumber() const
G4ThreeVector GetMomentum() const
const G4ThreeVector & GetVertexPosition() const
const G4VTouchable * GetOriginTouchable() const
G4double GetLocalTime() const
G4VAuxiliaryTrackInformation * GetAuxiliaryTrackInformation(G4int id) const
Definition G4Track.cc:229
void SetOriginTouchableHandle(const G4TouchableHandle &apValue)
G4Material * GetMaterial() const
G4ParticleDefinition * GetDefinition() const
G4String GetParentResonanceName() const
void SetLocalTime(const G4double aValue)
void AddTrackLength(const G4double aValue)
const G4DynamicParticle * GetDynamicParticle() const
G4bool operator!=(const G4Track &)
const G4TouchableHandle & GetTouchableHandle() const
const G4ThreeVector & GetMomentumDirection() const
G4VUserTrackInformation * GetUserInformation() const
G4double GetKineticEnergy() const
G4double CalculateVelocity() const
const G4ThreeVector & GetPolarization() const
G4int GetParentResonanceID() const
G4Track()
Definition G4Track.cc:62
void SetUserInformation(G4VUserTrackInformation *aValue) const
void SetBelowThresholdFlag(G4bool value=true)
G4int GetCreatorModelIndex() const
void SetCreatorModelID(const G4int id)
G4double GetStepLength() const
void UseGivenVelocity(G4bool val)
void SetVertexKineticEnergy(const G4double aValue)
void SetTrackID(const G4int aValue)
G4int GetParentID() const
void IncrementCurrentStepNumber()
void SetKineticEnergy(const G4double aValue)
G4bool IsBelowThreshold() const
void SetMomentumDirection(const G4ThreeVector &aValue)
const G4TouchableHandle & GetOriginTouchableHandle() const
void SetLogicalVolumeAtVertex(const G4LogicalVolume *)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
const G4VTouchable * GetNextTouchable() const
G4bool operator==(const G4Track &)
G4double GetParentResonanceMass() const
const G4VTouchable * GetTouchable() const
G4double GetTotalEnergy() const
const G4Step * GetStep() const
G4bool IsGoodForTracking() const
G4Track & operator=(const G4Track &)
Definition G4Track.cc:83
void SetProperTime(const G4double aValue)
void SetGlobalTime(const G4double aValue)
~G4Track()
Definition G4Track.cc:75
void SetParentResonanceDef(const G4ParticleDefinition *parent)
void SetParentID(const G4int aValue)
void SetGoodForTrackingFlag(G4bool value=true)
G4bool HasParentResonance() const
void SetCreatorProcess(const G4VProcess *aValue)