Geant4  10.00.p01
G4ChipsKaonPlusInelasticXS.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: G4QKaonPlusNuclearCrossSection 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 from the CHIPS package for
36 // kaon(minus)-nuclear interactions. Author: M. Kossov
37 // -------------------------------------------------------------------------------------
38 //
39 
41 #include "G4SystemOfUnits.hh"
42 #include "G4DynamicParticle.hh"
43 #include "G4ParticleDefinition.hh"
44 #include "G4KaonPlus.hh"
45 #include "G4Proton.hh"
46 #include "G4PionPlus.hh"
47 #include "G4AutoLock.hh"
48 
49 // factory
50 #include "G4CrossSectionFactory.hh"
51 //
53 
54 namespace {
55  const G4double THmin=27.; // default minimum Momentum (MeV/c) Threshold
56  const G4double THmiG=THmin*.001; // minimum Momentum (GeV/c) Threshold
57  const G4double dP=10.; // step for the LEN (Low ENergy) table MeV/c
58  const G4double dPG=dP*.001; // step for the LEN (Low ENergy) table GeV/c
59  const G4int nL=105; // A#of LEN points in E (step 10 MeV/c)
60  const G4double Pmin=THmin+(nL-1)*dP; // minP for the HighE part with safety
61  const G4double Pmax=227000.; // maxP for the HEN (High ENergy) part 227 GeV
62  const G4int nH=224; // A#of HEN points in lnE
63  const G4double milP=std::log(Pmin);// Low logarithm energy for the HEN part
64  const G4double malP=std::log(Pmax);// High logarithm energy (each 2.75 percent)
65  const G4double dlP=(malP-milP)/(nH-1); // Step in log energy in the HEN part
66  const G4double milPG=std::log(.001*Pmin);// Low logarithmEnergy for HEN part GeV/c
67  const G4double third=1./3.;
69  G4double prM;// = G4Proton::Proton()->GetPDGMass(); // Proton mass in MeV
70  G4double piM;// = G4PionPlus::PionPlus()->GetPDGMass()+.1; // Pion mass in MeV+Safety (WP)??
71  G4double pM;// = G4KaonPlus::KaonPlus()->GetPDGMass(); // Projectile mass in MeV
72  G4double tpM;//= pM+pM; // Doubled projectile mass (MeV)
73 }
74 
76 {
77  G4AutoLock l(&initM);
78  prM = G4Proton::Proton()->GetPDGMass(); // Proton mass in MeV
79  piM = G4PionPlus::PionPlus()->GetPDGMass()+.1; // Pion mass in MeV+Safety (WP)??
80  pM = G4KaonPlus::KaonPlus()->GetPDGMass(); // Projectile mass in MeV
81  tpM = pM+pM; // Doubled projectile mass (MeV)
82  l.unlock();
83  // Initialization of the
84  lastLEN=0; // Pointer to the lastArray of LowEn CS
85  lastHEN=0; // Pointer to the lastArray of HighEn CS
86  lastN=0; // The last N of calculated nucleus
87  lastZ=0; // The last Z of calculated nucleus
88  lastP=0.; // Last used in cross section Momentum
89  lastTH=0.; // Last threshold momentum
90  lastCS=0.; // Last value of the Cross Section
91  lastI=0; // The last position in the DAMDB
92  LEN = new std::vector<G4double*>;
93  HEN = new std::vector<G4double*>;
94 }
95 
97 {
98  G4int lens=LEN->size();
99  for(G4int i=0; i<lens; ++i) delete[] (*LEN)[i];
100  delete LEN;
101 
102  G4int hens=HEN->size();
103  for(G4int i=0; i<hens; ++i) delete[] (*HEN)[i];
104  delete HEN;
105 }
106 
107 
109  const G4Element*,
110  const G4Material*)
111 {
112  G4ParticleDefinition* particle = Pt->GetDefinition();
113  if (particle == G4KaonPlus::KaonPlus() ) return true;
114  return false;
115 }
116 
117 
118 // The main member function giving the collision cross section (P is in IU, CS is in mb)
119 // Make pMom in independent units ! (Now it is MeV)
121  const G4Isotope*,
122  const G4Element*,
123  const G4Material*)
124 {
125  G4double pMom=Pt->GetTotalMomentum();
126  G4int tgN = A - tgZ;
127 
128  return GetChipsCrossSection(pMom, tgZ, tgN, 321);
129 }
130 
132 {
133  static G4ThreadLocal G4int j; // A#0f Z/N-records already tested in AMDB
134  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)
135  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)
136  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
137  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
138  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
139  // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
140 
141  G4bool in=false; // By default the isotope must be found in the AMDB
142  if(tgN!=lastN || tgZ!=lastZ) // The nucleus was not the last used isotope
143  {
144  in = false; // By default the isotope haven't be found in AMDB
145  lastP = 0.; // New momentum history (nothing to compare with)
146  lastN = tgN; // The last N of the calculated nucleus
147  lastZ = tgZ; // The last Z of the calculated nucleus
148  lastI = colN.size(); // Size of the Associative Memory DB in the heap
149  j = 0; // A#0f records found in DB for this projectile
150 
151  if(lastI) for(G4int i=0; i<lastI; i++) // AMDB exists, try to find the (Z,N) isotope
152  {
153  if(colN[i]==tgN && colZ[i]==tgZ) // Try the record "i" in the AMDB
154  {
155  lastI=i; // Remember the index for future fast/last use
156  lastTH =colTH[i]; // The last THreshold (A-dependent)
157 
158  if(pMom<=lastTH)
159  {
160  return 0.; // Energy is below the Threshold value
161  }
162  lastP =colP [i]; // Last Momentum (A-dependent)
163  lastCS =colCS[i]; // Last CrossSect (A-dependent)
164  in = true; // This is the case when the isotop is found in DB
165  // Momentum pMom is in IU ! @@ Units
166  lastCS=CalculateCrossSection(-1,j,321,lastZ,lastN,pMom); // read & update
167 
168  if(lastCS<=0. && pMom>lastTH) // Correct the threshold (@@ No intermediate Zeros)
169  {
170  lastCS=0.;
171  lastTH=pMom;
172  }
173  break; // Go out of the LOOP
174  }
175  j++; // Increment a#0f records found in DB
176  }
177  if(!in) // This isotope has not been calculated previously
178  {
180  lastCS=CalculateCrossSection(0,j,321,lastZ,lastN,pMom); //calculate & create
181 
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,321,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  G4double sigma=0.;
217  if(F&&I) sigma=0.; // @@ *!* Fake line *!* to use F & I !!!Temporary!!!
218  G4double A=targN+targZ; // A of the target
219 
220  if(F<=0) // This isotope was not the last used isotop
221  {
222  if(F<0) // This isotope was found in DAMDB =-----=> RETRIEVE
223  {
224  G4int sync=LEN->size();
225  if(sync<=I) G4cerr<<"*!*G4ChipsKPlusNuclCS::CalcCrosSect:Sync="<<sync<<"<="<<I<<G4endl;
226  lastLEN=(*LEN)[I]; // Pointer to prepared LowEnergy cross sections
227  lastHEN=(*HEN)[I]; // Pointer to prepared High Energy cross sections
228  }
229  else // This isotope wasn't calculated before => CREATE
230  {
231  lastLEN = new G4double[nL]; // Allocate memory for the new LEN cross sections
232  lastHEN = new G4double[nH]; // Allocate memory for the new HEN cross sections
233  // --- Instead of making a separate function ---
234  G4double P=THmiG; // Table threshold in GeV/c
235  for(G4int k=0; k<nL; k++)
236  {
237  lastLEN[k] = CrossSectionLin(targZ, targN, P);
238  P+=dPG;
239  }
240  G4double lP=milPG;
241  for(G4int n=0; n<nH; n++)
242  {
243  lastHEN[n] = CrossSectionLog(targZ, targN, lP);
244  lP+=dlP;
245  }
246  // --- End of possible separate function
247  // *** The synchronization check ***
248  G4int sync=LEN->size();
249  if(sync!=I)
250  {
251  G4cerr<<"***G4ChipsKPlusNuclCS::CalcCrossSect: Sinc="<<sync<<"#"<<I<<", Z=" <<targZ
252  <<", N="<<targN<<", F="<<F<<G4endl;
253  //G4Exception("G4PiMinusNuclearCS::CalculateCS:","39",FatalException,"DBoverflow");
254  }
255  LEN->push_back(lastLEN); // remember the Low Energy Table
256  HEN->push_back(lastHEN); // remember the High Energy Table
257  } // End of creation of the new set of parameters
258  } // End of parameters udate
259  // =--------------------------= NOW the Magic Formula =---------------------------------=
260 
261  if (Momentum<lastTH) return 0.; // It must be already checked in the interface class
262  else if (Momentum<Pmin) // Low Energy region
263  {
264  if(A<=1. && Momentum < 600.) sigma=0.; // Approximation tot/el uncertainty
265  else sigma=EquLinearFit(Momentum,nL,THmin,dP,lastLEN);
266  }
267  else if (Momentum<Pmax) // High Energy region
268  {
269  G4double lP=std::log(Momentum);
270  sigma=EquLinearFit(lP,nH,milP,dlP,lastHEN);
271  }
272  else // UHE region (calculation, not frequent)
273  {
274  G4double P=0.001*Momentum; // Approximation formula is for P in GeV/c
275  sigma=CrossSectionFormula(targZ, targN, P, std::log(P));
276  }
277  if(sigma<0.) return 0.;
278  return sigma;
279 }
280 
281 // Electromagnetic momentum-threshold (in MeV/c)
283 {
284  G4double tA=tZ+tN;
285  if(tZ<.99 || tN<0.) return 0.;
286  G4double tM=931.5*tA;
287  G4double dE=piM; // At least one Pi0 must be created
288  if(tZ==1 && tN==0) tM=prM; // A threshold on the free proton
289  else dE=tZ/(1.+std::pow(tA,third)); // Safety for diffused edge of the nucleus (QE)
290  //G4double dE=1.263*tZ/(1.+std::pow(tA,third));
291  G4double T=dE+dE*(dE/2+pM)/tM;
292  return std::sqrt(T*(tpM+T));
293 }
294 
295 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) (P in GeV/c)
297 {
298  G4double lP=std::log(P);
299  return CrossSectionFormula(tZ, tN, P, lP);
300 }
301 
302 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
304 {
305  G4double P=std::exp(lP);
306  return CrossSectionFormula(tZ, tN, P, lP);
307 }
308 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
310  G4double P, G4double lP)
311 {
312  G4double sigma=0.;
313  if(tZ==1 && !tN) // KPlus-Proton interaction from G4QuasiElRatios
314  {
315  G4double ld=lP-3.5;
316  G4double ld2=ld*ld;
317  G4double sp=std::sqrt(P);
318  G4double p2=P*P;
319  G4double p4=p2*p2;
320  G4double lm=P-1.;
321  G4double md=lm*lm+.372;
322  G4double El=(.0557*ld2+2.23)/(1.-.7/sp+.1/p4);
323  G4double To=(.3*ld2+19.5)/(1.+.46/sp+1.6/p4);
324  sigma=(To-El)+.6/md;
325  }
326  else if(tZ<97 && tN<152) // General solution
327  {
328  G4double p2=P*P;
329  G4double p4=p2*p2;
330  G4double a=tN+tZ; // A of the target
331  G4double al=std::log(a);
332  G4double sa=std::sqrt(a);
333  G4double asa=a*sa;
334  G4double a2=a*a;
335  G4double a3=a2*a;
336  G4double a4=a2*a2;
337  G4double a8=a4*a4;
338  G4double a12=a8*a4;
339  G4double f=.6; // Default values for deutrons
340  G4double r=.5;
341  G4double gg=3.7;
342  G4double c=36.;
343  G4double ss=3.5;
344  G4double t=3.;
345  G4double u=.44;
346  G4double v=5.E-9;
347  if(tZ>1 && tN>1) // More than deuteron
348  {
349  f=1.;
350  r=1./(1.+.007*a2);
351  gg=4.2;
352  c=52.*std::exp(al*.6)*(1.+95./a2)/(1.+9./a)/(1.+46./a2);
353  ss=(40.+.14*a)/(1.+12./a);
354  G4double y=std::exp(al*1.7);
355  t=.185*y/(1.+.00012*y);
356  u=(1.+80./asa)/(1.+200./asa);
357  v=(1.+3.E-6*a4*(1.+6.E-7*a3+4.E10/a12))/a3/20000.;
358  }
359  G4double d=lP-gg;
360  G4double w=P-1.;
361  G4double rD=ss/(w*w+.36);
362  G4double h=P-.44;
363  G4double rR=t/(h*h+u*u);
364  sigma=(f*d*d+c)/(1.+r/std::sqrt(P)+1./p4)+(rD+rR)/(1+v/p4/p4);
365  }
366  else
367  {
368  G4cerr<<"-Warning-G4ChipsKaonPlusNuclearCroSect::CSForm:Bad A, Z="<<tZ<<", N="<<tN<<G4endl;
369  sigma=0.;
370  }
371  if(sigma<0.) return 0.;
372  return sigma;
373 }
374 
376 {
377  if(DX<=0. || N<2)
378  {
379  G4cerr<<"***G4ChipsKaonPlusInelasticXS::EquLinearFit: DX="<<DX<<", N="<<N<<G4endl;
380  return Y[0];
381  }
382 
383  G4int N2=N-2;
384  G4double d=(X-X0)/DX;
385  G4int j=static_cast<int>(d);
386  if (j<0) j=0;
387  else if(j>N2) j=N2;
388  d-=j; // excess
389  G4double yi=Y[j];
390  G4double sigma=yi+(Y[j+1]-yi)*d;
391 
392  return sigma;
393 }
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
static const G4int nH
G4double CrossSectionLog(G4int targZ, G4int targN, G4double lP)
G4double ThresholdMomentum(G4int targZ, G4int targN)
static const G4double a4
G4ParticleDefinition * GetDefinition() const
G4double a
Definition: TRTMaterials.hh:39
static const G4int nL
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int tgZ, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
std::vector< G4double * > * HEN
#define G4ThreadLocal
Definition: tls.hh:52
int G4int
Definition: G4Types.hh:78
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:158
std::vector< G4double * > * LEN
G4double GetTotalMomentum() const
static const G4double dE
G4double CalculateCrossSection(G4int F, G4int I, G4int PDG, G4int Z, G4int N, G4double Momentum)
G4_DECLARE_XS_FACTORY(G4ChipsKaonPlusInelasticXS)
bool G4bool
Definition: G4Types.hh:79
G4double CrossSectionLin(G4int targZ, G4int targN, G4double P)
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
const G4int n
G4double EquLinearFit(G4double X, G4int N, G4double X0, G4double DX, G4double *Y)
static const G4double A[nN]
static const G4double a3
G4int G4Mutex
Definition: G4Threading.hh:156
G4double GetPDGMass() const
G4double CrossSectionFormula(G4int targZ, G4int targN, G4double P, G4double lP)
static const double millibarn
Definition: G4SIunits.hh:96
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
static const G4double THmin
static const G4double a2
virtual G4bool IsIsoApplicable(const G4DynamicParticle *Pt, G4int Z, G4int A, const G4Element *elm, const G4Material *mat)
G4GLOB_DLL std::ostream G4cerr