Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/cross_sections/src/G4ChipsProtonInelasticXS.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

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