Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ComponentGGNuclNuclXsc.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 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
34// G4GlauberGribovCrossSection
35//
36// 04.09.18 V. Ivantchenko Major revision of interfaces and implementation
37// 27.05.19 V. Ivantchenko Removed obsolete methods and members
38//
39
40#ifndef G4ComponentGGNuclNuclXsc_h
41#define G4ComponentGGNuclNuclXsc_h
42
43#include "globals.hh"
45#include "G4DynamicParticle.hh"
46
50class G4Material;
51
53{
54public:
55
58
59 // virtual interface methods
61 G4double kinEnergy,
62 G4int Z, G4double A) final;
63
65 G4double kinEnergy,
66 G4int Z, G4int A) final;
67
69 G4double kinEnergy,
70 G4int Z, G4double A) final;
71
73 G4double kinEnergy,
74 G4int Z, G4int A) final;
75
77 G4double kinEnergy,
78 G4int Z, G4double A) final;
79
81 G4double kinEnergy,
82 G4int Z, G4int A) final;
83
85 G4double kinEnergy,
86 G4int Z, G4int A) final;
87
88 void BuildPhysicsTable(const G4ParticleDefinition&) final;
89
90 void DumpPhysicsTable(const G4ParticleDefinition&) final;
91
92 void Description(std::ostream&) const final;
93
94 // Extra methods
95 // inline G4double GetElementCrossSection(const G4DynamicParticle*,
96 // G4int Z, const G4Material*);
97
99 G4int Z, G4int A);
100
103 G4double pR, G4double tR);
104
106 G4double kinEnergy, G4int Z, G4int A,
107 G4double pR, G4double tR);
108
111
114
115 inline G4double GetTotalGlauberGribovXsc() const { return fTotalXsc; };
116 inline G4double GetElasticGlauberGribovXsc() const { return fElasticXsc; };
117 inline G4double GetInelasticGlauberGribovXsc() const { return fInelasticXsc; };
118 inline G4double GetProductionGlauberGribovXsc() const { return fProductionXsc; };
119 inline G4double GetDiffractionGlauberGribovXsc() const { return fDiffractionXsc; };
120
121private:
122
123 // Glauber-Gribov cross section
124 void ComputeCrossSections(const G4ParticleDefinition* aParticle,
125 G4double kinEnergy, G4int Z, G4int A);
126
127 G4double fTotalXsc, fElasticXsc, fInelasticXsc;
128 G4double fProductionXsc, fDiffractionXsc;
129 // Cache
130 G4double fEnergy;
131
132 const G4ParticleDefinition* theProton;
133 const G4ParticleDefinition* theNeutron;
134 const G4ParticleDefinition* theLambda;
135
137 G4HadronNucleonXsc* fHNXsc;
138
139 // Cache
140 const G4ParticleDefinition* fParticle;
141 G4int fZ, fA;
142};
143
144inline G4double
146 G4int Z, G4int A)
147{
148 ComputeCrossSections(dp->GetDefinition(), dp->GetKineticEnergy(), Z, A);
149 return fElasticXsc;
150}
151
152inline G4double
154 G4int Z, G4int A)
155{
156 ComputeCrossSections(dp->GetDefinition(), dp->GetKineticEnergy(), Z, A);
157 return fInelasticXsc;
158}
159
160/*
161inline G4double
162G4ComponentGGNuclNuclXsc::GetElementCrossSection(const G4DynamicParticle* dp,
163 G4int Z, const G4Material*)
164{
165 G4int A = G4lrint(fNist->GetAtomicMassAmu(Z));
166 ComputeCrossSections(dp->GetDefinition(), dp->GetKineticEnergy(), Z, A);
167 return fInelasticXsc;
168}
169*/
170inline G4double
172 G4int Z, G4int A)
173{
174 ComputeCrossSections(dp->GetDefinition(), dp->GetKineticEnergy(), Z, A);
175 return fInelasticXsc;
176}
177
178inline G4double
181 G4double pR, G4double tR)
182{
184 G4lrint(Z), G4lrint(A), pR, tR);
185}
186
187#endif
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
G4double GetTotalElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A) final
G4double GetDiffractionGlauberGribovXsc() const
G4double ComputeQuasiElasticRatio(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A) final
G4double GetRatioQE(const G4DynamicParticle *, G4double At, G4double Zt)
G4double GetElasticElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A) final
G4double GetRatioSD(const G4DynamicParticle *, G4double At, G4double Zt)
void DumpPhysicsTable(const G4ParticleDefinition &) final
G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A) final
G4double ComputeCoulombBarier(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A, G4double pR, G4double tR)
G4double GetTotalIsotopeCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A) final
void Description(std::ostream &) const final
G4double GetZandACrossSection(const G4DynamicParticle *, G4int Z, G4int A)
G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A) final
G4double GetElasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
G4double GetCoulombBarier(const G4DynamicParticle *, G4double Z, G4double A, G4double pR, G4double tR)
G4double GetProductionGlauberGribovXsc() const
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
G4double GetInelasticGlauberGribovXsc() const
G4double GetInelasticElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A) final
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
int G4lrint(double ad)
Definition: templates.hh:134