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