File indexing completed on 2025-01-18 09:59:05
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
0064 class G4SmoothTrajectory : public G4VTrajectory
0065 {
0066 using G4TrajectoryPointContainer = std::vector<G4VTrajectoryPoint*>;
0067
0068 public:
0069
0070
0071 G4SmoothTrajectory() = default;
0072 G4SmoothTrajectory(const G4Track* aTrack);
0073 ~G4SmoothTrajectory() override;
0074 G4SmoothTrajectory(G4SmoothTrajectory&);
0075 G4SmoothTrajectory& operator=(const G4SmoothTrajectory&) = delete;
0076
0077
0078
0079 inline G4bool operator==(const G4SmoothTrajectory& r) const;
0080 inline void* operator new(size_t);
0081 inline void operator delete(void*);
0082
0083
0084
0085 inline G4int GetTrackID() const override { return fTrackID; }
0086 inline G4int GetParentID() const override { return fParentID; }
0087 inline G4String GetParticleName() const override { return ParticleName; }
0088 inline G4double GetCharge() const override { return PDGCharge; }
0089 inline G4int GetPDGEncoding() const override { return PDGEncoding; }
0090 inline G4double GetInitialKineticEnergy() const { return initialKineticEnergy; }
0091 inline G4ThreeVector GetInitialMomentum() const override { return initialMomentum; }
0092
0093
0094
0095 void ShowTrajectory(std::ostream& os = G4cout) const override;
0096 void DrawTrajectory() const override;
0097 void AppendStep(const G4Step* aStep) override;
0098 G4int GetPointEntries() const override { return (G4int)positionRecord->size(); }
0099 G4VTrajectoryPoint* GetPoint(G4int i) const override { return (*positionRecord)[i]; }
0100 void MergeTrajectory(G4VTrajectory* secondTrajectory) override;
0101
0102 G4ParticleDefinition* GetParticleDefinition();
0103
0104
0105
0106 const std::map<G4String, G4AttDef>* GetAttDefs() const override;
0107 std::vector<G4AttValue>* CreateAttValues() const override;
0108
0109 private:
0110 G4TrajectoryPointContainer* positionRecord = nullptr;
0111 G4int fTrackID = 0;
0112 G4int fParentID = 0;
0113 G4int PDGEncoding = 0;
0114 G4double PDGCharge = 0.0;
0115 G4String ParticleName = "";
0116 G4double initialKineticEnergy = 0.0;
0117 G4ThreeVector initialMomentum;
0118 };
0119
0120 extern G4TRACKING_DLL G4Allocator<G4SmoothTrajectory>*& aSmoothTrajectoryAllocator();
0121
0122 inline void* G4SmoothTrajectory::operator new(size_t)
0123 {
0124 if (aSmoothTrajectoryAllocator() == nullptr) {
0125 aSmoothTrajectoryAllocator() = new G4Allocator<G4SmoothTrajectory>;
0126 }
0127 return (void*)aSmoothTrajectoryAllocator()->MallocSingle();
0128 }
0129
0130 inline void G4SmoothTrajectory::operator delete(void* aTrajectory)
0131 {
0132 aSmoothTrajectoryAllocator()->FreeSingle((G4SmoothTrajectory*)aTrajectory);
0133 }
0134
0135 inline G4bool G4SmoothTrajectory::operator==(const G4SmoothTrajectory& r) const
0136 {
0137 return (this == &r);
0138 }
0139
0140 #endif