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