Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 // G4Trajectory 27 // 28 // Class description: 29 // 30 // This class represents the trajectory of a particle being tracked. 31 // It includes information of: 32 // 1) List of trajectory points which compose the trajectory; 33 // 2) Static information of the particle which generated the 34 // trajectory; 35 // 3) Track ID and parent particle ID of the trajectory. 36 37 // Contact: 38 // Questions and comments to this code should be sent to 39 // Katsuya Amako (e-mail: Katsuya.Amako@kek.jp) 40 // Makoto Asai (e-mail: asai@slac.stanford.edu) 41 // Takashi Sasaki (e-mail: Takashi.Sasaki@kek.jp) 42 // -------------------------------------------------------------------- 43 #ifndef G4Trajectory_hh 44 #define G4Trajectory_hh 1 45 46 #include "G4Allocator.hh" 47 #include "G4ParticleDefinition.hh" // Include from 'particle+matter' 48 #include "G4Step.hh" 49 #include "G4Track.hh" 50 #include "G4TrajectoryPoint.hh" // Include from 'tracking' 51 #include "G4VTrajectory.hh" 52 #include "G4ios.hh" // Include from 'system' 53 #include "globals.hh" // Include from 'global' 54 #include "G4Threading.hh" 55 56 #include "trkgdefs.hh" 57 #include <stdlib.h> // Include from 'system' 58 59 #include <vector> 60 61 class G4Polyline; 62 class G4ClonedTrajectory; 63 64 class G4Trajectory : public G4VTrajectory 65 { 66 using G4TrajectoryPointContainer = std::vector<G4VTrajectoryPoint*>; 67 68 friend class G4ClonedTrajectory; 69 70 public: 71 // Constructors/Destructor 72 73 G4Trajectory() = default; 74 G4Trajectory(const G4Track* aTrack); 75 G4Trajectory(G4Trajectory&); 76 ~G4Trajectory() override; 77 78 // Operators 79 80 inline void* operator new(size_t); 81 inline void operator delete(void*); 82 inline G4bool operator==(const G4Trajectory& r) const; 83 84 // cloning with the master thread allocator 85 G4VTrajectory* CloneForMaster() const override; 86 87 // Get/Set functions 88 89 inline G4int GetTrackID() const override { return fTrackID; } 90 inline G4int GetParentID() const override { return fParentID; } 91 inline G4String GetParticleName() const override { return ParticleName; } 92 inline G4double GetCharge() const override { return PDGCharge; } 93 inline G4int GetPDGEncoding() const override { return PDGEncoding; } 94 inline G4double GetInitialKineticEnergy() const { return initialKineticEnergy; } 95 inline G4ThreeVector GetInitialMomentum() const override { return initialMomentum; } 96 97 // Other member functions 98 99 void ShowTrajectory(std::ostream& os = G4cout) const override; 100 void DrawTrajectory() const override; 101 void AppendStep(const G4Step* aStep) override; 102 G4int GetPointEntries() const override { return G4int(positionRecord->size()); } 103 G4VTrajectoryPoint* GetPoint(G4int i) const override { return (*positionRecord)[i]; } 104 void MergeTrajectory(G4VTrajectory* secondTrajectory) override; 105 106 G4ParticleDefinition* GetParticleDefinition(); 107 108 const std::map<G4String, G4AttDef>* GetAttDefs() const override; 109 std::vector<G4AttValue>* CreateAttValues() const override; 110 111 private: 112 G4TrajectoryPointContainer* positionRecord = nullptr; 113 G4int fTrackID = 0; 114 G4int fParentID = 0; 115 G4int PDGEncoding = 0; 116 G4double PDGCharge = 0.0; 117 G4String ParticleName = "dummy"; 118 G4double initialKineticEnergy = 0.0; 119 G4ThreeVector initialMomentum; 120 }; 121 122 extern G4TRACKING_DLL G4Allocator<G4Trajectory>*& aTrajectoryAllocator(); 123 124 inline void* G4Trajectory::operator new(size_t) 125 { 126 if (aTrajectoryAllocator() == nullptr) { 127 aTrajectoryAllocator() = new G4Allocator<G4Trajectory>; 128 } 129 return (void*)aTrajectoryAllocator()->MallocSingle(); 130 } 131 132 inline void G4Trajectory::operator delete(void* aTrajectory) 133 { 134 aTrajectoryAllocator()->FreeSingle((G4Trajectory*)aTrajectory); 135 } 136 137 inline G4bool G4Trajectory::operator==(const G4Trajectory& r) const { return (this == &r); } 138 139 #endif 140