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

#include <G4QuasiElasticChannel.hh>

Inheritance diagram for G4QuasiElasticChannel:
G4HadronicInteraction

Public Member Functions

void ActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Material *aMaterial)
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
void DeActivateFor (const G4Element *anElement)
 
void DeActivateFor (const G4Material *aMaterial)
 
 G4QuasiElasticChannel ()
 
 G4QuasiElasticChannel (const G4QuasiElasticChannel &)=delete
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
G4double GetFraction (G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
 
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 G4QuasiElasticChannel &) const =delete
 
const G4QuasiElasticChanneloperator= (const G4QuasiElasticChannel &)=delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4bool operator== (const G4QuasiElasticChannel &) const =delete
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
G4KineticTrackVectorScatter (G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
 
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)
 
 ~G4QuasiElasticChannel () override
 

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 Attributes

std::pair< G4double, G4doubleepCheckLevels
 
G4double recoilEnergyThreshold
 
G4HadronicInteractionRegistryregistry
 
G4int secID
 
G4V3DNucleusthe3DNucleus
 
std::vector< const G4Material * > theBlockedList
 
std::vector< const G4Element * > theBlockedListElements
 
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
 
G4QuasiElRatiostheQuasiElastic
 

Detailed Description

Definition at line 51 of file G4QuasiElasticChannel.hh.

Constructor & Destructor Documentation

◆ G4QuasiElasticChannel() [1/2]

G4QuasiElasticChannel::G4QuasiElasticChannel ( )
explicit

Definition at line 62 of file G4QuasiElasticChannel.cc.

63 : G4HadronicInteraction("QuasiElastic"),
66 secID(-1) {
67 secID = G4PhysicsModelCatalog::GetModelID( "model_QuasiElastic" );
68}
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static G4int GetModelID(const G4int modelIndex)
G4QuasiElRatios * theQuasiElastic

References G4PhysicsModelCatalog::GetModelID(), and secID.

◆ ~G4QuasiElasticChannel()

G4QuasiElasticChannel::~G4QuasiElasticChannel ( )
override

Definition at line 70 of file G4QuasiElasticChannel.cc.

71{
72 delete the3DNucleus;
73 delete theQuasiElastic;
74}

References the3DNucleus, and theQuasiElastic.

◆ G4QuasiElasticChannel() [2/2]

G4QuasiElasticChannel::G4QuasiElasticChannel ( const G4QuasiElasticChannel )
delete

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 * G4HadronicInteraction::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
)
virtualinherited

Reimplemented in G4WilsonAbrasionModel, G4EMDissociation, G4VLongitudinalStringDecay, G4FissLib, G4LENDorBERTModel, G4LENDCapture, G4LENDCombinedModel, G4LENDElastic, G4LENDFission, G4LENDGammaModel, G4LENDInelastic, G4LENDModel, G4ElectroVDNuclearModel, G4ParticleHPCapture, G4ParticleHPElastic, G4ParticleHPFission, G4ParticleHPInelastic, G4ParticleHPThermalScattering, G4GeneratorPrecompoundInterface, G4NeutrinoElectronNcModel, G4NeutronElectronElModel, G4LFission, G4ANuElNucleusCcModel, G4ANuElNucleusNcModel, G4ANuMuNucleusCcModel, G4ANuMuNucleusNcModel, G4MuonVDNuclearModel, G4NeutrinoElectronCcModel, G4NuElNucleusCcModel, G4NuElNucleusNcModel, G4NuMuNucleusCcModel, G4NuMuNucleusNcModel, G4QMDReaction, G4EmCaptureCascade, G4MuMinusCapturePrecompound, G4MuonMinusBoundDecay, G4NeutronRadCapture, G4LowEGammaNuclearModel, G4ChargeExchange, G4HadronElastic, G4LEHadronProtonElastic, G4LEnp, G4LEpp, G4NeutrinoNucleusModel, G4BinaryCascade, G4BinaryLightIonReaction, G4CascadeInterface, G4INCLXXInterface, G4LMsdGenerator, G4PreCompoundModel, G4LowEIonFragmentation, G4TheoFSGenerator, and G4AblaInterface.

Definition at line 63 of file G4HadronicInteraction.cc.

64{
65 return nullptr;
66}

Referenced by G4LENDorBERTModel::ApplyYourself(), G4BinaryCascade::ApplyYourself(), G4INCLXXInterface::ApplyYourself(), G4HadronStoppingProcess::AtRestDoIt(), G4MuonMinusAtomicCapture::AtRestDoIt(), G4MuonicAtomDecay::DecayIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4HadronicProcess::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), and G4MuNeutrinoNucleusProcess::PostStepDoIt().

◆ Block()

void G4HadronicInteraction::Block ( )
inlineprotectedinherited

◆ BuildPhysicsTable()

void G4HadronicInteraction::BuildPhysicsTable ( const G4ParticleDefinition )
virtualinherited

◆ 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
static constexpr double GeV
Definition: G4SIunits.hh:203

References GeV, and perCent.

Referenced by G4HadronicProcess::CheckResult().

◆ GetFraction()

G4double G4QuasiElasticChannel::GetFraction ( G4Nucleus theNucleus,
const G4DynamicParticle thePrimary 
)

Definition at line 76 of file G4QuasiElasticChannel.cc.

78{
79 #ifdef debug_scatter
80 G4cout << "G4QuasiElasticChannel:: P=" << thePrimary.GetTotalMomentum()
81 << ", pPDG=" << thePrimary.GetDefinition()->GetPDGEncoding()
82 << ", Z = " << theNucleus.GetZ_asInt())
83 << ", N = " << theNucleus.GetN_asInt())
84 << ", A = " << theNucleus.GetA_asInt() << G4endl;
85 #endif
86
87 std::pair<G4double,G4double> ratios;
88 ratios=theQuasiElastic->GetRatios(thePrimary.GetTotalMomentum(),
89 thePrimary.GetDefinition()->GetPDGEncoding(),
90 theNucleus.GetZ_asInt(),
91 theNucleus.GetN_asInt());
92 #ifdef debug_scatter
93 G4cout << "G4QuasiElasticChannel::ratios " << ratios.first << " x " <<ratios.second
94 << " = " << ratios.first*ratios.second << G4endl;
95 #endif
96
97 return ratios.first*ratios.second;
98}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ParticleDefinition * GetDefinition() const
G4double GetTotalMomentum() const
G4int GetA_asInt() const
Definition: G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105
G4int GetN_asInt() const
Definition: G4Nucleus.hh:102
std::pair< G4double, G4double > GetRatios(G4double pIU, G4int prPDG, G4int tgZ, G4int tgN)

References G4cout, G4endl, G4Nucleus::GetA_asInt(), G4DynamicParticle::GetDefinition(), G4Nucleus::GetN_asInt(), G4ParticleDefinition::GetPDGEncoding(), G4QuasiElRatios::GetRatios(), G4DynamicParticle::GetTotalMomentum(), G4Nucleus::GetZ_asInt(), and theQuasiElastic.

Referenced by G4TheoFSGenerator::ApplyYourself().

◆ 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::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 G4HadronicInteraction::ModelDescription ( std::ostream &  outFile) const
virtualinherited

◆ operator!=() [1/2]

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

◆ operator!=() [2/2]

G4bool G4QuasiElasticChannel::operator!= ( const G4QuasiElasticChannel ) const
delete

◆ operator=()

const G4QuasiElasticChannel & G4QuasiElasticChannel::operator= ( const G4QuasiElasticChannel )
delete

◆ operator==() [1/2]

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

◆ operator==() [2/2]

G4bool G4QuasiElasticChannel::operator== ( const G4QuasiElasticChannel ) const
delete

◆ SampleInvariantT()

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

◆ Scatter()

G4KineticTrackVector * G4QuasiElasticChannel::Scatter ( G4Nucleus theNucleus,
const G4DynamicParticle thePrimary 
)

Definition at line 100 of file G4QuasiElasticChannel.cc.

102{
103 G4int A=theNucleus.GetA_asInt();
104 G4int Z=theNucleus.GetZ_asInt();
105 // build Nucleus and choose random nucleon to scatter with
106 the3DNucleus->Init(theNucleus.GetA_asInt(),theNucleus.GetZ_asInt());
107 const std::vector<G4Nucleon>& nucleons=the3DNucleus->GetNucleons();
108 G4double targetNucleusMass=the3DNucleus->GetMass();
109 G4LorentzVector targetNucleus4Mom(0.,0.,0.,targetNucleusMass);
110 G4int index;
111 do {
112 index=G4lrint((A-1)*G4UniformRand());
113 } while (index < 0 || index >= static_cast<G4int>(nucleons.size())); /* Loop checking, 07.08.2015, A.Ribon */
114
115 const G4ParticleDefinition * pDef= nucleons[index].GetDefinition();
116
117 G4int resA=A - 1;
118 G4int resZ=Z - static_cast<int>(pDef->GetPDGCharge());
119 const G4ParticleDefinition* resDef;
120 G4double residualNucleusMass;
121 if(resZ)
122 {
123 resDef=G4ParticleTable::GetParticleTable()->GetIonTable()->GetIon(resZ,resA,0);
124 residualNucleusMass=resDef->GetPDGMass();
125 }
126 else {
127 resDef=G4Neutron::Neutron();
128 residualNucleusMass=resA * G4Neutron::Neutron()->GetPDGMass();
129 }
130 #ifdef debug_scatter
131 G4cout<<"G4QElChan::Scatter: neutron - proton? A ="<<A<<", Z="<<Z<<", projName="
132 <<pDef->GetParticleName()<<G4endl;
133 #endif
134
135 G4LorentzVector pNucleon=nucleons[index].Get4Momentum();
136 G4double residualNucleusEnergy=std::sqrt(sqr(residualNucleusMass) +
137 pNucleon.vect().mag2());
138 pNucleon.setE(targetNucleusMass-residualNucleusEnergy);
139 G4LorentzVector residualNucleus4Mom=targetNucleus4Mom-pNucleon;
140
141 std::pair<G4LorentzVector,G4LorentzVector> result;
142
143 result=theQuasiElastic->Scatter(pDef->GetPDGEncoding(),pNucleon,
144 thePrimary.GetDefinition()->GetPDGEncoding(),
145 thePrimary.Get4Momentum());
146 G4LorentzVector scatteredHadron4Mom;
147 if (result.first.e() > 0.)
148 scatteredHadron4Mom=result.second;
149 else { //scatter failed
150 //G4cout << "Warning - G4QuasiElasticChannel::Scatter no scattering" << G4endl;
151 //return 0; //no scatter
152 scatteredHadron4Mom=thePrimary.Get4Momentum();
153 residualNucleus4Mom=G4LorentzVector(0.,0.,0.,targetNucleusMass);
155 }
156
157#ifdef debug_scatter
158 G4LorentzVector EpConservation=pNucleon+thePrimary.Get4Momentum()
159 - result.first - result.second;
160 if ( (EpConservation.vect().mag2() > .01*MeV*MeV )
161 || (std::abs(EpConservation.e()) > 0.1 * MeV ) )
162 {
163 G4cout << "Warning - G4QuasiElasticChannel::Scatter E-p non conservation : "
164 << EpConservation << G4endl;
165 }
166#endif
167
169 G4KineticTrack * sPrim=new G4KineticTrack(thePrimary.GetDefinition(),
170 0.,G4ThreeVector(0), scatteredHadron4Mom);
171 sPrim->SetCreatorModelID( secID );
172 ktv->push_back(sPrim);
173 if (result.first.e() > 0.)
174 {
175 G4KineticTrack * sNuc=new G4KineticTrack(pDef, 0.,G4ThreeVector(0), result.first);
176 sNuc->SetCreatorModelID( secID );
177 ktv->push_back(sNuc);
178 }
179 if(resZ || resA==1) // For the only neutron or for tnuclei with Z>0
180 {
181 G4KineticTrack * rNuc=new G4KineticTrack(resDef,
182 0.,G4ThreeVector(0), residualNucleus4Mom);
183 rNuc->SetCreatorModelID( secID );
184 ktv->push_back(rNuc);
185 }
186 else // The residual nucleus consists of only neutrons
187 {
188 residualNucleus4Mom/=resA; // Split 4-mom of A*n system equally
189 for(G4int in=0; in<resA; in++) // Loop over neutrons in A*n system.
190 {
191 G4KineticTrack* rNuc=new G4KineticTrack(resDef,
192 0.,G4ThreeVector(0), residualNucleus4Mom);
193 rNuc->SetCreatorModelID( secID );
194 ktv->push_back(rNuc);
195 }
196 }
197#ifdef debug_scatter
198 G4cout<<"G4QElC::Scat: Nucleon: "<<result.first <<" mass "<<result.first.mag() << G4endl;
199 G4cout<<"G4QElC::Scat: Project: "<<result.second<<" mass "<<result.second.mag()<< G4endl;
200#endif
201 return ktv;
202}
CLHEP::HepLorentzVector G4LorentzVector
static constexpr double MeV
Definition: G4SIunits.hh:200
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4UniformRand()
Definition: Randomize.hh:52
double mag2() const
Hep3Vector vect() const
G4LorentzVector Get4Momentum() const
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
void SetCreatorModelID(G4int id)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
G4double GetPDGCharge() const
const G4String & GetParticleName() const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
std::pair< G4LorentzVector, G4LorentzVector > Scatter(G4int NPDG, G4LorentzVector N4M, G4int pPDG, G4LorentzVector p4M)
virtual G4double GetMass()=0
virtual void Init(G4int theA, G4int theZ, G4int numberOfLambdas=0)=0
virtual const std::vector< G4Nucleon > & GetNucleons()=0
int G4lrint(double ad)
Definition: templates.hh:134
T sqr(const T &x)
Definition: templates.hh:128

References A, CLHEP::HepLorentzVector::e(), G4cout, G4endl, G4lrint(), G4UniformRand, G4DynamicParticle::Get4Momentum(), G4Nucleus::GetA_asInt(), G4DynamicParticle::GetDefinition(), G4IonTable::GetIon(), G4ParticleTable::GetIonTable(), G4V3DNucleus::GetMass(), G4V3DNucleus::GetNucleons(), G4ParticleDefinition::GetParticleName(), G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGEncoding(), G4ParticleDefinition::GetPDGMass(), G4Nucleus::GetZ_asInt(), G4V3DNucleus::Init(), CLHEP::Hep3Vector::mag2(), MeV, G4Neutron::Neutron(), G4QuasiElRatios::Scatter(), secID, G4KineticTrack::SetCreatorModelID(), CLHEP::HepLorentzVector::setE(), sqr(), the3DNucleus, theQuasiElastic, CLHEP::HepLorentzVector::vect(), and Z.

Referenced by G4TheoFSGenerator::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().

◆ 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::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::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

◆ epCheckLevels

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

◆ isBlocked

G4bool G4HadronicInteraction::isBlocked
protectedinherited

◆ recoilEnergyThreshold

G4double G4HadronicInteraction::recoilEnergyThreshold
privateinherited

◆ registry

G4HadronicInteractionRegistry* G4HadronicInteraction::registry
privateinherited

◆ secID

G4int G4QuasiElasticChannel::secID
private

Definition at line 71 of file G4QuasiElasticChannel.hh.

Referenced by G4QuasiElasticChannel(), and Scatter().

◆ the3DNucleus

G4V3DNucleus* G4QuasiElasticChannel::the3DNucleus
private

Definition at line 70 of file G4QuasiElasticChannel.hh.

Referenced by Scatter(), and ~G4QuasiElasticChannel().

◆ theBlockedList

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

◆ theBlockedListElements

std::vector<const G4Element *> G4HadronicInteraction::theBlockedListElements
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

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

◆ theQuasiElastic

G4QuasiElRatios* G4QuasiElasticChannel::theQuasiElastic
private

Definition at line 69 of file G4QuasiElasticChannel.hh.

Referenced by GetFraction(), Scatter(), and ~G4QuasiElasticChannel().

◆ 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: