Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PreCompoundFragment Class Reference

#include <G4PreCompoundFragment.hh>

Inheritance diagram for G4PreCompoundFragment:

Public Member Functions

 G4PreCompoundFragment (const G4ParticleDefinition *, G4VCoulombBarrier *aCoulombBarrier)
 ~G4PreCompoundFragment () override=default
G4double CrossSection (G4double ekin)
G4double RecentXS () const
 G4PreCompoundFragment (const G4PreCompoundFragment &right)=delete
const G4PreCompoundFragmentoperator= (const G4PreCompoundFragment &right)=delete
G4bool operator== (const G4PreCompoundFragment &right) const =delete
G4bool operator!= (const G4PreCompoundFragment &right) const =delete
Public Member Functions inherited from G4VPreCompoundFragment
 G4VPreCompoundFragment (const G4ParticleDefinition *, G4VCoulombBarrier *)
 ~G4VPreCompoundFragment () override
G4bool Initialize (const G4Fragment &aFragment)
virtual G4double CalcEmissionProbability (const G4Fragment &)
virtual G4double SampleKineticEnergy (const G4Fragment &)
G4double ProbabilityDensityFunction (G4double energy) override
G4ReactionProductGetReactionProduct () const
G4int GetA () const
G4int GetZ () const
G4int GetRestA () const
G4int GetRestZ () const
G4double GetBindingEnergy () const
G4double GetEnergyThreshold () const
G4double GetEmissionProbability () const
G4double GetNuclearMass () const
G4double GetRestNuclearMass () const
const G4LorentzVectorGetMomentum () const
void SetMomentum (const G4LorentzVector &lv)
void SetOPTxs (G4int)
void UseSICB (G4bool use)
 G4VPreCompoundFragment (const G4VPreCompoundFragment &right)=delete
const G4VPreCompoundFragmentoperator= (const G4VPreCompoundFragment &right)=delete
G4bool operator== (const G4VPreCompoundFragment &right) const =delete
G4bool operator!= (const G4VPreCompoundFragment &right) const =delete

Additional Inherited Members

Protected Member Functions inherited from G4VPreCompoundFragment
virtual G4double ProbabilityDistributionFunction (G4double, const G4Fragment &)
virtual G4double GetAlpha () const =0
virtual G4double GetBeta () const
Protected Attributes inherited from G4VPreCompoundFragment
G4NuclearLevelDatafNucData
G4DeexPrecoParameterstheParameters
G4Powg4calc
G4InterfaceToXSfXSection {nullptr}
const G4FragmentpFragment {nullptr}
G4int theA
G4int theZ
G4int theResA {0}
G4int theResZ {0}
G4int theFragA {0}
G4int theFragZ {0}
G4int OPTxs
G4int index {0}
G4double theResA13 {0.0}
G4double theBindingEnergy {0.0}
G4double theMinKinEnergy {0.0}
G4double theMaxKinEnergy {0.0}
G4double theResMass {0.0}
G4double theReducedMass {0.0}
G4double theMass
G4double theEmissionProbability {0.0}
G4double theCoulombBarrier {0.0}
G4bool useSICB {true}

Detailed Description

Definition at line 43 of file G4PreCompoundFragment.hh.

Constructor & Destructor Documentation

◆ G4PreCompoundFragment() [1/2]

G4PreCompoundFragment::G4PreCompoundFragment ( const G4ParticleDefinition * p,
G4VCoulombBarrier * aCoulombBarrier )

Definition at line 45 of file G4PreCompoundFragment.cc.

47 : G4VPreCompoundFragment(p, aCoulBarrier)
48{}
G4VPreCompoundFragment(const G4ParticleDefinition *, G4VCoulombBarrier *)

Referenced by G4PreCompoundFragment(), G4PreCompoundIon::G4PreCompoundIon(), G4PreCompoundNucleon::G4PreCompoundNucleon(), operator!=(), operator=(), and operator==().

◆ ~G4PreCompoundFragment()

G4PreCompoundFragment::~G4PreCompoundFragment ( )
overridedefault

◆ G4PreCompoundFragment() [2/2]

G4PreCompoundFragment::G4PreCompoundFragment ( const G4PreCompoundFragment & right)
delete

Member Function Documentation

◆ CrossSection()

G4double G4PreCompoundFragment::CrossSection ( G4double ekin)

Definition at line 50 of file G4PreCompoundFragment.cc.

51{
52 // compute power once
53 if (OPTxs > 1 && 0 < index && theResA != lastA) {
54 lastA = theResA;
56 }
57 if (OPTxs == 0) {
58 recentXS = GetOpt0(ekin);
59 } else if (OPTxs == 1) {
60 G4int Z = std::min(theResZ, ZMAXNUCLEARDATA);
61 const G4double lim = 2*CLHEP::MeV;
62 const G4double Kmin = 20*CLHEP::keV;
63 G4double e = lowEnergyLimitMeV[Z];
64 if (e == 0.0) { e = lim; }
65 G4double K = std::max(ekin, Kmin);
66 e = std::max(e, K);
67 recentXS = fXSection->GetElementCrossSection(e, Z)/CLHEP::millibarn;
68 recentXS *= GetAlpha()*(1.0 + GetBeta()/K);
69
70 } else if (OPTxs == 2) {
73 theResA13, muu,
74 index, theZ, theResA);
75
76 } else {
78 theResA13, muu, index,
79 theZ, theA, theResA);
80 }
81 return recentXS;
82}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
static G4double ComputeCrossSection(G4double K, G4double cb, G4double resA13, G4double amu1, G4int idx, G4int Z, G4int resA)
static G4double ComputePowerParameter(G4int resA, G4int idx)
static G4double ComputeCrossSection(G4double K, G4double cb, G4double resA13, G4double amu1, G4int idx, G4int Z, G4int A, G4int resA)
virtual G4double GetAlpha() const =0
virtual G4double GetBeta() const

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

◆ operator!=()

G4bool G4PreCompoundFragment::operator!= ( const G4PreCompoundFragment & right) const
delete

◆ operator=()

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

◆ operator==()

G4bool G4PreCompoundFragment::operator== ( const G4PreCompoundFragment & right) const
delete

◆ RecentXS()

G4double G4PreCompoundFragment::RecentXS ( ) const
inline

Definition at line 56 of file G4PreCompoundFragment.hh.

56{ return recentXS; };

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