Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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//
32//
33//
34//
35// 25.04.12 V. Grichine - first implementation based on G4GlauberGribovCrossSection old interface
36//
37//
38
39#ifndef G4ComponentGGHadronNucleusXsc_h
40#define G4ComponentGGHadronNucleusXsc_h 1
41
42#include "globals.hh"
43#include "G4Proton.hh"
44#include "G4Nucleus.hh"
45
47
50
52{
53public:
54
57
58
59 // virtual interface methods
60
61virtual
63 G4double kinEnergy,
64 G4int Z, G4int A);
65
66virtual
68 G4double kinEnergy,
69 G4int Z, G4double A);
70
71virtual
73 G4double kinEnergy,
74 G4int Z, G4int A);
75
76virtual
78 G4double kinEnergy,
79 G4int Z, G4double A);
80
81virtual
83 G4double kinEnergy,
84 G4int Z, G4double A);
85
86virtual
88 G4double kinEnergy,
89 G4int Z, G4int A);
90
91virtual
93 G4double kinEnergy,
94 G4int Z, G4int A);
95
96
97
98 // virtual
100 const G4Element* elm = 0,
101 const G4Material* mat = 0);
102
103 // virtual
105 const G4Isotope* iso = 0,
106 const G4Element* elm = 0,
107 const G4Material* mat = 0);
108
111
114
120
124
125 G4double CalculateEcmValue ( const G4double , const G4double , const G4double );
126
127 G4double CalcMandelstamS( const G4double , const G4double , const G4double );
128
131
132 virtual void CrossSectionDescription(std::ostream&) const;
133
136
137 inline G4double GetTotalGlauberGribovXsc() { return fTotalXsc; };
138 inline G4double GetElasticGlauberGribovXsc() { return fElasticXsc; };
139 inline G4double GetInelasticGlauberGribovXsc(){ return fInelasticXsc; };
140 inline G4double GetProductionGlauberGribovXsc(){ return fProductionXsc; };
141 inline G4double GetDiffractionGlauberGribovXsc(){ return fDiffractionXsc; };
142 inline G4double GetRadiusConst() { return fRadiusConst; };
143
144 inline G4double GetParticleBarCorTot(const G4ParticleDefinition* theParticle, G4int Z);
145 inline G4double GetParticleBarCorIn(const G4ParticleDefinition* theParticle, G4int Z);
146
147 inline void SetEnergyLowerLimit(G4double E ){fLowerLimit=E;};
148
149private:
150
151 const G4double fUpperLimit;
152 G4double fLowerLimit;
153 const G4double fRadiusConst;
154
155 static const G4double fNeutronBarCorrectionTot[93];
156 static const G4double fNeutronBarCorrectionIn[93];
157
158 static const G4double fProtonBarCorrectionTot[93];
159 static const G4double fProtonBarCorrectionIn[93];
160
161 static const G4double fPionPlusBarCorrectionTot[93];
162 static const G4double fPionPlusBarCorrectionIn[93];
163
164 static const G4double fPionMinusBarCorrectionTot[93];
165 static const G4double fPionMinusBarCorrectionIn[93];
166
167 G4double fTotalXsc, fElasticXsc, fInelasticXsc, fProductionXsc, fDiffractionXsc;
168 G4double fHadronNucleonXsc;
169
170 G4ParticleDefinition* theGamma;
171 G4ParticleDefinition* theProton;
172 G4ParticleDefinition* theNeutron;
173 G4ParticleDefinition* theAProton;
174 G4ParticleDefinition* theANeutron;
175 G4ParticleDefinition* thePiPlus;
176 G4ParticleDefinition* thePiMinus;
177 G4ParticleDefinition* thePiZero;
178 G4ParticleDefinition* theKPlus;
179 G4ParticleDefinition* theKMinus;
180 G4ParticleDefinition* theK0S;
181 G4ParticleDefinition* theK0L;
183 G4ParticleDefinition* theAntiL;
184 G4ParticleDefinition* theSPlus;
185 G4ParticleDefinition* theASPlus;
186 G4ParticleDefinition* theSMinus;
187 G4ParticleDefinition* theASMinus;
189 G4ParticleDefinition* theAS0;
190 G4ParticleDefinition* theXiMinus;
191 G4ParticleDefinition* theXi0;
192 G4ParticleDefinition* theAXiMinus;
193 G4ParticleDefinition* theAXi0;
194 G4ParticleDefinition* theOmega;
195 G4ParticleDefinition* theAOmega;
199 G4ParticleDefinition* theHe3;
200
201 G4HadronNucleonXsc* hnXsc;
202
203};
204
205////////////////////////////////////////////////////////////////
206//
207// Inlines
208
209inline
212 G4int Z, G4int A)
213{
214 GetIsoCrossSection(dp, Z, A);
215 return fElasticXsc;
216}
217
218/////////////////////////////////////////////////////////////////
219
220inline
223 G4int Z, G4int A)
224{
225 GetIsoCrossSection(dp, Z, A);
226 return fInelasticXsc;
227}
228
229/////////////////////////////////////////////////////////////////////
230//
231// return correction at Tkin = 90*GeV GG -> Barashenkov tot xsc, when it
232// is available, else return 1.0
233
234
236 const G4ParticleDefinition* theParticle, G4int Z)
237{
238 if(Z >= 2 && Z <= 92)
239 {
240 if( theParticle == theProton ) return fProtonBarCorrectionTot[Z];
241 else if( theParticle == theNeutron) return fNeutronBarCorrectionTot[Z];
242 else if( theParticle == thePiPlus ) return fPionPlusBarCorrectionTot[Z];
243 else if( theParticle == thePiMinus) return fPionMinusBarCorrectionTot[Z];
244 else return 1.0;
245 }
246 else return 1.0;
247}
248
249/////////////////////////////////////////////////////////////////////
250//
251// return correction at Tkin = 90*GeV GG -> Barashenkov in xsc, when it
252// is available, else return 1.0
253
254
256 const G4ParticleDefinition* theParticle, G4int Z)
257{
258 if(Z >= 2 && Z <= 92)
259 {
260 if( theParticle == theProton ) return fProtonBarCorrectionIn[Z];
261 else if( theParticle == theNeutron) return fNeutronBarCorrectionIn[Z];
262 else if( theParticle == thePiPlus ) return fPionPlusBarCorrectionIn[Z];
263 else if( theParticle == thePiMinus) return fPionMinusBarCorrectionIn[Z];
264 else return 1.0;
265 }
266 else return 1.0;
267}
268
269#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4double GetNucleusRadius(const G4DynamicParticle *, const G4Element *)
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4Element *)
virtual G4double GetInelasticElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4Element *)
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)
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4double CalcMandelstamS(const G4double, const G4double, const G4double)
virtual G4double GetTotalIsotopeCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
virtual G4double ComputeQuasiElasticRatio(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
virtual G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
G4double GetRatioQE(const G4DynamicParticle *, G4int At, G4int Zt)
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
virtual G4double GetTotalElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
G4double GetKaonNucleonXscVector(const G4DynamicParticle *, G4int At, G4int Zt)
virtual G4double GetElasticElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
G4double GetRatioSD(const G4DynamicParticle *, G4int At, G4int Zt)
virtual void CrossSectionDescription(std::ostream &) const
G4double GetParticleBarCorIn(const G4ParticleDefinition *theParticle, G4int Z)
G4double GetHNinelasticXsc(const G4DynamicParticle *, const G4Element *)
virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
G4double CalculateEcmValue(const G4double, const G4double, const G4double)
G4bool IsIsoApplicable(const G4DynamicParticle *aDP, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)