Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | Static Private Attributes
G4INCL::NuclearPotential::NuclearPotentialEnergyIsospinSmooth Class Reference

#include <G4INCLNuclearPotentialEnergyIsospinSmooth.hh>

Inheritance diagram for G4INCL::NuclearPotential::NuclearPotentialEnergyIsospinSmooth:
G4INCL::NuclearPotential::NuclearPotentialIsospin G4INCL::NuclearPotential::INuclearPotential

Public Member Functions

virtual G4double computePotentialEnergy (const Particle *const p) const
 
G4double getFermiEnergy (const Particle *const p) const
 Return the Fermi energy for a particle. More...
 
G4double getFermiEnergy (const ParticleType t) const
 Return the Fermi energy for a particle type. More...
 
G4double getFermiMomentum (const Particle *const p) const
 Return the Fermi momentum for a particle. More...
 
G4double getFermiMomentum (const ParticleType t) const
 Return the Fermi momentum for a particle type. More...
 
G4double getSeparationEnergy (const Particle *const p) const
 Return the separation energy for a particle. More...
 
G4double getSeparationEnergy (const ParticleType t) const
 Return the separation energy for a particle type. More...
 
G4bool hasPionPotential () const
 Do we have a pion potential? More...
 
 NuclearPotentialEnergyIsospinSmooth (const G4int A, const G4int Z, const G4bool pionPotential)
 
virtual ~NuclearPotentialEnergyIsospinSmooth ()
 

Protected Member Functions

G4double computeKaonPotentialEnergy (const Particle *const p) const
 Compute the potential energy for the given kaon. More...
 
G4double computePionPotentialEnergy (const Particle *const p) const
 Compute the potential energy for the given pion. More...
 
G4double computePionResonancePotentialEnergy (const Particle *const p) const
 Compute the potential energy for the given pion resonances (Eta, Omega and EtaPrime and Gamma also). More...
 

Protected Attributes

std::map< ParticleType, G4doublefermiEnergy
 
std::map< ParticleType, G4doublefermiMomentum
 
std::map< ParticleType, G4doubleseparationEnergy
 
const G4int theA
 The mass number of the nucleus. More...
 
const G4int theZ
 The charge number of the nucleus. More...
 

Private Member Functions

void initialize ()
 

Private Attributes

const G4bool pionPotential
 
G4double vDeltaMinus
 
G4double vDeltaPlus
 
G4double vDeltaPlusPlus
 
G4double vDeltaZero
 
G4double vKMinus
 
G4double vKPlus
 
G4double vKZero
 
G4double vKZeroBar
 
G4double vLambda
 
G4double vNeutron
 
G4double vPiMinus
 
G4double vPiPlus
 
G4double vPiZero
 
G4double vProton
 
G4double vSigmaMinus
 
G4double vSigmaPlus
 
G4double vSigmaZero
 

Static Private Attributes

static const G4double alpha = 0.223
 Slope of the V(T) curve. More...
 
static const G4double deltaE = 25.
 Distance from the cusp where the exponential kicks in. More...
 
static const G4double vKMinusDefault = 60.
 
static const G4double vKPlusDefault = -25.
 
static const G4double vPionDefault = 30.6
 

Detailed Description

Definition at line 57 of file G4INCLNuclearPotentialEnergyIsospinSmooth.hh.

Constructor & Destructor Documentation

◆ NuclearPotentialEnergyIsospinSmooth()

G4INCL::NuclearPotential::NuclearPotentialEnergyIsospinSmooth::NuclearPotentialEnergyIsospinSmooth ( const G4int  A,
const G4int  Z,
const G4bool  pionPotential 
)

Definition at line 60 of file G4INCLNuclearPotentialEnergyIsospinSmooth.cc.

61 : NuclearPotentialIsospin(A,Z,aPionPotential)
62 {}
const G4int Z[17]
const G4double A[17]
NuclearPotentialIsospin(const G4int A, const G4int Z, const G4bool pionPotential)

◆ ~NuclearPotentialEnergyIsospinSmooth()

G4INCL::NuclearPotential::NuclearPotentialEnergyIsospinSmooth::~NuclearPotentialEnergyIsospinSmooth ( )
virtual

Definition at line 65 of file G4INCLNuclearPotentialEnergyIsospinSmooth.cc.

65{}

Member Function Documentation

◆ computeKaonPotentialEnergy()

G4double G4INCL::NuclearPotential::INuclearPotential::computeKaonPotentialEnergy ( const Particle *const  p) const
inlineprotectedinherited

Compute the potential energy for the given kaon.

Definition at line 197 of file G4INCLINuclearPotential.hh.

197 {
198// assert(p->getType()==KPlus || p->getType()==KZero || p->getType()==KZeroBar || p->getType()==KMinus|| p->getType()==KShort|| p->getType()==KLong);
199 if(pionPotential && !p->isOutOfWell()) { // if pionPotental false -> kaonPotential false
200 switch( p->getType() ) {
201 case KPlus:
202 return vKPlus;
203 break;
204 case KZero:
205 return vKZero;
206 break;
207 case KZeroBar:
208 return vKZeroBar;
209 break;
210 case KShort:
211 case KLong:
212 return 0.0; // Should never be in the nucleus
213 break;
214 case KMinus:
215 return vKMinus;
216 break;
217 default:
218 return 0.0;
219 break;
220 }
221 }
222 else
223 return 0.0;
224 }

References G4INCL::Particle::getType(), G4INCL::Particle::isOutOfWell(), G4INCL::KLong, G4INCL::KMinus, G4INCL::KPlus, G4INCL::KShort, G4INCL::KZero, G4INCL::KZeroBar, G4INCL::NuclearPotential::INuclearPotential::pionPotential, G4INCL::NuclearPotential::INuclearPotential::vKMinus, G4INCL::NuclearPotential::INuclearPotential::vKPlus, G4INCL::NuclearPotential::INuclearPotential::vKZero, and G4INCL::NuclearPotential::INuclearPotential::vKZeroBar.

Referenced by G4INCL::NuclearPotential::NuclearPotentialConstant::computePotentialEnergy(), and G4INCL::NuclearPotential::NuclearPotentialIsospin::computePotentialEnergy().

◆ computePionPotentialEnergy()

G4double G4INCL::NuclearPotential::INuclearPotential::computePionPotentialEnergy ( const Particle *const  p) const
inlineprotectedinherited

Compute the potential energy for the given pion.

Definition at line 173 of file G4INCLINuclearPotential.hh.

173 {
174// assert(p->getType()==PiPlus || p->getType()==PiZero || p->getType()==PiMinus);
175 if(pionPotential && !p->isOutOfWell()) {
176 switch( p->getType() ) {
177 case PiPlus:
178 return vPiPlus;
179 break;
180 case PiZero:
181 return vPiZero;
182 break;
183 case PiMinus:
184 return vPiMinus;
185 break;
186 default: // Pion potential is defined and non-zero only for pions
187 return 0.0;
188 break;
189 }
190 }
191 else
192 return 0.0;
193 }

References G4INCL::Particle::getType(), G4INCL::Particle::isOutOfWell(), G4INCL::PiMinus, G4INCL::NuclearPotential::INuclearPotential::pionPotential, G4INCL::PiPlus, G4INCL::PiZero, G4INCL::NuclearPotential::INuclearPotential::vPiMinus, G4INCL::NuclearPotential::INuclearPotential::vPiPlus, and G4INCL::NuclearPotential::INuclearPotential::vPiZero.

Referenced by G4INCL::NuclearPotential::NuclearPotentialConstant::computePotentialEnergy(), and G4INCL::NuclearPotential::NuclearPotentialIsospin::computePotentialEnergy().

◆ computePionResonancePotentialEnergy()

G4double G4INCL::NuclearPotential::INuclearPotential::computePionResonancePotentialEnergy ( const Particle *const  p) const
inlineprotectedinherited

Compute the potential energy for the given pion resonances (Eta, Omega and EtaPrime and Gamma also).

Definition at line 228 of file G4INCLINuclearPotential.hh.

228 {
229// assert(p->getType()==Eta || p->getType()==Omega || p->getType()==EtaPrime || p->getType()==Photon);
230 if(pionPotential && !p->isOutOfWell()) {
231 switch( p->getType() ) {
232 case Eta:
233//jcd return vPiZero;
234//jcd return vPiZero*1.5;
235 return 0.0; // (JCD: seems to give better results)
236 break;
237 case Omega:
238 return 15.0; // S.Friedrich et al., Physics Letters B736(2014)26-32. (V. Metag in Hyperfine Interact (2015) 234:25-31 gives 29 MeV)
239 break;
240 case EtaPrime:
241 return 37.0; // V. Metag in Hyperfine Interact (2015) 234:25-31
242 break;
243 case Photon:
244 return 0.0;
245 break;
246 default:
247 return 0.0;
248 break;
249 }
250 }
251 else
252 return 0.0;
253 }

References G4INCL::Eta, G4INCL::EtaPrime, G4INCL::Particle::getType(), G4INCL::Particle::isOutOfWell(), G4INCL::Omega, G4INCL::Photon, and G4INCL::NuclearPotential::INuclearPotential::pionPotential.

Referenced by G4INCL::NuclearPotential::NuclearPotentialConstant::computePotentialEnergy(), and G4INCL::NuclearPotential::NuclearPotentialIsospin::computePotentialEnergy().

◆ computePotentialEnergy()

G4double G4INCL::NuclearPotential::NuclearPotentialEnergyIsospinSmooth::computePotentialEnergy ( const Particle *const  p) const
virtual

Reimplemented from G4INCL::NuclearPotential::NuclearPotentialIsospin.

Definition at line 67 of file G4INCLNuclearPotentialEnergyIsospinSmooth.cc.

67 {
68
70
71 if(particle->isNucleon()) {
72 const G4double t = particle->getKineticEnergy();
73 const G4double tf = getFermiEnergy(particle);
74 // Constant potential for T<Tf
75 if(t < tf)
76 return v0;
77
78 // Linear function for Tf<T<T0, exponential function for T>T0
79 const G4double t0 = tf + v0*(1.-alpha)/alpha - deltaE; // deltaE before the linear potential vanishes
80 G4double v;
81 if(t<t0) {
82 v = v0 - alpha*(t-tf)/(1.-alpha);
83 } else {
84 const G4double v_at_t0 = v0 - alpha*(t0-tf)/(1.-alpha);
85 const G4double kappa = alpha / (v_at_t0 * (1.-alpha));
86 v = v_at_t0 * std::exp(kappa * (t0-t));
87 }
88 return (v>0.0) ? v : 0.0;
89 } else
90 return v0;
91 }
double G4double
Definition: G4Types.hh:83
G4double getFermiEnergy(const Particle *const p) const
Return the Fermi energy for a particle.
static const G4double deltaE
Distance from the cusp where the exponential kicks in.
virtual G4double computePotentialEnergy(const Particle *const p) const

References alpha, G4INCL::NuclearPotential::NuclearPotentialIsospin::computePotentialEnergy(), deltaE, G4INCL::NuclearPotential::INuclearPotential::getFermiEnergy(), G4INCL::Particle::getKineticEnergy(), and G4INCL::Particle::isNucleon().

◆ getFermiEnergy() [1/2]

G4double G4INCL::NuclearPotential::INuclearPotential::getFermiEnergy ( const Particle *const  p) const
inlineinherited

Return the Fermi energy for a particle.

Parameters
ppointer to a Particle
Returns
Fermi energy for that particle type

Definition at line 105 of file G4INCLINuclearPotential.hh.

105 {
106 std::map<ParticleType, G4double>::const_iterator i = fermiEnergy.find(p->getType());
107// assert(i!=fermiEnergy.end());
108 return i->second;
109 }
std::map< ParticleType, G4double > fermiEnergy

References G4INCL::NuclearPotential::INuclearPotential::fermiEnergy, and G4INCL::Particle::getType().

Referenced by G4INCL::NuclearPotential::NuclearPotentialEnergyIsospin::computePotentialEnergy(), computePotentialEnergy(), G4INCL::ParticleEntryChannel::fillFinalState(), G4INCL::SurfaceAvatar::getChannel(), G4INCL::NuclearPotential::INuclearPotential::getFermiMomentum(), and G4INCL::CDPP::processOneParticle().

◆ getFermiEnergy() [2/2]

G4double G4INCL::NuclearPotential::INuclearPotential::getFermiEnergy ( const ParticleType  t) const
inlineinherited

Return the Fermi energy for a particle type.

Parameters
tparticle type
Returns
Fermi energy for that particle type

Definition at line 116 of file G4INCLINuclearPotential.hh.

116 {
117 std::map<ParticleType, G4double>::const_iterator i = fermiEnergy.find(t);
118// assert(i!=fermiEnergy.end());
119 return i->second;
120 }

References G4INCL::NuclearPotential::INuclearPotential::fermiEnergy.

◆ getFermiMomentum() [1/2]

G4double G4INCL::NuclearPotential::INuclearPotential::getFermiMomentum ( const Particle *const  p) const
inlineinherited

Return the Fermi momentum for a particle.

Parameters
ppointer to a Particle
Returns
Fermi momentum for that particle type

Definition at line 149 of file G4INCLINuclearPotential.hh.

149 {
150 if(p->isDelta()) {
151 const G4double Tf = getFermiEnergy(p), mass = p->getMass();
152 return std::sqrt(Tf*(Tf+2.*mass));
153 } else {
154 std::map<ParticleType, G4double>::const_iterator i = fermiMomentum.find(p->getType());
155// assert(i!=fermiMomentum.end());
156 return i->second;
157 }
158 }
std::map< ParticleType, G4double > fermiMomentum

References G4INCL::NuclearPotential::INuclearPotential::fermiMomentum, G4INCL::NuclearPotential::INuclearPotential::getFermiEnergy(), G4INCL::Particle::getMass(), G4INCL::Particle::getType(), and G4INCL::Particle::isDelta().

Referenced by G4INCL::PauliStandard::getBlockingProbability(), G4INCL::Nucleus::getSurfaceRadius(), G4INCL::ParticleSampler::sampleOneParticleWithFuzzyRPCorrelation(), and G4INCL::ParticleSampler::sampleOneParticleWithRPCorrelation().

◆ getFermiMomentum() [2/2]

G4double G4INCL::NuclearPotential::INuclearPotential::getFermiMomentum ( const ParticleType  t) const
inlineinherited

Return the Fermi momentum for a particle type.

Parameters
tparticle type
Returns
Fermi momentum for that particle type

Definition at line 165 of file G4INCLINuclearPotential.hh.

165 {
166// assert(t!=DeltaPlusPlus && t!=DeltaPlus && t!=DeltaZero && t!=DeltaMinus);
167 std::map<ParticleType, G4double>::const_iterator i = fermiMomentum.find(t);
168 return i->second;
169 }

References G4INCL::NuclearPotential::INuclearPotential::fermiMomentum.

◆ getSeparationEnergy() [1/2]

G4double G4INCL::NuclearPotential::INuclearPotential::getSeparationEnergy ( const Particle *const  p) const
inlineinherited

Return the separation energy for a particle.

Parameters
ppointer to a Particle
Returns
separation energy for that particle type

Definition at line 127 of file G4INCLINuclearPotential.hh.

127 {
128 std::map<ParticleType, G4double>::const_iterator i = separationEnergy.find(p->getType());
129// assert(i!=separationEnergy.end());
130 return i->second;
131 }
std::map< ParticleType, G4double > separationEnergy

References G4INCL::Particle::getType(), and G4INCL::NuclearPotential::INuclearPotential::separationEnergy.

Referenced by G4INCL::Nucleus::computeSeparationEnergyBalance(), G4INCL::Nucleus::Nucleus(), and G4INCL::CDPP::processOneParticle().

◆ getSeparationEnergy() [2/2]

G4double G4INCL::NuclearPotential::INuclearPotential::getSeparationEnergy ( const ParticleType  t) const
inlineinherited

Return the separation energy for a particle type.

Parameters
tparticle type
Returns
separation energy for that particle type

Definition at line 138 of file G4INCLINuclearPotential.hh.

138 {
139 std::map<ParticleType, G4double>::const_iterator i = separationEnergy.find(t);
140// assert(i!=separationEnergy.end());
141 return i->second;
142 }

References G4INCL::NuclearPotential::INuclearPotential::separationEnergy.

◆ hasPionPotential()

G4bool G4INCL::NuclearPotential::INuclearPotential::hasPionPotential ( ) const
inlineinherited

Do we have a pion potential?

Definition at line 96 of file G4INCLINuclearPotential.hh.

96{ return pionPotential; }

References G4INCL::NuclearPotential::INuclearPotential::pionPotential.

Referenced by G4INCL::Nucleus::decayInsideDeltas().

◆ initialize()

void G4INCL::NuclearPotential::NuclearPotentialIsospin::initialize ( )
privateinherited

Definition at line 66 of file G4INCLNuclearPotentialIsospin.cc.

66 {
67 const G4double ZOverA = ((G4double) theZ) / ((G4double) theA);
68
72
73 const G4double theFermiMomentum = ParticleTable::getFermiMomentum(theA,theZ);
74
75 fermiMomentum[Proton] = theFermiMomentum * Math::pow13(2.*ZOverA);
76 const G4double theProtonFermiEnergy = std::sqrt(fermiMomentum[Proton]*fermiMomentum[Proton] + mp*mp) - mp;
77 fermiEnergy[Proton] = theProtonFermiEnergy;
78 // Use separation energies from the ParticleTable
79 const G4double theProtonSeparationEnergy = ParticleTable::getSeparationEnergy(Proton,theA,theZ);
80 separationEnergy[Proton] = theProtonSeparationEnergy;
81 vProton = theProtonFermiEnergy + theProtonSeparationEnergy;
82
83 fermiMomentum[Neutron] = theFermiMomentum * Math::pow13(2.*(1.-ZOverA));
84 const G4double theNeutronFermiEnergy = std::sqrt(fermiMomentum[Neutron]*fermiMomentum[Neutron] + mn*mn) - mn;
85 fermiEnergy[Neutron] = theNeutronFermiEnergy;
86 // Use separation energies from the ParticleTable
87 const G4double theNeutronSeparationEnergy = ParticleTable::getSeparationEnergy(Neutron,theA,theZ);
88 separationEnergy[Neutron] = theNeutronSeparationEnergy;
89 vNeutron = theNeutronFermiEnergy + theNeutronSeparationEnergy;
90
91 const G4double separationEnergyDeltaPlusPlus = 2.*theProtonSeparationEnergy - theNeutronSeparationEnergy;
92 separationEnergy[DeltaPlusPlus] = separationEnergyDeltaPlusPlus;
93 separationEnergy[DeltaPlus] = theProtonSeparationEnergy;
94 separationEnergy[DeltaZero] = theNeutronSeparationEnergy;
95 const G4double separationEnergyDeltaMinus = 2.*theNeutronSeparationEnergy - theProtonSeparationEnergy;
96 separationEnergy[DeltaMinus] = separationEnergyDeltaMinus;
97
98 const G4double tinyMargin = 1E-7;
101 vDeltaPlusPlus = std::max(separationEnergyDeltaPlusPlus + tinyMargin, 2.*vDeltaPlus - vDeltaZero);
102 vDeltaMinus = std::max(separationEnergyDeltaMinus + tinyMargin, 2.*vDeltaZero - vDeltaPlus);
103
104 vSigmaMinus = -16.; // Repulsive potential, from Eur. Phys.J.A. (2016) 52:21
105 vSigmaZero = -16.; // hypothesis: same potential for each sigma
106 vSigmaPlus = -16.;
107
108 vLambda = 30.;
109 const G4double asy = (theA - 2.*theZ)/theA;
110 // Jose Luis Rodriguez-Sanchez et al., Rapid Communication PRC 98, 021602 (2018)
111 if (asy > 0.236) vLambda = 40.91;
112 else if (asy > 0.133) vLambda = 56.549 - 678.73*asy + 4905.35*asy*asy - 9789.1*asy*asy*asy;
113
114 const G4double theLambdaSeparationEnergy = ParticleTable::getSeparationEnergy(Lambda,theA,theZ);
115
116 separationEnergy[PiPlus] = theProtonSeparationEnergy - theNeutronSeparationEnergy;
118 separationEnergy[PiMinus] = theNeutronSeparationEnergy - theProtonSeparationEnergy;
119
120 separationEnergy[Eta] = 0.;
124
125 separationEnergy[Lambda] = theLambdaSeparationEnergy;
126 separationEnergy[SigmaPlus] = theProtonSeparationEnergy + theLambdaSeparationEnergy - theNeutronSeparationEnergy;
127 separationEnergy[SigmaZero] = theLambdaSeparationEnergy;
128 separationEnergy[SigmaMinus] = theNeutronSeparationEnergy + theLambdaSeparationEnergy - theProtonSeparationEnergy;
129
130 separationEnergy[KPlus] = theProtonSeparationEnergy - theLambdaSeparationEnergy;
131 separationEnergy[KZero] = (theNeutronSeparationEnergy - theLambdaSeparationEnergy);
132 separationEnergy[KZeroBar] = (theLambdaSeparationEnergy - theNeutronSeparationEnergy);
133 separationEnergy[KMinus] = 2.*theNeutronSeparationEnergy - theProtonSeparationEnergy-theLambdaSeparationEnergy;
134
135 separationEnergy[KShort] = (theNeutronSeparationEnergy - theLambdaSeparationEnergy);
136 separationEnergy[KLong] = (theNeutronSeparationEnergy - theLambdaSeparationEnergy);
137
142
144 if (fermiEnergy[Lambda] <= 0.)
146 else
147 fermiMomentum[Lambda]=std::sqrt(std::pow(fermiEnergy[Lambda]+ml,2.0)-ml*ml);
148
152
153 INCL_DEBUG("Table of separation energies [MeV] for A=" << theA << ", Z=" << theZ << ":" << '\n'
154 << " proton: " << separationEnergy[Proton] << '\n'
155 << " neutron: " << separationEnergy[Neutron] << '\n'
156 << " delta++: " << separationEnergy[DeltaPlusPlus] << '\n'
157 << " delta+: " << separationEnergy[DeltaPlus] << '\n'
158 << " delta0: " << separationEnergy[DeltaZero] << '\n'
159 << " delta-: " << separationEnergy[DeltaMinus] << '\n'
160 << " pi+: " << separationEnergy[PiPlus] << '\n'
161 << " pi0: " << separationEnergy[PiZero] << '\n'
162 << " pi-: " << separationEnergy[PiMinus] << '\n'
163 << " eta: " << separationEnergy[Eta] << '\n'
164 << " omega: " << separationEnergy[Omega] << '\n'
165 << " etaprime:" << separationEnergy[EtaPrime] << '\n'
166 << " photon: " << separationEnergy[Photon] << '\n'
167 << " lambda: " << separationEnergy[Lambda] << '\n'
168 << " sigmaplus: " << separationEnergy[SigmaPlus] << '\n'
169 << " sigmazero: " << separationEnergy[SigmaZero] << '\n'
170 << " sigmaminus: " << separationEnergy[SigmaMinus] << '\n'
171 << " kplus: " << separationEnergy[KPlus] << '\n'
172 << " kzero: " << separationEnergy[KZero] << '\n'
173 << " kzerobar: " << separationEnergy[KZeroBar] << '\n'
174 << " kminus: " << separationEnergy[KMinus] << '\n'
175 << " kshort: " << separationEnergy[KShort] << '\n'
176 << " klong: " << separationEnergy[KLong] << '\n'
177 );
178
179 INCL_DEBUG("Table of Fermi energies [MeV] for A=" << theA << ", Z=" << theZ << ":" << '\n'
180 << " proton: " << fermiEnergy[Proton] << '\n'
181 << " neutron: " << fermiEnergy[Neutron] << '\n'
182 << " delta++: " << fermiEnergy[DeltaPlusPlus] << '\n'
183 << " delta+: " << fermiEnergy[DeltaPlus] << '\n'
184 << " delta0: " << fermiEnergy[DeltaZero] << '\n'
185 << " delta-: " << fermiEnergy[DeltaMinus] << '\n'
186 << " lambda: " << fermiEnergy[Lambda] << '\n'
187 << " sigma+: " << fermiEnergy[SigmaPlus] << '\n'
188 << " sigma0: " << fermiEnergy[SigmaZero] << '\n'
189 << " sigma-: " << fermiEnergy[SigmaMinus] << '\n'
190 );
191
192 INCL_DEBUG("Table of Fermi momenta [MeV/c] for A=" << theA << ", Z=" << theZ << ":" << '\n'
193 << " proton: " << fermiMomentum[Proton] << '\n'
194 << " neutron: " << fermiMomentum[Neutron] << '\n'
195 );
196 }
#define INCL_DEBUG(x)
const G4int theA
The mass number of the nucleus.
const G4int theZ
The charge number of the nucleus.
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double pow13(G4double x)
G4ThreadLocal FermiMomentumFn getFermiMomentum
G4ThreadLocal SeparationEnergyFn getSeparationEnergy
Static pointer to the separation-energy function.
G4double getINCLMass(const G4int A, const G4int Z, const G4int S)
Get INCL nuclear mass (in MeV/c^2)

References G4INCL::DeltaMinus, G4INCL::DeltaPlus, G4INCL::DeltaPlusPlus, G4INCL::DeltaZero, G4INCL::Eta, G4INCL::EtaPrime, G4INCL::NuclearPotential::INuclearPotential::fermiEnergy, G4INCL::NuclearPotential::INuclearPotential::fermiMomentum, G4INCL::ParticleTable::getFermiMomentum, G4INCL::ParticleTable::getINCLMass(), G4INCL::ParticleTable::getSeparationEnergy, INCL_DEBUG, G4INCL::KLong, G4INCL::KMinus, G4INCL::KPlus, G4INCL::KShort, G4INCL::KZero, G4INCL::KZeroBar, G4INCL::Lambda, G4INCL::Math::max(), G4INCL::Neutron, G4INCL::Omega, G4INCL::Photon, G4INCL::PiMinus, G4INCL::PiPlus, G4INCL::PiZero, G4INCL::Math::pow13(), G4INCL::Proton, G4INCL::NuclearPotential::INuclearPotential::separationEnergy, G4INCL::SigmaMinus, G4INCL::SigmaPlus, G4INCL::SigmaZero, G4INCL::NuclearPotential::INuclearPotential::theA, G4INCL::NuclearPotential::INuclearPotential::theZ, G4INCL::NuclearPotential::NuclearPotentialIsospin::vDeltaMinus, G4INCL::NuclearPotential::NuclearPotentialIsospin::vDeltaPlus, G4INCL::NuclearPotential::NuclearPotentialIsospin::vDeltaPlusPlus, G4INCL::NuclearPotential::NuclearPotentialIsospin::vDeltaZero, G4INCL::NuclearPotential::NuclearPotentialIsospin::vLambda, G4INCL::NuclearPotential::NuclearPotentialIsospin::vNeutron, G4INCL::NuclearPotential::NuclearPotentialIsospin::vProton, G4INCL::NuclearPotential::NuclearPotentialIsospin::vSigmaMinus, G4INCL::NuclearPotential::NuclearPotentialIsospin::vSigmaPlus, and G4INCL::NuclearPotential::NuclearPotentialIsospin::vSigmaZero.

Referenced by G4INCL::NuclearPotential::NuclearPotentialIsospin::NuclearPotentialIsospin().

Field Documentation

◆ alpha

const G4double G4INCL::NuclearPotential::NuclearPotentialEnergyIsospinSmooth::alpha = 0.223
staticprivate

Slope of the V(T) curve.

Definition at line 67 of file G4INCLNuclearPotentialEnergyIsospinSmooth.hh.

Referenced by computePotentialEnergy().

◆ deltaE

const G4double G4INCL::NuclearPotential::NuclearPotentialEnergyIsospinSmooth::deltaE = 25.
staticprivate

Distance from the cusp where the exponential kicks in.

Definition at line 70 of file G4INCLNuclearPotentialEnergyIsospinSmooth.hh.

Referenced by computePotentialEnergy().

◆ fermiEnergy

std::map<ParticleType,G4double> G4INCL::NuclearPotential::INuclearPotential::fermiEnergy
protectedinherited

◆ fermiMomentum

std::map<ParticleType,G4double> G4INCL::NuclearPotential::INuclearPotential::fermiMomentum
protectedinherited

◆ pionPotential

const G4bool G4INCL::NuclearPotential::INuclearPotential::pionPotential
privateinherited

◆ separationEnergy

std::map<ParticleType,G4double> G4INCL::NuclearPotential::INuclearPotential::separationEnergy
protectedinherited

◆ theA

const G4int G4INCL::NuclearPotential::INuclearPotential::theA
protectedinherited

◆ theZ

const G4int G4INCL::NuclearPotential::INuclearPotential::theZ
protectedinherited

◆ vDeltaMinus

G4double G4INCL::NuclearPotential::NuclearPotentialIsospin::vDeltaMinus
privateinherited

◆ vDeltaPlus

G4double G4INCL::NuclearPotential::NuclearPotentialIsospin::vDeltaPlus
privateinherited

◆ vDeltaPlusPlus

G4double G4INCL::NuclearPotential::NuclearPotentialIsospin::vDeltaPlusPlus
privateinherited

◆ vDeltaZero

G4double G4INCL::NuclearPotential::NuclearPotentialIsospin::vDeltaZero
privateinherited

◆ vKMinus

G4double G4INCL::NuclearPotential::INuclearPotential::vKMinus
privateinherited

◆ vKMinusDefault

const G4double G4INCL::NuclearPotential::INuclearPotential::vKMinusDefault = 60.
staticprivateinherited

◆ vKPlus

G4double G4INCL::NuclearPotential::INuclearPotential::vKPlus
privateinherited

◆ vKPlusDefault

const G4double G4INCL::NuclearPotential::INuclearPotential::vKPlusDefault = -25.
staticprivateinherited

◆ vKZero

G4double G4INCL::NuclearPotential::INuclearPotential::vKZero
privateinherited

◆ vKZeroBar

G4double G4INCL::NuclearPotential::INuclearPotential::vKZeroBar
privateinherited

◆ vLambda

G4double G4INCL::NuclearPotential::NuclearPotentialIsospin::vLambda
privateinherited

◆ vNeutron

G4double G4INCL::NuclearPotential::NuclearPotentialIsospin::vNeutron
privateinherited

◆ vPiMinus

G4double G4INCL::NuclearPotential::INuclearPotential::vPiMinus
privateinherited

◆ vPionDefault

const G4double G4INCL::NuclearPotential::INuclearPotential::vPionDefault = 30.6
staticprivateinherited

◆ vPiPlus

G4double G4INCL::NuclearPotential::INuclearPotential::vPiPlus
privateinherited

◆ vPiZero

G4double G4INCL::NuclearPotential::INuclearPotential::vPiZero
privateinherited

◆ vProton

G4double G4INCL::NuclearPotential::NuclearPotentialIsospin::vProton
privateinherited

◆ vSigmaMinus

G4double G4INCL::NuclearPotential::NuclearPotentialIsospin::vSigmaMinus
privateinherited

◆ vSigmaPlus

G4double G4INCL::NuclearPotential::NuclearPotentialIsospin::vSigmaPlus
privateinherited

◆ vSigmaZero

G4double G4INCL::NuclearPotential::NuclearPotentialIsospin::vSigmaZero
privateinherited

The documentation for this class was generated from the following files: