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