File indexing completed on 2025-01-18 09:59:13
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040 #ifndef G4Track_hh
0041 #define G4Track_hh 1
0042
0043 #include <cmath> // Include from 'system'
0044 #include <map>
0045 #include <CLHEP/Units/PhysicalConstants.h>
0046
0047 #include "globals.hh" // Include from 'global'
0048 #include "trkdefs.hh" // Include DLL defs...
0049 #include "G4ThreeVector.hh" // Include from 'geometry'
0050 #include "G4LogicalVolume.hh" // Include from 'geometry'
0051 #include "G4VPhysicalVolume.hh" // Include from 'geometry'
0052 #include "G4Allocator.hh" // Include from 'particle+matter'
0053 #include "G4DynamicParticle.hh" // Include from 'particle+matter'
0054 #include "G4TrackStatus.hh" // Include from 'tracking'
0055 #include "G4TouchableHandle.hh" // Include from 'geometry'
0056 #include "G4VUserTrackInformation.hh"
0057 #include "G4PhysicsModelCatalog.hh"
0058 #include "G4Material.hh"
0059 #include "G4Profiler.hh"
0060
0061 class G4Step;
0062 class G4MaterialCutsCouple;
0063 class G4VAuxiliaryTrackInformation;
0064 class G4VProcess;
0065
0066 class G4Track
0067 {
0068 public:
0069
0070 using ProfilerConfig = G4ProfilerConfig<G4ProfileType::Track>;
0071
0072 G4Track();
0073
0074 G4Track(G4DynamicParticle* apValueDynamicParticle,
0075 G4double aValueTime,
0076 const G4ThreeVector& aValuePosition);
0077
0078
0079 G4Track(const G4Track&);
0080
0081
0082 ~G4Track();
0083
0084
0085 inline void* operator new(std::size_t);
0086
0087 inline void operator delete(void* aTrack);
0088
0089
0090 G4Track& operator=(const G4Track&);
0091
0092
0093 inline G4bool operator==(const G4Track&);
0094 inline G4bool operator!=(const G4Track&);
0095
0096
0097 void CopyTrackInfo(const G4Track&);
0098
0099
0100 G4int GetTrackID() const;
0101 void SetTrackID(const G4int aValue);
0102
0103
0104 G4int GetParentID() const;
0105 void SetParentID(const G4int aValue);
0106
0107 const G4DynamicParticle* GetDynamicParticle() const;
0108
0109
0110 const G4ParticleDefinition* GetParticleDefinition() const;
0111
0112 G4ParticleDefinition* GetDefinition() const;
0113
0114
0115 const G4ThreeVector& GetPosition() const;
0116 void SetPosition(const G4ThreeVector& aValue);
0117
0118
0119 G4double GetGlobalTime() const;
0120 void SetGlobalTime(const G4double aValue);
0121
0122
0123 G4double GetLocalTime() const;
0124 void SetLocalTime(const G4double aValue);
0125
0126
0127 G4double GetProperTime() const;
0128 void SetProperTime(const G4double aValue);
0129
0130
0131 G4VPhysicalVolume* GetVolume() const;
0132 G4VPhysicalVolume* GetNextVolume() const;
0133
0134
0135 G4Material* GetMaterial() const;
0136 G4Material* GetNextMaterial() const;
0137
0138 const G4MaterialCutsCouple* GetMaterialCutsCouple() const;
0139 const G4MaterialCutsCouple* GetNextMaterialCutsCouple() const;
0140
0141 const G4VTouchable* GetTouchable() const;
0142 const G4TouchableHandle& GetTouchableHandle() const;
0143 void SetTouchableHandle(const G4TouchableHandle& apValue);
0144
0145 const G4VTouchable* GetNextTouchable() const;
0146 const G4TouchableHandle& GetNextTouchableHandle() const;
0147 void SetNextTouchableHandle(const G4TouchableHandle& apValue);
0148
0149 const G4VTouchable* GetOriginTouchable() const;
0150 const G4TouchableHandle& GetOriginTouchableHandle() const;
0151 void SetOriginTouchableHandle(const G4TouchableHandle& apValue);
0152
0153 G4double GetKineticEnergy() const;
0154 void SetKineticEnergy(const G4double aValue);
0155
0156
0157 G4double GetTotalEnergy() const;
0158
0159 const G4ThreeVector& GetMomentumDirection() const;
0160 void SetMomentumDirection(const G4ThreeVector& aValue);
0161
0162
0163 G4ThreeVector GetMomentum() const;
0164
0165 G4double GetVelocity() const;
0166 void SetVelocity(G4double val);
0167
0168
0169 G4double CalculateVelocity() const;
0170 G4double CalculateVelocityForOpticalPhoton() const;
0171
0172 G4bool UseGivenVelocity() const;
0173 void UseGivenVelocity(G4bool val);
0174
0175 const G4ThreeVector& GetPolarization() const;
0176 void SetPolarization(const G4ThreeVector& aValue);
0177
0178
0179 G4TrackStatus GetTrackStatus() const;
0180 void SetTrackStatus(const G4TrackStatus aTrackStatus);
0181
0182
0183 G4bool IsBelowThreshold() const;
0184 void SetBelowThresholdFlag(G4bool value = true);
0185
0186
0187
0188
0189 G4bool IsGoodForTracking() const;
0190 void SetGoodForTrackingFlag(G4bool value = true);
0191
0192
0193
0194
0195 G4double GetTrackLength() const;
0196 void AddTrackLength(const G4double aValue);
0197
0198
0199 const G4Step* GetStep() const;
0200 void SetStep(const G4Step* aValue);
0201
0202
0203 G4int GetCurrentStepNumber() const;
0204 void IncrementCurrentStepNumber();
0205
0206 G4double GetStepLength() const;
0207 void SetStepLength(G4double value);
0208
0209
0210
0211
0212
0213 const G4ThreeVector& GetVertexPosition() const;
0214 void SetVertexPosition(const G4ThreeVector& aValue);
0215
0216
0217 const G4ThreeVector& GetVertexMomentumDirection() const;
0218 void SetVertexMomentumDirection(const G4ThreeVector& aValue);
0219
0220 G4double GetVertexKineticEnergy() const;
0221 void SetVertexKineticEnergy(const G4double aValue);
0222
0223 const G4LogicalVolume* GetLogicalVolumeAtVertex() const;
0224 void SetLogicalVolumeAtVertex(const G4LogicalVolume*);
0225
0226 const G4VProcess* GetCreatorProcess() const;
0227 void SetCreatorProcess(const G4VProcess* aValue);
0228
0229 inline void SetCreatorModelID(const G4int id);
0230 inline G4int GetCreatorModelID() const;
0231 inline G4int GetCreatorModelIndex() const;
0232 inline const G4String GetCreatorModelName() const;
0233
0234
0235
0236
0237
0238 inline const G4ParticleDefinition* GetParentResonanceDef() const;
0239 inline void SetParentResonanceDef(const G4ParticleDefinition* parent);
0240 inline G4int GetParentResonanceID() const;
0241 inline void SetParentResonanceID(const G4int parentID );
0242 inline G4bool HasParentResonance() const;
0243 inline G4int GetParentResonancePDGEncoding() const;
0244 inline G4String GetParentResonanceName() const;
0245 inline G4double GetParentResonanceMass() const;
0246
0247
0248
0249
0250
0251
0252
0253
0254
0255
0256
0257
0258
0259
0260
0261 G4double GetWeight() const;
0262 void SetWeight(G4double aValue);
0263
0264
0265 G4VUserTrackInformation* GetUserInformation() const;
0266 void SetUserInformation(G4VUserTrackInformation* aValue) const;
0267
0268
0269 void SetAuxiliaryTrackInformation(G4int id,
0270 G4VAuxiliaryTrackInformation* info) const;
0271 G4VAuxiliaryTrackInformation* GetAuxiliaryTrackInformation(G4int id) const;
0272 inline std::map<G4int, G4VAuxiliaryTrackInformation*>*
0273 GetAuxiliaryTrackInformationMap() const;
0274
0275 void RemoveAuxiliaryTrackInformation(G4int id);
0276 void RemoveAuxiliaryTrackInformation(G4String& name);
0277
0278
0279 private:
0280
0281 void ClearAuxiliaryTrackInformation();
0282
0283
0284
0285 G4ThreeVector fPosition;
0286
0287 G4double fGlobalTime = 0.0;
0288
0289 G4double fLocalTime = 0.0;
0290
0291 G4double fTrackLength = 0.0;
0292
0293
0294 G4double fVelocity = 0.0;
0295
0296 G4TouchableHandle fpTouchable;
0297 G4TouchableHandle fpNextTouchable;
0298 G4TouchableHandle fpOriginTouchable;
0299
0300
0301 G4DynamicParticle* fpDynamicParticle = nullptr;
0302 mutable G4TrackStatus fTrackStatus = fAlive;
0303
0304 G4double fStepLength = 0.0;
0305
0306
0307
0308
0309
0310 G4double fWeight = 1.0;
0311
0312
0313 const G4Step* fpStep = nullptr;
0314
0315 G4ThreeVector fVtxPosition;
0316
0317 G4ThreeVector fVtxMomentumDirection;
0318
0319 G4double fVtxKineticEnergy = 0.0;
0320
0321 const G4LogicalVolume* fpLVAtVertex = nullptr;
0322
0323 const G4VProcess* fpCreatorProcess = nullptr;
0324
0325
0326 mutable G4VUserTrackInformation* fpUserInformation = nullptr;
0327
0328 mutable G4Material* prev_mat = nullptr;
0329 mutable G4MaterialPropertyVector* groupvel = nullptr;
0330 mutable G4double prev_velocity = 0.0;
0331 mutable G4double prev_momentum = 0.0;
0332
0333
0334 mutable std::map<G4int, G4VAuxiliaryTrackInformation*>*
0335 fpAuxiliaryTrackInformationMap = nullptr;
0336
0337 G4int fCurrentStepNumber = 0;
0338
0339
0340 G4int fCreatorModelID = -1;
0341
0342
0343 const G4ParticleDefinition* fParentResonanceDef = nullptr;
0344
0345
0346
0347 G4int fParentResonanceID = 0;
0348
0349
0350
0351 G4int fParentID = 0;
0352 G4int fTrackID = 0;
0353
0354 G4bool fBelowThreshold = false;
0355
0356
0357 G4bool fGoodForTracking = false;
0358
0359
0360
0361 G4bool is_OpticalPhoton = false;
0362
0363 G4bool useGivenVelocity = false;
0364
0365
0366 };
0367
0368 #include "G4Track.icc"
0369
0370 #endif