Geant4  10.00.p01
G4GlauberGribovCrossSection.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 // 17.07.06 V. Grichine - first implementation
36 // 22.01.07 V.Ivanchenko - add interface with Z and A
37 // 05.03.07 V.Ivanchenko - add IfZAApplicable
38 // 06.03.07 V.Ivanchenko - add GetElasticGlauberGribov and GetElasticGlauberGribov
39 // for combined dataset
40 //
41 //
42 
43 #ifndef G4GlauberGribovCrossSection_h
44 #define G4GlauberGribovCrossSection_h 1
45 
46 #include "globals.hh"
47 #include "G4Proton.hh"
48 #include "G4Nucleus.hh"
49 
51 
53 class G4HadronNucleonXsc;
54 
56 {
57 public:
58 
61 
62  static const char* Default_Name() {return "Glauber-Gribov";}
63 
64  virtual
66  const G4Element* elm = 0,
67  const G4Material* mat = 0);
68 
69  virtual
71  const G4Isotope* iso = 0,
72  const G4Element* elm = 0,
73  const G4Material* mat = 0);
74 
77 
80 
86 
90 
91  G4double CalculateEcmValue ( const G4double , const G4double , const G4double );
92 
93  G4double CalcMandelstamS( const G4double , const G4double , const G4double );
94 
97 
98  virtual void CrossSectionDescription(std::ostream&) const;
99 
102 
108  inline G4double GetRadiusConst() { return fRadiusConst; };
109 
110  inline G4double GetParticleBarCorTot(const G4ParticleDefinition* theParticle, G4int Z);
111  inline G4double GetParticleBarCorIn(const G4ParticleDefinition* theParticle, G4int Z);
112 
114 
115 private:
116 
117 // const G4double fUpperLimit;
120 
123 
126 
129 
132 
134 // G4double fHadronNucleonXsc;
135 
166 
168 
169 };
170 
172 //
173 // Inlines
174 
175 inline
176 G4double
178  G4int Z, G4int A)
179 {
180  GetIsoCrossSection(dp, Z, A);
181  return fElasticXsc;
182 }
183 
185 
186 inline
187 G4double
189  G4int Z, G4int A)
190 {
191  GetIsoCrossSection(dp, Z, A);
192  return fInelasticXsc;
193 }
194 
196 //
197 // return correction at Tkin = 90*GeV GG -> Barashenkov tot xsc, when it
198 // is available, else return 1.0
199 
200 
202  const G4ParticleDefinition* theParticle, G4int Z)
203 {
204  if(Z >= 2 && Z <= 92)
205  {
206  if( theParticle == theProton ) return fProtonBarCorrectionTot[Z];
207  else if( theParticle == theNeutron) return fNeutronBarCorrectionTot[Z];
208  else if( theParticle == thePiPlus ) return fPionPlusBarCorrectionTot[Z];
209  else if( theParticle == thePiMinus) return fPionMinusBarCorrectionTot[Z];
210  else return 1.0;
211  }
212  else return 1.0;
213 }
214 
216 //
217 // return correction at Tkin = 90*GeV GG -> Barashenkov in xsc, when it
218 // is available, else return 1.0
219 
220 
222  const G4ParticleDefinition* theParticle, G4int Z)
223 {
224  if(Z >= 2 && Z <= 92)
225  {
226  if( theParticle == theProton ) return fProtonBarCorrectionIn[Z];
227  else if( theParticle == theNeutron) return fNeutronBarCorrectionIn[Z];
228  else if( theParticle == thePiPlus ) return fPionPlusBarCorrectionIn[Z];
229  else if( theParticle == thePiMinus) return fPionMinusBarCorrectionIn[Z];
230  else return 1.0;
231  }
232  else return 1.0;
233 }
234 
235 #endif
G4double GetHNinelasticXscVU(const G4DynamicParticle *, G4int At, G4int Zt)
G4double GetRatioSD(const G4DynamicParticle *, G4int At, G4int Zt)
G4double GetParticleBarCorTot(const G4ParticleDefinition *theParticle, G4int Z)
static const G4double fNeutronBarCorrectionTot[93]
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4Element *)
static const G4double fPionMinusBarCorrectionIn[93]
int G4int
Definition: G4Types.hh:78
G4double GetHadronNucleonXsc(const G4DynamicParticle *, const G4Element *)
static const G4double fPionMinusBarCorrectionTot[93]
static const G4double fPionPlusBarCorrectionTot[93]
G4double GetKaonNucleonXscVector(const G4DynamicParticle *, G4int At, G4int Zt)
bool G4bool
Definition: G4Types.hh:79
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4Element *)
static const G4double fProtonBarCorrectionTot[93]
G4double GetParticleBarCorIn(const G4ParticleDefinition *theParticle, G4int Z)
static const G4double A[nN]
G4double CalcMandelstamS(const G4double, const G4double, const G4double)
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
G4double GetNucleusRadius(const G4DynamicParticle *, const G4Element *)
static const G4double fNeutronBarCorrectionIn[93]
G4double CalculateEcmValue(const G4double, const G4double, const G4double)
virtual void CrossSectionDescription(std::ostream &) const
G4double GetRatioQE(const G4DynamicParticle *, G4int At, G4int Zt)
static const G4double fPionPlusBarCorrectionIn[93]
double G4double
Definition: G4Types.hh:76
G4double GetHNinelasticXsc(const G4DynamicParticle *, const G4Element *)
virtual G4bool IsIsoApplicable(const G4DynamicParticle *aDP, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
G4double GetElasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
static const G4double fProtonBarCorrectionIn[93]