Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
G4ElectroVDNuclearModel Class Reference

#include <G4ElectroVDNuclearModel.hh>

Inheritance diagram for G4ElectroVDNuclearModel:
G4HadronicInteraction

Public Member Functions

void ActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Material *aMaterial)
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
void DeActivateFor (const G4Element *anElement)
 
void DeActivateFor (const G4Material *aMaterial)
 
 G4ElectroVDNuclearModel ()
 
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
 
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
 
virtual void ModelDescription (std::ostream &outFile) const
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
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)
 
 ~G4ElectroVDNuclearModel ()
 

Protected Member Functions

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

Protected Attributes

G4bool isBlocked
 
G4double theMaxEnergy
 
G4double theMinEnergy
 
G4HadFinalState theParticleChange
 
G4int verboseLevel
 

Private Member Functions

G4DynamicParticleCalculateEMVertex (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
 
void CalculateHadronicVertex (G4DynamicParticle *incident, G4Nucleus &target)
 

Private Attributes

G4CascadeInterfacebert
 
G4ElectroNuclearCrossSectionelectroXS
 
std::pair< G4double, G4doubleepCheckLevels
 
G4TheoFSGeneratorftfp
 
G4PhotoNuclearCrossSectiongammaXS
 
G4double leptonKE
 
G4double photonEnergy
 
G4double photonQ2
 
G4double recoilEnergyThreshold
 
G4HadronicInteractionRegistryregistry
 
G4int secID
 
std::vector< const G4Material * > theBlockedList
 
std::vector< const G4Element * > theBlockedListElements
 
G4LundStringFragmentationtheFragmentation
 
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
 
G4ExcitedStringDecaytheStringDecay
 

Detailed Description

Definition at line 49 of file G4ElectroVDNuclearModel.hh.

Constructor & Destructor Documentation

◆ G4ElectroVDNuclearModel()

G4ElectroVDNuclearModel::G4ElectroVDNuclearModel ( )

Definition at line 63 of file G4ElectroVDNuclearModel.cc.

64 : G4HadronicInteraction("G4ElectroVDNuclearModel"),
65 leptonKE(0.0), photonEnergy(0.0), photonQ2(0.0), secID(-1)
66{
67 SetMinEnergy(0.0);
69 electroXS =
71 GetCrossSectionDataSet(G4ElectroNuclearCrossSection::Default_Name());
72 gammaXS =
74 GetCrossSectionDataSet(G4PhotoNuclearCrossSection::Default_Name());
75
76 // reuse existing pre-compound model
77 G4GeneratorPrecompoundInterface* precoInterface
81 G4VPreCompoundModel* pre = static_cast<G4VPreCompoundModel*>(p);
82 if(!pre) { pre = new G4PreCompoundModel(); }
83 precoInterface->SetDeExcitation(pre);
84
85 // string model
86 ftfp = new G4TheoFSGenerator();
87 ftfp->SetTransport(precoInterface);
90 G4FTFModel* theStringModel = new G4FTFModel();
92 ftfp->SetHighEnergyGenerator(theStringModel);
93
94 // Build Bertini model
96
97 // Creator model ID
99}
static constexpr double PeV
Definition: G4SIunits.hh:205
static G4CrossSectionDataSetRegistry * Instance()
G4ExcitedStringDecay * theStringDecay
G4LundStringFragmentation * theFragmentation
G4PhotoNuclearCrossSection * gammaXS
G4ElectroNuclearCrossSection * electroXS
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
void SetMinEnergy(G4double anEnergy)
G4HadronicInteraction(const G4String &modelName="HadronicModel")
const G4String & GetModelName() const
void SetMaxEnergy(const G4double anEnergy)
static G4int GetModelID(const G4int modelIndex)
void SetTransport(G4VIntraNuclearTransportModel *const value)
void SetHighEnergyGenerator(G4VHighEnergyGenerator *const value)
void SetDeExcitation(G4VPreCompoundModel *ptr)
void SetFragmentationModel(G4VStringFragmentation *aModel)

References bert, G4ElectroNuclearCrossSection::Default_Name(), G4PhotoNuclearCrossSection::Default_Name(), electroXS, G4HadronicInteractionRegistry::FindModel(), ftfp, gammaXS, G4PhysicsModelCatalog::GetModelID(), G4HadronicInteraction::GetModelName(), G4CrossSectionDataSetRegistry::Instance(), G4HadronicInteractionRegistry::Instance(), PeV, secID, G4VIntraNuclearTransportModel::SetDeExcitation(), G4VPartonStringModel::SetFragmentationModel(), G4TheoFSGenerator::SetHighEnergyGenerator(), G4HadronicInteraction::SetMaxEnergy(), G4HadronicInteraction::SetMinEnergy(), G4TheoFSGenerator::SetTransport(), theFragmentation, and theStringDecay.

◆ ~G4ElectroVDNuclearModel()

G4ElectroVDNuclearModel::~G4ElectroVDNuclearModel ( )

Definition at line 101 of file G4ElectroVDNuclearModel.cc.

102{
103 delete theFragmentation;
104 delete theStringDecay;
105}

References theFragmentation, and theStringDecay.

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 * G4ElectroVDNuclearModel::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus aTargetNucleus 
)
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 123 of file G4ElectroVDNuclearModel.cc.

125{
126 // Set up default particle change (just returns initial state)
129 leptonKE = aTrack.GetKineticEnergy();
132
133 // Set up sanity checks for real photon production
134 G4DynamicParticle lepton(aTrack.GetDefinition(), aTrack.Get4Momentum() );
135
136 // Need to call GetElementCrossSection before calling GetEquivalentPhotonEnergy.
137 G4Material* mat = 0;
138 G4int targZ = targetNucleus.GetZ_asInt();
139 electroXS->GetElementCrossSection(&lepton, targZ, mat);
140
142 // Photon energy cannot exceed lepton energy
143 if (photonEnergy < leptonKE) {
146 // Photon
147 if (photonEnergy > photonQ2/dM) {
148 // Produce recoil lepton and transferred photon
149 G4DynamicParticle* transferredPhoton = CalculateEMVertex(aTrack, targetNucleus);
150 // Interact gamma with nucleus
151 if (transferredPhoton) CalculateHadronicVertex(transferredPhoton, targetNucleus);
152 }
153 }
154 return &theParticleChange;
155}
static const G4double dM
@ isAlive
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
Hep3Vector unit() const
Hep3Vector vect() const
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat)
void CalculateHadronicVertex(G4DynamicParticle *incident, G4Nucleus &target)
G4DynamicParticle * CalculateEMVertex(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
void SetStatusChange(G4HadFinalStateStatus aS)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4Proton * Proton()
Definition: G4Proton.cc:92

References CalculateEMVertex(), CalculateHadronicVertex(), G4HadFinalState::Clear(), dM, electroXS, G4HadProjectile::Get4Momentum(), G4HadProjectile::GetDefinition(), G4ElectroNuclearCrossSection::GetElementCrossSection(), G4ElectroNuclearCrossSection::GetEquivalentPhotonEnergy(), G4ElectroNuclearCrossSection::GetEquivalentPhotonQ2(), G4HadProjectile::GetKineticEnergy(), G4ParticleDefinition::GetPDGMass(), G4Nucleus::GetZ_asInt(), isAlive, leptonKE, G4Neutron::Neutron(), photonEnergy, photonQ2, G4Proton::Proton(), G4HadFinalState::SetEnergyChange(), G4HadFinalState::SetMomentumChange(), G4HadFinalState::SetStatusChange(), G4HadronicInteraction::theParticleChange, CLHEP::Hep3Vector::unit(), and CLHEP::HepLorentzVector::vect().

◆ Block()

void G4HadronicInteraction::Block ( )
inlineprotectedinherited

◆ BuildPhysicsTable()

void G4HadronicInteraction::BuildPhysicsTable ( const G4ParticleDefinition )
virtualinherited

◆ CalculateEMVertex()

G4DynamicParticle * G4ElectroVDNuclearModel::CalculateEMVertex ( const G4HadProjectile aTrack,
G4Nucleus aTargetNucleus 
)
private

Definition at line 159 of file G4ElectroVDNuclearModel.cc.

161{
163 G4ThreeVector(0.,0.,1.) );
164
165 // Get gamma cross section at Q**2 = 0 (real gamma)
166 G4int targZ = targetNucleus.GetZ_asInt();
167 G4Material* mat = 0;
168 G4double sigNu =
170
171 // Change real gamma energy to equivalent energy and get cross section at that energy
173 photon.SetKineticEnergy(photonEnergy - photonQ2/dM);
174 G4double sigK =
177
178 // No gamma produced, return null ptr
179 if (sigNu*G4UniformRand() > sigK*rndFraction) return 0;
180
181 // Scatter the lepton
182 G4double mProj = aTrack.GetDefinition()->GetPDGMass();
183 G4double mProj2 = mProj*mProj;
184 G4double iniE = leptonKE + mProj; // Total energy of incident lepton
185 G4double finE = iniE - photonEnergy; // Total energy of scattered lepton
187 G4double iniP = std::sqrt(iniE*iniE-mProj2); // Incident lepton momentum
188 G4double finP = std::sqrt(finE*finE-mProj2); // Scattered lepton momentum
189 G4double cost = (iniE*finE - mProj2 - photonQ2/2.)/iniP/finP; // cos(theta) from Q**2
190 if (cost > 1.) cost= 1.;
191 if (cost < -1.) cost=-1.;
192 G4double sint = std::sqrt(1.-cost*cost);
193
194 G4ThreeVector dir = aTrack.Get4Momentum().vect().unit();
195 G4ThreeVector ortx = dir.orthogonal().unit(); // Ortho-normal to scattering plane
196 G4ThreeVector orty = dir.cross(ortx); // Third unit vector
198 G4double sinx = sint*std::sin(phi);
199 G4double siny = sint*std::cos(phi);
200 G4ThreeVector findir = cost*dir+sinx*ortx+siny*orty;
201 theParticleChange.SetMomentumChange(findir); // change lepton direction
202
203 // Create a gamma with momentum equal to momentum transfer
204 G4ThreeVector photonMomentum = iniP*dir - finP*findir;
206 photonEnergy, photonMomentum);
207 return gamma;
208}
static constexpr double twopi
Definition: G4SIunits.hh:56
CLHEP::Hep3Vector G4ThreeVector
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
G4double GetVirtualFactor(G4double nu, G4double Q2)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *)

References CLHEP::Hep3Vector::cross(), dM, electroXS, G4UniformRand, G4Gamma::Gamma(), gammaXS, G4HadProjectile::Get4Momentum(), G4HadProjectile::GetDefinition(), G4PhotoNuclearCrossSection::GetElementCrossSection(), G4ParticleDefinition::GetPDGMass(), G4ElectroNuclearCrossSection::GetVirtualFactor(), G4Nucleus::GetZ_asInt(), leptonKE, G4Neutron::Neutron(), CLHEP::Hep3Vector::orthogonal(), G4InuclParticleNames::photon, photonEnergy, photonQ2, G4Proton::Proton(), G4HadFinalState::SetEnergyChange(), G4HadFinalState::SetMomentumChange(), G4HadronicInteraction::theParticleChange, twopi, CLHEP::Hep3Vector::unit(), and CLHEP::HepLorentzVector::vect().

Referenced by ApplyYourself().

◆ CalculateHadronicVertex()

void G4ElectroVDNuclearModel::CalculateHadronicVertex ( G4DynamicParticle incident,
G4Nucleus target 
)
private

Definition at line 212 of file G4ElectroVDNuclearModel.cc.

214{
215 G4HadFinalState* hfs = 0;
216 G4double gammaE = incident->GetTotalEnergy();
217
218 if (gammaE < 10*GeV) {
219 G4HadProjectile projectile(*incident);
220 hfs = bert->ApplyYourself(projectile, target);
221 } else {
222 // At high energies convert incident gamma to a pion
224 G4double piMom = std::sqrt(gammaE*gammaE - piMass*piMass);
225 G4ThreeVector piMomentum(incident->GetMomentumDirection() );
226 piMomentum *= piMom;
227 G4DynamicParticle theHadron(G4PionZero::PionZero(), piMomentum);
228 G4HadProjectile projectile(theHadron);
229 hfs = ftfp->ApplyYourself(projectile, target);
230 }
231
232 delete incident;
233
234 // Assign the creator model ID to the secondaries
235 for ( size_t i = 0; i < hfs->GetNumberOfSecondaries(); ++i ) {
237 }
238
239 // Copy secondaries from sub-model to model
241}
static constexpr double GeV
Definition: G4SIunits.hh:203
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
const G4ThreeVector & GetMomentumDirection() const
G4double GetTotalEnergy() const
void AddSecondaries(const std::vector< G4HadSecondary > &addSecs)
std::size_t GetNumberOfSecondaries() const
G4HadSecondary * GetSecondary(size_t i)
void SetCreatorModelID(G4int id)
static G4PionZero * PionZero()
Definition: G4PionZero.cc:107
G4HadFinalState * ApplyYourself(const G4HadProjectile &thePrimary, G4Nucleus &theNucleus) override

References G4HadFinalState::AddSecondaries(), G4CascadeInterface::ApplyYourself(), G4TheoFSGenerator::ApplyYourself(), bert, ftfp, G4DynamicParticle::GetMomentumDirection(), G4HadFinalState::GetNumberOfSecondaries(), G4ParticleDefinition::GetPDGMass(), G4HadFinalState::GetSecondary(), G4DynamicParticle::GetTotalEnergy(), GeV, G4PionZero::PionZero(), secID, G4HadSecondary::SetCreatorModelID(), and G4HadronicInteraction::theParticleChange.

Referenced by ApplyYourself().

◆ 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(), G4EMDissociation::G4EMDissociation(), G4ExcitedStringDecay::G4ExcitedStringDecay(), G4LEHadronProtonElastic::G4LEHadronProtonElastic(), G4LENDModel::G4LENDModel(), G4LENDorBERTModel::G4LENDorBERTModel(), G4LEnp::G4LEnp(), G4LEpp::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

◆ GetVerboseLevel()

G4int G4HadronicInteraction::GetVerboseLevel ( ) const
inlineinherited

Definition at line 109 of file G4HadronicInteraction.hh.

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.

◆ ModelDescription()

void G4ElectroVDNuclearModel::ModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 107 of file G4ElectroVDNuclearModel.cc.

108{
109 outFile << "G4ElectroVDNuclearModel handles the inelastic scattering\n"
110 << "of e- and e+ from nuclei using the equivalent photon\n"
111 << "approximation in which the incoming lepton generates a\n"
112 << "virtual photon at the electromagnetic vertex, and the\n"
113 << "virtual photon is converted to a real photon. At low\n"
114 << "energies, the photon interacts directly with the nucleus\n"
115 << "using the Bertini cascade. At high energies the photon\n"
116 << "is converted to a pi0 which interacts using the FTFP\n"
117 << "model. The electro- and gamma-nuclear cross sections of\n"
118 << "M. Kossov are used to generate the virtual photon spectrum\n";
119}

◆ operator!=()

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

◆ operator==()

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

◆ SampleInvariantT()

G4double G4HadronicInteraction::SampleInvariantT ( const G4ParticleDefinition p,
G4double  plab,
G4int  Z,
G4int  A 
)
virtualinherited

◆ 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().

◆ 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(), 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::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(), 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::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

◆ bert

G4CascadeInterface* G4ElectroVDNuclearModel::bert
private

Definition at line 79 of file G4ElectroVDNuclearModel.hh.

Referenced by CalculateHadronicVertex(), and G4ElectroVDNuclearModel().

◆ electroXS

G4ElectroNuclearCrossSection* G4ElectroVDNuclearModel::electroXS
private

◆ epCheckLevels

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

◆ ftfp

G4TheoFSGenerator* G4ElectroVDNuclearModel::ftfp
private

Definition at line 76 of file G4ElectroVDNuclearModel.hh.

Referenced by CalculateHadronicVertex(), and G4ElectroVDNuclearModel().

◆ gammaXS

G4PhotoNuclearCrossSection* G4ElectroVDNuclearModel::gammaXS
private

Definition at line 74 of file G4ElectroVDNuclearModel.hh.

Referenced by CalculateEMVertex(), and G4ElectroVDNuclearModel().

◆ isBlocked

G4bool G4HadronicInteraction::isBlocked
protectedinherited

◆ leptonKE

G4double G4ElectroVDNuclearModel::leptonKE
private

Definition at line 69 of file G4ElectroVDNuclearModel.hh.

Referenced by ApplyYourself(), and CalculateEMVertex().

◆ photonEnergy

G4double G4ElectroVDNuclearModel::photonEnergy
private

Definition at line 70 of file G4ElectroVDNuclearModel.hh.

Referenced by ApplyYourself(), and CalculateEMVertex().

◆ photonQ2

G4double G4ElectroVDNuclearModel::photonQ2
private

Definition at line 71 of file G4ElectroVDNuclearModel.hh.

Referenced by ApplyYourself(), and CalculateEMVertex().

◆ recoilEnergyThreshold

G4double G4HadronicInteraction::recoilEnergyThreshold
privateinherited

◆ registry

G4HadronicInteractionRegistry* G4HadronicInteraction::registry
privateinherited

◆ secID

G4int G4ElectroVDNuclearModel::secID
private

Definition at line 81 of file G4ElectroVDNuclearModel.hh.

Referenced by CalculateHadronicVertex(), and G4ElectroVDNuclearModel().

◆ theBlockedList

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

◆ theBlockedListElements

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

◆ theFragmentation

G4LundStringFragmentation* G4ElectroVDNuclearModel::theFragmentation
private

Definition at line 77 of file G4ElectroVDNuclearModel.hh.

Referenced by G4ElectroVDNuclearModel(), and ~G4ElectroVDNuclearModel().

◆ 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

◆ 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(), 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(), G4LEpp::ApplyYourself(), G4BinaryCascade::ApplyYourself(), G4CascadeInterface::ApplyYourself(), G4LMsdGenerator::ApplyYourself(), CalculateEMVertex(), G4MuonVDNuclearModel::CalculateEMVertex(), 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().

◆ theStringDecay

G4ExcitedStringDecay* G4ElectroVDNuclearModel::theStringDecay
private

Definition at line 78 of file G4ElectroVDNuclearModel.hh.

Referenced by G4ElectroVDNuclearModel(), and ~G4ElectroVDNuclearModel().

◆ 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(), G4LEpp::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: