Geant4 Cross Reference |
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> 45 #include <CLHEP/Units/PhysicalConstants.h> 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' 56 #include "G4VUserTrackInformation.hh" 57 #include "G4PhysicsModelCatalog.hh" 58 #include "G4Material.hh" 59 60 class G4Step; // Forward declaration 61 class G4MaterialCutsCouple; 62 class G4VAuxiliaryTrackInformation; 63 class G4VProcess; 64 65 class G4Track 66 { 67 public: 68 G4Track(); 69 // Default constructor 70 G4Track(G4DynamicParticle* apValueDynamicParticle, 71 G4double aValueTime, 72 const G4ThreeVector& aValuePosition); 73 // Constructor - aValueTime is a global time 74 75 G4Track(const G4Track&,G4bool copyTouchables = true); 76 // Copy Constructor - copies members other than tracking information 77 78 ~G4Track(); 79 // Destructor 80 81 inline void* operator new(std::size_t); 82 // Override "new" for "G4Allocator". 83 inline void operator delete(void* aTrack); 84 // Override "delete" for "G4Allocator". 85 86 G4Track& operator=(const G4Track&); 87 // Assignment operator 88 89 inline G4bool operator==(const G4Track&); 90 inline G4bool operator!=(const G4Track&); 91 // Equality operators 92 93 void CopyTrackInfo(const G4Track&, G4bool copyTouchables = true); 94 // Copy information of the track (w/o tracking information) 95 96 inline G4int GetTrackID() const; 97 inline void SetTrackID(const G4int aValue); 98 // Get/Set functions track ID 99 100 inline G4int GetParentID() const; 101 inline void SetParentID(const G4int aValue); 102 103 inline const G4DynamicParticle* GetDynamicParticle() const; 104 // Dynamic particle 105 106 inline const G4ParticleDefinition* GetParticleDefinition() const; 107 // Particle definition 108 inline G4ParticleDefinition* GetDefinition() const; 109 // Obsolete, for backwards compatibility... 110 111 inline const G4ThreeVector& GetPosition() const; 112 inline void SetPosition(const G4ThreeVector& aValue); 113 // Position, time 114 115 inline G4double GetGlobalTime() const; 116 inline void SetGlobalTime(const G4double aValue); 117 // Time since the event in which the track belongs is created 118 119 inline G4double GetLocalTime() const; 120 inline void SetLocalTime(const G4double aValue); 121 // Time since the current track is created 122 123 inline G4double GetProperTime() const; 124 inline void SetProperTime(const G4double aValue); 125 // Proper time of the current track 126 127 inline G4VPhysicalVolume* GetVolume() const; 128 inline G4VPhysicalVolume* GetNextVolume() const; 129 // Volume, material, touchable 130 131 inline G4Material* GetMaterial() const; 132 inline G4Material* GetNextMaterial() const; 133 134 inline const G4MaterialCutsCouple* GetMaterialCutsCouple() const; 135 inline const G4MaterialCutsCouple* GetNextMaterialCutsCouple() const; 136 137 inline const G4VTouchable* GetTouchable() const; 138 inline const G4TouchableHandle& GetTouchableHandle() const; 139 inline void SetTouchableHandle(const G4TouchableHandle& apValue); 140 141 inline const G4VTouchable* GetNextTouchable() const; 142 inline const G4TouchableHandle& GetNextTouchableHandle() const; 143 inline void SetNextTouchableHandle(const G4TouchableHandle& apValue); 144 145 inline const G4VTouchable* GetOriginTouchable() const; 146 inline const G4TouchableHandle& GetOriginTouchableHandle() const; 147 inline void SetOriginTouchableHandle(const G4TouchableHandle& apValue); 148 149 inline G4double GetKineticEnergy() const; 150 inline G4double GetTotalEnergy() const; 151 inline void SetKineticEnergy(const G4double aValue); 152 // Energy 153 154 inline G4ThreeVector GetMomentum() const; 155 inline const G4ThreeVector& GetMomentumDirection() const; 156 inline void SetMomentumDirection(const G4ThreeVector& aValue); 157 // Momentum 158 159 inline G4double GetVelocity() const; 160 inline void SetVelocity(G4double val); 161 // Velocity 162 163 inline G4double CalculateVelocity() const; 164 G4double CalculateVelocityForOpticalPhoton() const; 165 166 inline G4bool UseGivenVelocity() const; 167 inline void UseGivenVelocity(G4bool val); 168 169 inline const G4ThreeVector& GetPolarization() const; 170 inline void SetPolarization(const G4ThreeVector& aValue); 171 // Polarization 172 173 inline G4TrackStatus GetTrackStatus() const; 174 inline void SetTrackStatus(const G4TrackStatus aTrackStatus); 175 // Track status, flags for tracking 176 177 inline G4bool IsBelowThreshold() const; 178 inline void SetBelowThresholdFlag(G4bool value = true); 179 // The flag of "BelowThreshold" is set to true 180 // If this track energy is below threshold energy 181 // in this material is determined by the range cut value 182 183 inline G4bool IsGoodForTracking() const; 184 inline void SetGoodForTrackingFlag(G4bool value = true); 185 // The flag of "GoodForTracking" is set by processes 186 // if this track should be tracked 187 // even if the energy is below threshold 188 189 inline G4double GetTrackLength() const; 190 inline void AddTrackLength(const G4double aValue); 191 // Accumulated track length 192 193 inline const G4Step* GetStep() const; 194 inline void SetStep(const G4Step* aValue); 195 // Step information 196 197 inline G4int GetCurrentStepNumber() const; 198 inline void IncrementCurrentStepNumber(); 199 200 inline G4double GetStepLength() const; 201 inline void SetStepLength(G4double value); 202 // Before the end of the AlongStepDoIt() loop, StepLength keeps 203 // the initial value which is determined by the shortest geometrical Step 204 // proposed by a physics process. After finishing the AlongStepDoIt(), 205 // it will be set equal to 'StepLength' in G4Step 206 207 inline const G4ThreeVector& GetVertexPosition() const; 208 inline void SetVertexPosition(const G4ThreeVector& aValue); 209 // Vertex (where this track was created) information 210 211 inline const G4ThreeVector& GetVertexMomentumDirection() const; 212 inline void SetVertexMomentumDirection(const G4ThreeVector& aValue); 213 214 inline G4double GetVertexKineticEnergy() const; 215 inline void SetVertexKineticEnergy(const G4double aValue); 216 217 inline const G4LogicalVolume* GetLogicalVolumeAtVertex() const; 218 inline void SetLogicalVolumeAtVertex(const G4LogicalVolume*); 219 220 inline const G4VProcess* GetCreatorProcess() const; 221 inline void SetCreatorProcess(const G4VProcess* aValue); 222 223 inline void SetCreatorModelID(const G4int id); 224 inline G4int GetCreatorModelID() const; 225 inline G4int GetCreatorModelIndex() const; 226 inline const G4String GetCreatorModelName() const; 227 // Identification of the physics model that created the track: 228 // each of the three information (ID, index, name) is unique 229 // (the model ID and its name are supposed to be used in Geant4 230 // code, whereas the index is meant for plotting in user code) 231 232 inline const G4ParticleDefinition* GetParentResonanceDef() const; 233 inline void SetParentResonanceDef(const G4ParticleDefinition* parent); 234 inline G4int GetParentResonanceID() const; 235 inline void SetParentResonanceID(const G4int parentID ); 236 inline G4bool HasParentResonance() const; 237 inline G4int GetParentResonancePDGEncoding() const; 238 inline G4String GetParentResonanceName() const; 239 inline G4double GetParentResonanceMass() const; 240 // Because short-lived resonances (e.g. omega, phi, rho, delta, etc.) 241 // do not have corresponding track objects, if the track is produced 242 // by a resonance parent, these methods allow to get/set information 243 // regarding this short-lived parent. 244 // The ID is a unique (integer) identifier for each resonance (which 245 // corresponds to the rounded integer of the mass of the resonance 246 // in keV), which allows to know if two (or more) tracks originated 247 // from the same parent resonance: this should not be confused with 248 // the parent-track-ID (fParentID) which corresponds to its closest 249 // ancestor which is not a short-lived resonance (and therefore has 250 // a corresponding track object). 251 // In the case of a track non originating from a resonance parent, 252 // the above "Get" methods return, respectively: nullptr, 0, false, 253 // 0, "", 0. 254 255 inline G4double GetWeight() const; 256 inline void SetWeight(G4double aValue); 257 // Track weight; methods for manipulating a weight for this track 258 259 inline G4VUserTrackInformation* GetUserInformation() const; 260 inline void SetUserInformation(G4VUserTrackInformation* aValue) const; 261 // User information 262 263 void SetAuxiliaryTrackInformation(G4int id, 264 G4VAuxiliaryTrackInformation* info) const; 265 G4VAuxiliaryTrackInformation* GetAuxiliaryTrackInformation(G4int id) const; 266 inline std::map<G4int, G4VAuxiliaryTrackInformation*>* 267 GetAuxiliaryTrackInformationMap() const; 268 269 void RemoveAuxiliaryTrackInformation(G4int id); 270 void RemoveAuxiliaryTrackInformation(G4String& name); 271 // Note: G4VAuxiliaryTrackInformation object itself is *NOT* deleted 272 273 private: 274 275 void ClearAuxiliaryTrackInformation(); 276 277 // Member data 278 279 G4ThreeVector fPosition; 280 // Current positon 281 G4double fGlobalTime = 0.0; 282 // Time since the event is created 283 G4double fLocalTime = 0.0; 284 // Time since the track is created 285 G4double fTrackLength = 0.0; 286 // Accumulated track length 287 288 G4double fVelocity = 0.0; 289 290 G4TouchableHandle fpTouchable; 291 G4TouchableHandle fpNextTouchable; 292 G4TouchableHandle fpOriginTouchable; 293 // Touchable Handle 294 295 G4DynamicParticle* fpDynamicParticle = nullptr; 296 mutable G4TrackStatus fTrackStatus = fAlive; 297 298 G4double fStepLength = 0.0; 299 // Before the end of the AlongStepDoIt loop, this keeps the initial 300 // Step length which is determined by the shortest geometrical Step 301 // proposed by a physics process. After finishing the AlongStepDoIt, 302 // this will be set equal to 'StepLength' in G4Step. 303 304 G4double fWeight = 1.0; 305 // This is a weight for this track 306 307 const G4Step* fpStep = nullptr; 308 309 G4ThreeVector fVtxPosition; 310 // (x,y,z) of the vertex 311 G4ThreeVector fVtxMomentumDirection; 312 // Momentum direction at the vertex 313 G4double fVtxKineticEnergy = 0.0; 314 // Kinetic energy at the vertex 315 const G4LogicalVolume* fpLVAtVertex = nullptr; 316 // Logical Volume at the vertex 317 const G4VProcess* fpCreatorProcess = nullptr; 318 // Process which created the track 319 320 mutable G4VUserTrackInformation* fpUserInformation = nullptr; 321 322 mutable G4Material* prev_mat = nullptr; 323 mutable G4MaterialPropertyVector* groupvel = nullptr; 324 mutable G4double prev_velocity = 0.0; 325 mutable G4double prev_momentum = 0.0; 326 // cached values for CalculateVelocity 327 328 mutable std::map<G4int, G4VAuxiliaryTrackInformation*>* 329 fpAuxiliaryTrackInformationMap = nullptr; 330 331 G4int fCurrentStepNumber = 0; 332 // Total steps number up to now 333 334 G4int fCreatorModelID = -1; 335 // ID of the physics model which created the track 336 337 const G4ParticleDefinition* fParentResonanceDef = nullptr; 338 // Pointer to the particle definition of a short-lived resonance, 339 // in the case that the track is produced by a resonance parent 340 // (which does not have a corresponding track object) 341 G4int fParentResonanceID = 0; 342 // Unique ID for the parent resonance, in the case that the track 343 // is produced by a resonance parent, else 0 344 345 G4int fParentID = 0; 346 G4int fTrackID = 0; 347 348 G4bool fBelowThreshold = false; 349 // This flag is set to true if this track energy is below 350 // threshold energy in this material determined by the range cut value 351 G4bool fGoodForTracking = false; 352 // This flag is set by processes if this track should be tracked 353 // even if the energy is below threshold 354 355 G4bool is_OpticalPhoton = false; 356 357 G4bool useGivenVelocity = false; 358 // do not calculate velocity and just use current fVelocity 359 // if this flag is set 360 361 G4bool fCopyTouchables = true; 362 }; 363 364 #include "G4Track.icc" 365 366 #endif 367