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