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