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 | // |
---|
27 | // |
---|
28 | // Hadronic Process: Nuclear De-excitations |
---|
29 | // by V. Lara (Sept 2001) |
---|
30 | // |
---|
31 | |
---|
32 | |
---|
33 | #include "G4GEMProbability.hh" |
---|
34 | #include "G4PairingCorrection.hh" |
---|
35 | |
---|
36 | |
---|
37 | |
---|
38 | G4GEMProbability::G4GEMProbability(const G4GEMProbability &) : G4VEmissionProbability() |
---|
39 | { |
---|
40 | throw G4HadronicException(__FILE__, __LINE__, "G4GEMProbability::copy_constructor meant to not be accessable"); |
---|
41 | } |
---|
42 | |
---|
43 | |
---|
44 | |
---|
45 | |
---|
46 | const G4GEMProbability & G4GEMProbability:: |
---|
47 | operator=(const G4GEMProbability &) |
---|
48 | { |
---|
49 | throw G4HadronicException(__FILE__, __LINE__, "G4GEMProbability::operator= meant to not be accessable"); |
---|
50 | return *this; |
---|
51 | } |
---|
52 | |
---|
53 | |
---|
54 | G4bool G4GEMProbability::operator==(const G4GEMProbability &) const |
---|
55 | { |
---|
56 | return false; |
---|
57 | } |
---|
58 | |
---|
59 | G4bool G4GEMProbability::operator!=(const G4GEMProbability &) const |
---|
60 | { |
---|
61 | return true; |
---|
62 | } |
---|
63 | |
---|
64 | G4double G4GEMProbability::EmissionProbability(const G4Fragment & fragment, |
---|
65 | const G4double MaximalKineticEnergy) |
---|
66 | { |
---|
67 | G4double EmissionProbability = 0.0; |
---|
68 | |
---|
69 | if (MaximalKineticEnergy > 0.0 && fragment.GetExcitationEnergy() > 0.0) { |
---|
70 | G4double CoulombBarrier = GetCoulombBarrier(fragment); |
---|
71 | |
---|
72 | EmissionProbability = CalcProbability(fragment,MaximalKineticEnergy,CoulombBarrier); |
---|
73 | Normalization = EmissionProbability; |
---|
74 | // Next there is a loop over excited states for this channel summing probabilities |
---|
75 | if (ExcitationEnergies && ExcitationSpins && ExcitationLifetimes) { |
---|
76 | G4double SavedSpin = Spin; |
---|
77 | for (unsigned int i = 0; i < ExcitationEnergies->size(); i++) { |
---|
78 | Spin = ExcitationSpins->operator[](i); |
---|
79 | // substract excitation energies |
---|
80 | G4double Tmax = MaximalKineticEnergy - ExcitationEnergies->operator[](i); |
---|
81 | if (Tmax > 0.0) { |
---|
82 | G4double width = CalcProbability(fragment,Tmax,CoulombBarrier); |
---|
83 | // update probability |
---|
84 | if (hbar_Planck*std::log(2.0)/width < ExcitationLifetimes->operator[](i)) { |
---|
85 | EmissionProbability += width; |
---|
86 | } |
---|
87 | } |
---|
88 | } |
---|
89 | // Restore Spin |
---|
90 | Spin = SavedSpin; |
---|
91 | } |
---|
92 | } |
---|
93 | return EmissionProbability; |
---|
94 | } |
---|
95 | |
---|
96 | G4double G4GEMProbability::CalcProbability(const G4Fragment & fragment, |
---|
97 | const G4double MaximalKineticEnergy, |
---|
98 | const G4double V) |
---|
99 | // Calculate integrated probability (width) for evaporation channel |
---|
100 | { |
---|
101 | G4double ResidualA = static_cast<G4double>(fragment.GetA() - theA); |
---|
102 | // G4double ResidualZ = static_cast<G4double>(fragment.GetZ() - theZ); |
---|
103 | G4double U = fragment.GetExcitationEnergy(); |
---|
104 | |
---|
105 | G4double NuclearMass = G4ParticleTable::GetParticleTable()-> |
---|
106 | GetIonTable()->GetNucleusMass(theZ,theA); |
---|
107 | |
---|
108 | G4double a = theEvapLDPptr->LevelDensityParameter(static_cast<G4int>(fragment.GetA()), |
---|
109 | static_cast<G4int>(fragment.GetZ()),U); |
---|
110 | //090115 |
---|
111 | //G4double delta0 = G4PairingCorrection::GetInstance()->GetPairingCorrection(static_cast<G4int>(fragment.GetA()), |
---|
112 | // static_cast<G4int>(fragment.GetZ())); |
---|
113 | G4double ResidualZ = static_cast<G4double>(fragment.GetZ() - theZ); |
---|
114 | G4double delta0 = G4PairingCorrection::GetInstance()->GetPairingCorrection( static_cast<G4int>( ResidualA ) , static_cast<G4int>( ResidualZ ) ); |
---|
115 | //090115 |
---|
116 | |
---|
117 | G4double Alpha = CalcAlphaParam(fragment); |
---|
118 | G4double Beta = CalcBetaParam(fragment); |
---|
119 | |
---|
120 | |
---|
121 | G4double Ux = (2.5 + 150.0/ResidualA)*MeV; |
---|
122 | G4double Ex = Ux + delta0; |
---|
123 | G4double T = 1.0/(std::sqrt(a/Ux) - 1.5/Ux); |
---|
124 | G4double E0 = Ex - T*(std::log(T*MeV) - std::log(a/MeV)/4.0 - 1.25*std::log(Ux*MeV) + 2.0*std::sqrt(a*Ux)); |
---|
125 | |
---|
126 | G4double Width; |
---|
127 | G4double InitialLevelDensity; |
---|
128 | if ( MaximalKineticEnergy < Ex ) { |
---|
129 | G4double t = MaximalKineticEnergy/T; |
---|
130 | Width = (I1(t,t) + (Beta+V)*I0(t))/std::exp(E0/T); |
---|
131 | InitialLevelDensity = (pi/12.0)*std::exp((U-E0)/T)/T; |
---|
132 | } else { |
---|
133 | G4double t = MaximalKineticEnergy/T; |
---|
134 | G4double tx = Ex/T; |
---|
135 | G4double s = 2.0*std::sqrt(a*(MaximalKineticEnergy-delta0)); |
---|
136 | G4double sx = 2.0*std::sqrt(a*(Ex-delta0)); |
---|
137 | Width = I1(t,tx)/std::exp(E0/T) + I3(s,sx)*std::exp(s)/(std::sqrt(2.0)*a); |
---|
138 | // For charged particles (Beta+V) = 0 beacuse Beta = -V |
---|
139 | if (theZ == 0) { |
---|
140 | Width += (Beta+V)*(I0(tx)/std::exp(E0/T) + 2.0*std::sqrt(2.0)*I2(s,sx)*std::exp(s)); |
---|
141 | } |
---|
142 | InitialLevelDensity = (pi/12.0)*std::exp(2*std::sqrt(a*(U-delta0)))/std::pow(a*std::pow(U-delta0,5.0),1.0/4.0); |
---|
143 | } |
---|
144 | |
---|
145 | |
---|
146 | G4double g = (2.0*Spin+1.0)*NuclearMass/(pi2* hbar_Planck*hbar_Planck); |
---|
147 | |
---|
148 | G4double RN = 0.0; |
---|
149 | if (theA > 4) |
---|
150 | { |
---|
151 | G4double R1 = std::pow(ResidualA,1.0/3.0); |
---|
152 | G4double R2 = std::pow(G4double(theA),1.0/3.0); |
---|
153 | RN = 1.12*(R1 + R2) - 0.86*((R1+R2)/(R1*R2)); |
---|
154 | RN *= fermi; |
---|
155 | } |
---|
156 | else |
---|
157 | { |
---|
158 | RN = 1.5*fermi; |
---|
159 | } |
---|
160 | G4double GeometricalXS = pi*RN*RN*std::pow(ResidualA,2./3.); |
---|
161 | |
---|
162 | |
---|
163 | Width *= std::sqrt(pi)*g*GeometricalXS*Alpha/(12.0*InitialLevelDensity); |
---|
164 | return Width; |
---|
165 | } |
---|
166 | |
---|
167 | |
---|
168 | |
---|
169 | |
---|
170 | G4double G4GEMProbability::I0(const G4double t) |
---|
171 | { |
---|
172 | G4double result = (std::exp(t) - 1.0); |
---|
173 | return result; |
---|
174 | } |
---|
175 | |
---|
176 | G4double G4GEMProbability::I1(const G4double t, const G4double tx) |
---|
177 | { |
---|
178 | G4double result = t - tx + 1.0; |
---|
179 | result *= std::exp(tx); |
---|
180 | result -= (t + 1.0); |
---|
181 | return result; |
---|
182 | } |
---|
183 | |
---|
184 | |
---|
185 | G4double G4GEMProbability::I2(const G4double s, const G4double sx) |
---|
186 | { |
---|
187 | G4double S = 1.0/std::sqrt(s); |
---|
188 | G4double Sx = 1.0/std::sqrt(sx); |
---|
189 | |
---|
190 | G4double p1 = S*S*S*( 1.0 + S*S*( 1.5 + 3.75*S*S) ); |
---|
191 | G4double p2 = Sx*Sx*Sx*( 1.0 + Sx*Sx*( 1.5 + 3.75*Sx*Sx) )*std::exp(sx-s); |
---|
192 | |
---|
193 | return p1-p2; |
---|
194 | } |
---|
195 | |
---|
196 | G4double G4GEMProbability::I3(const G4double s, const G4double sx) |
---|
197 | { |
---|
198 | G4double s2 = s*s; |
---|
199 | G4double sx2 = sx*sx; |
---|
200 | G4double S = 1.0/std::sqrt(s); |
---|
201 | G4double S2 = S*S; |
---|
202 | G4double Sx = 1.0/std::sqrt(sx); |
---|
203 | G4double Sx2 = Sx*Sx; |
---|
204 | |
---|
205 | G4double p1 = S *(2.0 + S2 *( 4.0 + S2 *( 13.5 + S2 *( 60.0 + S2 * 325.125 )))); |
---|
206 | G4double p2 = Sx*Sx2 *( |
---|
207 | (s2-sx2) + Sx2 *( |
---|
208 | (1.5*s2+0.5*sx2) + Sx2 *( |
---|
209 | (3.75*s2+0.25*sx2) + Sx2 *( |
---|
210 | (12.875*s2+0.625*sx2) + Sx2 *( |
---|
211 | (59.0625*s2+0.9375*sx2) + Sx2 *(324.8*s2+3.28*sx2)))))); |
---|
212 | |
---|
213 | p2 *= std::exp(sx-s); |
---|
214 | |
---|
215 | return p1-p2; |
---|
216 | } |
---|