Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CrossSectionPairGG Class Reference

#include <G4CrossSectionPairGG.hh>

+ Inheritance 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 IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
virtual const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
virtual G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
bool ForAllAtomsAndEnergies () const
 
void SetForAllAtomsAndEnergies (G4bool val)
 
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 49 of file G4CrossSectionPairGG.hh.

Constructor & Destructor Documentation

◆ G4CrossSectionPairGG()

G4CrossSectionPairGG::G4CrossSectionPairGG ( G4VCrossSectionDataSet low,
G4double  Etransit 
)

Definition at line 47 of file G4CrossSectionPairGG.cc.

48 :
49 G4VCrossSectionDataSet("G4CrossSectionPairGG"), theLowX(low), ETransition(
50 Etransit) {
51 NistMan = G4NistManager::Instance();
52 theHighX = new G4ComponentGGHadronNucleusXsc();
53 verboseLevel = 0;
54}
static G4NistManager * Instance()

Referenced by G4CrossSectionPairGG().

◆ ~G4CrossSectionPairGG()

G4CrossSectionPairGG::~G4CrossSectionPairGG ( )
virtual

Definition at line 56 of file G4CrossSectionPairGG.cc.

56 {
57}

Member Function Documentation

◆ BuildPhysicsTable()

void G4CrossSectionPairGG::BuildPhysicsTable ( const G4ParticleDefinition pDef)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 123 of file G4CrossSectionPairGG.cc.

123 {
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}
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)
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

Referenced by GetElementCrossSection().

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 59 of file G4CrossSectionPairGG.cc.

60 {
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}

◆ DumpPhysicsTable()

void G4CrossSectionPairGG::DumpPhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 197 of file G4CrossSectionPairGG.cc.

197 {
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}

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 84 of file G4CrossSectionPairGG.cc.

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}
double G4double
Definition: G4Types.hh:83
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const

Referenced by GetElementCrossSection().

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 72 of file G4CrossSectionPairGG.cc.

73 {
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}

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