Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ChipsAntiBaryonInelasticXS.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: G4ChipsAntiBaryonInelasticXS 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 // anti-baryoninteractions. Original author: M. Kossov
38 // -------------------------------------------------------------------------------------
39 //
40 
42 #include "G4SystemOfUnits.hh"
43 #include "G4DynamicParticle.hh"
44 #include "G4ParticleDefinition.hh"
45 #include "G4AntiNeutron.hh"
46 #include "G4AntiProton.hh"
47 #include "G4AntiLambda.hh"
48 #include "G4AntiSigmaPlus.hh"
49 #include "G4AntiSigmaMinus.hh"
50 #include "G4AntiSigmaZero.hh"
51 #include "G4AntiXiMinus.hh"
52 #include "G4AntiXiZero.hh"
53 #include "G4AntiOmegaMinus.hh"
54 
55 // factory
56 #include "G4CrossSectionFactory.hh"
57 //
59 
61 {
62  lastLEN=0; // Pointer to lastArray of LowEn CS
63  lastHEN=0; // Pointer to lastArray of HighEn CS
64  lastN=0; // The last N of calculated nucleus
65  lastZ=0; // The last Z of calculated nucleus
66  lastP=0.; // Last used Cross Section Momentum
67  lastTH=0.; // Last threshold momentum
68  lastCS=0.; // Last value of the Cross Section
69  lastI=0; // The last position in the DAMDB
70  LEN = new std::vector<G4double*>;
71  HEN = new std::vector<G4double*>;
72 }
73 
75 {
76  G4int lens=LEN->size();
77  for(G4int i=0; i<lens; ++i) delete[] (*LEN)[i];
78  delete LEN;
79  G4int hens=HEN->size();
80  for(G4int i=0; i<hens; ++i) delete[] (*HEN)[i];
81  delete HEN;
82 }
83 
85  const G4Element*,
86  const G4Material*)
87 {
88  G4ParticleDefinition* particle = Pt->GetDefinition();
89 
90  if(particle == G4AntiNeutron::AntiNeutron())
91  {
92  return true;
93  }
94  else if(particle == G4AntiProton::AntiProton())
95  {
96  return true;
97  }
98  else if(particle == G4AntiLambda::AntiLambda())
99  {
100  return true;
101  }
102  else if(particle == G4AntiSigmaPlus::AntiSigmaPlus())
103  {
104  return true;
105  }
106  else if(particle == G4AntiSigmaMinus::AntiSigmaMinus())
107  {
108  return true;
109  }
110  else if(particle == G4AntiSigmaZero::AntiSigmaZero())
111  {
112  return true;
113  }
114  else if(particle == G4AntiXiMinus::AntiXiMinus())
115  {
116  return true;
117  }
118  else if(particle == G4AntiXiZero::AntiXiZero())
119  {
120  return true;
121  }
122  else if(particle == G4AntiOmegaMinus::AntiOmegaMinus())
123  {
124  return true;
125  }
126  return false;
127 }
128 
129 // The main member function giving the collision cross section (P is in IU, CS is in mb)
130 // Make pMom in independent units ! (Now it is MeV)
132  const G4Isotope*,
133  const G4Element*,
134  const G4Material*)
135 {
136  G4double pMom=Pt->GetTotalMomentum();
137  G4int tgN = A - tgZ;
138  G4int pdg = Pt->GetDefinition()->GetPDGEncoding();
139 
140  return GetChipsCrossSection(pMom, tgZ, tgN, pdg);
141 }
142 
144 {
145  static G4int j; // A#0f Z/N-records already tested in AMDB
146  static std::vector <G4int> colN; // Vector of N for calculated nuclei (isotops)
147  static std::vector <G4int> colZ; // Vector of Z for calculated nuclei (isotops)
148  static std::vector <G4double> colP; // Vector of last momenta for the reaction
149  static std::vector <G4double> colTH; // Vector of energy thresholds for the reaction
150  static std::vector <G4double> colCS; // Vector of last cross sections for the reaction
151  // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
152 
153  G4bool in=false; // By default the isotope must be found in the AMDB
154  if(tgN!=lastN || tgZ!=lastZ) // The nucleus was not the last used isotope
155  {
156  in = false; // By default the isotope haven't be found in AMDB
157  lastP = 0.; // New momentum history (nothing to compare with)
158  lastN = tgN; // The last N of the calculated nucleus
159  lastZ = tgZ; // The last Z of the calculated nucleus
160  lastI = colN.size(); // Size of the Associative Memory DB in the heap
161  j = 0; // A#0f records found in DB for this projectile
162  if(lastI) for(G4int i=0; i<lastI; i++) // AMDB exists, try to find the (Z,N) isotope
163  {
164  if(colN[i]==tgN && colZ[i]==tgZ) // Try the record "i" in the AMDB
165  {
166  lastI=i; // Remember the index for future fast/last use
167  lastTH =colTH[i]; // The last THreshold (A-dependent)
168  if(pMom<=lastTH)
169  {
170  return 0.; // Energy is below the Threshold value
171  }
172  lastP =colP [i]; // Last Momentum (A-dependent)
173  lastCS =colCS[i]; // Last CrossSect (A-dependent)
174  in = true; // This is the case when the isotop is found in DB
175  // Momentum pMom is in IU ! @@ Units
176  lastCS=CalculateCrossSection(-1,j,cPDG,lastZ,lastN,pMom); // read & update
177  if(lastCS<=0. && pMom>lastTH) // Correct the threshold (@@ No intermediate Zeros)
178  {
179  lastCS=0.;
180  lastTH=pMom;
181  }
182  break; // Go out of the LOOP
183  }
184  j++; // Increment a#0f records found in DB
185  }
186  if(!in) // This isotope has not been calculated previously
187  {
189  lastCS=CalculateCrossSection(0,j,cPDG,lastZ,lastN,pMom); //calculate & create
190  //if(lastCS>0.) // It means that the AMBD was initialized
191  //{
192 
193  lastTH = 0; //ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
194  colN.push_back(tgN);
195  colZ.push_back(tgZ);
196  colP.push_back(pMom);
197  colTH.push_back(lastTH);
198  colCS.push_back(lastCS);
199  //} // M.K. Presence of H1 with high threshold breaks the syncronization
200  return lastCS*millibarn;
201  } // End of creation of the new set of parameters
202  else
203  {
204  colP[lastI]=pMom;
205  colCS[lastI]=lastCS;
206  }
207  } // End of parameters udate
208  else if(pMom<=lastTH)
209  {
210  return 0.; // Momentum is below the Threshold Value -> CS=0
211  }
212  else // It is the last used -> use the current tables
213  {
214  lastCS=CalculateCrossSection(1,j,cPDG,lastZ,lastN,pMom); // Only read and UpdateDB
215  lastP=pMom;
216  }
217  return lastCS*millibarn;
218 }
219 
220 // The main member function giving the gamma-A cross section (E in GeV, CS in mb)
221 G4double G4ChipsAntiBaryonInelasticXS::CalculateCrossSection(G4int F, G4int I,
222  G4int, G4int targZ, G4int targN, G4double Momentum)
223 {
224  static const G4double THmin=27.; // default minimum Momentum (MeV/c) Threshold
225  static const G4double THmiG=THmin*.001; // minimum Momentum (GeV/c) Threshold
226  static const G4double dP=10.; // step for the LEN (Low ENergy) table MeV/c
227  static const G4double dPG=dP*.001; // step for the LEN (Low ENergy) table GeV/c
228  static const G4int nL=105; // A#of LEN points in E (step 10 MeV/c)
229  static const G4double Pmin=THmin+(nL-1)*dP; // minP for the HighE part with safety
230  static const G4double Pmax=227000.; // maxP for the HEN (High ENergy) part 227 GeV
231  static const G4int nH=224; // A#of HEN points in lnE
232  static const G4double milP=std::log(Pmin);// Low logarithm energy for the HEN part
233  static const G4double malP=std::log(Pmax);// High logarithm energy (each 2.75 percent)
234  static const G4double dlP=(malP-milP)/(nH-1); // Step in log energy in the HEN part
235  static const G4double milPG=std::log(.001*Pmin);// Low logarithmEnergy for HEN part GeV/c
236  G4double sigma=0.;
237  if(F&&I) sigma=0.; // @@ *!* Fake line *!* to use F & I !!!Temporary!!!
238  //G4double A=targN+targZ; // A of the target
239  if(F<=0) // This isotope was not the last used isotop
240  {
241  if(F<0) // This isotope was found in DAMDB =-----=> RETRIEVE
242  {
243  G4int sync=LEN->size();
244  if(sync<=I) G4cerr<<"*!*G4QPiMinusNuclCS::CalcCrosSect:Sync="<<sync<<"<="<<I<<G4endl;
245  lastLEN=(*LEN)[I]; // Pointer to prepared LowEnergy cross sections
246  lastHEN=(*HEN)[I]; // Pointer to prepared High Energy cross sections
247  }
248  else // This isotope wasn't calculated before => CREATE
249  {
250  lastLEN = new G4double[nL]; // Allocate memory for the new LEN cross sections
251  lastHEN = new G4double[nH]; // Allocate memory for the new HEN cross sections
252  // --- Instead of making a separate function ---
253  G4double P=THmiG; // Table threshold in GeV/c
254  for(G4int k=0; k<nL; k++)
255  {
256  lastLEN[k] = CrossSectionLin(targZ, targN, P);
257  P+=dPG;
258  }
259  G4double lP=milPG;
260  for(G4int n=0; n<nH; n++)
261  {
262  lastHEN[n] = CrossSectionLog(targZ, targN, lP);
263  lP+=dlP;
264  }
265  // --- End of possible separate function
266  // *** The synchronization check ***
267  G4int sync=LEN->size();
268  if(sync!=I)
269  {
270  G4cerr<<"***G4QPiMinusNuclCS::CalcCrossSect: Sinc="<<sync<<"#"<<I<<", Z=" <<targZ
271  <<", N="<<targN<<", F="<<F<<G4endl;
272  //G4Exception("G4PiMinusNuclearCS::CalculateCS:","39",FatalException,"DBoverflow");
273  }
274  LEN->push_back(lastLEN); // remember the Low Energy Table
275  HEN->push_back(lastHEN); // remember the High Energy Table
276  } // End of creation of the new set of parameters
277  } // End of parameters udate
278  // =-------------------= NOW the Magic Formula =--------------------=
279  if (Momentum<lastTH) return 0.; // It must be already checked in the interface class
280  else if (Momentum<Pmin) // High Energy region
281  {
282  sigma=EquLinearFit(Momentum,nL,THmin,dP,lastLEN);
283  }
284  else if (Momentum<Pmax) // High Energy region
285  {
286  G4double lP=std::log(Momentum);
287  sigma=EquLinearFit(lP,nH,milP,dlP,lastHEN);
288  }
289  else // UHE region (calculation, not frequent)
290  {
291  G4double P=0.001*Momentum; // Approximation formula is for P in GeV/c
292  sigma=CrossSectionFormula(targZ, targN, P, std::log(P));
293  }
294  if(sigma<0.) return 0.;
295  return sigma;
296 }
297 
298 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) (P in GeV/c)
299 G4double G4ChipsAntiBaryonInelasticXS::CrossSectionLin(G4int tZ, G4int tN, G4double P)
300 {
301  G4double lP=std::log(P);
302  return CrossSectionFormula(tZ, tN, P, lP);
303 }
304 
305 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
306 G4double G4ChipsAntiBaryonInelasticXS::CrossSectionLog(G4int tZ, G4int tN, G4double lP)
307 {
308  G4double P=std::exp(lP);
309  return CrossSectionFormula(tZ, tN, P, lP);
310 }
311 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
312 G4double G4ChipsAntiBaryonInelasticXS::CrossSectionFormula(G4int tZ, G4int tN,
313  G4double P, G4double lP)
314 {
315  G4double sigma=0.;
316  if(tZ==1 && !tN) // AntiBar-Prot interaction from G4QuasiElRatios
317  {
318  G4double ld=lP-3.5;
319  G4double ld2=ld*ld;
320  G4double ye=std::exp(lP*1.25);
321  G4double yt=std::exp(lP*0.35);
322  G4double El=80./(ye+1.);
323  G4double To=(80./yt+.3)/yt;
324  sigma=(To-El)+.2443*ld2+31.48;
325  }
326  else if(tZ==1 && tN==1)
327  {
328  G4double r=lP-3.7;
329  sigma=0.6*r*r+67.+90.*std::exp(-lP*.666);
330  }
331  else if(tZ<97 && tN<152) // General solution
332  {
333  G4double d=lP-4.2;
334  G4double sp=std::sqrt(P);
335  G4double a=tN+tZ; // A of the target
336  G4double sa=std::sqrt(a);
337  G4double a2=a*a;
338  G4double a3=a2*a;
339  G4double a2s=a2*sa;
340  G4double c=(170.+3600./a2s)/(1.+65./a2s)+40.*std::pow(a,0.712)/(1.+12.2/a)/(1.+34./a2);
341  G4double r=(170.+0.01*a3)/(1.+a3/28000.);
342  sigma=c+d*d+r/sp;
343  }
344  else
345  {
346  G4cerr<<"-Warning-G4QAntiBarNuclearCroSect::CSForm:*Bad A* Z="<<tZ<<", N="<<tN<<G4endl;
347  sigma=0.;
348  }
349  if(sigma<0.) return 0.;
350  return sigma;
351 }
352 
353 G4double G4ChipsAntiBaryonInelasticXS::EquLinearFit(G4double X, G4int N, G4double X0, G4double DX, G4double* Y)
354 {
355  if(DX<=0. || N<2)
356  {
357  G4cerr<<"***G4ChipsAntiBaryonInelasticXS::EquLinearFit: DX="<<DX<<", N="<<N<<G4endl;
358  return Y[0];
359  }
360 
361  G4int N2=N-2;
362  G4double d=(X-X0)/DX;
363  G4int j=static_cast<int>(d);
364  if (j<0) j=0;
365  else if(j>N2) j=N2;
366  d-=j; // excess
367  G4double yi=Y[j];
368  G4double sigma=yi+(Y[j+1]-yi)*d;
369 
370  return sigma;
371 }