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 // G4RichTrajectory 27 // 28 // Class description: 29 // 30 // This class extends G4Trajectory, which includes the following: 31 // 1) List of trajectory points which compose the trajectory, 32 // 2) static information of particle which generated the 33 // trajectory, 34 // 3) trackID and parent particle ID of the trajectory. 35 // The extended information, only publicly accessible through AttValues, 36 // includes: 37 // 4) physical volume and next physical volume; 38 // 5) creator process; 39 // ...and much more. 40 41 // Contact: 42 // Questions and comments on G4Trajectory should be sent to 43 // Katsuya Amako (e-mail: Katsuya.Amako@kek.jp) 44 // Makoto Asai (e-mail: asai@slac.stanford.edu) 45 // Takashi Sasaki (e-mail: Takashi.Sasaki@kek.jp) 46 // and on the extended code to: 47 // John Allison (e-mail: John.Allison@manchester.ac.uk) 48 // Joseph Perl (e-mail: perl@slac.stanford.edu) 49 // -------------------------------------------------------------------- 50 #ifndef G4RICHTRAJECTORY_HH 51 #define G4RICHTRAJECTORY_HH 52 53 #include "G4TouchableHandle.hh" 54 #include "G4VTrajectory.hh" 55 #include "G4ParticleDefinition.hh" // Include from 'particle+matter' 56 #include "G4Step.hh" 57 #include "G4Track.hh" 58 #include "G4RichTrajectoryPoint.hh" // Include from 'tracking' 59 #include "G4ios.hh" // Include from 'system' 60 #include "globals.hh" // Include from 'global' 61 62 #include "trkgdefs.hh" 63 #include <stdlib.h> // Include from 'system' 64 65 #include <vector> 66 67 class G4Polyline; 68 class G4ClonedRichTrajectory; 69 70 class G4RichTrajectory : public G4VTrajectory 71 { 72 using G4TrajectoryPointContainer = std::vector<G4VTrajectoryPoint*>; 73 74 friend class G4ClonedRichTrajectory; 75 76 public: 77 // Constructors/destructor 78 // 79 G4RichTrajectory() = default; 80 G4RichTrajectory(const G4Track* aTrack); 81 ~G4RichTrajectory() override; 82 G4RichTrajectory(G4RichTrajectory&); 83 G4RichTrajectory& operator=(const G4RichTrajectory&) = delete; 84 85 // Operators 86 // 87 G4bool operator==(const G4RichTrajectory& r) const { return (this == &r); } 88 89 inline void* operator new(size_t); 90 inline void operator delete(void*); 91 92 // cloning with the master thread allocator 93 G4VTrajectory* CloneForMaster() const override; 94 95 // Get/Set functions 96 97 inline G4int GetTrackID() const override { return fTrackID; } 98 inline G4int GetParentID() const override { return fParentID; } 99 inline G4String GetParticleName() const override { return ParticleName; } 100 inline G4double GetCharge() const override { return PDGCharge; } 101 inline G4int GetPDGEncoding() const override { return PDGEncoding; } 102 inline G4double GetInitialKineticEnergy() const { return initialKineticEnergy; } 103 inline G4ThreeVector GetInitialMomentum() const override { return initialMomentum; } 104 105 // Other (virtual) member functions 106 // 107 void ShowTrajectory(std::ostream& os = G4cout) const override; 108 void DrawTrajectory() const override; 109 void AppendStep(const G4Step* aStep) override; 110 void MergeTrajectory(G4VTrajectory* secondTrajectory) override; 111 inline G4int GetPointEntries() const override; 112 inline G4VTrajectoryPoint* GetPoint(G4int i) const override; 113 114 G4ParticleDefinition* GetParticleDefinition(); 115 116 // Get methods for HepRep style attributes 117 // 118 const std::map<G4String, G4AttDef>* GetAttDefs() const override; 119 std::vector<G4AttValue>* CreateAttValues() const override; 120 121 private: 122 G4TrajectoryPointContainer* positionRecord = nullptr; 123 G4int fTrackID = 0; 124 G4int fParentID = 0; 125 G4int PDGEncoding = 0; 126 G4double PDGCharge = 0.0; 127 G4String ParticleName = "dummy"; 128 G4double initialKineticEnergy = 0.0; 129 G4ThreeVector initialMomentum; 130 131 // Extended information (only publicly accessible through AttValues)... 132 // 133 G4TrajectoryPointContainer* fpRichPointContainer = nullptr; 134 G4TouchableHandle fpInitialVolume; 135 G4TouchableHandle fpInitialNextVolume; 136 const G4VProcess* fpCreatorProcess = nullptr; 137 G4int fCreatorModelID = 0; 138 G4TouchableHandle fpFinalVolume; 139 G4TouchableHandle fpFinalNextVolume; 140 const G4VProcess* fpEndingProcess = nullptr; 141 G4double fFinalKineticEnergy = 0.0; 142 }; 143 144 extern G4TRACKING_DLL G4Allocator<G4RichTrajectory>*& aRichTrajectoryAllocator(); 145 146 inline void* G4RichTrajectory::operator new(size_t) 147 { 148 if (aRichTrajectoryAllocator() == nullptr) { 149 aRichTrajectoryAllocator() = new G4Allocator<G4RichTrajectory>; 150 } 151 return (void*)aRichTrajectoryAllocator()->MallocSingle(); 152 } 153 154 inline void G4RichTrajectory::operator delete(void* aRichTrajectory) 155 { 156 aRichTrajectoryAllocator()->FreeSingle((G4RichTrajectory*)aRichTrajectory); 157 } 158 159 inline G4int G4RichTrajectory::GetPointEntries() const 160 { 161 return G4int(fpRichPointContainer->size()); 162 } 163 164 inline G4VTrajectoryPoint* G4RichTrajectory::GetPoint(G4int i) const 165 { 166 return (*fpRichPointContainer)[i]; 167 } 168 169 #endif 170