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