#include <G4CrossSectionPairGG.hh>
|
| 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 &) |
|
| 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 G4Isotope * | SelectIsotope (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 G4String & | GetName () const |
|
Definition at line 51 of file G4CrossSectionPairGG.hh.
Definition at line 50 of file G4CrossSectionPairGG.cc.
G4VCrossSectionDataSet(const G4String &nam="")
static G4NistManager * Instance()
G4CrossSectionPairGG::~G4CrossSectionPairGG |
( |
| ) |
|
|
virtual |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 128 of file G4CrossSectionPairGG.cc.
133 G4cout <<
"G4CrossSectionPairGG::BuildPhysicsTable "
138 std::vector<ParticleXScale>::iterator iter;
139 iter = scale_factors.begin();
140 while (iter != scale_factors.end() && (*iter).first != myDef)
149 if (iter == scale_factors.end()) {
150 XS_factors factors(93);
155 G4cout <<
"G4CrossSectionPairGG::BuildPhysicsTable for particle "
158 for (
G4int aZ = 1; aZ < 93; ++aZ) {
170 G4cout <<
"Z=" << aZ <<
", A=" << AA <<
", scale="
187 ParticleXScale forPart(myDef, factors);
188 scale_factors.push_back(forPart);
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
const G4String & GetName() const
const G4String & GetParticleName() const
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
G4GLOB_DLL std::ostream G4cout
G4double GetAtomicMassAmu(const G4String &symb) const
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
const G4String & GetName() const
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
void G4CrossSectionPairGG::CrossSectionDescription |
( |
std::ostream & |
outFile | ) |
const |
|
virtual |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 64 of file G4CrossSectionPairGG.cc.
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";
Reimplemented from G4VCrossSectionDataSet.
Definition at line 202 of file G4CrossSectionPairGG.cc.
203 G4cout << std::setw(24) <<
" " <<
" G4CrossSectionPairGG: "
205 G4cout << std::setw(27) <<
" " <<
"below " << ETransition /
GeV
206 <<
" GeV, Glauber-Gribov above " <<
G4endl;
const G4String & GetName() const
G4GLOB_DLL std::ostream G4cout
static constexpr double GeV
Reimplemented from G4VCrossSectionDataSet.
Definition at line 89 of file G4CrossSectionPairGG.cc.
99 std::vector<ParticleXScale>::iterator iter = scale_factors.begin();
101 while (iter != scale_factors.end() && (*iter).first != pDef)
105 if (iter != scale_factors.end() )
109 * (*iter).second[ZZ];
112 G4cout <<
" scaling .." << ZZ <<
" " << AA <<
" "
113 << (*iter).second[ZZ] <<
" "
G4double GetKineticEnergy() const
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4ParticleDefinition * GetDefinition() const
G4GLOB_DLL std::ostream G4cout
G4double GetAtomicMassAmu(const G4String &symb) const
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
Reimplemented from G4VCrossSectionDataSet.
Definition at line 77 of file G4CrossSectionPairGG.cc.
79 G4bool isApplicable(
false);
81 if (Ekin <= ETransition) {
G4double GetKineticEnergy() const
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
The documentation for this class was generated from the following files: