Geant4  10.01.p03
G4ChipsKaonMinusInelasticXS.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: G4ChipsKaonMinusInelasticXS 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 // Short description: Cross-sections extracted from the CHIPS package for
35 // kaon(minus)-nuclear interactions. Author: M. Kossov
36 // -------------------------------------------------------------------------------------
37 //
38 
40 #include "G4SystemOfUnits.hh"
41 #include "G4DynamicParticle.hh"
42 #include "G4ParticleDefinition.hh"
43 #include "G4KaonMinus.hh"
44 
45 // factory
46 #include "G4CrossSectionFactory.hh"
47 //
49 
50 namespace {
51  const G4double THmin=27.; // default minimum Momentum (MeV/c) Threshold
52  const G4double THmiG=THmin*.001; // minimum Momentum (GeV/c) Threshold
53  const G4double dP=10.; // step for the LEN (Low ENergy) table MeV/c
54  const G4double dPG=dP*.001; // step for the LEN (Low ENergy) table GeV/c
55  const G4int nL=105; // A#of LEN points in E (step 10 MeV/c)
56  const G4double Pmin=THmin+(nL-1)*dP; // minP for the HighE part with safety
57  const G4double Pmax=227000.; // maxP for the HEN (High ENergy) part 227 GeV
58  const G4int nH=224; // A#of HEN points in lnE
59  const G4double milP=std::log(Pmin);// Low logarithm energy for the HEN part
60  const G4double malP=std::log(Pmax);// High logarithm energy (each 2.75 percent)
61  const G4double dlP=(malP-milP)/(nH-1); // Step in log energy in the HEN part
62  const G4double milPG=std::log(.001*Pmin);// Low logarithmEnergy for HEN part GeV/c
63 }
64 // Initialization of the
65 
67 {
68  lastLEN=0; // Pointer to lastArray of LowEn CS
69  lastHEN=0; // Pointer to lastArray of HighEn CS
70  lastN=0; // The last N of calculated nucleus
71  lastZ=0; // The last Z of calculated nucleus
72  lastP=0.; // Last used in CrossSection Momentum
73  lastTH=0.; // Last threshold momentum
74  lastCS=0.; // Last value of the Cross Section
75  lastI=0; // The last position in the DAMDB
76  LEN = new std::vector<G4double*>;
77  HEN = new std::vector<G4double*>;
78 }
79 
81 {
82  G4int lens=LEN->size();
83  for(G4int i=0; i<lens; ++i) delete[] (*LEN)[i];
84  delete LEN;
85 
86  G4int hens=HEN->size();
87  for(G4int i=0; i<hens; ++i) delete[] (*HEN)[i];
88  delete HEN;
89 }
90 
92  const G4Element*,
93  const G4Material*)
94 {
95  return true;
96 }
97 
98 
99 // The main member function giving the collision cross section (P is in IU, CS is in mb)
100 // Make pMom in independent units ! (Now it is MeV)
102  const G4Isotope*,
103  const G4Element*,
104  const G4Material*)
105 {
106  G4double pMom=Pt->GetTotalMomentum();
107  G4int tgN = A - tgZ;
108 
109  return GetChipsCrossSection(pMom, tgZ, tgN, -321);
110 }
111 
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 be 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,-321,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,-321,lastZ,lastN,pMom); //calculate & create
151  //if(lastCS>0.) // It means that the AMBD was initialized
152  //{
153 
154  // lastTH = ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
155 
156  lastTH = 0; // WP - to be checked!!!
157  colN.push_back(tgN);
158  colZ.push_back(tgZ);
159  colP.push_back(pMom);
160  colTH.push_back(lastTH);
161  colCS.push_back(lastCS);
162  //} // M.K. Presence of H1 with high threshold breaks the syncronization
163  return lastCS*millibarn;
164  } // End of creation of the new set of parameters
165  else
166  {
167  colP[lastI]=pMom;
168  colCS[lastI]=lastCS;
169  }
170  } // End of parameters udate
171  else if(pMom<=lastTH)
172  {
173  return 0.; // Momentum is below the Threshold Value -> CS=0
174  }
175  else // It is the last used -> use the current tables
176  {
177  lastCS=CalculateCrossSection(1,j,-321,lastZ,lastN,pMom); // Only read and UpdateDB
178  lastP=pMom;
179  }
180  return lastCS*millibarn;
181 }
182 
183 // The main member function giving the gamma-A cross section (E in GeV, CS in mb)
185  G4int, G4int targZ, G4int targN, G4double Momentum)
186 {
187  G4double sigma=0.;
188  if(F&&I) sigma=0.; // @@ *!* Fake line *!* to use F & I !!!Temporary!!!
189  //G4double A=targN+targZ; // A of the target
190  if(F<=0) // This isotope was not the last used isotop
191  {
192  if(F<0) // This isotope was found in DAMDB =-----=> RETRIEVE
193  {
194  G4int sync=LEN->size();
195  if(sync<=I) G4cerr<<"*!*G4QPiMinusNuclCS::CalcCrosSect:Sync="<<sync<<"<="<<I<<G4endl;
196  lastLEN=(*LEN)[I]; // Pointer to prepared LowEnergy cross sections
197  lastHEN=(*HEN)[I]; // Pointer to prepared High Energy cross sections
198  }
199  else // This isotope wasn't calculated before => CREATE
200  {
201  lastLEN = new G4double[nL]; // Allocate memory for the new LEN cross sections
202  lastHEN = new G4double[nH]; // Allocate memory for the new HEN cross sections
203  // --- Instead of making a separate function ---
204  G4double P=THmiG; // Table threshold in GeV/c
205  for(G4int k=0; k<nL; k++)
206  {
207  lastLEN[k] = CrossSectionLin(targZ, targN, P);
208  P+=dPG;
209  }
210  G4double lP=milPG;
211  for(G4int n=0; n<nH; n++)
212  {
213  lastHEN[n] = CrossSectionLog(targZ, targN, lP);
214  lP+=dlP;
215  }
216  // --- End of possible separate function
217  // *** The synchronization check ***
218  G4int sync=LEN->size();
219  if(sync!=I)
220  {
221  G4cerr<<"***G4ChipsKaonMinusCS::CalcCrossSect: Sinc="<<sync<<"#"<<I<<", Z=" <<targZ
222  <<", N="<<targN<<", F="<<F<<G4endl;
223  //G4Exception("G4PiMinusNuclearCS::CalculateCS:","39",FatalException,"DBoverflow");
224  }
225  LEN->push_back(lastLEN); // remember the Low Energy Table
226  HEN->push_back(lastHEN); // remember the High Energy Table
227  } // End of creation of the new set of parameters
228  } // End of parameters udate
229  // =------------------= NOW the Magic Formula =--------------------------=
230  if (Momentum<lastTH) return 0.; // It must be already checked in the interface class
231  else if (Momentum<Pmin) // High Energy region
232  {
233  sigma=EquLinearFit(Momentum,nL,THmin,dP,lastLEN);
234  }
235  else if (Momentum<Pmax) // High Energy region
236  {
237  G4double lP=std::log(Momentum);
238  sigma=EquLinearFit(lP,nH,milP,dlP,lastHEN);
239  }
240  else // UHE region (calculation, not frequent)
241  {
242  G4double P=0.001*Momentum; // Approximation formula is for P in GeV/c
243  sigma=CrossSectionFormula(targZ, targN, P, std::log(P));
244  }
245  if(sigma<0.) return 0.;
246  return sigma;
247 }
248 
249 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) (P in GeV/c)
251 {
252  G4double lP=std::log(P);
253  return CrossSectionFormula(tZ, tN, P, lP);
254 }
255 
256 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
258 {
259  G4double P=std::exp(lP);
260  return CrossSectionFormula(tZ, tN, P, lP);
261 }
262 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
264  G4double P, G4double lP)
265 {
266  G4double sigma=0.;
267  if(tZ==1 && !tN) // PiMin-Proton interaction from G4QuasiElRatios
268  {
269  G4double ld=lP-3.5;
270  G4double ld2=ld*ld;
271  G4double p2=P*P;
272  G4double p4=p2*p2;
273  G4double sp=std::sqrt(P);
274  G4double psp=P*sp;
275  G4double lm=P-.39;
276  G4double md=lm*lm+.000156;
277  G4double lh=P-1.;
278  G4double hd=lh*lh+.0156;
279  G4double El=(.0557*ld2+2.23)/(1.-.7/sp+.075/p4);
280  G4double To=(.3*ld2+19.5)/(1.-.21/sp+.52/p4);
281  sigma=8.8/psp+(To-El)+.002/md+.15/hd;
282  }
283  else if(tZ==1 && tN==1) // kmp_tot
284  {
285  G4double p2=P*P;
286  G4double dX=lP-3.7;
287  G4double dR=P-.94;
288  G4double sp=std::sqrt(P);
289  sigma=(.6*dX*dX+36.)/(1.-.11/sp+.52/p2/p2)+.7/(dR*dR+.0256)+18./P/sp;
290  }
291  else if(tZ<97 && tN<152) // General solution
292  {
293  G4double d=lP-4.2;
294  G4double sp=std::sqrt(P);
295  G4double p2=P*P;
296  G4double a=tN+tZ; // A of the target
297  G4double sa=std::sqrt(a);
298  G4double al=std::log(a);
299  G4double a2=a*a;
300  G4double c=52.*std::exp(al*0.6)*(1.+97./a2)/(1.+9.8/a)/(1.+47./a2);
301  G4double gg=-.2-.003*a;
302  G4double h=.5+.07*a;
303  G4double v=P-1.;
304  G4double f=.6*a*sa/(1.+.00002*a2);
305  G4double u=.125+.127*al;
306  sigma=(c+d*d)/(1.+gg/sp+h/p2/p2)+f/(v*v+u*u)+20.*sa/P/sp;
307  }
308  else
309  {
310  G4cerr<<"-Warning-G4ChipsKMinusNuclearCroSect::CSForm:*Bad A* Z="<<tZ<<", N="<<tN<<G4endl;
311  sigma=0.;
312  }
313  if(sigma<0.) return 0.;
314  return sigma;
315 }
316 
318 {
319  if(DX<=0. || N<2)
320  {
321  G4cerr<<"***G4ChipsKaonMinusInelasticXS::EquLinearFit: DX="<<DX<<", N="<<N<<G4endl;
322  return Y[0];
323  }
324 
325  G4int N2=N-2;
326  G4double d=(X-X0)/DX;
327  G4int jj=static_cast<int>(d);
328  if (jj<0) jj=0;
329  else if(jj>N2) jj=N2;
330  d-=jj; // excess
331  G4double yi=Y[jj];
332  G4double sigma=yi+(Y[jj+1]-yi)*d;
333 
334  return sigma;
335 }
static const G4int nH
G4double a
Definition: TRTMaterials.hh:39
static const G4int nL
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)
G4double CrossSectionLin(G4int targZ, G4int targN, G4double P)
bool G4bool
Definition: G4Types.hh:79
const G4double p2
G4double CalculateCrossSection(G4int F, G4int I, G4int PDG, G4int Z, G4int N, G4double Momentum)
const G4int n
static const G4double A[nN]
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
G4double CrossSectionFormula(G4int targZ, G4int targN, G4double P, G4double lP)
G4double CrossSectionLog(G4int targZ, G4int targN, G4double lP)
G4_DECLARE_XS_FACTORY(G4ChipsKaonMinusInelasticXS)
static const double millibarn
Definition: G4SIunits.hh:96
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int tgZ, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
static const G4double THmin
static const G4double a2
G4GLOB_DLL std::ostream G4cerr
G4double EquLinearFit(G4double X, G4int N, G4double X0, G4double DX, G4double *Y)