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