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