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