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