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 // G4MCTSimVertex 27 28 // Author: Youhei Morita, 12.09.2001 29 // -------------------------------------------------------------------- 30 #ifndef G4MCTSIMVERTEX_HH 31 #define G4MCTSIMVERTEX_HH 1 32 33 #include <iostream> 34 #include <vector> 35 36 #include "G4Types.hh" 37 #include "G4ThreeVector.hh" 38 #include "G4MCTSimParticle.hh" 39 40 class G4MCTSimVertex 41 { 42 public: 43 44 G4MCTSimVertex(); 45 G4MCTSimVertex(const G4ThreeVector& x, G4double t); 46 G4MCTSimVertex(const G4ThreeVector& x, G4double t, 47 const G4String& vname, G4int ncopy, 48 const G4String& pname); 49 ~G4MCTSimVertex(); 50 51 inline G4MCTSimVertex(const G4MCTSimVertex& right); 52 inline G4MCTSimVertex& operator=(const G4MCTSimVertex& right); 53 // copy constructor and assignment operator 54 55 inline void SetID(G4int i); 56 inline G4int GetID() const; 57 58 inline void SetPosition(const G4ThreeVector& x); 59 inline const G4ThreeVector& GetPosition() const; 60 61 inline void SetTime(G4double t); 62 inline G4double GetTime() const; 63 64 inline void SetVolumeName(const G4String& vname); 65 inline const G4String& GetVolumeName() const; 66 67 inline void SetVolumeNumber(G4int n); 68 inline G4int GetVolumeNumber() const; 69 70 inline void SetCreatorProcessName(const G4String& pname); 71 inline const G4String& GetCreatorProcessName() const; 72 73 inline void SetStoreFlag(G4bool q); 74 inline G4bool GetStoreFlag() const; 75 76 inline void SetInParticle(const G4MCTSimParticle* in); 77 inline void SetInParticle(G4int in); 78 inline G4int GetInParticleTrackID() const; 79 80 inline G4int GetNofOutParticles() const; 81 inline G4int AddOutParticle(const G4MCTSimParticle* out); 82 inline G4int AddOutParticle(G4int out); 83 inline G4int GetOutParticleTrackID(G4int i) const; 84 85 void Print(std::ostream& ostr = std::cout) const; 86 87 private: 88 89 G4int inParticleTrackID = 0; 90 std::vector<G4int> outParticleTrackIDList; 91 92 G4String volumeName = ""; 93 G4String creatorProcessName = "none"; 94 G4ThreeVector position; 95 G4double time = 0.0; 96 G4int id = -1; // assigned independently from G4 97 G4int volumeNumber = -1; 98 G4bool storeFlag = false; 99 }; 100 101 // ==================================================================== 102 // inline methods 103 // ==================================================================== 104 105 inline G4MCTSimVertex::G4MCTSimVertex(const G4MCTSimVertex& right) 106 { 107 *this = right; 108 } 109 110 inline G4MCTSimVertex& G4MCTSimVertex::operator=( 111 const G4MCTSimVertex& right) 112 { 113 inParticleTrackID = right.inParticleTrackID; 114 outParticleTrackIDList = right.outParticleTrackIDList; 115 116 id = right.id; 117 position = right.position; 118 time = right.time; 119 volumeName = right.volumeName; 120 volumeNumber = right.volumeNumber; 121 creatorProcessName = right.creatorProcessName; 122 123 return *this; 124 } 125 126 inline void G4MCTSimVertex::SetID(G4int i) 127 { 128 id = i; 129 } 130 131 inline G4int G4MCTSimVertex::GetID() const 132 { 133 return id; 134 } 135 136 inline void G4MCTSimVertex::SetPosition(const G4ThreeVector& x) 137 { 138 position = x; 139 } 140 141 inline const G4ThreeVector& G4MCTSimVertex::GetPosition() const 142 { 143 return position; 144 } 145 146 inline void G4MCTSimVertex::SetTime(G4double t) 147 { 148 time = t; 149 } 150 151 inline G4double G4MCTSimVertex::GetTime() const 152 { 153 return time; 154 } 155 156 inline void G4MCTSimVertex::SetVolumeName(const G4String& vname) 157 { 158 volumeName = vname; 159 } 160 161 inline const G4String& G4MCTSimVertex::GetVolumeName() const 162 { 163 return volumeName; 164 } 165 166 inline void G4MCTSimVertex::SetVolumeNumber(G4int n) 167 { 168 volumeNumber = n; 169 } 170 171 inline G4int G4MCTSimVertex::GetVolumeNumber() const 172 { 173 return volumeNumber; 174 } 175 176 inline void G4MCTSimVertex::SetCreatorProcessName(const G4String& pname) 177 { 178 creatorProcessName = pname; 179 } 180 181 inline const G4String& G4MCTSimVertex::GetCreatorProcessName() const 182 { 183 return creatorProcessName; 184 } 185 186 inline void G4MCTSimVertex::SetStoreFlag(G4bool q) 187 { 188 storeFlag = q; 189 } 190 191 inline G4bool G4MCTSimVertex::GetStoreFlag() const 192 { 193 return storeFlag; 194 } 195 196 inline void G4MCTSimVertex::SetInParticle(const G4MCTSimParticle* in) 197 { 198 inParticleTrackID = in->GetTrackID(); 199 } 200 201 inline void G4MCTSimVertex::SetInParticle(G4int in) 202 { 203 inParticleTrackID = in; 204 } 205 206 inline G4int G4MCTSimVertex::GetInParticleTrackID() const 207 { 208 return inParticleTrackID; 209 } 210 211 inline G4int G4MCTSimVertex::GetNofOutParticles() const 212 { 213 return (G4int)outParticleTrackIDList.size(); 214 } 215 216 inline G4int G4MCTSimVertex::AddOutParticle(const G4MCTSimParticle* out) 217 { 218 outParticleTrackIDList.push_back(out->GetTrackID()); 219 return (G4int)outParticleTrackIDList.size(); 220 } 221 222 inline G4int G4MCTSimVertex::AddOutParticle(G4int out) 223 { 224 outParticleTrackIDList.push_back(out); 225 return (G4int)outParticleTrackIDList.size(); 226 } 227 228 inline G4int G4MCTSimVertex::GetOutParticleTrackID(G4int i) const 229 { 230 G4int size = (G4int)outParticleTrackIDList.size(); 231 if(i >= 0 && i < size) 232 return outParticleTrackIDList[i]; 233 else 234 return 0; 235 } 236 237 #endif 238