Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4E1SingleProbability1.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: G4E1SingleProbability1.cc 68724 2013-04-05 09:26:32Z gcosmo $
27 //
28 // Class G4E1SingleProbability1.cc
29 //
30 
33 #include "Randomize.hh"
34 #include "G4Pow.hh"
35 #include "G4PhysicalConstants.hh"
36 #include "G4SystemOfUnits.hh"
37 
38 // Constructors and operators
39 //
40 
42 {}
43 
45 {}
46 
47 // Calculate the emission probability
48 //
49 
51  G4double exciteE)
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  static 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 }
120 
122  G4double exciteE)
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  G4double Uexcite = frag.GetExcitationEnergy();
135  G4double Uafter = Uexcite - exciteE;
136 
137  static const G4double normC = 3.0;
138 
139  G4double upperLim = Uexcite;
140  G4double lowerLim = Uafter;
141  static 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 }
154 
155 G4double G4E1SingleProbability1::EmissionIntegration(const G4Fragment& frag,
156  G4double ,
157  G4double lowLim, G4double upLim,
158  G4int numIters)
159 
160 {
161 
162  // Simple Gaussian quadrature integration
163 
164  G4double x;
165  static const G4double root3 = 1.0/std::sqrt(3.0);
166 
167  G4double Step = (upLim-lowLim)/(2.0*numIters);
168  G4double Delta = Step*root3;
169 
170  G4double mean = 0.0;
171 
172  G4double theInt = 0.0;
173 
174  for(G4int i = 0; i < numIters; i++) {
175 
176  x = (2*i + 1)/Step;
177  G4double E1ProbDensityA = EmissionProbDensity(frag,x+Delta);
178  G4double E1ProbDensityB = EmissionProbDensity(frag,x-Delta);
179 
180  mean += E1ProbDensityA + E1ProbDensityB;
181 
182  }
183 
184  if(mean*Step > 0.0) theInt = mean*Step;
185 
186  return theInt;
187 
188 }
189 
190 
191 
static G4Pow * GetInstance()
Definition: G4Pow.cc:53
int G4int
Definition: G4Types.hh:78
int millibarn
Definition: hepunit.py:40
G4int GetA_asInt() const
Definition: G4Fragment.hh:238
G4double EmissionProbability(const G4Fragment &frag, G4double excite)
G4double EmissionProbDensity(const G4Fragment &frag, G4double ePhoton)
G4int GetZ_asInt() const
Definition: G4Fragment.hh:243
Definition: Step.hh:41
G4double powZ(G4int Z, G4double y) const
Definition: G4Pow.hh:258
G4double LevelDensityParameter(const G4int A, const G4int, const G4double) const
double G4double
Definition: G4Types.hh:76
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:255