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

#include <G4LMsdGenerator.hh>

Inheritance diagram for G4LMsdGenerator:
G4HadronicInteraction

Public Member Functions

void ActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Material *aMaterial)
 
G4HadFinalStateApplyYourself (const G4HadProjectile &thePrimary, G4Nucleus &theNucleus)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
void DeActivateFor (const G4Element *anElement)
 
void DeActivateFor (const G4Material *aMaterial)
 
 G4LMsdGenerator (const G4String &name="LMsdGenerator")
 
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 ()
 
G4bool IsApplicable (const G4HadProjectile &thePrimary, G4Nucleus &theNucleus)
 
G4bool IsBlocked (const G4Element *anElement) const
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
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)
 
G4double SampleMx (const G4HadProjectile *aParticle)
 
G4double SampleT (const G4HadProjectile *aParticle, G4double Mx)
 
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)
 
 ~G4LMsdGenerator ()
 

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

 G4LMsdGenerator (const G4LMsdGenerator &right)
 
int operator!= (const G4LMsdGenerator &right) const
 
const G4LMsdGeneratoroperator= (const G4LMsdGenerator &right)
 
int operator== (const G4LMsdGenerator &right) const
 

Private Attributes

std::pair< G4double, G4doubleepCheckLevels
 
G4int fPDGencoding
 
G4double recoilEnergyThreshold
 
G4HadronicInteractionRegistryregistry
 
G4int secID
 
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
 

Static Private Attributes

static const G4double fMxBdata [23][2]
 
static const G4double fProbMx [60][2]
 

Detailed Description

Definition at line 59 of file G4LMsdGenerator.hh.

Constructor & Destructor Documentation

◆ G4LMsdGenerator() [1/2]

G4LMsdGenerator::G4LMsdGenerator ( const G4String name = "LMsdGenerator")

Definition at line 46 of file G4LMsdGenerator.cc.

48
49{
50 fPDGencoding = 0;
51 secID = G4PhysicsModelCatalog::GetModelID( "model_LMsdGenerator" );
52
53 // theParticleChange = new G4HadFinalState;
54}
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static G4int GetModelID(const G4int modelIndex)
const char * name(G4int ptype)

References fPDGencoding, G4PhysicsModelCatalog::GetModelID(), and secID.

◆ ~G4LMsdGenerator()

G4LMsdGenerator::~G4LMsdGenerator ( )

Definition at line 56 of file G4LMsdGenerator.cc.

57{
58 // delete theParticleChange;
59}

◆ G4LMsdGenerator() [2/2]

G4LMsdGenerator::G4LMsdGenerator ( const G4LMsdGenerator right)
private

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 * G4LMsdGenerator::ApplyYourself ( const G4HadProjectile thePrimary,
G4Nucleus theNucleus 
)
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 114 of file G4LMsdGenerator.cc.

116{
118
119 const G4HadProjectile* aParticle = &aTrack;
120 G4double eTkin = aParticle->GetKineticEnergy();
121
122 if( eTkin <= 1.*CLHEP::GeV && aTrack.GetDefinition() != G4Proton::Proton())
123 {
125 theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
126 return &theParticleChange;
127 }
128
129 G4int A = targetNucleus.GetA_asInt();
130 G4int Z = targetNucleus.GetZ_asInt();
131
132 G4double plab = aParticle->GetTotalMomentum();
133 G4double plab2 = plab*plab;
134
135 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
136 G4double partMass = theParticle->GetPDGMass();
137
138 G4double oldE = partMass + eTkin;
139
141 G4double targMass2 = targMass*targMass;
142
143 G4LorentzVector partLV = aParticle->Get4Momentum();
144
145 G4double sumE = oldE + targMass;
146 G4double sumE2 = sumE*sumE;
147
148 G4ThreeVector p1 = partLV.vect();
149 // G4cout<<"p1 = "<<p1<<G4endl;
150 G4ParticleMomentum p1unit = p1.unit();
151
152 G4double Mx = SampleMx(aParticle); // in GeV
153 G4double t = SampleT( aParticle, Mx); // in GeV
154
155 Mx *= CLHEP::GeV;
156
157 G4double Mx2 = Mx*Mx;
158
159 // equation for q|| based on sum-E-P and new invariant mass
160
161 G4double B = sumE2 + targMass2 - Mx2 - plab2;
162
163 G4double a = 4*(plab2 - sumE2);
164 G4double b = 4*plab*B;
165 G4double c = B*B - 4*sumE2*targMass2;
166 G4double det2 = b*b - 4*a*c;
167 G4double qLong, det, eRetard; // , x2, x3, e2;
168
169 if( det2 >= 0.)
170 {
171 det = std::sqrt(det2);
172 qLong = (-b - det)/2./a;
173 eRetard = std::sqrt((plab-qLong)*(plab-qLong)+Mx2);
174 }
175 else
176 {
178 theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
179 return &theParticleChange;
180 }
182
183 plab -= qLong;
184
185 G4ThreeVector pRetard = plab*p1unit;
186
187 G4ThreeVector pTarg = p1 - pRetard;
188
189 G4double eTarg = std::sqrt( targMass2 + pTarg.mag2()); // std::sqrt( targMass*targMass + pTarg.mag2() );
190
191 G4LorentzVector lvRetard(pRetard, eRetard);
192 G4LorentzVector lvTarg(pTarg, eTarg);
193
194 lvTarg += lvRetard; // sum LV
195
196 G4ThreeVector bst = lvTarg.boostVector();
197
198 lvRetard.boost(-bst); // to CNS
199
200 G4ThreeVector pCMS = lvRetard.vect();
201 G4double momentumCMS = pCMS.mag();
202 G4double tMax = 4.0*momentumCMS*momentumCMS;
203
204 if( t > tMax ) t = tMax*G4UniformRand();
205
206 G4double cost = 1. - 2.0*t/tMax;
207
208
210 G4double sint;
211
212 if( cost > 1.0 || cost < -1.0 ) //
213 {
214 cost = 1.0;
215 sint = 0.0;
216 }
217 else // normal situation
218 {
219 sint = std::sqrt( (1.0-cost)*(1.0+cost) );
220 }
221 G4ThreeVector v1( sint*std::cos(phi), sint*std::sin(phi), cost);
222
223 v1 *= momentumCMS;
224
225 G4LorentzVector lvRes( v1.x(),v1.y(),v1.z(), std::sqrt( momentumCMS*momentumCMS + Mx2));
226
227 lvRes.boost(bst); // to LS
228
229 lvTarg -= lvRes;
230
231 G4double eRecoil = lvTarg.e() - targMass;
232
233 if( eRecoil > 100.*CLHEP::MeV ) // add recoil nucleus
234 {
235 G4ParticleDefinition * recoilDef = 0;
236
237 if ( Z == 1 && A == 1 ) { recoilDef = G4Proton::Proton(); }
238 else if ( Z == 1 && A == 2 ) { recoilDef = G4Deuteron::Deuteron(); }
239 else if ( Z == 1 && A == 3 ) { recoilDef = G4Triton::Triton(); }
240 else if ( Z == 2 && A == 3 ) { recoilDef = G4He3::He3(); }
241 else if ( Z == 2 && A == 4 ) { recoilDef = G4Alpha::Alpha(); }
242 else
243 {
244 recoilDef =
246 }
247 G4DynamicParticle * aSec = new G4DynamicParticle( recoilDef, lvTarg);
249 }
250 else if( eRecoil > 0.0 )
251 {
253 }
254
256 FindParticle(fPDGencoding);
257
258 // G4cout<<fPDGencoding<<", "<<ddPart->GetParticleName()<<", "<<ddPart->GetPDGMass()<<" MeV; lvRes = "<<lvRes<<G4endl;
259
260 // G4DynamicParticle * aRes = new G4DynamicParticle( ddPart, lvRes);
261 // theParticleChange.AddSecondary(aRes); // simply return resonance
262
263
264
265 // Recursive decay using methods of G4KineticTrack
266
267 G4KineticTrack ddkt( ddPart, 0., G4ThreeVector(0.,0.,0.), lvRes);
268 G4KineticTrackVector* ddktv = ddkt.Decay();
269
271
272 for( unsigned int i = 0; i < ddktv->size(); i++ ) // add products to partchange
273 {
274 G4DynamicParticle * aNew =
275 new G4DynamicParticle( ddktv->operator[](i)->GetDefinition(),
276 ddktv->operator[](i)->Get4Momentum());
277
278 // G4cout<<" "<<i<<", "<<aNew->GetDefinition()->GetParticleName()<<", "<<aNew->Get4Momentum()<<G4endl;
279
281 delete ddktv->operator[](i);
282 }
283 delete ddktv;
284
285 return &theParticleChange;
286}
G4double B(G4double temperature)
@ stopAndKill
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
Hep3Vector unit() const
double mag2() const
double mag() const
Hep3Vector vect() const
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:93
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
void SetLocalEnergyDeposit(G4double aE)
G4double GetTotalMomentum() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
static G4He3 * He3()
Definition: G4He3.cc:93
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
G4double SampleMx(const G4HadProjectile *aParticle)
G4double SampleT(const G4HadProjectile *aParticle, G4double Mx)
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
static G4Proton * Proton()
Definition: G4Proton.cc:92
static G4Triton * Triton()
Definition: G4Triton.cc:93
static constexpr double GeV
static constexpr double twopi
Definition: SystemOfUnits.h:56
static constexpr double MeV
ParticleList decay(Cluster *const c)
Carries out a cluster decay.

References A, G4HadFinalState::AddSecondary(), G4Alpha::Alpha(), B(), CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::boostVector(), G4HadFinalState::Clear(), G4KineticTrack::Decay(), G4INCL::ClusterDecay::decay(), G4Deuteron::Deuteron(), CLHEP::HepLorentzVector::e(), fPDGencoding, G4UniformRand, G4HadProjectile::Get4Momentum(), G4Nucleus::GetA_asInt(), G4HadProjectile::GetDefinition(), G4IonTable::GetIon(), G4ParticleTable::GetIonTable(), G4HadProjectile::GetKineticEnergy(), G4NucleiProperties::GetNuclearMass(), G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetPDGMass(), G4HadProjectile::GetTotalMomentum(), G4Nucleus::GetZ_asInt(), CLHEP::GeV, G4He3::He3(), CLHEP::Hep3Vector::mag(), CLHEP::Hep3Vector::mag2(), CLHEP::MeV, G4Proton::Proton(), SampleMx(), SampleT(), secID, G4HadFinalState::SetEnergyChange(), G4HadFinalState::SetLocalEnergyDeposit(), G4HadFinalState::SetMomentumChange(), G4HadFinalState::SetStatusChange(), stopAndKill, G4HadronicInteraction::theParticleChange, G4Triton::Triton(), CLHEP::twopi, CLHEP::Hep3Vector::unit(), CLHEP::HepLorentzVector::vect(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), CLHEP::Hep3Vector::z(), and Z.

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

◆ 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(), 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 G4LMsdGenerator::IsApplicable ( const G4HadProjectile thePrimary,
G4Nucleus theNucleus 
)
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 78 of file G4LMsdGenerator.cc.

80{
81 G4bool applied = false;
82
83 if( ( aTrack.GetDefinition() == G4Proton::Proton() ||
84 aTrack.GetDefinition() == G4Neutron::Neutron() ) &&
85 targetNucleus.GetA_asInt() >= 1 &&
86 // aTrack.GetKineticEnergy() > 1800*CLHEP::MeV
87 aTrack.GetKineticEnergy() > 300*CLHEP::MeV
88 ) // 750*CLHEP::MeV )
89 {
90 applied = true;
91 }
92 else if( ( aTrack.GetDefinition() == G4PionPlus::PionPlus() ||
93 aTrack.GetDefinition() == G4PionMinus::PionMinus() ) &&
94 targetNucleus.GetA_asInt() >= 1 &&
95 aTrack.GetKineticEnergy() > 2340*CLHEP::MeV )
96 {
97 applied = true;
98 }
99 else if( ( aTrack.GetDefinition() == G4KaonPlus::KaonPlus() ||
100 aTrack.GetDefinition() == G4KaonMinus::KaonMinus() ) &&
101 targetNucleus.GetA_asInt() >= 1 &&
102 aTrack.GetKineticEnergy() > 1980*CLHEP::MeV )
103 {
104 applied = true;
105 }
106 return applied;
107}
bool G4bool
Definition: G4Types.hh:86
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:112
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:112
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:97
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:97

References G4Nucleus::GetA_asInt(), G4HadProjectile::GetDefinition(), G4HadProjectile::GetKineticEnergy(), G4KaonMinus::KaonMinus(), G4KaonPlus::KaonPlus(), CLHEP::MeV, G4Neutron::Neutron(), G4PionMinus::PionMinus(), G4PionPlus::PionPlus(), and G4Proton::Proton().

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

Reimplemented from G4HadronicInteraction.

Definition at line 61 of file G4LMsdGenerator.cc.

62{
63 outFile << GetModelName() <<" consists of a "
64 << " string model and a stage to de-excite the excited nuclear fragment."
65 << "\n<p>"
66 << "The string model simulates the interaction of\n"
67 << "an incident hadron with a nucleus, forming \n"
68 << "excited strings, decays these strings into hadrons,\n"
69 << "and leaves an excited nucleus. \n"
70 << "<p>The string model:\n";
71}
const G4String & GetModelName() const

References G4HadronicInteraction::GetModelName().

◆ operator!=() [1/2]

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

◆ operator!=() [2/2]

int G4LMsdGenerator::operator!= ( const G4LMsdGenerator right) const
private

◆ operator=()

const G4LMsdGenerator & G4LMsdGenerator::operator= ( const G4LMsdGenerator right)
private

◆ operator==() [1/2]

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

◆ operator==() [2/2]

int G4LMsdGenerator::operator== ( const G4LMsdGenerator right) const
private

◆ SampleInvariantT()

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

◆ SampleMx()

G4double G4LMsdGenerator::SampleMx ( const G4HadProjectile aParticle)

Definition at line 292 of file G4LMsdGenerator.cc.

293{
294 G4double Mx = 0.;
295 G4int i;
296 G4double rand = G4UniformRand();
297
298 for( i = 0; i < 60; i++)
299 {
300 if( rand >= fProbMx[i][1] ) break;
301 }
302 if(i <= 0) Mx = fProbMx[0][0];
303 else if(i >= 59) Mx = fProbMx[59][0];
304 else Mx = fProbMx[i][0];
305
306 fPDGencoding = 0;
307
308 if ( Mx <= 1.45 )
309 {
310 if( aParticle->GetDefinition() == G4Proton::Proton() )
311 {
312 Mx = 1.44;
313 // fPDGencoding = 12212;
314 fPDGencoding = 2214;
315 }
316 else if( aParticle->GetDefinition() == G4Neutron::Neutron() )
317 {
318 Mx = 1.44;
319 fPDGencoding = 12112;
320 }
321 else if( aParticle->GetDefinition() == G4PionPlus::PionPlus() )
322 {
323 // Mx = 1.3;
324 // fPDGencoding = 100211;
325 Mx = 1.26;
326 fPDGencoding = 20213; // a1(1260)+
327 }
328 else if( aParticle->GetDefinition() == G4PionMinus::PionMinus() )
329 {
330 // Mx = 1.3;
331 // fPDGencoding = -100211;
332 Mx = 1.26;
333 fPDGencoding = -20213; // a1(1260)-
334 }
335 else if( aParticle->GetDefinition() == G4KaonPlus::KaonPlus() )
336 {
337 Mx = 1.27;
338 fPDGencoding = 10323;
339 }
340 else if( aParticle->GetDefinition() == G4KaonMinus::KaonMinus() )
341 {
342 Mx = 1.27;
343 fPDGencoding = -10323;
344 }
345 }
346 else if ( Mx <= 1.55 )
347 {
348 if( aParticle->GetDefinition() == G4Proton::Proton() )
349 {
350 Mx = 1.52;
351 // fPDGencoding = 2124;
352 fPDGencoding = 2214;
353 }
354 else if( aParticle->GetDefinition() == G4Neutron::Neutron() )
355 {
356 Mx = 1.52;
357 fPDGencoding = 1214;
358 }
359 else if( aParticle->GetDefinition() == G4PionPlus::PionPlus() )
360 {
361 // Mx = 1.45;
362 // fPDGencoding = 10211;
363 Mx = 1.32;
364 fPDGencoding = 215; // a2(1320)+
365 }
366 else if( aParticle->GetDefinition() == G4PionMinus::PionMinus() )
367 {
368 // Mx = 1.45;
369 // fPDGencoding = -10211;
370 Mx = 1.32;
371 fPDGencoding = -215; // a2(1320)-
372 }
373 else if( aParticle->GetDefinition() == G4KaonPlus::KaonPlus() )
374 {
375 Mx = 1.46;
376 fPDGencoding = 100321;
377 }
378 else if( aParticle->GetDefinition() == G4KaonMinus::KaonMinus() )
379 {
380 Mx = 1.46;
381 fPDGencoding = -100321;
382 }
383 }
384 else
385 {
386 if( aParticle->GetDefinition() == G4Proton::Proton() )
387 {
388 Mx = 1.68;
389 // fPDGencoding = 12216;
390 fPDGencoding = 2214;
391 }
392 else if( aParticle->GetDefinition() == G4Neutron::Neutron() )
393 {
394 Mx = 1.68;
395 fPDGencoding = 12116;
396 }
397 else if( aParticle->GetDefinition() == G4PionPlus::PionPlus() )
398 {
399 Mx = 1.67;
400 fPDGencoding = 10215; // pi2(1670)+
401 // Mx = 1.45;
402 // fPDGencoding = 10211;
403 }
404 else if( aParticle->GetDefinition() == G4PionMinus::PionMinus() )
405 {
406 Mx = 1.67; // f0 problems->4pi vmg 20.11.14
407 fPDGencoding = -10215; // pi2(1670)-
408 // Mx = 1.45;
409 // fPDGencoding = -10211;
410 }
411 else if( aParticle->GetDefinition() == G4KaonPlus::KaonPlus() )
412 {
413 Mx = 1.68;
414 fPDGencoding = 30323;
415 }
416 else if( aParticle->GetDefinition() == G4KaonMinus::KaonMinus() )
417 {
418 Mx = 1.68;
419 fPDGencoding = -30323;
420 }
421 }
422 if(fPDGencoding == 0)
423 {
424 Mx = 1.44;
425 // fPDGencoding = 12212;
426 fPDGencoding = 2214;
427 }
428 G4ParticleDefinition* myResonance =
430
431 if ( myResonance ) Mx = myResonance->GetPDGMass();
432
433 // G4cout<<"PDG-ID = "<<fPDGencoding<<"; with mass = "<<Mx/CLHEP::GeV<<" GeV"<<G4endl;
434
435 return Mx/CLHEP::GeV;
436}
static const G4double fProbMx[60][2]
G4ParticleDefinition * FindParticle(G4int PDGEncoding)

References G4ParticleTable::FindParticle(), fPDGencoding, fProbMx, G4UniformRand, G4HadProjectile::GetDefinition(), G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetPDGMass(), CLHEP::GeV, G4KaonMinus::KaonMinus(), G4KaonPlus::KaonPlus(), G4Neutron::Neutron(), G4PionMinus::PionMinus(), G4PionPlus::PionPlus(), and G4Proton::Proton().

Referenced by ApplyYourself().

◆ SampleT()

G4double G4LMsdGenerator::SampleT ( const G4HadProjectile aParticle,
G4double  Mx 
)

Definition at line 442 of file G4LMsdGenerator.cc.

444{
445 G4double t=0., b=0., rTkin = 50.*CLHEP::GeV, eTkin = aParticle->GetKineticEnergy();
446 G4int i;
447
448 for( i = 0; i < 23; ++i)
449 {
450 if( Mx <= fMxBdata[i][0] ) break;
451 }
452 if( i <= 0 ) b = fMxBdata[0][1];
453 else if( i >= 22 ) b = fMxBdata[22][1];
454 else b = fMxBdata[i][1];
455
456 if( eTkin > rTkin ) b *= 1. + G4Log(eTkin/rTkin);
457
458 G4double rand = G4UniformRand();
459
460 t = -G4Log(rand)/b;
461
462 t *= (CLHEP::GeV*CLHEP::GeV); // in G4 internal units
463
464 return t;
465}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
static const G4double fMxBdata[23][2]

References fMxBdata, G4Log(), G4UniformRand, G4HadProjectile::GetKineticEnergy(), and CLHEP::GeV.

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

◆ 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

◆ fMxBdata

const G4double G4LMsdGenerator::fMxBdata
staticprivate
Initial value:
=
{
{1.09014, 17.8620},
{1.12590, 19.2831},
{1.18549, 17.6907},
{1.21693, 16.4760},
{1.25194, 15.3867},
{1.26932, 14.4236},
{1.29019, 13.2931},
{1.30755, 12.2882},
{1.31790, 11.4509},
{1.33888, 10.6969},
{1.34911, 9.44130},
{1.37711, 8.56148},
{1.39101, 7.76593},
{1.42608, 6.88582},
{1.48593, 6.13019},
{1.53179, 5.87723},
{1.58111, 5.37308},
{1.64105, 4.95217},
{1.69037, 4.44803},
{1.81742, 3.89879},
{1.88096, 3.68693},
{1.95509, 3.43278},
{2.02219, 3.30445}
}

Definition at line 95 of file G4LMsdGenerator.hh.

Referenced by SampleT().

◆ fPDGencoding

G4int G4LMsdGenerator::fPDGencoding
private

Definition at line 92 of file G4LMsdGenerator.hh.

Referenced by ApplyYourself(), G4LMsdGenerator(), and SampleMx().

◆ fProbMx

const G4double G4LMsdGenerator::fProbMx
staticprivate

Definition at line 96 of file G4LMsdGenerator.hh.

Referenced by SampleMx().

◆ isBlocked

G4bool G4HadronicInteraction::isBlocked
protectedinherited

◆ recoilEnergyThreshold

G4double G4HadronicInteraction::recoilEnergyThreshold
privateinherited

◆ registry

G4HadronicInteractionRegistry* G4HadronicInteraction::registry
privateinherited

◆ secID

G4int G4LMsdGenerator::secID
private

Definition at line 93 of file G4LMsdGenerator.hh.

Referenced by ApplyYourself(), and G4LMsdGenerator().

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

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