Geant4  10.01
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  const 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 
134  G4bool in=false; // By default the isotope must be found in the AMDB
135  if(tgN!=lastN || tgZ!=lastZ) // The nucleus was not the last used isotope
136  {
137  in = false; // By default the isotope haven't be found in AMDB
138  lastP = 0.; // New momentum history (nothing to compare with)
139  lastN = tgN; // The last N of the calculated nucleus
140  lastZ = tgZ; // The last Z of the calculated nucleus
141  lastI = colN.size(); // Size of the Associative Memory DB in the heap
142  j = 0; // A#0f records found in DB for this projectile
143 
144  if(lastI) for(G4int i=0; i<lastI; i++) // AMDB exists, try to find the (Z,N) isotope
145  {
146  if(colN[i]==tgN && colZ[i]==tgZ) // Try the record "i" in the AMDB
147  {
148  lastI=i; // Remember the index for future fast/last use
149  lastTH =colTH[i]; // The last THreshold (A-dependent)
150 
151  if(pMom<=lastTH)
152  {
153  return 0.; // Energy is below the Threshold value
154  }
155  lastP =colP [i]; // Last Momentum (A-dependent)
156  lastCS =colCS[i]; // Last CrossSect (A-dependent)
157  in = true; // This is the case when the isotop is found in DB
158  // Momentum pMom is in IU ! @@ Units
159  lastCS=CalculateCrossSection(-1,j,321,lastZ,lastN,pMom); // read & update
160 
161  if(lastCS<=0. && pMom>lastTH) // Correct the threshold (@@ No intermediate Zeros)
162  {
163  lastCS=0.;
164  lastTH=pMom;
165  }
166  break; // Go out of the LOOP
167  }
168  j++; // Increment a#0f records found in DB
169  }
170  if(!in) // This isotope has not been calculated previously
171  {
173  lastCS=CalculateCrossSection(0,j,321,lastZ,lastN,pMom); //calculate & create
174 
175  //if(lastCS>0.) // It means that the AMBD was initialized
176  //{
177 
178  lastTH = 0; //ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
179  colN.push_back(tgN);
180  colZ.push_back(tgZ);
181  colP.push_back(pMom);
182  colTH.push_back(lastTH);
183  colCS.push_back(lastCS);
184  //} // M.K. Presence of H1 with high threshold breaks the syncronization
185  return lastCS*millibarn;
186  } // End of creation of the new set of parameters
187  else
188  {
189  colP[lastI]=pMom;
190  colCS[lastI]=lastCS;
191  }
192  } // End of parameters udate
193  else if(pMom<=lastTH)
194  {
195  return 0.; // Momentum is below the Threshold Value -> CS=0
196  }
197  else // It is the last used -> use the current tables
198  {
199  lastCS=CalculateCrossSection(1,j,321,lastZ,lastN,pMom); // Only read and UpdateDB
200  lastP=pMom;
201  }
202  return lastCS*millibarn;
203 }
204 
205 // The main member function giving the gamma-A cross section (E in GeV, CS in mb)
207  G4int, G4int targZ, G4int targN, G4double Momentum)
208 {
209  G4double sigma=0.;
210  if(F&&I) sigma=0.; // @@ *!* Fake line *!* to use F & I !!!Temporary!!!
211  G4double A=targN+targZ; // A of the target
212 
213  if(F<=0) // This isotope was not the last used isotop
214  {
215  if(F<0) // This isotope was found in DAMDB =-----=> RETRIEVE
216  {
217  G4int sync=LEN->size();
218  if(sync<=I) G4cerr<<"*!*G4ChipsKPlusNuclCS::CalcCrosSect:Sync="<<sync<<"<="<<I<<G4endl;
219  lastLEN=(*LEN)[I]; // Pointer to prepared LowEnergy cross sections
220  lastHEN=(*HEN)[I]; // Pointer to prepared High Energy cross sections
221  }
222  else // This isotope wasn't calculated before => CREATE
223  {
224  lastLEN = new G4double[nL]; // Allocate memory for the new LEN cross sections
225  lastHEN = new G4double[nH]; // Allocate memory for the new HEN cross sections
226  // --- Instead of making a separate function ---
227  G4double P=THmiG; // Table threshold in GeV/c
228  for(G4int k=0; k<nL; k++)
229  {
230  lastLEN[k] = CrossSectionLin(targZ, targN, P);
231  P+=dPG;
232  }
233  G4double lP=milPG;
234  for(G4int n=0; n<nH; n++)
235  {
236  lastHEN[n] = CrossSectionLog(targZ, targN, lP);
237  lP+=dlP;
238  }
239  // --- End of possible separate function
240  // *** The synchronization check ***
241  G4int sync=LEN->size();
242  if(sync!=I)
243  {
244  G4cerr<<"***G4ChipsKPlusNuclCS::CalcCrossSect: Sinc="<<sync<<"#"<<I<<", Z=" <<targZ
245  <<", N="<<targN<<", F="<<F<<G4endl;
246  //G4Exception("G4PiMinusNuclearCS::CalculateCS:","39",FatalException,"DBoverflow");
247  }
248  LEN->push_back(lastLEN); // remember the Low Energy Table
249  HEN->push_back(lastHEN); // remember the High Energy Table
250  } // End of creation of the new set of parameters
251  } // End of parameters udate
252  // =--------------------------= NOW the Magic Formula =---------------------------------=
253 
254  if (Momentum<lastTH) return 0.; // It must be already checked in the interface class
255  else if (Momentum<Pmin) // Low Energy region
256  {
257  if(A<=1. && Momentum < 600.) sigma=0.; // Approximation tot/el uncertainty
258  else sigma=EquLinearFit(Momentum,nL,THmin,dP,lastLEN);
259  }
260  else if (Momentum<Pmax) // High Energy region
261  {
262  G4double lP=std::log(Momentum);
263  sigma=EquLinearFit(lP,nH,milP,dlP,lastHEN);
264  }
265  else // UHE region (calculation, not frequent)
266  {
267  G4double P=0.001*Momentum; // Approximation formula is for P in GeV/c
268  sigma=CrossSectionFormula(targZ, targN, P, std::log(P));
269  }
270  if(sigma<0.) return 0.;
271  return sigma;
272 }
273 
274 // Electromagnetic momentum-threshold (in MeV/c)
276 {
277  G4double tA=tZ+tN;
278  if(tZ<.99 || tN<0.) return 0.;
279  G4double tM=931.5*tA;
280  G4double dE=piM; // At least one Pi0 must be created
281  if(tZ==1 && tN==0) tM=prM; // A threshold on the free proton
282  else dE=tZ/(1.+std::pow(tA,third)); // Safety for diffused edge of the nucleus (QE)
283  //G4double dE=1.263*tZ/(1.+std::pow(tA,third));
284  G4double T=dE+dE*(dE/2+pM)/tM;
285  return std::sqrt(T*(tpM+T));
286 }
287 
288 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) (P in GeV/c)
290 {
291  G4double lP=std::log(P);
292  return CrossSectionFormula(tZ, tN, P, lP);
293 }
294 
295 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
297 {
298  G4double P=std::exp(lP);
299  return CrossSectionFormula(tZ, tN, P, lP);
300 }
301 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
303  G4double P, G4double lP)
304 {
305  G4double sigma=0.;
306  if(tZ==1 && !tN) // KPlus-Proton interaction from G4QuasiElRatios
307  {
308  G4double ld=lP-3.5;
309  G4double ld2=ld*ld;
310  G4double sp=std::sqrt(P);
311  G4double p2=P*P;
312  G4double p4=p2*p2;
313  G4double lm=P-1.;
314  G4double md=lm*lm+.372;
315  G4double El=(.0557*ld2+2.23)/(1.-.7/sp+.1/p4);
316  G4double To=(.3*ld2+19.5)/(1.+.46/sp+1.6/p4);
317  sigma=(To-El)+.6/md;
318  }
319  else if(tZ<97 && tN<152) // General solution
320  {
321  G4double p2=P*P;
322  G4double p4=p2*p2;
323  G4double a=tN+tZ; // A of the target
324  G4double al=std::log(a);
325  G4double sa=std::sqrt(a);
326  G4double asa=a*sa;
327  G4double a2=a*a;
328  G4double a3=a2*a;
329  G4double a4=a2*a2;
330  G4double a8=a4*a4;
331  G4double a12=a8*a4;
332  G4double f=.6; // Default values for deutrons
333  G4double r=.5;
334  G4double gg=3.7;
335  G4double c=36.;
336  G4double ss=3.5;
337  G4double t=3.;
338  G4double u=.44;
339  G4double v=5.E-9;
340  if(tZ>1 && tN>1) // More than deuteron
341  {
342  f=1.;
343  r=1./(1.+.007*a2);
344  gg=4.2;
345  c=52.*std::exp(al*.6)*(1.+95./a2)/(1.+9./a)/(1.+46./a2);
346  ss=(40.+.14*a)/(1.+12./a);
347  G4double y=std::exp(al*1.7);
348  t=.185*y/(1.+.00012*y);
349  u=(1.+80./asa)/(1.+200./asa);
350  v=(1.+3.E-6*a4*(1.+6.E-7*a3+4.E10/a12))/a3/20000.;
351  }
352  G4double d=lP-gg;
353  G4double w=P-1.;
354  G4double rD=ss/(w*w+.36);
355  G4double h=P-.44;
356  G4double rR=t/(h*h+u*u);
357  sigma=(f*d*d+c)/(1.+r/std::sqrt(P)+1./p4)+(rD+rR)/(1+v/p4/p4);
358  }
359  else
360  {
361  G4cerr<<"-Warning-G4ChipsKaonPlusNuclearCroSect::CSForm:Bad A, Z="<<tZ<<", N="<<tN<<G4endl;
362  sigma=0.;
363  }
364  if(sigma<0.) return 0.;
365  return sigma;
366 }
367 
369 {
370  if(DX<=0. || N<2)
371  {
372  G4cerr<<"***G4ChipsKaonPlusInelasticXS::EquLinearFit: DX="<<DX<<", N="<<N<<G4endl;
373  return Y[0];
374  }
375 
376  G4int N2=N-2;
377  G4double d=(X-X0)/DX;
378  G4int jj=static_cast<int>(d);
379  if (jj<0) jj=0;
380  else if(jj>N2) jj=N2;
381  d-=jj; // excess
382  G4double yi=Y[jj];
383  G4double sigma=yi+(Y[jj+1]-yi)*d;
384 
385  return sigma;
386 }
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
int G4int
Definition: G4Types.hh:78
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:163
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:161
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