Geant4  10.01.p02
G4ChipsProtonInelasticXS.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 //
27 // The lust update: M.V. Kossov, CERN/ITEP(Moscow) 17-June-02
28 //
29 //
30 // G4 Physics class: G4ChipsProtonInelasticXS for gamma+A cross sections
31 // Created: M.V. Kossov, CERN/ITEP(Moscow), 20-Dec-03
32 // The last update: M.V. Kossov, CERN/ITEP (Moscow) 15-Feb-04
33 //
34 //
35 // ****************************************************************************************
36 // Short description: Cross-sections extracted (by W.Pokorski) from the CHIPS package for
37 // proton-nuclear interactions. Original author: M. Kossov
38 // -------------------------------------------------------------------------------------
39 //
40 
41 
43 #include "G4SystemOfUnits.hh"
44 #include "G4DynamicParticle.hh"
45 #include "G4ParticleDefinition.hh"
46 #include "G4Proton.hh"
47 
48 // factory
49 #include "G4CrossSectionFactory.hh"
50 //
52 
54 {
55  // Initialization of the
56  lastLEN=0; // Pointer to the lastArray of LowEn CS
57  lastHEN=0; // Pointer to the lastArray of HighEn CS
58  lastN=0; // The last N of calculated nucleus
59  lastZ=0; // The last Z of calculated nucleus
60  lastP=0.; // Last used in cross section Momentum
61  lastTH=0.; // Last threshold momentum
62  lastCS=0.; // Last value of the Cross Section
63  lastI=0; // The last position in the DAMDB
64 
65  LEN = new std::vector<G4double*>;
66  HEN = new std::vector<G4double*>;
67 }
68 
70 {
71  G4int lens=LEN->size();
72  for(G4int i=0; i<lens; ++i) delete[] (*LEN)[i];
73  delete LEN;
74  G4int hens=HEN->size();
75  for(G4int i=0; i<hens; ++i) delete[] (*HEN)[i];
76  delete HEN;
77 }
78 
80  const G4Element*,
81  const G4Material*)
82 {
83  return true;
84 }
85 
86 
87 // The main member function giving the collision cross section (P is in IU, CS is in mb)
88 // Make pMom in independent units ! (Now it is MeV)
90  const G4Isotope*,
91  const G4Element*,
92  const G4Material*)
93 {
94  G4double pMom=Pt->GetTotalMomentum();
95  G4int tgN = A - tgZ;
96 
97  return GetChipsCrossSection(pMom, tgZ, tgN, 2212);
98 }
99 
101 {
102 
103  G4bool in=false; // By default the isotope must be found in the AMDB
104  if(tgN!=lastN || tgZ!=lastZ) // The nucleus was not the last used isotope
105  {
106  in = false; // By default the isotope haven't been found in AMDB
107  lastP = 0.; // New momentum history (nothing to compare with)
108  lastN = tgN; // The last N of the calculated nucleus
109  lastZ = tgZ; // The last Z of the calculated nucleus
110  lastI = colN.size(); // Size of the Associative Memory DB in the heap
111  j = 0; // A#0f records found in DB for this projectile
112  if(lastI) for(G4int i=0; i<lastI; i++) // AMDB exists, try to find the (Z,N) isotope
113  {
114  if(colN[i]==tgN && colZ[i]==tgZ) // Try the record "i" in the AMDB
115  {
116  lastI=i; // Remember the index for future fast/last use
117  lastTH =colTH[i]; // The last THreshold (A-dependent)
118  if(pMom<=lastTH)
119  {
120  return 0.; // Energy is below the Threshold value
121  }
122  lastP =colP [i]; // Last Momentum (A-dependent)
123  lastCS =colCS[i]; // Last CrossSect (A-dependent)
124  in = true; // This is the case when the isotop is found in DB
125  // Momentum pMom is in IU ! @@ Units
126  lastCS=CalculateCrossSection(-1,j,2212,lastZ,lastN,pMom); // read & update
127  if(lastCS<=0. && pMom>lastTH) // Correct the threshold (@@ No intermediate Zeros)
128  {
129  lastCS=0.;
130  lastTH=pMom;
131  }
132  break; // Go out of the LOOP
133  }
134  j++; // Increment a#0f records found in DB
135  }
136  if(!in) // This isotope has not been calculated previously
137  {
139  lastCS=CalculateCrossSection(0,j,2212,lastZ,lastN,pMom); //calculate & create
140  //if(lastCS>0.) // It means that the AMBD was initialized
141  //{
142 
143  lastTH = 0; //ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
144  colN.push_back(tgN);
145  colZ.push_back(tgZ);
146  colP.push_back(pMom);
147  colTH.push_back(lastTH);
148  colCS.push_back(lastCS);
149  //} // M.K. Presence of H1 with high threshold breaks the syncronization
150  return lastCS*millibarn;
151  } // End of creation of the new set of parameters
152  else
153  {
154  colP[lastI]=pMom;
155  colCS[lastI]=lastCS;
156  }
157  } // End of parameters udate
158  else if(pMom<=lastTH)
159  {
160  return 0.; // Momentum is below the Threshold Value -> CS=0
161  }
162  else // It is the last used -> use the current tables
163  {
164  lastCS=CalculateCrossSection(1,j,2212,lastZ,lastN,pMom); // Only read and UpdateDB
165  lastP=pMom;
166  }
167  return lastCS*millibarn;
168 }
169 
170 // The main member function giving the gamma-A cross section (E in GeV, CS in mb)
172  G4int, G4int targZ, G4int targN, G4double Momentum)
173 {
174  static const G4double THmin=27.; // default minimum Momentum (MeV/c) Threshold
175  static const G4double THmiG=THmin*.001; // minimum Momentum (GeV/c) Threshold
176  static const G4double dP=10.; // step for the LEN (Low ENergy) table MeV/c
177  static const G4double dPG=dP*.001; // step for the LEN (Low ENergy) table GeV/c
178  static const G4int nL=105; // A#of LEN points in E (step 10 MeV/c)
179  static const G4double Pmin=THmin+(nL-1)*dP; // minP for the HighE part with safety
180  static const G4double Pmax=227000.; // maxP for the HEN (High ENergy) part 227 GeV
181  static const G4int nH=224; // A#of HEN points in lnE
182  static const G4double milP=std::log(Pmin);// Low logarithm energy for the HEN part
183  static const G4double malP=std::log(Pmax);// High logarithm energy (each 2.75 percent)
184  static const G4double dlP=(malP-milP)/(nH-1); // Step in log energy in the HEN part
185  static const G4double milPG=std::log(.001*Pmin);// Low logarithmEnergy for HEN part GeV/c
186  G4double sigma=0.;
187  if(F&&I) sigma=0.; // @@ *!* Fake line *!* to use F & I !!!Temporary!!!
188  //G4double A=targN+targZ; // A of the target
189  if(F<=0) // This isotope was not the last used isotop
190  {
191  if(F<0) // This isotope was found in DAMDB =-----=> RETRIEVE
192  {
193  G4int sync=LEN->size();
194  if(sync<=I) G4cout<<"*!*G4QProtonNuclCS::CalcCrossSect:Sync="<<sync<<"<="<<I<<G4endl;
195  lastLEN=(*LEN)[I]; // Pointer to prepared LowEnergy cross sections
196  lastHEN=(*HEN)[I]; // Pointer to prepared High Energy cross sections
197  }
198  else // This isotope wasn't calculated before => CREATE
199  {
200  lastLEN = new G4double[nL]; // Allocate memory for the new LEN cross sections
201  lastHEN = new G4double[nH]; // Allocate memory for the new HEN cross sections
202  // --- Instead of making a separate function ---
203  G4double P=THmiG; // Table threshold in GeV/c
204  for(G4int k=0; k<nL; k++)
205  {
206  lastLEN[k] = CrossSectionLin(targZ, targN, P);
207  P+=dPG;
208  }
209  G4double lP=milPG;
210  for(G4int n=0; n<nH; n++)
211  {
212  lastHEN[n] = CrossSectionLog(targZ, targN, lP);
213  lP+=dlP;
214  }
215  // --- End of possible separate function
216  // *** The synchronization check ***
217  G4int sync=LEN->size();
218  if(sync!=I)
219  {
220  G4cout<<"***G4ChipsProtonNuclCS::CalcCrossSect: Sinc="<<sync<<"#"<<I<<", Z=" <<targZ
221  <<", N="<<targN<<", F="<<F<<G4endl;
222  //G4Exception("G4ProtonNuclearCS::CalculateCS:","39",FatalException,"overflow DB");
223  }
224  LEN->push_back(lastLEN); // remember the Low Energy Table
225  HEN->push_back(lastHEN); // remember the High Energy Table
226  } // End of creation of the new set of parameters
227  } // End of parameters udate
228  // =------------------= NOW the Magic Formula =-----------------------=
229  if (Momentum<lastTH) return 0.; // It must be already checked in the interface class
230  else if (Momentum<Pmin) // High Energy region
231  {
232  sigma=EquLinearFit(Momentum,nL,THmin,dP,lastLEN);
233  }
234  else if (Momentum<Pmax) // High Energy region
235  {
236  G4double lP=std::log(Momentum);
237  sigma=EquLinearFit(lP,nH,milP,dlP,lastHEN);
238  }
239  else // UHE region (calculation, not frequent)
240  {
241  G4double P=0.001*Momentum; // Approximation formula is for P in GeV/c
242  sigma=CrossSectionFormula(targZ, targN, P, std::log(P));
243  }
244  if(sigma<0.) return 0.;
245  return sigma;
246 }
247 
248 // Electromagnetic momentum-threshold (in MeV/c)
250 {
251  static const G4double third=1./3.;
252  static const G4double pM = G4Proton::Proton()->Definition()->GetPDGMass(); // Projectile mass in MeV
253  static const G4double tpM= pM+pM; // Doubled projectile mass (MeV)
254 
255  G4double tA=tZ+tN;
256  if(tZ<.99 || tN<0.) return 0.;
257  else if(tZ==1 && tN==0) return 800.; // A threshold on the free proton
258  //G4double dE=1.263*tZ/(1.+std::pow(tA,third));
259  G4double dE=tZ/(1.+std::pow(tA,third)); // Safety for diffused edge of the nucleus (QE)
260  G4double tM=931.5*tA;
261  G4double T=dE+dE*(dE/2+pM)/tM;
262  return std::sqrt(T*(tpM+T));
263 }
264 
265 // Calculation formula for proton-nuclear inelastic cross-section (mb) (P in GeV/c)
267 {
268  G4double sigma=0.;
269  if(P<ThresholdMomentum(tZ,tN)*.001) return sigma;
270  G4double lP=std::log(P);
271  if(tZ==1&&!tN){if(P>.35) sigma=CrossSectionFormula(tZ,tN,P,lP);}// s(pp)=0 below 350Mev/c
272  else if(tZ<97 && tN<152) // General solution
273  {
274  G4double pex=0.;
275  G4double pos=0.;
276  G4double wid=1.;
277  if(tZ==13 && tN==14) // Excited metastable states
278  {
279  pex=230.;
280  pos=.13;
281  wid=8.e-5;
282  }
283  else if(tZ<7)
284  {
285  if(tZ==6 && tN==6)
286  {
287  pex=320.;
288  pos=.14;
289  wid=7.e-6;
290  }
291  else if(tZ==5 && tN==6)
292  {
293  pex=270.;
294  pos=.17;
295  wid=.002;
296  }
297  else if(tZ==4 && tN==5)
298  {
299  pex=600.;
300  pos=.132;
301  wid=.005;
302  }
303  else if(tZ==3 && tN==4)
304  {
305  pex=280.;
306  pos=.19;
307  wid=.0025;
308  }
309  else if(tZ==3 && tN==3)
310  {
311  pex=370.;
312  pos=.171;
313  wid=.006;
314  }
315  else if(tZ==2 && tN==1)
316  {
317  pex=30.;
318  pos=.22;
319  wid=.0005;
320  }
321  }
322  sigma=CrossSectionFormula(tZ,tN,P,lP);
323  if(pex>0.)
324  {
325  G4double dp=P-pos;
326  sigma+=pex*std::exp(-dp*dp/wid);
327  }
328  }
329  else
330  {
331  G4cerr<<"-Warning-G4ChipsProtonNuclearXS::CSLin:*Bad A* Z="<<tZ<<", N="<<tN<<G4endl;
332  sigma=0.;
333  }
334  if(sigma<0.) return 0.;
335  return sigma;
336 }
337 
338 // Calculation formula for proton-nuclear inelastic cross-section (mb) log(P in GeV/c)
340 {
341  G4double P=std::exp(lP);
342  return CrossSectionFormula(tZ, tN, P, lP);
343 }
344 // Calculation formula for proton-nuclear inelastic cross-section (mb) log(P in GeV/c)
346  G4double P, G4double lP)
347 {
348  G4double sigma=0.;
349  if(tZ==1 && !tN) // pp interaction (from G4QuasiElasticRatios)
350  {
351  G4double El(0.),To(0.); // Uzhi
352  if(P<0.1) // Copied from G4QuasiElasticRatios Uzhi / start
353  {
354  G4double p2=P*P;
355  El=1./(0.00012+p2*0.2);
356  To=El;
357  }
358  else if(P>1000.)
359  {
360  G4double lp=std::log(P)-3.5;
361  G4double lp2=lp*lp;
362  El=0.0557*lp2+6.72;
363  To=0.3*lp2+38.2;
364  }
365  else
366  {
367  G4double p2=P*P;
368  G4double LE=1./(0.00012+p2*0.2);
369  G4double lp=std::log(P)-3.5;
370  G4double lp2=lp*lp;
371  G4double rp2=1./p2;
372  El=LE+(0.0557*lp2+6.72+32.6/P)/(1.+rp2/P);
373  To=LE+(0.3 *lp2+38.2+52.7*rp2)/(1.+2.72*rp2*rp2);
374  } // Copied from G4QuasiElasticRatios Uzhi / end
375 
376 /* // Uzhi 4.03.2013
377  G4double p2=P*P;
378  G4double lp=lP-3.5;
379  G4double lp2=lp*lp;
380  G4double rp2=1./p2;
381  G4double El=(.0557*lp2+6.72+30./P)/(1.+.49*rp2/P);
382  G4double To=(.3*lp2+38.2)/(1.+.54*rp2*rp2);
383 */ // Uzhi 4.03.2013
384 
385  sigma=To-El;
386  }
387  else if(tZ<97 && tN<152) // General solution
388  {
389  //G4double lP=std::log(P); // Already calculated
390  G4double d=lP-4.2;
391  G4double p2=P*P;
392  G4double p4=p2*p2;
393  G4double a=tN+tZ; // A of the target
394  G4double al=std::log(a);
395  G4double sa=std::sqrt(a);
396  G4double a2=a*a;
397  G4double a2s=a2*sa;
398  G4double a4=a2*a2;
399  G4double a8=a4*a4;
400  G4double a12=a8*a4;
401  G4double a16=a8*a8;
402  G4double c=(170.+3600./a2s)/(1.+65./a2s);
403  G4double dl=al-3.;
404  G4double dl2=dl*dl;
405  G4double r=.21+.62*dl2/(1.+.5*dl2);
406  G4double gg=40.*std::exp(al*0.712)/(1.+12.2/a)/(1.+34./a2);
407  G4double e=318.+a4/(1.+.0015*a4/std::exp(al*0.09))/(1.+4.e-28*a12)+
408  8.e-18/(1./a16+1.3e-20)/(1.+1.e-21*a12);
409  G4double ss=3.57+.009*a2/(1.+.0001*a2*a);
410  G4double h=(.01/a4+2.5e-6/a)*(1.+6.e-6*a2*a)/(1.+6.e7/a12/a2);
411  sigma=(c+d*d)/(1.+r/p4)+(gg+e*std::exp(-ss*P))/(1.+h/p4/p4);
412  }
413  else
414  {
415  G4cerr<<"-Warning-G4QProtonNuclearCroSect::CSForm:*Bad A* Z="<<tZ<<", N="<<tN<<G4endl;
416  sigma=0.;
417  }
418  if(sigma<0.) return 0.;
419  return sigma;
420 }
421 
423 {
424  if(DX<=0. || N<2)
425  {
426  G4cerr<<"***G4ChipsProtonInelasticXS::EquLinearFit: DX="<<DX<<", N="<<N<<G4endl;
427  return Y[0];
428  }
429 
430  G4int N2=N-2;
431  G4double d=(X-X0)/DX;
432  G4int jj=static_cast<int>(d);
433  if (jj<0) jj=0;
434  else if(jj>N2) jj=N2;
435  d-=jj; // excess
436  G4double yi=Y[jj];
437  G4double sigma=yi+(Y[jj+1]-yi)*d;
438 
439  return sigma;
440 }
G4double CrossSectionLin(G4int targZ, G4int targN, G4double P)
static const G4int nH
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
G4double ThresholdMomentum(G4int targZ, G4int targN)
static const G4double a4
G4double a
Definition: TRTMaterials.hh:39
static const G4int nL
G4double CrossSectionFormula(G4int targZ, G4int targN, G4double P, G4double lP)
int G4int
Definition: G4Types.hh:78
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int tgZ, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4double GetTotalMomentum() const
G4_DECLARE_XS_FACTORY(G4ChipsProtonInelasticXS)
static G4Proton * Definition()
Definition: G4Proton.cc:49
static const G4double dE
G4double EquLinearFit(G4double X, G4int N, G4double X0, G4double DX, G4double *Y)
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
const G4double p2
static G4Proton * Proton()
Definition: G4Proton.cc:93
virtual G4bool IsIsoApplicable(const G4DynamicParticle *Pt, G4int Z, G4int A, const G4Element *elm, const G4Material *mat)
const G4int n
static const G4double A[nN]
std::vector< G4double * > * LEN
Definition: Evaluator.cc:66
std::vector< G4double > colP
G4double GetPDGMass() const
G4double CrossSectionLog(G4int targZ, G4int targN, G4double lP)
static const double millibarn
Definition: G4SIunits.hh:96
std::vector< G4double > colCS
#define G4endl
Definition: G4ios.hh:61
std::vector< G4double * > * HEN
std::vector< G4double > colTH
double G4double
Definition: G4Types.hh:76
static const G4double pos
static const G4double THmin
G4double CalculateCrossSection(G4int F, G4int I, G4int PDG, G4int Z, G4int N, G4double Momentum)
static const G4double a2
G4GLOB_DLL std::ostream G4cerr