Geant4 Cross Reference |
1 // 1 // 2 // ******************************************* 2 // ******************************************************************** 3 // * License and Disclaimer 3 // * License and Disclaimer * 4 // * 4 // * * 5 // * The Geant4 software is copyright of th 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/ 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. 9 // * include a list of copyright holders. * 10 // * 10 // * * 11 // * Neither the authors of this software syst 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing fin 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warran 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assum 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitatio 16 // * for the full disclaimer and the limitation of liability. * 17 // * 17 // * * 18 // * This code implementation is the result 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboratio 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distri 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you ag 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publicati 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Sof 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************* 24 // ******************************************************************** 25 // 25 // 26 // Calculation of the total, elastic and inela 26 // Calculation of the total, elastic and inelastic cross-sections 27 // based on parametrisations of (proton, pion, 27 // based on parametrisations of (proton, pion, kaon, photon) nucleon 28 // cross-sections and the hadron-nucleous cros 28 // cross-sections and the hadron-nucleous cross-section model in 29 // the framework of Glauber-Gribov approach 29 // the framework of Glauber-Gribov approach 30 // 30 // 31 // 25.04.12 V. Grichine - first implementation 31 // 25.04.12 V. Grichine - first implementation based on 32 // G4GlauberGribovCross 32 // G4GlauberGribovCrossSection old interface 33 // 33 // 34 // 04.09.18 V. Ivantchenko Major revision of i 34 // 04.09.18 V. Ivantchenko Major revision of interfaces and implementation 35 // 01.10.18 V. Grichine strange hyperon xsc 35 // 01.10.18 V. Grichine strange hyperon xsc 36 // 27.05.19 V. Ivantchenko Removed obsolete me 36 // 27.05.19 V. Ivantchenko Removed obsolete methods and members 37 37 38 #ifndef G4ComponentGGHadronNucleusXsc_h 38 #ifndef G4ComponentGGHadronNucleusXsc_h 39 #define G4ComponentGGHadronNucleusXsc_h 1 39 #define G4ComponentGGHadronNucleusXsc_h 1 40 40 41 #include "globals.hh" 41 #include "globals.hh" 42 #include "G4Proton.hh" 42 #include "G4Proton.hh" 43 #include "G4Nucleus.hh" 43 #include "G4Nucleus.hh" 44 44 45 #include "G4VComponentCrossSection.hh" 45 #include "G4VComponentCrossSection.hh" 46 46 47 class G4ParticleDefinition; 47 class G4ParticleDefinition; 48 class G4HadronNucleonXsc; 48 class G4HadronNucleonXsc; 49 class G4Pow; 49 class G4Pow; 50 50 51 class G4ComponentGGHadronNucleusXsc final : pu 51 class G4ComponentGGHadronNucleusXsc final : public G4VComponentCrossSection 52 { 52 { 53 public: 53 public: 54 54 55 explicit G4ComponentGGHadronNucleusXsc(); 55 explicit G4ComponentGGHadronNucleusXsc(); 56 ~G4ComponentGGHadronNucleusXsc() final; 56 ~G4ComponentGGHadronNucleusXsc() final; 57 57 58 static const char* Default_Name() { return " 58 static const char* Default_Name() { return "Glauber-Gribov"; } 59 59 60 // virtual interface methods 60 // virtual interface methods 61 G4double GetTotalElementCrossSection(const G 61 G4double GetTotalElementCrossSection(const G4ParticleDefinition* aParticle, 62 G4double kinEnergy, 62 G4double kinEnergy, 63 G4int Z, G4double A) final; 63 G4int Z, G4double A) final; 64 64 65 G4double GetTotalIsotopeCrossSection(const G 65 G4double GetTotalIsotopeCrossSection(const G4ParticleDefinition* aParticle, 66 G4double kinEnergy, 66 G4double kinEnergy, 67 G4int Z, G4int A) final; 67 G4int Z, G4int A) final; 68 68 69 G4double GetInelasticElementCrossSection(con 69 G4double GetInelasticElementCrossSection(const G4ParticleDefinition* aParticle, 70 G4double kinEnergy, 70 G4double kinEnergy, 71 G4int Z, G4double A) final; 71 G4int Z, G4double A) final; 72 72 73 G4double GetInelasticIsotopeCrossSection(con 73 G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition* aParticle, 74 G4double kinEnergy, 74 G4double kinEnergy, 75 G4int Z, G4int A) final; 75 G4int Z, G4int A) final; 76 76 77 G4double GetElasticElementCrossSection(const 77 G4double GetElasticElementCrossSection(const G4ParticleDefinition* aParticle, 78 G4double kinEnergy, 78 G4double kinEnergy, 79 G4int Z, G4double A) final; 79 G4int Z, G4double A) final; 80 80 81 G4double GetElasticIsotopeCrossSection(const 81 G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition* aParticle, 82 G4double kinEnergy, 82 G4double kinEnergy, 83 G4int Z, G4int A) final; 83 G4int Z, G4int A) final; 84 84 85 G4double ComputeQuasiElasticRatio(const G4Pa 85 G4double ComputeQuasiElasticRatio(const G4ParticleDefinition* aParticle, 86 G4double kinEnergy, 86 G4double kinEnergy, 87 G4int Z, G4int A) final; 87 G4int Z, G4int A) final; 88 88 89 // Glauber-Gribov cross section 89 // Glauber-Gribov cross section 90 void ComputeCrossSections(const G4ParticleDe 90 void ComputeCrossSections(const G4ParticleDefinition* aParticle, 91 G4double kinEnergy, G4int Z, G4int A 91 G4double kinEnergy, G4int Z, G4int A, G4int nL = 0); 92 92 93 // additional public methods 93 // additional public methods 94 G4double GetProductionElementCrossSection(co 94 G4double GetProductionElementCrossSection(const G4ParticleDefinition* aParticle, 95 G4double kinEnergy, 95 G4double kinEnergy, 96 G4int Z, G4double A); 96 G4int Z, G4double A); 97 97 98 G4double GetProductionIsotopeCrossSection(co 98 G4double GetProductionIsotopeCrossSection(const G4ParticleDefinition* aParticle, 99 G4double kinEnergy, 99 G4double kinEnergy, 100 G4int Z, G4int A); 100 G4int Z, G4int A); 101 101 102 G4double GetRatioSD(const G4DynamicParticle* 102 G4double GetRatioSD(const G4DynamicParticle*, G4int At, G4int Zt); 103 G4double GetRatioQE(const G4DynamicParticle* 103 G4double GetRatioQE(const G4DynamicParticle*, G4int At, G4int Zt); 104 104 105 G4double GetHadronNucleonXsc(const G4Dynamic 105 G4double GetHadronNucleonXsc(const G4DynamicParticle*, const G4Element*); 106 G4double GetHadronNucleonXsc(const G4Dynamic 106 G4double GetHadronNucleonXsc(const G4DynamicParticle*, G4int At, G4int Zt); 107 107 108 G4double GetHadronNucleonXscPDG(const G4Dyna 108 G4double GetHadronNucleonXscPDG(const G4DynamicParticle*, const G4Element*); 109 G4double GetHadronNucleonXscPDG(const G4Dyna 109 G4double GetHadronNucleonXscPDG(const G4DynamicParticle*, G4int At, G4int Zt); 110 G4double GetHadronNucleonXscNS(const G4Dynam 110 G4double GetHadronNucleonXscNS(const G4DynamicParticle*, const G4Element*); 111 G4double GetHadronNucleonXscNS(const G4Dynam 111 G4double GetHadronNucleonXscNS(const G4DynamicParticle*, G4int At, G4int Zt); 112 112 113 G4double GetHNinelasticXsc(const G4DynamicPa 113 G4double GetHNinelasticXsc(const G4DynamicParticle*, const G4Element*); 114 G4double GetHNinelasticXsc(const G4DynamicPa 114 G4double GetHNinelasticXsc(const G4DynamicParticle*, G4int At, G4int Zt); 115 G4double GetHNinelasticXscVU(const G4Dynamic 115 G4double GetHNinelasticXscVU(const G4DynamicParticle*, G4int At, G4int Zt); 116 116 117 void Description(std::ostream&) const final; 117 void Description(std::ostream&) const final; 118 118 119 inline G4double GetIsoCrossSection(const G4D 119 inline G4double GetIsoCrossSection(const G4DynamicParticle*, G4int Z, G4int A, 120 const G4Isotope* iso = null 120 const G4Isotope* iso = nullptr, 121 const G4Element* elm = null 121 const G4Element* elm = nullptr, 122 const G4Material* mat = nul 122 const G4Material* mat = nullptr); 123 123 124 inline G4double GetElasticGlauberGribov(cons 124 inline G4double GetElasticGlauberGribov(const G4DynamicParticle*, G4int Z, G4int A); 125 inline G4double GetInelasticGlauberGribov(co 125 inline G4double GetInelasticGlauberGribov(const G4DynamicParticle*, G4int Z, G4int A); 126 126 127 inline G4double GetAxsc2piR2() const 127 inline G4double GetAxsc2piR2() const { return fAxsc2piR2; }; 128 inline G4double GetModelInLog() const 128 inline G4double GetModelInLog() const { return fModelInLog; }; 129 inline G4double GetTotalGlauberGribovXsc() c 129 inline G4double GetTotalGlauberGribovXsc() const { return fTotalXsc; }; 130 inline G4double GetElasticGlauberGribovXsc() 130 inline G4double GetElasticGlauberGribovXsc() const { return fElasticXsc; }; 131 inline G4double GetInelasticGlauberGribovXsc 131 inline G4double GetInelasticGlauberGribovXsc() const { return fInelasticXsc; }; 132 inline G4double GetProductionGlauberGribovXs 132 inline G4double GetProductionGlauberGribovXsc() const { return fProductionXsc; }; 133 inline G4double GetDiffractionGlauberGribovX 133 inline G4double GetDiffractionGlauberGribovXsc() const { return fDiffractionXsc; }; 134 134 135 inline G4double GetParticleBarCorTot(const G 135 inline G4double GetParticleBarCorTot(const G4ParticleDefinition* theParticle, G4int Z); 136 inline G4double GetParticleBarCorIn(const G4 136 inline G4double GetParticleBarCorIn(const G4ParticleDefinition* theParticle, G4int Z); 137 137 138 private: 138 private: 139 139 140 static const G4double fNeutronBarCorrectionT 140 static const G4double fNeutronBarCorrectionTot[93]; 141 static const G4double fNeutronBarCorrectionI 141 static const G4double fNeutronBarCorrectionIn[93]; 142 142 143 static const G4double fProtonBarCorrectionTo 143 static const G4double fProtonBarCorrectionTot[93]; 144 static const G4double fProtonBarCorrectionIn 144 static const G4double fProtonBarCorrectionIn[93]; 145 145 146 static const G4double fPionPlusBarCorrection 146 static const G4double fPionPlusBarCorrectionTot[93]; 147 static const G4double fPionPlusBarCorrection 147 static const G4double fPionPlusBarCorrectionIn[93]; 148 148 149 static const G4double fPionMinusBarCorrectio 149 static const G4double fPionMinusBarCorrectionTot[93]; 150 static const G4double fPionMinusBarCorrectio 150 static const G4double fPionMinusBarCorrectionIn[93]; 151 151 152 G4double fTotalXsc, fElasticXsc, fInelasticX 152 G4double fTotalXsc, fElasticXsc, fInelasticXsc, fProductionXsc, fDiffractionXsc; 153 G4double fAxsc2piR2, fModelInLog; 153 G4double fAxsc2piR2, fModelInLog; 154 G4double fEnergy; //Cache 154 G4double fEnergy; //Cache 155 155 156 const G4ParticleDefinition* theGamma; 156 const G4ParticleDefinition* theGamma; 157 const G4ParticleDefinition* theProton; 157 const G4ParticleDefinition* theProton; 158 const G4ParticleDefinition* theNeutron; 158 const G4ParticleDefinition* theNeutron; 159 const G4ParticleDefinition* theAProton; 159 const G4ParticleDefinition* theAProton; 160 const G4ParticleDefinition* theANeutron; 160 const G4ParticleDefinition* theANeutron; 161 const G4ParticleDefinition* thePiPlus; 161 const G4ParticleDefinition* thePiPlus; 162 const G4ParticleDefinition* thePiMinus; 162 const G4ParticleDefinition* thePiMinus; 163 const G4ParticleDefinition* theKPlus; 163 const G4ParticleDefinition* theKPlus; 164 const G4ParticleDefinition* theKMinus; 164 const G4ParticleDefinition* theKMinus; 165 const G4ParticleDefinition* theK0S; 165 const G4ParticleDefinition* theK0S; 166 const G4ParticleDefinition* theK0L; 166 const G4ParticleDefinition* theK0L; 167 const G4ParticleDefinition* theLambda; 167 const G4ParticleDefinition* theLambda; 168 168 169 G4HadronNucleonXsc* hnXsc; 169 G4HadronNucleonXsc* hnXsc; 170 170 171 // Cache 171 // Cache 172 const G4ParticleDefinition* fParticle; 172 const G4ParticleDefinition* fParticle; 173 G4int fZ, fA, fL; 173 G4int fZ, fA, fL; 174 174 175 }; 175 }; 176 176 177 ////////////////////////////////////////////// 177 //////////////////////////////////////////////////////////////// 178 // 178 // 179 // Inlines 179 // Inlines 180 180 181 inline G4double 181 inline G4double 182 G4ComponentGGHadronNucleusXsc::GetIsoCrossSect 182 G4ComponentGGHadronNucleusXsc::GetIsoCrossSection(const G4DynamicParticle* dp, 183 G4int Z, G4int A, 183 G4int Z, G4int A, 184 const G4Isotope*, 184 const G4Isotope*, 185 const G4Element*, 185 const G4Element*, 186 const G4Material*) 186 const G4Material*) 187 { 187 { 188 ComputeCrossSections(dp->GetDefinition(), dp 188 ComputeCrossSections(dp->GetDefinition(), dp->GetKineticEnergy(), Z, A); 189 return fTotalXsc; 189 return fTotalXsc; 190 } 190 } 191 191 192 inline G4double 192 inline G4double 193 G4ComponentGGHadronNucleusXsc::GetElasticGlaub 193 G4ComponentGGHadronNucleusXsc::GetElasticGlauberGribov( 194 const G4DynamicParticle* dp, G4int 194 const G4DynamicParticle* dp, G4int Z, G4int A) 195 { 195 { 196 ComputeCrossSections(dp->GetDefinition(), dp 196 ComputeCrossSections(dp->GetDefinition(), dp->GetKineticEnergy(), Z, A); 197 return fElasticXsc; 197 return fElasticXsc; 198 } 198 } 199 199 200 ////////////////////////////////////////////// 200 ///////////////////////////////////////////////////////////////// 201 201 202 inline G4double 202 inline G4double 203 G4ComponentGGHadronNucleusXsc::GetInelasticGla 203 G4ComponentGGHadronNucleusXsc::GetInelasticGlauberGribov( 204 const G4DynamicParticle* dp, G4int 204 const G4DynamicParticle* dp, G4int Z, G4int A) 205 { 205 { 206 ComputeCrossSections(dp->GetDefinition(), dp 206 ComputeCrossSections(dp->GetDefinition(), dp->GetKineticEnergy(), Z, A); 207 return fInelasticXsc; 207 return fInelasticXsc; 208 } 208 } 209 209 210 ////////////////////////////////////////////// 210 ///////////////////////////////////////////////////////////////////// 211 // 211 // 212 // return correction at Tkin = 90*GeV GG -> Ba 212 // return correction at Tkin = 90*GeV GG -> Barashenkov tot xsc, when it 213 // is available, else return 1.0 213 // is available, else return 1.0 214 214 215 inline G4double G4ComponentGGHadronNucleusXsc: 215 inline G4double G4ComponentGGHadronNucleusXsc::GetParticleBarCorTot( 216 const G4ParticleDef 216 const G4ParticleDefinition* theParticle, G4int ZZ) 217 { 217 { 218 G4double cor = 1.0; 218 G4double cor = 1.0; 219 G4int z = std::min(92, std::max(ZZ, 1)); 219 G4int z = std::min(92, std::max(ZZ, 1)); 220 if( theParticle == theProton ) cor = fP 220 if( theParticle == theProton ) cor = fProtonBarCorrectionTot[z]; 221 else if( theParticle == theNeutron) cor = fN 221 else if( theParticle == theNeutron) cor = fNeutronBarCorrectionTot[z]; 222 else if( theParticle == thePiPlus ) cor = fP 222 else if( theParticle == thePiPlus ) cor = fPionPlusBarCorrectionTot[z]; 223 else if( theParticle == thePiMinus) cor = fP 223 else if( theParticle == thePiMinus) cor = fPionMinusBarCorrectionTot[z]; 224 return cor; 224 return cor; 225 } 225 } 226 226 227 ////////////////////////////////////////////// 227 ///////////////////////////////////////////////////////////////////// 228 // 228 // 229 // return correction at Tkin = 90*GeV GG -> Ba 229 // return correction at Tkin = 90*GeV GG -> Barashenkov in xsc, when it 230 // is available, else return 1.0 230 // is available, else return 1.0 231 231 232 232 233 inline G4double G4ComponentGGHadronNucleusXsc: 233 inline G4double G4ComponentGGHadronNucleusXsc::GetParticleBarCorIn( 234 const G4ParticleDef 234 const G4ParticleDefinition* theParticle, G4int ZZ) 235 { 235 { 236 G4double cor = 1.0; 236 G4double cor = 1.0; 237 G4int z = std::min(92, std::max(ZZ, 1)); 237 G4int z = std::min(92, std::max(ZZ, 1)); 238 if( theParticle == theProton ) cor = fP 238 if( theParticle == theProton ) cor = fProtonBarCorrectionIn[z]; 239 else if( theParticle == theNeutron) cor = fN 239 else if( theParticle == theNeutron) cor = fNeutronBarCorrectionIn[z]; 240 else if( theParticle == thePiPlus ) cor = fP 240 else if( theParticle == thePiPlus ) cor = fPionPlusBarCorrectionIn[z]; 241 else if( theParticle == thePiMinus) cor = fP 241 else if( theParticle == thePiMinus) cor = fPionMinusBarCorrectionIn[z]; 242 return cor; 242 return cor; 243 } 243 } 244 244 245 #endif 245 #endif 246 246