G4ChipsAntiBaryonInelasticXS.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: G4ChipsAntiBaryonInelasticXS 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 // Short description: Cross-sections extracted (by W.Pokorski) from the CHIPS package for 
00037 // anti-baryoninteractions. Original author: M. Kossov
00038 // -------------------------------------------------------------------------------------
00039 //
00040 
00041 #include "G4ChipsAntiBaryonInelasticXS.hh"
00042 #include "G4SystemOfUnits.hh"
00043 #include "G4DynamicParticle.hh"
00044 #include "G4ParticleDefinition.hh"
00045 #include "G4AntiNeutron.hh"
00046 #include "G4AntiProton.hh"
00047 #include "G4AntiLambda.hh"
00048 #include "G4AntiSigmaPlus.hh"
00049 #include "G4AntiSigmaMinus.hh"
00050 #include "G4AntiSigmaZero.hh"
00051 #include "G4AntiXiMinus.hh"
00052 #include "G4AntiXiZero.hh"
00053 #include "G4AntiOmegaMinus.hh"
00054 
00055 // factory
00056 #include "G4CrossSectionFactory.hh"
00057 //
00058 G4_DECLARE_XS_FACTORY(G4ChipsAntiBaryonInelasticXS);
00059 
00060 G4ChipsAntiBaryonInelasticXS::G4ChipsAntiBaryonInelasticXS():G4VCrossSectionDataSet(Default_Name())
00061 {
00062   lastLEN=0; // Pointer to lastArray of LowEn CS
00063   lastHEN=0; // Pointer to lastArray of HighEn CS
00064   lastN=0;   // The last N of calculated nucleus
00065   lastZ=0;   // The last Z of calculated nucleus
00066   lastP=0.;  // Last used Cross Section Momentum
00067   lastTH=0.; // Last threshold momentum
00068   lastCS=0.; // Last value of the Cross Section
00069   lastI=0;   // The last position in the DAMDB
00070   LEN = new std::vector<G4double*>;
00071   HEN = new std::vector<G4double*>;
00072 }
00073 
00074 G4ChipsAntiBaryonInelasticXS::~G4ChipsAntiBaryonInelasticXS()
00075 {
00076   G4int lens=LEN->size();
00077   for(G4int i=0; i<lens; ++i) delete[] (*LEN)[i];
00078   delete LEN;
00079   G4int hens=HEN->size();
00080   for(G4int i=0; i<hens; ++i) delete[] (*HEN)[i];
00081   delete HEN;
00082 }
00083 
00084 G4bool G4ChipsAntiBaryonInelasticXS::IsIsoApplicable(const G4DynamicParticle* Pt, G4int, G4int,    
00085                                  const G4Element*,
00086                                  const G4Material*)
00087 {
00088   G4ParticleDefinition* particle = Pt->GetDefinition();
00089 
00090   if(particle == G4AntiNeutron::AntiNeutron())
00091   {
00092     return true;
00093   }
00094   else if(particle == G4AntiProton::AntiProton())
00095   {
00096     return true;
00097   }
00098   else if(particle == G4AntiLambda::AntiLambda())
00099   {
00100     return true;
00101   }
00102   else if(particle == G4AntiSigmaPlus::AntiSigmaPlus())
00103   {
00104     return true;
00105   }
00106   else if(particle == G4AntiSigmaMinus::AntiSigmaMinus())
00107   {
00108     return true;
00109   }
00110   else if(particle == G4AntiSigmaZero::AntiSigmaZero())
00111   {
00112     return true;
00113   }
00114   else if(particle == G4AntiXiMinus::AntiXiMinus())
00115   {
00116     return true;
00117   }
00118   else if(particle == G4AntiXiZero::AntiXiZero())
00119   {
00120     return true;
00121   }
00122   else if(particle == G4AntiOmegaMinus::AntiOmegaMinus())
00123   {
00124     return true;
00125   }
00126   return false;
00127 }
00128 
00129 // The main member function giving the collision cross section (P is in IU, CS is in mb)
00130 // Make pMom in independent units ! (Now it is MeV)
00131 G4double G4ChipsAntiBaryonInelasticXS::GetIsoCrossSection(const G4DynamicParticle* Pt, G4int tgZ, G4int A,  
00132                                                         const G4Isotope*,
00133                                                         const G4Element*,
00134                                                         const G4Material*)
00135 {
00136   G4double pMom=Pt->GetTotalMomentum();
00137   G4int tgN = A - tgZ;
00138   G4int pdg = Pt->GetDefinition()->GetPDGEncoding();
00139   
00140   return GetChipsCrossSection(pMom, tgZ, tgN, pdg);
00141 }
00142 
00143 G4double G4ChipsAntiBaryonInelasticXS::GetChipsCrossSection(G4double pMom, G4int tgZ, G4int tgN, G4int cPDG)
00144 {
00145   static G4int j;                      // A#0f Z/N-records already tested in AMDB
00146   static std::vector <G4int>    colN;  // Vector of N for calculated nuclei (isotops)
00147   static std::vector <G4int>    colZ;  // Vector of Z for calculated nuclei (isotops)
00148   static std::vector <G4double> colP;  // Vector of last momenta for the reaction
00149   static std::vector <G4double> colTH; // Vector of energy thresholds for the reaction
00150   static std::vector <G4double> colCS; // Vector of last cross sections for the reaction
00151   // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
00152 
00153   G4bool in=false;                     // By default the isotope must be found in the AMDB
00154   if(tgN!=lastN || tgZ!=lastZ)         // The nucleus was not the last used isotope
00155   {
00156     in = false;                        // By default the isotope haven't be found in AMDB  
00157     lastP   = 0.;                      // New momentum history (nothing to compare with)
00158     lastN   = tgN;                     // The last N of the calculated nucleus
00159     lastZ   = tgZ;                     // The last Z of the calculated nucleus
00160     lastI   = colN.size();             // Size of the Associative Memory DB in the heap
00161     j  = 0;                            // A#0f records found in DB for this projectile
00162     if(lastI) for(G4int i=0; i<lastI; i++) // AMDB exists, try to find the (Z,N) isotope
00163     {
00164       if(colN[i]==tgN && colZ[i]==tgZ) // Try the record "i" in the AMDB
00165       {
00166         lastI=i;                       // Remember the index for future fast/last use
00167         lastTH =colTH[i];              // The last THreshold (A-dependent)
00168         if(pMom<=lastTH)
00169         {
00170           return 0.;                   // Energy is below the Threshold value
00171         }
00172         lastP  =colP [i];              // Last Momentum  (A-dependent)
00173         lastCS =colCS[i];              // Last CrossSect (A-dependent)
00174         in = true;                     // This is the case when the isotop is found in DB
00175         // Momentum pMom is in IU ! @@ Units
00176         lastCS=CalculateCrossSection(-1,j,cPDG,lastZ,lastN,pMom); // read & update
00177         if(lastCS<=0. && pMom>lastTH)  // Correct the threshold (@@ No intermediate Zeros)
00178         {
00179           lastCS=0.;
00180           lastTH=pMom;
00181         }
00182         break;                         // Go out of the LOOP
00183       }
00184       j++;                             // Increment a#0f records found in DB
00185     }
00186     if(!in)                            // This isotope has not been calculated previously
00187     {
00189       lastCS=CalculateCrossSection(0,j,cPDG,lastZ,lastN,pMom); //calculate & create
00190       //if(lastCS>0.)                   // It means that the AMBD was initialized
00191       //{
00192 
00193       lastTH = 0; //ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
00194         colN.push_back(tgN);
00195         colZ.push_back(tgZ);
00196         colP.push_back(pMom);
00197         colTH.push_back(lastTH);
00198         colCS.push_back(lastCS);
00199       //} // M.K. Presence of H1 with high threshold breaks the syncronization
00200       return lastCS*millibarn;
00201     } // End of creation of the new set of parameters
00202     else
00203     {
00204       colP[lastI]=pMom;
00205       colCS[lastI]=lastCS;
00206     }
00207   } // End of parameters udate
00208   else if(pMom<=lastTH)
00209   {
00210     return 0.;                         // Momentum is below the Threshold Value -> CS=0
00211   }
00212   else                                 // It is the last used -> use the current tables
00213   {
00214     lastCS=CalculateCrossSection(1,j,cPDG,lastZ,lastN,pMom); // Only read and UpdateDB
00215     lastP=pMom;
00216   }
00217   return lastCS*millibarn;
00218 }
00219 
00220 // The main member function giving the gamma-A cross section (E in GeV, CS in mb)
00221 G4double G4ChipsAntiBaryonInelasticXS::CalculateCrossSection(G4int F, G4int I,
00222                                         G4int, G4int targZ, G4int targN, G4double Momentum)
00223 {
00224   static const G4double THmin=27.;     // default minimum Momentum (MeV/c) Threshold
00225   static const G4double THmiG=THmin*.001; // minimum Momentum (GeV/c) Threshold
00226   static const G4double dP=10.;        // step for the LEN (Low ENergy) table MeV/c
00227   static const G4double dPG=dP*.001;   // step for the LEN (Low ENergy) table GeV/c
00228   static const G4int    nL=105;        // A#of LEN points in E (step 10 MeV/c)
00229   static const G4double Pmin=THmin+(nL-1)*dP; // minP for the HighE part with safety
00230   static const G4double Pmax=227000.;  // maxP for the HEN (High ENergy) part 227 GeV
00231   static const G4int    nH=224;        // A#of HEN points in lnE
00232   static const G4double milP=std::log(Pmin);// Low logarithm energy for the HEN part
00233   static const G4double malP=std::log(Pmax);// High logarithm energy (each 2.75 percent)
00234   static const G4double dlP=(malP-milP)/(nH-1); // Step in log energy in the HEN part
00235   static const G4double milPG=std::log(.001*Pmin);// Low logarithmEnergy for HEN part GeV/c
00236   G4double sigma=0.;
00237   if(F&&I) sigma=0.;                   // @@ *!* Fake line *!* to use F & I !!!Temporary!!!
00238   //G4double A=targN+targZ;              // A of the target
00239   if(F<=0)                             // This isotope was not the last used isotop
00240   {
00241     if(F<0)                            // This isotope was found in DAMDB =-----=> RETRIEVE
00242     {
00243       G4int sync=LEN->size();
00244       if(sync<=I) G4cerr<<"*!*G4QPiMinusNuclCS::CalcCrosSect:Sync="<<sync<<"<="<<I<<G4endl;
00245       lastLEN=(*LEN)[I];               // Pointer to prepared LowEnergy cross sections
00246       lastHEN=(*HEN)[I];               // Pointer to prepared High Energy cross sections
00247     }
00248     else                               // This isotope wasn't calculated before => CREATE
00249     {
00250       lastLEN = new G4double[nL];      // Allocate memory for the new LEN cross sections
00251       lastHEN = new G4double[nH];      // Allocate memory for the new HEN cross sections
00252       // --- Instead of making a separate function ---
00253       G4double P=THmiG;                // Table threshold in GeV/c
00254       for(G4int k=0; k<nL; k++)
00255       {
00256         lastLEN[k] = CrossSectionLin(targZ, targN, P);
00257         P+=dPG;
00258       }
00259       G4double lP=milPG;
00260       for(G4int n=0; n<nH; n++)
00261       {
00262         lastHEN[n] = CrossSectionLog(targZ, targN, lP);
00263         lP+=dlP;
00264       }
00265       // --- End of possible separate function
00266       // *** The synchronization check ***
00267       G4int sync=LEN->size();
00268       if(sync!=I)
00269       {
00270         G4cerr<<"***G4QPiMinusNuclCS::CalcCrossSect: Sinc="<<sync<<"#"<<I<<", Z=" <<targZ
00271               <<", N="<<targN<<", F="<<F<<G4endl;
00272         //G4Exception("G4PiMinusNuclearCS::CalculateCS:","39",FatalException,"DBoverflow");
00273       }
00274       LEN->push_back(lastLEN);         // remember the Low Energy Table
00275       HEN->push_back(lastHEN);         // remember the High Energy Table
00276     } // End of creation of the new set of parameters
00277   } // End of parameters udate
00278   // =-------------------= NOW the Magic Formula =--------------------=
00279   if (Momentum<lastTH) return 0.;      // It must be already checked in the interface class
00280   else if (Momentum<Pmin)              // High Energy region
00281   {
00282     sigma=EquLinearFit(Momentum,nL,THmin,dP,lastLEN);
00283   }
00284   else if (Momentum<Pmax)              // High Energy region
00285   {
00286     G4double lP=std::log(Momentum);
00287     sigma=EquLinearFit(lP,nH,milP,dlP,lastHEN);
00288   }
00289   else                                 // UHE region (calculation, not frequent)
00290   {
00291     G4double P=0.001*Momentum;         // Approximation formula is for P in GeV/c
00292     sigma=CrossSectionFormula(targZ, targN, P, std::log(P));
00293   }
00294   if(sigma<0.) return 0.;
00295   return sigma;
00296 }
00297 
00298 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) (P in GeV/c)
00299 G4double G4ChipsAntiBaryonInelasticXS::CrossSectionLin(G4int tZ, G4int tN, G4double P)
00300 {
00301   G4double lP=std::log(P);
00302   return CrossSectionFormula(tZ, tN, P, lP);
00303 }
00304 
00305 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
00306 G4double G4ChipsAntiBaryonInelasticXS::CrossSectionLog(G4int tZ, G4int tN, G4double lP)
00307 {
00308   G4double P=std::exp(lP);
00309   return CrossSectionFormula(tZ, tN, P, lP);
00310 }
00311 // Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
00312 G4double G4ChipsAntiBaryonInelasticXS::CrossSectionFormula(G4int tZ, G4int tN,
00313                                                               G4double P, G4double lP)
00314 {
00315   G4double sigma=0.;
00316   if(tZ==1 && !tN)                        // AntiBar-Prot interaction from G4QuasiElRatios
00317   {
00318     G4double ld=lP-3.5;
00319     G4double ld2=ld*ld;
00320     G4double ye=std::exp(lP*1.25);
00321     G4double yt=std::exp(lP*0.35);
00322     G4double El=80./(ye+1.);
00323     G4double To=(80./yt+.3)/yt;
00324     sigma=(To-El)+.2443*ld2+31.48;
00325   }
00326   else if(tZ==1 && tN==1)
00327   {
00328     G4double r=lP-3.7;
00329     sigma=0.6*r*r+67.+90.*std::exp(-lP*.666);
00330   }
00331   else if(tZ<97 && tN<152)                // General solution
00332   {
00333     G4double d=lP-4.2;
00334     G4double sp=std::sqrt(P);
00335     G4double a=tN+tZ;                      // A of the target
00336     G4double sa=std::sqrt(a);
00337     G4double a2=a*a;
00338     G4double a3=a2*a;
00339     G4double a2s=a2*sa;
00340     G4double c=(170.+3600./a2s)/(1.+65./a2s)+40.*std::pow(a,0.712)/(1.+12.2/a)/(1.+34./a2);
00341     G4double r=(170.+0.01*a3)/(1.+a3/28000.);
00342     sigma=c+d*d+r/sp;
00343   }
00344   else
00345   {
00346     G4cerr<<"-Warning-G4QAntiBarNuclearCroSect::CSForm:*Bad A* Z="<<tZ<<", N="<<tN<<G4endl;
00347     sigma=0.;
00348   }
00349   if(sigma<0.) return 0.;
00350   return sigma;  
00351 }
00352 
00353 G4double G4ChipsAntiBaryonInelasticXS::EquLinearFit(G4double X, G4int N, G4double X0, G4double DX, G4double* Y)
00354 {
00355   if(DX<=0. || N<2)
00356     {
00357       G4cerr<<"***G4ChipsAntiBaryonInelasticXS::EquLinearFit: DX="<<DX<<", N="<<N<<G4endl;
00358       return Y[0];
00359     }
00360   
00361   G4int    N2=N-2;
00362   G4double d=(X-X0)/DX;
00363   G4int         j=static_cast<int>(d);
00364   if     (j<0)  j=0;
00365   else if(j>N2) j=N2;
00366   d-=j; // excess
00367   G4double yi=Y[j];
00368   G4double sigma=yi+(Y[j+1]-yi)*d;
00369   
00370   return sigma;
00371 }

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