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 | //J.M. Quesada (August2008). Based on: |
---|
27 | // |
---|
28 | // Hadronic Process: Nuclear De-excitations |
---|
29 | // by V. Lara (Oct 1998) |
---|
30 | // |
---|
31 | // Modif (03 September 2008) by J. M. Quesada for external choice of inverse |
---|
32 | // cross section option |
---|
33 | // JMQ (06 September 2008) Also external choices have been added for |
---|
34 | // superimposed Coulomb barrier (if useSICB is set true, by default is false) |
---|
35 | // |
---|
36 | // JMQ (14 february 2009) bug fixed in emission width: hbarc instead of hbar_Planck in the denominator |
---|
37 | // |
---|
38 | #include <iostream> |
---|
39 | using namespace std; |
---|
40 | |
---|
41 | #include "G4EvaporationProbability.hh" |
---|
42 | #include "G4PairingCorrection.hh" |
---|
43 | #include "G4ParticleTable.hh" |
---|
44 | #include "G4IonTable.hh" |
---|
45 | |
---|
46 | G4EvaporationProbability::G4EvaporationProbability(G4int anA, G4int aZ, |
---|
47 | G4double aGamma, |
---|
48 | G4VCoulombBarrier * aCoulombBarrier) |
---|
49 | : theA(anA), |
---|
50 | theZ(aZ), |
---|
51 | Gamma(aGamma), |
---|
52 | theCoulombBarrierptr(aCoulombBarrier) |
---|
53 | {} |
---|
54 | |
---|
55 | G4EvaporationProbability::G4EvaporationProbability() |
---|
56 | : theA(0), |
---|
57 | theZ(0), |
---|
58 | Gamma(0.0), |
---|
59 | theCoulombBarrierptr(0) |
---|
60 | {} |
---|
61 | |
---|
62 | G4EvaporationProbability::~G4EvaporationProbability() |
---|
63 | {} |
---|
64 | |
---|
65 | G4double |
---|
66 | G4EvaporationProbability::EmissionProbability(const G4Fragment & fragment, G4double anEnergy) |
---|
67 | { |
---|
68 | G4double EmissionProbability = 0.0; |
---|
69 | G4double MaximalKineticEnergy = anEnergy; |
---|
70 | |
---|
71 | if (MaximalKineticEnergy > 0.0 && fragment.GetExcitationEnergy() > 0.0) { |
---|
72 | EmissionProbability = CalculateProbability(fragment, MaximalKineticEnergy); |
---|
73 | |
---|
74 | } |
---|
75 | return EmissionProbability; |
---|
76 | } |
---|
77 | |
---|
78 | //////////////////////////////////// |
---|
79 | |
---|
80 | // Computes the integrated probability of evaporation channel |
---|
81 | G4double |
---|
82 | G4EvaporationProbability::CalculateProbability(const G4Fragment & fragment, |
---|
83 | G4double MaximalKineticEnergy) |
---|
84 | { |
---|
85 | G4int ResidualA = fragment.GetA_asInt() - theA; |
---|
86 | G4int ResidualZ = fragment.GetZ_asInt() - theZ; |
---|
87 | G4double U = fragment.GetExcitationEnergy(); |
---|
88 | |
---|
89 | if (OPTxs==0) { |
---|
90 | |
---|
91 | G4double NuclearMass = fragment.ComputeGroundStateMass(theZ,theA); |
---|
92 | |
---|
93 | G4double delta0 = fPairCorr->GetPairingCorrection(fragment.GetA_asInt(), |
---|
94 | fragment.GetZ_asInt()); |
---|
95 | |
---|
96 | G4double SystemEntropy = 2.0*std::sqrt( |
---|
97 | theEvapLDPptr->LevelDensityParameter(fragment.GetA_asInt(),fragment.GetZ_asInt(),U)* |
---|
98 | (U-delta0)); |
---|
99 | |
---|
100 | const G4double RN = 1.5*fermi; |
---|
101 | |
---|
102 | G4double Alpha = CalcAlphaParam(fragment); |
---|
103 | G4double Beta = CalcBetaParam(fragment); |
---|
104 | |
---|
105 | G4double Rmax = MaximalKineticEnergy; |
---|
106 | G4double a = theEvapLDPptr->LevelDensityParameter(ResidualA,ResidualZ,Rmax); |
---|
107 | G4double GlobalFactor = Gamma * Alpha/(a*a) * |
---|
108 | (NuclearMass*RN*RN*fG4pow->Z23(ResidualA))/ |
---|
109 | (twopi* hbar_Planck*hbar_Planck); |
---|
110 | G4double Term1 = (2.0*Beta*a-3.0)/2.0 + Rmax*a; |
---|
111 | G4double Term2 = (2.0*Beta*a-3.0)*std::sqrt(Rmax*a) + 2.0*a*Rmax; |
---|
112 | |
---|
113 | G4double ExpTerm1 = 0.0; |
---|
114 | if (SystemEntropy <= 600.0) { ExpTerm1 = std::exp(-SystemEntropy); } |
---|
115 | |
---|
116 | G4double ExpTerm2 = 2.*std::sqrt(a*Rmax) - SystemEntropy; |
---|
117 | if (ExpTerm2 > 700.0) { ExpTerm2 = 700.0; } |
---|
118 | ExpTerm2 = std::exp(ExpTerm2); |
---|
119 | |
---|
120 | G4double Width = GlobalFactor*(Term1*ExpTerm1 + Term2*ExpTerm2); |
---|
121 | |
---|
122 | return Width; |
---|
123 | |
---|
124 | } else if (OPTxs==1 || OPTxs==2 ||OPTxs==3 || OPTxs==4) { |
---|
125 | |
---|
126 | G4double EvaporatedMass = fragment.ComputeGroundStateMass(theZ,theA); |
---|
127 | G4double ResidulalMass = fragment.ComputeGroundStateMass(ResidualZ,ResidualA); |
---|
128 | G4double limit = std::max(0.0,fragment.GetGroundStateMass()-EvaporatedMass-ResidulalMass); |
---|
129 | if (useSICB) { |
---|
130 | limit = std::max(limit,theCoulombBarrierptr->GetCoulombBarrier(ResidualA,ResidualZ,U)); |
---|
131 | } |
---|
132 | |
---|
133 | if (MaximalKineticEnergy <= limit) { return 0.0; } |
---|
134 | |
---|
135 | // if Coulomb barrier cutoff is superimposed for all cross sections |
---|
136 | // then the limit is the Coulomb Barrier |
---|
137 | G4double LowerLimit= limit; |
---|
138 | |
---|
139 | //MaximalKineticEnergy: asimptotic value (already accounted for in G4EvaporationChannel) |
---|
140 | |
---|
141 | G4double UpperLimit = MaximalKineticEnergy; |
---|
142 | |
---|
143 | G4double Width = IntegrateEmissionProbability(fragment,LowerLimit,UpperLimit); |
---|
144 | |
---|
145 | return Width; |
---|
146 | } else { |
---|
147 | std::ostringstream errOs; |
---|
148 | errOs << "Bad option for cross sections at evaporation" <<G4endl; |
---|
149 | throw G4HadronicException(__FILE__, __LINE__, errOs.str()); |
---|
150 | } |
---|
151 | |
---|
152 | } |
---|
153 | |
---|
154 | ///////////////////////////////////////////////////////////////////// |
---|
155 | |
---|
156 | G4double G4EvaporationProbability:: |
---|
157 | IntegrateEmissionProbability(const G4Fragment & fragment, |
---|
158 | const G4double & Low, const G4double & Up ) |
---|
159 | { |
---|
160 | static const G4int N = 10; |
---|
161 | // 10-Points Gauss-Legendre abcisas and weights |
---|
162 | static const G4double w[N] = { |
---|
163 | 0.0666713443086881, |
---|
164 | 0.149451349150581, |
---|
165 | 0.219086362515982, |
---|
166 | 0.269266719309996, |
---|
167 | 0.295524224714753, |
---|
168 | 0.295524224714753, |
---|
169 | 0.269266719309996, |
---|
170 | 0.219086362515982, |
---|
171 | 0.149451349150581, |
---|
172 | 0.0666713443086881 |
---|
173 | }; |
---|
174 | static const G4double x[N] = { |
---|
175 | -0.973906528517172, |
---|
176 | -0.865063366688985, |
---|
177 | -0.679409568299024, |
---|
178 | -0.433395394129247, |
---|
179 | -0.148874338981631, |
---|
180 | 0.148874338981631, |
---|
181 | 0.433395394129247, |
---|
182 | 0.679409568299024, |
---|
183 | 0.865063366688985, |
---|
184 | 0.973906528517172 |
---|
185 | }; |
---|
186 | |
---|
187 | G4double Total = 0.0; |
---|
188 | |
---|
189 | |
---|
190 | for (G4int i = 0; i < N; i++) |
---|
191 | { |
---|
192 | |
---|
193 | G4double KineticE = ((Up-Low)*x[i]+(Up+Low))/2.0; |
---|
194 | |
---|
195 | Total += w[i]*ProbabilityDistributionFunction(fragment, KineticE); |
---|
196 | |
---|
197 | } |
---|
198 | Total *= (Up-Low)/2.0; |
---|
199 | return Total; |
---|
200 | } |
---|
201 | |
---|
202 | |
---|
203 | ///////////////////////////////////////////////////////// |
---|
204 | //New method (OPT=1,2,3,4) |
---|
205 | |
---|
206 | G4double |
---|
207 | G4EvaporationProbability::ProbabilityDistributionFunction( const G4Fragment & fragment, |
---|
208 | G4double K) |
---|
209 | { |
---|
210 | G4int ResidualA = fragment.GetA_asInt() - theA; |
---|
211 | G4int ResidualZ = fragment.GetZ_asInt() - theZ; |
---|
212 | G4double U = fragment.GetExcitationEnergy(); |
---|
213 | //G4cout << "### G4EvaporationProbability::ProbabilityDistributionFunction" << G4endl; |
---|
214 | //G4cout << "FragZ= " << fragment.GetZ_asInt() << " FragA= " << fragment.GetA_asInt() |
---|
215 | // << " Z= " << theZ << " A= " << theA << G4endl; |
---|
216 | //G4cout << "PC " << fPairCorr << " DP " << theEvapLDPptr << G4endl; |
---|
217 | |
---|
218 | // if(K <= theCoulombBarrierptr->GetCoulombBarrier(ResidualA,ResidualZ,U)) return 0.0; |
---|
219 | |
---|
220 | G4double delta1 = fPairCorr->GetPairingCorrection(ResidualA,ResidualZ); |
---|
221 | |
---|
222 | G4double delta0 = fPairCorr->GetPairingCorrection(fragment.GetA_asInt(), |
---|
223 | fragment.GetZ_asInt()); |
---|
224 | |
---|
225 | |
---|
226 | G4double ParticleMass = fragment.ComputeGroundStateMass(theZ,theA); |
---|
227 | G4double ResidualMass = fragment.ComputeGroundStateMass(ResidualZ,ResidualA); |
---|
228 | |
---|
229 | G4double theSeparationEnergy = ParticleMass + ResidualMass |
---|
230 | - fragment.GetGroundStateMass(); |
---|
231 | |
---|
232 | G4double a0 = theEvapLDPptr->LevelDensityParameter(fragment.GetA_asInt(), |
---|
233 | fragment.GetZ_asInt(), |
---|
234 | U - delta0); |
---|
235 | |
---|
236 | G4double a1 = theEvapLDPptr->LevelDensityParameter(ResidualA, ResidualZ, |
---|
237 | U - theSeparationEnergy - delta1); |
---|
238 | |
---|
239 | |
---|
240 | G4double E0 = U - delta0; |
---|
241 | |
---|
242 | G4double E1 = U - theSeparationEnergy - delta1 - K; |
---|
243 | |
---|
244 | if (E1<0.) { return 0.; } |
---|
245 | |
---|
246 | //JMQ 14/02/09 BUG fixed: hbarc should be in the denominator instead of hbar_Planck |
---|
247 | //Without 1/hbar_Panck remains as a width |
---|
248 | |
---|
249 | G4double Prob=Gamma*ParticleMass/((pi*hbarc)*(pi*hbarc)*std::exp(2*std::sqrt(a0*E0))) |
---|
250 | *K*CrossSection(fragment,K)*std::exp(2*std::sqrt(a1*E1))*millibarn; |
---|
251 | |
---|
252 | return Prob; |
---|
253 | } |
---|
254 | |
---|
255 | |
---|