Geant4_10
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  G4ParticleDefinition* particle = Pt->GetDefinition();
96  if (particle == G4KaonMinus::KaonMinus() ) return true;
97  return false;
98 }
99 
100 
101 // The main member function giving the collision cross section (P is in IU, CS is in mb)
102 // Make pMom in independent units ! (Now it is MeV)
104  const G4Isotope*,
105  const G4Element*,
106  const G4Material*)
107 {
108  G4double pMom=Pt->GetTotalMomentum();
109  G4int tgN = A - tgZ;
110 
111  return GetChipsCrossSection(pMom, tgZ, tgN, -321);
112 }
113 
115 {
116  static G4ThreadLocal G4int j; // A#0f Z/N-records already tested in AMDB
117  static G4ThreadLocal std::vector <G4int> *colN_G4MT_TLS_ = 0 ; if (!colN_G4MT_TLS_) colN_G4MT_TLS_ = new std::vector <G4int> ; std::vector <G4int> &colN = *colN_G4MT_TLS_; // Vector of N for calculated nuclei (isotops)
118  static G4ThreadLocal std::vector <G4int> *colZ_G4MT_TLS_ = 0 ; if (!colZ_G4MT_TLS_) colZ_G4MT_TLS_ = new std::vector <G4int> ; std::vector <G4int> &colZ = *colZ_G4MT_TLS_; // Vector of Z for calculated nuclei (isotops)
119  static G4ThreadLocal std::vector <G4double> *colP_G4MT_TLS_ = 0 ; if (!colP_G4MT_TLS_) colP_G4MT_TLS_ = new std::vector <G4double> ; std::vector <G4double> &colP = *colP_G4MT_TLS_; // Vector of last momenta for the reaction
120  static G4ThreadLocal std::vector <G4double> *colTH_G4MT_TLS_ = 0 ; if (!colTH_G4MT_TLS_) colTH_G4MT_TLS_ = new std::vector <G4double> ; std::vector <G4double> &colTH = *colTH_G4MT_TLS_; // Vector of energy thresholds for the reaction
121  static G4ThreadLocal std::vector <G4double> *colCS_G4MT_TLS_ = 0 ; if (!colCS_G4MT_TLS_) colCS_G4MT_TLS_ = new std::vector <G4double> ; std::vector <G4double> &colCS = *colCS_G4MT_TLS_; // Vector of last cross sections for the reaction
122  // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
123 
124  G4bool in=false; // By default the isotope must be found in the AMDB
125  if(tgN!=lastN || tgZ!=lastZ) // The nucleus was not the last used isotope
126  {
127  in = false; // By default the isotope haven't be found in AMDB
128  lastP = 0.; // New momentum history (nothing to compare with)
129  lastN = tgN; // The last N of the calculated nucleus
130  lastZ = tgZ; // The last Z of the calculated nucleus
131  lastI = colN.size(); // Size of the Associative Memory DB in the heap
132  j = 0; // A#0f records found in DB for this projectile
133  if(lastI) for(G4int i=0; i<lastI; i++) // AMDB exists, try to find the (Z,N) isotope
134  {
135  if(colN[i]==tgN && colZ[i]==tgZ) // Try the record "i" in the AMDB
136  {
137  lastI=i; // Remember the index for future fast/last use
138  lastTH =colTH[i]; // The last THreshold (A-dependent)
139  if(pMom<=lastTH)
140  {
141  return 0.; // Energy is below the Threshold value
142  }
143  lastP =colP [i]; // Last Momentum (A-dependent)
144  lastCS =colCS[i]; // Last CrossSect (A-dependent)
145  in = true; // This is the case when the isotop is found in DB
146  // Momentum pMom is in IU ! @@ Units
147  lastCS=CalculateCrossSection(-1,j,-321,lastZ,lastN,pMom); // read & update
148  if(lastCS<=0. && pMom>lastTH) // Correct the threshold (@@ No intermediate Zeros)
149  {
150  lastCS=0.;
151  lastTH=pMom;
152  }
153  break; // Go out of the LOOP
154  }
155  j++; // Increment a#0f records found in DB
156  }
157  if(!in) // This isotope has not been calculated previously
158  {
160  lastCS=CalculateCrossSection(0,j,-321,lastZ,lastN,pMom); //calculate & create
161  //if(lastCS>0.) // It means that the AMBD was initialized
162  //{
163 
164  // lastTH = ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
165 
166  lastTH = 0; // WP - to be checked!!!
167  colN.push_back(tgN);
168  colZ.push_back(tgZ);
169  colP.push_back(pMom);
170  colTH.push_back(lastTH);
171  colCS.push_back(lastCS);
172  //} // M.K. Presence of H1 with high threshold breaks the syncronization
173  return lastCS*millibarn;
174  } // End of creation of the new set of parameters
175  else
176  {
177  colP[lastI]=pMom;
178  colCS[lastI]=lastCS;
179  }
180  } // End of parameters udate
181  else if(pMom<=lastTH)
182  {
183  return 0.; // Momentum is below the Threshold Value -> CS=0
184  }
185  else // It is the last used -> use the current tables
186  {
187  lastCS=CalculateCrossSection(1,j,-321,lastZ,lastN,pMom); // Only read and UpdateDB
188  lastP=pMom;
189  }
190  return lastCS*millibarn;
191 }
192 
193 // The main member function giving the gamma-A cross section (E in GeV, CS in mb)
194 G4double G4ChipsKaonMinusInelasticXS::CalculateCrossSection(G4int F, G4int I,
195  G4int, G4int targZ, G4int targN, G4double Momentum)
196 {
197  G4double sigma=0.;
198  if(F&&I) sigma=0.; // @@ *!* Fake line *!* to use F & I !!!Temporary!!!
199  //G4double A=targN+targZ; // A of the target
200  if(F<=0) // This isotope was not the last used isotop
201  {
202  if(F<0) // This isotope was found in DAMDB =-----=> RETRIEVE
203  {
204  G4int sync=LEN->size();
205  if(sync<=I) G4cerr<<"*!*G4QPiMinusNuclCS::CalcCrosSect:Sync="<<sync<<"<="<<I<<G4endl;
206  lastLEN=(*LEN)[I]; // Pointer to prepared LowEnergy cross sections
207  lastHEN=(*HEN)[I]; // Pointer to prepared High Energy cross sections
208  }
209  else // This isotope wasn't calculated before => CREATE
210  {
211  lastLEN = new G4double[nL]; // Allocate memory for the new LEN cross sections
212  lastHEN = new G4double[nH]; // Allocate memory for the new HEN cross sections
213  // --- Instead of making a separate function ---
214  G4double P=THmiG; // Table threshold in GeV/c
215  for(G4int k=0; k<nL; k++)
216  {
217  lastLEN[k] = CrossSectionLin(targZ, targN, P);
218  P+=dPG;
219  }
220  G4double lP=milPG;
221  for(G4int n=0; n<nH; n++)
222  {
223  lastHEN[n] = CrossSectionLog(targZ, targN, lP);
224  lP+=dlP;
225  }
226  // --- End of possible separate function
227  // *** The synchronization check ***
228  G4int sync=LEN->size();
229  if(sync!=I)
230  {
231  G4cerr<<"***G4ChipsKaonMinusCS::CalcCrossSect: Sinc="<<sync<<"#"<<I<<", Z=" <<targZ
232  <<", N="<<targN<<", F="<<F<<G4endl;
233  //G4Exception("G4PiMinusNuclearCS::CalculateCS:","39",FatalException,"DBoverflow");
234  }
235  LEN->push_back(lastLEN); // remember the Low Energy Table
236  HEN->push_back(lastHEN); // remember the High Energy Table
237  } // End of creation of the new set of parameters
238  } // End of parameters udate
239  // =------------------= NOW the Magic Formula =--------------------------=
240  if (Momentum<lastTH) return 0.; // It must be already checked in the interface class
241  else if (Momentum<Pmin) // High Energy region
242  {
243  sigma=EquLinearFit(Momentum,nL,THmin,dP,lastLEN);
244  }
245  else if (Momentum<Pmax) // High Energy region
246  {
247  G4double lP=std::log(Momentum);
248  sigma=EquLinearFit(lP,nH,milP,dlP,lastHEN);
249  }
250  else // UHE region (calculation, not frequent)
251  {
252  G4double P=0.001*Momentum; // Approximation formula is for P in GeV/c
253  sigma=CrossSectionFormula(targZ, targN, P, std::log(P));
254  }
255  if(sigma<0.) return 0.;
256  return sigma;
257 }
258 
259 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) (P in GeV/c)
260 G4double G4ChipsKaonMinusInelasticXS::CrossSectionLin(G4int tZ, G4int tN, G4double P)
261 {
262  G4double lP=std::log(P);
263  return CrossSectionFormula(tZ, tN, P, lP);
264 }
265 
266 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
267 G4double G4ChipsKaonMinusInelasticXS::CrossSectionLog(G4int tZ, G4int tN, G4double lP)
268 {
269  G4double P=std::exp(lP);
270  return CrossSectionFormula(tZ, tN, P, lP);
271 }
272 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
273 G4double G4ChipsKaonMinusInelasticXS::CrossSectionFormula(G4int tZ, G4int tN,
274  G4double P, G4double lP)
275 {
276  G4double sigma=0.;
277  if(tZ==1 && !tN) // PiMin-Proton interaction from G4QuasiElRatios
278  {
279  G4double ld=lP-3.5;
280  G4double ld2=ld*ld;
281  G4double p2=P*P;
282  G4double p4=p2*p2;
283  G4double sp=std::sqrt(P);
284  G4double psp=P*sp;
285  G4double lm=P-.39;
286  G4double md=lm*lm+.000156;
287  G4double lh=P-1.;
288  G4double hd=lh*lh+.0156;
289  G4double El=(.0557*ld2+2.23)/(1.-.7/sp+.075/p4);
290  G4double To=(.3*ld2+19.5)/(1.-.21/sp+.52/p4);
291  sigma=8.8/psp+(To-El)+.002/md+.15/hd;
292  }
293  else if(tZ==1 && tN==1) // kmp_tot
294  {
295  G4double p2=P*P;
296  G4double dX=lP-3.7;
297  G4double dR=P-.94;
298  G4double sp=std::sqrt(P);
299  sigma=(.6*dX*dX+36.)/(1.-.11/sp+.52/p2/p2)+.7/(dR*dR+.0256)+18./P/sp;
300  }
301  else if(tZ<97 && tN<152) // General solution
302  {
303  G4double d=lP-4.2;
304  G4double sp=std::sqrt(P);
305  G4double p2=P*P;
306  G4double a=tN+tZ; // A of the target
307  G4double sa=std::sqrt(a);
308  G4double al=std::log(a);
309  G4double a2=a*a;
310  G4double c=52.*std::exp(al*0.6)*(1.+97./a2)/(1.+9.8/a)/(1.+47./a2);
311  G4double gg=-.2-.003*a;
312  G4double h=.5+.07*a;
313  G4double v=P-1.;
314  G4double f=.6*a*sa/(1.+.00002*a2);
315  G4double u=.125+.127*al;
316  sigma=(c+d*d)/(1.+gg/sp+h/p2/p2)+f/(v*v+u*u)+20.*sa/P/sp;
317  }
318  else
319  {
320  G4cerr<<"-Warning-G4ChipsKMinusNuclearCroSect::CSForm:*Bad A* Z="<<tZ<<", N="<<tN<<G4endl;
321  sigma=0.;
322  }
323  if(sigma<0.) return 0.;
324  return sigma;
325 }
326 
327 G4double G4ChipsKaonMinusInelasticXS::EquLinearFit(G4double X, G4int N, G4double X0, G4double DX, G4double* Y)
328 {
329  if(DX<=0. || N<2)
330  {
331  G4cerr<<"***G4ChipsKaonMinusInelasticXS::EquLinearFit: DX="<<DX<<", N="<<N<<G4endl;
332  return Y[0];
333  }
334 
335  G4int N2=N-2;
336  G4double d=(X-X0)/DX;
337  G4int j=static_cast<int>(d);
338  if (j<0) j=0;
339  else if(j>N2) j=N2;
340  d-=j; // excess
341  G4double yi=Y[j];
342  G4double sigma=yi+(Y[j+1]-yi)*d;
343 
344  return sigma;
345 }
tuple a
Definition: test.py:11
Float_t d
Definition: plot.C:237
ifstream in
Definition: comparison.C:7
Float_t Y
Definition: plot.C:39
G4ParticleDefinition * GetDefinition() const
#define G4ThreadLocal
Definition: tls.hh:52
int G4int
Definition: G4Types.hh:78
int millibarn
Definition: hepunit.py:40
TFile f
Definition: plotHisto.C:6
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
Float_t X
Definition: plot.C:39
G4double GetTotalMomentum() const
virtual G4bool IsIsoApplicable(const G4DynamicParticle *Pt, G4int Z, G4int A, const G4Element *elm, const G4Material *mat)
Char_t n[5]
bool G4bool
Definition: G4Types.hh:79
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
tuple v
Definition: test.py:18
#define G4_DECLARE_XS_FACTORY(cross_section)
#define G4endl
Definition: G4ios.hh:61
**D E S C R I P T I O N
Definition: HEPEvtcom.cc:77
double G4double
Definition: G4Types.hh:76
tuple c
Definition: test.py:13
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int tgZ, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4GLOB_DLL std::ostream G4cerr