Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4E1Probability.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id$
27//
28//---------------------------------------------------------------------
29//
30// Geant4 class G4E1Probability
31//
32// by V. Lara (May 2003)
33//
34// Modifications:
35// 18.05.2010 V.Ivanchenko trying to speedup the most slow method
36// by usage of G4Pow, integer A and introduction of const members
37// 17.11.2010 V.Ivanchenko perform general cleanup and simplification
38// of integration method; low-limit of integration is defined
39// by gamma energy or is zero (was always zero before)
40//
41
42#include "G4E1Probability.hh"
43#include "Randomize.hh"
44#include "G4Pow.hh"
45#include "G4SystemOfUnits.hh"
46
47// Constructors and operators
48//
49
51{
52 G4double x = CLHEP::pi*CLHEP::hbarc;
53 normC = 1.0 / (x*x);
54 theLevelDensityParameter = 0.125/MeV;
55 fG4pow = G4Pow::GetInstance();
56}
57
59{}
60
61// Calculate the emission probability
62//
63
65 G4double gammaE)
66{
67
68 // Calculate the probability density here
69
70 // From nuclear fragment properties and the excitation energy, calculate
71 // the probability density for photon evaporation from U to U - gammaE
72 // (U = nucleus excitation energy, gammaE = total evaporated photon
73 // energy). Fragment = nuclear fragment BEFORE de-excitation
74
75 G4double theProb = 0.0;
76
77 G4int Afrag = frag.GetA_asInt();
78 G4double Uexcite = frag.GetExcitationEnergy();
79 G4double U = std::max(0.0,Uexcite-gammaE);
80
81 if(gammaE < 0.0) { return theProb; }
82
83 // Need a level density parameter.
84 // For now, just use the constant approximation (not reliable near magic
85 // nuclei).
86
87 G4double aLevelDensityParam = Afrag*theLevelDensityParameter;
88
89 // G4double levelDensBef = std::exp(2*std::sqrt(aLevelDensityParam*Uexcite));
90 // G4double levelDensAft = std::exp(2*std::sqrt(aLevelDensityParam*(Uexcite-gammaE)));
91 // VI reduce number of calls to exp
92 G4double levelDens =
93 std::exp(2*(std::sqrt(aLevelDensityParam*U)-std::sqrt(aLevelDensityParam*Uexcite)));
94 // Now form the probability density
95
96 // Define constants for the photoabsorption cross-section (the reverse
97 // process of our de-excitation)
98
99 G4double sigma0 = 2.5 * Afrag * millibarn; // millibarns
100
101 G4double Egdp = (40.3 / fG4pow->powZ(Afrag,0.2) )*MeV;
102 G4double GammaR = 0.30 * Egdp;
103
104 // CD
105 //cout<<" PROB TESTS "<<G4endl;
106 //cout<<" hbarc = "<<hbarc<<G4endl;
107 //cout<<" pi = "<<pi<<G4endl;
108 //cout<<" Uexcite, gammaE = "<<Uexcite<<" "<<gammaE<<G4endl;
109 //cout<<" Uexcite, gammaE = "<<Uexcite*MeV<<" "<<gammaE*MeV<<G4endl;
110 //cout<<" lev density param = "<<aLevelDensityParam<<G4endl;
111 //cout<<" level densities = "<<levelDensBef<<" "<<levelDensAft<<G4endl;
112 //cout<<" sigma0 = "<<sigma0<<G4endl;
113 //cout<<" Egdp, GammaR = "<<Egdp<<" "<<GammaR<<G4endl;
114 //cout<<" normC = "<<normC<<G4endl;
115
116 // VI implementation 18.05.2010
117 G4double gammaE2 = gammaE*gammaE;
118 G4double gammaR2 = gammaE2*GammaR*GammaR;
119 G4double egdp2 = gammaE2 - Egdp*Egdp;
120 G4double sigmaAbs = sigma0*gammaR2/(egdp2*egdp2 + gammaR2);
121 theProb = normC * sigmaAbs * gammaE2 * levelDens;
122
123 // old implementation
124 // G4double numerator = sigma0 * gammaE*gammaE * GammaR*GammaR;
125 // G4double denominator = (gammaE*gammaE - Egdp*Egdp)*
126 // (gammaE*gammaE - Egdp*Egdp) + GammaR*GammaR*gammaE*gammaE;
127
128 //G4double sigmaAbs = numerator/denominator;
129 //theProb = normC * sigmaAbs * gammaE2 * levelDensAft/levelDensBef;
130
131 // CD
132 //cout<<" sigmaAbs = "<<sigmaAbs<<G4endl;
133 //cout<<" Probability = "<<theProb<<G4endl;
134
135 return theProb;
136
137}
138
140 G4double gammaE)
141{
142 // From nuclear fragment properties and the excitation energy, calculate
143 // the probability for photon evaporation down to last ground level.
144 // fragment = nuclear fragment BEFORE de-excitation
145
146 G4double upperLim = gammaE;
147 G4double lowerLim = 0.0;
148
149 //G4cout << "G4E1Probability::EmissionProbability: Emin= " << lowerLim
150 // << " Emax= " << upperLim << G4endl;
151 if( upperLim - lowerLim <= CLHEP::keV ) { return 0.0; }
152
153 // Need to integrate EmissionProbDensity from lowerLim to upperLim
154 // and multiply by factor 3 (?!)
155
156 G4double integ = 3.0 * EmissionIntegration(frag,lowerLim,upperLim);
157
158 return integ;
159
160}
161
162G4double G4E1Probability::EmissionIntegration(const G4Fragment& frag,
163 G4double lowLim, G4double upLim)
164
165{
166 // Simple integration
167 // VI replace by direct integration over 100 point
168
169 const G4int numIters = 100;
170 G4double Step = (upLim-lowLim)/G4double(numIters);
171
172 G4double res = 0.0;
173 G4double x = lowLim - 0.5*Step;
174
175 for(G4int i = 0; i < numIters; ++i) {
176 x += Step;
177 res += EmissionProbDensity(frag, x);
178 }
179
180 if(res > 0.0) { res /= G4double(numIters); }
181 else { res = 0.0; }
182
183 return res;
184
185}
186
187
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4double EmissionProbDensity(const G4Fragment &frag, G4double ePhoton)
G4double EmissionProbability(const G4Fragment &frag, G4double excite)
virtual ~G4E1Probability()
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:235
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