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 | |
---|
44 | |
---|
45 | |
---|
46 | G4EvaporationProbability::G4EvaporationProbability(const G4EvaporationProbability &) : G4VEmissionProbability() |
---|
47 | { |
---|
48 | throw G4HadronicException(__FILE__, __LINE__, "G4EvaporationProbability::copy_constructor meant to not be accessable"); |
---|
49 | } |
---|
50 | |
---|
51 | |
---|
52 | |
---|
53 | |
---|
54 | const G4EvaporationProbability & G4EvaporationProbability:: |
---|
55 | operator=(const G4EvaporationProbability &) |
---|
56 | { |
---|
57 | throw G4HadronicException(__FILE__, __LINE__, "G4EvaporationProbability::operator= meant to not be accessable"); |
---|
58 | return *this; |
---|
59 | } |
---|
60 | |
---|
61 | |
---|
62 | G4bool G4EvaporationProbability::operator==(const G4EvaporationProbability &) const |
---|
63 | { |
---|
64 | return false; |
---|
65 | } |
---|
66 | |
---|
67 | G4bool G4EvaporationProbability::operator!=(const G4EvaporationProbability &) const |
---|
68 | { |
---|
69 | return true; |
---|
70 | } |
---|
71 | |
---|
72 | G4double G4EvaporationProbability::EmissionProbability(const G4Fragment & fragment, const G4double anEnergy) |
---|
73 | { |
---|
74 | G4double EmissionProbability = 0.0; |
---|
75 | G4double MaximalKineticEnergy = anEnergy; |
---|
76 | |
---|
77 | if (MaximalKineticEnergy > 0.0 && fragment.GetExcitationEnergy() > 0.0) { |
---|
78 | EmissionProbability = CalculateProbability(fragment, MaximalKineticEnergy); |
---|
79 | |
---|
80 | } |
---|
81 | return EmissionProbability; |
---|
82 | } |
---|
83 | |
---|
84 | //////////////////////////////////// |
---|
85 | |
---|
86 | // Computes the integrated probability of evaporation channel |
---|
87 | G4double G4EvaporationProbability::CalculateProbability(const G4Fragment & fragment, const G4double MaximalKineticEnergy) |
---|
88 | { |
---|
89 | G4double ResidualA = fragment.GetA() - theA; |
---|
90 | G4double ResidualZ = fragment.GetZ() - theZ; |
---|
91 | G4double U = fragment.GetExcitationEnergy(); |
---|
92 | |
---|
93 | if (OPTxs==0) { |
---|
94 | |
---|
95 | |
---|
96 | G4double NuclearMass = G4ParticleTable::GetParticleTable()->GetIonTable()->GetNucleusMass(theZ,theA); |
---|
97 | |
---|
98 | |
---|
99 | G4double delta0 = G4PairingCorrection::GetInstance()->GetPairingCorrection(static_cast<G4int>(fragment.GetA()), |
---|
100 | static_cast<G4int>(fragment.GetZ())); |
---|
101 | |
---|
102 | G4double SystemEntropy = 2.0*std::sqrt(theEvapLDPptr->LevelDensityParameter(static_cast<G4int>(fragment.GetA()), |
---|
103 | static_cast<G4int>(fragment.GetZ()),U)* |
---|
104 | (U-delta0)); |
---|
105 | |
---|
106 | |
---|
107 | G4double RN = 1.5*fermi; |
---|
108 | |
---|
109 | G4double Alpha = CalcAlphaParam(fragment); |
---|
110 | G4double Beta = CalcBetaParam(fragment); |
---|
111 | |
---|
112 | G4double Rmax = MaximalKineticEnergy; |
---|
113 | G4double a = theEvapLDPptr->LevelDensityParameter(static_cast<G4int>(ResidualA), |
---|
114 | static_cast<G4int>(ResidualZ), |
---|
115 | Rmax); |
---|
116 | G4double GlobalFactor = static_cast<G4double>(Gamma) * (Alpha/(a*a)) * |
---|
117 | (NuclearMass*RN*RN*std::pow(ResidualA,2./3.))/ |
---|
118 | (2.*pi* hbar_Planck*hbar_Planck); |
---|
119 | G4double Term1 = (2.0*Beta*a-3.0)/2.0 + Rmax*a; |
---|
120 | G4double Term2 = (2.0*Beta*a-3.0)*std::sqrt(Rmax*a) + 2.0*a*Rmax; |
---|
121 | |
---|
122 | G4double ExpTerm1 = 0.0; |
---|
123 | if (SystemEntropy <= 600.0) ExpTerm1 = std::exp(-SystemEntropy); |
---|
124 | |
---|
125 | G4double ExpTerm2 = 2.*std::sqrt(a*Rmax) - SystemEntropy; |
---|
126 | if (ExpTerm2 > 700.0) ExpTerm2 = 700.0; |
---|
127 | ExpTerm2 = std::exp(ExpTerm2); |
---|
128 | |
---|
129 | G4double Width = GlobalFactor*(Term1*ExpTerm1 + Term2*ExpTerm2); |
---|
130 | |
---|
131 | return Width; |
---|
132 | |
---|
133 | } else if (OPTxs==1 || OPTxs==2 ||OPTxs==3 || OPTxs==4) { |
---|
134 | |
---|
135 | G4double limit; |
---|
136 | if (useSICB) limit=theCoulombBarrierptr->GetCoulombBarrier(G4lrint(ResidualA),G4lrint(ResidualZ),U); |
---|
137 | else limit=0.; |
---|
138 | |
---|
139 | if (MaximalKineticEnergy <= limit) return 0.0; |
---|
140 | |
---|
141 | |
---|
142 | // if Coulomb barrier cutoff is superimposed for all cross sections the limit is the Coulomb Barrier |
---|
143 | G4double LowerLimit= limit; |
---|
144 | |
---|
145 | // MaximalKineticEnergy: asimptotic value (already accounted for in G4EvaporationChannel) |
---|
146 | |
---|
147 | G4double UpperLimit = MaximalKineticEnergy; |
---|
148 | |
---|
149 | |
---|
150 | G4double Width = IntegrateEmissionProbability(fragment,LowerLimit,UpperLimit); |
---|
151 | |
---|
152 | return Width; |
---|
153 | } else{ |
---|
154 | std::ostringstream errOs; |
---|
155 | errOs << "Bad option for cross sections at evaporation" <<G4endl; |
---|
156 | throw G4HadronicException(__FILE__, __LINE__, errOs.str()); |
---|
157 | } |
---|
158 | |
---|
159 | } |
---|
160 | |
---|
161 | ///////////////////////////////////////////////////////////////////// |
---|
162 | |
---|
163 | G4double G4EvaporationProbability:: |
---|
164 | IntegrateEmissionProbability(const G4Fragment & fragment, const G4double & Low, const G4double & Up ) |
---|
165 | { |
---|
166 | |
---|
167 | static const G4int N = 10; |
---|
168 | // 10-Points Gauss-Legendre abcisas and weights |
---|
169 | static const G4double w[N] = { |
---|
170 | 0.0666713443086881, |
---|
171 | 0.149451349150581, |
---|
172 | 0.219086362515982, |
---|
173 | 0.269266719309996, |
---|
174 | 0.295524224714753, |
---|
175 | 0.295524224714753, |
---|
176 | 0.269266719309996, |
---|
177 | 0.219086362515982, |
---|
178 | 0.149451349150581, |
---|
179 | 0.0666713443086881 |
---|
180 | }; |
---|
181 | static const G4double x[N] = { |
---|
182 | -0.973906528517172, |
---|
183 | -0.865063366688985, |
---|
184 | -0.679409568299024, |
---|
185 | -0.433395394129247, |
---|
186 | -0.148874338981631, |
---|
187 | 0.148874338981631, |
---|
188 | 0.433395394129247, |
---|
189 | 0.679409568299024, |
---|
190 | 0.865063366688985, |
---|
191 | 0.973906528517172 |
---|
192 | }; |
---|
193 | |
---|
194 | G4double Total = 0.0; |
---|
195 | |
---|
196 | |
---|
197 | for (G4int i = 0; i < N; i++) |
---|
198 | { |
---|
199 | |
---|
200 | G4double KineticE = ((Up-Low)*x[i]+(Up+Low))/2.0; |
---|
201 | |
---|
202 | Total += w[i]*ProbabilityDistributionFunction(fragment, KineticE); |
---|
203 | |
---|
204 | } |
---|
205 | Total *= (Up-Low)/2.0; |
---|
206 | return Total; |
---|
207 | } |
---|
208 | |
---|
209 | |
---|
210 | ///////////////////////////////////////////////////////// |
---|
211 | //New method (OPT=1,2,3,4) |
---|
212 | |
---|
213 | G4double G4EvaporationProbability::ProbabilityDistributionFunction( const G4Fragment & fragment, const G4double K) |
---|
214 | { |
---|
215 | |
---|
216 | |
---|
217 | |
---|
218 | |
---|
219 | G4double ResidualA = fragment.GetA() - theA; |
---|
220 | G4double ResidualZ = fragment.GetZ() - theZ; |
---|
221 | G4double U = fragment.GetExcitationEnergy(); |
---|
222 | |
---|
223 | // if(K <= theCoulombBarrierptr->GetCoulombBarrier(G4lrint(ResidualA),G4lrint(ResidualZ),U)) return 0.0; |
---|
224 | |
---|
225 | G4double delta1 = G4PairingCorrection::GetInstance()->GetPairingCorrection(static_cast<G4int>(ResidualA),static_cast<G4int>(ResidualZ)); |
---|
226 | |
---|
227 | |
---|
228 | G4double delta0 = G4PairingCorrection::GetInstance()->GetPairingCorrection(static_cast<G4int>(fragment.GetA()),static_cast<G4int>(fragment.GetZ())); |
---|
229 | |
---|
230 | |
---|
231 | G4double ParticleMass = G4ParticleTable::GetParticleTable()->GetIonTable()->GetNucleusMass(theZ,theA); |
---|
232 | |
---|
233 | G4double theSeparationEnergy= G4NucleiProperties::GetMassExcess(static_cast<G4int>(theA),static_cast<G4int>(theZ)) + |
---|
234 | G4NucleiProperties::GetMassExcess(static_cast<G4int>(ResidualA),static_cast<G4int>(ResidualZ)) - |
---|
235 | G4NucleiProperties::GetMassExcess(static_cast<G4int>(fragment.GetA()),static_cast<G4int>(fragment.GetZ())); |
---|
236 | |
---|
237 | G4double a0 = theEvapLDPptr->LevelDensityParameter(static_cast<G4int>(fragment.GetA()),static_cast<G4int>(fragment.GetZ()),U - delta0); |
---|
238 | |
---|
239 | G4double a1 = theEvapLDPptr->LevelDensityParameter(static_cast<G4int>(ResidualA),static_cast<G4int>(ResidualZ),U - theSeparationEnergy - delta1); |
---|
240 | |
---|
241 | |
---|
242 | G4double E0=U-delta0; |
---|
243 | |
---|
244 | G4double E1=U-theSeparationEnergy-delta1-K; |
---|
245 | |
---|
246 | if (E1<0.) return 0.; |
---|
247 | |
---|
248 | //JMQ 14/02/09 BUG fixed: hbarc should be in the denominator instead of hbar_Planck |
---|
249 | //Without 1/hbar_Panck remains as a width |
---|
250 | // G4double Prob=Gamma*ParticleMass/((pi*hbar_Planck)*(pi*hbar_Planck)* |
---|
251 | //std::exp(2*std::sqrt(a0*E0)))*K*CrossSection(fragment,K)*std::exp(2*std::sqrt(a1*E1))*millibarn; |
---|
252 | |
---|
253 | G4double Prob=Gamma*ParticleMass/((pi*hbarc)*(pi*hbarc)*std::exp(2*std::sqrt(a0*E0))) |
---|
254 | *K*CrossSection(fragment,K)*std::exp(2*std::sqrt(a1*E1))*millibarn; |
---|
255 | |
---|
256 | return Prob; |
---|
257 | } |
---|
258 | |
---|
259 | |
---|