Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ChipsPionPlusInelasticXS.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: G4ChipsPionPlusInelasticXS 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 // pion interactions. Original author: M. Kossov
37 // -------------------------------------------------------------------------------------
38 //
39 
41 #include "G4SystemOfUnits.hh"
42 #include "G4DynamicParticle.hh"
43 #include "G4ParticleDefinition.hh"
44 #include "G4PionPlus.hh"
45 
46 // factory
47 #include "G4CrossSectionFactory.hh"
48 //
50 
52 {
53  // Initialization of the
54  lastLEN=0; // Pointer to lastArray of LowEn CS
55  lastHEN=0; // Pointer to lastArray of HighEn CS
56  lastN=0; // The last N of calculated nucleus
57  lastZ=0; // The last Z of calculated nucleus
58  lastP=0.; // Last used in cross section Momentum
59  lastTH=0.; // Last threshold momentum
60  lastCS=0.; // Last value of the Cross Section
61  lastI=0; // The last position in the DAMDB
62  LEN = new std::vector<G4double*>;
63  HEN = new std::vector<G4double*>;
64 }
65 
66 
68 {
69  G4int lens=LEN->size();
70  for(G4int i=0; i<lens; ++i) delete[] (*LEN)[i];
71  delete LEN;
72  G4int hens=HEN->size();
73  for(G4int i=0; i<hens; ++i) delete[] (*HEN)[i];
74  delete HEN;
75 }
76 
78  const G4Element*,
79  const G4Material*)
80 {
81  G4ParticleDefinition* particle = Pt->GetDefinition();
82  if (particle == G4PionPlus::PionPlus() ) return true;
83  return false;
84 }
85 
86 // The main member function giving the collision cross section (P is in IU, CS is in mb)
87 // Make pMom in independent units ! (Now it is MeV)
89  const G4Isotope*,
90  const G4Element*,
91  const G4Material*)
92 {
93  G4double pMom=Pt->GetTotalMomentum();
94  G4int tgN = A - tgZ;
95 
96  return GetChipsCrossSection(pMom, tgZ, tgN, 211);
97 }
98 
99 
101 {
102  static G4ThreadLocal G4int j; // A#0f Z/N-records already tested in AMDB
103  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)
104  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)
105  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
106  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
107  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
108  // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
109 
110  G4bool in=false; // By default the isotope must be found in the AMDB
111  if(tgN!=lastN || tgZ!=lastZ) // The nucleus was not the last used isotope
112  {
113  in = false; // By default the isotope haven't be found in AMDB
114  lastP = 0.; // New momentum history (nothing to compare with)
115  lastN = tgN; // The last N of the calculated nucleus
116  lastZ = tgZ; // The last Z of the calculated nucleus
117  lastI = colN.size(); // Size of the Associative Memory DB in the heap
118  j = 0; // A#0f records found in DB for this projectile
119  if(lastI) for(G4int i=0; i<lastI; i++) // AMDB exists, try to find the (Z,N) isotope
120  {
121  if(colN[i]==tgN && colZ[i]==tgZ) // Try the record "i" in the AMDB
122  {
123  lastI=i; // Remember the index for future fast/last use
124  lastTH =colTH[i]; // The last THreshold (A-dependent)
125  if(pMom<=lastTH)
126  {
127  return 0.; // Energy is below the Threshold value
128  }
129  lastP =colP [i]; // Last Momentum (A-dependent)
130  lastCS =colCS[i]; // Last CrossSect (A-dependent)
131  in = true; // This is the case when the isotop is found in DB
132  // Momentum pMom is in IU ! @@ Units
133  lastCS=CalculateCrossSection(-1,j,211,lastZ,lastN,pMom); // read & update
134  if(lastCS<=0. && pMom>lastTH) // Correct the threshold (@@ No intermediate Zeros)
135  {
136  lastCS=0.;
137  lastTH=pMom;
138  }
139  break; // Go out of the LOOP
140  }
141  j++; // Increment a#0f records found in DB
142  }
143  if(!in) // This isotope has not been calculated previously
144  {
145  //!!The slave functions must provide cross-sections in millibarns (mb) !! (not in IU)
146  lastCS=CalculateCrossSection(0,j,211,lastZ,lastN,pMom); //calculate & create
147  //if(lastCS>0.) // It means that the AMBD was initialized
148  //{
149 
150  lastTH = 0; //ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
151  colN.push_back(tgN);
152  colZ.push_back(tgZ);
153  colP.push_back(pMom);
154  colTH.push_back(lastTH);
155  colCS.push_back(lastCS);
156  //} // M.K. Presence of H1 with high threshold breaks the syncronization
157  return lastCS*millibarn;
158  } // End of creation of the new set of parameters
159  else
160  {
161  colP[lastI]=pMom;
162  colCS[lastI]=lastCS;
163  }
164  } // End of parameters udate
165  else if(pMom<=lastTH)
166  {
167  return 0.; // Momentum is below the Threshold Value -> CS=0
168  }
169  else // It is the last used -> use the current tables
170  {
171  lastCS=CalculateCrossSection(1,j,211,lastZ,lastN,pMom); // Only read and UpdateDB
172  lastP=pMom;
173  }
174  return lastCS*millibarn;
175 }
176 
177 // The main member function giving the gamma-A cross section (E in GeV, CS in mb)
178 G4double G4ChipsPionPlusInelasticXS::CalculateCrossSection(G4int F, G4int I,
179  G4int, G4int targZ, G4int targN, G4double Momentum)
180 {
181  static const G4double THmin=27.; // default minimum Momentum (MeV/c) Threshold
182  static const G4double THmiG=THmin*.001; // minimum Momentum (GeV/c) Threshold
183  static const G4double dP=10.; // step for the LEN (Low ENergy) table MeV/c
184  static const G4double dPG=dP*.001; // step for the LEN (Low ENergy) table GeV/c
185  static const G4int nL=105; // A#of LEN points in E (step 10 MeV/c)
186  static const G4double Pmin=THmin+(nL-1)*dP; // minP for the HighE part with safety
187  static const G4double Pmax=227000.; // maxP for the HEN (High ENergy) part 227 GeV
188  static const G4int nH=224; // A#of HEN points in lnE
189  static const G4double milP=std::log(Pmin);// Low logarithm energy for the HEN part
190  static const G4double malP=std::log(Pmax);// High logarithm energy (each 2.75 percent)
191  static const G4double dlP=(malP-milP)/(nH-1); // Step in log energy in the HEN part
192  static const G4double milPG=std::log(.001*Pmin);// Low logarithmEnergy for HEN part GeV/c
193  G4double sigma=0.;
194  if(F&&I) sigma=0.; // @@ *!* Fake line *!* to use F & I !!!Temporary!!!
195  //G4double A=targN+targZ; // A of the target
196  if(F<=0) // This isotope was not the last used isotop
197  {
198  if(F<0) // This isotope was found in DAMDB =-----=> RETRIEVE
199  {
200  G4int sync=LEN->size();
201  if(sync<=I) G4cerr<<"*!*G4ChipsPiMinusNuclCS::CalcCrosSect:Sync="<<sync<<"<="<<I<<G4endl;
202  lastLEN=(*LEN)[I]; // Pointer to prepared LowEnergy cross sections
203  lastHEN=(*HEN)[I]; // Pointer to prepared High Energy cross sections
204  }
205  else // This isotope wasn't calculated before => CREATE
206  {
207  lastLEN = new G4double[nL]; // Allocate memory for the new LEN cross sections
208  lastHEN = new G4double[nH]; // Allocate memory for the new HEN cross sections
209  // --- Instead of making a separate function ---
210  G4double P=THmiG; // Table threshold in GeV/c
211  for(G4int k=0; k<nL; k++)
212  {
213  lastLEN[k] = CrossSectionLin(targZ, targN, P);
214  P+=dPG;
215  }
216  G4double lP=milPG;
217  for(G4int n=0; n<nH; n++)
218  {
219  lastHEN[n] = CrossSectionLog(targZ, targN, lP);
220  lP+=dlP;
221  }
222  // --- End of possible separate function
223  // *** The synchronization check ***
224  G4int sync=LEN->size();
225  if(sync!=I)
226  {
227  G4cerr<<"***G4ChipsPiMinusNuclCS::CalcCrossSect: Sinc="<<sync<<"#"<<I<<", Z=" <<targZ
228  <<", N="<<targN<<", F="<<F<<G4endl;
229  //G4Exception("G4PiMinusNuclearCS::CalculateCS:","39",FatalException,"DBoverflow");
230  }
231  LEN->push_back(lastLEN); // remember the Low Energy Table
232  HEN->push_back(lastHEN); // remember the High Energy Table
233  } // End of creation of the new set of parameters
234  } // End of parameters udate
235  // =-----------------= NOW the Magic Formula =-------------------------=
236  if (Momentum<lastTH) return 0.; // It must be already checked in the interface class
237  else if (Momentum<Pmin) // High Energy region
238  {
239  sigma=EquLinearFit(Momentum,nL,THmin,dP,lastLEN);
240  }
241  else if (Momentum<Pmax) // High Energy region
242  {
243  G4double lP=std::log(Momentum);
244  sigma=EquLinearFit(lP,nH,milP,dlP,lastHEN);
245  }
246  else // UHE region (calculation, not frequent)
247  {
248  G4double P=0.001*Momentum; // Approximation formula is for P in GeV/c
249  sigma=CrossSectionFormula(targZ, targN, P, std::log(P));
250  }
251  if(sigma<0.) return 0.;
252  return sigma;
253 }
254 
255 // Electromagnetic momentum-threshold (in MeV/c)
256 G4double G4ChipsPionPlusInelasticXS::ThresholdMomentum(G4int tZ, G4int tN)
257 {
258  static const G4double third=1./3.;
259  static const G4double pM = G4PionPlus::PionPlus()->Definition()->GetPDGMass(); // Projectile mass in MeV
260  static const G4double tpM= pM+pM; // Doubled projectile mass (MeV)
261  G4double tA=tZ+tN;
262  if(tZ<.99 || tN<0.) return 0.;
263  else if(tZ==1 && tN==0) return 300.; // A threshold on the free proton
264  //G4double dE=1.263*tZ/(1.+std::pow(tA,third));
265  G4double dE=tZ/(1.+std::pow(tA,third)); // Safety for diffused edge of the nucleus (QE)
266  G4double tM=931.5*tA;
267  G4double T=dE+dE*(dE/2+pM)/tM;
268  return std::sqrt(T*(tpM+T));
269 }
270 
271 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) (P in GeV/c)
272 G4double G4ChipsPionPlusInelasticXS::CrossSectionLin(G4int tZ, G4int tN, G4double P)
273 {
274  G4double lP=std::log(P);
275  return CrossSectionFormula(tZ, tN, P, lP);
276 }
277 
278 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
279 G4double G4ChipsPionPlusInelasticXS::CrossSectionLog(G4int tZ, G4int tN, G4double lP)
280 {
281  G4double P=std::exp(lP);
282  return CrossSectionFormula(tZ, tN, P, lP);
283 }
284 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
285 G4double G4ChipsPionPlusInelasticXS::CrossSectionFormula(G4int tZ, G4int tN,
286  G4double P, G4double lP)
287 {
288  G4double sigma=0.;
289  if(tZ==1 && !tN) // PiPlus-Proton interaction from G4QuasiElRatios
290  {
291  G4double ld=lP-3.5;
292  G4double ld2=ld*ld;
293  G4double p2=P*P;
294  G4double p4=p2*p2;
295  G4double sp=std::sqrt(P);
296  G4double lm=lP-.32;
297  G4double md=lm*lm+.04;
298  G4double El=(.0557*ld2+2.4+6./sp)/(1.+3./p4);
299  G4double To=(.3*ld2+22.3+5./sp)/(1.+1./p4);
300  sigma=(To-El)+.1/md;
301  }
302  else if(tZ==1 && tN==1) // pimp_tot
303  {
304  G4double p2=P*P;
305  G4double d=lP-2.7;
306  G4double f=lP+1.25;
307  G4double gg=lP-.017;
308  sigma=(.55*d*d+38.+23./std::sqrt(P))/(1.+.3/p2/p2)+18./(f*f+.1089)+.02/(gg*gg+.0025);
309  }
310  else if(tZ<97 && tN<152) // General solution
311  {
312  G4double d=lP-4.2;
313  G4double p2=P*P;
314  G4double p4=p2*p2;
315  G4double a=tN+tZ; // A of the target
316  G4double al=std::log(a);
317  G4double sa=std::sqrt(a);
318  G4double ssa=std::sqrt(sa);
319  G4double a2=a*a;
320  G4double c=41.*std::exp(al*.68)*(1.+44./a2)/(1.+8./a)/(1.+200./a2/a2);
321  G4double f=290.*ssa/(1.+34./a/ssa);
322  G4double gg=-1.32-al*.043;
323  G4double u=lP-gg;
324  G4double h=al*(.4-.055*al);
325  G4double r=.01+a2*5.E-8;
326  sigma=(c+d*d)/(1.+(.2-.009*sa)/p4)+f/(u*u+h*h)/(1.+r/p2);
327  }
328  else
329  {
330  G4cerr<<"-Warning-G4ChipsPiPlusNuclearCroSect::CSForm:*Bad A* Z="<<tZ<<", N="<<tN<<G4endl;
331  sigma=0.;
332  }
333  if(sigma<0.) return 0.;
334  return sigma;
335 }
336 
337 G4double G4ChipsPionPlusInelasticXS::EquLinearFit(G4double X, G4int N, G4double X0, G4double DX, G4double* Y)
338 {
339  if(DX<=0. || N<2)
340  {
341  G4cerr<<"***G4ChipsPionPlusInelasticXS::EquLinearFit: DX="<<DX<<", N="<<N<<G4endl;
342  return Y[0];
343  }
344 
345  G4int N2=N-2;
346  G4double d=(X-X0)/DX;
347  G4int j=static_cast<int>(d);
348  if (j<0) j=0;
349  else if(j>N2) j=N2;
350  d-=j; // excess
351  G4double yi=Y[j];
352  G4double sigma=yi+(Y[j+1]-yi)*d;
353 
354  return sigma;
355 }
static G4PionPlus * Definition()
Definition: G4PionPlus.cc:52
G4ParticleDefinition * GetDefinition() const
#define G4ThreadLocal
Definition: tls.hh:52
int G4int
Definition: G4Types.hh:78
int millibarn
Definition: hepunit.py:40
virtual G4bool IsIsoApplicable(const G4DynamicParticle *Pt, G4int Z, G4int A, const G4Element *elm, const G4Material *mat)
G4_DECLARE_XS_FACTORY(G4ChipsPionPlusInelasticXS)
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
G4double GetTotalMomentum() const
bool G4bool
Definition: G4Types.hh:79
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
const G4int n
G4double GetPDGMass() const
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int tgZ, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
#define G4endl
Definition: G4ios.hh:61
**D E S C R I P T I O N
double G4double
Definition: G4Types.hh:76
G4GLOB_DLL std::ostream G4cerr