Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4E1SingleProbability1 Class Reference

#include <G4E1SingleProbability1.hh>

+ Inheritance diagram for G4E1SingleProbability1:

Public Member Functions

 G4E1SingleProbability1 ()
 
virtual ~G4E1SingleProbability1 ()
 
G4double EmissionProbability (const G4Fragment &frag, G4double excite)
 
G4double EmissionProbDensity (const G4Fragment &frag, G4double ePhoton)
 
- Public Member Functions inherited from G4VEmissionProbability
 G4VEmissionProbability ()
 
virtual ~G4VEmissionProbability ()
 
virtual G4double EmissionProbability (const G4Fragment &fragment, const G4double anEnergy)=0
 
void SetOPTxs (G4int opt)
 
void UseSICB (G4bool use)
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmissionProbability
G4int OPTxs
 
G4bool useSICB
 
G4PowfG4pow
 
G4PairingCorrectionfPairCorr
 
G4EvaporationLevelDensityParametertheEvapLDPptr
 

Detailed Description

Definition at line 38 of file G4E1SingleProbability1.hh.

Constructor & Destructor Documentation

◆ G4E1SingleProbability1()

G4E1SingleProbability1::G4E1SingleProbability1 ( )

Definition at line 41 of file G4E1SingleProbability1.cc.

42{}

◆ ~G4E1SingleProbability1()

G4E1SingleProbability1::~G4E1SingleProbability1 ( )
virtual

Definition at line 44 of file G4E1SingleProbability1.cc.

45{}

Member Function Documentation

◆ EmissionProbability()

G4double G4E1SingleProbability1::EmissionProbability ( const G4Fragment frag,
G4double  excite 
)
virtual

Implements G4VEmissionProbability.

Definition at line 121 of file G4E1SingleProbability1.cc.

123{
124
125 // From nuclear fragment properties and the excitation energy, calculate
126 // the probability for photon evaporation down to the level
127 // Uexcite-exciteE.
128 // fragment = nuclear fragment BEFORE de-excitation
129
130 G4double theProb = 0.0;
131
132 G4double ScaleFactor = 1.0; // playing with scale factors
133
134 const G4double Uexcite = frag.GetExcitationEnergy();
135 G4double Uafter = Uexcite - exciteE;
136
137 G4double normC = 3.0;
138
139 const G4double upperLim = Uexcite;
140 const G4double lowerLim = Uafter;
141 const G4int numIters = 25;
142
143 // Need to integrate EmissionProbDensity from lowerLim to upperLim
144 // and multiply by normC
145
146 G4double integ = normC *
147 EmissionIntegration(frag,exciteE,lowerLim,upperLim,numIters);
148
149 if(integ > 0.0) theProb = integ;
150
151 return theProb * ScaleFactor;
152
153}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:235

◆ EmissionProbDensity()

G4double G4E1SingleProbability1::EmissionProbDensity ( const G4Fragment frag,
G4double  ePhoton 
)

Definition at line 50 of file G4E1SingleProbability1.cc.

52{
53
54 // Calculate the probability density here
55
56 // From nuclear fragment properties and the excitation energy, calculate
57 // the probability density for photon evaporation from U to U - exciteE
58 // (U = nucleus excitation energy, exciteE = total evaporated photon
59 // energy).
60 // fragment = nuclear fragment BEFORE de-excitation
61
62 G4double theProb = 0.0;
63
64 G4int Afrag = frag.GetA_asInt();
65 G4int Zfrag = frag.GetZ_asInt();
66 G4double Uexcite = frag.GetExcitationEnergy();
67
68 if( (Uexcite-exciteE) < 0.0 || exciteE < 0 || Uexcite <= 0) return theProb;
69
70 // Need a level density parameter.
71 // For now, just use the constant approximation (not reliable near magic
72 // nuclei).
73
75 G4double aLevelDensityParam = a.LevelDensityParameter(Afrag,Zfrag,Uexcite);
76
77 G4double levelDensBef = std::exp(2.0*std::sqrt(aLevelDensityParam*Uexcite));
78 G4double levelDensAft = std::exp(2.0*std::sqrt(aLevelDensityParam*(Uexcite-exciteE)));
79
80 // Now form the probability density
81
82 // Define constants for the photoabsorption cross-section (the reverse
83 // process of our de-excitation)
84
85 G4double sigma0 = 2.5 * Afrag * millibarn; // millibarns
86
87 G4double Egdp = (40.3 / G4Pow::GetInstance()->powZ(Afrag,0.2) )*MeV;
88 G4double GammaR = 0.30 * Egdp;
89
90 const G4double normC = 1.0 / ((pi * hbarc)*(pi * hbarc));
91
92 // CD
93 //cout<<" PROB TESTS "<<G4endl;
94 //cout<<" hbarc = "<<hbarc<<G4endl;
95 //cout<<" pi = "<<pi<<G4endl;
96 //cout<<" Uexcite, exciteE = "<<Uexcite<<" "<<exciteE<<G4endl;
97 //cout<<" Uexcite, exciteE = "<<Uexcite*MeV<<" "<<exciteE*MeV<<G4endl;
98 //cout<<" lev density param = "<<aLevelDensityParam<<G4endl;
99 //cout<<" level densities = "<<levelDensBef<<" "<<levelDensAft<<G4endl;
100 //cout<<" sigma0 = "<<sigma0<<G4endl;
101 //cout<<" Egdp, GammaR = "<<Egdp<<" "<<GammaR<<G4endl;
102 //cout<<" normC = "<<normC<<G4endl;
103
104 G4double numerator = sigma0 * exciteE*exciteE * GammaR*GammaR;
105 G4double denominator = (exciteE*exciteE - Egdp*Egdp)*
106 (exciteE*exciteE - Egdp*Egdp) + GammaR*GammaR*exciteE*exciteE;
107
108 G4double sigmaAbs = numerator/denominator;
109
110 theProb = normC * sigmaAbs * exciteE*exciteE *
111 levelDensAft/levelDensBef;
112
113 // CD
114 //cout<<" sigmaAbs = "<<sigmaAbs<<G4endl;
115 //cout<<" Probability = "<<theProb<<G4endl;
116
117 return theProb;
118
119}
G4double LevelDensityParameter(const G4int A, const G4int, const G4double) const
G4int GetZ_asInt() const
Definition: G4Fragment.hh:223
G4int GetA_asInt() const
Definition: G4Fragment.hh:218
static G4Pow * GetInstance()
Definition: G4Pow.cc:50
G4double powZ(G4int Z, G4double y)
Definition: G4Pow.hh:180
const G4double pi

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