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

#include <G4BGGPionElasticXS.hh>

+ Inheritance diagram for G4BGGPionElasticXS:

Public Member Functions

 G4BGGPionElasticXS (const G4ParticleDefinition *)
 
 ~G4BGGPionElasticXS () final
 
G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *) final
 
G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm, const G4Material *mat) final
 
G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat) final
 
G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
 
void BuildPhysicsTable (const G4ParticleDefinition &) final
 
void CrossSectionDescription (std::ostream &) const final
 
- 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 64 of file G4BGGPionElasticXS.hh.

Constructor & Destructor Documentation

◆ G4BGGPionElasticXS()

G4BGGPionElasticXS::G4BGGPionElasticXS ( const G4ParticleDefinition p)
explicit

Definition at line 68 of file G4BGGPionElasticXS.cc.

69 : G4VCrossSectionDataSet("BarashenkovGlauberGribov")
70{
71 verboseLevel = 0;
72 fGlauberEnergy = 91.*GeV;
73 fLowEnergy = 20.*MeV;
74 fLowestEnergy = 1.*MeV;
75 SetMinKinEnergy(0.0);
77
78 fPion = nullptr;
79 fGlauber = nullptr;
80 fHadron = nullptr;
81
82 fG4pow = G4Pow::GetInstance();
83
84 theProton= G4Proton::Proton();
85 thePiPlus= G4PionPlus::PionPlus();
86 isPiplus = (p == thePiPlus);
87 isMaster = false;
89}
static G4HadronicParameters * Instance()
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:97
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
static G4Proton * Proton()
Definition: G4Proton.cc:92
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
void SetForAllAtomsAndEnergies(G4bool val)

◆ ~G4BGGPionElasticXS()

G4BGGPionElasticXS::~G4BGGPionElasticXS ( )
final

Definition at line 93 of file G4BGGPionElasticXS.cc.

94{
95 delete fHadron;
96}

Member Function Documentation

◆ BuildPhysicsTable()

void G4BGGPionElasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 175 of file G4BGGPionElasticXS.cc.

176{
177 if(fPion) { return; }
178 if(verboseLevel > 1) {
179 G4cout << "G4BGGPionElasticXS::BuildPhysicsTable for "
180 << p.GetParticleName() << G4endl;
181 }
182 if(&p == G4PionPlus::PionPlus() || &p == G4PionMinus::PionMinus()) {
183 isPiplus = (&p == G4PionPlus::PionPlus());
184 } else {
186 ed << "This BGG cross section is applicable only to pions and not to "
187 << p.GetParticleName() << G4endl;
188 G4Exception("G4BGGPionElasticXS::BuildPhysicsTable", "had001",
189 FatalException, ed);
190 return;
191 }
192
193 fPion = new G4UPiNuclearCrossSection();
194 fGlauber = new G4ComponentGGHadronNucleusXsc();
195 fHadron = new G4HadronNucleonXsc();
196
197 fPion->BuildPhysicsTable(p);
198
199 if(0 == theA[0]) {
200#ifdef G4MULTITHREADED
201 G4MUTEXLOCK(&pionElasticXSMutex);
202 if(0 == theA[0]) {
203#endif
204 isMaster = true;
205#ifdef G4MULTITHREADED
206 }
207 G4MUTEXUNLOCK(&pionElasticXSMutex);
208#endif
209 } else {
210 return;
211 }
212
213 if(isMaster && 0 == theA[0]) {
214
215 theA[0] = theA[1] = 1;
216 G4ThreeVector mom(0.0,0.0,1.0);
217 G4DynamicParticle dp(thePiPlus, mom, fGlauberEnergy);
218
220
221 G4double csup, csdn;
222 for(G4int iz=2; iz<93; ++iz) {
223
224 G4int A = G4lrint(nist->GetAtomicMassAmu(iz));
225 theA[iz] = A;
226
227 csup = fGlauber->GetElasticGlauberGribov(&dp, iz, A);
228 csdn = fPion->GetElasticCrossSection(&dp, iz, A);
229 theGlauberFacPiPlus[iz] = csdn/csup;
230 }
231
232 dp.SetDefinition(G4PionMinus::PionMinus());
233 for(G4int iz=2; iz<93; ++iz) {
234 csup = fGlauber->GetElasticGlauberGribov(&dp, iz, theA[iz]);
235 csdn = fPion->GetElasticCrossSection(&dp, iz, theA[iz]);
236 theGlauberFacPiMinus[iz] = csdn/csup;
237
238 if(verboseLevel > 0) {
239 G4cout << "Z= " << iz << " A= " << theA[iz]
240 << " factorPiPlus= " << theGlauberFacPiPlus[iz]
241 << " factorPiMinus= " << theGlauberFacPiMinus[iz]
242 << G4endl;
243 }
244 }
245 theCoulombFacPiPlus[1] = 1.0;
246 theCoulombFacPiMinus[1]= 1.0;
247 dp.SetKineticEnergy(fLowEnergy);
248 dp.SetDefinition(thePiPlus);
249 for(G4int iz=2; iz<93; ++iz) {
250 theCoulombFacPiPlus[iz] = fPion->GetElasticCrossSection(&dp, iz, theA[iz])
251 /CoulombFactorPiPlus(fLowEnergy, iz);
252 }
253 dp.SetDefinition(G4PionMinus::PionMinus());
254 for(G4int iz=2; iz<93; ++iz) {
255 theCoulombFacPiMinus[iz] = fPion->GetElasticCrossSection(&dp, iz, theA[iz])
256 /FactorPiMinus(fLowEnergy);
257
258 if(verboseLevel > 0) {
259 G4cout << "Z= " << iz << " A= " << theA[iz]
260 << " CoulombFactorPiPlus= " << theCoulombFacPiPlus[iz]
261 << " CoulombFactorPiMinus= " << theCoulombFacPiMinus[iz]
262 << G4endl;
263 }
264 }
265 }
266}
double A(double temperature)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetElasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
const G4String & GetParticleName() const
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:97
G4double GetElasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A) const
void BuildPhysicsTable(const G4ParticleDefinition &) final
int G4lrint(double ad)
Definition: templates.hh:134

Referenced by G4QMDReaction::G4QMDReaction().

◆ CrossSectionDescription()

void G4BGGPionElasticXS::CrossSectionDescription ( std::ostream &  outFile) const
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 286 of file G4BGGPionElasticXS.cc.

287{
288 outFile << "The Barashenkov-Glauber-Gribov cross section handles elastic\n"
289 << "scattering of pions from nuclei at all energies. The\n"
290 << "Barashenkov parameterization is used below 91 GeV and the\n"
291 << "Glauber-Gribov parameterization is used above 91 GeV.\n";
292}

◆ GetElementCrossSection()

G4double G4BGGPionElasticXS::GetElementCrossSection ( const G4DynamicParticle dp,
G4int  Z,
const G4Material mat 
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 119 of file G4BGGPionElasticXS.cc.

121{
122 // this method should be called only for Z > 1
123 G4double cross = 0.0;
124 G4double ekin = std::max(dp->GetKineticEnergy(), fLowestEnergy);
125 G4int Z = std::min(ZZ, 92);
126 if(1 == Z) {
127 cross = 1.0115*GetIsoCrossSection(dp,1,1);
128 } else {
129 if(ekin <= fLowEnergy) {
130 cross = (isPiplus) ? theCoulombFacPiPlus[Z]*CoulombFactorPiPlus(ekin, Z)
131 : theCoulombFacPiMinus[Z]*FactorPiMinus(ekin);
132 } else if(ekin > fGlauberEnergy) {
133 cross = (isPiplus) ? theGlauberFacPiPlus[Z] : theGlauberFacPiMinus[Z];
134 cross *= fGlauber->GetElasticGlauberGribov(dp, Z, theA[Z]);
135 } else {
136 cross = fPion->GetElasticCrossSection(dp, Z, theA[Z]);
137 }
138 }
139 if(verboseLevel > 1) {
140 G4cout << "G4BGGPionElasticXS::GetElementCrossSection for "
142 << " Ekin(GeV)= " << dp->GetKineticEnergy()
143 << " in nucleus Z= " << Z << " A= " << theA[Z]
144 << " XS(b)= " << cross/barn
145 << G4endl;
146 }
147 return cross;
148}
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const

Referenced by G4QMDReaction::ApplyYourself().

◆ GetIsoCrossSection()

G4double G4BGGPionElasticXS::GetIsoCrossSection ( const G4DynamicParticle dp,
G4int  Z,
G4int  A,
const G4Isotope iso = nullptr,
const G4Element elm = nullptr,
const G4Material mat = nullptr 
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 151 of file G4BGGPionElasticXS.cc.

156{
157 // this method should be called only for Z = 1
158 fHadron->HadronNucleonXscNS(dp->GetDefinition(), theProton,
159 dp->GetKineticEnergy());
160 G4double cross = A*fHadron->GetElasticHadronNucleonXsc();
161
162 if(verboseLevel > 1) {
163 G4cout << "G4BGGPionElasticXS::GetIsoCrossSection for "
165 << " Ekin(GeV)= " << dp->GetKineticEnergy()
166 << " in nucleus Z= " << Z << " A= " << A
167 << " XS(b)= " << cross/barn
168 << G4endl;
169 }
170 return cross;
171}
G4double GetElasticHadronNucleonXsc() const
G4double HadronNucleonXscNS(const G4ParticleDefinition *theParticle, const G4ParticleDefinition *nucleon, G4double ekin)

Referenced by GetElementCrossSection().

◆ IsElementApplicable()

G4bool G4BGGPionElasticXS::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material  
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 101 of file G4BGGPionElasticXS.cc.

103{
104 return true;
105}

◆ IsIsoApplicable()

G4bool G4BGGPionElasticXS::IsIsoApplicable ( const G4DynamicParticle ,
G4int  Z,
G4int  A,
const G4Element elm,
const G4Material mat 
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 109 of file G4BGGPionElasticXS.cc.

112{
113 return (1 == Z);
114}

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