Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Types | Private Attributes | Static Private Attributes
G4LEpp Class Reference

#include <G4LEpp.hh>

Inheritance diagram for G4LEpp:
G4HadronElastic G4HadronicInteraction

Public Member Functions

void ActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Material *aMaterial)
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) override
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4double ComputeMomentumCMS (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
void DeActivateFor (const G4Element *anElement)
 
void DeActivateFor (const G4Material *aMaterial)
 
 G4LEpp ()
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
const G4StringGetModelName () const
 
G4double GetRecoilEnergyThreshold () const
 
G4double GetSlopeCof (const G4int pdg)
 
G4int GetVerboseLevel () const
 
virtual void InitialiseModel ()
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4bool IsBlocked (const G4Element *anElement) const
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4double LowestEnergyLimit () const
 
void ModelDescription (std::ostream &) const override
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A) override
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
void SetLowestEnergyLimit (G4double value)
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
void SetRecoilEnergyThreshold (G4double val)
 
void SetVerboseLevel (G4int value)
 
 ~G4LEpp () override
 

Protected Member Functions

void Block ()
 
G4bool IsBlocked () const
 
void SetModelName (const G4String &nam)
 

Protected Attributes

G4bool isBlocked
 
G4double pLocalTmax
 
G4int secID
 
G4double theMaxEnergy
 
G4double theMinEnergy
 
G4HadFinalState theParticleChange
 
G4int verboseLevel
 

Private Types

enum  { NENERGY =40 , NANGLE =180 }
 

Private Attributes

std::pair< G4double, G4doubleepCheckLevels
 
G4double lowestEnergyLimit
 
G4int nwarn
 
G4double recoilEnergyThreshold
 
G4HadronicInteractionRegistryregistry
 
G4ParticleDefinitiontheAlpha
 
std::vector< const G4Material * > theBlockedList
 
std::vector< const G4Element * > theBlockedListElements
 
G4ParticleDefinitiontheDeuteron
 
std::vector< std::pair< G4double, const G4Material * > > theMaxEnergyList
 
std::vector< std::pair< G4double, const G4Element * > > theMaxEnergyListElements
 
std::vector< std::pair< G4double, const G4Material * > > theMinEnergyList
 
std::vector< std::pair< G4double, const G4Element * > > theMinEnergyListElements
 
G4String theModelName
 
G4ParticleDefinitiontheNeutron
 
G4ParticleDefinitiontheProton
 

Static Private Attributes

static const G4float dSigmax [NENERGY]
 
static const G4float elab [NENERGY]
 
static const G4float Sig [NENERGY][NANGLE]
 
static const G4float Sigtot [NENERGY]
 

Detailed Description

Definition at line 41 of file G4LEpp.hh.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
private
Enumerator
NENERGY 
NANGLE 

Definition at line 45 of file G4LEpp.hh.

45{ NENERGY=40, NANGLE=180 };
@ NENERGY
Definition: G4LEpp.hh:45
@ NANGLE
Definition: G4LEpp.hh:45

Constructor & Destructor Documentation

◆ G4LEpp()

G4LEpp::G4LEpp ( )
explicit

Definition at line 44 of file G4LEpp.cc.

44 :G4HadronElastic("G4LEpp")
45{
47 SetMinEnergy(0.);
48 SetMaxEnergy(5.*GeV);
49}
static constexpr double GeV
Definition: G4SIunits.hh:203
G4HadronElastic(const G4String &name="hElasticLHEP")
void SetMinEnergy(G4double anEnergy)
const G4String & GetModelName() const
void SetMaxEnergy(const G4double anEnergy)
static G4int GetModelID(const G4int modelIndex)

References G4PhysicsModelCatalog::GetModelID(), G4HadronicInteraction::GetModelName(), GeV, G4HadronElastic::secID, G4HadronicInteraction::SetMaxEnergy(), and G4HadronicInteraction::SetMinEnergy().

◆ ~G4LEpp()

G4LEpp::~G4LEpp ( )
override

Definition at line 51 of file G4LEpp.cc.

52{}

Member Function Documentation

◆ ActivateFor() [1/2]

void G4HadronicInteraction::ActivateFor ( const G4Element anElement)
inlineinherited

◆ ActivateFor() [2/2]

void G4HadronicInteraction::ActivateFor ( const G4Material aMaterial)
inlineinherited

◆ ApplyYourself()

G4HadFinalState * G4LEpp::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
)
overridevirtual

Reimplemented from G4HadronElastic.

Definition at line 55 of file G4LEpp.cc.

56{
58 const G4HadProjectile* aParticle = &aTrack;
59
60 G4double P = aParticle->GetTotalMomentum();
61 G4double Px = aParticle->Get4Momentum().x();
62 G4double Py = aParticle->Get4Momentum().y();
63 G4double Pz = aParticle->Get4Momentum().z();
64 G4double E = aParticle->GetTotalEnergy();
65 G4ThreeVector theInitial = aParticle->Get4Momentum().vect().unit();
66
67 if (verboseLevel > 1) {
68 G4double ek = aParticle->GetKineticEnergy();
69 G4double E0 = aParticle->GetDefinition()->GetPDGMass();
70 G4double Q = aParticle->GetDefinition()->GetPDGCharge();
71 G4int A = targetNucleus.GetA_asInt();
72 G4int Z = targetNucleus.GetZ_asInt();
73 G4cout << "G4LEpp:ApplyYourself: incident particle: "
74 << aParticle->GetDefinition()->GetParticleName() << G4endl;
75 G4cout << "P = " << P/GeV << " GeV/c"
76 << ", Px = " << Px/GeV << " GeV/c"
77 << ", Py = " << Py/GeV << " GeV/c"
78 << ", Pz = " << Pz/GeV << " GeV/c" << G4endl;
79 G4cout << "E = " << E/GeV << " GeV"
80 << ", kinetic energy = " << ek/GeV << " GeV"
81 << ", mass = " << E0/GeV << " GeV"
82 << ", charge = " << Q << G4endl;
83 G4cout << "G4LEpp:ApplyYourself: material:" << G4endl;
84 G4cout << "A = " << A
85 << ", Z = " << Z
86 << ", atomic mass "
87 << G4Proton::Proton()->GetPDGMass()/GeV << "GeV"
88 << G4endl;
89 //
90 // GHEISHA ADD operation to get total energy, mass, charge
91 //
92 E += proton_mass_c2;
93 G4double E02 = E*E - P*P;
94 E0 = std::sqrt(std::fabs(E02));
95 if (E02 < 0)E0 *= -1;
96 Q += Z;
97 G4cout << "G4LEpp:ApplyYourself: total:" << G4endl;
98 G4cout << "E = " << E/GeV << " GeV"
99 << ", mass = " << E0/GeV << " GeV"
100 << ", charge = " << Q << G4endl;
101 }
102 G4double t = SampleInvariantT(aParticle->GetDefinition(), P, 0, 0);
103 G4double cost = 1.0 - 2*t/(P*P);
104 if(cost > 1.0) { cost = 1.0; }
105 if(cost <-1.0) { cost =-1.0; }
106 G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
108 // Get the target particle
109 G4DynamicParticle* targetParticle = targetNucleus.ReturnTargetParticle();
110
111 G4double E1 = aParticle->GetTotalEnergy();
112 G4double M1 = aParticle->GetDefinition()->GetPDGMass();
113 G4double E2 = targetParticle->GetTotalEnergy();
114 G4double M2 = targetParticle->GetDefinition()->GetPDGMass();
115 G4double totalEnergy = E1 + E2;
116 G4double pseudoMass = std::sqrt(totalEnergy*totalEnergy - P*P);
117
118 // Transform into centre of mass system
119
120 G4double px = (M2/pseudoMass)*Px;
121 G4double py = (M2/pseudoMass)*Py;
122 G4double pz = (M2/pseudoMass)*Pz;
123 G4double p = std::sqrt(px*px + py*py + pz*pz);
124
125 if (verboseLevel > 1) {
126 G4cout << " E1, M1 (GeV) " << E1/GeV << " " << M1/GeV << G4endl;
127 G4cout << " E2, M2 (GeV) " << E2/GeV << " " << M2/GeV << G4endl;
128 G4cout << " particle 1 momentum in CM " << px/GeV
129 << " " << py/GeV << " "
130 << pz/GeV << " " << p/GeV << G4endl;
131 }
132
133 // First scatter w.r.t. Z axis
134 G4double pxnew = p*sint*std::cos(phi);
135 G4double pynew = p*sint*std::sin(phi);
136 G4double pznew = p*cost;
137
138 // Rotate according to the direction of the incident particle
139 if (px*px + py*py > 0) {
140 G4double ph, cosp, sinp;
141 cost = pz/p;
142 sint = (std::sqrt((1-cost)*(1+cost)) + std::sqrt(px*px+py*py)/p)/2;
143 py < 0 ? ph = 3*halfpi : ph = halfpi;
144 if (std::fabs(px) > 0.000001*GeV) ph = std::atan2(py,px);
145 cosp = std::cos(ph);
146 sinp = std::sin(ph);
147 px = (cost*cosp*pxnew - sinp*pynew + sint*cosp*pznew);
148 py = (cost*sinp*pxnew + cosp*pynew + sint*sinp*pznew);
149 pz = (-sint*pxnew + cost*pznew);
150 }
151 else {
152 px = pxnew;
153 py = pynew;
154 pz = pznew;
155 }
156
157 if (verboseLevel > 1) {
158 G4cout << " AFTER SCATTER..." << G4endl;
159 G4cout << " particle 1 momentum in CM " << px/GeV << " " << py/GeV << " "
160 << pz/GeV << " " << p/GeV << G4endl;
161 }
162
163 // Transform to lab system
164
165 G4double E1pM2 = E1 + M2;
166 G4double betaCM = P/E1pM2;
167 G4double betaCMx = Px/E1pM2;
168 G4double betaCMy = Py/E1pM2;
169 G4double betaCMz = Pz/E1pM2;
170 G4double gammaCM = E1pM2/std::sqrt(E1pM2*E1pM2 - P*P);
171
172 if (verboseLevel > 1) {
173 G4cout << " betaCM " << betaCMx << " " << betaCMy << " "
174 << betaCMz << " " << betaCM << G4endl;
175 G4cout << " gammaCM " << gammaCM << G4endl;
176 }
177
178 // Now following GLOREN...
179
180 G4double BETA[5], PA[5], PB[5];
181 BETA[1] = -betaCMx;
182 BETA[2] = -betaCMy;
183 BETA[3] = -betaCMz;
184 BETA[4] = gammaCM;
185
186 //The incident particle...
187
188 PA[1] = px;
189 PA[2] = py;
190 PA[3] = pz;
191 PA[4] = std::sqrt(M1*M1 + p*p);
192
193 G4double BETPA = BETA[1]*PA[1] + BETA[2]*PA[2] + BETA[3]*PA[3];
194 G4double BPGAM = (BETPA * BETA[4]/(BETA[4] + 1.) - PA[4]) * BETA[4];
195
196 PB[1] = PA[1] + BPGAM * BETA[1];
197 PB[2] = PA[2] + BPGAM * BETA[2];
198 PB[3] = PA[3] + BPGAM * BETA[3];
199 PB[4] = (PA[4] - BETPA) * BETA[4];
200
202 newP->SetDefinition(aParticle->GetDefinition());
203 newP->SetMomentum(G4ThreeVector(PB[1], PB[2], PB[3]));
204
205 //The target particle...
206
207 PA[1] = -px;
208 PA[2] = -py;
209 PA[3] = -pz;
210 PA[4] = std::sqrt(M2*M2 + p*p);
211
212 BETPA = BETA[1]*PA[1] + BETA[2]*PA[2] + BETA[3]*PA[3];
213 BPGAM = (BETPA * BETA[4]/(BETA[4] + 1.) - PA[4]) * BETA[4];
214
215 PB[1] = PA[1] + BPGAM * BETA[1];
216 PB[2] = PA[2] + BPGAM * BETA[2];
217 PB[3] = PA[3] + BPGAM * BETA[3];
218 PB[4] = (PA[4] - BETPA) * BETA[4];
219
220 targetParticle->SetMomentum(G4ThreeVector(PB[1], PB[2], PB[3]));
221
222 if (verboseLevel > 1) {
223 G4cout << " particle 1 momentum in LAB "
224 << newP->GetMomentum()/GeV
225 << " " << newP->GetTotalMomentum()/GeV << G4endl;
226 G4cout << " particle 2 momentum in LAB "
227 << targetParticle->GetMomentum()/GeV
228 << " " << targetParticle->GetTotalMomentum()/GeV << G4endl;
229 G4cout << " TOTAL momentum in LAB "
230 << (newP->GetMomentum()+targetParticle->GetMomentum())/GeV
231 << " "
232 << (newP->GetMomentum()+targetParticle->GetMomentum()).mag()/GeV
233 << G4endl;
234 }
235
238 delete newP;
239
240 // Recoil particle
241 theParticleChange.AddSecondary(targetParticle, secID);
242 return &theParticleChange;
243}
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double halfpi
Definition: G4SIunits.hh:57
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector vect() const
const G4ThreeVector & GetMomentumDirection() const
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
void SetMomentum(const G4ThreeVector &momentum)
G4ThreeVector GetMomentum() const
G4double GetTotalMomentum() const
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
G4double GetTotalMomentum() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetTotalEnergy() const
G4double SampleInvariantT(const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A) override
Definition: G4LEpp.cc:249
G4int GetA_asInt() const
Definition: G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105
G4DynamicParticle * ReturnTargetParticle() const
Definition: G4Nucleus.cc:340
G4double GetPDGCharge() const
const G4String & GetParticleName() const
static G4Proton * Proton()
Definition: G4Proton.cc:92
float proton_mass_c2
Definition: hepunit.py:274
static double Q[]
static double P[]

References A, G4HadFinalState::AddSecondary(), G4HadFinalState::Clear(), G4cout, G4endl, G4UniformRand, G4HadProjectile::Get4Momentum(), G4Nucleus::GetA_asInt(), G4DynamicParticle::GetDefinition(), G4HadProjectile::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4HadProjectile::GetKineticEnergy(), G4DynamicParticle::GetMomentum(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), G4DynamicParticle::GetTotalEnergy(), G4HadProjectile::GetTotalEnergy(), G4DynamicParticle::GetTotalMomentum(), G4HadProjectile::GetTotalMomentum(), G4Nucleus::GetZ_asInt(), GeV, halfpi, P, G4Proton::Proton(), source.hepunit::proton_mass_c2, Q, G4Nucleus::ReturnTargetParticle(), SampleInvariantT(), G4HadronElastic::secID, G4DynamicParticle::SetDefinition(), G4HadFinalState::SetEnergyChange(), G4DynamicParticle::SetMomentum(), G4HadFinalState::SetMomentumChange(), G4HadronicInteraction::theParticleChange, twopi, CLHEP::Hep3Vector::unit(), CLHEP::HepLorentzVector::vect(), G4HadronicInteraction::verboseLevel, CLHEP::HepLorentzVector::x(), CLHEP::HepLorentzVector::y(), CLHEP::HepLorentzVector::z(), and Z.

◆ Block()

void G4HadronicInteraction::Block ( )
inlineprotectedinherited

◆ BuildPhysicsTable()

void G4HadronicInteraction::BuildPhysicsTable ( const G4ParticleDefinition )
virtualinherited

◆ ComputeMomentumCMS()

G4double G4HadronElastic::ComputeMomentumCMS ( const G4ParticleDefinition p,
G4double  plab,
G4int  Z,
G4int  A 
)
inlineinherited

Definition at line 102 of file G4HadronElastic.hh.

104{
105 G4double m1 = p->GetPDGMass();
106 G4double m12= m1*m1;
108 return plab*mass2/std::sqrt(m12 + mass2*mass2 + 2.*mass2*std::sqrt(m12 + plab*plab));
109}
static G4double GetNuclearMass(const G4double A, const G4double Z)

References A, G4NucleiProperties::GetNuclearMass(), G4ParticleDefinition::GetPDGMass(), and Z.

◆ DeActivateFor() [1/2]

void G4HadronicInteraction::DeActivateFor ( const G4Element anElement)
inherited

Definition at line 186 of file G4HadronicInteraction.cc.

187{
188 Block();
189 theBlockedListElements.push_back(anElement);
190}
std::vector< const G4Element * > theBlockedListElements

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theBlockedListElements.

◆ DeActivateFor() [2/2]

void G4HadronicInteraction::DeActivateFor ( const G4Material aMaterial)
inherited

Definition at line 180 of file G4HadronicInteraction.cc.

181{
182 Block();
183 theBlockedList.push_back(aMaterial);
184}
std::vector< const G4Material * > theBlockedList

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theBlockedList.

Referenced by G4HadronHElasticPhysics::ConstructProcess().

◆ GetEnergyMomentumCheckLevels()

std::pair< G4double, G4double > G4HadronicInteraction::GetEnergyMomentumCheckLevels ( ) const
virtualinherited

◆ GetFatalEnergyCheckLevels()

const std::pair< G4double, G4double > G4HadronicInteraction::GetFatalEnergyCheckLevels ( ) const
virtualinherited

Reimplemented in G4FissLib, G4LFission, G4LENDFission, G4ParticleHPCapture, G4ParticleHPElastic, G4ParticleHPFission, G4ParticleHPInelastic, and G4ParticleHPThermalScattering.

Definition at line 210 of file G4HadronicInteraction.cc.

211{
212 // default level of Check
213 return std::pair<G4double, G4double>(2.*perCent, 1. * GeV);
214}
static constexpr double perCent
Definition: G4SIunits.hh:325

References GeV, and perCent.

Referenced by G4HadronicProcess::CheckResult().

◆ GetMaxEnergy() [1/2]

G4double G4HadronicInteraction::GetMaxEnergy ( ) const
inlineinherited

◆ GetMaxEnergy() [2/2]

G4double G4HadronicInteraction::GetMaxEnergy ( const G4Material aMaterial,
const G4Element anElement 
) const
inherited

Definition at line 131 of file G4HadronicInteraction.cc.

133{
134 if(!IsBlocked()) { return theMaxEnergy; }
135 if( IsBlocked(aMaterial) || IsBlocked(anElement) ) { return 0.0; }
136 if(!theMaxEnergyListElements.empty()) {
137 for(auto const& elmlist : theMaxEnergyListElements) {
138 if( anElement == elmlist.second )
139 { return elmlist.first; }
140 }
141 }
142 if(!theMaxEnergyList.empty()) {
143 for(auto const& matlist : theMaxEnergyList) {
144 if( aMaterial == matlist.second )
145 { return matlist.first; }
146 }
147 }
148 return theMaxEnergy;
149}
std::vector< std::pair< G4double, const G4Material * > > theMaxEnergyList
std::vector< std::pair< G4double, const G4Element * > > theMaxEnergyListElements

References G4HadronicInteraction::IsBlocked(), G4HadronicInteraction::theMaxEnergy, G4HadronicInteraction::theMaxEnergyList, and G4HadronicInteraction::theMaxEnergyListElements.

◆ GetMinEnergy() [1/2]

G4double G4HadronicInteraction::GetMinEnergy ( ) const
inlineinherited

◆ GetMinEnergy() [2/2]

G4double G4HadronicInteraction::GetMinEnergy ( const G4Material aMaterial,
const G4Element anElement 
) const
inherited

Definition at line 81 of file G4HadronicInteraction.cc.

83{
84 if(!IsBlocked()) { return theMinEnergy; }
85 if( IsBlocked(aMaterial) || IsBlocked(anElement) ) { return DBL_MAX; }
86 if(!theMinEnergyListElements.empty()) {
87 for(auto const& elmlist : theMinEnergyListElements) {
88 if( anElement == elmlist.second )
89 { return elmlist.first; }
90 }
91 }
92 if(!theMinEnergyList.empty()) {
93 for(auto const & matlist : theMinEnergyList) {
94 if( aMaterial == matlist.second )
95 { return matlist.first; }
96 }
97 }
98 return theMinEnergy;
99}
std::vector< std::pair< G4double, const G4Element * > > theMinEnergyListElements
std::vector< std::pair< G4double, const G4Material * > > theMinEnergyList
#define DBL_MAX
Definition: templates.hh:62

References DBL_MAX, G4HadronicInteraction::IsBlocked(), G4HadronicInteraction::theMinEnergy, G4HadronicInteraction::theMinEnergyList, and G4HadronicInteraction::theMinEnergyListElements.

◆ GetModelName()

const G4String & G4HadronicInteraction::GetModelName ( ) const
inlineinherited

Definition at line 115 of file G4HadronicInteraction.hh.

116 { return theModelName; }

References G4HadronicInteraction::theModelName.

Referenced by G4MuMinusCapturePrecompound::ApplyYourself(), G4HadronElastic::ApplyYourself(), G4INCLXXInterface::ApplyYourself(), G4TheoFSGenerator::ApplyYourself(), G4HadronStoppingProcess::AtRestDoIt(), G4VHadronPhysics::BuildModel(), G4HadronicProcess::CheckEnergyMomentumConservation(), G4HadronicProcess::CheckResult(), G4ChargeExchangePhysics::ConstructProcess(), G4MuonicAtomDecay::DecayIt(), G4LENDModel::DumpLENDTargetInfo(), G4AblaInterface::G4AblaInterface(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4EMDissociation::G4EMDissociation(), G4ExcitedStringDecay::G4ExcitedStringDecay(), G4LEHadronProtonElastic::G4LEHadronProtonElastic(), G4LENDModel::G4LENDModel(), G4LENDorBERTModel::G4LENDorBERTModel(), G4LEnp::G4LEnp(), G4LEpp(), G4LFission::G4LFission(), G4LowEGammaNuclearModel::G4LowEGammaNuclearModel(), G4LowEIonFragmentation::G4LowEIonFragmentation(), G4MuonVDNuclearModel::G4MuonVDNuclearModel(), G4NeutrinoElectronCcModel::G4NeutrinoElectronCcModel(), G4NeutrinoNucleusModel::G4NeutrinoNucleusModel(), G4WilsonAbrasionModel::G4WilsonAbrasionModel(), G4INCLXXInterface::GetDeExcitationModelName(), G4EnergyRangeManager::GetHadronicInteraction(), G4VHighEnergyGenerator::GetProjectileNucleus(), G4NeutronRadCapture::InitialiseModel(), G4BinaryCascade::ModelDescription(), G4LMsdGenerator::ModelDescription(), G4VPartonStringModel::ModelDescription(), G4TheoFSGenerator::ModelDescription(), G4VHadronPhysics::NewModel(), G4NeutrinoElectronProcess::PostStepDoIt(), G4HadronicProcess::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4HadronicProcessStore::PrintModelHtml(), G4BinaryCascade::PropagateModelDescription(), G4HadronicProcessStore::RegisterInteraction(), and G4LENDModel::returnUnchanged().

◆ GetRecoilEnergyThreshold()

G4double G4HadronicInteraction::GetRecoilEnergyThreshold ( ) const
inlineinherited

◆ GetSlopeCof()

G4double G4HadronElastic::GetSlopeCof ( const G4int  pdg)
inherited

Definition at line 277 of file G4HadronElastic.cc.

278{
279 // The input parameter "pdg" should be the absolute value of the PDG code
280 // (i.e. the same value for a particle and its antiparticle).
281
282 G4double coeff = 1.0;
283
284 // heavy barions
285
286 static const G4double lBarCof1S = 0.88;
287 static const G4double lBarCof2S = 0.76;
288 static const G4double lBarCof3S = 0.64;
289 static const G4double lBarCof1C = 0.784378;
290 static const G4double lBarCofSC = 0.664378;
291 static const G4double lBarCof2SC = 0.544378;
292 static const G4double lBarCof1B = 0.740659;
293 static const G4double lBarCofSB = 0.620659;
294 static const G4double lBarCof2SB = 0.500659;
295
296 if( pdg == 3122 || pdg == 3222 || pdg == 3112 || pdg == 3212 )
297 {
298 coeff = lBarCof1S; // Lambda, Sigma+, Sigma-, Sigma0
299
300 } else if( pdg == 3322 || pdg == 3312 )
301 {
302 coeff = lBarCof2S; // Xi-, Xi0
303 }
304 else if( pdg == 3324)
305 {
306 coeff = lBarCof3S; // Omega
307 }
308 else if( pdg == 4122 || pdg == 4212 || pdg == 4222 || pdg == 4112 )
309 {
310 coeff = lBarCof1C; // LambdaC+, SigmaC+, SigmaC++, SigmaC0
311 }
312 else if( pdg == 4332 )
313 {
314 coeff = lBarCof2SC; // OmegaC
315 }
316 else if( pdg == 4232 || pdg == 4132 )
317 {
318 coeff = lBarCofSC; // XiC+, XiC0
319 }
320 else if( pdg == 5122 || pdg == 5222 || pdg == 5112 || pdg == 5212 )
321 {
322 coeff = lBarCof1B; // LambdaB, SigmaB+, SigmaB-, SigmaB0
323 }
324 else if( pdg == 5332 )
325 {
326 coeff = lBarCof2SB; // OmegaB-
327 }
328 else if( pdg == 5132 || pdg == 5232 ) // XiB-, XiB0
329 {
330 coeff = lBarCofSB;
331 }
332 // heavy mesons Kaons?
333 static const G4double lMesCof1S = 0.82; // Kp/piP kaons?
334 static const G4double llMesCof1C = 0.676568;
335 static const G4double llMesCof1B = 0.610989;
336 static const G4double llMesCof2C = 0.353135;
337 static const G4double llMesCof2B = 0.221978;
338 static const G4double llMesCofSC = 0.496568;
339 static const G4double llMesCofSB = 0.430989;
340 static const G4double llMesCofCB = 0.287557;
341 static const G4double llMesCofEtaP = 0.88;
342 static const G4double llMesCofEta = 0.76;
343
344 if( pdg == 321 || pdg == 311 || pdg == 310 )
345 {
346 coeff = lMesCof1S; //K+-0
347 }
348 else if( pdg == 511 || pdg == 521 )
349 {
350 coeff = llMesCof1B; // BMeson0, BMeson+
351 }
352 else if(pdg == 421 || pdg == 411 )
353 {
354 coeff = llMesCof1C; // DMeson+, DMeson0
355 }
356 else if( pdg == 531 )
357 {
358 coeff = llMesCofSB; // BSMeson0
359 }
360 else if( pdg == 541 )
361 {
362 coeff = llMesCofCB; // BCMeson+-
363 }
364 else if(pdg == 431 )
365 {
366 coeff = llMesCofSC; // DSMeson+-
367 }
368 else if(pdg == 441 || pdg == 443 )
369 {
370 coeff = llMesCof2C; // Etac, JPsi
371 }
372 else if(pdg == 553 )
373 {
374 coeff = llMesCof2B; // Upsilon
375 }
376 else if(pdg == 221 )
377 {
378 coeff = llMesCofEta; // Eta
379 }
380 else if(pdg == 331 )
381 {
382 coeff = llMesCofEtaP; // Eta'
383 }
384 return coeff;
385}

◆ GetVerboseLevel()

G4int G4HadronicInteraction::GetVerboseLevel ( ) const
inlineinherited

Definition at line 109 of file G4HadronicInteraction.hh.

110 { return verboseLevel; }

References G4HadronicInteraction::verboseLevel.

◆ InitialiseModel()

void G4HadronicInteraction::InitialiseModel ( )
virtualinherited

◆ IsApplicable()

G4bool G4HadronicInteraction::IsApplicable ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
)
virtualinherited

◆ IsBlocked() [1/3]

G4bool G4HadronicInteraction::IsBlocked ( ) const
inlineprotectedinherited

◆ IsBlocked() [2/3]

G4bool G4HadronicInteraction::IsBlocked ( const G4Element anElement) const
inherited

Definition at line 202 of file G4HadronicInteraction.cc.

203{
204 for (auto const& elm : theBlockedListElements) {
205 if (anElement == elm) return true;
206 }
207 return false;
208}

References G4HadronicInteraction::theBlockedListElements.

◆ IsBlocked() [3/3]

G4bool G4HadronicInteraction::IsBlocked ( const G4Material aMaterial) const
inherited

Definition at line 193 of file G4HadronicInteraction.cc.

194{
195 for (auto const& mat : theBlockedList) {
196 if (aMaterial == mat) return true;
197 }
198 return false;
199}

References G4HadronicInteraction::theBlockedList.

◆ LowestEnergyLimit()

G4double G4HadronElastic::LowestEnergyLimit ( ) const
inlineinherited

◆ ModelDescription()

void G4HadronElastic::ModelDescription ( std::ostream &  outFile) const
overridevirtualinherited

Reimplemented from G4HadronicInteraction.

Reimplemented in G4NeutrinoElectronNcModel, and G4NeutronElectronElModel.

Definition at line 70 of file G4HadronElastic.cc.

71{
72 outFile << "G4HadronElastic is the base class for all hadron-nucleus\n"
73 << "elastic scattering models except HP.\n"
74 << "By default it uses the Gheisha two-exponential momentum\n"
75 << "transfer parameterization. The model is fully relativistic\n"
76 << "as opposed to the original Gheisha model which was not.\n"
77 << "This model may be used for all long-lived hadrons at all\n"
78 << "incident energies but fit the data only for relativistic scattering.\n";
79}

◆ operator!=()

G4bool G4HadronicInteraction::operator!= ( const G4HadronicInteraction right) const
deleteinherited

◆ operator==()

G4bool G4HadronicInteraction::operator== ( const G4HadronicInteraction right) const
deleteinherited

◆ SampleInvariantT()

G4double G4LEpp::SampleInvariantT ( const G4ParticleDefinition p,
G4double  plab,
G4int  Z,
G4int  A 
)
overridevirtual

Reimplemented from G4HadronElastic.

Definition at line 249 of file G4LEpp.cc.

251{
252 G4double nMass = p->GetPDGMass(); // 939.565346*MeV;
253 G4double ek = std::sqrt(plab*plab+nMass*nMass) - nMass;
254
255 // Find energy bin
256
257 G4int je1 = 0;
258 G4int je2 = NENERGY - 1;
259 ek /= GeV;
260
261 do
262 {
263 G4int midBin = (je1 + je2)/2;
264
265 if (ek < elab[midBin]) je2 = midBin;
266 else je1 = midBin;
267 }
268 while (je2 - je1 > 1); /* Loop checking, 10.08.2015, A.Ribon */
269
270 G4double delab = elab[je2] - elab[je1];
271
272 // Sample the angle
273
274 G4double sample = G4UniformRand();
275 G4int ke1 = 0;
276 G4int ke2 = NANGLE - 1;
277 G4double dsig, b, rc;
278
279 dsig = Sig[je2][0] - Sig[je1][0];
280 rc = dsig/delab;
281 b = Sig[je1][0] - rc*elab[je1];
282
283 G4double sigint1 = rc*ek + b;
284 G4double sigint2 = 0.;
285
286 do
287 {
288 G4int midBin = (ke1 + ke2)/2;
289 dsig = Sig[je2][midBin] - Sig[je1][midBin];
290 rc = dsig/delab;
291 b = Sig[je1][midBin] - rc*elab[je1];
292 G4double sigint = rc*ek + b;
293
294 if (sample < sigint)
295 {
296 ke2 = midBin;
297 sigint2 = sigint;
298 }
299 else
300 {
301 ke1 = midBin;
302 sigint1 = sigint;
303 }
304 }
305 while (ke2 - ke1 > 1); /* Loop checking, 10.08.2015, A.Ribon */
306
307 dsig = sigint2 - sigint1;
308 rc = 1./dsig;
309 b = ke1 - rc*sigint1;
310
311 G4double kint = rc*sample + b;
312 G4double theta = (0.5 + kint)*pi/180.;
313 G4double t = 0.5*plab*plab*(1 - std::cos(theta));
314
315 return t;
316}
static constexpr double pi
Definition: G4SIunits.hh:55
static const G4float Sig[NENERGY][NANGLE]
Definition: G4LEpp.hh:65
static const G4float elab[NENERGY]
Definition: G4LEpp.hh:66

References elab, G4UniformRand, G4ParticleDefinition::GetPDGMass(), GeV, NANGLE, NENERGY, pi, and Sig.

Referenced by ApplyYourself().

◆ SetEnergyMomentumCheckLevels()

void G4HadronicInteraction::SetEnergyMomentumCheckLevels ( G4double  relativeLevel,
G4double  absoluteLevel 
)
inlineinherited

Definition at line 149 of file G4HadronicInteraction.hh.

150 { epCheckLevels.first = relativeLevel;
151 epCheckLevels.second = absoluteLevel; }

References G4HadronicInteraction::epCheckLevels.

Referenced by G4BinaryCascade::G4BinaryCascade(), G4CascadeInterface::G4CascadeInterface(), and G4FTFModel::G4FTFModel().

◆ SetLowestEnergyLimit()

void G4HadronElastic::SetLowestEnergyLimit ( G4double  value)
inlineinherited

◆ SetMaxEnergy() [1/3]

void G4HadronicInteraction::SetMaxEnergy ( const G4double  anEnergy)
inlineinherited

Definition at line 102 of file G4HadronicInteraction.hh.

103 { theMaxEnergy = anEnergy; }

References G4HadronicInteraction::theMaxEnergy.

Referenced by G4HadronicInteraction::ActivateFor(), G4IonINCLXXPhysics::AddProcess(), G4BertiniElectroNuclearBuilder::Build(), G4LENDBertiniGammaElectroNuclearBuilder::Build(), G4NeutronLENDBuilder::Build(), G4NeutronPHPBuilder::Build(), G4AlphaPHPBuilder::Build(), G4BertiniKaonBuilder::Build(), G4BertiniNeutronBuilder::Build(), G4BertiniPiKBuilder::Build(), G4BertiniPionBuilder::Build(), G4BertiniProtonBuilder::Build(), G4BinaryAlphaBuilder::Build(), G4BinaryDeuteronBuilder::Build(), G4BinaryHe3Builder::Build(), G4BinaryNeutronBuilder::Build(), G4BinaryPiKBuilder::Build(), G4BinaryPionBuilder::Build(), G4BinaryProtonBuilder::Build(), G4BinaryTritonBuilder::Build(), G4DeuteronPHPBuilder::Build(), G4FTFBinaryKaonBuilder::Build(), G4FTFBinaryNeutronBuilder::Build(), G4FTFBinaryPionBuilder::Build(), G4FTFBinaryProtonBuilder::Build(), G4FTFPAntiBarionBuilder::Build(), G4FTFPKaonBuilder::Build(), G4FTFPNeutronBuilder::Build(), G4FTFPPiKBuilder::Build(), G4FTFPPionBuilder::Build(), G4FTFPProtonBuilder::Build(), G4He3PHPBuilder::Build(), G4HyperonFTFPBuilder::Build(), G4HyperonQGSPBuilder::Build(), G4INCLXXNeutronBuilder::Build(), G4INCLXXPionBuilder::Build(), G4INCLXXProtonBuilder::Build(), G4PrecoNeutronBuilder::Build(), G4PrecoProtonBuilder::Build(), G4ProtonPHPBuilder::Build(), G4QGSBinaryKaonBuilder::Build(), G4QGSBinaryNeutronBuilder::Build(), G4QGSBinaryPiKBuilder::Build(), G4QGSBinaryPionBuilder::Build(), G4QGSBinaryProtonBuilder::Build(), G4QGSPAntiBarionBuilder::Build(), G4QGSPKaonBuilder::Build(), G4QGSPLundStrFragmProtonBuilder::Build(), G4QGSPNeutronBuilder::Build(), G4QGSPPiKBuilder::Build(), G4QGSPPionBuilder::Build(), G4TritonPHPBuilder::Build(), G4QGSPProtonBuilder::Build(), G4HadronicBuilder::BuildFTFP_BERT(), G4HadronicBuilder::BuildFTFQGSP_BERT(), G4QGSBuilder::BuildModel(), G4VHadronPhysics::BuildModel(), G4HadronicBuilder::BuildQGSP_FTFP_BERT(), G4EmExtraPhysics::ConstructGammaElectroNuclear(), LBE::ConstructHad(), G4EmExtraPhysics::ConstructLENDGammaNuclear(), G4HadronDElasticPhysics::ConstructProcess(), G4HadronElasticPhysics::ConstructProcess(), G4HadronHElasticPhysics::ConstructProcess(), G4IonINCLXXPhysics::ConstructProcess(), G4IonPhysics::ConstructProcess(), G4IonPhysicsPHP::ConstructProcess(), G4IonQMDPhysics::ConstructProcess(), G4ANuElNucleusNcModel::G4ANuElNucleusNcModel(), G4ANuMuNucleusNcModel::G4ANuMuNucleusNcModel(), G4BertiniKaonBuilder::G4BertiniKaonBuilder(), G4BertiniPiKBuilder::G4BertiniPiKBuilder(), G4BertiniPionBuilder::G4BertiniPionBuilder(), G4BinaryCascade::G4BinaryCascade(), G4BinaryPiKBuilder::G4BinaryPiKBuilder(), G4BinaryPionBuilder::G4BinaryPionBuilder(), G4ChargeExchange::G4ChargeExchange(), G4DiffuseElastic::G4DiffuseElastic(), G4DiffuseElasticV2::G4DiffuseElasticV2(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4EMDissociation::G4EMDissociation(), G4FissLib::G4FissLib(), G4FTFBinaryKaonBuilder::G4FTFBinaryKaonBuilder(), G4FTFBinaryNeutronBuilder::G4FTFBinaryNeutronBuilder(), G4FTFBinaryPiKBuilder::G4FTFBinaryPiKBuilder(), G4FTFBinaryPionBuilder::G4FTFBinaryPionBuilder(), G4FTFBinaryProtonBuilder::G4FTFBinaryProtonBuilder(), G4FTFPAntiBarionBuilder::G4FTFPAntiBarionBuilder(), G4FTFPKaonBuilder::G4FTFPKaonBuilder(), G4FTFPNeutronBuilder::G4FTFPNeutronBuilder(), G4FTFPPiKBuilder::G4FTFPPiKBuilder(), G4FTFPPionBuilder::G4FTFPPionBuilder(), G4FTFPProtonBuilder::G4FTFPProtonBuilder(), G4HadronElastic::G4HadronElastic(), G4HadronicAbsorptionFritiof::G4HadronicAbsorptionFritiof(), G4HadronicAbsorptionFritiofWithBinaryCascade::G4HadronicAbsorptionFritiofWithBinaryCascade(), G4hhElastic::G4hhElastic(), G4HyperonFTFPBuilder::G4HyperonFTFPBuilder(), G4HyperonQGSPBuilder::G4HyperonQGSPBuilder(), G4INCLXXPionBuilder::G4INCLXXPionBuilder(), G4LEHadronProtonElastic::G4LEHadronProtonElastic(), G4LENDModel::G4LENDModel(), G4LEnp::G4LEnp(), G4LEpp(), G4LFission::G4LFission(), G4LowEGammaNuclearModel::G4LowEGammaNuclearModel(), G4MuonVDNuclearModel::G4MuonVDNuclearModel(), G4NeutrinoElectronCcModel::G4NeutrinoElectronCcModel(), G4NeutrinoElectronNcModel::G4NeutrinoElectronNcModel(), G4NeutrinoNucleusModel::G4NeutrinoNucleusModel(), G4NeutronElectronElModel::G4NeutronElectronElModel(), G4NeutronRadCapture::G4NeutronRadCapture(), G4NuclNuclDiffuseElastic::G4NuclNuclDiffuseElastic(), G4NuElNucleusNcModel::G4NuElNucleusNcModel(), G4NuMuNucleusNcModel::G4NuMuNucleusNcModel(), G4ParticleHPCapture::G4ParticleHPCapture(), G4ParticleHPElastic::G4ParticleHPElastic(), G4ParticleHPFission::G4ParticleHPFission(), G4ParticleHPInelastic::G4ParticleHPInelastic(), G4ParticleHPThermalScattering::G4ParticleHPThermalScattering(), G4QGSPAntiBarionBuilder::G4QGSPAntiBarionBuilder(), G4WilsonAbrasionModel::G4WilsonAbrasionModel(), G4HadronPhysicsFTFP_BERT_HP::Neutron(), G4HadronPhysicsINCLXX::Neutron(), G4HadronPhysicsQGSP_BERT_HP::Neutron(), G4HadronPhysicsQGSP_BIC_HP::Neutron(), G4HadronPhysicsShielding::Neutron(), and G4VHadronPhysics::NewModel().

◆ SetMaxEnergy() [2/3]

void G4HadronicInteraction::SetMaxEnergy ( G4double  anEnergy,
const G4Element anElement 
)
inherited

Definition at line 151 of file G4HadronicInteraction.cc.

153{
154 Block();
155 if(!theMaxEnergyListElements.empty()) {
156 for(auto & elmlist : theMaxEnergyListElements) {
157 if( anElement == elmlist.second ) {
158 elmlist.first = anEnergy;
159 return;
160 }
161 }
162 }
163 theMaxEnergyListElements.push_back(std::pair<G4double, const G4Element *>(anEnergy, anElement));
164}

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theMaxEnergyListElements.

◆ SetMaxEnergy() [3/3]

void G4HadronicInteraction::SetMaxEnergy ( G4double  anEnergy,
const G4Material aMaterial 
)
inherited

Definition at line 166 of file G4HadronicInteraction.cc.

167{
168 Block();
169 if(!theMaxEnergyList.empty()) {
170 for(auto & matlist: theMaxEnergyList) {
171 if( aMaterial == matlist.second ) {
172 matlist.first = anEnergy;
173 return;
174 }
175 }
176 }
177 theMaxEnergyList.push_back(std::pair<G4double, const G4Material *>(anEnergy, aMaterial));
178}

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theMaxEnergyList.

◆ SetMinEnergy() [1/3]

void G4HadronicInteraction::SetMinEnergy ( G4double  anEnergy)
inlineinherited

Definition at line 89 of file G4HadronicInteraction.hh.

90 { theMinEnergy = anEnergy; }

References G4HadronicInteraction::theMinEnergy.

Referenced by G4HadronicInteraction::ActivateFor(), G4BertiniElectroNuclearBuilder::Build(), G4LENDBertiniGammaElectroNuclearBuilder::Build(), G4NeutronLENDBuilder::Build(), G4NeutronPHPBuilder::Build(), G4AlphaPHPBuilder::Build(), G4BertiniKaonBuilder::Build(), G4BertiniNeutronBuilder::Build(), G4BertiniPiKBuilder::Build(), G4BertiniPionBuilder::Build(), G4BertiniProtonBuilder::Build(), G4BinaryAlphaBuilder::Build(), G4BinaryDeuteronBuilder::Build(), G4BinaryHe3Builder::Build(), G4BinaryNeutronBuilder::Build(), G4BinaryPiKBuilder::Build(), G4BinaryPionBuilder::Build(), G4BinaryProtonBuilder::Build(), G4BinaryTritonBuilder::Build(), G4DeuteronPHPBuilder::Build(), G4FTFBinaryKaonBuilder::Build(), G4FTFBinaryNeutronBuilder::Build(), G4FTFBinaryPiKBuilder::Build(), G4FTFBinaryPionBuilder::Build(), G4FTFBinaryProtonBuilder::Build(), G4FTFPAntiBarionBuilder::Build(), G4FTFPKaonBuilder::Build(), G4FTFPNeutronBuilder::Build(), G4FTFPPiKBuilder::Build(), G4FTFPPionBuilder::Build(), G4FTFPProtonBuilder::Build(), G4He3PHPBuilder::Build(), G4HyperonFTFPBuilder::Build(), G4HyperonQGSPBuilder::Build(), G4INCLXXNeutronBuilder::Build(), G4INCLXXPionBuilder::Build(), G4INCLXXProtonBuilder::Build(), G4PrecoNeutronBuilder::Build(), G4PrecoProtonBuilder::Build(), G4ProtonPHPBuilder::Build(), G4QGSBinaryKaonBuilder::Build(), G4QGSBinaryNeutronBuilder::Build(), G4QGSBinaryPiKBuilder::Build(), G4QGSBinaryPionBuilder::Build(), G4QGSBinaryProtonBuilder::Build(), G4QGSPAntiBarionBuilder::Build(), G4QGSPKaonBuilder::Build(), G4QGSPLundStrFragmProtonBuilder::Build(), G4QGSPNeutronBuilder::Build(), G4QGSPPiKBuilder::Build(), G4QGSPPionBuilder::Build(), G4TritonPHPBuilder::Build(), G4QGSPProtonBuilder::Build(), G4QGSBuilder::BuildModel(), G4VHadronPhysics::BuildModel(), G4EmExtraPhysics::ConstructGammaElectroNuclear(), LBE::ConstructHad(), G4EmExtraPhysics::ConstructLENDGammaNuclear(), G4HadronElasticPhysicsHP::ConstructProcess(), G4HadronElasticPhysicsLEND::ConstructProcess(), G4HadronElasticPhysicsPHP::ConstructProcess(), G4HadronDElasticPhysics::ConstructProcess(), G4HadronElasticPhysics::ConstructProcess(), G4HadronHElasticPhysics::ConstructProcess(), G4IonElasticPhysics::ConstructProcess(), G4IonINCLXXPhysics::ConstructProcess(), G4IonPhysics::ConstructProcess(), G4IonPhysicsPHP::ConstructProcess(), G4IonQMDPhysics::ConstructProcess(), G4ANuElNucleusNcModel::G4ANuElNucleusNcModel(), G4ANuMuNucleusNcModel::G4ANuMuNucleusNcModel(), G4BertiniKaonBuilder::G4BertiniKaonBuilder(), G4BertiniPiKBuilder::G4BertiniPiKBuilder(), G4BertiniPionBuilder::G4BertiniPionBuilder(), G4BinaryCascade::G4BinaryCascade(), G4BinaryPiKBuilder::G4BinaryPiKBuilder(), G4BinaryPionBuilder::G4BinaryPionBuilder(), G4ChargeExchange::G4ChargeExchange(), G4DiffuseElastic::G4DiffuseElastic(), G4DiffuseElasticV2::G4DiffuseElasticV2(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4EMDissociation::G4EMDissociation(), G4FissLib::G4FissLib(), G4FTFBinaryKaonBuilder::G4FTFBinaryKaonBuilder(), G4FTFBinaryNeutronBuilder::G4FTFBinaryNeutronBuilder(), G4FTFBinaryPiKBuilder::G4FTFBinaryPiKBuilder(), G4FTFBinaryPionBuilder::G4FTFBinaryPionBuilder(), G4FTFBinaryProtonBuilder::G4FTFBinaryProtonBuilder(), G4FTFPAntiBarionBuilder::G4FTFPAntiBarionBuilder(), G4FTFPKaonBuilder::G4FTFPKaonBuilder(), G4FTFPNeutronBuilder::G4FTFPNeutronBuilder(), G4FTFPPiKBuilder::G4FTFPPiKBuilder(), G4FTFPPionBuilder::G4FTFPPionBuilder(), G4FTFPProtonBuilder::G4FTFPProtonBuilder(), G4HadronElastic::G4HadronElastic(), G4HadronicAbsorptionBertini::G4HadronicAbsorptionBertini(), G4HadronicAbsorptionFritiof::G4HadronicAbsorptionFritiof(), G4HadronicAbsorptionFritiofWithBinaryCascade::G4HadronicAbsorptionFritiofWithBinaryCascade(), G4hhElastic::G4hhElastic(), G4HyperonFTFPBuilder::G4HyperonFTFPBuilder(), G4HyperonQGSPBuilder::G4HyperonQGSPBuilder(), G4INCLXXPionBuilder::G4INCLXXPionBuilder(), G4LEHadronProtonElastic::G4LEHadronProtonElastic(), G4LENDModel::G4LENDModel(), G4LEnp::G4LEnp(), G4LEpp(), G4LFission::G4LFission(), G4LowEGammaNuclearModel::G4LowEGammaNuclearModel(), G4MuonVDNuclearModel::G4MuonVDNuclearModel(), G4NeutrinoElectronCcModel::G4NeutrinoElectronCcModel(), G4NeutrinoElectronNcModel::G4NeutrinoElectronNcModel(), G4NeutrinoNucleusModel::G4NeutrinoNucleusModel(), G4NeutronElectronElModel::G4NeutronElectronElModel(), G4NeutronRadCapture::G4NeutronRadCapture(), G4NuclNuclDiffuseElastic::G4NuclNuclDiffuseElastic(), G4NuElNucleusNcModel::G4NuElNucleusNcModel(), G4NuMuNucleusNcModel::G4NuMuNucleusNcModel(), G4ParticleHPCapture::G4ParticleHPCapture(), G4ParticleHPElastic::G4ParticleHPElastic(), G4ParticleHPFission::G4ParticleHPFission(), G4ParticleHPInelastic::G4ParticleHPInelastic(), G4ParticleHPThermalScattering::G4ParticleHPThermalScattering(), G4QGSPAntiBarionBuilder::G4QGSPAntiBarionBuilder(), G4WilsonAbrasionModel::G4WilsonAbrasionModel(), G4NeutrinoElectronCcModel::IsApplicable(), G4HadronPhysicsFTFP_BERT_HP::Neutron(), G4HadronPhysicsINCLXX::Neutron(), G4HadronPhysicsQGSP_BERT_HP::Neutron(), G4HadronPhysicsQGSP_BIC_HP::Neutron(), G4HadronPhysicsShielding::Neutron(), and G4VHadronPhysics::NewModel().

◆ SetMinEnergy() [2/3]

void G4HadronicInteraction::SetMinEnergy ( G4double  anEnergy,
const G4Element anElement 
)
inherited

Definition at line 101 of file G4HadronicInteraction.cc.

103{
104 Block();
105 if(!theMinEnergyListElements.empty()) {
106 for(auto & elmlist : theMinEnergyListElements) {
107 if( anElement == elmlist.second ) {
108 elmlist.first = anEnergy;
109 return;
110 }
111 }
112 }
113 theMinEnergyListElements.push_back(std::pair<G4double, const G4Element *>(anEnergy, anElement));
114}

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theMinEnergyListElements.

◆ SetMinEnergy() [3/3]

void G4HadronicInteraction::SetMinEnergy ( G4double  anEnergy,
const G4Material aMaterial 
)
inherited

Definition at line 116 of file G4HadronicInteraction.cc.

118{
119 Block();
120 if(!theMinEnergyList.empty()) {
121 for(auto & matlist : theMinEnergyList) {
122 if( aMaterial == matlist.second ) {
123 matlist.first = anEnergy;
124 return;
125 }
126 }
127 }
128 theMinEnergyList.push_back(std::pair<G4double, const G4Material *>(anEnergy, aMaterial));
129}

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theMinEnergyList.

◆ SetModelName()

void G4HadronicInteraction::SetModelName ( const G4String nam)
inlineprotectedinherited

◆ SetRecoilEnergyThreshold()

void G4HadronicInteraction::SetRecoilEnergyThreshold ( G4double  val)
inlineinherited

◆ SetVerboseLevel()

void G4HadronicInteraction::SetVerboseLevel ( G4int  value)
inlineinherited

Field Documentation

◆ dSigmax

const G4float G4LEpp::dSigmax
staticprivate
Initial value:
= {
63.9f, 29.9f, 14.5f, 8.97f, 7.23f, 6.66f, 6.30f,
5.89f, 5.68f, 5.57f, 5.56f, 5.57f, 5.74f, 6.38f,
7.66f, 9.65f, 12.4f, 15.6f, 18.0f, 20.1f, 22.0f, 23.8f,
26.710f, 28.910f, 31.000f, 32.990f, 34.900f, 36.730f, 38.520f,
40.270f, 43.680f, 46.950f, 50.100f, 53.160f, 55.880f, 58.820f,
57.420f, 57.820f, 58.590f, 57.800f
}

Definition at line 67 of file G4LEpp.hh.

◆ elab

const G4float G4LEpp::elab
staticprivate
Initial value:
= {
0.100E-01f, 0.200E-01f, 0.400E-01f, 0.700E-01f, 0.100f, 0.120f, 0.140f,
0.180f, 0.220f, 0.260f, 0.280f, 0.300f, 0.340f, 0.420f,
0.520f, 0.620f, 0.700f, 0.800f, 0.900f, 1.00f, 1.10f, 1.20f,
1.30f, 1.40f, 1.50f, 1.60f, 1.70f, 1.80f, 1.90f,
2.00f, 2.20f, 2.40f, 2.60f, 2.80f, 3.00f, 3.40f,
3.80f, 4.20f, 4.60f, 5.00f
}

Definition at line 66 of file G4LEpp.hh.

Referenced by SampleInvariantT().

◆ epCheckLevels

std::pair<G4double, G4double> G4HadronicInteraction::epCheckLevels
privateinherited

◆ isBlocked

G4bool G4HadronicInteraction::isBlocked
protectedinherited

◆ lowestEnergyLimit

G4double G4HadronElastic::lowestEnergyLimit
privateinherited

◆ nwarn

G4int G4HadronElastic::nwarn
privateinherited

◆ pLocalTmax

G4double G4HadronElastic::pLocalTmax
protectedinherited

◆ recoilEnergyThreshold

G4double G4HadronicInteraction::recoilEnergyThreshold
privateinherited

◆ registry

G4HadronicInteractionRegistry* G4HadronicInteraction::registry
privateinherited

◆ secID

G4int G4HadronElastic::secID
protectedinherited

◆ Sig

const G4float G4LEpp::Sig
staticprivate

Definition at line 65 of file G4LEpp.hh.

Referenced by SampleInvariantT().

◆ Sigtot

const G4float G4LEpp::Sigtot
staticprivate
Initial value:
= {
396.f, 179.f, 80.7f, 44.5f, 32.8f, 29.0f, 26.9f,
24.9f, 24.4f, 24.3f, 24.4f, 24.4f, 24.7f, 25.5f,
26.4f, 26.1f, 25.6f, 25.1f, 24.4f, 23.4f, 22.4f, 21.5f,
22.765f, 23.275f, 23.744f, 24.188f, 24.620f, 25.058f, 25.509f,
25.981f, 26.984f, 28.055f, 29.242f, 30.608f, 32.026f, 34.457f,
34.394f, 33.885f, 34.122f, 33.603f
}

Definition at line 68 of file G4LEpp.hh.

◆ theAlpha

G4ParticleDefinition* G4HadronElastic::theAlpha
privateinherited

◆ theBlockedList

std::vector<const G4Material *> G4HadronicInteraction::theBlockedList
privateinherited

◆ theBlockedListElements

std::vector<const G4Element *> G4HadronicInteraction::theBlockedListElements
privateinherited

◆ theDeuteron

G4ParticleDefinition* G4HadronElastic::theDeuteron
privateinherited

◆ theMaxEnergy

G4double G4HadronicInteraction::theMaxEnergy
protectedinherited

◆ theMaxEnergyList

std::vector<std::pair<G4double, const G4Material *> > G4HadronicInteraction::theMaxEnergyList
privateinherited

◆ theMaxEnergyListElements

std::vector<std::pair<G4double, const G4Element *> > G4HadronicInteraction::theMaxEnergyListElements
privateinherited

◆ theMinEnergy

G4double G4HadronicInteraction::theMinEnergy
protectedinherited

◆ theMinEnergyList

std::vector<std::pair<G4double, const G4Material *> > G4HadronicInteraction::theMinEnergyList
privateinherited

◆ theMinEnergyListElements

std::vector<std::pair<G4double, const G4Element *> > G4HadronicInteraction::theMinEnergyListElements
privateinherited

◆ theModelName

G4String G4HadronicInteraction::theModelName
privateinherited

◆ theNeutron

G4ParticleDefinition* G4HadronElastic::theNeutron
privateinherited

Definition at line 83 of file G4HadronElastic.hh.

Referenced by G4HadronElastic::G4HadronElastic().

◆ theParticleChange

G4HadFinalState G4HadronicInteraction::theParticleChange
protectedinherited

Definition at line 172 of file G4HadronicInteraction.hh.

Referenced by G4WilsonAbrasionModel::ApplyYourself(), G4EMDissociation::ApplyYourself(), G4LENDCapture::ApplyYourself(), G4LENDElastic::ApplyYourself(), G4LENDFission::ApplyYourself(), G4LENDInelastic::ApplyYourself(), G4ElectroVDNuclearModel::ApplyYourself(), G4ParticleHPThermalScattering::ApplyYourself(), G4NeutrinoElectronNcModel::ApplyYourself(), G4NeutronElectronElModel::ApplyYourself(), G4LFission::ApplyYourself(), G4ANuElNucleusCcModel::ApplyYourself(), G4ANuElNucleusNcModel::ApplyYourself(), G4ANuMuNucleusCcModel::ApplyYourself(), G4ANuMuNucleusNcModel::ApplyYourself(), G4MuonVDNuclearModel::ApplyYourself(), G4NeutrinoElectronCcModel::ApplyYourself(), G4NuElNucleusCcModel::ApplyYourself(), G4NuElNucleusNcModel::ApplyYourself(), G4NuMuNucleusCcModel::ApplyYourself(), G4NuMuNucleusNcModel::ApplyYourself(), G4QMDReaction::ApplyYourself(), G4NeutronRadCapture::ApplyYourself(), G4LowEGammaNuclearModel::ApplyYourself(), G4ChargeExchange::ApplyYourself(), G4HadronElastic::ApplyYourself(), G4LEHadronProtonElastic::ApplyYourself(), G4LEnp::ApplyYourself(), ApplyYourself(), G4BinaryCascade::ApplyYourself(), G4CascadeInterface::ApplyYourself(), G4LMsdGenerator::ApplyYourself(), G4ElectroVDNuclearModel::CalculateEMVertex(), G4MuonVDNuclearModel::CalculateEMVertex(), G4ElectroVDNuclearModel::CalculateHadronicVertex(), G4MuonVDNuclearModel::CalculateHadronicVertex(), G4NeutrinoNucleusModel::CoherentPion(), G4CascadeInterface::copyOutputToHadronicResult(), G4BinaryCascade::DebugEpConservation(), G4BinaryCascade::DebugFinalEpConservation(), G4NeutrinoNucleusModel::FinalBarion(), G4NeutrinoNucleusModel::FinalMeson(), G4WilsonAbrasionModel::GetAbradedNucleons(), G4CascadeInterface::NoInteraction(), G4CascadeInterface::Propagate(), G4NeutrinoNucleusModel::RecoilDeexcitation(), G4LEHadronProtonElastic::~G4LEHadronProtonElastic(), G4LEnp::~G4LEnp(), and G4LFission::~G4LFission().

◆ theProton

G4ParticleDefinition* G4HadronElastic::theProton
privateinherited

◆ verboseLevel

G4int G4HadronicInteraction::verboseLevel
protectedinherited

Definition at line 177 of file G4HadronicInteraction.hh.

Referenced by G4WilsonAbrasionModel::ApplyYourself(), G4EMDissociation::ApplyYourself(), G4LFission::ApplyYourself(), G4MuMinusCapturePrecompound::ApplyYourself(), G4NeutronRadCapture::ApplyYourself(), G4LowEGammaNuclearModel::ApplyYourself(), G4ChargeExchange::ApplyYourself(), G4HadronElastic::ApplyYourself(), G4LEHadronProtonElastic::ApplyYourself(), G4LEnp::ApplyYourself(), ApplyYourself(), G4CascadeInterface::ApplyYourself(), G4CascadeInterface::checkFinalResult(), G4CascadeInterface::copyOutputToHadronicResult(), G4CascadeInterface::copyOutputToReactionProducts(), G4LENDModel::create_used_target_map(), G4CascadeInterface::createBullet(), G4CascadeInterface::createTarget(), G4ElasticHadrNucleusHE::DefineHadronValues(), G4ElasticHadrNucleusHE::FillData(), G4ElasticHadrNucleusHE::FillFq2(), G4DiffuseElastic::G4DiffuseElastic(), G4DiffuseElasticV2::G4DiffuseElasticV2(), G4ElasticHadrNucleusHE::G4ElasticHadrNucleusHE(), G4EMDissociation::G4EMDissociation(), G4hhElastic::G4hhElastic(), G4NuclNuclDiffuseElastic::G4NuclNuclDiffuseElastic(), G4WilsonAbrasionModel::G4WilsonAbrasionModel(), G4ElasticHadrNucleusHE::GetFt(), G4ElasticHadrNucleusHE::GetLightFq2(), G4ElasticHadrNucleusHE::GetQ2_2(), G4HadronicInteraction::GetVerboseLevel(), G4ElasticHadrNucleusHE::HadronNucleusQ2_2(), G4ElasticHadrNucleusHE::HadronProtonQ2(), G4LFission::init(), G4DiffuseElastic::Initialise(), G4DiffuseElasticV2::Initialise(), G4NuclNuclDiffuseElastic::Initialise(), G4DiffuseElastic::InitialiseOnFly(), G4DiffuseElasticV2::InitialiseOnFly(), G4NuclNuclDiffuseElastic::InitialiseOnFly(), G4CascadeInterface::makeDynamicParticle(), G4CascadeInterface::NoInteraction(), G4CascadeInterface::Propagate(), G4ElasticHadrNucleusHE::SampleInvariantT(), G4AntiNuclElastic::SampleThetaCMS(), G4DiffuseElastic::SampleThetaLab(), G4NuclNuclDiffuseElastic::SampleThetaLab(), G4AntiNuclElastic::SampleThetaLab(), G4WilsonAbrasionModel::SetUseAblation(), G4HadronicInteraction::SetVerboseLevel(), G4WilsonAbrasionModel::SetVerboseLevel(), G4DiffuseElastic::ThetaCMStoThetaLab(), G4DiffuseElasticV2::ThetaCMStoThetaLab(), G4NuclNuclDiffuseElastic::ThetaCMStoThetaLab(), G4DiffuseElastic::ThetaLabToThetaCMS(), G4DiffuseElasticV2::ThetaLabToThetaCMS(), and G4NuclNuclDiffuseElastic::ThetaLabToThetaCMS().


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