Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ChipsProtonInelasticXS.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: G4ChipsProtonInelasticXS 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// ****************************************************************************************
37// Short description: Cross-sections extracted (by W.Pokorski) from the CHIPS package for
38// proton-nuclear interactions. Original author: M. Kossov
39// -------------------------------------------------------------------------------------
40//
41
42
44#include "G4SystemOfUnits.hh"
45#include "G4DynamicParticle.hh"
47#include "G4Proton.hh"
48
49// factory
51//
53
55{
56 // Initialization of the
57 lastLEN=0; // Pointer to the lastArray of LowEn CS
58 lastHEN=0; // Pointer to the lastArray of HighEn CS
59 lastN=0; // The last N of calculated nucleus
60 lastZ=0; // The last Z of calculated nucleus
61 lastP=0.; // Last used in cross section Momentum
62 lastTH=0.; // Last threshold momentum
63 lastCS=0.; // Last value of the Cross Section
64 lastI=0; // The last position in the DAMDB
65
66 LEN = new std::vector<G4double*>;
67 HEN = new std::vector<G4double*>;
68}
69
71{
72 /*
73 G4int lens=LEN->size();
74 for(G4int i=0; i<lens; ++i) delete[] (*LEN)[i];
75 delete LEN;
76 G4int hens=HEN->size();
77 for(G4int i=0; i<hens; ++i) delete[] (*HEN)[i];
78 delete HEN;
79 */
80}
81
83 const G4Element*,
84 const G4Material*)
85{
86 G4ParticleDefinition* particle = Pt->GetDefinition();
87 if (particle == G4Proton::Proton() ) return true;
88 return false;
89}
90
91
92// The main member function giving the collision cross section (P is in IU, CS is in mb)
93// Make pMom in independent units ! (Now it is MeV)
95 const G4Isotope*,
96 const G4Element*,
97 const G4Material*)
98{
99 G4double pMom=Pt->GetTotalMomentum();
100 G4int tgN = A - tgZ;
101
102 return GetChipsCrossSection(pMom, tgZ, tgN, 2212);
103}
104
106{
107 static G4int j; // A#0f Z/N-records already tested in AMDB
108 static std::vector <G4int> colN; // Vector of N for calculated nuclei (isotops)
109 static std::vector <G4int> colZ; // Vector of Z for calculated nuclei (isotops)
110 static std::vector <G4double> colP; // Vector of last momenta for the reaction
111 static std::vector <G4double> colTH; // Vector of energy thresholds for the reaction
112 static std::vector <G4double> colCS; // Vector of last cross sections for the reaction
113 // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
114
115 G4bool in=false; // By default the isotope must be found in the AMDB
116 if(tgN!=lastN || tgZ!=lastZ) // The nucleus was not the last used isotope
117 {
118 in = false; // By default the isotope haven't been found in AMDB
119 lastP = 0.; // New momentum history (nothing to compare with)
120 lastN = tgN; // The last N of the calculated nucleus
121 lastZ = tgZ; // The last Z of the calculated nucleus
122 lastI = colN.size(); // Size of the Associative Memory DB in the heap
123 j = 0; // A#0f records found in DB for this projectile
124 if(lastI) for(G4int i=0; i<lastI; i++) // AMDB exists, try to find the (Z,N) isotope
125 {
126 if(colN[i]==tgN && colZ[i]==tgZ) // Try the record "i" in the AMDB
127 {
128 lastI=i; // Remember the index for future fast/last use
129 lastTH =colTH[i]; // The last THreshold (A-dependent)
130 if(pMom<=lastTH)
131 {
132 return 0.; // Energy is below the Threshold value
133 }
134 lastP =colP [i]; // Last Momentum (A-dependent)
135 lastCS =colCS[i]; // Last CrossSect (A-dependent)
136 in = true; // This is the case when the isotop is found in DB
137 // Momentum pMom is in IU ! @@ Units
138 lastCS=CalculateCrossSection(-1,j,2212,lastZ,lastN,pMom); // read & update
139 if(lastCS<=0. && pMom>lastTH) // Correct the threshold (@@ No intermediate Zeros)
140 {
141 lastCS=0.;
142 lastTH=pMom;
143 }
144 break; // Go out of the LOOP
145 }
146 j++; // Increment a#0f records found in DB
147 }
148 if(!in) // This isotope has not been calculated previously
149 {
150 //!!The slave functions must provide cross-sections in millibarns (mb) !! (not in IU)
151 lastCS=CalculateCrossSection(0,j,2212,lastZ,lastN,pMom); //calculate & create
152 //if(lastCS>0.) // It means that the AMBD was initialized
153 //{
154
155 lastTH = 0; //ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
156 colN.push_back(tgN);
157 colZ.push_back(tgZ);
158 colP.push_back(pMom);
159 colTH.push_back(lastTH);
160 colCS.push_back(lastCS);
161 //} // M.K. Presence of H1 with high threshold breaks the syncronization
162 return lastCS*millibarn;
163 } // End of creation of the new set of parameters
164 else
165 {
166 colP[lastI]=pMom;
167 colCS[lastI]=lastCS;
168 }
169 } // End of parameters udate
170 else if(pMom<=lastTH)
171 {
172 return 0.; // Momentum is below the Threshold Value -> CS=0
173 }
174 else // It is the last used -> use the current tables
175 {
176 lastCS=CalculateCrossSection(1,j,2212,lastZ,lastN,pMom); // Only read and UpdateDB
177 lastP=pMom;
178 }
179 return lastCS*millibarn;
180}
181
182// The main member function giving the gamma-A cross section (E in GeV, CS in mb)
183G4double G4ChipsProtonInelasticXS::CalculateCrossSection(G4int F, G4int I,
184 G4int, G4int targZ, G4int targN, G4double Momentum)
185{
186 static const G4double THmin=27.; // default minimum Momentum (MeV/c) Threshold
187 static const G4double THmiG=THmin*.001; // minimum Momentum (GeV/c) Threshold
188 static const G4double dP=10.; // step for the LEN (Low ENergy) table MeV/c
189 static const G4double dPG=dP*.001; // step for the LEN (Low ENergy) table GeV/c
190 static const G4int nL=105; // A#of LEN points in E (step 10 MeV/c)
191 static const G4double Pmin=THmin+(nL-1)*dP; // minP for the HighE part with safety
192 static const G4double Pmax=227000.; // maxP for the HEN (High ENergy) part 227 GeV
193 static const G4int nH=224; // A#of HEN points in lnE
194 static const G4double milP=std::log(Pmin);// Low logarithm energy for the HEN part
195 static const G4double malP=std::log(Pmax);// High logarithm energy (each 2.75 percent)
196 static const G4double dlP=(malP-milP)/(nH-1); // Step in log energy in the HEN part
197 static const G4double milPG=std::log(.001*Pmin);// Low logarithmEnergy for HEN part GeV/c
198 G4double sigma=0.;
199 if(F&&I) sigma=0.; // @@ *!* Fake line *!* to use F & I !!!Temporary!!!
200 //G4double A=targN+targZ; // A of the target
201 if(F<=0) // This isotope was not the last used isotop
202 {
203 if(F<0) // This isotope was found in DAMDB =-----=> RETRIEVE
204 {
205 G4int sync=LEN->size();
206 if(sync<=I) G4cout<<"*!*G4QProtonNuclCS::CalcCrossSect:Sync="<<sync<<"<="<<I<<G4endl;
207 lastLEN=(*LEN)[I]; // Pointer to prepared LowEnergy cross sections
208 lastHEN=(*HEN)[I]; // Pointer to prepared High Energy cross sections
209 }
210 else // This isotope wasn't calculated before => CREATE
211 {
212 lastLEN = new G4double[nL]; // Allocate memory for the new LEN cross sections
213 lastHEN = new G4double[nH]; // Allocate memory for the new HEN cross sections
214 // --- Instead of making a separate function ---
215 G4double P=THmiG; // Table threshold in GeV/c
216 for(G4int k=0; k<nL; k++)
217 {
218 lastLEN[k] = CrossSectionLin(targZ, targN, P);
219 P+=dPG;
220 }
221 G4double lP=milPG;
222 for(G4int n=0; n<nH; n++)
223 {
224 lastHEN[n] = CrossSectionLog(targZ, targN, lP);
225 lP+=dlP;
226 }
227 // --- End of possible separate function
228 // *** The synchronization check ***
229 G4int sync=LEN->size();
230 if(sync!=I)
231 {
232 G4cout<<"***G4ChipsProtonNuclCS::CalcCrossSect: Sinc="<<sync<<"#"<<I<<", Z=" <<targZ
233 <<", N="<<targN<<", F="<<F<<G4endl;
234 //G4Exception("G4ProtonNuclearCS::CalculateCS:","39",FatalException,"overflow DB");
235 }
236 LEN->push_back(lastLEN); // remember the Low Energy Table
237 HEN->push_back(lastHEN); // remember the High Energy Table
238 } // End of creation of the new set of parameters
239 } // End of parameters udate
240 // =------------------= NOW the Magic Formula =-----------------------=
241 if (Momentum<lastTH) return 0.; // It must be already checked in the interface class
242 else if (Momentum<Pmin) // High Energy region
243 {
244 sigma=EquLinearFit(Momentum,nL,THmin,dP,lastLEN);
245 }
246 else if (Momentum<Pmax) // High Energy region
247 {
248 G4double lP=std::log(Momentum);
249 sigma=EquLinearFit(lP,nH,milP,dlP,lastHEN);
250 }
251 else // UHE region (calculation, not frequent)
252 {
253 G4double P=0.001*Momentum; // Approximation formula is for P in GeV/c
254 sigma=CrossSectionFormula(targZ, targN, P, std::log(P));
255 }
256 if(sigma<0.) return 0.;
257 return sigma;
258}
259
260// Electromagnetic momentum-threshold (in MeV/c)
261G4double G4ChipsProtonInelasticXS::ThresholdMomentum(G4int tZ, G4int tN)
262{
263 static const G4double third=1./3.;
264 static const G4double pM = G4Proton::Proton()->Definition()->GetPDGMass(); // Projectile mass in MeV
265 static const G4double tpM= pM+pM; // Doubled projectile mass (MeV)
266
267 G4double tA=tZ+tN;
268 if(tZ<.99 || tN<0.) return 0.;
269 else if(tZ==1 && tN==0) return 800.; // A threshold on the free proton
270 //G4double dE=1.263*tZ/(1.+std::pow(tA,third));
271 G4double dE=tZ/(1.+std::pow(tA,third)); // Safety for diffused edge of the nucleus (QE)
272 G4double tM=931.5*tA;
273 G4double T=dE+dE*(dE/2+pM)/tM;
274 return std::sqrt(T*(tpM+T));
275}
276
277// Calculation formula for proton-nuclear inelastic cross-section (mb) (P in GeV/c)
278G4double G4ChipsProtonInelasticXS::CrossSectionLin(G4int tZ, G4int tN, G4double P)
279{
280 G4double sigma=0.;
281 if(P<ThresholdMomentum(tZ,tN)*.001) return sigma;
282 G4double lP=std::log(P);
283 if(tZ==1&&!tN){if(P>.35) sigma=CrossSectionFormula(tZ,tN,P,lP);}// s(pp)=0 below 350Mev/c
284 else if(tZ<97 && tN<152) // General solution
285 {
286 G4double pex=0.;
287 G4double pos=0.;
288 G4double wid=1.;
289 if(tZ==13 && tN==14) // Excited metastable states
290 {
291 pex=230.;
292 pos=.13;
293 wid=8.e-5;
294 }
295 else if(tZ<7)
296 {
297 if(tZ==6 && tN==6)
298 {
299 pex=320.;
300 pos=.14;
301 wid=7.e-6;
302 }
303 else if(tZ==5 && tN==6)
304 {
305 pex=270.;
306 pos=.17;
307 wid=.002;
308 }
309 else if(tZ==4 && tN==5)
310 {
311 pex=600.;
312 pos=.132;
313 wid=.005;
314 }
315 else if(tZ==3 && tN==4)
316 {
317 pex=280.;
318 pos=.19;
319 wid=.0025;
320 }
321 else if(tZ==3 && tN==3)
322 {
323 pex=370.;
324 pos=.171;
325 wid=.006;
326 }
327 else if(tZ==2 && tN==1)
328 {
329 pex=30.;
330 pos=.22;
331 wid=.0005;
332 }
333 }
334 sigma=CrossSectionFormula(tZ,tN,P,lP);
335 if(pex>0.)
336 {
337 G4double dp=P-pos;
338 sigma+=pex*std::exp(-dp*dp/wid);
339 }
340 }
341 else
342 {
343 G4cerr<<"-Warning-G4ChipsProtonNuclearXS::CSLin:*Bad A* Z="<<tZ<<", N="<<tN<<G4endl;
344 sigma=0.;
345 }
346 if(sigma<0.) return 0.;
347 return sigma;
348}
349
350// Calculation formula for proton-nuclear inelastic cross-section (mb) log(P in GeV/c)
351G4double G4ChipsProtonInelasticXS::CrossSectionLog(G4int tZ, G4int tN, G4double lP)
352{
353 G4double P=std::exp(lP);
354 return CrossSectionFormula(tZ, tN, P, lP);
355}
356// Calculation formula for proton-nuclear inelastic cross-section (mb) log(P in GeV/c)
357G4double G4ChipsProtonInelasticXS::CrossSectionFormula(G4int tZ, G4int tN,
358 G4double P, G4double lP)
359{
360 G4double sigma=0.;
361 if(tZ==1 && !tN) // pp interaction (from G4QuasiElasticRatios)
362 {
363 G4double p2=P*P;
364 G4double lp=lP-3.5;
365 G4double lp2=lp*lp;
366 G4double rp2=1./p2;
367 G4double El=(.0557*lp2+6.72+30./P)/(1.+.49*rp2/P);
368 G4double To=(.3*lp2+38.2)/(1.+.54*rp2*rp2);
369 sigma=To-El;
370 }
371 else if(tZ<97 && tN<152) // General solution
372 {
373 //G4double lP=std::log(P); // Already calculated
374 G4double d=lP-4.2;
375 G4double p2=P*P;
376 G4double p4=p2*p2;
377 G4double a=tN+tZ; // A of the target
378 G4double al=std::log(a);
379 G4double sa=std::sqrt(a);
380 G4double a2=a*a;
381 G4double a2s=a2*sa;
382 G4double a4=a2*a2;
383 G4double a8=a4*a4;
384 G4double a12=a8*a4;
385 G4double a16=a8*a8;
386 G4double c=(170.+3600./a2s)/(1.+65./a2s);
387 G4double dl=al-3.;
388 G4double dl2=dl*dl;
389 G4double r=.21+.62*dl2/(1.+.5*dl2);
390 G4double gg=40.*std::exp(al*0.712)/(1.+12.2/a)/(1.+34./a2);
391 G4double e=318.+a4/(1.+.0015*a4/std::exp(al*0.09))/(1.+4.e-28*a12)+
392 8.e-18/(1./a16+1.3e-20)/(1.+1.e-21*a12);
393 G4double ss=3.57+.009*a2/(1.+.0001*a2*a);
394 G4double h=(.01/a4+2.5e-6/a)*(1.+6.e-6*a2*a)/(1.+6.e7/a12/a2);
395 sigma=(c+d*d)/(1.+r/p4)+(gg+e*std::exp(-ss*P))/(1.+h/p4/p4);
396 }
397 else
398 {
399 G4cerr<<"-Warning-G4QProtonNuclearCroSect::CSForm:*Bad A* Z="<<tZ<<", N="<<tN<<G4endl;
400 sigma=0.;
401 }
402 if(sigma<0.) return 0.;
403 return sigma;
404}
405
406G4double G4ChipsProtonInelasticXS::EquLinearFit(G4double X, G4int N, G4double X0, G4double DX, G4double* Y)
407{
408 if(DX<=0. || N<2)
409 {
410 G4cerr<<"***G4ChipsProtonInelasticXS::EquLinearFit: DX="<<DX<<", N="<<N<<G4endl;
411 return Y[0];
412 }
413
414 G4int N2=N-2;
415 G4double d=(X-X0)/DX;
416 G4int j=static_cast<int>(d);
417 if (j<0) j=0;
418 else if(j>N2) j=N2;
419 d-=j; // excess
420 G4double yi=Y[j];
421 G4double sigma=yi+(Y[j+1]-yi)*d;
422
423 return sigma;
424}
#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 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)
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
G4ParticleDefinition * GetDefinition() const
G4double GetTotalMomentum() const
static G4Proton * Definition()
Definition: G4Proton.cc:49
static G4Proton * Proton()
Definition: G4Proton.cc:93