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

#include <G4PreCompoundIon.hh>

+ Inheritance diagram for G4PreCompoundIon:

Public Member Functions

 G4PreCompoundIon (const G4ParticleDefinition *, G4VCoulombBarrier *aCoulombBarrier)
 
 ~G4PreCompoundIon () override=default
 
 G4PreCompoundIon (const G4PreCompoundIon &right)=delete
 
const G4PreCompoundIonoperator= (const G4PreCompoundIon &right)=delete
 
G4bool operator== (const G4PreCompoundIon &right) const =delete
 
G4bool operator!= (const G4PreCompoundIon &right) const =delete
 
- Public Member Functions inherited from G4PreCompoundFragment
 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 ProbabilityDistributionFunction (G4double eKin, const G4Fragment &) override
 
virtual G4double GetRj (G4int NumberParticles, G4int NumberCharged) const =0
 
virtual G4double FactorialFactor (G4int N, G4int P) const =0
 
virtual G4double CoalescenceFactor (G4int A) const =0
 
- Protected Member Functions inherited from G4PreCompoundFragment
G4double CrossSection (G4double ekin) const
 
- 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 41 of file G4PreCompoundIon.hh.

Constructor & Destructor Documentation

◆ G4PreCompoundIon() [1/2]

G4PreCompoundIon::G4PreCompoundIon ( const G4ParticleDefinition * part,
G4VCoulombBarrier * aCoulombBarrier )

Definition at line 47 of file G4PreCompoundIon.cc.

50 : G4PreCompoundFragment(part,aCoulombBarrier)
51{
53 fact = 0.75*CLHEP::millibarn/(CLHEP::pi*r0*r0*r0);
54}
double G4double
Definition G4Types.hh:83
G4PreCompoundFragment(const G4ParticleDefinition *, G4VCoulombBarrier *aCoulombBarrier)
G4DeexPrecoParameters * theParameters

◆ ~G4PreCompoundIon()

G4PreCompoundIon::~G4PreCompoundIon ( )
overridedefault

◆ G4PreCompoundIon() [2/2]

G4PreCompoundIon::G4PreCompoundIon ( const G4PreCompoundIon & right)
delete

Member Function Documentation

◆ CoalescenceFactor()

virtual G4double G4PreCompoundIon::CoalescenceFactor ( G4int A) const
protectedpure virtual

◆ FactorialFactor()

virtual G4double G4PreCompoundIon::FactorialFactor ( G4int N,
G4int P ) const
protectedpure virtual

◆ GetRj()

virtual G4double G4PreCompoundIon::GetRj ( G4int NumberParticles,
G4int NumberCharged ) const
protectedpure virtual

◆ operator!=()

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

◆ operator=()

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

◆ operator==()

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

◆ ProbabilityDistributionFunction()

G4double G4PreCompoundIon::ProbabilityDistributionFunction ( G4double eKin,
const G4Fragment & aFragment )
overrideprotectedvirtual

Implements G4PreCompoundFragment.

Definition at line 56 of file G4PreCompoundIon.cc.

59{
60 G4double efinal = eKin + theBindingEnergy;
61 if(efinal <= 0.0 ) { return 0.0; }
62
63 G4double U = aFragment.GetExcitationEnergy();
64 G4int P = aFragment.GetNumberOfParticles();
65 G4int H = aFragment.GetNumberOfHoles();
66 G4int A = GetA();
67 G4int N = P + H;
68
69 static const G4double sixoverpi2 = 6.0/CLHEP::pi2;
71 G4double g1 = sixoverpi2*fNucData->GetLevelDensity(theResZ, theResA, 0.0);
72
73 G4double gj = g1;
74
75 G4double A0 = (P*P+H*H+P-3*H)/(4.0*g0);
76 G4double A1 = std::max(0.0,(A0*g0 + A*(A-2*P-1)*0.25)/g1);
77
78 G4double E0 = U - A0;
79 if (E0 <= 0.0) { return 0.0; }
80
81 G4double E1 = std::max(0.0,theMaxKinEnergy - eKin - A1);
82
83 G4double Aj = A*(A+1)/(4.0*gj);
84 G4double Ej = std::max(0.0,efinal - Aj);
85
86 G4double rj = GetRj(P, aFragment.GetNumberOfCharged());
87 G4double xs = CrossSection(eKin);
88
89 G4double pA = fact*eKin*xs*rj
91 * std::sqrt(2.0/(theReducedMass*efinal))
92 * g4calc->powN(g1*E1/(g0*E0), N-A-1)
93 * g4calc->powN(gj*Ej/(g0*E0), A-1)*gj*g1/(g0*g0*E0*theResA);
94
95 return pA;
96}
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4int GetNumberOfParticles() const
G4int GetNumberOfHoles() const
G4double GetExcitationEnergy() const
G4int GetNumberOfCharged() const
G4double GetLevelDensity(G4int Z, G4int A, G4double U)
G4double powN(G4double x, G4int n) const
Definition G4Pow.cc:162
G4double CrossSection(G4double ekin) const
virtual G4double FactorialFactor(G4int N, G4int P) const =0
virtual G4double CoalescenceFactor(G4int A) const =0
virtual G4double GetRj(G4int NumberParticles, G4int NumberCharged) const =0
#define N
Definition crc32.c:57

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