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

#include <G4PreCompoundFragment.hh>

+ Inheritance diagram for G4PreCompoundFragment:

Public Member Functions

 G4PreCompoundFragment (const G4ParticleDefinition *, G4VCoulombBarrier *aCoulombBarrier)
 
 ~G4PreCompoundFragment () override=default
 
G4double CalcEmissionProbability (const G4Fragment &aFragment) override
 
G4double SampleKineticEnergy (const G4Fragment &aFragment) override
 
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 *)
 
virtual ~G4VPreCompoundFragment ()
 
G4bool Initialize (const G4Fragment &aFragment)
 
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 ekin, const G4Fragment &aFragment)=0
 
- Protected Member Functions inherited from G4VPreCompoundFragment
virtual G4double GetAlpha () const =0
 
virtual G4double GetBeta () const
 

Additional Inherited Members

- Protected Attributes inherited from G4VPreCompoundFragment
G4NuclearLevelDatafNucData
 
G4DeexPrecoParameterstheParameters
 
G4Powg4calc
 
G4InterfaceToXSfXSection {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 44 of file G4PreCompoundFragment.cc.

46 : G4VPreCompoundFragment(p, aCoulBarrier)
47{}
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

◆ CalcEmissionProbability()

G4double G4PreCompoundFragment::CalcEmissionProbability ( const G4Fragment & aFragment)
overridevirtual

Implements G4VPreCompoundFragment.

Definition at line 49 of file G4PreCompoundFragment.cc.

50{
52 IntegrateEmissionProbability(theMinKinEnergy, theMaxKinEnergy, fr) : 0.0;
53 /*
54 G4cout << "## G4PreCompoundFragment::CalcEmisProb "
55 << "Zf= " << fr.GetZ_asInt()
56 << " Af= " << fr.GetA_asInt()
57 << " Elow= " << theMinKinEnergy
58 << " Eup= " << theMaxKinEnergy
59 << " prob= " << theEmissionProbability
60 << " index=" << index << " Z=" << theZ << " A=" << theA
61 << G4endl;
62 */
64}
G4bool Initialize(const G4Fragment &aFragment)

◆ CrossSection()

G4double G4PreCompoundFragment::CrossSection ( G4double ekin)

Definition at line 94 of file G4PreCompoundFragment.cc.

95{
96 /*
97 G4cout << "G4PreCompoundFragment::CrossSection OPTxs=" << OPTxs << " E=" << ekin
98 << " resZ=" << theResZ << " resA=" << theResA << " index=" << index
99 << " fXSection:" << fXSection << G4endl;
100 */
101 // compute power once
102 if (OPTxs > 1 && 0 < index && theResA != lastA) {
103 lastA = theResA;
105 }
106 if (OPTxs == 0) {
107 recentXS = GetOpt0(ekin);
108 } else if (OPTxs == 1) {
109 G4int Z = std::min(theResZ, ZMAXNUCLEARDATA);
110 //G4double e = std::max(ekin, lowEnergyLimitMeV[Z]);
111 recentXS = fXSection->GetElementCrossSection(ekin, Z)/CLHEP::millibarn;
112
113 } else if (OPTxs == 2) {
116 theResA13, muu,
117 index, theZ, theResA);
118
119 } else {
121 theResA13, muu, index,
122 theZ, theA, theResA);
123 }
124 return recentXS;
125}
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)

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

◆ ProbabilityDistributionFunction()

virtual G4double G4PreCompoundFragment::ProbabilityDistributionFunction ( G4double ekin,
const G4Fragment & aFragment )
protectedpure virtual

Implemented in G4PreCompoundIon, and G4PreCompoundNucleon.

Referenced by SampleKineticEnergy().

◆ RecentXS()

G4double G4PreCompoundFragment::RecentXS ( ) const
inline

Definition at line 64 of file G4PreCompoundFragment.hh.

64{ return recentXS; };

◆ SampleKineticEnergy()

G4double G4PreCompoundFragment::SampleKineticEnergy ( const G4Fragment & aFragment)
overridevirtual

Implements G4VPreCompoundFragment.

Definition at line 136 of file G4PreCompoundFragment.cc.

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

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