Geant4 11.2.2
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
 
 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 ()
 
void Initialize (const G4Fragment &aFragment)
 
G4bool IsItPossible (const G4Fragment &aFragment) const
 
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 opt)
 
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

G4double CrossSection (G4double ekin) const
 
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
 
G4int theA
 
G4int theZ
 
G4int theResA {0}
 
G4int theResZ {0}
 
G4int theFragA {0}
 
G4int theFragZ {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}
 
G4int OPTxs {3}
 
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 42 of file G4PreCompoundFragment.cc.

44 : G4VPreCompoundFragment(p, aCoulBarrier)
45{
46 muu = probmax = 0.0;
47 if(0 == theZ) { index = 0; }
48 else if(1 == theZ) { index = theA; }
49 else { index = theA + 1; }
50}
G4VPreCompoundFragment(const G4ParticleDefinition *, G4VCoulombBarrier *)

◆ ~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 52 of file G4PreCompoundFragment.cc.

53{
54 //G4cout << theCoulombBarrier << " " << GetMaximalKineticEnergy() << G4endl;
55 // If theCoulombBarrier effect is included in the emission probabilities
56 // Coulomb barrier is the lower limit of integration over kinetic energy
57
59
60 if (theMaxKinEnergy <= theMinKinEnergy) { return 0.0; }
61
62 // compute power once
63 if(0 < index) {
65 }
66
68 IntegrateEmissionProbability(theMinKinEnergy, theMaxKinEnergy, fr);
69 /*
70 G4cout << "## G4PreCompoundFragment::CalcEmisProb "
71 << "Z= " << fr.GetZ_asInt()
72 << " A= " << fr.GetA_asInt()
73 << " Elow= " << LowerLimit/MeV
74 << " Eup= " << UpperLimit/MeV
75 << " prob= " << theEmissionProbability
76 << G4endl;
77 */
79}
static G4double ComputePowerParameter(G4int resA, G4int idx)

◆ CrossSection()

G4double G4PreCompoundFragment::CrossSection ( G4double ekin) const
protected

Definition at line 109 of file G4PreCompoundFragment.cc.

110{
111 G4double res;
112 if(OPTxs == 0 || (OPTxs == 4 && theMaxKinEnergy < 10.)) {
113 res = GetOpt0(ekin);
114
115 } else if(OPTxs <= 2) {
118 theResA13, muu,
119 index, theZ, theResA);
120
121 } else {
123 theResA13, muu, index,
124 theZ, theA, theResA);
125 }
126 return res;
127}
double G4double
Definition G4Types.hh:83
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)

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().

◆ SampleKineticEnergy()

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

Implements G4VPreCompoundFragment.

Definition at line 138 of file G4PreCompoundFragment.cc.

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

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