Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:58:03

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 // Calculation of the total, elastic and inelastic cross-sections
0027 // based on parametrisations of (proton, pion, kaon, photon) nucleon
0028 // cross-sections and the hadron-nucleous cross-section model in 
0029 // the framework of Glauber-Gribov approach
0030 //
0031 // 25.04.12 V. Grichine - first implementation based on 
0032 //                        G4GlauberGribovCrossSection old interface
0033 //
0034 // 04.09.18 V. Ivantchenko Major revision of interfaces and implementation
0035 // 01.10.18 V. Grichine strange hyperon xsc
0036 // 27.05.19 V. Ivantchenko Removed obsolete methods and members 
0037 
0038 #ifndef G4ComponentGGHadronNucleusXsc_h
0039 #define G4ComponentGGHadronNucleusXsc_h 1
0040 
0041 #include "globals.hh"
0042 #include "G4Proton.hh"
0043 #include "G4Nucleus.hh"
0044 
0045 #include "G4VComponentCrossSection.hh"
0046 
0047 class G4ParticleDefinition;
0048 class G4HadronNucleonXsc;
0049 class G4Pow;
0050 
0051 class G4ComponentGGHadronNucleusXsc final : public G4VComponentCrossSection
0052 {
0053 public:
0054 
0055   explicit G4ComponentGGHadronNucleusXsc();
0056   ~G4ComponentGGHadronNucleusXsc() final;
0057 
0058   static const char* Default_Name() { return "Glauber-Gribov"; }
0059 
0060   // virtual interface methods
0061   G4double GetTotalElementCrossSection(const G4ParticleDefinition* aParticle,
0062                        G4double kinEnergy, 
0063                        G4int Z, G4double A) final;
0064 
0065   G4double GetTotalIsotopeCrossSection(const G4ParticleDefinition* aParticle,
0066                        G4double kinEnergy,
0067                        G4int Z, G4int A) final;
0068 
0069   G4double GetInelasticElementCrossSection(const G4ParticleDefinition* aParticle,
0070                        G4double kinEnergy, 
0071                        G4int Z, G4double A) final;
0072 
0073   G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition* aParticle,
0074                        G4double kinEnergy, 
0075                        G4int Z, G4int A) final;
0076 
0077   G4double GetElasticElementCrossSection(const G4ParticleDefinition* aParticle,
0078                      G4double kinEnergy, 
0079                      G4int Z, G4double A) final;
0080 
0081   G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition* aParticle,
0082                      G4double kinEnergy, 
0083                      G4int Z, G4int A) final;
0084  
0085   G4double ComputeQuasiElasticRatio(const G4ParticleDefinition* aParticle,
0086                     G4double kinEnergy, 
0087                     G4int Z, G4int A) final;
0088 
0089   // Glauber-Gribov cross section
0090   void ComputeCrossSections(const G4ParticleDefinition* aParticle,
0091                 G4double kinEnergy, G4int Z, G4int A, G4int nL = 0);
0092 
0093   // additional public methods
0094   G4double GetProductionElementCrossSection(const G4ParticleDefinition* aParticle,
0095                         G4double kinEnergy, 
0096                         G4int Z, G4double A);
0097 
0098   G4double GetProductionIsotopeCrossSection(const G4ParticleDefinition* aParticle,
0099                         G4double kinEnergy, 
0100                         G4int Z, G4int A);
0101 
0102   G4double GetRatioSD(const G4DynamicParticle*, G4int At, G4int Zt);
0103   G4double GetRatioQE(const G4DynamicParticle*, G4int At, G4int Zt);
0104 
0105   G4double GetHadronNucleonXsc(const G4DynamicParticle*, const G4Element*);
0106   G4double GetHadronNucleonXsc(const G4DynamicParticle*, G4int At, G4int Zt);
0107 
0108   G4double GetHadronNucleonXscPDG(const G4DynamicParticle*, const G4Element*);
0109   G4double GetHadronNucleonXscPDG(const G4DynamicParticle*, G4int At, G4int Zt);
0110   G4double GetHadronNucleonXscNS(const G4DynamicParticle*, const G4Element*);
0111   G4double GetHadronNucleonXscNS(const G4DynamicParticle*, G4int At, G4int Zt);
0112 
0113   G4double GetHNinelasticXsc(const G4DynamicParticle*, const G4Element*);
0114   G4double GetHNinelasticXsc(const G4DynamicParticle*, G4int At, G4int Zt);
0115   G4double GetHNinelasticXscVU(const G4DynamicParticle*, G4int At, G4int Zt);
0116 
0117   void Description(std::ostream&) const final;
0118 
0119   inline G4double GetIsoCrossSection(const G4DynamicParticle*, G4int Z, G4int A,  
0120                          const G4Isotope* iso = nullptr,
0121                          const G4Element* elm = nullptr,
0122                          const G4Material* mat = nullptr);
0123 
0124   inline G4double GetElasticGlauberGribov(const G4DynamicParticle*, G4int Z, G4int A);
0125   inline G4double GetInelasticGlauberGribov(const G4DynamicParticle*, G4int Z, G4int A);
0126 
0127   inline G4double GetAxsc2piR2() const                   { return fAxsc2piR2; };
0128   inline G4double GetModelInLog() const                  { return fModelInLog; };
0129   inline G4double GetTotalGlauberGribovXsc() const       { return fTotalXsc; }; 
0130   inline G4double GetElasticGlauberGribovXsc() const     { return fElasticXsc; }; 
0131   inline G4double GetInelasticGlauberGribovXsc() const   { return fInelasticXsc; }; 
0132   inline G4double GetProductionGlauberGribovXsc() const  { return fProductionXsc; }; 
0133   inline G4double GetDiffractionGlauberGribovXsc() const { return fDiffractionXsc; }; 
0134 
0135   inline G4double GetParticleBarCorTot(const G4ParticleDefinition* theParticle, G4int Z);
0136   inline G4double GetParticleBarCorIn(const G4ParticleDefinition* theParticle, G4int Z);
0137 
0138 private:
0139 
0140   static const G4double fNeutronBarCorrectionTot[93];
0141   static const G4double fNeutronBarCorrectionIn[93];
0142 
0143   static const G4double fProtonBarCorrectionTot[93];
0144   static const G4double fProtonBarCorrectionIn[93];
0145 
0146   static const G4double fPionPlusBarCorrectionTot[93];
0147   static const G4double fPionPlusBarCorrectionIn[93];
0148 
0149   static const G4double fPionMinusBarCorrectionTot[93];
0150   static const G4double fPionMinusBarCorrectionIn[93];
0151 
0152   G4double fTotalXsc, fElasticXsc, fInelasticXsc, fProductionXsc, fDiffractionXsc;
0153   G4double fAxsc2piR2, fModelInLog;
0154   G4double fEnergy; //Cache
0155  
0156   const G4ParticleDefinition* theGamma;
0157   const G4ParticleDefinition* theProton;
0158   const G4ParticleDefinition* theNeutron;
0159   const G4ParticleDefinition* theAProton;
0160   const G4ParticleDefinition* theANeutron;
0161   const G4ParticleDefinition* thePiPlus;
0162   const G4ParticleDefinition* thePiMinus;
0163   const G4ParticleDefinition* theKPlus;
0164   const G4ParticleDefinition* theKMinus;
0165   const G4ParticleDefinition* theK0S;
0166   const G4ParticleDefinition* theK0L;
0167   const G4ParticleDefinition* theLambda;
0168 
0169   G4HadronNucleonXsc* hnXsc;
0170 
0171   // Cache
0172   const G4ParticleDefinition* fParticle;
0173   G4int fZ, fA, fL;
0174 
0175 };
0176 
0177 ////////////////////////////////////////////////////////////////
0178 //
0179 // Inlines
0180 
0181 inline G4double 
0182 G4ComponentGGHadronNucleusXsc::GetIsoCrossSection(const G4DynamicParticle* dp, 
0183                           G4int Z, G4int A,  
0184                           const G4Isotope*,
0185                           const G4Element*,
0186                           const G4Material*)
0187 {
0188   ComputeCrossSections(dp->GetDefinition(), dp->GetKineticEnergy(), Z, A);
0189   return fTotalXsc;
0190 }
0191 
0192 inline G4double
0193 G4ComponentGGHadronNucleusXsc::GetElasticGlauberGribov(
0194            const G4DynamicParticle* dp, G4int Z, G4int A)
0195 {
0196   ComputeCrossSections(dp->GetDefinition(), dp->GetKineticEnergy(), Z, A);
0197   return fElasticXsc;
0198 }
0199 
0200 /////////////////////////////////////////////////////////////////
0201 
0202 inline G4double
0203 G4ComponentGGHadronNucleusXsc::GetInelasticGlauberGribov(
0204            const G4DynamicParticle* dp, G4int Z, G4int A)
0205 {
0206   ComputeCrossSections(dp->GetDefinition(), dp->GetKineticEnergy(), Z, A);
0207   return fInelasticXsc;
0208 }
0209 
0210 /////////////////////////////////////////////////////////////////////
0211 //
0212 // return correction at Tkin = 90*GeV GG -> Barashenkov tot xsc, when it 
0213 // is available, else return 1.0
0214 
0215 inline G4double G4ComponentGGHadronNucleusXsc::GetParticleBarCorTot( 
0216                            const G4ParticleDefinition* theParticle, G4int ZZ)
0217 {
0218   G4double cor = 1.0;
0219   G4int z = std::min(92, std::max(ZZ, 1));
0220   if(      theParticle == theProton ) cor = fProtonBarCorrectionTot[z]; 
0221   else if( theParticle == theNeutron) cor = fNeutronBarCorrectionTot[z]; 
0222   else if( theParticle == thePiPlus ) cor = fPionPlusBarCorrectionTot[z];
0223   else if( theParticle == thePiMinus) cor = fPionMinusBarCorrectionTot[z];
0224   return cor;
0225 }
0226 
0227 /////////////////////////////////////////////////////////////////////
0228 //
0229 // return correction at Tkin = 90*GeV GG -> Barashenkov in xsc, when it 
0230 // is available, else return 1.0
0231 
0232 
0233 inline G4double G4ComponentGGHadronNucleusXsc::GetParticleBarCorIn( 
0234                            const G4ParticleDefinition* theParticle, G4int ZZ)
0235 {
0236   G4double cor = 1.0;
0237   G4int z = std::min(92, std::max(ZZ, 1));
0238   if(      theParticle == theProton ) cor = fProtonBarCorrectionIn[z]; 
0239   else if( theParticle == theNeutron) cor = fNeutronBarCorrectionIn[z]; 
0240   else if( theParticle == thePiPlus ) cor = fPionPlusBarCorrectionIn[z];
0241   else if( theParticle == thePiMinus) cor = fPionMinusBarCorrectionIn[z];
0242   return cor;
0243 }
0244 
0245 #endif