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 // 27 //--------------------------------------------------------------------- 28 // 29 // Geant4 header G4Fragment 30 // 31 // by V. Lara (May 1998) 32 // 33 // Modifications: 34 // 03.05.2010 V.Ivanchenko General cleanup of inline functions: objects 35 // are accessed by reference; remove double return 36 // tolerance of excitation energy at modent it is computed; 37 // safe computation of excitation for exotic fragments 38 // 18.05.2010 V.Ivanchenko added member theGroundStateMass and inline 39 // method which allowing to compute this value once and use 40 // many times 41 // 26.09.2010 V.Ivanchenko added number of protons, neutrons, proton holes 42 // and neutron holes as members of the class and Get/Set methods; 43 // removed not needed 'const'; removed old debug staff and unused 44 // private methods; add comments and reorder methods for 45 // better reading 46 // 27.10.2021 A.Ribon extension for hypernuclei. 47 48 #ifndef G4Fragment_h 49 #define G4Fragment_h 1 50 51 #include "globals.hh" 52 #include "G4Allocator.hh" 53 #include "G4LorentzVector.hh" 54 #include "G4ThreeVector.hh" 55 #include "G4NuclearPolarization.hh" 56 #include "G4NucleiProperties.hh" 57 #include "G4HyperNucleiProperties.hh" 58 #include "G4Proton.hh" 59 #include "G4Neutron.hh" 60 #include <vector> 61 62 class G4ParticleDefinition; 63 64 class G4Fragment; 65 typedef std::vector<G4Fragment*> G4FragmentVector; 66 67 class G4Fragment 68 { 69 public: 70 71 // ============= CONSTRUCTORS ================== 72 73 // Default constructor - obsolete 74 G4Fragment(); 75 76 // Destructor 77 ~G4Fragment(); 78 79 // Copy constructor 80 G4Fragment(const G4Fragment &right); 81 82 // A,Z and 4-momentum - main constructor for fragment 83 G4Fragment(G4int A, G4int Z, const G4LorentzVector& aMomentum); 84 85 // A,Z,numberOfLambdas and 4-momentum 86 G4Fragment(G4int A, G4int Z, G4int numberOfLambdas, 87 const G4LorentzVector& aMomentum); 88 89 // 4-momentum and pointer to G4particleDefinition (for gammas, e-) 90 G4Fragment(const G4LorentzVector& aMomentum, 91 const G4ParticleDefinition* aParticleDefinition); 92 93 // ============= OPERATORS ================== 94 95 G4Fragment & operator=(const G4Fragment &right); 96 G4bool operator==(const G4Fragment &right) const; 97 G4bool operator!=(const G4Fragment &right) const; 98 99 friend std::ostream& operator<<(std::ostream&, const G4Fragment&); 100 101 // new/delete operators are overloded to use G4Allocator 102 inline void *operator new(size_t); 103 inline void operator delete(void *aFragment); 104 105 // ============= GENERAL METHODS ================== 106 107 inline G4int GetZ_asInt() const; 108 inline G4int GetA_asInt() const; 109 110 // update number of nucleons without check on input 111 // ground state mass is not recomputed 112 inline void SetZandA_asInt(G4int Znew, G4int Anew, G4int Lnew=0); 113 114 // non-negative number of lambdas/anti-lambdas 115 // ground state mass is not recomputed 116 void SetNumberOfLambdas(G4int numberOfLambdas); 117 118 inline G4int GetNumberOfLambdas() const; 119 120 inline G4double GetExcitationEnergy() const; 121 122 inline G4double GetGroundStateMass() const; 123 124 inline const G4LorentzVector& GetMomentum() const; 125 126 // ground state mass and excitation energy are recomputed 127 inline G4double RecomputeGroundStateMass(); 128 129 // update main fragment parameters full check on input 130 // ground state mass and excitation energy are recomputed 131 inline void SetMomentum(const G4LorentzVector& value); 132 inline void SetZAandMomentum(const G4LorentzVector&, 133 G4int Z, G4int A, 134 G4int nLambdas = 0); 135 136 // ground state mass is not recomputed 137 void SetExcEnergyAndMomentum(G4double eexc, const G4LorentzVector&); 138 G4double GetBindingEnergy() const; 139 140 // computation of mass for any imput Z, A and number of Lambdas 141 // no check on input values 142 inline G4double ComputeGroundStateMass(G4int Z, G4int A, 143 G4int nLambdas = 0) const; 144 145 // extra methods 146 inline G4double GetSpin() const; 147 inline void SetSpin(G4double value); 148 149 inline G4int GetCreatorModelID() const; 150 inline void SetCreatorModelID(G4int value); 151 152 inline G4bool IsLongLived() const; 153 inline void SetLongLived(G4bool value); 154 155 // obsolete methods 156 inline G4double GetZ() const; 157 inline G4double GetA() const; 158 inline void SetZ(G4double value); 159 inline void SetA(G4double value); 160 161 // ============= METHODS FOR PRE-COMPOUND MODEL =============== 162 163 inline G4int GetNumberOfExcitons() const; 164 165 inline G4int GetNumberOfParticles() const; 166 inline G4int GetNumberOfCharged() const; 167 inline void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP); 168 169 inline G4int GetNumberOfHoles() const; 170 inline G4int GetNumberOfChargedHoles() const; 171 inline void SetNumberOfHoles(G4int valueTot, G4int valueP=0); 172 173 // these methods will be removed in future 174 inline void SetNumberOfParticles(G4int value); 175 inline void SetNumberOfCharged(G4int value); 176 177 // ============= METHODS FOR PHOTON EVAPORATION =============== 178 179 inline G4int GetNumberOfElectrons() const; 180 inline void SetNumberOfElectrons(G4int value); 181 182 inline G4int GetFloatingLevelNumber() const; 183 inline void SetFloatingLevelNumber(G4int value); 184 185 inline const G4ParticleDefinition * GetParticleDefinition() const; 186 inline void SetParticleDefinition(const G4ParticleDefinition * p); 187 188 inline G4double GetCreationTime() const; 189 inline void SetCreationTime(G4double time); 190 191 // G4Fragment class is not responsible for creation and delition of 192 // G4NuclearPolarization object 193 inline G4NuclearPolarization* NuclearPolarization(); 194 inline G4NuclearPolarization* GetNuclearPolarization() const; 195 inline void SetNuclearPolarization(G4NuclearPolarization*); 196 197 void SetAngularMomentum(const G4ThreeVector&); 198 G4ThreeVector GetAngularMomentum() const; 199 200 // ============= PRIVATE METHODS ============================== 201 202 private: 203 204 void CalculateMassAndExcitationEnergy(); 205 206 void ExcitationEnergyWarning(); 207 208 void NumberOfExitationWarning(const G4String&); 209 210 // ============= DATA MEMBERS ================== 211 212 G4int theA; 213 214 G4int theZ; 215 216 // Non-negative number of lambdas/anti-lambdas inside the nucleus/anti-nucleus 217 G4int theL; 218 219 G4double theExcitationEnergy; 220 221 G4double theGroundStateMass; 222 223 G4LorentzVector theMomentum; 224 225 // Nuclear polarisation by default is nullptr 226 G4NuclearPolarization* thePolarization; 227 228 // creator model type 229 G4int creatorModel; 230 231 // Exciton model data members 232 G4int numberOfParticles; 233 G4int numberOfCharged; 234 G4int numberOfHoles; 235 G4int numberOfChargedHoles; 236 237 // Gamma evaporation data members 238 G4int numberOfShellElectrons; 239 G4int xLevel; 240 241 const G4ParticleDefinition* theParticleDefinition; 242 243 G4double spin; 244 G4double theCreationTime; 245 246 G4bool isLongLived = false; 247 }; 248 249 // ============= INLINE METHOD IMPLEMENTATIONS =================== 250 251 #if defined G4HADRONIC_ALLOC_EXPORT 252 extern G4DLLEXPORT G4Allocator<G4Fragment>*& pFragmentAllocator(); 253 #else 254 extern G4DLLIMPORT G4Allocator<G4Fragment>*& pFragmentAllocator(); 255 #endif 256 257 inline void * G4Fragment::operator new(size_t) 258 { 259 if (!pFragmentAllocator()) { 260 pFragmentAllocator() = new G4Allocator<G4Fragment>; 261 } 262 return (void*) pFragmentAllocator()->MallocSingle(); 263 } 264 265 inline void G4Fragment::operator delete(void * aFragment) 266 { 267 pFragmentAllocator()->FreeSingle((G4Fragment *) aFragment); 268 } 269 270 inline G4double 271 G4Fragment::ComputeGroundStateMass(G4int Z, G4int A, G4int nLambdas) const 272 { 273 return ( nLambdas <= 0 ) 274 ? G4NucleiProperties::GetNuclearMass(A, Z) 275 : G4HyperNucleiProperties::GetNuclearMass(A, Z, nLambdas); 276 } 277 278 inline G4double G4Fragment::RecomputeGroundStateMass() 279 { 280 CalculateMassAndExcitationEnergy(); 281 return theGroundStateMass; 282 } 283 284 inline G4int G4Fragment::GetA_asInt() const 285 { 286 return theA; 287 } 288 289 inline G4int G4Fragment::GetZ_asInt() const 290 { 291 return theZ; 292 } 293 294 inline void 295 G4Fragment::SetZandA_asInt(G4int Znew, G4int Anew, G4int Lnew) 296 { 297 theZ = Znew; 298 theA = Anew; 299 theL = std::max(Lnew, 0); 300 } 301 302 inline void G4Fragment::SetNumberOfLambdas(G4int Lnew) 303 { 304 theL = std::max(Lnew, 0); 305 } 306 307 inline G4int G4Fragment::GetNumberOfLambdas() const 308 { 309 return theL; 310 } 311 312 inline G4double G4Fragment::GetExcitationEnergy() const 313 { 314 return theExcitationEnergy; 315 } 316 317 inline G4double G4Fragment::GetGroundStateMass() const 318 { 319 return theGroundStateMass; 320 } 321 322 inline const G4LorentzVector& G4Fragment::GetMomentum() const 323 { 324 return theMomentum; 325 } 326 327 inline void G4Fragment::SetMomentum(const G4LorentzVector& value) 328 { 329 theMomentum = value; 330 CalculateMassAndExcitationEnergy(); 331 } 332 333 inline void 334 G4Fragment::SetZAandMomentum(const G4LorentzVector& v, 335 G4int Z, G4int A, G4int nLambdas) 336 { 337 SetZandA_asInt(Z, A, nLambdas); 338 SetMomentum(v); 339 } 340 341 inline G4double G4Fragment::GetZ() const 342 { 343 return static_cast<G4double>(theZ); 344 } 345 346 inline G4double G4Fragment::GetA() const 347 { 348 return static_cast<G4double>(theA); 349 } 350 351 inline void G4Fragment::SetZ(const G4double value) 352 { 353 theZ = G4lrint(value); 354 } 355 356 inline void G4Fragment::SetA(const G4double value) 357 { 358 theA = G4lrint(value); 359 } 360 361 inline G4int G4Fragment::GetNumberOfExcitons() const 362 { 363 return numberOfParticles + numberOfHoles; 364 } 365 366 inline G4int G4Fragment::GetNumberOfParticles() const 367 { 368 return numberOfParticles; 369 } 370 371 inline G4int G4Fragment::GetNumberOfCharged() const 372 { 373 return numberOfCharged; 374 } 375 376 inline 377 void G4Fragment::SetNumberOfExcitedParticle(G4int valueTot, G4int valueP) 378 { 379 numberOfParticles = valueTot; 380 numberOfCharged = valueP; 381 if(valueTot < valueP) { 382 NumberOfExitationWarning("SetNumberOfExcitedParticle"); 383 } 384 } 385 386 inline G4int G4Fragment::GetNumberOfHoles() const 387 { 388 return numberOfHoles; 389 } 390 391 inline G4int G4Fragment::GetNumberOfChargedHoles() const 392 { 393 return numberOfChargedHoles; 394 } 395 396 inline void G4Fragment::SetNumberOfHoles(G4int valueTot, G4int valueP) 397 { 398 numberOfHoles = valueTot; 399 numberOfChargedHoles = valueP; 400 if(valueTot < valueP) { 401 NumberOfExitationWarning("SetNumberOfHoles"); 402 } 403 } 404 405 inline void G4Fragment::SetNumberOfParticles(G4int value) 406 { 407 numberOfParticles = value; 408 } 409 410 inline void G4Fragment::SetNumberOfCharged(G4int value) 411 { 412 numberOfCharged = value; 413 if(value > numberOfParticles) { 414 NumberOfExitationWarning("SetNumberOfCharged"); 415 } 416 } 417 418 inline G4int G4Fragment::GetNumberOfElectrons() const 419 { 420 return numberOfShellElectrons; 421 } 422 423 inline void G4Fragment::SetNumberOfElectrons(G4int value) 424 { 425 numberOfShellElectrons = value; 426 } 427 428 inline G4int G4Fragment::GetCreatorModelID() const 429 { 430 return creatorModel; 431 } 432 433 inline void G4Fragment::SetCreatorModelID(G4int value) 434 { 435 creatorModel = value; 436 } 437 438 inline G4double G4Fragment::GetSpin() const 439 { 440 return spin; 441 } 442 443 inline void G4Fragment::SetSpin(G4double value) 444 { 445 spin = value; 446 } 447 448 inline G4bool G4Fragment::IsLongLived() const 449 { 450 return isLongLived; 451 } 452 453 inline void G4Fragment::SetLongLived(G4bool value) 454 { 455 isLongLived = value; 456 } 457 458 inline G4int G4Fragment::GetFloatingLevelNumber() const 459 { 460 return xLevel; 461 } 462 463 inline void G4Fragment::SetFloatingLevelNumber(G4int value) 464 { 465 xLevel = value; 466 } 467 468 inline 469 const G4ParticleDefinition* G4Fragment::GetParticleDefinition(void) const 470 { 471 return theParticleDefinition; 472 } 473 474 inline void G4Fragment::SetParticleDefinition(const G4ParticleDefinition * p) 475 { 476 theParticleDefinition = p; 477 } 478 479 inline G4double G4Fragment::GetCreationTime() const 480 { 481 return theCreationTime; 482 } 483 484 inline void G4Fragment::SetCreationTime(G4double time) 485 { 486 theCreationTime = time; 487 } 488 489 inline G4NuclearPolarization* G4Fragment::NuclearPolarization() 490 { 491 return thePolarization; 492 } 493 494 inline G4NuclearPolarization* G4Fragment::GetNuclearPolarization() const 495 { 496 return thePolarization; 497 } 498 499 inline void G4Fragment::SetNuclearPolarization(G4NuclearPolarization* p) 500 { 501 thePolarization = p; 502 } 503 504 #endif 505