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 // G4PrimaryParticle 27 // 28 // Class description: 29 // 30 // This class represents a primary particle. 31 // G4PrimaryParticle is a completely different object from G4Track or 32 // G4DynamicParticle. G4PrimaryParticle is designed to take into account 33 // the possibility of making the object persistent, i.e. kept with G4Event 34 // class object within an ODBMS. This class is therefore almost independent 35 // from any other Geant4 object. The only exception is a pointer to 36 // G4ParticleDefinition, which can be rebuilt by its PDGcode. 37 // 38 // Primary particles are stored in G4PrimaryVertex through a form of linked 39 // list. A G4PrimaryParticle object can have one or more objects of this 40 // class as its daughters as a form of linked list. 41 // A primary particle represented by this class object needs not to be 42 // a particle type which Geant4 can simulate: 43 // Case a) mother particle is not a particle Geant4 can simulate 44 // daughters associated to the mother will be examined. 45 // Case b) mother particle is a perticle Geant4 can simulate 46 // daughters associated to the mother will be converted to 47 // G4DynamicParticle and be set to the mother G4Track. 48 // For this case, daugthers are used as the "pre-fixed" 49 // decay channel. 50 51 // Authors: G.Cosmo, 2 December 1995 - Design, based on object model 52 // M.Asai, 29 January 1996 - First implementation 53 // -------------------------------------------------------------------- 54 #ifndef G4PrimaryParticle_hh 55 #define G4PrimaryParticle_hh 1 56 57 #include "G4Allocator.hh" 58 #include "G4ThreeVector.hh" 59 #include "globals.hh" 60 61 #include "pwdefs.hh" 62 63 class G4ParticleDefinition; 64 class G4VUserPrimaryParticleInformation; 65 66 class G4PrimaryParticle 67 { 68 public: 69 // Constructors 70 G4PrimaryParticle(); 71 G4PrimaryParticle(G4int Pcode); 72 G4PrimaryParticle(G4int Pcode, G4double px, G4double py, G4double pz); 73 G4PrimaryParticle(G4int Pcode, G4double px, G4double py, G4double pz, G4double E); 74 G4PrimaryParticle(const G4ParticleDefinition* Gcode); 75 G4PrimaryParticle(const G4ParticleDefinition* Gcode, G4double px, G4double py, G4double pz); 76 G4PrimaryParticle(const G4ParticleDefinition* Gcode, G4double px, G4double py, G4double pz, 77 G4double E); 78 79 // Destructor 80 virtual ~G4PrimaryParticle(); 81 82 // Copy constructor and assignment operator 83 // NOTE: nextParticle and daughterParticle are copied by object 84 // (i.e. deep copy); userInfo will not be copied 85 G4PrimaryParticle(const G4PrimaryParticle& right); 86 G4PrimaryParticle& operator=(const G4PrimaryParticle& right); 87 88 // Equality operator returns 'true' only if same object 89 // (i.e. comparison by pointer value) 90 G4bool operator==(const G4PrimaryParticle& right) const; 91 G4bool operator!=(const G4PrimaryParticle& right) const; 92 93 // Overloaded new/delete operators 94 inline void* operator new(std::size_t); 95 inline void operator delete(void* aPrimaryParticle); 96 97 // Print the properties of the particle 98 void Print() const; 99 100 // Followings are the available accessors/modifiers. 101 // "trackID" will be set if this particle is sent to G4EventManager 102 // and converted to G4Track. Otherwise = -1. 103 // The mass and charge in G4ParticleDefinition will be used by default. 104 // "SetMass" and "SetCharge" methods are used to set dynamical mass and 105 // charge of G4DynamicParticle. 106 // "GetMass" and "GetCharge" methods will return those in 107 // G4ParticleDefinition if users do not set dynamical ones 108 109 inline G4int GetPDGcode() const; 110 void SetPDGcode(G4int Pcode); 111 inline G4ParticleDefinition* GetG4code() const; 112 inline void SetG4code(const G4ParticleDefinition* Gcode); 113 inline const G4ParticleDefinition* GetParticleDefinition() const; 114 void SetParticleDefinition(const G4ParticleDefinition* pdef); 115 inline G4double GetMass() const; 116 inline void SetMass(G4double mas); 117 inline G4double GetCharge() const; 118 inline void SetCharge(G4double chg); 119 inline G4double GetKineticEnergy() const; 120 inline void SetKineticEnergy(G4double eKin); 121 inline const G4ThreeVector& GetMomentumDirection() const; 122 inline void SetMomentumDirection(const G4ThreeVector& p); 123 inline G4double GetTotalMomentum() const; 124 void Set4Momentum(G4double px, G4double py, G4double pz, G4double E); 125 inline G4double GetTotalEnergy() const; 126 inline void SetTotalEnergy(G4double eTot); 127 inline G4ThreeVector GetMomentum() const; 128 void SetMomentum(G4double px, G4double py, G4double pz); 129 inline G4double GetPx() const; 130 inline G4double GetPy() const; 131 inline G4double GetPz() const; 132 inline G4PrimaryParticle* GetNext() const; 133 inline void SetNext(G4PrimaryParticle* np); 134 inline void ClearNext(); 135 inline G4PrimaryParticle* GetDaughter() const; 136 inline void SetDaughter(G4PrimaryParticle* np); 137 inline G4int GetTrackID() const; 138 inline void SetTrackID(G4int id); 139 inline G4ThreeVector GetPolarization() const; 140 inline void SetPolarization(const G4ThreeVector& pol); 141 inline void SetPolarization(G4double px, G4double py, G4double pz); 142 inline G4double GetPolX() const; 143 inline G4double GetPolY() const; 144 inline G4double GetPolZ() const; 145 inline G4double GetWeight() const; 146 inline void SetWeight(G4double w); 147 inline G4double GetProperTime() const; 148 inline void SetProperTime(G4double t); 149 inline G4VUserPrimaryParticleInformation* GetUserInformation() const; 150 inline void SetUserInformation(G4VUserPrimaryParticleInformation* anInfo); 151 152 private: 153 const G4ParticleDefinition* G4code = nullptr; 154 155 G4ThreeVector direction; 156 G4double kinE = 0.0; 157 158 G4PrimaryParticle* nextParticle = nullptr; 159 G4PrimaryParticle* daughterParticle = nullptr; 160 161 G4double mass = -1.0; 162 G4double charge = 0.0; 163 G4double polX = 0.0; 164 G4double polY = 0.0; 165 G4double polZ = 0.0; 166 G4double Weight0 = 1.0; 167 G4double properTime = -1.0; 168 G4VUserPrimaryParticleInformation* userInfo = nullptr; 169 170 G4int PDGcode = 0; 171 G4int trackID = -1; // This will be set if this particle is 172 // sent to G4EventManager and converted to 173 // G4Track. Otherwise = -1 174 }; 175 176 extern G4PART_DLL G4Allocator<G4PrimaryParticle>*& aPrimaryParticleAllocator(); 177 178 // ------------------------ 179 // Inline methods 180 // ------------------------ 181 182 inline void* G4PrimaryParticle::operator new(std::size_t) 183 { 184 if (aPrimaryParticleAllocator() == nullptr) { 185 aPrimaryParticleAllocator() = new G4Allocator<G4PrimaryParticle>; 186 } 187 return (void*)aPrimaryParticleAllocator()->MallocSingle(); 188 } 189 190 inline void G4PrimaryParticle::operator delete(void* aPrimaryParticle) 191 { 192 aPrimaryParticleAllocator()->FreeSingle((G4PrimaryParticle*)aPrimaryParticle); 193 } 194 195 inline G4double G4PrimaryParticle::GetMass() const 196 { 197 return mass; 198 } 199 200 inline G4double G4PrimaryParticle::GetCharge() const 201 { 202 return charge; 203 } 204 205 inline G4int G4PrimaryParticle::GetPDGcode() const 206 { 207 return PDGcode; 208 } 209 210 inline G4ParticleDefinition* G4PrimaryParticle::GetG4code() const 211 { 212 return const_cast<G4ParticleDefinition*>(G4code); 213 } 214 215 inline const G4ParticleDefinition* G4PrimaryParticle::GetParticleDefinition() const 216 { 217 return G4code; 218 } 219 220 inline G4double G4PrimaryParticle::GetTotalMomentum() const 221 { 222 if (mass < 0.) return kinE; 223 return std::sqrt(kinE * (kinE + 2. * mass)); 224 } 225 226 inline G4ThreeVector G4PrimaryParticle::GetMomentum() const 227 { 228 return GetTotalMomentum() * direction; 229 } 230 231 inline const G4ThreeVector& G4PrimaryParticle::GetMomentumDirection() const 232 { 233 return direction; 234 } 235 236 inline void G4PrimaryParticle::SetMomentumDirection(const G4ThreeVector& p) 237 { 238 direction = p; 239 } 240 241 inline G4double G4PrimaryParticle::GetPx() const 242 { 243 return GetTotalMomentum() * direction.x(); 244 } 245 246 inline G4double G4PrimaryParticle::GetPy() const 247 { 248 return GetTotalMomentum() * direction.y(); 249 } 250 251 inline G4double G4PrimaryParticle::GetPz() const 252 { 253 return GetTotalMomentum() * direction.z(); 254 } 255 256 inline G4double G4PrimaryParticle::GetTotalEnergy() const 257 { 258 if (mass < 0.) return kinE; 259 return kinE + mass; 260 } 261 262 inline void G4PrimaryParticle::SetTotalEnergy(G4double eTot) 263 { 264 if (mass < 0.) 265 kinE = eTot; 266 else 267 kinE = eTot - mass; 268 } 269 270 inline G4double G4PrimaryParticle::GetKineticEnergy() const 271 { 272 return kinE; 273 } 274 275 inline void G4PrimaryParticle::SetKineticEnergy(G4double eKin) 276 { 277 kinE = eKin; 278 } 279 280 inline G4PrimaryParticle* G4PrimaryParticle::GetNext() const 281 { 282 return nextParticle; 283 } 284 285 inline G4PrimaryParticle* G4PrimaryParticle::GetDaughter() const 286 { 287 return daughterParticle; 288 } 289 290 inline G4int G4PrimaryParticle::GetTrackID() const 291 { 292 return trackID; 293 } 294 295 inline G4ThreeVector G4PrimaryParticle::GetPolarization() const 296 { 297 return G4ThreeVector(polX, polY, polZ); 298 } 299 300 inline G4double G4PrimaryParticle::GetPolX() const 301 { 302 return polX; 303 } 304 305 inline G4double G4PrimaryParticle::GetPolY() const 306 { 307 return polY; 308 } 309 310 inline G4double G4PrimaryParticle::GetPolZ() const 311 { 312 return polZ; 313 } 314 315 inline G4double G4PrimaryParticle::GetWeight() const 316 { 317 return Weight0; 318 } 319 320 inline void G4PrimaryParticle::SetWeight(G4double w) 321 { 322 Weight0 = w; 323 } 324 325 inline void G4PrimaryParticle::SetProperTime(G4double t) 326 { 327 properTime = t; 328 } 329 330 inline G4double G4PrimaryParticle::GetProperTime() const 331 { 332 return properTime; 333 } 334 335 inline void G4PrimaryParticle::SetUserInformation(G4VUserPrimaryParticleInformation* anInfo) 336 { 337 userInfo = anInfo; 338 } 339 340 inline G4VUserPrimaryParticleInformation* G4PrimaryParticle::GetUserInformation() const 341 { 342 return userInfo; 343 } 344 345 inline void G4PrimaryParticle::SetG4code(const G4ParticleDefinition* Gcode) 346 { 347 SetParticleDefinition(Gcode); 348 } 349 350 inline void G4PrimaryParticle::SetNext(G4PrimaryParticle* np) 351 { 352 if (nextParticle == nullptr) { 353 nextParticle = np; 354 } 355 else { 356 nextParticle->SetNext(np); 357 } 358 } 359 360 inline void G4PrimaryParticle::ClearNext() 361 { 362 nextParticle = nullptr; 363 } 364 365 inline void G4PrimaryParticle::SetDaughter(G4PrimaryParticle* np) 366 { 367 if (daughterParticle == nullptr) { 368 daughterParticle = np; 369 } 370 else { 371 daughterParticle->SetNext(np); 372 } 373 } 374 375 inline void G4PrimaryParticle::SetTrackID(G4int id) 376 { 377 trackID = id; 378 } 379 380 inline void G4PrimaryParticle::SetMass(G4double mas) 381 { 382 mass = mas; 383 } 384 385 inline void G4PrimaryParticle::SetCharge(G4double chg) 386 { 387 charge = chg; 388 } 389 390 inline void G4PrimaryParticle::SetPolarization(G4double px, G4double py, G4double pz) 391 { 392 polX = px; 393 polY = py; 394 polZ = pz; 395 } 396 397 inline void G4PrimaryParticle::SetPolarization(const G4ThreeVector& pol) 398 { 399 polX = pol.x(); 400 polY = pol.y(); 401 polZ = pol.z(); 402 } 403 404 #endif 405