source: trunk/source/processes/hadronic/cross_sections/include/G4GGNuclNuclCrossSection.hh @ 990

Last change on this file since 990 was 966, checked in by garnier, 15 years ago

fichier ajoutes

File size: 5.9 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 nucleus-nucleus total, inelastic, production,
27// elastic and quasi-elastic  cross-sections
28// based on parametrisations of nucleon-nucleon
29// cross-sections  in
30// the framework of simplified Glauber-Gribov approach
31//
32//
33//
34//
35//
36// 24.11.08 V. Grichine - first implementation based on G4GlauberGribovCrossSection
37//
38//
39
40#ifndef G4GGNuclNuclCrossSection_h
41#define G4GGNuclNuclCrossSection_h
42
43#include "globals.hh"
44#include "G4Proton.hh"
45#include "G4Nucleus.hh"
46
47#include "G4VCrossSectionDataSet.hh"
48
49class G4ParticleDefinition;
50
51class G4GGNuclNuclCrossSection : public G4VCrossSectionDataSet
52{
53public:
54
55  G4GGNuclNuclCrossSection ();
56  virtual ~G4GGNuclNuclCrossSection ();
57   
58  virtual
59  G4bool IsApplicable(const G4DynamicParticle* aDP, const G4Element*);
60
61  virtual
62  G4bool IsZAApplicable(const G4DynamicParticle* aDP, G4double Z, G4double A);
63
64  virtual
65  G4double GetCrossSection(const G4DynamicParticle*, 
66                           const G4Element*, 
67                           G4double aTemperature = 0.0);
68
69  virtual
70  G4double GetIsoZACrossSection(const G4DynamicParticle*, 
71                                G4double Z, G4double A, 
72                                G4double aTemperature = 0.0);
73
74  virtual
75  void BuildPhysicsTable(const G4ParticleDefinition&)
76  {}
77
78  virtual
79  void DumpPhysicsTable(const G4ParticleDefinition&) 
80  {G4cout << "G4NuclNuclCrossSection: uses Glauber-Gribov formula"<<G4endl;}
81
82  G4double GetRatioSD(const G4DynamicParticle*, G4double At, G4double Zt);
83  G4double GetRatioQE(const G4DynamicParticle*, G4double At, G4double Zt);
84
85  G4double GetHadronNucleonXsc(const G4DynamicParticle*, const G4Element*);
86  G4double GetHadronNucleonXsc(const G4DynamicParticle*, G4double At, G4double Zt);
87
88  G4double GetHadronNucleonXscPDG(const G4DynamicParticle*, const G4Element*);
89  G4double GetHadronNucleonXscPDG(const G4DynamicParticle*, G4double At, G4double Zt);
90
91
92
93
94  // G4double GetHadronNucleonXscNS(const G4DynamicParticle*, const G4Element*);
95  // G4double GetHadronNucleonXscNS(const G4DynamicParticle*,G4double At, G4double Zt);
96
97
98
99
100  G4double GetHadronNucleonXscNS(G4ParticleDefinition*,G4double pTkin, G4ParticleDefinition*);
101
102
103
104
105
106
107
108  // G4double GetHNinelasticXsc(const G4DynamicParticle*, const G4Element*);
109  // G4double GetHNinelasticXsc(const G4DynamicParticle*, G4double At, G4double Zt);
110
111  G4double GetHNinelasticXscVU(const G4DynamicParticle*, G4double At, G4double Zt);
112
113  G4double GetHadronNucleonXscMK(G4ParticleDefinition* pParticle, G4double pTkin, 
114                                 G4ParticleDefinition* nucleon  ); 
115
116  G4double CalculateEcmValue ( const G4double , const G4double , const G4double ); 
117
118  G4double CalcMandelstamS( const G4double , const G4double , const G4double );
119
120  G4double GetElasticGlauberGribov(const G4DynamicParticle*,G4double Z, G4double A);
121  G4double GetInelasticGlauberGribov(const G4DynamicParticle*,G4double Z, G4double A);
122
123  G4double GetTotalGlauberGribovXsc()    { return fTotalXsc;     }; 
124  G4double GetElasticGlauberGribovXsc()  { return fElasticXsc;   }; 
125  G4double GetInelasticGlauberGribovXsc(){ return fInelasticXsc; }; 
126  G4double GetProductionGlauberGribovXsc(){ return fProductionXsc; }; 
127  G4double GetDiffractionGlauberGribovXsc(){ return fDiffractionXsc; }; 
128  G4double GetRadiusConst()              { return fRadiusConst;  }; 
129
130  G4double GetNucleusRadius(const G4DynamicParticle*, const G4Element*);
131
132
133  G4double GetNucleusRadius(G4double At);
134  G4double GetNucleusRadiusGG(G4double At);
135  G4double GetNucleusRadiusDE(G4double At);
136
137
138  inline void SetEnergyLowerLimit(G4double E ){fLowerLimit=E;};
139
140private:
141
142  const G4double fUpperLimit;
143  G4double fLowerLimit; 
144  const G4double fRadiusConst;
145 
146  G4double fTotalXsc, fElasticXsc, fInelasticXsc, fProductionXsc, fDiffractionXsc;
147  G4double fHadronNucleonXsc;
148 
149  G4ParticleDefinition* theProton;
150  G4ParticleDefinition* theNeutron;
151
152};
153
154////////////////////////////////////////////////////////////////
155//
156// Inlines
157
158inline
159G4double G4GGNuclNuclCrossSection::GetElasticGlauberGribov(
160         const G4DynamicParticle* dp, G4double Z, G4double A)
161{
162  GetIsoZACrossSection(dp, Z, A);
163  return fElasticXsc;
164}
165
166/////////////////////////////////////////////////////////////////
167
168inline
169G4double G4GGNuclNuclCrossSection::GetInelasticGlauberGribov(
170         const G4DynamicParticle* dp, G4double Z, G4double A)
171{
172  GetIsoZACrossSection(dp, Z, A);
173  return fInelasticXsc;
174}
175
176
177
178#endif
Note: See TracBrowser for help on using the repository browser.