G4ComponentGGHadronNucleusXsc.hh

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 // Calculation of the total, elastic and inelastic cross-sections
00027 // based on parametrisations of (proton, pion, kaon, photon) nucleon
00028 // cross-sections and the hadron-nucleous cross-section model in 
00029 // the framework of Glauber-Gribov approach
00030 //
00031 //
00032 //
00033 //
00034 //
00035 // 25.04.12 V. Grichine - first implementation based on G4GlauberGribovCrossSection old interface
00036 //
00037 //
00038 
00039 #ifndef G4ComponentGGHadronNucleusXsc_h
00040 #define G4ComponentGGHadronNucleusXsc_h 1
00041 
00042 #include "globals.hh"
00043 #include "G4Proton.hh"
00044 #include "G4Nucleus.hh"
00045 
00046 #include "G4VComponentCrossSection.hh"
00047 
00048 class G4ParticleDefinition;
00049 class G4HadronNucleonXsc;
00050 
00051 class G4ComponentGGHadronNucleusXsc : public G4VComponentCrossSection
00052 {
00053 public:
00054 
00055   G4ComponentGGHadronNucleusXsc ();
00056   virtual ~G4ComponentGGHadronNucleusXsc ();
00057 
00058 
00059   // virtual interface methods
00060 
00061 virtual
00062   G4double GetTotalIsotopeCrossSection(const G4ParticleDefinition* aParticle,
00063                                        G4double kinEnergy,
00064                                        G4int Z, G4int A);
00065 
00066 virtual
00067   G4double GetTotalElementCrossSection(const G4ParticleDefinition* aParticle,
00068                                        G4double kinEnergy, 
00069                                        G4int Z, G4double A);
00070 
00071 virtual
00072   G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition* aParticle,
00073                                            G4double kinEnergy, 
00074                                            G4int Z, G4int A);
00075 
00076 virtual
00077   G4double GetInelasticElementCrossSection(const G4ParticleDefinition* aParticle,
00078                                            G4double kinEnergy, 
00079                                            G4int Z, G4double A);
00080 
00081 virtual
00082   G4double GetElasticElementCrossSection(const G4ParticleDefinition* aParticle,
00083                                          G4double kinEnergy, 
00084                                          G4int Z, G4double A);
00085 
00086 virtual
00087   G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition* aParticle,
00088                                          G4double kinEnergy, 
00089                                          G4int Z, G4int A);
00090  
00091 virtual
00092   G4double ComputeQuasiElasticRatio(const G4ParticleDefinition* aParticle,
00093                                          G4double kinEnergy, 
00094                                          G4int Z, G4int A);
00095  
00096 
00097    
00098   //  virtual
00099   G4bool IsIsoApplicable(const G4DynamicParticle* aDP, G4int Z, G4int A, 
00100                          const G4Element* elm = 0,
00101                          const G4Material* mat = 0);
00102 
00103   //  virtual
00104   G4double GetIsoCrossSection(const G4DynamicParticle*, G4int Z, G4int A,  
00105                               const G4Isotope* iso = 0,
00106                               const G4Element* elm = 0,
00107                               const G4Material* mat = 0);
00108 
00109   G4double GetRatioSD(const G4DynamicParticle*, G4int At, G4int Zt);
00110   G4double GetRatioQE(const G4DynamicParticle*, G4int At, G4int Zt);
00111 
00112   G4double GetHadronNucleonXsc(const G4DynamicParticle*, const G4Element*);
00113   G4double GetHadronNucleonXsc(const G4DynamicParticle*, G4int At, G4int Zt);
00114 
00115   G4double GetHadronNucleonXscPDG(const G4DynamicParticle*, const G4Element*);
00116   G4double GetHadronNucleonXscPDG(const G4DynamicParticle*, G4int At, G4int Zt);
00117   G4double GetHadronNucleonXscNS(const G4DynamicParticle*, const G4Element*);
00118   G4double GetHadronNucleonXscNS(const G4DynamicParticle*, G4int At, G4int Zt);
00119   G4double GetKaonNucleonXscVector(const G4DynamicParticle*, G4int At, G4int Zt);
00120 
00121   G4double GetHNinelasticXsc(const G4DynamicParticle*, const G4Element*);
00122   G4double GetHNinelasticXsc(const G4DynamicParticle*, G4int At, G4int Zt);
00123   G4double GetHNinelasticXscVU(const G4DynamicParticle*, G4int At, G4int Zt);
00124 
00125   G4double CalculateEcmValue ( const G4double , const G4double , const G4double ); 
00126 
00127   G4double CalcMandelstamS( const G4double , const G4double , const G4double );
00128 
00129   G4double GetNucleusRadius(const G4DynamicParticle*, const G4Element*);
00130   G4double GetNucleusRadius(G4int At);
00131 
00132   virtual void CrossSectionDescription(std::ostream&) const;
00133 
00134   inline G4double GetElasticGlauberGribov(const G4DynamicParticle*, G4int Z, G4int A);
00135   inline G4double GetInelasticGlauberGribov(const G4DynamicParticle*, G4int Z, G4int A);
00136 
00137   inline G4double GetTotalGlauberGribovXsc()    { return fTotalXsc;     }; 
00138   inline G4double GetElasticGlauberGribovXsc()  { return fElasticXsc;   }; 
00139   inline G4double GetInelasticGlauberGribovXsc(){ return fInelasticXsc; }; 
00140   inline G4double GetProductionGlauberGribovXsc(){ return fProductionXsc; }; 
00141   inline G4double GetDiffractionGlauberGribovXsc(){ return fDiffractionXsc; }; 
00142   inline G4double GetRadiusConst()              { return fRadiusConst;  }; 
00143 
00144   inline G4double GetParticleBarCorTot(const G4ParticleDefinition* theParticle, G4int Z);
00145   inline G4double GetParticleBarCorIn(const G4ParticleDefinition* theParticle, G4int Z);
00146 
00147   inline void SetEnergyLowerLimit(G4double E ){fLowerLimit=E;};
00148 
00149 private:
00150 
00151   const G4double fUpperLimit;
00152   G4double fLowerLimit; 
00153   const G4double fRadiusConst;
00154 
00155   static const G4double fNeutronBarCorrectionTot[93];
00156   static const G4double fNeutronBarCorrectionIn[93];
00157 
00158   static const G4double fProtonBarCorrectionTot[93];
00159   static const G4double fProtonBarCorrectionIn[93];
00160 
00161   static const G4double fPionPlusBarCorrectionTot[93];
00162   static const G4double fPionPlusBarCorrectionIn[93];
00163 
00164   static const G4double fPionMinusBarCorrectionTot[93];
00165   static const G4double fPionMinusBarCorrectionIn[93];
00166 
00167   G4double fTotalXsc, fElasticXsc, fInelasticXsc, fProductionXsc, fDiffractionXsc;
00168   G4double fHadronNucleonXsc;
00169  
00170   G4ParticleDefinition* theGamma;
00171   G4ParticleDefinition* theProton;
00172   G4ParticleDefinition* theNeutron;
00173   G4ParticleDefinition* theAProton;
00174   G4ParticleDefinition* theANeutron;
00175   G4ParticleDefinition* thePiPlus;
00176   G4ParticleDefinition* thePiMinus;
00177   G4ParticleDefinition* thePiZero;
00178   G4ParticleDefinition* theKPlus;
00179   G4ParticleDefinition* theKMinus;
00180   G4ParticleDefinition* theK0S;
00181   G4ParticleDefinition* theK0L;
00182   G4ParticleDefinition* theL;
00183   G4ParticleDefinition* theAntiL;
00184   G4ParticleDefinition* theSPlus;
00185   G4ParticleDefinition* theASPlus;
00186   G4ParticleDefinition* theSMinus;
00187   G4ParticleDefinition* theASMinus;
00188   G4ParticleDefinition* theS0;
00189   G4ParticleDefinition* theAS0;
00190   G4ParticleDefinition* theXiMinus;
00191   G4ParticleDefinition* theXi0;
00192   G4ParticleDefinition* theAXiMinus;
00193   G4ParticleDefinition* theAXi0;
00194   G4ParticleDefinition* theOmega;
00195   G4ParticleDefinition* theAOmega;
00196   G4ParticleDefinition* theD;
00197   G4ParticleDefinition* theT;
00198   G4ParticleDefinition* theA;
00199   G4ParticleDefinition* theHe3;
00200 
00201   G4HadronNucleonXsc* hnXsc;
00202 
00203 };
00204 
00206 //
00207 // Inlines
00208 
00209 inline
00210 G4double
00211 G4ComponentGGHadronNucleusXsc::GetElasticGlauberGribov(const G4DynamicParticle* dp,
00212                                                      G4int Z, G4int A)
00213 {
00214   GetIsoCrossSection(dp, Z, A);
00215   return fElasticXsc;
00216 }
00217 
00219 
00220 inline
00221 G4double
00222 G4ComponentGGHadronNucleusXsc::GetInelasticGlauberGribov(const G4DynamicParticle* dp,
00223                                                        G4int Z, G4int A)
00224 {
00225   GetIsoCrossSection(dp, Z, A);
00226   return fInelasticXsc;
00227 }
00228 
00230 //
00231 // return correction at Tkin = 90*GeV GG -> Barashenkov tot xsc, when it 
00232 // is available, else return 1.0
00233 
00234 
00235 inline G4double G4ComponentGGHadronNucleusXsc::GetParticleBarCorTot( 
00236                           const G4ParticleDefinition* theParticle, G4int Z)
00237 {
00238   if(Z >= 2 && Z <= 92)
00239   {
00240     if(      theParticle == theProton ) return fProtonBarCorrectionTot[Z]; 
00241     else if( theParticle == theNeutron) return fNeutronBarCorrectionTot[Z]; 
00242     else if( theParticle == thePiPlus ) return fPionPlusBarCorrectionTot[Z];
00243     else if( theParticle == thePiMinus) return fPionMinusBarCorrectionTot[Z];
00244     else return 1.0;
00245   }
00246   else return 1.0;
00247 }
00248 
00250 //
00251 // return correction at Tkin = 90*GeV GG -> Barashenkov in xsc, when it 
00252 // is available, else return 1.0
00253 
00254 
00255 inline G4double G4ComponentGGHadronNucleusXsc::GetParticleBarCorIn( 
00256                           const G4ParticleDefinition* theParticle, G4int Z)
00257 {
00258   if(Z >= 2 && Z <= 92)
00259   {
00260     if(      theParticle == theProton ) return fProtonBarCorrectionIn[Z]; 
00261     else if( theParticle == theNeutron) return fNeutronBarCorrectionIn[Z]; 
00262     else if( theParticle == thePiPlus ) return fPionPlusBarCorrectionIn[Z];
00263     else if( theParticle == thePiMinus) return fPionMinusBarCorrectionIn[Z];
00264     else return 1.0;
00265   }
00266   else return 1.0;
00267 }
00268 
00269 #endif

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