File indexing completed on 2025-09-15 08:59:32
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 #ifndef G4Trajectory_hh
0044 #define G4Trajectory_hh 1
0045
0046 #include "G4Allocator.hh"
0047 #include "G4ParticleDefinition.hh" // Include from 'particle+matter'
0048 #include "G4Step.hh"
0049 #include "G4Track.hh"
0050 #include "G4TrajectoryPoint.hh" // Include from 'tracking'
0051 #include "G4VTrajectory.hh"
0052 #include "G4ios.hh" // Include from 'system'
0053 #include "globals.hh" // Include from 'global'
0054 #include "G4Threading.hh"
0055
0056 #include "trkgdefs.hh"
0057 #include <stdlib.h> // Include from 'system'
0058
0059 #include <vector>
0060
0061 class G4Polyline;
0062 class G4ClonedTrajectory;
0063
0064 class G4Trajectory : public G4VTrajectory
0065 {
0066 using G4TrajectoryPointContainer = std::vector<G4VTrajectoryPoint*>;
0067
0068 friend class G4ClonedTrajectory;
0069
0070 public:
0071
0072
0073 G4Trajectory() = default;
0074 G4Trajectory(const G4Track* aTrack);
0075 G4Trajectory(G4Trajectory&);
0076 ~G4Trajectory() override;
0077
0078
0079
0080 inline void* operator new(size_t);
0081 inline void operator delete(void*);
0082 inline G4bool operator==(const G4Trajectory& r) const;
0083
0084
0085 G4VTrajectory* CloneForMaster() const override;
0086
0087
0088
0089 inline G4int GetTrackID() const override { return fTrackID; }
0090 inline G4int GetParentID() const override { return fParentID; }
0091 inline G4String GetParticleName() const override { return ParticleName; }
0092 inline G4double GetCharge() const override { return PDGCharge; }
0093 inline G4int GetPDGEncoding() const override { return PDGEncoding; }
0094 inline G4double GetInitialKineticEnergy() const { return initialKineticEnergy; }
0095 inline G4ThreeVector GetInitialMomentum() const override { return initialMomentum; }
0096
0097
0098
0099 void ShowTrajectory(std::ostream& os = G4cout) const override;
0100 void DrawTrajectory() const override;
0101 void AppendStep(const G4Step* aStep) override;
0102 G4int GetPointEntries() const override { return G4int(positionRecord->size()); }
0103 G4VTrajectoryPoint* GetPoint(G4int i) const override { return (*positionRecord)[i]; }
0104 void MergeTrajectory(G4VTrajectory* secondTrajectory) override;
0105
0106 G4ParticleDefinition* GetParticleDefinition();
0107
0108 const std::map<G4String, G4AttDef>* GetAttDefs() const override;
0109 std::vector<G4AttValue>* CreateAttValues() const override;
0110
0111 private:
0112 G4TrajectoryPointContainer* positionRecord = nullptr;
0113 G4int fTrackID = 0;
0114 G4int fParentID = 0;
0115 G4int PDGEncoding = 0;
0116 G4double PDGCharge = 0.0;
0117 G4String ParticleName = "dummy";
0118 G4double initialKineticEnergy = 0.0;
0119 G4ThreeVector initialMomentum;
0120 };
0121
0122 extern G4TRACKING_DLL G4Allocator<G4Trajectory>*& aTrajectoryAllocator();
0123
0124 inline void* G4Trajectory::operator new(size_t)
0125 {
0126 if (aTrajectoryAllocator() == nullptr) {
0127 aTrajectoryAllocator() = new G4Allocator<G4Trajectory>;
0128 }
0129 return (void*)aTrajectoryAllocator()->MallocSingle();
0130 }
0131
0132 inline void G4Trajectory::operator delete(void* aTrajectory)
0133 {
0134 aTrajectoryAllocator()->FreeSingle((G4Trajectory*)aTrajectory);
0135 }
0136
0137 inline G4bool G4Trajectory::operator==(const G4Trajectory& r) const { return (this == &r); }
0138
0139 #endif