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

#include <G4VPreCompoundFragment.hh>

Inheritance diagram for G4VPreCompoundFragment:

Public Member Functions

 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

Protected Member Functions

virtual G4double ProbabilityDistributionFunction (G4double, const G4Fragment &)
virtual G4double GetAlpha () const =0
virtual G4double GetBeta () const

Protected Attributes

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}

Friends

std::ostream & operator<< (std::ostream &out, const G4VPreCompoundFragment *theFragment)
std::ostream & operator<< (std::ostream &out, const G4VPreCompoundFragment &theFragment)

Detailed Description

Definition at line 57 of file G4VPreCompoundFragment.hh.

Constructor & Destructor Documentation

◆ G4VPreCompoundFragment() [1/2]

G4VPreCompoundFragment::G4VPreCompoundFragment ( const G4ParticleDefinition * part,
G4VCoulombBarrier * aCoulombBarrier )
explicit

Definition at line 41 of file G4VPreCompoundFragment.cc.

43 : theA(part->GetBaryonNumber()),
44 theZ(G4lrint(part->GetPDGCharge()/CLHEP::eplus)),
45 particle(part),
46 theCoulombBarrierPtr(aCoulombBarrier)
47{
48 theMass = particle->GetPDGMass();
50 theParameters = fNucData->GetParameters();
51 OPTxs = theParameters->GetDeexModelType();
53
54 if (1 == theZ && 1 == theA) { index = 1; }
55 else if (1 == theZ && 2 == theA) { index = 2; }
56 else if (1 == theZ && 3 == theA) { index = 3; }
57 else if (2 == theZ && 3 == theA) { index = 4; }
58 else if (2 == theZ && 4 == theA) { index = 5; }
59
60 if (OPTxs == 1) {
61 fXSection = new G4InterfaceToXS(particle, index);
62 }
63 InitialiseIntegrator(0.005, 0.25, 1.10, 0.5*CLHEP::MeV,
64 0.2*CLHEP::MeV, 5*CLHEP::MeV);
65}
static G4NuclearLevelData * GetInstance()
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4DeexPrecoParameters * theParameters
void InitialiseIntegrator(G4double accuracy, G4double fact1, G4double fact2, G4double de, G4double dmin, G4double dmax)
int G4lrint(double ad)
Definition templates.hh:134

Referenced by G4HETCFragment::G4HETCFragment(), G4PreCompoundFragment::G4PreCompoundFragment(), G4VPreCompoundFragment(), operator!=(), operator<<, operator<<, operator=(), and operator==().

◆ ~G4VPreCompoundFragment()

G4VPreCompoundFragment::~G4VPreCompoundFragment ( )
override

Definition at line 67 of file G4VPreCompoundFragment.cc.

68{
69 delete theCoulombBarrierPtr;
70 delete fXSection;
71}

◆ G4VPreCompoundFragment() [2/2]

G4VPreCompoundFragment::G4VPreCompoundFragment ( const G4VPreCompoundFragment & right)
delete

Member Function Documentation

◆ CalcEmissionProbability()

G4double G4VPreCompoundFragment::CalcEmissionProbability ( const G4Fragment & fr)
virtual

Reimplemented in G4HETCFragment.

Definition at line 139 of file G4VPreCompoundFragment.cc.

140{
141 G4bool ok = Initialize(fr);
143 if (ok) {
145 }
147}
bool G4bool
Definition G4Types.hh:86
G4bool Initialize(const G4Fragment &aFragment)
G4double ComputeIntegral(const G4double emin, const G4double emax)

◆ GetA()

◆ GetAlpha()

◆ GetBeta()

virtual G4double G4VPreCompoundFragment::GetBeta ( ) const
inlineprotectedvirtual

Reimplemented in G4HETCNeutron, and G4PreCompoundNeutron.

Definition at line 129 of file G4VPreCompoundFragment.hh.

Referenced by G4PreCompoundFragment::CrossSection().

◆ GetBindingEnergy()

G4double G4VPreCompoundFragment::GetBindingEnergy ( ) const
inline

Definition at line 88 of file G4VPreCompoundFragment.hh.

◆ GetEmissionProbability()

G4double G4VPreCompoundFragment::GetEmissionProbability ( ) const
inline

Definition at line 95 of file G4VPreCompoundFragment.hh.

95{ return theEmissionProbability; }

◆ GetEnergyThreshold()

G4double G4VPreCompoundFragment::GetEnergyThreshold ( ) const
inline

Definition at line 90 of file G4VPreCompoundFragment.hh.

91 {
93 }

Referenced by G4HETCFragment::CalcEmissionProbability().

◆ GetMomentum()

const G4LorentzVector & G4VPreCompoundFragment::GetMomentum ( ) const
inline

Definition at line 101 of file G4VPreCompoundFragment.hh.

101{ return theMomentum; }

Referenced by GetReactionProduct().

◆ GetNuclearMass()

G4double G4VPreCompoundFragment::GetNuclearMass ( ) const
inline

◆ GetReactionProduct()

G4ReactionProduct * G4VPreCompoundFragment::GetReactionProduct ( ) const
inline

Definition at line 168 of file G4VPreCompoundFragment.hh.

169{
170 G4ReactionProduct* theReactionProduct = new G4ReactionProduct(particle);
171 theReactionProduct->SetMomentum(GetMomentum().vect());
172 theReactionProduct->SetTotalEnergy(GetMomentum().e());
173 return theReactionProduct;
174}
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetTotalEnergy(const G4double en)
const G4LorentzVector & GetMomentum() const

Referenced by G4PreCompoundEmission::PerformEmission(), and G4PreCompoundEmissionInt::PerformEmission().

◆ GetRestA()

G4int G4VPreCompoundFragment::GetRestA ( ) const
inline

◆ GetRestNuclearMass()

G4double G4VPreCompoundFragment::GetRestNuclearMass ( ) const
inline

Definition at line 99 of file G4VPreCompoundFragment.hh.

◆ GetRestZ()

G4int G4VPreCompoundFragment::GetRestZ ( ) const
inline

◆ GetZ()

G4int G4VPreCompoundFragment::GetZ ( ) const
inline

◆ Initialize()

G4bool G4VPreCompoundFragment::Initialize ( const G4Fragment & aFragment)

Definition at line 91 of file G4VPreCompoundFragment.cc.

92{
93 theFragA = aFragment.GetA_asInt();
94 theFragZ = aFragment.GetZ_asInt();
97
99 if ((theResA < theResZ) || (theResA < theA) || (theResZ < theZ)
100 || (theResA == theA && theResZ < theZ)
101 || ((theResA > 1) && (theResA == theResZ || theResZ == 0))) {
102 return false;
103 }
104 pFragment = &aFragment;
106 G4double Ecm = aFragment.GetMomentum().m();
107 if (Ecm <= theResMass + theMass) { return 0.0; }
108
109 theResA13 = g4calc->Z13(theResA);
110 G4int nex = aFragment.GetNumberOfExcitons();
111
112 G4double elim = 0.0;
113 if (0 < theZ) {
114 theCoulombBarrier = theCoulombBarrierPtr->
115 GetCoulombBarrier(theResA + nex, theResZ, aFragment.GetExcitationEnergy());
116 elim = (0 < OPTxs) ? theCoulombBarrier*0.5 : theCoulombBarrier;
117 }
118
119 // Compute Maximal Kinetic Energy which can be carried by fragments
120 // after separation - the true assimptotic value
122 0.5*((Ecm - theResMass)*(Ecm + theResMass) + theMass*theMass)/Ecm - theMass;
123 G4double resM = Ecm - theMass - elim;
124 if (resM < theResMass) { return false; }
126 0.5*((Ecm - resM)*(Ecm + resM) + theMass*theMass)/Ecm - theMass;
127 theMinKinEnergy = std::max(theMinKinEnergy, 0.0);
128
129 if (theMinKinEnergy >= theMaxKinEnergy) { return false; }
130 // Calculate masses
132
133 // Compute Binding Energies for fragments
134 // needed to separate a fragment from the nucleus
136 return true;
137}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4double GetGroundStateMass() const
G4double GetExcitationEnergy() const
const G4LorentzVector & GetMomentum() const
G4int GetZ_asInt() const
G4int GetNumberOfExcitons() const
G4int GetA_asInt() const
static G4double GetNuclearMass(const G4double A, const G4double Z)

Referenced by CalcEmissionProbability().

◆ operator!=()

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

◆ operator=()

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

◆ operator==()

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

◆ ProbabilityDensityFunction()

G4double G4VPreCompoundFragment::ProbabilityDensityFunction ( G4double energy)
overridevirtual

Implements G4VSIntegration.

Definition at line 155 of file G4VPreCompoundFragment.cc.

156{
157 G4double e = std::max(ekin, 0.02);
159}
virtual G4double ProbabilityDistributionFunction(G4double, const G4Fragment &)

◆ ProbabilityDistributionFunction()

virtual G4double G4VPreCompoundFragment::ProbabilityDistributionFunction ( G4double ,
const G4Fragment &  )
inlineprotectedvirtual

Reimplemented in G4PreCompoundIon, and G4PreCompoundNucleon.

Definition at line 124 of file G4VPreCompoundFragment.hh.

125 { return 0.0; };

Referenced by ProbabilityDensityFunction().

◆ SampleKineticEnergy()

G4double G4VPreCompoundFragment::SampleKineticEnergy ( const G4Fragment & )
virtual

Reimplemented in G4HETCChargedFragment, and G4HETCNeutron.

Definition at line 149 of file G4VPreCompoundFragment.cc.

150{
151 G4double ekin = SampleValue();
152 return ekin;
153}

Referenced by G4PreCompoundEmission::PerformEmission(), and G4PreCompoundEmissionInt::PerformEmission().

◆ SetMomentum()

void G4VPreCompoundFragment::SetMomentum ( const G4LorentzVector & lv)
inline

Definition at line 103 of file G4VPreCompoundFragment.hh.

103{ theMomentum = lv; }

Referenced by G4PreCompoundEmission::PerformEmission(), and G4PreCompoundEmissionInt::PerformEmission().

◆ SetOPTxs()

void G4VPreCompoundFragment::SetOPTxs ( G4int )
inline

Definition at line 106 of file G4VPreCompoundFragment.hh.

106{}

◆ UseSICB()

void G4VPreCompoundFragment::UseSICB ( G4bool use)
inline

Definition at line 108 of file G4VPreCompoundFragment.hh.

◆ operator<< [1/2]

std::ostream & operator<< ( std::ostream & out,
const G4VPreCompoundFragment & theFragment )
friend

Definition at line 73 of file G4VPreCompoundFragment.cc.

75{
76 out << &theFragment;
77 return out;
78}

◆ operator<< [2/2]

std::ostream & operator<< ( std::ostream & out,
const G4VPreCompoundFragment * theFragment )
friend

Definition at line 80 of file G4VPreCompoundFragment.cc.

82{
83 out
84 << "PreCompoundModel Emitted Fragment: Z= " << theFragment->GetZ()
85 << " A= " << theFragment->GetA()
86 << " Mass(GeV)= " << theFragment->GetNuclearMass()/CLHEP::GeV;
87 return out;
88}

Member Data Documentation

◆ fNucData

◆ fXSection

G4InterfaceToXS* G4VPreCompoundFragment::fXSection {nullptr}
protected

◆ g4calc

◆ index

G4int G4VPreCompoundFragment::index {0}
protected

Definition at line 145 of file G4VPreCompoundFragment.hh.

145{0};

Referenced by G4PreCompoundFragment::CrossSection(), and G4VPreCompoundFragment().

◆ OPTxs

G4int G4VPreCompoundFragment::OPTxs
protected

◆ pFragment

const G4Fragment* G4VPreCompoundFragment::pFragment {nullptr}
protected

◆ theA

◆ theBindingEnergy

G4double G4VPreCompoundFragment::theBindingEnergy {0.0}
protected

◆ theCoulombBarrier

◆ theEmissionProbability

G4double G4VPreCompoundFragment::theEmissionProbability {0.0}
protected

◆ theFragA

◆ theFragZ

◆ theMass

G4double G4VPreCompoundFragment::theMass
protected

Definition at line 153 of file G4VPreCompoundFragment.hh.

Referenced by G4VPreCompoundFragment(), GetNuclearMass(), and Initialize().

◆ theMaxKinEnergy

◆ theMinKinEnergy

G4double G4VPreCompoundFragment::theMinKinEnergy {0.0}
protected

Definition at line 149 of file G4VPreCompoundFragment.hh.

149{0.0};

Referenced by CalcEmissionProbability(), and Initialize().

◆ theParameters

G4DeexPrecoParameters* G4VPreCompoundFragment::theParameters
protected

◆ theReducedMass

G4double G4VPreCompoundFragment::theReducedMass {0.0}
protected

◆ theResA

◆ theResA13

◆ theResMass

G4double G4VPreCompoundFragment::theResMass {0.0}
protected

Definition at line 151 of file G4VPreCompoundFragment.hh.

151{0.0};

Referenced by GetRestNuclearMass(), and Initialize().

◆ theResZ

◆ theZ

◆ useSICB

G4bool G4VPreCompoundFragment::useSICB {true}
protected

Definition at line 159 of file G4VPreCompoundFragment.hh.

159{true};

Referenced by UseSICB().


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