Geant4 10.7.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)
 
virtual ~G4PreCompoundFragment ()
 
G4double CalcEmissionProbability (const G4Fragment &aFragment)
 
G4double SampleKineticEnergy (const G4Fragment &aFragment)
 
- Public Member Functions inherited from G4VPreCompoundFragment
 G4VPreCompoundFragment (const G4ParticleDefinition *, G4VCoulombBarrier *aCoulombBarrier)
 
virtual ~G4VPreCompoundFragment ()
 
void Initialize (const G4Fragment &aFragment)
 
virtual G4double CalcEmissionProbability (const G4Fragment &aFragment)=0
 
virtual G4double SampleKineticEnergy (const G4Fragment &aFragment)=0
 
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 &value)
 
void SetOPTxs (G4int)
 
void UseSICB (G4bool)
 

Protected Member Functions

virtual G4double GetAlpha () const =0
 
virtual G4double GetBeta () const =0
 
G4double CrossSection (G4double ekin) const
 
virtual G4double ProbabilityDistributionFunction (G4double K, const G4Fragment &aFragment)=0
 
virtual G4double GetAlpha () const =0
 
virtual G4double GetBeta () const =0
 

Additional Inherited Members

- Protected Attributes inherited from G4VPreCompoundFragment
G4NuclearLevelDatafNucData
 
G4DeexPrecoParameterstheParameters
 
G4Powg4calc
 
G4int theA
 
G4int theZ
 
G4int theResA
 
G4int theResZ
 
G4int theFragA
 
G4int theFragZ
 
G4double theResA13
 
G4double theBindingEnergy
 
G4double theMinKinEnergy
 
G4double theMaxKinEnergy
 
G4double theResMass
 
G4double theReducedMass
 
G4double theMass
 
G4double theEmissionProbability
 
G4double theCoulombBarrier
 
G4int OPTxs
 
G4bool useSICB
 

Detailed Description

Definition at line 43 of file G4PreCompoundFragment.hh.

Constructor & Destructor Documentation

◆ G4PreCompoundFragment()

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}

◆ ~G4PreCompoundFragment()

G4PreCompoundFragment::~G4PreCompoundFragment ( )
virtual

Definition at line 52 of file G4PreCompoundFragment.cc.

53{}

Member Function Documentation

◆ CalcEmissionProbability()

G4double G4PreCompoundFragment::CalcEmissionProbability ( const G4Fragment aFragment)
virtual

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
72 IntegrateEmissionProbability(theMinKinEnergy,theMaxKinEnergy,aFragment);
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)

◆ CrossSection()

G4double G4PreCompoundFragment::CrossSection ( G4double  ekin) const
protected

Definition at line 112 of file G4PreCompoundFragment.cc.

113{
114 G4double res;
115 if(OPTxs == 0 || (OPTxs == 4 && theMaxKinEnergy < 10.)) {
116 res = GetOpt0(ekin);
117
118 } else if(OPTxs <= 2) {
120 theResA13, muu,
121 index, theZ, theResA);
122
123 } else {
125 theResA13, muu,
126 index, theZ, theA, theResA);
127 }
128 return res;
129}
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().

◆ GetAlpha()

virtual G4double G4PreCompoundFragment::GetAlpha ( ) const
protectedpure virtual

◆ GetBeta()

virtual G4double G4PreCompoundFragment::GetBeta ( ) const
protectedpure virtual

◆ ProbabilityDistributionFunction()

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

Implemented in G4PreCompoundIon, and G4PreCompoundNucleon.

Referenced by SampleKineticEnergy().

◆ SampleKineticEnergy()

G4double G4PreCompoundFragment::SampleKineticEnergy ( const G4Fragment aFragment)
virtual

Implements G4VPreCompoundFragment.

Definition at line 140 of file G4PreCompoundFragment.cc.

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

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