Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FermiNucleiProperties.cc
Go to the documentation of this file.
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// G4FermiBreakUpAN alternative FermiBreakUp model
28// by A. Novikov (January 2025)
29//
30
32
34#include <G4NucleiProperties.hh>
36
37namespace
38{
39std::size_t GetSlot(G4FermiAtomicMass atomicMass, G4FermiChargeNumber chargeNumber)
40{
41 const auto mass = static_cast<std::uint32_t>(atomicMass);
42 const auto charge = static_cast<std::uint32_t>(chargeNumber);
43 return (mass * (mass + 1)) / 2 + charge;
44}
45} // namespace
46
47G4FermiNucleiProperties::G4FermiNucleiProperties()
48{
49 for (auto a = 1; a < MAX_A; ++a) {
50 for (auto z = 0; z <= a; ++z) {
51 const auto atomicMass = G4FermiAtomicMass(a);
52 const auto chargeNumber = G4FermiChargeNumber(z);
53
54 const auto mass = G4NucleiProperties::GetNuclearMass(a, z);
55 if (mass > 0.) {
56 InsertNuclei(atomicMass, chargeNumber, mass, G4NucleiProperties::IsInStableTable(a, z));
57 }
58 }
59 }
60}
61
62G4double G4FermiNucleiProperties::GetNuclearMassImpl(G4FermiAtomicMass atomicMass,
63 G4FermiChargeNumber chargeNumber) const
64{
65 FERMI_ASSERT_MSG(static_cast<std::uint32_t>(atomicMass)
66 >= static_cast<std::uint32_t>(chargeNumber),
67 "invalid nuclei A = " << atomicMass << ", Z = " << chargeNumber);
68
69 const auto slot = GetSlot(atomicMass, chargeNumber);
70 if (slot < nucleiMasses_.size() && nucleiMasses_[slot].isCached) {
71 return nucleiMasses_[slot].mass;
72 }
73
74 return G4NucleiProperties::GetNuclearMass(G4int(atomicMass), G4int(chargeNumber));
75}
76
77G4bool G4FermiNucleiProperties::IsStableImpl(G4FermiAtomicMass atomicMass,
78 G4FermiChargeNumber chargeNumber) const
79{
80 if (atomicMass < 1_m || chargeNumber < 0_c
81 || static_cast<std::uint32_t>(chargeNumber)
82 > static_cast<std::uint32_t>(atomicMass))
83 {
84 return false;
85 }
86
87 const auto slot = GetSlot(atomicMass, chargeNumber);
88
89 return slot < nucleiMasses_.size() && nucleiMasses_[slot].isStable;
90}
91
93 G4FermiChargeNumber chargeNumber, G4double mass,
94 G4bool isStable)
95{
96 const auto slot = GetSlot(atomicMass, chargeNumber);
97 if (slot >= nucleiMasses_.size()) {
98 nucleiMasses_.resize(slot + static_cast<std::uint32_t>(atomicMass));
99 }
100
101 nucleiMasses_[slot] = G4FermiMassData{
102 mass, // mass
103 isStable, // isStable
104 true, // isCached
105 };
106}
#define FERMI_ASSERT_MSG(COND, MSG)
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void InsertNuclei(G4FermiAtomicMass atomicMass, G4FermiChargeNumber chargeNumber, G4double mass, G4bool isStable=true)
static G4bool IsInStableTable(const G4double A, const G4double Z)
static G4double GetNuclearMass(const G4double A, const G4double Z)