Geant4  10.01.p02
G4ChipsPionPlusInelasticXS.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: G4ChipsPionPlusInelasticXS 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 // Short description: Cross-sections extracted (by W.Pokorski) from the CHIPS package for
36 // pion interactions. Original author: M. Kossov
37 // -------------------------------------------------------------------------------------
38 //
39 
41 #include "G4SystemOfUnits.hh"
42 #include "G4DynamicParticle.hh"
43 #include "G4ParticleDefinition.hh"
44 #include "G4PionPlus.hh"
45 
46 // factory
47 #include "G4CrossSectionFactory.hh"
48 //
50 
52 {
53  // Initialization of the
54  lastLEN=0; // Pointer to lastArray of LowEn CS
55  lastHEN=0; // Pointer to lastArray of HighEn CS
56  lastN=0; // The last N of calculated nucleus
57  lastZ=0; // The last Z of calculated nucleus
58  lastP=0.; // Last used in cross section Momentum
59  lastTH=0.; // Last threshold momentum
60  lastCS=0.; // Last value of the Cross Section
61  lastI=0; // The last position in the DAMDB
62  LEN = new std::vector<G4double*>;
63  HEN = new std::vector<G4double*>;
64 }
65 
66 
68 {
69  G4int lens=LEN->size();
70  for(G4int i=0; i<lens; ++i) delete[] (*LEN)[i];
71  delete LEN;
72  G4int hens=HEN->size();
73  for(G4int i=0; i<hens; ++i) delete[] (*HEN)[i];
74  delete HEN;
75 }
76 
78  const G4Element*,
79  const G4Material*)
80 {
81  return true;
82 }
83 
84 // The main member function giving the collision cross section (P is in IU, CS is in mb)
85 // Make pMom in independent units ! (Now it is MeV)
87  const G4Isotope*,
88  const G4Element*,
89  const G4Material*)
90 {
91  G4double pMom=Pt->GetTotalMomentum();
92  G4int tgN = A - tgZ;
93 
94  return GetChipsCrossSection(pMom, tgZ, tgN, 211);
95 }
96 
97 
99 {
100 
101  G4bool in=false; // By default the isotope must be found in the AMDB
102  if(tgN!=lastN || tgZ!=lastZ) // The nucleus was not the last used isotope
103  {
104  in = false; // By default the isotope haven't be found in AMDB
105  lastP = 0.; // New momentum history (nothing to compare with)
106  lastN = tgN; // The last N of the calculated nucleus
107  lastZ = tgZ; // The last Z of the calculated nucleus
108  lastI = colN.size(); // Size of the Associative Memory DB in the heap
109  j = 0; // A#0f records found in DB for this projectile
110  if(lastI) for(G4int i=0; i<lastI; i++) // AMDB exists, try to find the (Z,N) isotope
111  {
112  if(colN[i]==tgN && colZ[i]==tgZ) // Try the record "i" in the AMDB
113  {
114  lastI=i; // Remember the index for future fast/last use
115  lastTH =colTH[i]; // The last THreshold (A-dependent)
116  if(pMom<=lastTH)
117  {
118  return 0.; // Energy is below the Threshold value
119  }
120  lastP =colP [i]; // Last Momentum (A-dependent)
121  lastCS =colCS[i]; // Last CrossSect (A-dependent)
122  in = true; // This is the case when the isotop is found in DB
123  // Momentum pMom is in IU ! @@ Units
124  lastCS=CalculateCrossSection(-1,j,211,lastZ,lastN,pMom); // read & update
125  if(lastCS<=0. && pMom>lastTH) // Correct the threshold (@@ No intermediate Zeros)
126  {
127  lastCS=0.;
128  lastTH=pMom;
129  }
130  break; // Go out of the LOOP
131  }
132  j++; // Increment a#0f records found in DB
133  }
134  if(!in) // This isotope has not been calculated previously
135  {
137  lastCS=CalculateCrossSection(0,j,211,lastZ,lastN,pMom); //calculate & create
138  //if(lastCS>0.) // It means that the AMBD was initialized
139  //{
140 
141  lastTH = 0; //ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
142  colN.push_back(tgN);
143  colZ.push_back(tgZ);
144  colP.push_back(pMom);
145  colTH.push_back(lastTH);
146  colCS.push_back(lastCS);
147  //} // M.K. Presence of H1 with high threshold breaks the syncronization
148  return lastCS*millibarn;
149  } // End of creation of the new set of parameters
150  else
151  {
152  colP[lastI]=pMom;
153  colCS[lastI]=lastCS;
154  }
155  } // End of parameters udate
156  else if(pMom<=lastTH)
157  {
158  return 0.; // Momentum is below the Threshold Value -> CS=0
159  }
160  else // It is the last used -> use the current tables
161  {
162  lastCS=CalculateCrossSection(1,j,211,lastZ,lastN,pMom); // Only read and UpdateDB
163  lastP=pMom;
164  }
165  return lastCS*millibarn;
166 }
167 
168 // The main member function giving the gamma-A cross section (E in GeV, CS in mb)
170  G4int, G4int targZ, G4int targN, G4double Momentum)
171 {
172  static const G4double THmin=27.; // default minimum Momentum (MeV/c) Threshold
173  static const G4double THmiG=THmin*.001; // minimum Momentum (GeV/c) Threshold
174  static const G4double dP=10.; // step for the LEN (Low ENergy) table MeV/c
175  static const G4double dPG=dP*.001; // step for the LEN (Low ENergy) table GeV/c
176  static const G4int nL=105; // A#of LEN points in E (step 10 MeV/c)
177  static const G4double Pmin=THmin+(nL-1)*dP; // minP for the HighE part with safety
178  static const G4double Pmax=227000.; // maxP for the HEN (High ENergy) part 227 GeV
179  static const G4int nH=224; // A#of HEN points in lnE
180  static const G4double milP=std::log(Pmin);// Low logarithm energy for the HEN part
181  static const G4double malP=std::log(Pmax);// High logarithm energy (each 2.75 percent)
182  static const G4double dlP=(malP-milP)/(nH-1); // Step in log energy in the HEN part
183  static const G4double milPG=std::log(.001*Pmin);// Low logarithmEnergy for HEN part GeV/c
184  G4double sigma=0.;
185  if(F&&I) sigma=0.; // @@ *!* Fake line *!* to use F & I !!!Temporary!!!
186  //G4double A=targN+targZ; // A of the target
187  if(F<=0) // This isotope was not the last used isotop
188  {
189  if(F<0) // This isotope was found in DAMDB =-----=> RETRIEVE
190  {
191  G4int sync=LEN->size();
192  if(sync<=I) G4cerr<<"*!*G4ChipsPiMinusNuclCS::CalcCrosSect:Sync="<<sync<<"<="<<I<<G4endl;
193  lastLEN=(*LEN)[I]; // Pointer to prepared LowEnergy cross sections
194  lastHEN=(*HEN)[I]; // Pointer to prepared High Energy cross sections
195  }
196  else // This isotope wasn't calculated before => CREATE
197  {
198  lastLEN = new G4double[nL]; // Allocate memory for the new LEN cross sections
199  lastHEN = new G4double[nH]; // Allocate memory for the new HEN cross sections
200  // --- Instead of making a separate function ---
201  G4double P=THmiG; // Table threshold in GeV/c
202  for(G4int k=0; k<nL; k++)
203  {
204  lastLEN[k] = CrossSectionLin(targZ, targN, P);
205  P+=dPG;
206  }
207  G4double lP=milPG;
208  for(G4int n=0; n<nH; n++)
209  {
210  lastHEN[n] = CrossSectionLog(targZ, targN, lP);
211  lP+=dlP;
212  }
213  // --- End of possible separate function
214  // *** The synchronization check ***
215  G4int sync=LEN->size();
216  if(sync!=I)
217  {
218  G4cerr<<"***G4ChipsPiMinusNuclCS::CalcCrossSect: Sinc="<<sync<<"#"<<I<<", Z=" <<targZ
219  <<", N="<<targN<<", F="<<F<<G4endl;
220  //G4Exception("G4PiMinusNuclearCS::CalculateCS:","39",FatalException,"DBoverflow");
221  }
222  LEN->push_back(lastLEN); // remember the Low Energy Table
223  HEN->push_back(lastHEN); // remember the High Energy Table
224  } // End of creation of the new set of parameters
225  } // End of parameters udate
226  // =-----------------= NOW the Magic Formula =-------------------------=
227  if (Momentum<lastTH) return 0.; // It must be already checked in the interface class
228  else if (Momentum<Pmin) // High Energy region
229  {
230  sigma=EquLinearFit(Momentum,nL,THmin,dP,lastLEN);
231  }
232  else if (Momentum<Pmax) // High Energy region
233  {
234  G4double lP=std::log(Momentum);
235  sigma=EquLinearFit(lP,nH,milP,dlP,lastHEN);
236  }
237  else // UHE region (calculation, not frequent)
238  {
239  G4double P=0.001*Momentum; // Approximation formula is for P in GeV/c
240  sigma=CrossSectionFormula(targZ, targN, P, std::log(P));
241  }
242  if(sigma<0.) return 0.;
243  return sigma;
244 }
245 
246 // Electromagnetic momentum-threshold (in MeV/c)
248 {
249  static const G4double third=1./3.;
250  static const G4double pM = G4PionPlus::PionPlus()->Definition()->GetPDGMass(); // Projectile mass in MeV
251  static const G4double tpM= pM+pM; // Doubled projectile mass (MeV)
252  G4double tA=tZ+tN;
253  if(tZ<.99 || tN<0.) return 0.;
254  else if(tZ==1 && tN==0) return 300.; // A threshold on the free proton
255  //G4double dE=1.263*tZ/(1.+std::pow(tA,third));
256  G4double dE=tZ/(1.+std::pow(tA,third)); // Safety for diffused edge of the nucleus (QE)
257  G4double tM=931.5*tA;
258  G4double T=dE+dE*(dE/2+pM)/tM;
259  return std::sqrt(T*(tpM+T));
260 }
261 
262 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) (P in GeV/c)
264 {
265  G4double lP=std::log(P);
266  return CrossSectionFormula(tZ, tN, P, lP);
267 }
268 
269 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
271 {
272  G4double P=std::exp(lP);
273  return CrossSectionFormula(tZ, tN, P, lP);
274 }
275 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
277  G4double P, G4double lP)
278 {
279  G4double sigma=0.;
280  if(tZ==1 && !tN) // PiPlus-Proton interaction from G4QuasiElRatios
281  {
282  G4double ld=lP-3.5;
283  G4double ld2=ld*ld;
284  G4double p2=P*P;
285  G4double p4=p2*p2;
286  G4double sp=std::sqrt(P);
287  G4double lm=lP-.32;
288  G4double md=lm*lm+.04;
289  G4double El=(.0557*ld2+2.4+6./sp)/(1.+3./p4);
290  G4double To=(.3*ld2+22.3+5./sp)/(1.+1./p4);
291  sigma=(To-El)+.1/md;
292  }
293  else if(tZ==1 && tN==1) // pimp_tot
294  {
295  G4double p2=P*P;
296  G4double d=lP-2.7;
297  G4double f=lP+1.25;
298  G4double gg=lP-.017;
299  sigma=(.55*d*d+38.+23./std::sqrt(P))/(1.+.3/p2/p2)+18./(f*f+.1089)+.02/(gg*gg+.0025);
300  }
301  else if(tZ<97 && tN<152) // General solution
302  {
303  G4double d=lP-4.2;
304  G4double p2=P*P;
305  G4double p4=p2*p2;
306  G4double a=tN+tZ; // A of the target
307  G4double al=std::log(a);
308  G4double sa=std::sqrt(a);
309  G4double ssa=std::sqrt(sa);
310  G4double a2=a*a;
311  G4double c=41.*std::exp(al*.68)*(1.+44./a2)/(1.+8./a)/(1.+200./a2/a2);
312  G4double f=290.*ssa/(1.+34./a/ssa);
313  G4double gg=-1.32-al*.043;
314  G4double u=lP-gg;
315  G4double h=al*(.4-.055*al);
316  G4double r=.01+a2*5.E-8;
317  sigma=(c+d*d)/(1.+(.2-.009*sa)/p4)+f/(u*u+h*h)/(1.+r/p2);
318  }
319  else
320  {
321  G4cerr<<"-Warning-G4ChipsPiPlusNuclearCroSect::CSForm:*Bad A* Z="<<tZ<<", N="<<tN<<G4endl;
322  sigma=0.;
323  }
324  if(sigma<0.) return 0.;
325  return sigma;
326 }
327 
329 {
330  if(DX<=0. || N<2)
331  {
332  G4cerr<<"***G4ChipsPionPlusInelasticXS::EquLinearFit: DX="<<DX<<", N="<<N<<G4endl;
333  return Y[0];
334  }
335 
336  G4int N2=N-2;
337  G4double d=(X-X0)/DX;
338  G4int jj=static_cast<int>(d);
339  if (jj<0) jj=0;
340  else if(jj>N2) jj=N2;
341  d-=jj; // excess
342  G4double yi=Y[jj];
343  G4double sigma=yi+(Y[jj+1]-yi)*d;
344 
345  return sigma;
346 }
static G4PionPlus * Definition()
Definition: G4PionPlus.cc:52
G4double CrossSectionFormula(G4int targZ, G4int targN, G4double P, G4double lP)
static const G4int nH
G4double a
Definition: TRTMaterials.hh:39
static const G4int nL
G4double EquLinearFit(G4double X, G4int N, G4double X0, G4double DX, G4double *Y)
int G4int
Definition: G4Types.hh:78
virtual G4bool IsIsoApplicable(const G4DynamicParticle *Pt, G4int Z, G4int A, const G4Element *elm, const G4Material *mat)
G4_DECLARE_XS_FACTORY(G4ChipsPionPlusInelasticXS)
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
G4double GetTotalMomentum() const
static const G4double dE
G4double CalculateCrossSection(G4int F, G4int I, G4int PDG, G4int Z, G4int N, G4double Momentum)
bool G4bool
Definition: G4Types.hh:79
G4double CrossSectionLin(G4int targZ, G4int targN, G4double P)
const G4double p2
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
const G4int n
static const G4double A[nN]
G4double ThresholdMomentum(G4int targZ, G4int targN)
std::vector< G4double * > * LEN
std::vector< G4double * > * HEN
G4double GetPDGMass() const
static const double millibarn
Definition: G4SIunits.hh:96
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int tgZ, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double CrossSectionLog(G4int targZ, G4int targN, G4double lP)
static const G4double THmin
static const G4double a2
G4GLOB_DLL std::ostream G4cerr