Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
G4PreCompoundHe3 Class Reference

#include <G4PreCompoundHe3.hh>

Inheritance diagram for G4PreCompoundHe3:
G4PreCompoundIon G4PreCompoundFragment G4VPreCompoundFragment

Public Member Functions

G4double CalcEmissionProbability (const G4Fragment &aFragment) override
 
 G4PreCompoundHe3 ()
 
 G4PreCompoundHe3 (const G4PreCompoundHe3 &right)=delete
 
G4int GetA () const
 
G4double GetBindingEnergy () const
 
G4double GetEmissionProbability () const
 
G4double GetEnergyThreshold () const
 
const G4LorentzVectorGetMomentum () const
 
G4double GetNuclearMass () const
 
G4ReactionProductGetReactionProduct () const
 
G4int GetRestA () const
 
G4double GetRestNuclearMass () const
 
G4int GetRestZ () const
 
G4int GetZ () const
 
void Initialize (const G4Fragment &aFragment)
 
G4bool IsItPossible (const G4Fragment &aFragment) const
 
G4bool operator!= (const G4PreCompoundFragment &right) const =delete
 
G4bool operator!= (const G4PreCompoundHe3 &right) const =delete
 
G4bool operator!= (const G4PreCompoundIon &right) const =delete
 
G4bool operator!= (const G4VPreCompoundFragment &right) const =delete
 
const G4PreCompoundHe3operator= (const G4PreCompoundHe3 &right)=delete
 
G4bool operator== (const G4PreCompoundFragment &right) const =delete
 
G4bool operator== (const G4PreCompoundHe3 &right) const =delete
 
G4bool operator== (const G4PreCompoundIon &right) const =delete
 
G4bool operator== (const G4VPreCompoundFragment &right) const =delete
 
G4double SampleKineticEnergy (const G4Fragment &aFragment) override
 
void SetMomentum (const G4LorentzVector &value)
 
void SetOPTxs (G4int)
 
void UseSICB (G4bool)
 
virtual ~G4PreCompoundHe3 ()
 

Protected Member Functions

G4double CoalescenceFactor (G4int A) const override
 
G4double CrossSection (G4double ekin) const
 
G4double FactorialFactor (G4int N, G4int P) const override
 
G4double GetAlpha () const override
 
G4double GetBeta () const override
 
G4double GetRj (G4int NumberParticles, G4int NumberCharged) const override
 
G4double ProbabilityDistributionFunction (G4double eKin, const G4Fragment &) override
 

Protected Attributes

G4NuclearLevelDatafNucData
 
G4Powg4calc
 
G4int OPTxs
 
G4int theA
 
G4double theBindingEnergy
 
G4double theCoulombBarrier
 
G4double theEmissionProbability
 
G4int theFragA
 
G4int theFragZ
 
G4double theMass
 
G4double theMaxKinEnergy
 
G4double theMinKinEnergy
 
G4DeexPrecoParameterstheParameters
 
G4double theReducedMass
 
G4int theResA
 
G4double theResA13
 
G4double theResMass
 
G4int theResZ
 
G4int theZ
 
G4bool useSICB
 

Private Member Functions

G4double GetOpt0 (G4double ekin) const
 
G4double IntegrateEmissionProbability (G4double low, G4double up, const G4Fragment &aFragment)
 

Private Attributes

G4double fact
 
G4int index
 
G4double muu
 
const G4ParticleDefinitionparticle
 
G4double probmax
 
G4VCoulombBarriertheCoulombBarrierPtr
 
G4He3CoulombBarrier theHe3CoulombBarrier
 
G4LorentzVector theMomentum
 

Detailed Description

Definition at line 41 of file G4PreCompoundHe3.hh.

Constructor & Destructor Documentation

◆ G4PreCompoundHe3() [1/2]

G4PreCompoundHe3::G4PreCompoundHe3 ( )

Definition at line 47 of file G4PreCompoundHe3.cc.

49{}
static G4He3 * He3()
Definition: G4He3.cc:93
G4He3CoulombBarrier theHe3CoulombBarrier
G4PreCompoundIon(const G4ParticleDefinition *, G4VCoulombBarrier *aCoulombBarrier)

◆ ~G4PreCompoundHe3()

G4PreCompoundHe3::~G4PreCompoundHe3 ( )
virtual

Definition at line 51 of file G4PreCompoundHe3.cc.

52{}

◆ G4PreCompoundHe3() [2/2]

G4PreCompoundHe3::G4PreCompoundHe3 ( const G4PreCompoundHe3 right)
delete

Member Function Documentation

◆ CalcEmissionProbability()

G4double G4PreCompoundFragment::CalcEmissionProbability ( const G4Fragment aFragment)
overridevirtualinherited

Implements G4VPreCompoundFragment.

Definition at line 55 of file G4PreCompoundFragment.cc.

57{
58 //G4cout << theCoulombBarrier << " " << GetMaximalKineticEnergy() << G4endl;
59 // If theCoulombBarrier effect is included in the emission probabilities
60 // Coulomb barrier is the lower limit of integration over kinetic energy
61
63
64 if (theMaxKinEnergy <= theMinKinEnergy) { return 0.0; }
65
66 // compute power once
67 if(0 < index) {
69 }
70
73 /*
74 G4cout << "## G4PreCompoundFragment::CalcEmisProb "
75 << "Z= " << aFragment.GetZ_asInt()
76 << " A= " << aFragment.GetA_asInt()
77 << " Elow= " << LowerLimit/MeV
78 << " Eup= " << UpperLimit/MeV
79 << " prob= " << theEmissionProbability
80 << G4endl;
81 */
83}
static G4double ComputePowerParameter(G4int resA, G4int idx)
G4double IntegrateEmissionProbability(G4double low, G4double up, const G4Fragment &aFragment)

References G4KalbachCrossSection::ComputePowerParameter(), G4PreCompoundFragment::index, G4PreCompoundFragment::IntegrateEmissionProbability(), G4PreCompoundFragment::muu, G4VPreCompoundFragment::theEmissionProbability, G4VPreCompoundFragment::theMaxKinEnergy, G4VPreCompoundFragment::theMinKinEnergy, and G4VPreCompoundFragment::theResA.

◆ CoalescenceFactor()

G4double G4PreCompoundHe3::CoalescenceFactor ( G4int  A) const
overrideprotectedvirtual

Implements G4PreCompoundIon.

Definition at line 59 of file G4PreCompoundHe3.cc.

60{
61 return 243.0/static_cast<G4double>(A*A);
62}
double G4double
Definition: G4Types.hh:83
const G4double A[17]

References A.

◆ CrossSection()

G4double G4PreCompoundFragment::CrossSection ( G4double  ekin) const
protectedinherited

Definition at line 113 of file G4PreCompoundFragment.cc.

114{
115 G4double res;
116 if(OPTxs == 0 || (OPTxs == 4 && theMaxKinEnergy < 10.)) {
117 res = GetOpt0(ekin);
118
119 } else if(OPTxs <= 2) {
122 theResA13, muu,
123 index, theZ, theResA);
124
125 } else {
128 theZ, theA, theResA);
129 }
130 return res;
131}
static G4double ComputeCrossSection(G4double K, G4double cb, G4double resA13, G4double amu1, G4int idx, G4int Z, G4int resA)
static G4double ComputeCrossSection(G4double K, G4double cb, G4double resA13, G4double amu1, G4int idx, G4int Z, G4int A, G4int resA)
G4double GetOpt0(G4double ekin) const

References G4KalbachCrossSection::ComputeCrossSection(), G4ChatterjeeCrossSection::ComputeCrossSection(), G4PreCompoundFragment::GetOpt0(), G4PreCompoundFragment::index, G4PreCompoundFragment::muu, G4VPreCompoundFragment::OPTxs, G4VPreCompoundFragment::theA, G4VPreCompoundFragment::theCoulombBarrier, G4VPreCompoundFragment::theMaxKinEnergy, G4VPreCompoundFragment::theResA, G4VPreCompoundFragment::theResA13, and G4VPreCompoundFragment::theZ.

Referenced by G4PreCompoundIon::ProbabilityDistributionFunction(), and G4PreCompoundNucleon::ProbabilityDistributionFunction().

◆ FactorialFactor()

G4double G4PreCompoundHe3::FactorialFactor ( G4int  N,
G4int  P 
) const
overrideprotectedvirtual

Implements G4PreCompoundIon.

Definition at line 54 of file G4PreCompoundHe3.cc.

55{
56 return static_cast<G4double>(((N-3)*(P-2)*(N-2))*((P-1)*(N-1)*P))/12.0;
57}
static double P[]

References P.

◆ GetA()

G4int G4VPreCompoundFragment::GetA ( ) const
inlineinherited

◆ GetAlpha()

G4double G4PreCompoundHe3::GetAlpha ( ) const
overrideprotectedvirtual

Implements G4VPreCompoundFragment.

Definition at line 74 of file G4PreCompoundHe3.cc.

75{
76 G4double C = 0.0;
77 if (theFragZ <= 30)
78 {
79 C = 0.10;
80 }
81 else if (theFragZ <= 50)
82 {
83 C = 0.1 - (theFragZ - 30)*0.001;
84 }
85 else if (theFragZ < 70)
86 {
87 C = 0.08 - (theFragZ - 50)*0.001;
88 }
89 else
90 {
91 C = 0.06;
92 }
93 return 1.0 + C*(4.0/3.0);
94}
G4double C(G4double temp)

References C(), and G4VPreCompoundFragment::theFragZ.

◆ GetBeta()

G4double G4PreCompoundIon::GetBeta ( ) const
overrideprotectedvirtualinherited

Implements G4VPreCompoundFragment.

Definition at line 101 of file G4PreCompoundIon.cc.

102{
103 return -theCoulombBarrier;
104}

References G4VPreCompoundFragment::theCoulombBarrier.

◆ GetBindingEnergy()

G4double G4VPreCompoundFragment::GetBindingEnergy ( ) const
inlineinherited

◆ GetEmissionProbability()

G4double G4VPreCompoundFragment::GetEmissionProbability ( ) const
inlineinherited

◆ GetEnergyThreshold()

G4double G4VPreCompoundFragment::GetEnergyThreshold ( ) const
inlineinherited

◆ GetMomentum()

const G4LorentzVector & G4VPreCompoundFragment::GetMomentum ( ) const
inlineinherited

◆ GetNuclearMass()

G4double G4VPreCompoundFragment::GetNuclearMass ( ) const
inlineinherited

◆ GetOpt0()

G4double G4PreCompoundFragment::GetOpt0 ( G4double  ekin) const
privateinherited

Definition at line 133 of file G4PreCompoundFragment.cc.

135{
137 // cross section is now given in mb (r0 is in mm) for the sake of consistency
138 //with the rest of the options
139 return 1.e+25*CLHEP::pi*r0*r0*theResA13*GetAlpha()*(1.0 + GetBeta()/ekin);
140}
G4DeexPrecoParameters * theParameters
virtual G4double GetBeta() const =0
virtual G4double GetAlpha() const =0
static constexpr double pi
Definition: SystemOfUnits.h:55

References CLHEP::pi.

Referenced by G4PreCompoundFragment::CrossSection().

◆ GetReactionProduct()

G4ReactionProduct * G4VPreCompoundFragment::GetReactionProduct ( ) const
inlineinherited

◆ GetRestA()

G4int G4VPreCompoundFragment::GetRestA ( ) const
inlineinherited

◆ GetRestNuclearMass()

G4double G4VPreCompoundFragment::GetRestNuclearMass ( ) const
inlineinherited

◆ GetRestZ()

G4int G4VPreCompoundFragment::GetRestZ ( ) const
inlineinherited

◆ GetRj()

G4double G4PreCompoundHe3::GetRj ( G4int  NumberParticles,
G4int  NumberCharged 
) const
overrideprotectedvirtual

Implements G4PreCompoundIon.

Definition at line 64 of file G4PreCompoundHe3.cc.

65{
66 G4double rj = 0.0;
67 if(nCharged >=2 && (nParticles-nCharged) >= 1) {
68 G4double denominator = (nParticles*(nParticles-1)*(nParticles-2));
69 rj = (3*nCharged*(nCharged-1)*(nParticles-nCharged))/denominator;
70 }
71 return rj;
72}

◆ GetZ()

G4int G4VPreCompoundFragment::GetZ ( ) const
inlineinherited

◆ Initialize()

void G4VPreCompoundFragment::Initialize ( const G4Fragment aFragment)
inherited

Definition at line 79 of file G4VPreCompoundFragment.cc.

80{
81 theFragA = aFragment.GetA_asInt();
82 theFragZ = aFragment.GetZ_asInt();
85
87 if ((theResA < theResZ) || (theResA < theA) || (theResZ < theZ)) {
88 return;
89 }
90
93 GetCoulombBarrier(theResA,theResZ,aFragment.GetExcitationEnergy());
94
96
97 // Calculate masses
100
101 // Compute Binding Energies for fragments
102 // needed to separate a fragment from the nucleus
104
105 // Compute Maximal Kinetic Energy which can be carried by fragments
106 // after separation - the true assimptotic value
107 G4double Ecm = aFragment.GetMomentum().m();
108 G4double twoEcm = Ecm + Ecm;
110 /twoEcm - theMass,0.0);
111 theMinKinEnergy = (elim == 0.0) ? 0.0 :
112 std::max(((theMass+elim)*(twoEcm-theMass-elim) +
113 theMass*theMass)/twoEcm - theMass,0.0);
114}
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:304
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:299
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:323
G4int GetZ_asInt() const
Definition: G4Fragment.hh:276
G4int GetA_asInt() const
Definition: G4Fragment.hh:271
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double Z13(G4int Z) const
Definition: G4Pow.hh:123
G4VCoulombBarrier * theCoulombBarrierPtr
T max(const T t1, const T t2)
brief Return the largest of the two arguments

References G4VPreCompoundFragment::g4calc, G4Fragment::GetA_asInt(), G4Fragment::GetExcitationEnergy(), G4Fragment::GetGroundStateMass(), G4Fragment::GetMomentum(), G4NucleiProperties::GetNuclearMass(), G4Fragment::GetZ_asInt(), CLHEP::HepLorentzVector::m(), G4INCL::Math::max(), G4VPreCompoundFragment::OPTxs, G4VPreCompoundFragment::theA, G4VPreCompoundFragment::theBindingEnergy, G4VPreCompoundFragment::theCoulombBarrier, G4VPreCompoundFragment::theCoulombBarrierPtr, G4VPreCompoundFragment::theFragA, G4VPreCompoundFragment::theFragZ, G4VPreCompoundFragment::theMass, G4VPreCompoundFragment::theMaxKinEnergy, G4VPreCompoundFragment::theMinKinEnergy, G4VPreCompoundFragment::theReducedMass, G4VPreCompoundFragment::theResA, G4VPreCompoundFragment::theResA13, G4VPreCompoundFragment::theResMass, G4VPreCompoundFragment::theResZ, G4VPreCompoundFragment::theZ, and G4Pow::Z13().

◆ IntegrateEmissionProbability()

G4double G4PreCompoundFragment::IntegrateEmissionProbability ( G4double  low,
G4double  up,
const G4Fragment aFragment 
)
privateinherited

Definition at line 85 of file G4PreCompoundFragment.cc.

88{
89 static const G4double den = 1.0/CLHEP::MeV;
90 G4double del = (up - low);
91 G4int nbins = del*den;
92 nbins = std::max(nbins, 4);
93 del /= static_cast<G4double>(nbins);
94 G4double e = low + 0.5*del;
96 //G4cout << " 0. e= " << e << " y= " << probmax << G4endl;
97
98 G4double sum = probmax;
99 for (G4int i=1; i<nbins; ++i) {
100 e += del;
101
102 G4double y = ProbabilityDistributionFunction(e, aFragment);
104 sum += y;
105 if(y < sum*0.01) { break; }
106 //G4cout <<" "<<i<<". e= "<<e<<" y= "<<y<<" sum= "<<sum<< G4endl;
107 }
108 sum *= del;
109 //G4cout << "Evap prob: " << sum << " probmax= " << probmax << G4endl;
110 return sum;
111}
int G4int
Definition: G4Types.hh:85
virtual G4double ProbabilityDistributionFunction(G4double ekin, const G4Fragment &aFragment)=0
static constexpr double MeV

References G4INCL::Math::max(), CLHEP::MeV, G4PreCompoundFragment::ProbabilityDistributionFunction(), and G4PreCompoundFragment::probmax.

Referenced by G4PreCompoundFragment::CalcEmissionProbability().

◆ IsItPossible()

G4bool G4VPreCompoundFragment::IsItPossible ( const G4Fragment aFragment) const
inlineinherited

◆ operator!=() [1/4]

G4bool G4PreCompoundFragment::operator!= ( const G4PreCompoundFragment right) const
deleteinherited

◆ operator!=() [2/4]

G4bool G4PreCompoundHe3::operator!= ( const G4PreCompoundHe3 right) const
delete

◆ operator!=() [3/4]

G4bool G4PreCompoundIon::operator!= ( const G4PreCompoundIon right) const
deleteinherited

◆ operator!=() [4/4]

G4bool G4VPreCompoundFragment::operator!= ( const G4VPreCompoundFragment right) const
deleteinherited

◆ operator=()

const G4PreCompoundHe3 & G4PreCompoundHe3::operator= ( const G4PreCompoundHe3 right)
delete

◆ operator==() [1/4]

G4bool G4PreCompoundFragment::operator== ( const G4PreCompoundFragment right) const
deleteinherited

◆ operator==() [2/4]

G4bool G4PreCompoundHe3::operator== ( const G4PreCompoundHe3 right) const
delete

◆ operator==() [3/4]

G4bool G4PreCompoundIon::operator== ( const G4PreCompoundIon right) const
deleteinherited

◆ operator==() [4/4]

G4bool G4VPreCompoundFragment::operator== ( const G4VPreCompoundFragment right) const
deleteinherited

◆ ProbabilityDistributionFunction()

G4double G4PreCompoundIon::ProbabilityDistributionFunction ( G4double  eKin,
const G4Fragment aFragment 
)
overrideprotectedvirtualinherited

Implements G4PreCompoundFragment.

Definition at line 59 of file G4PreCompoundIon.cc.

62{
63 G4double efinal = eKin + theBindingEnergy;
64 if(efinal <= 0.0 ) { return 0.0; }
65
66 G4double U = aFragment.GetExcitationEnergy();
67 G4int P = aFragment.GetNumberOfParticles();
68 G4int H = aFragment.GetNumberOfHoles();
69 G4int A = GetA();
70 G4int N = P + H;
71
72 static const G4double sixoverpi2 = 6.0/CLHEP::pi2;
74 G4double g1 = sixoverpi2*fNucData->GetLevelDensity(theResZ, theResA, 0.0);
75
76 G4double gj = g1;
77
78 G4double A0 = (P*P+H*H+P-3*H)/(4.0*g0);
79 G4double A1 = std::max(0.0,(A0*g0 + A*(A-2*P-1)*0.25)/g1);
80
81 G4double E0 = U - A0;
82 if (E0 <= 0.0) { return 0.0; }
83
84 G4double E1 = std::max(0.0,theMaxKinEnergy - eKin - A1);
85
86 G4double Aj = A*(A+1)/(4.0*gj);
87 G4double Ej = std::max(0.0,efinal - Aj);
88
89 G4double rj = GetRj(P, aFragment.GetNumberOfCharged());
90 G4double xs = CrossSection(eKin);
91
92 G4double pA = fact*eKin*xs*rj
94 * std::sqrt(2.0/(theReducedMass*efinal))
95 * g4calc->powN(g1*E1/(g0*E0), N-A-1)
96 * g4calc->powN(gj*Ej/(g0*E0), A-1)*gj*g1/(g0*g0*E0*theResA);
97
98 return pA;
99}
G4int GetNumberOfParticles() const
Definition: G4Fragment.hh:361
G4int GetNumberOfHoles() const
Definition: G4Fragment.hh:381
G4int GetNumberOfCharged() const
Definition: G4Fragment.hh:366
G4double GetLevelDensity(G4int Z, G4int A, G4double U)
G4double powN(G4double x, G4int n) const
Definition: G4Pow.cc:166
G4double CrossSection(G4double ekin) const
virtual G4double FactorialFactor(G4int N, G4int P) const =0
virtual G4double CoalescenceFactor(G4int A) const =0
virtual G4double GetRj(G4int NumberParticles, G4int NumberCharged) const =0
G4NuclearLevelData * fNucData
static constexpr double pi2
Definition: SystemOfUnits.h:58

References A, G4PreCompoundIon::CoalescenceFactor(), G4PreCompoundFragment::CrossSection(), G4PreCompoundIon::fact, G4PreCompoundIon::FactorialFactor(), G4VPreCompoundFragment::fNucData, G4VPreCompoundFragment::g4calc, G4VPreCompoundFragment::GetA(), G4Fragment::GetExcitationEnergy(), G4NuclearLevelData::GetLevelDensity(), G4Fragment::GetNumberOfCharged(), G4Fragment::GetNumberOfHoles(), G4Fragment::GetNumberOfParticles(), G4PreCompoundIon::GetRj(), G4INCL::Math::max(), P, CLHEP::pi2, G4Pow::powN(), G4VPreCompoundFragment::theBindingEnergy, G4VPreCompoundFragment::theFragA, G4VPreCompoundFragment::theFragZ, G4VPreCompoundFragment::theMaxKinEnergy, G4VPreCompoundFragment::theReducedMass, G4VPreCompoundFragment::theResA, and G4VPreCompoundFragment::theResZ.

◆ SampleKineticEnergy()

G4double G4PreCompoundFragment::SampleKineticEnergy ( const G4Fragment aFragment)
overridevirtualinherited

Implements G4VPreCompoundFragment.

Definition at line 142 of file G4PreCompoundFragment.cc.

143{
145 static const G4double toler = 1.25;
146 probmax *= toler;
147 G4double prob, T(0.0);
148 CLHEP::HepRandomEngine* rndm = G4Random::getTheEngine();
149 G4int i;
150 for(i=0; i<100; ++i) {
151 T = theMinKinEnergy + delta*rndm->flat();
152 prob = ProbabilityDistributionFunction(T, fragment);
153 /*
154 if(prob > probmax) {
155 G4cout << "G4PreCompoundFragment WARNING: prob= " << prob
156 << " probmax= " << probmax << G4endl;
157 G4cout << "i= " << i << " Z= " << theZ << " A= " << theA
158 << " resZ= " << theResZ << " resA= " << theResA << "\n"
159 << " T= " << T << " Tmax= " << theMaxKinEnergy
160 << " Tmin= " << limit
161 << G4endl;
162 for(G4int i=0; i<N; ++i) { G4cout << " " << probability[i]; }
163 G4cout << G4endl;
164 }
165 */
166 // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
167 if(probmax*rndm->flat() <= prob) { break; }
168 }
169 /*
170 G4cout << "G4PreCompoundFragment: i= " << i << " Z= " << theZ
171 << " A= " << theA <<" T(MeV)= " << T << " Emin(MeV)= "
172 << theMinKinEnergy << " Emax= " << theMaxKinEnergy << G4endl;
173 */
174 return T;
175}
virtual double flat()=0

References CLHEP::HepRandomEngine::flat(), G4PreCompoundFragment::ProbabilityDistributionFunction(), G4PreCompoundFragment::probmax, G4VPreCompoundFragment::theMaxKinEnergy, G4VPreCompoundFragment::theMinKinEnergy, and anonymous_namespace{G4QuasiElRatios.cc}::toler.

◆ SetMomentum()

void G4VPreCompoundFragment::SetMomentum ( const G4LorentzVector value)
inlineinherited

◆ SetOPTxs()

void G4VPreCompoundFragment::SetOPTxs ( G4int  )
inlineinherited

◆ UseSICB()

void G4VPreCompoundFragment::UseSICB ( G4bool  )
inlineinherited

Field Documentation

◆ fact

G4double G4PreCompoundIon::fact
privateinherited

◆ fNucData

G4NuclearLevelData* G4VPreCompoundFragment::fNucData
protectedinherited

◆ g4calc

G4Pow* G4VPreCompoundFragment::g4calc
protectedinherited

◆ index

G4int G4PreCompoundFragment::index
privateinherited

◆ muu

G4double G4PreCompoundFragment::muu
privateinherited

◆ OPTxs

G4int G4VPreCompoundFragment::OPTxs
protectedinherited

◆ particle

const G4ParticleDefinition* G4VPreCompoundFragment::particle
privateinherited

◆ probmax

G4double G4PreCompoundFragment::probmax
privateinherited

◆ theA

G4int G4VPreCompoundFragment::theA
protectedinherited

◆ theBindingEnergy

G4double G4VPreCompoundFragment::theBindingEnergy
protectedinherited

◆ theCoulombBarrier

G4double G4VPreCompoundFragment::theCoulombBarrier
protectedinherited

◆ theCoulombBarrierPtr

G4VCoulombBarrier* G4VPreCompoundFragment::theCoulombBarrierPtr
privateinherited

Definition at line 132 of file G4VPreCompoundFragment.hh.

Referenced by G4VPreCompoundFragment::Initialize().

◆ theEmissionProbability

G4double G4VPreCompoundFragment::theEmissionProbability
protectedinherited

◆ theFragA

G4int G4VPreCompoundFragment::theFragA
protectedinherited

◆ theFragZ

G4int G4VPreCompoundFragment::theFragZ
protectedinherited

◆ theHe3CoulombBarrier

G4He3CoulombBarrier G4PreCompoundHe3::theHe3CoulombBarrier
private

Definition at line 67 of file G4PreCompoundHe3.hh.

◆ theMass

G4double G4VPreCompoundFragment::theMass
protectedinherited

◆ theMaxKinEnergy

G4double G4VPreCompoundFragment::theMaxKinEnergy
protectedinherited

◆ theMinKinEnergy

G4double G4VPreCompoundFragment::theMinKinEnergy
protectedinherited

◆ theMomentum

G4LorentzVector G4VPreCompoundFragment::theMomentum
privateinherited

Definition at line 134 of file G4VPreCompoundFragment.hh.

◆ theParameters

G4DeexPrecoParameters* G4VPreCompoundFragment::theParameters
protectedinherited

◆ theReducedMass

G4double G4VPreCompoundFragment::theReducedMass
protectedinherited

◆ theResA

G4int G4VPreCompoundFragment::theResA
protectedinherited

◆ theResA13

G4double G4VPreCompoundFragment::theResA13
protectedinherited

◆ theResMass

G4double G4VPreCompoundFragment::theResMass
protectedinherited

Definition at line 153 of file G4VPreCompoundFragment.hh.

Referenced by G4VPreCompoundFragment::Initialize().

◆ theResZ

G4int G4VPreCompoundFragment::theResZ
protectedinherited

◆ theZ

G4int G4VPreCompoundFragment::theZ
protectedinherited

◆ useSICB

G4bool G4VPreCompoundFragment::useSICB
protectedinherited

Definition at line 163 of file G4VPreCompoundFragment.hh.


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