Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4CrossSectionPairGG Class Reference

#include <G4CrossSectionPairGG.hh>

Inheritance diagram for G4CrossSectionPairGG:
Collaboration diagram for G4CrossSectionPairGG:

Public Member Functions

 G4CrossSectionPairGG (G4VCrossSectionDataSet *low, G4double Etransit)
 
virtual ~G4CrossSectionPairGG ()
 
virtual void CrossSectionDescription (std::ostream &) const
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
virtual G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
virtual G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
const G4StringGetName () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VCrossSectionDataSet
void SetName (const G4String &)
 
- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 

Detailed Description

Definition at line 51 of file G4CrossSectionPairGG.hh.

Constructor & Destructor Documentation

G4CrossSectionPairGG::G4CrossSectionPairGG ( G4VCrossSectionDataSet low,
G4double  Etransit 
)

Definition at line 50 of file G4CrossSectionPairGG.cc.

51  :
52  G4VCrossSectionDataSet("G4CrossSectionPairGG"), theLowX(low), ETransition(
53  Etransit) {
54  NistMan = G4NistManager::Instance();
55  theHighX = new G4ComponentGGHadronNucleusXsc();
56  verboseLevel = 0;
57 }
G4VCrossSectionDataSet(const G4String &nam="")
static G4NistManager * Instance()

Here is the call graph for this function:

G4CrossSectionPairGG::~G4CrossSectionPairGG ( )
virtual

Definition at line 59 of file G4CrossSectionPairGG.cc.

59  {
60  delete theHighX;
61  // The cross section registry will delete theLowX
62 }

Member Function Documentation

void G4CrossSectionPairGG::BuildPhysicsTable ( const G4ParticleDefinition pDef)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 128 of file G4CrossSectionPairGG.cc.

128  {
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 }
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
const G4String & GetName() const
int G4int
Definition: G4Types.hh:78
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
int G4lrint(double ad)
Definition: templates.hh:163
G4double GetAtomicMassAmu(const G4String &symb) const
#define G4endl
Definition: G4ios.hh:61
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
const G4String & GetName() const
virtual void BuildPhysicsTable(const G4ParticleDefinition &)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4CrossSectionPairGG::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 64 of file G4CrossSectionPairGG.cc.

65  {
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 }
void G4CrossSectionPairGG::DumpPhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 202 of file G4CrossSectionPairGG.cc.

202  {
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 }
const G4String & GetName() const
G4GLOB_DLL std::ostream G4cout
static constexpr double GeV
Definition: G4SIunits.hh:217
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

G4double G4CrossSectionPairGG::GetElementCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 89 of file G4CrossSectionPairGG.cc.

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 }
G4double GetKineticEnergy() const
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
int G4lrint(double ad)
Definition: templates.hh:163
G4double GetAtomicMassAmu(const G4String &symb) const
#define G4endl
Definition: G4ios.hh:61
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
double G4double
Definition: G4Types.hh:76
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)

Here is the call graph for this function:

G4bool G4CrossSectionPairGG::IsElementApplicable ( const G4DynamicParticle aParticle,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 77 of file G4CrossSectionPairGG.cc.

78  {
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 }
G4double GetKineticEnergy() const
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
bool G4bool
Definition: G4Types.hh:79
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:


The documentation for this class was generated from the following files: