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