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

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