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

Last change on this file since 1347 was 1340, checked in by garnier, 15 years ago

update ti head

File size: 5.5 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// 24.11.08 V. Grichine - first implementation based on G4GlauberGribovCrossSection
34//
35//
36
37#ifndef G4GGNuclNuclCrossSection_h
38#define G4GGNuclNuclCrossSection_h
39
40#include "globals.hh"
41#include "G4Proton.hh"
42#include "G4Nucleus.hh"
43
44#include "G4VCrossSectionDataSet.hh"
45
46class G4ParticleDefinition;
47
48class G4GGNuclNuclCrossSection : public G4VCrossSectionDataSet
49{
50public:
51
52 G4GGNuclNuclCrossSection ();
53 virtual ~G4GGNuclNuclCrossSection ();
54
55 virtual
56 G4bool IsApplicable(const G4DynamicParticle* aDP, const G4Element*);
57
58 virtual
59 G4bool IsIsoApplicable(const G4DynamicParticle* aDP, G4int Z, G4int A);
60
61 virtual
62 G4double GetCrossSection(const G4DynamicParticle*,
63 const G4Element*,
64 G4double aTemperature = 0.0);
65
66 virtual
67 G4double GetZandACrossSection(const G4DynamicParticle*,
68 G4int Z, G4int A,
69 G4double aTemperature = 0.0);
70
71 G4double GetCoulombBarier(const G4DynamicParticle*,
72 G4double Z, G4double A, G4double pR, G4double tR);
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*, G4int At, G4int Zt);
87
88 G4double GetHadronNucleonXscPDG(const G4DynamicParticle*, const G4Element*);
89 G4double GetHadronNucleonXscPDG(const G4DynamicParticle*, G4int At, G4int Zt);
90 G4double GetHadronNucleonXscNS(G4ParticleDefinition*,G4double pTkin, G4ParticleDefinition*);
91
92 G4double GetHNinelasticXscVU(const G4DynamicParticle*, G4int At, G4int Zt);
93 G4double CalculateEcmValue(const G4double, const G4double, const G4double);
94 G4double CalcMandelstamS( const G4double , const G4double , const G4double );
95
96 G4double GetElasticGlauberGribov(const G4DynamicParticle*,G4int Z, G4int A);
97 G4double GetInelasticGlauberGribov(const G4DynamicParticle*,G4int Z, G4int A);
98
99 G4double GetTotalGlauberGribovXsc() { return fTotalXsc; };
100 G4double GetElasticGlauberGribovXsc() { return fElasticXsc; };
101 G4double GetInelasticGlauberGribovXsc(){ return fInelasticXsc; };
102 G4double GetProductionGlauberGribovXsc(){ return fProductionXsc; };
103 G4double GetDiffractionGlauberGribovXsc(){ return fDiffractionXsc; };
104 G4double GetRadiusConst() { return fRadiusConst; };
105
106 G4double GetNucleusRadius(const G4DynamicParticle*, const G4Element*);
107
108 G4double GetNucleusRadius(G4double At);
109 G4double GetNucleusRadiusGG(G4double At);
110 G4double GetNucleusRadiusDE(G4double At);
111
112 inline void SetEnergyLowerLimit(G4double E ){fLowerLimit=E;};
113
114private:
115
116 const G4double fUpperLimit;
117 G4double fLowerLimit;
118 const G4double fRadiusConst;
119
120 G4double fTotalXsc, fElasticXsc, fInelasticXsc, fProductionXsc, fDiffractionXsc;
121 G4double fHadronNucleonXsc;
122
123 G4ParticleDefinition* theProton;
124 G4ParticleDefinition* theNeutron;
125
126};
127
128////////////////////////////////////////////////////////////////
129//
130// Inlines
131
132inline G4double
133G4GGNuclNuclCrossSection::GetElasticGlauberGribov(const G4DynamicParticle* dp,
134 G4int Z, G4int A)
135{
136 GetZandACrossSection(dp, Z, A);
137 return fElasticXsc;
138}
139
140/////////////////////////////////////////////////////////////////
141
142inline G4double
143G4GGNuclNuclCrossSection::GetInelasticGlauberGribov(const G4DynamicParticle* dp,
144 G4int Z, G4int A)
145{
146 GetZandACrossSection(dp, Z, A);
147 return fInelasticXsc;
148}
149
150#endif
Note: See TracBrowser for help on using the repository browser.