Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CrossSectionPairGG.cc
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// Class G4CrossSectionPairGG
28//
29// smoothly join two cross section sets by scaling the second at a given
30// transition energy to match the first.
31//
32// Author: Gunter Folger
33// November 2009
34//
35
37
38#include "globals.hh"
40#include "G4SystemOfUnits.hh"
41#include "G4HadTmpUtil.hh"
42#include "G4NistManager.hh"
43#include "G4ThreeVector.hh"
44#include "G4NistManager.hh"
46
48 G4double Etransit) :
49 G4VCrossSectionDataSet("G4CrossSectionPairGG"), theLowX(low), ETransition(
50 Etransit) {
51 NistMan = G4NistManager::Instance();
52 theHighX = new G4ComponentGGHadronNucleusXsc();
53 verboseLevel = 0;
54}
55
57}
58
60 std::ostream& outFile) const {
61 outFile << "G4CrossSectionPairGG is used to add the relativistic rise to\n"
62 << "hadronic cross section data sets above a given energy. In this\n"
63 << "case, the Glauber-Gribov cross section is used above 91 GeV.\n"
64 << "At this energy the low energy cross section is smoothly joined\n"
65 << "to the high energy cross section. Below 91 GeV the Barashenkov\n"
66 << "cross section is used for pions (G4PiNuclearCrossSection), the\n"
67 << "Axen-Wellisch cross section is used for protons\n"
68 << "(G4ProtonInelasticCrossSection), and the Wellisch-Laidlaw cross\n"
69 << "section is used for neutrons (G4NeutronInelasticCrossSection).\n";
70}
71
73 const G4DynamicParticle* aParticle, G4int Z, const G4Material* mat) {
74 G4bool isApplicable(false);
75 G4double Ekin = aParticle->GetKineticEnergy();
76 if (Ekin <= ETransition) {
77 isApplicable = theLowX->IsElementApplicable(aParticle, Z, mat);
78 } else if (Z > 1) {
79 isApplicable = true;
80 }
81 return isApplicable;
82}
83
85 const G4DynamicParticle* aParticle, G4int ZZ, const G4Material* mat)
86{
87 G4double Xsec(0.);
88
89 if (aParticle->GetKineticEnergy() < ETransition)
90 {
91 Xsec = theLowX->GetElementCrossSection(aParticle, ZZ, mat);
92 } else {
93
94 std::vector<ParticleXScale>::iterator iter = scale_factors.begin();
95 const G4ParticleDefinition * pDef = aParticle->GetDefinition();
96 while (iter != scale_factors.end() && (*iter).first != pDef) /* Loop checking, 08.01.2016, W. Pokorski */
97 {
98 ++iter;
99 }
100 if (iter != scale_factors.end() )
101 {
102 G4int AA = G4lrint(NistMan->GetAtomicMassAmu(ZZ));
103 Xsec = theHighX->GetInelasticGlauberGribov(aParticle, ZZ, AA)
104 * (*iter).second[ZZ];
105 if (verboseLevel > 2)
106 {
107 G4cout << " scaling .." << ZZ << " " << AA << " "
108 << (*iter).second[ZZ] << " "
109 << theHighX->GetInelasticGlauberGribov(aParticle, ZZ, AA)
110 << " " << Xsec << G4endl;
111 }
112 } else {
113 // BuildPhysicsTable not done for pDef=aParticle->GetDefinition
114 // build table, and recurse
115 BuildPhysicsTable(*pDef);
116 Xsec=GetElementCrossSection(aParticle, ZZ, mat);
117 }
118 }
119
120 return Xsec;
121}
122
124 theLowX->BuildPhysicsTable(pDef);
125 theHighX->BuildPhysicsTable(pDef);
126
127 if (verboseLevel > 0) {
128 G4cout << "G4CrossSectionPairGG::BuildPhysicsTable "
129 << theLowX->GetName() << " " << theHighX->GetName() << G4endl;
130 }
131
132 const G4ParticleDefinition * myDef = &pDef;
133 std::vector<ParticleXScale>::iterator iter;
134 iter = scale_factors.begin();
135 while (iter != scale_factors.end() && (*iter).first != myDef) /* Loop checking, 08.01.2016, W. Pokorski */
136 {
137 ++iter;
138 }
139
140 // new particle, initialise
141
142 G4Material* mat = 0;
143
144 if (iter == scale_factors.end()) {
145 XS_factors factors(93);
146 G4ThreeVector mom(0.0, 0.0, 1.0);
147 G4DynamicParticle DynPart(myDef, mom, ETransition); // last is kinetic Energy
148
149 if (verboseLevel > 0) {
150 G4cout << "G4CrossSectionPairGG::BuildPhysicsTable for particle "
151 << pDef.GetParticleName() << G4endl;
152 }
153 for (G4int aZ = 1; aZ < 93; ++aZ) {
154 factors[aZ] = 1.; // default, to give reasonable value if only high is applicable
155 G4int AA = G4lrint(NistMan->GetAtomicMassAmu(aZ));
156 G4bool isApplicable = theLowX->IsElementApplicable(&DynPart, aZ,
157 mat) && (aZ > 1);
158
159 if (isApplicable) {
160 factors[aZ] = theLowX->GetElementCrossSection(&DynPart, aZ, mat)
161 / theHighX->GetInelasticGlauberGribov(&DynPart, aZ, AA);
162
163 }
164 if (verboseLevel > 0) {
165 G4cout << "Z=" << aZ << ", A=" << AA << ", scale="
166 << factors[aZ];
167 if (verboseLevel == 1) {
168 G4cout << G4endl;
169 } else {
170 if (isApplicable) {
171 G4cout << ", low / high "
172 << theLowX->GetElementCrossSection(&DynPart, aZ,
173 mat) << " "
174 << theHighX->GetInelasticGlauberGribov(&DynPart,
175 aZ, AA) << G4endl;
176 } else {
177 G4cout << ", N/A" << G4endl;
178 }
179 }
180 }
181 }
182 ParticleXScale forPart(myDef, factors);
183 scale_factors.push_back(forPart);
184 }
185}
186
187/*
188 void G4CrossSectionPairGG::DumpHtml(const G4ParticleDefinition&,
189 std::ofstream outFile)
190 {
191 outFile << " <li><b>"
192 << " G4CrossSectionPairGG: " << theLowX->GetName() << " cross sections \n";
193 outFile << "below " << ETransition/GeV << " GeV, Glauber-Gribov above \n";
194 }
195 */
196
198 G4cout << std::setw(24) << " " << " G4CrossSectionPairGG: "
199 << theLowX->GetName() << " cross sections " << G4endl;
200 G4cout << std::setw(27) << " " << "below " << ETransition / GeV
201 << " GeV, Glauber-Gribov above " << G4endl;
202}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
virtual void DumpPhysicsTable(const G4ParticleDefinition &)
G4CrossSectionPairGG(G4VCrossSectionDataSet *low, G4double Etransit)
virtual void CrossSectionDescription(std::ostream &) const
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
const G4String & GetParticleName() const
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
const G4String & GetName() const
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
const G4String & GetName() const
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
int G4lrint(double ad)
Definition: templates.hh:134