G4ChipsProtonInelasticXS.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // The lust update: M.V. Kossov, CERN/ITEP(Moscow) 17-June-02
00028 // GEANT4 tag $Name: not supported by cvs2svn $
00029 //
00030 //
00031 // G4 Physics class: G4ChipsProtonInelasticXS for gamma+A cross sections
00032 // Created: M.V. Kossov, CERN/ITEP(Moscow), 20-Dec-03
00033 // The last update: M.V. Kossov, CERN/ITEP (Moscow) 15-Feb-04
00034 //
00035 //
00036 // ****************************************************************************************
00037 // Short description: Cross-sections extracted (by W.Pokorski) from the CHIPS package for 
00038 // proton-nuclear  interactions. Original author: M. Kossov
00039 // -------------------------------------------------------------------------------------
00040 //
00041 
00042 
00043 #include "G4ChipsProtonInelasticXS.hh"
00044 #include "G4SystemOfUnits.hh"
00045 #include "G4DynamicParticle.hh"
00046 #include "G4ParticleDefinition.hh"
00047 #include "G4Proton.hh"
00048 
00049 // factory
00050 #include "G4CrossSectionFactory.hh"
00051 //
00052 G4_DECLARE_XS_FACTORY(G4ChipsProtonInelasticXS);
00053 
00054 G4ChipsProtonInelasticXS::G4ChipsProtonInelasticXS():G4VCrossSectionDataSet(Default_Name())
00055 {
00056   // Initialization of the
00057   lastLEN=0; // Pointer to the lastArray of LowEn CS
00058   lastHEN=0; // Pointer to the lastArray of HighEn CS
00059   lastN=0;   // The last N of calculated nucleus
00060   lastZ=0;   // The last Z of calculated nucleus
00061   lastP=0.;  // Last used in cross section Momentum
00062   lastTH=0.; // Last threshold momentum
00063   lastCS=0.; // Last value of the Cross Section
00064   lastI=0;   // The last position in the DAMDB
00065 
00066   LEN = new std::vector<G4double*>;
00067   HEN = new std::vector<G4double*>;
00068 }
00069 
00070 G4ChipsProtonInelasticXS::~G4ChipsProtonInelasticXS()
00071 {
00072   /*
00073   G4int lens=LEN->size();
00074   for(G4int i=0; i<lens; ++i) delete[] (*LEN)[i];
00075   delete LEN;
00076   G4int hens=HEN->size();
00077   for(G4int i=0; i<hens; ++i) delete[] (*HEN)[i];
00078   delete HEN;
00079   */
00080 }
00081 
00082 G4bool G4ChipsProtonInelasticXS::IsIsoApplicable(const G4DynamicParticle* Pt, G4int, G4int,    
00083                                  const G4Element*,
00084                                  const G4Material*)
00085 {
00086   G4ParticleDefinition* particle = Pt->GetDefinition();
00087   if (particle == G4Proton::Proton()      ) return true;
00088   return false;
00089 }
00090 
00091 
00092 // The main member function giving the collision cross section (P is in IU, CS is in mb)
00093 // Make pMom in independent units ! (Now it is MeV)
00094 G4double G4ChipsProtonInelasticXS::GetIsoCrossSection(const G4DynamicParticle* Pt, G4int tgZ, G4int A,  
00095                                                         const G4Isotope*,
00096                                                         const G4Element*,
00097                                                         const G4Material*)
00098 {
00099   G4double pMom=Pt->GetTotalMomentum();
00100   G4int tgN = A - tgZ;
00101 
00102   return GetChipsCrossSection(pMom, tgZ, tgN, 2212);
00103 }
00104 
00105 G4double G4ChipsProtonInelasticXS::GetChipsCrossSection(G4double pMom, G4int tgZ, G4int tgN, G4int)
00106 {
00107   static G4int j;                      // A#0f Z/N-records already tested in AMDB
00108   static std::vector <G4int>    colN;  // Vector of N for calculated nuclei (isotops)
00109   static std::vector <G4int>    colZ;  // Vector of Z for calculated nuclei (isotops)
00110   static std::vector <G4double> colP;  // Vector of last momenta for the reaction
00111   static std::vector <G4double> colTH; // Vector of energy thresholds for the reaction
00112   static std::vector <G4double> colCS; // Vector of last cross sections for the reaction
00113   // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
00114 
00115   G4bool in=false;                     // By default the isotope must be found in the AMDB
00116   if(tgN!=lastN || tgZ!=lastZ)         // The nucleus was not the last used isotope
00117   {
00118     in      = false;                   // By default the isotope haven't been found in AMDB
00119     lastP   = 0.;                      // New momentum history (nothing to compare with)
00120     lastN   = tgN;                     // The last N of the calculated nucleus
00121     lastZ   = tgZ;                     // The last Z of the calculated nucleus
00122     lastI   = colN.size();             // Size of the Associative Memory DB in the heap
00123     j  = 0;                            // A#0f records found in DB for this projectile
00124     if(lastI) for(G4int i=0; i<lastI; i++) // AMDB exists, try to find the (Z,N) isotope
00125     {
00126       if(colN[i]==tgN && colZ[i]==tgZ) // Try the record "i" in the AMDB
00127       {
00128         lastI=i;                       // Remember the index for future fast/last use
00129         lastTH =colTH[i];              // The last THreshold (A-dependent)
00130         if(pMom<=lastTH)
00131         {
00132           return 0.;                   // Energy is below the Threshold value
00133         }
00134         lastP  =colP [i];              // Last Momentum  (A-dependent)
00135         lastCS =colCS[i];              // Last CrossSect (A-dependent)
00136         in = true;                     // This is the case when the isotop is found in DB
00137         // Momentum pMom is in IU ! @@ Units
00138         lastCS=CalculateCrossSection(-1,j,2212,lastZ,lastN,pMom); // read & update
00139         if(lastCS<=0. && pMom>lastTH)  // Correct the threshold (@@ No intermediate Zeros)
00140         {
00141           lastCS=0.;
00142           lastTH=pMom;
00143         }
00144         break;                         // Go out of the LOOP
00145       }
00146       j++;                             // Increment a#0f records found in DB
00147     }
00148     if(!in)                            // This isotope has not been calculated previously
00149     {
00151       lastCS=CalculateCrossSection(0,j,2212,lastZ,lastN,pMom); //calculate & create
00152       //if(lastCS>0.)                   // It means that the AMBD was initialized
00153       //{
00154 
00155       lastTH = 0; //ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
00156         colN.push_back(tgN);
00157         colZ.push_back(tgZ);
00158         colP.push_back(pMom);
00159         colTH.push_back(lastTH);
00160         colCS.push_back(lastCS);
00161         //} // M.K. Presence of H1 with high threshold breaks the syncronization
00162       return lastCS*millibarn;
00163     } // End of creation of the new set of parameters
00164     else
00165     {
00166       colP[lastI]=pMom;
00167       colCS[lastI]=lastCS;
00168     }
00169   } // End of parameters udate
00170   else if(pMom<=lastTH)
00171   {
00172     return 0.;                         // Momentum is below the Threshold Value -> CS=0
00173   }
00174   else                                 // It is the last used -> use the current tables
00175   {
00176     lastCS=CalculateCrossSection(1,j,2212,lastZ,lastN,pMom); // Only read and UpdateDB
00177     lastP=pMom;
00178   }
00179   return lastCS*millibarn;
00180 }
00181 
00182 // The main member function giving the gamma-A cross section (E in GeV, CS in mb)
00183 G4double G4ChipsProtonInelasticXS::CalculateCrossSection(G4int F, G4int I,
00184                                         G4int, G4int targZ, G4int targN, G4double Momentum)
00185 {
00186   static const G4double THmin=27.;     // default minimum Momentum (MeV/c) Threshold
00187   static const G4double THmiG=THmin*.001; // minimum Momentum (GeV/c) Threshold
00188   static const G4double dP=10.;        // step for the LEN (Low ENergy) table MeV/c
00189   static const G4double dPG=dP*.001;   // step for the LEN (Low ENergy) table GeV/c
00190   static const G4int    nL=105;        // A#of LEN points in E (step 10 MeV/c)
00191   static const G4double Pmin=THmin+(nL-1)*dP; // minP for the HighE part with safety
00192   static const G4double Pmax=227000.;  // maxP for the HEN (High ENergy) part 227 GeV
00193   static const G4int    nH=224;        // A#of HEN points in lnE
00194   static const G4double milP=std::log(Pmin);// Low logarithm energy for the HEN part
00195   static const G4double malP=std::log(Pmax);// High logarithm energy (each 2.75 percent)
00196   static const G4double dlP=(malP-milP)/(nH-1); // Step in log energy in the HEN part
00197   static const G4double milPG=std::log(.001*Pmin);// Low logarithmEnergy for HEN part GeV/c
00198   G4double sigma=0.;
00199   if(F&&I) sigma=0.;                   // @@ *!* Fake line *!* to use F & I !!!Temporary!!!
00200   //G4double A=targN+targZ;              // A of the target
00201   if(F<=0)                             // This isotope was not the last used isotop
00202   {
00203     if(F<0)                            // This isotope was found in DAMDB =-----=> RETRIEVE
00204     {
00205       G4int sync=LEN->size();
00206       if(sync<=I) G4cout<<"*!*G4QProtonNuclCS::CalcCrossSect:Sync="<<sync<<"<="<<I<<G4endl;
00207       lastLEN=(*LEN)[I];               // Pointer to prepared LowEnergy cross sections
00208       lastHEN=(*HEN)[I];               // Pointer to prepared High Energy cross sections
00209     }
00210     else                               // This isotope wasn't calculated before => CREATE
00211     {
00212       lastLEN = new G4double[nL];      // Allocate memory for the new LEN cross sections
00213       lastHEN = new G4double[nH];      // Allocate memory for the new HEN cross sections
00214       // --- Instead of making a separate function ---
00215       G4double P=THmiG;                // Table threshold in GeV/c
00216       for(G4int k=0; k<nL; k++)
00217       {
00218         lastLEN[k] = CrossSectionLin(targZ, targN, P);
00219         P+=dPG;
00220       }
00221       G4double lP=milPG;
00222       for(G4int n=0; n<nH; n++)
00223       {
00224         lastHEN[n] = CrossSectionLog(targZ, targN, lP);
00225         lP+=dlP;
00226       }
00227       // --- End of possible separate function
00228       // *** The synchronization check ***
00229       G4int sync=LEN->size();
00230       if(sync!=I)
00231       {
00232         G4cout<<"***G4ChipsProtonNuclCS::CalcCrossSect: Sinc="<<sync<<"#"<<I<<", Z=" <<targZ
00233               <<", N="<<targN<<", F="<<F<<G4endl;
00234         //G4Exception("G4ProtonNuclearCS::CalculateCS:","39",FatalException,"overflow DB");
00235       }
00236       LEN->push_back(lastLEN);          // remember the Low Energy Table
00237       HEN->push_back(lastHEN);          // remember the High Energy Table
00238     } // End of creation of the new set of parameters
00239   } // End of parameters udate
00240   // =------------------= NOW the Magic Formula =-----------------------=
00241   if (Momentum<lastTH) return 0.;      // It must be already checked in the interface class
00242   else if (Momentum<Pmin)              // High Energy region
00243   {
00244     sigma=EquLinearFit(Momentum,nL,THmin,dP,lastLEN);
00245   }
00246   else if (Momentum<Pmax)              // High Energy region
00247   {
00248     G4double lP=std::log(Momentum);
00249     sigma=EquLinearFit(lP,nH,milP,dlP,lastHEN);
00250   }
00251   else                                 // UHE region (calculation, not frequent)
00252   {
00253     G4double P=0.001*Momentum;         // Approximation formula is for P in GeV/c
00254     sigma=CrossSectionFormula(targZ, targN, P, std::log(P));
00255   }
00256   if(sigma<0.) return 0.;
00257   return sigma;
00258 }
00259 
00260 // Electromagnetic momentum-threshold (in MeV/c) 
00261 G4double G4ChipsProtonInelasticXS::ThresholdMomentum(G4int tZ, G4int tN)
00262 {
00263   static const G4double third=1./3.;
00264   static const G4double pM = G4Proton::Proton()->Definition()->GetPDGMass(); // Projectile mass in MeV
00265   static const G4double tpM= pM+pM;       // Doubled projectile mass (MeV)
00266 
00267   G4double tA=tZ+tN;
00268   if(tZ<.99 || tN<0.) return 0.;
00269   else if(tZ==1 && tN==0) return 800.;    // A threshold on the free proton
00270   //G4double dE=1.263*tZ/(1.+std::pow(tA,third));
00271   G4double dE=tZ/(1.+std::pow(tA,third)); // Safety for diffused edge of the nucleus (QE)
00272   G4double tM=931.5*tA;
00273   G4double T=dE+dE*(dE/2+pM)/tM;
00274   return std::sqrt(T*(tpM+T));
00275 }
00276 
00277 // Calculation formula for proton-nuclear inelastic cross-section (mb) (P in GeV/c)
00278 G4double G4ChipsProtonInelasticXS::CrossSectionLin(G4int tZ, G4int tN, G4double P)
00279 {
00280   G4double sigma=0.;
00281   if(P<ThresholdMomentum(tZ,tN)*.001) return sigma;
00282   G4double lP=std::log(P);
00283   if(tZ==1&&!tN){if(P>.35) sigma=CrossSectionFormula(tZ,tN,P,lP);}// s(pp)=0 below 350Mev/c
00284   else if(tZ<97 && tN<152)                // General solution
00285   {
00286     G4double pex=0.;
00287     G4double pos=0.;
00288     G4double wid=1.;
00289     if(tZ==13 && tN==14)                  // Excited metastable states
00290     {
00291       pex=230.;
00292       pos=.13;
00293       wid=8.e-5;
00294     }
00295     else if(tZ<7)
00296     {
00297       if(tZ==6 && tN==6)
00298       {
00299         pex=320.;
00300         pos=.14;
00301         wid=7.e-6;
00302       }
00303       else if(tZ==5 && tN==6)
00304       {
00305         pex=270.;
00306         pos=.17;
00307         wid=.002;
00308       }
00309       else if(tZ==4 && tN==5)
00310       {
00311         pex=600.;
00312         pos=.132;
00313         wid=.005;
00314       }
00315       else if(tZ==3 && tN==4)
00316       {
00317         pex=280.;
00318         pos=.19;
00319         wid=.0025;
00320       }
00321       else if(tZ==3 && tN==3)
00322       {
00323         pex=370.;
00324         pos=.171;
00325         wid=.006;
00326       }
00327       else if(tZ==2 && tN==1)
00328       {
00329         pex=30.;
00330         pos=.22;
00331         wid=.0005;
00332       }
00333     }
00334     sigma=CrossSectionFormula(tZ,tN,P,lP);
00335     if(pex>0.)
00336     {
00337       G4double dp=P-pos;
00338       sigma+=pex*std::exp(-dp*dp/wid);
00339     }
00340   }
00341   else
00342   {
00343     G4cerr<<"-Warning-G4ChipsProtonNuclearXS::CSLin:*Bad A* Z="<<tZ<<", N="<<tN<<G4endl;
00344     sigma=0.;
00345   }
00346   if(sigma<0.) return 0.;
00347   return sigma;  
00348 }
00349 
00350 // Calculation formula for proton-nuclear inelastic cross-section (mb) log(P in GeV/c)
00351 G4double G4ChipsProtonInelasticXS::CrossSectionLog(G4int tZ, G4int tN, G4double lP)
00352 {
00353   G4double P=std::exp(lP);
00354   return CrossSectionFormula(tZ, tN, P, lP);
00355 }
00356 // Calculation formula for proton-nuclear inelastic cross-section (mb) log(P in GeV/c)
00357 G4double G4ChipsProtonInelasticXS::CrossSectionFormula(G4int tZ, G4int tN,
00358                                                            G4double P, G4double lP)
00359 {
00360   G4double sigma=0.;
00361   if(tZ==1 && !tN)                        // pp interaction (from G4QuasiElasticRatios)
00362   {
00363     G4double p2=P*P;
00364     G4double lp=lP-3.5;
00365     G4double lp2=lp*lp;
00366     G4double rp2=1./p2;
00367     G4double El=(.0557*lp2+6.72+30./P)/(1.+.49*rp2/P);
00368     G4double To=(.3*lp2+38.2)/(1.+.54*rp2*rp2);
00369     sigma=To-El;
00370   }
00371   else if(tZ<97 && tN<152)                // General solution
00372   {
00373     //G4double lP=std::log(P);            // Already calculated
00374     G4double d=lP-4.2;
00375     G4double p2=P*P;
00376     G4double p4=p2*p2;
00377     G4double a=tN+tZ;                       // A of the target
00378     G4double al=std::log(a);
00379     G4double sa=std::sqrt(a);
00380     G4double a2=a*a;
00381     G4double a2s=a2*sa;
00382     G4double a4=a2*a2;
00383     G4double a8=a4*a4;
00384     G4double a12=a8*a4;
00385     G4double a16=a8*a8;
00386     G4double c=(170.+3600./a2s)/(1.+65./a2s);
00387     G4double dl=al-3.;
00388     G4double dl2=dl*dl;
00389     G4double r=.21+.62*dl2/(1.+.5*dl2);
00390     G4double gg=40.*std::exp(al*0.712)/(1.+12.2/a)/(1.+34./a2);
00391     G4double e=318.+a4/(1.+.0015*a4/std::exp(al*0.09))/(1.+4.e-28*a12)+
00392                8.e-18/(1./a16+1.3e-20)/(1.+1.e-21*a12);
00393     G4double ss=3.57+.009*a2/(1.+.0001*a2*a);
00394     G4double h=(.01/a4+2.5e-6/a)*(1.+6.e-6*a2*a)/(1.+6.e7/a12/a2);
00395     sigma=(c+d*d)/(1.+r/p4)+(gg+e*std::exp(-ss*P))/(1.+h/p4/p4);
00396   }
00397   else
00398   {
00399     G4cerr<<"-Warning-G4QProtonNuclearCroSect::CSForm:*Bad A* Z="<<tZ<<", N="<<tN<<G4endl;
00400     sigma=0.;
00401   }
00402   if(sigma<0.) return 0.;
00403   return sigma;  
00404 }
00405 
00406 G4double G4ChipsProtonInelasticXS::EquLinearFit(G4double X, G4int N, G4double X0, G4double DX, G4double* Y)
00407 {
00408   if(DX<=0. || N<2)
00409     {
00410       G4cerr<<"***G4ChipsProtonInelasticXS::EquLinearFit: DX="<<DX<<", N="<<N<<G4endl;
00411       return Y[0];
00412     }
00413   
00414   G4int    N2=N-2;
00415   G4double d=(X-X0)/DX;
00416   G4int         j=static_cast<int>(d);
00417   if     (j<0)  j=0;
00418   else if(j>N2) j=N2;
00419   d-=j; // excess
00420   G4double yi=Y[j];
00421   G4double sigma=yi+(Y[j+1]-yi)*d;
00422   
00423   return sigma;
00424 }

Generated on Mon May 27 17:47:53 2013 for Geant4 by  doxygen 1.4.7