source: trunk/source/processes/hadronic/cross_sections/include/G4GlauberGribovCrossSection.hh @ 1350

Last change on this file since 1350 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 8.7 KB
Line 
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
45
46#include "globals.hh"
47#include "G4Proton.hh"
48#include "G4Nucleus.hh"
49#include "G4HadTmpUtil.hh"
50
51#include "G4VCrossSectionDataSet.hh"
52
53class G4ParticleDefinition;
54
55class G4GlauberGribovCrossSection : public G4VCrossSectionDataSet
56{
57public:
58
59  G4GlauberGribovCrossSection ();
60  virtual ~G4GlauberGribovCrossSection ();
61   
62  virtual
63  G4bool IsApplicable(const G4DynamicParticle* aDP, const G4Element*);
64
65  virtual
66  G4bool IsIsoApplicable(const G4DynamicParticle* aDP, G4int Z, G4int A);
67
68  virtual
69  G4double GetCrossSection(const G4DynamicParticle*, 
70                           const G4Element*, 
71                           G4double aTemperature = 0.0);
72
73  virtual
74  G4double GetZandACrossSection(const G4DynamicParticle*,
75                                G4int Z, G4int A,
76                                G4double aTemperature = 0.0);
77
78  virtual
79  void BuildPhysicsTable(const G4ParticleDefinition&)
80  {}
81
82  virtual
83  void DumpPhysicsTable(const G4ParticleDefinition&) 
84  {G4cout << "G4GlauberGribovCrossSection: uses Glauber-Gribov formula"<<G4endl;}
85
86  G4double GetRatioSD(const G4DynamicParticle*, G4int At, G4int Zt);
87  G4double GetRatioQE(const G4DynamicParticle*, G4int At, G4int Zt);
88
89  G4double GetHadronNucleonXsc(const G4DynamicParticle*, const G4Element*);
90  G4double GetHadronNucleonXsc(const G4DynamicParticle*, G4int At, G4int Zt);
91
92  G4double GetHadronNucleonXscPDG(const G4DynamicParticle*, const G4Element*);
93  G4double GetHadronNucleonXscPDG(const G4DynamicParticle*, G4int At, G4int Zt);
94
95  G4double GetHadronNucleonXscNS(const G4DynamicParticle*, const G4Element*);
96  G4double GetHadronNucleonXscNS(const G4DynamicParticle*, G4int At, G4int Zt);
97
98  G4double GetHNinelasticXsc(const G4DynamicParticle*, const G4Element*);
99  G4double GetHNinelasticXsc(const G4DynamicParticle*, G4int At, G4int Zt);
100  G4double GetHNinelasticXscVU(const G4DynamicParticle*, G4int At, G4int Zt);
101
102  G4double CalculateEcmValue ( const G4double , const G4double , const G4double ); 
103
104  G4double CalcMandelstamS( const G4double , const G4double , const G4double );
105
106  G4double GetElasticGlauberGribov(const G4DynamicParticle*, G4int Z, G4int A);
107  G4double GetInelasticGlauberGribov(const G4DynamicParticle*, G4int Z, G4int A);
108
109  G4double GetTotalGlauberGribovXsc()    { return fTotalXsc;     }; 
110  G4double GetElasticGlauberGribovXsc()  { return fElasticXsc;   }; 
111  G4double GetInelasticGlauberGribovXsc(){ return fInelasticXsc; }; 
112  G4double GetProductionGlauberGribovXsc(){ return fProductionXsc; }; 
113  G4double GetDiffractionGlauberGribovXsc(){ return fDiffractionXsc; }; 
114  G4double GetRadiusConst()              { return fRadiusConst;  }; 
115
116  G4double GetNucleusRadius(const G4DynamicParticle*, const G4Element*);
117  G4double GetNucleusRadius(G4int At);
118
119  inline G4double GetParticleBarCorTot(const G4ParticleDefinition* theParticle, G4int Z);
120  inline G4double GetParticleBarCorIn(const G4ParticleDefinition* theParticle, G4int Z);
121
122  inline void SetEnergyLowerLimit(G4double E ){fLowerLimit=E;};
123
124private:
125
126  const G4double fUpperLimit;
127  G4double fLowerLimit; 
128  const G4double fRadiusConst;
129
130  static const G4double fNeutronBarCorrectionTot[93];
131  static const G4double fNeutronBarCorrectionIn[93];
132
133  static const G4double fProtonBarCorrectionTot[93];
134  static const G4double fProtonBarCorrectionIn[93];
135
136  static const G4double fPionPlusBarCorrectionTot[93];
137  static const G4double fPionPlusBarCorrectionIn[93];
138
139  static const G4double fPionMinusBarCorrectionTot[93];
140  static const G4double fPionMinusBarCorrectionIn[93];
141
142  G4double fTotalXsc, fElasticXsc, fInelasticXsc, fProductionXsc, fDiffractionXsc;
143  G4double fHadronNucleonXsc;
144 
145  G4ParticleDefinition* theGamma;
146  G4ParticleDefinition* theProton;
147  G4ParticleDefinition* theNeutron;
148  G4ParticleDefinition* theAProton;
149  G4ParticleDefinition* theANeutron;
150  G4ParticleDefinition* thePiPlus;
151  G4ParticleDefinition* thePiMinus;
152  G4ParticleDefinition* thePiZero;
153  G4ParticleDefinition* theKPlus;
154  G4ParticleDefinition* theKMinus;
155  G4ParticleDefinition* theK0S;
156  G4ParticleDefinition* theK0L;
157  G4ParticleDefinition* theL;
158  G4ParticleDefinition* theAntiL;
159  G4ParticleDefinition* theSPlus;
160  G4ParticleDefinition* theASPlus;
161  G4ParticleDefinition* theSMinus;
162  G4ParticleDefinition* theASMinus;
163  G4ParticleDefinition* theS0;
164  G4ParticleDefinition* theAS0;
165  G4ParticleDefinition* theXiMinus;
166  G4ParticleDefinition* theXi0;
167  G4ParticleDefinition* theAXiMinus;
168  G4ParticleDefinition* theAXi0;
169  G4ParticleDefinition* theOmega;
170  G4ParticleDefinition* theAOmega;
171  G4ParticleDefinition* theD;
172  G4ParticleDefinition* theT;
173  G4ParticleDefinition* theA;
174  G4ParticleDefinition* theHe3;
175
176};
177
178////////////////////////////////////////////////////////////////
179//
180// Inlines
181
182inline
183G4double
184G4GlauberGribovCrossSection::GetElasticGlauberGribov(const G4DynamicParticle* dp,
185                                                     G4int Z, G4int A)
186{
187  GetZandACrossSection(dp, Z, A);
188  return fElasticXsc;
189}
190
191/////////////////////////////////////////////////////////////////
192
193inline
194G4double
195G4GlauberGribovCrossSection::GetInelasticGlauberGribov(const G4DynamicParticle* dp,
196                                                       G4int Z, G4int A)
197{
198  GetZandACrossSection(dp, Z, A);
199  return fInelasticXsc;
200}
201
202/////////////////////////////////////////////////////////////////////
203//
204// return correction at Tkin = 90*GeV GG -> Barashenkov tot xsc, when it
205// is available, else return 1.0
206
207
208inline G4double G4GlauberGribovCrossSection::GetParticleBarCorTot( 
209                          const G4ParticleDefinition* theParticle, G4int Z)
210{
211  if(Z >= 2 && Z <= 92)
212  {
213    if(      theParticle == theProton ) return fProtonBarCorrectionTot[Z]; 
214    else if( theParticle == theNeutron) return fNeutronBarCorrectionTot[Z]; 
215    else if( theParticle == thePiPlus ) return fPionPlusBarCorrectionTot[Z];
216    else if( theParticle == thePiMinus) return fPionMinusBarCorrectionTot[Z];
217    else return 1.0;
218  }
219  else return 1.0;
220}
221
222/////////////////////////////////////////////////////////////////////
223//
224// return correction at Tkin = 90*GeV GG -> Barashenkov in xsc, when it
225// is available, else return 1.0
226
227
228inline G4double G4GlauberGribovCrossSection::GetParticleBarCorIn( 
229                          const G4ParticleDefinition* theParticle, G4int Z)
230{
231  if(Z >= 2 && Z <= 92)
232  {
233    if(      theParticle == theProton ) return fProtonBarCorrectionIn[Z]; 
234    else if( theParticle == theNeutron) return fNeutronBarCorrectionIn[Z]; 
235    else if( theParticle == thePiPlus ) return fPionPlusBarCorrectionIn[Z];
236    else if( theParticle == thePiMinus) return fPionMinusBarCorrectionIn[Z];
237    else return 1.0;
238  }
239  else return 1.0;
240}
241
242#endif
Note: See TracBrowser for help on using the repository browser.