Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleInelasticXS.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// -------------------------------------------------------------------
27//
28// GEANT4 Class header file
29//
30//
31// File name: G4ParticleInelasticXS
32//
33// Author Ivantchenko, Geant4, 24 May 2018
34//
35
36// Class Description:
37// This is a base class for n,p,d,t,he3,he4 inelastic hadronic cross
38// section based on data files from G4PARTICLEXSDATA data set
39// and Glauber-Gribov model for high energy
40//
41
42#ifndef G4ParticleInelasticXS_h
43#define G4ParticleInelasticXS_h 1
44
46#include "globals.hh"
47#include "G4ElementData.hh"
48#include "G4PhysicsVector.hh"
49#include <vector>
50
53class G4Element;
55
57{
58public:
59
61
63
65 const G4Material* mat = nullptr) final;
66
68 const G4Element* elm = nullptr,
69 const G4Material* mat = nullptr) final;
70
72 const G4Material* mat = nullptr) final;
73
76 const G4Element*,
77 const G4Material*) final;
78
81 G4int Z, G4int A,
82 const G4Isotope* iso,
83 const G4Element* elm,
84 const G4Material* mat) final;
85
87 const G4Isotope* iso = nullptr,
88 const G4Element* elm = nullptr,
89 const G4Material* mat = nullptr) final;
90
91 const G4Isotope* SelectIsotope(const G4Element*,
92 G4double kinEnergy, G4double logE) final;
93
94 void BuildPhysicsTable(const G4ParticleDefinition&) final;
95
96 void CrossSectionDescription(std::ostream&) const final;
97
99
101
103 (const G4ParticleInelasticXS &right) = delete;
105
106private:
107
108 void Initialise(G4int Z);
109
110 void InitialiseOnFly(G4int Z);
111
112 const G4String& FindDirectoryPath();
113
114 inline const G4PhysicsVector* GetPhysicsVector(G4int Z);
115
116 G4PhysicsVector* RetrieveVector(std::ostringstream& in, G4bool warn);
117
118 G4VComponentCrossSection* highEnergyXsection = nullptr;
119 const G4ParticleDefinition* particle;
120
121 std::vector<G4double> temp;
122 G4double elimit;
123
124 G4int index = 0;
125 G4bool isMaster = false;
126
127 static const G4int MAXZINELP = 93;
128 static G4ElementData* data[5];
129 static G4double coeff[MAXZINELP][5];
130 static G4String gDataDirectory[5];
131};
132
133inline
134const G4PhysicsVector* G4ParticleInelasticXS::GetPhysicsVector(G4int Z)
135{
136 const G4PhysicsVector* pv = data[index]->GetElementData(Z);
137 if(pv == nullptr) {
138 InitialiseOnFly(Z);
139 pv = data[index]->GetElementData(Z);
140 }
141 return pv;
142}
143
144#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
G4PhysicsVector * GetElementData(G4int Z)
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr) final
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
G4double ComputeCrossSectionPerElement(G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *) final
G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr) final
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
G4ParticleInelasticXS(const G4ParticleInelasticXS &)=delete
G4double IsoCrossSection(G4double ekin, G4double logE, G4int Z, G4int A)
G4double ComputeIsoCrossSection(G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) final
const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy, G4double logE) final
G4double ElementCrossSection(G4double kinEnergy, G4double loge, G4int Z)
void CrossSectionDescription(std::ostream &) const final