File indexing completed on 2025-09-18 09:15:37
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
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050 #ifndef G4RICHTRAJECTORY_HH
0051 #define G4RICHTRAJECTORY_HH
0052
0053 #include "G4TouchableHandle.hh"
0054 #include "G4VTrajectory.hh"
0055 #include "G4ParticleDefinition.hh" // Include from 'particle+matter'
0056 #include "G4Step.hh"
0057 #include "G4Track.hh"
0058 #include "G4RichTrajectoryPoint.hh" // Include from 'tracking'
0059 #include "G4ios.hh" // Include from 'system'
0060 #include "globals.hh" // Include from 'global'
0061
0062 #include "trkgdefs.hh"
0063 #include <stdlib.h> // Include from 'system'
0064
0065 #include <vector>
0066
0067 class G4Polyline;
0068 class G4ClonedRichTrajectory;
0069
0070 class G4RichTrajectory : public G4VTrajectory
0071 {
0072 using G4TrajectoryPointContainer = std::vector<G4VTrajectoryPoint*>;
0073
0074 friend class G4ClonedRichTrajectory;
0075
0076 public:
0077
0078
0079 G4RichTrajectory() = default;
0080 G4RichTrajectory(const G4Track* aTrack);
0081 ~G4RichTrajectory() override;
0082 G4RichTrajectory(G4RichTrajectory&);
0083 G4RichTrajectory& operator=(const G4RichTrajectory&) = delete;
0084
0085
0086
0087 G4bool operator==(const G4RichTrajectory& r) const { return (this == &r); }
0088
0089 inline void* operator new(size_t);
0090 inline void operator delete(void*);
0091
0092
0093 G4VTrajectory* CloneForMaster() const override;
0094
0095
0096
0097 inline G4int GetTrackID() const override { return fTrackID; }
0098 inline G4int GetParentID() const override { return fParentID; }
0099 inline G4String GetParticleName() const override { return ParticleName; }
0100 inline G4double GetCharge() const override { return PDGCharge; }
0101 inline G4int GetPDGEncoding() const override { return PDGEncoding; }
0102 inline G4double GetInitialKineticEnergy() const { return initialKineticEnergy; }
0103 inline G4ThreeVector GetInitialMomentum() const override { return initialMomentum; }
0104
0105
0106
0107 void ShowTrajectory(std::ostream& os = G4cout) const override;
0108 void DrawTrajectory() const override;
0109 void AppendStep(const G4Step* aStep) override;
0110 void MergeTrajectory(G4VTrajectory* secondTrajectory) override;
0111 inline G4int GetPointEntries() const override;
0112 inline G4VTrajectoryPoint* GetPoint(G4int i) const override;
0113
0114 G4ParticleDefinition* GetParticleDefinition();
0115
0116
0117
0118 const std::map<G4String, G4AttDef>* GetAttDefs() const override;
0119 std::vector<G4AttValue>* CreateAttValues() const override;
0120
0121 private:
0122 G4TrajectoryPointContainer* positionRecord = nullptr;
0123 G4int fTrackID = 0;
0124 G4int fParentID = 0;
0125 G4int PDGEncoding = 0;
0126 G4double PDGCharge = 0.0;
0127 G4String ParticleName = "dummy";
0128 G4double initialKineticEnergy = 0.0;
0129 G4ThreeVector initialMomentum;
0130
0131
0132
0133 G4TrajectoryPointContainer* fpRichPointContainer = nullptr;
0134 G4TouchableHandle fpInitialVolume;
0135 G4TouchableHandle fpInitialNextVolume;
0136 const G4VProcess* fpCreatorProcess = nullptr;
0137 G4int fCreatorModelID = 0;
0138 G4TouchableHandle fpFinalVolume;
0139 G4TouchableHandle fpFinalNextVolume;
0140 const G4VProcess* fpEndingProcess = nullptr;
0141 G4double fFinalKineticEnergy = 0.0;
0142 };
0143
0144 extern G4TRACKING_DLL G4Allocator<G4RichTrajectory>*& aRichTrajectoryAllocator();
0145
0146 inline void* G4RichTrajectory::operator new(size_t)
0147 {
0148 if (aRichTrajectoryAllocator() == nullptr) {
0149 aRichTrajectoryAllocator() = new G4Allocator<G4RichTrajectory>;
0150 }
0151 return (void*)aRichTrajectoryAllocator()->MallocSingle();
0152 }
0153
0154 inline void G4RichTrajectory::operator delete(void* aRichTrajectory)
0155 {
0156 aRichTrajectoryAllocator()->FreeSingle((G4RichTrajectory*)aRichTrajectory);
0157 }
0158
0159 inline G4int G4RichTrajectory::GetPointEntries() const
0160 {
0161 return G4int(fpRichPointContainer->size());
0162 }
0163
0164 inline G4VTrajectoryPoint* G4RichTrajectory::GetPoint(G4int i) const
0165 {
0166 return (*fpRichPointContainer)[i];
0167 }
0168
0169 #endif