Geant4  10.01.p03
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: G4E1Probability.cc 86783 2014-11-18 08:43:58Z gcosmo $
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 "G4Exp.hh"
46 #include "G4SystemOfUnits.hh"
47 
48 // Constructors and operators
49 //
50 
52 {
53  G4double x = CLHEP::pi*CLHEP::hbarc;
54  normC = 1.0 / (x*x);
57 }
58 
60 {}
61 
62 
64  G4double gammaE)
65 {
66  // Calculate the probability density here
67 
68  // From nuclear fragment properties and the excitation energy, calculate
69  // the probability density for photon evaporation from U to U - gammaE
70  // (U = nucleus excitation energy, gammaE = total evaporated photon
71  // energy). Fragment = nuclear fragment BEFORE de-excitation
72 
73  G4double theProb = 0.0;
74 
75  G4int Afrag = frag.GetA_asInt();
76  G4double Uexcite = frag.GetExcitationEnergy();
77  G4double U = Uexcite - gammaE;
78 
79  if(U < 0.0) { return theProb; }
80 
81  // Need a level density parameter.
82  // For now, just use the constant approximation (not reliable near magic
83  // nuclei) - is equivalent to G4ConstantLevelDensityParameter class
84 
85  G4double aLevelDensityParam = Afrag*theLevelDensityParameter;
86 
87  // VI reduce number of calls to exp
88  G4double levelDens =
89  G4Exp(2*(std::sqrt(aLevelDensityParam*U)-std::sqrt(aLevelDensityParam*Uexcite)));
90  // Now form the probability density
91 
92  // Define constants for the photoabsorption cross-section (the reverse
93  // process of our de-excitation)
94 
95  G4double sigma0 = 2.5 * Afrag * millibarn; // millibarns
96 
97  G4double Egdp = (40.3 / fG4pow->powZ(Afrag,0.2) )*MeV;
98  G4double GammaR = 0.30 * Egdp;
99 
100  // CD
101  //cout<<" PROB TESTS "<<G4endl;
102  //cout<<" hbarc = "<<hbarc<<G4endl;
103  //cout<<" pi = "<<pi<<G4endl;
104  //cout<<" Uexcite, gammaE = "<<Uexcite<<" "<<gammaE<<G4endl;
105  //cout<<" Uexcite, gammaE = "<<Uexcite*MeV<<" "<<gammaE*MeV<<G4endl;
106  //cout<<" lev density param = "<<aLevelDensityParam<<G4endl;
107  //cout<<" level densities = "<<levelDensBef<<" "<<levelDensAft<<G4endl;
108  //cout<<" sigma0 = "<<sigma0<<G4endl;
109  //cout<<" Egdp, GammaR = "<<Egdp<<" "<<GammaR<<G4endl;
110  //cout<<" normC = "<<normC<<G4endl;
111 
112  // VI implementation 18.05.2010
113  G4double gammaE2 = gammaE*gammaE;
114  G4double gammaR2 = gammaE2*GammaR*GammaR;
115  G4double egdp2 = gammaE2 - Egdp*Egdp;
116  G4double sigmaAbs = sigma0*gammaR2/(egdp2*egdp2 + gammaR2);
117  theProb = normC * sigmaAbs * gammaE2 * levelDens;
118 
119  // CD
120  //cout<<" sigmaAbs = "<<sigmaAbs<<G4endl;
121  //cout<<" Probability = "<<theProb<<G4endl;
122 
123  return theProb;
124 
125 }
126 
128  G4double gammaE)
129 {
130  // From nuclear fragment properties and the excitation energy, calculate
131  // the probability for photon evaporation down to last ground level.
132  // fragment = nuclear fragment BEFORE de-excitation
133 
134  G4double upperLim = gammaE;
135  G4double lowerLim = 0.0;
136 
137  //G4cout << "G4E1Probability::EmissionProbability: Emin= " << lowerLim
138  // << " Emax= " << upperLim << G4endl;
139  if( upperLim - lowerLim <= CLHEP::keV ) { return 0.0; }
140 
141  // Need to integrate EmissionProbDensity from lowerLim to upperLim
142  // and multiply by factor 3 (?!)
143 
144  G4double integ = EmissionIntegration(frag,lowerLim,upperLim);
145 
146  return integ;
147 
148 }
149 
151  G4double lowLim, G4double upLim)
152 
153 {
154  // Simple integration
155  // VI replace by direct integration over 100 point
156 
157  static const G4int numIters = 100;
158  G4double Step = (upLim-lowLim)/G4double(numIters);
159 
160  G4double res = 0.0;
161  G4double x = lowLim - 0.5*Step;
162 
163  for(G4int i = 0; i < numIters; ++i) {
164  x += Step;
165  res += EmissionProbDensity(frag, x);
166  }
167 
168  if(res > 0.0) { res *= Step; }
169  else { res = 0.0; }
170 
171  return res;
172 
173 }
174 
175 
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double theLevelDensityParameter
static const double MeV
Definition: G4SIunits.hh:193
const G4double pi
G4double EmissionProbability(const G4Fragment &frag, G4double excite)
G4double EmissionProbDensity(const G4Fragment &frag, G4double ePhoton)
int G4int
Definition: G4Types.hh:78
G4int GetA_asInt() const
Definition: G4Fragment.hh:243
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
virtual ~G4E1Probability()
G4double EmissionIntegration(const G4Fragment &frag, G4double lowLim, G4double upLim)
static const double millibarn
Definition: G4SIunits.hh:96
Definition: Step.hh:41
static const double keV
Definition: G4SIunits.hh:195
G4double powZ(G4int Z, G4double y) const
Definition: G4Pow.hh:256
double G4double
Definition: G4Types.hh:76
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:260