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