Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ProtonEvaporationProbability.cc
Go to the documentation of this file.
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$
27 //
28 // J.M. Quesada (August2008). Based on:
29 //
30 // Hadronic Process: Nuclear De-excitations
31 // by V. Lara (Oct 1998)
32 //
33 // Modified:
34 // 03-09-2008 J.M. Quesada for external choice of inverse cross section option
35 // 17-11-2010 V.Ivanchenko integer Z and A
36 
38 #include "G4PhysicalConstants.hh"
39 #include "G4SystemOfUnits.hh"
40 
42  G4EvaporationProbability(1,1,2,&theCoulombBarrier) // A,Z,Gamma,&theCoulombBarrier
43 {
44  ResidualA = ResidualZ = theA = theZ = FragmentA = 0;
45  ResidualAthrd = FragmentAthrd = U = 0.0;
46 }
47 
49 {}
50 
51 G4double G4ProtonEvaporationProbability::CalcAlphaParam(const G4Fragment & fragment)
52  { return 1.0 + CCoeficient(fragment.GetZ_asInt()-GetZ());}
53 
54 G4double G4ProtonEvaporationProbability::CalcBetaParam(const G4Fragment & )
55  { return 0.0; }
56 
57 G4double G4ProtonEvaporationProbability::CCoeficient(G4int aZ)
58 {
59  // Data comes from
60  // Dostrovsky, Fraenkel and Friedlander
61  // Physical Review, vol 116, num. 3 1959
62  //
63  // const G4int size = 5;
64  // G4double Zlist[5] = { 10.0, 20.0, 30.0, 50.0, 70.0};
65  // G4double Cp[5] = { 0.50, 0.28, 0.20, 0.15, 0.10};
66  G4double C = 0.0;
67 
68  if (aZ >= 70) {
69  C = 0.10;
70  } else {
71  C = ((((0.15417e-06*aZ) - 0.29875e-04)*aZ + 0.21071e-02)*aZ - 0.66612e-01)*aZ + 0.98375;
72  }
73 
74  return C;
75 
76 }
77 
79 //J. M. Quesada (Dec 2007-June 2008): New inverse reaction cross sections for protons
80 //OPT=0 Dostrovski's parameterization
81 //OPT=1 Chatterjee's parameterization
82 //OPT=2,4 Wellisch's parameterization
83 //OPT=3 Kalbach's parameterization
84 //
85 G4double
86 G4ProtonEvaporationProbability::CrossSection(const G4Fragment & fragment, G4double K)
87 {
88  // G4cout<<" In G4ProtonEVaporationProbability OPTxs="<<OPTxs<<G4endl;
89  // G4cout<<" In G4ProtonEVaporationProbability useSICB="<<useSICB<<G4endl;
90 
91  theA=GetA();
92  theZ=GetZ();
93  ResidualA=fragment.GetA_asInt()-theA;
94  ResidualZ=fragment.GetZ_asInt()-theZ;
95 
96  ResidualAthrd=fG4pow->Z13(ResidualA);
97  FragmentA=fragment.GetA_asInt();
98  FragmentAthrd=fG4pow->Z13(FragmentA);
99 
100  U=fragment.GetExcitationEnergy();
101 
102  if (OPTxs==0) {std::ostringstream errOs;
103  errOs << "We should'n be here (OPT =0) at evaporation cross section calculation (protons)!!" <<G4endl;
104  throw G4HadronicException(__FILE__, __LINE__, errOs.str());
105  return 0.;}
106  else if( OPTxs==1 ) return GetOpt1( K);
107  else if( OPTxs==2 ||OPTxs==4) return GetOpt2( K);
108  else if (OPTxs==3 ) return GetOpt3( K);
109  else{
110  std::ostringstream errOs;
111  errOs << "BAD PROTON CROSS SECTION OPTION AT EVAPORATION!!" <<G4endl;
112  throw G4HadronicException(__FILE__, __LINE__, errOs.str());
113  return 0.;
114  }
115 }
116 
117 //********************* OPT=1 : Chatterjee's cross section *********************
118 //(fitting to cross section from Bechetti & Greenles OM potential)
119 
120 G4double G4ProtonEvaporationProbability::GetOpt1(G4double K)
121 {
122  G4double Kc=K;
123 
124  // JMQ xsec is set constat above limit of validity
125  if (K > 50*MeV) { Kc = 50*MeV; }
126 
127  G4double landa, landa0, landa1, mu, mum0, mu1,nu, nu0, nu1, nu2,xs;
128  G4double p, p0, p1, p2,Ec,delta,q,r,ji;
129 
130  p0 = 15.72;
131  p1 = 9.65;
132  p2 = -449.0;
133  landa0 = 0.00437;
134  landa1 = -16.58;
135  mum0 = 244.7;
136  mu1 = 0.503;
137  nu0 = 273.1;
138  nu1 = -182.4;
139  nu2 = -1.872;
140  delta=0.;
141 
142  Ec = 1.44*theZ*ResidualZ/(1.5*ResidualAthrd+delta);
143  p = p0 + p1/Ec + p2/(Ec*Ec);
144  landa = landa0*ResidualA + landa1;
145 
146  G4double resmu1 = fG4pow->powZ(ResidualA,mu1);
147  mu = mum0*resmu1;
148  nu = resmu1*(nu0 + nu1*Ec + nu2*(Ec*Ec));
149  q = landa - nu/(Ec*Ec) - 2*p*Ec;
150  r = mu + 2*nu/Ec + p*(Ec*Ec);
151 
152  ji=std::max(Kc,Ec);
153  if(Kc < Ec) { xs = p*Kc*Kc + q*Kc + r;}
154  else {xs = p*(Kc - ji)*(Kc - ji) + landa*Kc + mu + nu*(2 - Kc/ji)/ji ;}
155  if (xs <0.0) {xs=0.0;}
156 
157  return xs;
158 }
159 
160 //************* OPT=2 : Welisch's proton reaction cross section ***************
161 
162 G4double G4ProtonEvaporationProbability::GetOpt2(G4double K)
163 {
164 
165  G4double eekin,ekin,ff1,ff2,ff3,r0,fac,fac1,fac2,b0,xine_th(0);
166 
167  // This is redundant when the Coulomb barrier is overimposed to all
168  // cross sections
169  // It should be kept when Coulomb barrier only imposed at OPTxs=2
170 
171  if(!useSICB && K<=theCoulombBarrier.GetCoulombBarrier(ResidualA,ResidualZ,U))
172  { return 0.0; }
173 
174  eekin=K;
175  G4int rnneu=ResidualA-ResidualZ;
176  ekin=eekin/1000;
177  r0=1.36*1.e-15;
178  fac=pi*r0*r0;
179  b0=2.247-0.915*(1.-1./ResidualAthrd);
180  fac1=b0*(1.-1./ResidualAthrd);
181  fac2=1.;
182  if(rnneu > 1.5) { fac2 = fG4pow->logZ(rnneu); }
183  xine_th= 1.e+31*fac*fac2*(1.+ResidualAthrd-fac1);
184  xine_th=(1.-0.15*std::exp(-ekin))*xine_th/(1.00-0.0007*ResidualA);
185  ff1=0.70-0.0020*ResidualA;
186  ff2=1.00+1/G4double(ResidualA);
187  ff3=0.8+18/G4double(ResidualA)-0.002*ResidualA;
188  fac=1.-(1./(1.+std::exp(-8.*ff1*(std::log10(ekin)+1.37*ff2))));
189  xine_th=xine_th*(1.+ff3*fac);
190  ff1=1.-1/G4double(ResidualA)-0.001*ResidualA;
191  ff2=1.17-2.7/G4double(ResidualA)-0.0014*ResidualA;
192  fac=-8.*ff1*(std::log10(ekin)+2.0*ff2);
193  fac=1./(1.+std::exp(fac));
194  xine_th=xine_th*fac;
195  if (xine_th < 0.0){
196  std::ostringstream errOs;
197  G4cout<<"WARNING: negative Wellisch cross section "<<G4endl;
198  errOs << "RESIDUAL: A=" << ResidualA << " Z=" << ResidualZ <<G4endl;
199  errOs <<" xsec("<<ekin<<" MeV) ="<<xine_th <<G4endl;
200  throw G4HadronicException(__FILE__, __LINE__, errOs.str());
201  }
202  return xine_th;
203 }
204 
205 // *********** OPT=3 : Kalbach's cross sections (from PRECO code)*************
206 G4double G4ProtonEvaporationProbability::GetOpt3(const G4double K)
207 {
208  // ** p from becchetti and greenlees (but modified with sub-barrier
209  // ** correction function and xp2 changed from -449)
210 
211  G4double landa, landa0, landa1, mu, mum0, mu1,nu, nu0, nu1, nu2;
212  G4double p, p0, p1, p2;
213  p0 = 15.72;
214  p1 = 9.65;
215  p2 = -300.;
216  landa0 = 0.00437;
217  landa1 = -16.58;
218  mum0 = 244.7;
219  mu1 = 0.503;
220  nu0 = 273.1;
221  nu1 = -182.4;
222  nu2 = -1.872;
223 
224  // parameters for proton cross section refinement
225  /*
226  G4double afit,bfit,a2,b2;
227  afit=-0.0785656;
228  bfit=5.10789;
229  a2= -0.00089076;
230  b2= 0.0231597;
231  */
232  G4double ec,ecsq,xnulam,etest(0.),ra(0.),a,w,c,signor(1.),signor2,sig;
233  G4double b,ecut,cut,ecut2,geom,elab;
234 
235  G4double flow = 1.e-18;
236  G4double spill= 1.e+18;
237 
238  if (ResidualA <= 60) { signor = 0.92; }
239  else if (ResidualA < 100) { signor = 0.8 + ResidualA*0.002; }
240 
241  ec = 1.44 * theZ * ResidualZ / (1.5*ResidualAthrd+ra);
242  ecsq = ec * ec;
243  p = p0 + p1/ec + p2/ecsq;
244  landa = landa0*ResidualA + landa1;
245  a = fG4pow->powZ(ResidualA,mu1);
246  mu = mum0 * a;
247  nu = a* (nu0+nu1*ec+nu2*ecsq);
248 
249  c =std::min(3.15,ec*0.5);
250  w = 0.7 * c / 3.15;
251 
252  xnulam = nu / landa;
253  if (xnulam > spill) { xnulam=0.; }
254  if (xnulam >= flow) { etest =std::sqrt(xnulam) + 7.; }
255 
256  a = -2.*p*ec + landa - nu/ecsq;
257  b = p*ecsq + mu + 2.*nu/ec;
258  ecut = 0.;
259  cut = a*a - 4.*p*b;
260  if (cut > 0.) { ecut = std::sqrt(cut); }
261  ecut = (ecut-a) / (p+p);
262  ecut2 = ecut;
263  //JMQ 290310 for avoiding unphysical increase below minimum (at ecut)
264  // ecut<0 means that there is no cut with energy axis, i.e. xs is set
265  // to 0 bellow minimum
266  // if (cut < 0.) ecut2 = ecut - 2.;
267  if (cut < 0.) { ecut2 = ecut; }
268  elab = K * FragmentA /G4double(ResidualA);
269  sig = 0.;
270  if (elab <= ec) { //start for E<Ec
271  if (elab > ecut2) { sig = (p*elab*elab+a*elab+b) * signor; }
272 
273  signor2 = (ec-elab-c) / w;
274  signor2 = 1. + std::exp(signor2);
275  sig = sig / signor2;
276  } //end for E<=Ec
277  else{ //start for E>Ec
278  sig = (landa*elab+mu+nu/elab) * signor;
279  geom = 0.;
280 
281  if (xnulam < flow || elab < etest)
282  {
283  if (sig <0.0) {sig=0.0;}
284  return sig;
285  }
286  geom = std::sqrt(theA*K);
287  geom = 1.23*ResidualAthrd + ra + 4.573/geom;
288  geom = 31.416 * geom * geom;
289  sig = std::max(geom,sig);
290 
291  } //end for E>Ec
292  return sig;
293 }
294