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