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: G4EvaporationChannel.cc,v 1.19 2010/11/24 15:30:49 vnivanch Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-04-ref-00 $ |
---|
28 | // |
---|
29 | //J.M. Quesada (August2008). Based on: |
---|
30 | // |
---|
31 | // Hadronic Process: Nuclear De-excitations |
---|
32 | // by V. Lara (Oct 1998) |
---|
33 | // |
---|
34 | // Modified: |
---|
35 | // 03-09-2008 J.M. Quesada for external choice of inverse cross section option |
---|
36 | // 06-09-2008 J.M. Quesada Also external choices have been added for superimposed |
---|
37 | // Coulomb barrier (if useSICB is set true, by default is false) |
---|
38 | // 17-11-2010 V.Ivanchenko in constructor replace G4VEmissionProbability by |
---|
39 | // G4EvaporationProbability and do not new and delete probability |
---|
40 | // object at each call; use G4Pow |
---|
41 | |
---|
42 | #include "G4EvaporationChannel.hh" |
---|
43 | #include "G4PairingCorrection.hh" |
---|
44 | #include "G4NucleiProperties.hh" |
---|
45 | #include "G4Pow.hh" |
---|
46 | #include "G4EvaporationLevelDensityParameter.hh" |
---|
47 | #include "Randomize.hh" |
---|
48 | #include "G4Alpha.hh" |
---|
49 | |
---|
50 | G4EvaporationChannel::G4EvaporationChannel(G4int anA, G4int aZ, |
---|
51 | const G4String & aName, |
---|
52 | G4EvaporationProbability* aEmissionStrategy, |
---|
53 | G4VCoulombBarrier* aCoulombBarrier): |
---|
54 | G4VEvaporationChannel(aName), |
---|
55 | theA(anA), |
---|
56 | theZ(aZ), |
---|
57 | theEvaporationProbabilityPtr(aEmissionStrategy), |
---|
58 | theCoulombBarrierPtr(aCoulombBarrier), |
---|
59 | EmissionProbability(0.0), |
---|
60 | MaximalKineticEnergy(-1000.0) |
---|
61 | { |
---|
62 | ResidualA = 0; |
---|
63 | ResidualZ = 0; |
---|
64 | ResidualMass = CoulombBarrier=0.0; |
---|
65 | EvaporatedMass = G4NucleiProperties::GetNuclearMass(theA, theZ); |
---|
66 | theLevelDensityPtr = new G4EvaporationLevelDensityParameter; |
---|
67 | } |
---|
68 | |
---|
69 | G4EvaporationChannel::G4EvaporationChannel(): |
---|
70 | G4VEvaporationChannel(""), |
---|
71 | theA(0), |
---|
72 | theZ(0), |
---|
73 | theEvaporationProbabilityPtr(0), |
---|
74 | theCoulombBarrierPtr(0), |
---|
75 | EmissionProbability(0.0), |
---|
76 | MaximalKineticEnergy(-1000.0) |
---|
77 | { |
---|
78 | ResidualA = 0; |
---|
79 | ResidualZ = 0; |
---|
80 | EvaporatedMass = ResidualMass = CoulombBarrier = 0.0; |
---|
81 | theLevelDensityPtr = new G4EvaporationLevelDensityParameter; |
---|
82 | } |
---|
83 | |
---|
84 | G4EvaporationChannel::~G4EvaporationChannel() |
---|
85 | { |
---|
86 | delete theLevelDensityPtr; |
---|
87 | } |
---|
88 | |
---|
89 | void G4EvaporationChannel::Initialize(const G4Fragment & fragment) |
---|
90 | { |
---|
91 | //for inverse cross section choice |
---|
92 | theEvaporationProbabilityPtr->SetOPTxs(OPTxs); |
---|
93 | // for superimposed Coulomb Barrier for inverse cross sections |
---|
94 | theEvaporationProbabilityPtr->UseSICB(useSICB); |
---|
95 | |
---|
96 | G4int FragmentA = fragment.GetA_asInt(); |
---|
97 | G4int FragmentZ = fragment.GetZ_asInt(); |
---|
98 | ResidualA = FragmentA - theA; |
---|
99 | ResidualZ = FragmentZ - theZ; |
---|
100 | //G4cout << "G4EvaporationChannel::Initialize Z= " << theZ << " A= " << theA |
---|
101 | // << " FragZ= " << FragmentZ << " FragA= " << FragmentA << G4endl; |
---|
102 | |
---|
103 | //Effective excitation energy |
---|
104 | G4double ExEnergy = fragment.GetExcitationEnergy() - |
---|
105 | G4PairingCorrection::GetInstance()->GetPairingCorrection(FragmentA,FragmentZ); |
---|
106 | |
---|
107 | // Only channels which are physically allowed are taken into account |
---|
108 | if (ResidualA <= 0 || ResidualZ <= 0 || ResidualA < ResidualZ || |
---|
109 | (ResidualA == ResidualZ && ResidualA > 1) || ExEnergy <= 0.0) { |
---|
110 | CoulombBarrier = ResidualMass = 0.0; |
---|
111 | MaximalKineticEnergy = -1000.0*MeV; |
---|
112 | EmissionProbability = 0.0; |
---|
113 | } else { |
---|
114 | ResidualMass = G4NucleiProperties::GetNuclearMass(ResidualA, ResidualZ); |
---|
115 | G4double FragmentMass = fragment.GetGroundStateMass(); |
---|
116 | CoulombBarrier = theCoulombBarrierPtr->GetCoulombBarrier(ResidualA,ResidualZ,ExEnergy); |
---|
117 | // Maximal Kinetic Energy |
---|
118 | MaximalKineticEnergy = CalcMaximalKineticEnergy(FragmentMass + ExEnergy); |
---|
119 | //MaximalKineticEnergy = ExEnergy + fragment.GetGroundStateMass() |
---|
120 | // - EvaporatedMass - ResidualMass; |
---|
121 | |
---|
122 | // Emission probability |
---|
123 | // Protection for the case Tmax<V. If not set in this way we could end up in an |
---|
124 | // infinite loop in the method GetKineticEnergy if OPTxs!=0 && useSICB=true. |
---|
125 | // Of course for OPTxs=0 we have the Coulomb barrier |
---|
126 | |
---|
127 | G4double limit; |
---|
128 | if (OPTxs==0 || (OPTxs!=0 && useSICB)) |
---|
129 | limit= CoulombBarrier; |
---|
130 | else limit=0.; |
---|
131 | limit = std::max(limit, FragmentMass - ResidualMass - EvaporatedMass); |
---|
132 | |
---|
133 | // The threshold for charged particle emission must be set to 0 if Coulomb |
---|
134 | //cutoff is included in the cross sections |
---|
135 | if (MaximalKineticEnergy <= limit) EmissionProbability = 0.0; |
---|
136 | else { |
---|
137 | // Total emission probability for this channel |
---|
138 | EmissionProbability = theEvaporationProbabilityPtr-> |
---|
139 | EmissionProbability(fragment, MaximalKineticEnergy); |
---|
140 | } |
---|
141 | } |
---|
142 | //G4cout << "G4EvaporationChannel:: probability= " << EmissionProbability << G4endl; |
---|
143 | |
---|
144 | return; |
---|
145 | } |
---|
146 | |
---|
147 | G4FragmentVector * G4EvaporationChannel::BreakUp(const G4Fragment & theNucleus) |
---|
148 | { |
---|
149 | /* |
---|
150 | G4double Ecm = GetKineticEnergy(theNucleus) + ResidualMass + EvaporatedMass; |
---|
151 | |
---|
152 | G4double EvaporatedEnergy = |
---|
153 | ((Ecm-ResidualMass)*(Ecm+ResidualMass) + EvaporatedMass*EvaporatedMass)/(2*Ecm); |
---|
154 | */ |
---|
155 | G4double EvaporatedEnergy = GetKineticEnergy(theNucleus) + EvaporatedMass; |
---|
156 | |
---|
157 | G4ThreeVector momentum(IsotropicVector |
---|
158 | (std::sqrt((EvaporatedEnergy - EvaporatedMass)* |
---|
159 | (EvaporatedEnergy + EvaporatedMass)))); |
---|
160 | |
---|
161 | G4LorentzVector EvaporatedMomentum(momentum,EvaporatedEnergy); |
---|
162 | G4LorentzVector ResidualMomentum = theNucleus.GetMomentum(); |
---|
163 | EvaporatedMomentum.boost(ResidualMomentum.boostVector()); |
---|
164 | |
---|
165 | G4Fragment * EvaporatedFragment = new G4Fragment(theA,theZ,EvaporatedMomentum); |
---|
166 | ResidualMomentum -= EvaporatedMomentum; |
---|
167 | |
---|
168 | G4Fragment * ResidualFragment = new G4Fragment(ResidualA, ResidualZ, ResidualMomentum); |
---|
169 | |
---|
170 | G4FragmentVector * theResult = new G4FragmentVector; |
---|
171 | |
---|
172 | #ifdef debug |
---|
173 | G4double Efinal = ResidualMomentum.e() + EvaporatedMomentum.e(); |
---|
174 | G4ThreeVector Pfinal = ResidualMomentum.vect() + EvaporatedMomentum.vect(); |
---|
175 | if (std::abs(Efinal-theNucleus.GetMomentum().e()) > 1.0*keV) { |
---|
176 | G4cout << "@@@@@@@@@@@@@@@@@@@@@ G4Evaporation Chanel: ENERGY @@@@@@@@@@@@@@@@" << G4endl; |
---|
177 | G4cout << "Initial : " << theNucleus.GetMomentum().e()/MeV << " MeV Final :" |
---|
178 | <<Efinal/MeV << " MeV Delta: " << (Efinal-theNucleus.GetMomentum().e())/MeV |
---|
179 | << " MeV" << G4endl; |
---|
180 | } |
---|
181 | if (std::abs(Pfinal.x()-theNucleus.GetMomentum().x()) > 1.0*keV || |
---|
182 | std::abs(Pfinal.y()-theNucleus.GetMomentum().y()) > 1.0*keV || |
---|
183 | std::abs(Pfinal.z()-theNucleus.GetMomentum().z()) > 1.0*keV ) { |
---|
184 | G4cout << "@@@@@@@@@@@@@@@@@@@@@ G4Evaporation Chanel: MOMENTUM @@@@@@@@@@@@@@@@" << G4endl; |
---|
185 | G4cout << "Initial : " << theNucleus.GetMomentum().vect() << " MeV Final :" |
---|
186 | <<Pfinal/MeV << " MeV Delta: " << Pfinal-theNucleus.GetMomentum().vect() |
---|
187 | << " MeV" << G4endl; |
---|
188 | |
---|
189 | } |
---|
190 | #endif |
---|
191 | theResult->push_back(EvaporatedFragment); |
---|
192 | theResult->push_back(ResidualFragment); |
---|
193 | return theResult; |
---|
194 | } |
---|
195 | |
---|
196 | ///////////////////////////////////////// |
---|
197 | // Calculates the maximal kinetic energy that can be carried by fragment. |
---|
198 | G4double G4EvaporationChannel::CalcMaximalKineticEnergy(G4double NucleusTotalE) |
---|
199 | { |
---|
200 | // This is the "true" assimptotic kinetic energy (from energy conservation) |
---|
201 | G4double Tmax = |
---|
202 | ((NucleusTotalE-ResidualMass)*(NucleusTotalE+ResidualMass) + EvaporatedMass*EvaporatedMass) |
---|
203 | /(2.0*NucleusTotalE) - EvaporatedMass; |
---|
204 | |
---|
205 | //JMQ (13-09-08) bug fixed: in the original version the Tmax is calculated |
---|
206 | //at the Coulomb barrier |
---|
207 | //IMPORTANT: meaning of Tmax differs in OPTxs=0 and OPTxs!=0 |
---|
208 | //When OPTxs!=0 Tmax is the TRUE (assimptotic) maximal kinetic energy |
---|
209 | |
---|
210 | if(OPTxs==0) |
---|
211 | Tmax=Tmax- CoulombBarrier; |
---|
212 | |
---|
213 | return Tmax; |
---|
214 | } |
---|
215 | |
---|
216 | /////////////////////////////////////////// |
---|
217 | //JMQ: New method for MC sampling of kinetic energy. Substitutes old CalcKineticEnergy |
---|
218 | G4double G4EvaporationChannel::GetKineticEnergy(const G4Fragment & aFragment) |
---|
219 | { |
---|
220 | |
---|
221 | if (OPTxs==0) { |
---|
222 | // It uses Dostrovsky's approximation for the inverse reaction cross |
---|
223 | // in the probability for fragment emission |
---|
224 | // MaximalKineticEnergy energy in the original version (V.Lara) was calculated at |
---|
225 | //the Coulomb barrier. |
---|
226 | |
---|
227 | |
---|
228 | if (MaximalKineticEnergy < 0.0) { |
---|
229 | throw G4HadronicException(__FILE__, __LINE__, |
---|
230 | "G4EvaporationChannel::CalcKineticEnergy: maximal kinetic at the Coulomb barrier is less than 0"); |
---|
231 | } |
---|
232 | G4double Rb = 4.0*theLevelDensityPtr-> |
---|
233 | LevelDensityParameter(ResidualA+theA,ResidualZ+theZ,MaximalKineticEnergy)* |
---|
234 | MaximalKineticEnergy; |
---|
235 | G4double RbSqrt = std::sqrt(Rb); |
---|
236 | G4double PEX1 = 0.0; |
---|
237 | if (RbSqrt < 160.0) PEX1 = std::exp(-RbSqrt); |
---|
238 | G4double Rk = 0.0; |
---|
239 | G4double FRk = 0.0; |
---|
240 | do { |
---|
241 | G4double RandNumber = G4UniformRand(); |
---|
242 | Rk = 1.0 + (1./RbSqrt)*std::log(RandNumber + (1.0-RandNumber)*PEX1); |
---|
243 | G4double Q1 = 1.0; |
---|
244 | G4double Q2 = 1.0; |
---|
245 | if (theZ == 0) { // for emitted neutron |
---|
246 | G4double Beta = (2.12/G4Pow::GetInstance()->Z23(ResidualA) - 0.05)*MeV/ |
---|
247 | (0.76 + 2.2/G4Pow::GetInstance()->Z13(ResidualA)); |
---|
248 | Q1 = 1.0 + Beta/(MaximalKineticEnergy); |
---|
249 | Q2 = Q1*std::sqrt(Q1); |
---|
250 | } |
---|
251 | |
---|
252 | FRk = (3.0*std::sqrt(3.0)/2.0)/Q2 * Rk * (Q1 - Rk*Rk); |
---|
253 | |
---|
254 | } while (FRk < G4UniformRand()); |
---|
255 | |
---|
256 | G4double result = MaximalKineticEnergy * (1.0-Rk*Rk) + CoulombBarrier; |
---|
257 | return result; |
---|
258 | } else if (OPTxs==1 || OPTxs==2 || OPTxs==3 || OPTxs==4) { |
---|
259 | |
---|
260 | // Coulomb barrier is just included in the cross sections |
---|
261 | G4double V = 0; |
---|
262 | if(useSICB) { V= CoulombBarrier; } |
---|
263 | |
---|
264 | V = std::max(V, aFragment.GetGroundStateMass()-EvaporatedMass-ResidualMass); |
---|
265 | |
---|
266 | G4double Tmax=MaximalKineticEnergy; |
---|
267 | G4double T(0.0); |
---|
268 | G4double NormalizedProbability(1.0); |
---|
269 | |
---|
270 | // VI: This is very ineffective - create new objects at each call to the method |
---|
271 | /* |
---|
272 | // A pointer is created in order to access the distribution function. |
---|
273 | G4EvaporationProbability * G4EPtemp = 0; |
---|
274 | |
---|
275 | if (theA==1 && theZ==0) G4EPtemp=new G4NeutronEvaporationProbability(); |
---|
276 | else if (theA==1 && theZ==1) G4EPtemp=new G4ProtonEvaporationProbability(); |
---|
277 | else if (theA==2 && theZ==1 ) G4EPtemp=new G4DeuteronEvaporationProbability(); |
---|
278 | else if (theA==3 && theZ==1 ) G4EPtemp=new G4TritonEvaporationProbability(); |
---|
279 | else if (theA==3 && theZ==2 ) G4EPtemp=new G4He3EvaporationProbability(); |
---|
280 | else if (theA==4 && theZ==2) G4EPtemp=new G4AlphaEvaporationProbability(); |
---|
281 | else { |
---|
282 | std::ostringstream errOs; |
---|
283 | errOs << "ejected particle out of range in G4EvaporationChannel" << G4endl; |
---|
284 | throw G4HadronicException(__FILE__, __LINE__, errOs.str()); |
---|
285 | } |
---|
286 | |
---|
287 | //for cross section selection and superimposed Coulom Barrier for xs |
---|
288 | G4EPtemp->SetOPTxs(OPTxs); |
---|
289 | G4EPtemp->UseSICB(useSICB); |
---|
290 | */ |
---|
291 | |
---|
292 | // use local pointer and not create a new one |
---|
293 | do |
---|
294 | { |
---|
295 | T=V+G4UniformRand()*(Tmax-V); |
---|
296 | NormalizedProbability = |
---|
297 | theEvaporationProbabilityPtr->ProbabilityDistributionFunction(aFragment,T)/ |
---|
298 | GetEmissionProbability(); |
---|
299 | |
---|
300 | } |
---|
301 | while (G4UniformRand() > NormalizedProbability); |
---|
302 | // delete G4EPtemp; |
---|
303 | return T; |
---|
304 | } else{ |
---|
305 | std::ostringstream errOs; |
---|
306 | errOs << "Bad option for energy sampling in evaporation" <<G4endl; |
---|
307 | throw G4HadronicException(__FILE__, __LINE__, errOs.str()); |
---|
308 | } |
---|
309 | } |
---|
310 | |
---|
311 | G4ThreeVector G4EvaporationChannel::IsotropicVector(G4double Magnitude) |
---|
312 | // Samples a isotropic random vectorwith a magnitud given by Magnitude. |
---|
313 | // By default Magnitude = 1.0 |
---|
314 | { |
---|
315 | G4double CosTheta = 1.0 - 2.0*G4UniformRand(); |
---|
316 | G4double SinTheta = std::sqrt(1.0 - CosTheta*CosTheta); |
---|
317 | G4double Phi = twopi*G4UniformRand(); |
---|
318 | G4ThreeVector Vector(Magnitude*std::cos(Phi)*SinTheta, |
---|
319 | Magnitude*std::sin(Phi)*SinTheta, |
---|
320 | Magnitude*CosTheta); |
---|
321 | return Vector; |
---|
322 | } |
---|
323 | |
---|
324 | |
---|
325 | |
---|
326 | |
---|
327 | |
---|
328 | |
---|
329 | |
---|
330 | |
---|
331 | |
---|
332 | |
---|