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,v 1.5 2010/11/17 16:50:53 vnivanch Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-04-ref-00 $ |
---|
28 | // |
---|
29 | // Class G4E1SingleProbability1.cc |
---|
30 | // |
---|
31 | |
---|
32 | #include "G4E1SingleProbability1.hh" |
---|
33 | #include "G4ConstantLevelDensityParameter.hh" |
---|
34 | #include "Randomize.hh" |
---|
35 | #include "G4Pow.hh" |
---|
36 | |
---|
37 | // Constructors and operators |
---|
38 | // |
---|
39 | |
---|
40 | G4E1SingleProbability1::G4E1SingleProbability1() |
---|
41 | {} |
---|
42 | |
---|
43 | G4E1SingleProbability1::~G4E1SingleProbability1() |
---|
44 | {} |
---|
45 | |
---|
46 | // Calculate the emission probability |
---|
47 | // |
---|
48 | |
---|
49 | G4double G4E1SingleProbability1::EmissionProbDensity(const G4Fragment& frag, |
---|
50 | G4double exciteE) |
---|
51 | { |
---|
52 | |
---|
53 | // Calculate the probability density here |
---|
54 | |
---|
55 | // From nuclear fragment properties and the excitation energy, calculate |
---|
56 | // the probability density for photon evaporation from U to U - exciteE |
---|
57 | // (U = nucleus excitation energy, exciteE = total evaporated photon |
---|
58 | // energy). |
---|
59 | // fragment = nuclear fragment BEFORE de-excitation |
---|
60 | |
---|
61 | G4double theProb = 0.0; |
---|
62 | |
---|
63 | G4int Afrag = frag.GetA_asInt(); |
---|
64 | G4int Zfrag = frag.GetZ_asInt(); |
---|
65 | G4double Uexcite = frag.GetExcitationEnergy(); |
---|
66 | |
---|
67 | if( (Uexcite-exciteE) < 0.0 || exciteE < 0 || Uexcite <= 0) return theProb; |
---|
68 | |
---|
69 | // Need a level density parameter. |
---|
70 | // For now, just use the constant approximation (not reliable near magic |
---|
71 | // nuclei). |
---|
72 | |
---|
73 | G4ConstantLevelDensityParameter a; |
---|
74 | G4double aLevelDensityParam = a.LevelDensityParameter(Afrag,Zfrag,Uexcite); |
---|
75 | |
---|
76 | G4double levelDensBef = std::exp(2.0*std::sqrt(aLevelDensityParam*Uexcite)); |
---|
77 | G4double levelDensAft = std::exp(2.0*std::sqrt(aLevelDensityParam*(Uexcite-exciteE))); |
---|
78 | |
---|
79 | // Now form the probability density |
---|
80 | |
---|
81 | // Define constants for the photoabsorption cross-section (the reverse |
---|
82 | // process of our de-excitation) |
---|
83 | |
---|
84 | G4double sigma0 = 2.5 * Afrag * millibarn; // millibarns |
---|
85 | |
---|
86 | G4double Egdp = (40.3 / G4Pow::GetInstance()->powZ(Afrag,0.2) )*MeV; |
---|
87 | G4double GammaR = 0.30 * Egdp; |
---|
88 | |
---|
89 | const G4double normC = 1.0 / ((pi * hbarc)*(pi * hbarc)); |
---|
90 | |
---|
91 | // CD |
---|
92 | //cout<<" PROB TESTS "<<G4endl; |
---|
93 | //cout<<" hbarc = "<<hbarc<<G4endl; |
---|
94 | //cout<<" pi = "<<pi<<G4endl; |
---|
95 | //cout<<" Uexcite, exciteE = "<<Uexcite<<" "<<exciteE<<G4endl; |
---|
96 | //cout<<" Uexcite, exciteE = "<<Uexcite*MeV<<" "<<exciteE*MeV<<G4endl; |
---|
97 | //cout<<" lev density param = "<<aLevelDensityParam<<G4endl; |
---|
98 | //cout<<" level densities = "<<levelDensBef<<" "<<levelDensAft<<G4endl; |
---|
99 | //cout<<" sigma0 = "<<sigma0<<G4endl; |
---|
100 | //cout<<" Egdp, GammaR = "<<Egdp<<" "<<GammaR<<G4endl; |
---|
101 | //cout<<" normC = "<<normC<<G4endl; |
---|
102 | |
---|
103 | G4double numerator = sigma0 * exciteE*exciteE * GammaR*GammaR; |
---|
104 | G4double denominator = (exciteE*exciteE - Egdp*Egdp)* |
---|
105 | (exciteE*exciteE - Egdp*Egdp) + GammaR*GammaR*exciteE*exciteE; |
---|
106 | |
---|
107 | G4double sigmaAbs = numerator/denominator; |
---|
108 | |
---|
109 | theProb = normC * sigmaAbs * exciteE*exciteE * |
---|
110 | levelDensAft/levelDensBef; |
---|
111 | |
---|
112 | // CD |
---|
113 | //cout<<" sigmaAbs = "<<sigmaAbs<<G4endl; |
---|
114 | //cout<<" Probability = "<<theProb<<G4endl; |
---|
115 | |
---|
116 | return theProb; |
---|
117 | |
---|
118 | } |
---|
119 | |
---|
120 | G4double G4E1SingleProbability1::EmissionProbability(const G4Fragment& frag, |
---|
121 | G4double exciteE) |
---|
122 | { |
---|
123 | |
---|
124 | // From nuclear fragment properties and the excitation energy, calculate |
---|
125 | // the probability for photon evaporation down to the level |
---|
126 | // Uexcite-exciteE. |
---|
127 | // fragment = nuclear fragment BEFORE de-excitation |
---|
128 | |
---|
129 | G4double theProb = 0.0; |
---|
130 | |
---|
131 | G4double ScaleFactor = 1.0; // playing with scale factors |
---|
132 | |
---|
133 | const G4double Uexcite = frag.GetExcitationEnergy(); |
---|
134 | G4double Uafter = Uexcite - exciteE; |
---|
135 | |
---|
136 | G4double normC = 3.0; |
---|
137 | |
---|
138 | const G4double upperLim = Uexcite; |
---|
139 | const G4double lowerLim = Uafter; |
---|
140 | const G4int numIters = 25; |
---|
141 | |
---|
142 | // Need to integrate EmissionProbDensity from lowerLim to upperLim |
---|
143 | // and multiply by normC |
---|
144 | |
---|
145 | G4double integ = normC * |
---|
146 | EmissionIntegration(frag,exciteE,lowerLim,upperLim,numIters); |
---|
147 | |
---|
148 | if(integ > 0.0) theProb = integ; |
---|
149 | |
---|
150 | return theProb * ScaleFactor; |
---|
151 | |
---|
152 | } |
---|
153 | |
---|
154 | G4double G4E1SingleProbability1::EmissionIntegration(const G4Fragment& frag, |
---|
155 | G4double , |
---|
156 | G4double lowLim, G4double upLim, |
---|
157 | G4int numIters) |
---|
158 | |
---|
159 | { |
---|
160 | |
---|
161 | // Simple Gaussian quadrature integration |
---|
162 | |
---|
163 | G4double x; |
---|
164 | const G4double root3 = 1.0/std::sqrt(3.0); |
---|
165 | |
---|
166 | G4double Step = (upLim-lowLim)/(2.0*numIters); |
---|
167 | G4double Delta = Step*root3; |
---|
168 | |
---|
169 | G4double mean = 0.0; |
---|
170 | |
---|
171 | G4double theInt = 0.0; |
---|
172 | |
---|
173 | for(G4int i = 0; i < numIters; i++) { |
---|
174 | |
---|
175 | x = (2*i + 1)/Step; |
---|
176 | G4double E1ProbDensityA = EmissionProbDensity(frag,x+Delta); |
---|
177 | G4double E1ProbDensityB = EmissionProbDensity(frag,x-Delta); |
---|
178 | |
---|
179 | mean += E1ProbDensityA + E1ProbDensityB; |
---|
180 | |
---|
181 | } |
---|
182 | |
---|
183 | if(mean*Step > 0.0) theInt = mean*Step; |
---|
184 | |
---|
185 | return theInt; |
---|
186 | |
---|
187 | } |
---|
188 | |
---|
189 | |
---|
190 | |
---|