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

#include <G4BGGPionInelasticXS.hh>

+ Inheritance diagram for G4BGGPionInelasticXS:

Public Member Functions

 G4BGGPionInelasticXS (const G4ParticleDefinition *)
 
 ~G4BGGPionInelasticXS () final
 
G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat) 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 65 of file G4BGGPionInelasticXS.hh.

Constructor & Destructor Documentation

◆ G4BGGPionInelasticXS()

G4BGGPionInelasticXS::G4BGGPionInelasticXS ( const G4ParticleDefinition p)
explicit

Definition at line 66 of file G4BGGPionInelasticXS.cc.

67 : G4VCrossSectionDataSet("BarashenkovGlauberGribov")
68{
69 verboseLevel = 0;
70 fGlauberEnergy = 91.*CLHEP::GeV;
71 fLowEnergy = 20.*CLHEP::MeV;
72 fLowestEnergy = 1.*CLHEP::MeV;
73 SetMinKinEnergy(0.0);
75
76 fPion = nullptr;
77 fGlauber = nullptr;
78 fHadron = nullptr;
79
80 fG4pow = G4Pow::GetInstance();
81
82 theProton = G4Proton::Proton();
83 thePiPlus = G4PionPlus::PionPlus();
84 isPiplus = (p == thePiPlus);
85 isMaster = false;
87}
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)

◆ ~G4BGGPionInelasticXS()

G4BGGPionInelasticXS::~G4BGGPionInelasticXS ( )
final

Definition at line 91 of file G4BGGPionInelasticXS.cc.

92{
93 delete fHadron;
94}

Member Function Documentation

◆ BuildPhysicsTable()

void G4BGGPionInelasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 176 of file G4BGGPionInelasticXS.cc.

177{
178 if(fPion) { return; }
179 if(verboseLevel > 1) {
180 G4cout << "G4BGGPionInelasticXS::BuildPhysicsTable for "
181 << p.GetParticleName() << G4endl;
182 }
183 if(&p == G4PionPlus::PionPlus() || &p == G4PionMinus::PionMinus()) {
184 isPiplus = (&p == G4PionPlus::PionPlus());
185 } else {
187 ed << "This BGG cross section is applicable only to pions and not to "
188 << p.GetParticleName() << G4endl;
189 G4Exception("G4BGGPionInelasticXS::BuildPhysicsTable", "had001",
190 FatalException, ed);
191 return;
192 }
193
194 fPion = new G4UPiNuclearCrossSection();
195 fGlauber = new G4ComponentGGHadronNucleusXsc();
196 fHadron = new G4HadronNucleonXsc();
197
198 fPion->BuildPhysicsTable(p);
199
200 if(0 == theA[0]) {
201#ifdef G4MULTITHREADED
202 G4MUTEXLOCK(&pionInelasticXSMutex);
203 if(0 == theA[0]) {
204#endif
205 isMaster = true;
206#ifdef G4MULTITHREADED
207 }
208 G4MUTEXUNLOCK(&pionInelasticXSMutex);
209#endif
210 } else {
211 return;
212 }
213
214 if(isMaster && 0 == theA[0]) {
215
216 theA[0] = theA[1] = 1;
217 G4ThreeVector mom(0.0,0.0,1.0);
218 G4DynamicParticle dp(thePiPlus, mom, fGlauberEnergy);
219
221 G4double csup, csdn;
222
223 if(verboseLevel > 0) {
224 G4cout << "### G4BGGPionInelasticXS::Initialise for "
225 << p.GetParticleName()
226 << " isPiplus: " << isPiplus
227 << G4endl;
228 }
229 for(G4int iz=2; iz<93; ++iz) {
230 G4int A = G4lrint(nist->GetAtomicMassAmu(iz));
231 theA[iz] = A;
232
233 csup = fGlauber->GetInelasticGlauberGribov(&dp, iz, A);
234 csdn = fPion->GetInelasticCrossSection(&dp, iz, A);
235 theGlauberFacPiPlus[iz] = csdn/csup;
236 }
237
238 dp.SetDefinition(G4PionMinus::PionMinus());
239 for(G4int iz=2; iz<93; ++iz) {
240 csup = fGlauber->GetInelasticGlauberGribov(&dp, iz, theA[iz]);
241 csdn = fPion->GetInelasticCrossSection(&dp, iz, theA[iz]);
242 theGlauberFacPiMinus[iz] = csdn/csup;
243
244 if(verboseLevel > 0) {
245 G4cout << "Z= " << iz << " A= " << theA[iz]
246 << " factorPiPlus= " << theGlauberFacPiPlus[iz]
247 << " factorPiMinus= " << theGlauberFacPiMinus[iz]
248 << G4endl;
249 }
250 }
251
252 theLowEPiPlus[1] = theLowEPiMinus[1]= 1.0;
253 dp.SetDefinition(thePiPlus);
254 dp.SetKineticEnergy(fLowEnergy);
255 for(G4int iz=2; iz<93; ++iz) {
256 theLowEPiPlus[iz] = fPion->GetInelasticCrossSection(&dp, iz, theA[iz])
257 /CoulombFactorPiPlus(fLowEnergy, iz);
258 }
259
260 dp.SetDefinition(G4PionMinus::PionMinus());
261 for(G4int iz=2; iz<93; ++iz) {
262 theLowEPiMinus[iz] = fPion->GetInelasticCrossSection(&dp, iz, theA[iz])
263 /FactorPiMinus(fLowEnergy);
264
265 if(verboseLevel > 0) {
266 G4cout << "Z= " << iz << " A= " << theA[iz]
267 << " LowEtorPiPlus= " << theLowEPiPlus[iz]
268 << " LowEtorPiMinus= " << theLowEPiMinus[iz]
269 << G4endl;
270 }
271 }
272 }
273}
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 GetInelasticGlauberGribov(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
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4double GetInelasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A) const
int G4lrint(double ad)
Definition: templates.hh:134

Referenced by G4QMDReaction::G4QMDReaction().

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 293 of file G4BGGPionInelasticXS.cc.

294{
295 outFile << "The Barashenkov-Glauber-Gribov cross section handles inelastic\n"
296 << "pion scattering from nuclei at all energies. The Barashenkov\n"
297 << "parameterization is used below 91 GeV and the Glauber-Gribov\n"
298 << "parameterization is used above 91 GeV.\n";
299}

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 118 of file G4BGGPionInelasticXS.cc.

120{
121 // this method should be called only for Z > 1
122
123 G4double cross = 0.0;
124 G4double ekin = std::max(dp->GetKineticEnergy(), fLowestEnergy);
125 G4int Z = std::min(ZZ, 92);
126
127 if(1 == Z) {
128 cross = 1.0115*GetIsoCrossSection(dp,1,1);
129 } else if(ekin < fLowEnergy) {
130 cross = (isPiplus) ? theLowEPiPlus[Z]*CoulombFactorPiPlus(ekin, Z)
131 : theLowEPiMinus[Z]*FactorPiMinus(ekin);
132 } else if(ekin > fGlauberEnergy) {
133 cross = (isPiplus) ? theGlauberFacPiPlus[Z] : theGlauberFacPiMinus[Z];
134 cross *= fGlauber->GetInelasticGlauberGribov(dp, Z, theA[Z]);
135 } else {
136 cross = fPion->GetInelasticCrossSection(dp, Z, theA[Z]);
137 }
138 if(verboseLevel > 1) {
139 G4cout << "G4BGGPionInelasticXS::GetCrossSection for "
141 << " Ekin(GeV)= " << dp->GetKineticEnergy()
142 << " in nucleus Z= " << Z << " A= " << theA[Z]
143 << " XS(b)= " << cross/barn
144 << G4endl;
145 }
146 return cross;
147}
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 G4BGGPionInelasticXS::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 152 of file G4BGGPionInelasticXS.cc.

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

Referenced by GetElementCrossSection().

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 99 of file G4BGGPionInelasticXS.cc.

101{
102 return true;
103}

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 107 of file G4BGGPionInelasticXS.cc.

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

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