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

#include <G4PreCompoundModel.hh>

Inheritance diagram for G4PreCompoundModel:
G4VPreCompoundModel G4HadronicInteraction

Public Member Functions

void ActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Material *aMaterial)
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &thePrimary, G4Nucleus &theNucleus) final
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) final
 
void DeActivateFor (const G4Element *anElement)
 
void DeActivateFor (const G4Material *aMaterial)
 
virtual G4ReactionProductVectorDeExcite (G4Fragment &aFragment) final
 
virtual void DeExciteModelDescription (std::ostream &outFile) const final
 
 G4PreCompoundModel (G4ExcitationHandler *ptr=nullptr)
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
G4ExcitationHandlerGetExcitationHandler () 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 () final
 
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 final
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 
G4bool operator!= (const G4VPreCompoundModel &right) const =delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4bool operator== (const G4VPreCompoundModel &right) const =delete
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
void SetExcitationHandler (G4ExcitationHandler *ptr)
 
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)
 
virtual ~G4PreCompoundModel ()
 

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

 G4PreCompoundModel (const G4PreCompoundModel &)=delete
 
G4bool operator!= (const G4PreCompoundModel &right) const =delete
 
const G4PreCompoundModeloperator= (const G4PreCompoundModel &right)=delete
 
G4bool operator== (const G4PreCompoundModel &right) const =delete
 
void PerformEquilibriumEmission (const G4Fragment &aFragment, G4ReactionProductVector *result) const
 

Private Attributes

std::pair< G4double, G4doubleepCheckLevels
 
G4double fHighLimitExc = DBL_MAX
 
G4double fLowLimitExc = 0.0
 
G4NuclearLevelDatafNuclData = nullptr
 
G4bool isActive = true
 
G4bool isInitialised = false
 
G4int minA = 5
 
G4int minZ = 3
 
G4int modelID = -1
 
const G4ParticleDefinitionneutron
 
const G4ParticleDefinitionproton
 
G4double recoilEnergyThreshold
 
G4HadronicInteractionRegistryregistry
 
std::vector< const G4Material * > theBlockedList
 
std::vector< const G4Element * > theBlockedListElements
 
G4PreCompoundEmissiontheEmission = nullptr
 
G4ExcitationHandlertheExcitationHandler
 
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
 
G4HadFinalState theResult
 
G4VPreCompoundTransitionstheTransition = nullptr
 
G4bool useSCO = false
 

Detailed Description

Definition at line 62 of file G4PreCompoundModel.hh.

Constructor & Destructor Documentation

◆ G4PreCompoundModel() [1/2]

G4PreCompoundModel::G4PreCompoundModel ( G4ExcitationHandler ptr = nullptr)
explicit

Definition at line 70 of file G4PreCompoundModel.cc.

71 : G4VPreCompoundModel(ptr,"PRECO")
72{
73 //G4cout << "### NEW PrecompoundModel " << this << G4endl;
74 if(nullptr == ptr) { SetExcitationHandler(new G4ExcitationHandler()); }
75
80}
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4NuclearLevelData * GetInstance()
static G4int GetModelID(const G4int modelIndex)
const G4ParticleDefinition * neutron
G4NuclearLevelData * fNuclData
const G4ParticleDefinition * proton
static G4Proton * Proton()
Definition: G4Proton.cc:92
G4VPreCompoundModel(G4ExcitationHandler *ptr=nullptr, const G4String &modelName="PrecompoundModel")
void SetExcitationHandler(G4ExcitationHandler *ptr)

References fNuclData, G4NuclearLevelData::GetInstance(), G4PhysicsModelCatalog::GetModelID(), modelID, G4Neutron::Neutron(), neutron, G4Proton::Proton(), proton, and G4VPreCompoundModel::SetExcitationHandler().

◆ ~G4PreCompoundModel()

G4PreCompoundModel::~G4PreCompoundModel ( )
virtual

Definition at line 84 of file G4PreCompoundModel.cc.

85{
86 delete theEmission;
87 delete theTransition;
88 delete GetExcitationHandler();
90}
G4HadFinalState theResult
G4PreCompoundEmission * theEmission
G4VPreCompoundTransitions * theTransition
G4ExcitationHandler * GetExcitationHandler() const

References G4HadFinalState::Clear(), G4VPreCompoundModel::GetExcitationHandler(), theEmission, theResult, and theTransition.

◆ G4PreCompoundModel() [2/2]

G4PreCompoundModel::G4PreCompoundModel ( const G4PreCompoundModel )
privatedelete

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

Reimplemented from G4HadronicInteraction.

Definition at line 135 of file G4PreCompoundModel.cc.

137{
138 const G4ParticleDefinition* primary = thePrimary.GetDefinition();
139 if(primary != neutron && primary != proton) {
141 ed << "G4PreCompoundModel is used for ";
142 if(primary) { ed << primary->GetParticleName(); }
143 G4Exception("G4PreCompoundModel::ApplyYourself()","had0033",FatalException,
144 ed,"");
145 return nullptr;
146 }
147
148 G4int Zp = 0;
149 G4int Ap = 1;
150 if(primary == proton) { Zp = 1; }
151
152 G4double timePrimary=thePrimary.GetGlobalTime();
153
154 G4int A = theNucleus.GetA_asInt();
155 G4int Z = theNucleus.GetZ_asInt();
156
157 //G4cout << "### G4PreCompoundModel::ApplyYourself: A= " << A << " Z= " << Z
158 // << " Ap= " << Ap << " Zp= " << Zp << G4endl;
159 // 4-Momentum
160 G4LorentzVector p = thePrimary.Get4Momentum();
162 p += G4LorentzVector(0.0,0.0,0.0,mass);
163 //G4cout << "Primary 4-mom " << p << " mass= " << mass << G4endl;
164
165 // prepare fragment
166 G4Fragment anInitialState(A + Ap, Z + Zp, p);
167 anInitialState.SetNumberOfExcitedParticle(2, 1);
168 anInitialState.SetNumberOfHoles(1,0);
169 anInitialState.SetCreationTime(thePrimary.GetGlobalTime());
170 anInitialState.SetCreatorModelID(modelID);
171
172 // call excitation handler
173 G4ReactionProductVector* result = DeExcite(anInitialState);
174
175 // fill particle change
178 for(auto const & prod : *result) {
179 G4DynamicParticle * aNewDP = new G4DynamicParticle(prod->GetDefinition(),
180 prod->GetTotalEnergy(),
181 prod->GetMomentum());
182 G4HadSecondary aNew = G4HadSecondary(aNewDP);
183 G4double time = std::max(prod->GetFormationTime(), 0.0);
184 aNew.SetTime(timePrimary + time);
185 aNew.SetCreatorModelID(prod->GetCreatorModelID());
186 delete prod;
188 }
189 delete result;
190
191 //return the filled particle change
192 return &theResult;
193}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
@ stopAndKill
CLHEP::HepLorentzVector G4LorentzVector
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & Get4Momentum() const
G4double GetGlobalTime() const
void SetTime(G4double aT)
void SetCreatorModelID(G4int id)
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetA_asInt() const
Definition: G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105
const G4String & GetParticleName() const
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment) final
T max(const T t1, const T t2)
brief Return the largest of the two arguments

References A, G4HadFinalState::AddSecondary(), G4HadFinalState::Clear(), DeExcite(), FatalException, G4Exception(), G4HadProjectile::Get4Momentum(), G4Nucleus::GetA_asInt(), G4HadProjectile::GetDefinition(), G4HadProjectile::GetGlobalTime(), G4NucleiProperties::GetNuclearMass(), G4ParticleDefinition::GetParticleName(), G4Nucleus::GetZ_asInt(), G4INCL::Math::max(), modelID, neutron, proton, G4Fragment::SetCreationTime(), G4HadSecondary::SetCreatorModelID(), G4Fragment::SetCreatorModelID(), G4Fragment::SetNumberOfExcitedParticle(), G4Fragment::SetNumberOfHoles(), G4HadFinalState::SetStatusChange(), G4HadSecondary::SetTime(), stopAndKill, theResult, and Z.

◆ Block()

void G4HadronicInteraction::Block ( )
inlineprotectedinherited

◆ BuildPhysicsTable()

void G4PreCompoundModel::BuildPhysicsTable ( const G4ParticleDefinition )
finalvirtual

Reimplemented from G4HadronicInteraction.

Definition at line 94 of file G4PreCompoundModel.cc.

95{
97}
virtual void InitialiseModel() final

References InitialiseModel().

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

◆ DeExcite()

G4ReactionProductVector * G4PreCompoundModel::DeExcite ( G4Fragment aFragment)
finalvirtual

Implements G4VPreCompoundModel.

Definition at line 197 of file G4PreCompoundModel.cc.

198{
200
202 G4double U = aFragment.GetExcitationEnergy();
203 G4int Z = aFragment.GetZ_asInt();
204 G4int A = aFragment.GetA_asInt();
205
206 //G4cout << "### G4PreCompoundModel::DeExcite" << G4endl;
207 //G4cout << aFragment << G4endl;
208
209 // Perform Equilibrium Emission
210 if (!isActive || (Z < minZ && A < minA) ||
211 U < fLowLimitExc*A || U > A*fHighLimitExc) {
212 PerformEquilibriumEmission(aFragment, Result);
213 return Result;
214 }
215
216 // main loop
217 G4int count = 0;
218 const G4double ldfact = 12.0/CLHEP::pi2;
219 const G4int countmax = 1000;
220 for (;;) {
221 //G4cout << "### PreCompound loop over fragment" << G4endl;
222 //G4cout << aFragment << G4endl;
223 U = aFragment.GetExcitationEnergy();
224 Z = aFragment.GetZ_asInt();
225 A = aFragment.GetA_asInt();
226 G4int eqExcitonNumber =
227 G4lrint(std::sqrt(ldfact*U*fNuclData->GetLevelDensity(Z, A, U)));
228 //
229 // G4cout<<"Neq="<<EquilibriumExcitonNumber<<G4endl;
230 //
231 // J. M. Quesada (Jan. 08) equilibrium hole number could be used as preeq.
232 // evap. delimiter (IAEA report)
233
234 // Loop for transitions, it is performed while there are
235 // preequilibrium transitions.
236 G4bool isTransition = false;
237
238 // G4cout<<"----------------------------------------"<<G4endl;
239 // G4double NP=aFragment.GetNumberOfParticles();
240 // G4double NH=aFragment.GetNumberOfHoles();
241 // G4double NE=aFragment.GetNumberOfExcitons();
242 // G4cout<<" Ex. Energy="<<aFragment.GetExcitationEnergy()<<G4endl;
243 // G4cout<<"N. excitons="<<NE<<" N. Part="<<NP<<"N. Holes ="<<NH<<G4endl;
244 do {
245 ++count;
246 //G4cout<<"transition number .."<<count
247 // <<" n ="<<aFragment.GetNumberOfExcitons()<<G4endl;
248 // soft cutoff criterium as an "ad-hoc" solution to force
249 // increase in evaporation
250 G4int ne = aFragment.GetNumberOfExcitons();
251 G4bool go_ahead = (ne <= eqExcitonNumber);
252
253 //J. M. Quesada (Apr. 08): soft-cutoff switched off by default
254 if (useSCO && go_ahead) {
255 G4double x = (G4double)(ne - eqExcitonNumber)/(G4double)eqExcitonNumber;
256 if( G4UniformRand() < 1.0 - G4Exp(-x*x/0.32) ) { go_ahead = false; }
257 }
258
259 // JMQ: WARNING: CalculateProbability MUST be called prior to Get!!
260 // (O values would be returned otherwise)
261 G4double transProbability =
266 //G4cout<<"#0 P1="<<P1<<" P2="<<P2<<" P3="<<P3<<G4endl;
267
268 //J.M. Quesada (May 2008) Physical criterium (lamdas) PREVAILS over
269 // approximation (critical exciton number)
270 //V.Ivanchenko (May 2011) added check on number of nucleons
271 // to send a fragment to FermiBreakUp
272 // or check on limits of excitation
273 if(!go_ahead || P1 <= P2+P3 || Z < minZ || A < minA ||
274 U <= fLowLimitExc*A || U > A*fHighLimitExc ||
275 aFragment.GetNumberOfExcitons() <= 0) {
276 //G4cout<<"#4 EquilibriumEmission"<<G4endl;
277 PerformEquilibriumEmission(aFragment,Result);
278 return Result;
279 }
280 G4double emissionProbability =
282 //G4cout<<"#1 TotalEmissionProbability="<<TotalEmissionProbability
283 // <<" Nex= " <<aFragment.GetNumberOfExcitons()<<G4endl;
284 //J.M.Quesada (May 08) this has already been done in order to decide
285 // what to do (preeq-eq)
286 // Sum of all probabilities
287 G4double TotalProbability = emissionProbability + transProbability;
288
289 // Select subprocess
290 if (TotalProbability*G4UniformRand() > emissionProbability) {
291 //G4cout<<"#2 Transition"<<G4endl;
292 // It will be transition to state with a new number of excitons
293 isTransition = true;
294 // Perform the transition
296 } else {
297 //G4cout<<"#3 Emission"<<G4endl;
298 // It will be fragment emission
299 isTransition = false;
300 Result->push_back(theEmission->PerformEmission(aFragment));
301 }
302 // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
303 } while (isTransition); // end of do loop
304
305 // stop if too many iterations
306 if(count >= countmax) {
308 ed << "G4PreCompoundModel loop over " << countmax << " iterations; "
309 << "current G4Fragment: \n" << aFragment;
310 G4Exception("G4PreCompoundModel::DeExcite()","had0034",JustWarning,
311 ed,"");
312 PerformEquilibriumEmission(aFragment, Result);
313 return Result;
314 }
315 } // end of for (;;) loop
316 return Result;
317}
static const G4double * P1[nN]
static const G4double * P2[nN]
@ JustWarning
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
bool G4bool
Definition: G4Types.hh:86
#define G4UniformRand()
Definition: Randomize.hh:52
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:299
G4int GetZ_asInt() const
Definition: G4Fragment.hh:276
G4int GetNumberOfExcitons() const
Definition: G4Fragment.hh:356
G4int GetA_asInt() const
Definition: G4Fragment.hh:271
G4double GetLevelDensity(G4int Z, G4int A, G4double U)
G4ReactionProduct * PerformEmission(G4Fragment &aFragment)
G4double GetTotalProbability(const G4Fragment &aFragment)
void PerformEquilibriumEmission(const G4Fragment &aFragment, G4ReactionProductVector *result) const
virtual G4double CalculateProbability(const G4Fragment &aFragment)=0
virtual void PerformTransition(G4Fragment &aFragment)=0
static constexpr double pi2
Definition: SystemOfUnits.h:58
int G4lrint(double ad)
Definition: templates.hh:134

References A, G4VPreCompoundTransitions::CalculateProbability(), fHighLimitExc, fNuclData, G4Exception(), G4Exp(), G4lrint(), G4UniformRand, G4Fragment::GetA_asInt(), G4Fragment::GetExcitationEnergy(), G4NuclearLevelData::GetLevelDensity(), G4Fragment::GetNumberOfExcitons(), G4PreCompoundEmission::GetTotalProbability(), G4VPreCompoundTransitions::GetTransitionProb1(), G4VPreCompoundTransitions::GetTransitionProb2(), G4VPreCompoundTransitions::GetTransitionProb3(), G4Fragment::GetZ_asInt(), InitialiseModel(), isActive, isInitialised, JustWarning, minA, minZ, P1, P2, G4PreCompoundEmission::PerformEmission(), PerformEquilibriumEmission(), G4VPreCompoundTransitions::PerformTransition(), CLHEP::pi2, theEmission, theTransition, useSCO, and Z.

Referenced by G4LowEGammaNuclearModel::ApplyYourself(), ApplyYourself(), G4LowEIonFragmentation::ApplyYourself(), and G4NeutrinoNucleusModel::RecoilDeexcitation().

◆ DeExciteModelDescription()

void G4PreCompoundModel::DeExciteModelDescription ( std::ostream &  outFile) const
finalvirtual

Implements G4VPreCompoundModel.

Definition at line 349 of file G4PreCompoundModel.cc.

350{
351 outFile << "description of precompound model as used with DeExcite()" << "\n";
352}

◆ GetEnergyMomentumCheckLevels()

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

◆ GetExcitationHandler()

G4ExcitationHandler * G4VPreCompoundModel::GetExcitationHandler ( ) const
inlineinherited

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

Reimplemented from G4HadronicInteraction.

Definition at line 101 of file G4PreCompoundModel.cc.

102{
103 if(isInitialised) { return; }
104 isInitialised = true;
105
106 //G4cout << "G4PreCompoundModel::InitialiseModel() started" << G4endl;
107
109
112
113 useSCO = param->UseSoftCutoff();
114
115 minZ = param->GetMinZForPreco();
116 minA = param->GetMinAForPreco();
117
119 if(param->UseHETC()) { theEmission->SetHETCModel(); }
121
122 if(param->UseGNASH()) { theTransition = new G4GNASHTransitions; }
125 theTransition->UseCEMtr(param->UseCEM());
126
127 if(param->PrecoDummy()) { isActive = false; }
128
130}
G4double GetPrecoLowEnergy() const
G4double GetPrecoHighEnergy() const
G4DeexPrecoParameters * GetParameters()

References fHighLimitExc, fLowLimitExc, fNuclData, G4VPreCompoundModel::GetExcitationHandler(), G4DeexPrecoParameters::GetMinAForPreco(), G4DeexPrecoParameters::GetMinZForPreco(), G4NuclearLevelData::GetParameters(), G4DeexPrecoParameters::GetPrecoHighEnergy(), G4DeexPrecoParameters::GetPrecoLowEnergy(), G4DeexPrecoParameters::GetPrecoModelType(), G4ExcitationHandler::Initialise(), isActive, isInitialised, minA, minZ, G4DeexPrecoParameters::NeverGoBack(), G4DeexPrecoParameters::PrecoDummy(), G4PreCompoundEmission::SetHETCModel(), G4PreCompoundEmission::SetOPTxs(), theEmission, theTransition, G4DeexPrecoParameters::UseCEM(), G4VPreCompoundTransitions::UseCEMtr(), G4DeexPrecoParameters::UseGNASH(), G4DeexPrecoParameters::UseHETC(), G4VPreCompoundTransitions::UseNGB(), useSCO, and G4DeexPrecoParameters::UseSoftCutoff().

Referenced by BuildPhysicsTable(), and DeExcite().

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

Reimplemented from G4HadronicInteraction.

Definition at line 323 of file G4PreCompoundModel.cc.

324{
325 outFile
326 << "The GEANT4 precompound model is considered as an extension of the\n"
327 << "hadron kinetic model. It gives a possibility to extend the low energy range\n"
328 << "of the hadron kinetic model for nucleon-nucleus inelastic collision and it \n"
329 << "provides a ”smooth” transition from kinetic stage of reaction described by the\n"
330 << "hadron kinetic model to the equilibrium stage of reaction described by the\n"
331 << "equilibrium deexcitation models.\n"
332 << "The initial information for calculation of pre-compound nuclear stage\n"
333 << "consists of the atomic mass number A, charge Z of residual nucleus, its\n"
334 << "four momentum P0 , excitation energy U and number of excitons n, which equals\n"
335 << "the sum of the number of particles p (from them p_Z are charged) and the number of\n"
336 << "holes h.\n"
337 << "At the preequilibrium stage of reaction, we follow the exciton model approach in ref. [1],\n"
338 << "taking into account the competition among all possible nuclear transitions\n"
339 << "with ∆n = +2, −2, 0 (which are defined by their associated transition probabilities) and\n"
340 << "the emission of neutrons, protons, deuterons, thritium and helium nuclei (also defined by\n"
341 << "their associated emission probabilities according to exciton model)\n"
342 << "\n"
343 << "[1] K.K. Gudima, S.G. Mashnik, V.D. Toneev, Nucl. Phys. A401 329 (1983)\n"
344 << "\n";
345}

◆ operator!=() [1/3]

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

◆ operator!=() [2/3]

G4bool G4PreCompoundModel::operator!= ( const G4PreCompoundModel right) const
privatedelete

◆ operator!=() [3/3]

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

◆ operator=()

const G4PreCompoundModel & G4PreCompoundModel::operator= ( const G4PreCompoundModel right)
privatedelete

◆ operator==() [1/3]

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

◆ operator==() [2/3]

G4bool G4PreCompoundModel::operator== ( const G4PreCompoundModel right) const
privatedelete

◆ operator==() [3/3]

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

◆ PerformEquilibriumEmission()

void G4PreCompoundModel::PerformEquilibriumEmission ( const G4Fragment aFragment,
G4ReactionProductVector result 
) const
inlineprivate

Definition at line 114 of file G4PreCompoundModel.hh.

117{
118 G4ReactionProductVector* deexResult =
119 GetExcitationHandler()->BreakItUp(aFragment);
120 result->insert(result->end(),deexResult->begin(), deexResult->end());
121 delete deexResult;
122}
G4ReactionProductVector * BreakItUp(const G4Fragment &theInitialState)

References G4ExcitationHandler::BreakItUp(), and G4VPreCompoundModel::GetExcitationHandler().

Referenced by DeExcite().

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

◆ SetExcitationHandler()

void G4VPreCompoundModel::SetExcitationHandler ( G4ExcitationHandler ptr)
inlineinherited

◆ SetMaxEnergy() [1/3]

void G4HadronicInteraction::SetMaxEnergy ( const G4double  anEnergy)
inlineinherited

Definition at line 102 of file G4HadronicInteraction.hh.

103 { theMaxEnergy = anEnergy; }

References G4HadronicInteraction::theMaxEnergy.

Referenced by G4HadronicInteraction::ActivateFor(), G4IonINCLXXPhysics::AddProcess(), G4BertiniElectroNuclearBuilder::Build(), G4LENDBertiniGammaElectroNuclearBuilder::Build(), G4NeutronLENDBuilder::Build(), G4NeutronPHPBuilder::Build(), G4AlphaPHPBuilder::Build(), G4BertiniKaonBuilder::Build(), G4BertiniNeutronBuilder::Build(), G4BertiniPiKBuilder::Build(), G4BertiniPionBuilder::Build(), G4BertiniProtonBuilder::Build(), G4BinaryAlphaBuilder::Build(), G4BinaryDeuteronBuilder::Build(), G4BinaryHe3Builder::Build(), G4BinaryNeutronBuilder::Build(), G4BinaryPiKBuilder::Build(), G4BinaryPionBuilder::Build(), G4BinaryProtonBuilder::Build(), G4BinaryTritonBuilder::Build(), G4DeuteronPHPBuilder::Build(), G4FTFBinaryKaonBuilder::Build(), G4FTFBinaryNeutronBuilder::Build(), G4FTFBinaryPionBuilder::Build(), G4FTFBinaryProtonBuilder::Build(), G4FTFPAntiBarionBuilder::Build(), G4FTFPKaonBuilder::Build(), G4FTFPNeutronBuilder::Build(), G4FTFPPiKBuilder::Build(), G4FTFPPionBuilder::Build(), G4FTFPProtonBuilder::Build(), G4He3PHPBuilder::Build(), G4HyperonFTFPBuilder::Build(), G4HyperonQGSPBuilder::Build(), G4INCLXXNeutronBuilder::Build(), G4INCLXXPionBuilder::Build(), G4INCLXXProtonBuilder::Build(), G4PrecoNeutronBuilder::Build(), G4PrecoProtonBuilder::Build(), G4ProtonPHPBuilder::Build(), G4QGSBinaryKaonBuilder::Build(), G4QGSBinaryNeutronBuilder::Build(), G4QGSBinaryPiKBuilder::Build(), G4QGSBinaryPionBuilder::Build(), G4QGSBinaryProtonBuilder::Build(), G4QGSPAntiBarionBuilder::Build(), G4QGSPKaonBuilder::Build(), G4QGSPLundStrFragmProtonBuilder::Build(), G4QGSPNeutronBuilder::Build(), G4QGSPPiKBuilder::Build(), G4QGSPPionBuilder::Build(), G4TritonPHPBuilder::Build(), G4QGSPProtonBuilder::Build(), G4HadronicBuilder::BuildFTFP_BERT(), G4HadronicBuilder::BuildFTFQGSP_BERT(), G4QGSBuilder::BuildModel(), G4VHadronPhysics::BuildModel(), G4HadronicBuilder::BuildQGSP_FTFP_BERT(), G4EmExtraPhysics::ConstructGammaElectroNuclear(), LBE::ConstructHad(), G4EmExtraPhysics::ConstructLENDGammaNuclear(), G4HadronDElasticPhysics::ConstructProcess(), G4HadronElasticPhysics::ConstructProcess(), G4HadronHElasticPhysics::ConstructProcess(), G4IonINCLXXPhysics::ConstructProcess(), G4IonPhysics::ConstructProcess(), G4IonPhysicsPHP::ConstructProcess(), G4IonQMDPhysics::ConstructProcess(), G4ANuElNucleusNcModel::G4ANuElNucleusNcModel(), G4ANuMuNucleusNcModel::G4ANuMuNucleusNcModel(), G4BertiniKaonBuilder::G4BertiniKaonBuilder(), G4BertiniPiKBuilder::G4BertiniPiKBuilder(), G4BertiniPionBuilder::G4BertiniPionBuilder(), G4BinaryCascade::G4BinaryCascade(), G4BinaryPiKBuilder::G4BinaryPiKBuilder(), G4BinaryPionBuilder::G4BinaryPionBuilder(), G4ChargeExchange::G4ChargeExchange(), G4DiffuseElastic::G4DiffuseElastic(), G4DiffuseElasticV2::G4DiffuseElasticV2(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4EMDissociation::G4EMDissociation(), G4FissLib::G4FissLib(), G4FTFBinaryKaonBuilder::G4FTFBinaryKaonBuilder(), G4FTFBinaryNeutronBuilder::G4FTFBinaryNeutronBuilder(), G4FTFBinaryPiKBuilder::G4FTFBinaryPiKBuilder(), G4FTFBinaryPionBuilder::G4FTFBinaryPionBuilder(), G4FTFBinaryProtonBuilder::G4FTFBinaryProtonBuilder(), G4FTFPAntiBarionBuilder::G4FTFPAntiBarionBuilder(), G4FTFPKaonBuilder::G4FTFPKaonBuilder(), G4FTFPNeutronBuilder::G4FTFPNeutronBuilder(), G4FTFPPiKBuilder::G4FTFPPiKBuilder(), G4FTFPPionBuilder::G4FTFPPionBuilder(), G4FTFPProtonBuilder::G4FTFPProtonBuilder(), G4HadronElastic::G4HadronElastic(), G4HadronicAbsorptionFritiof::G4HadronicAbsorptionFritiof(), G4HadronicAbsorptionFritiofWithBinaryCascade::G4HadronicAbsorptionFritiofWithBinaryCascade(), G4hhElastic::G4hhElastic(), G4HyperonFTFPBuilder::G4HyperonFTFPBuilder(), G4HyperonQGSPBuilder::G4HyperonQGSPBuilder(), G4INCLXXPionBuilder::G4INCLXXPionBuilder(), G4LEHadronProtonElastic::G4LEHadronProtonElastic(), G4LENDModel::G4LENDModel(), G4LEnp::G4LEnp(), G4LEpp::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

◆ fHighLimitExc

G4double G4PreCompoundModel::fHighLimitExc = DBL_MAX
private

Definition at line 101 of file G4PreCompoundModel.hh.

Referenced by DeExcite(), and InitialiseModel().

◆ fLowLimitExc

G4double G4PreCompoundModel::fLowLimitExc = 0.0
private

Definition at line 100 of file G4PreCompoundModel.hh.

Referenced by InitialiseModel().

◆ fNuclData

G4NuclearLevelData* G4PreCompoundModel::fNuclData = nullptr
private

Definition at line 95 of file G4PreCompoundModel.hh.

Referenced by DeExcite(), G4PreCompoundModel(), and InitialiseModel().

◆ isActive

G4bool G4PreCompoundModel::isActive = true
private

Definition at line 105 of file G4PreCompoundModel.hh.

Referenced by DeExcite(), and InitialiseModel().

◆ isBlocked

G4bool G4HadronicInteraction::isBlocked
protectedinherited

◆ isInitialised

G4bool G4PreCompoundModel::isInitialised = false
private

Definition at line 104 of file G4PreCompoundModel.hh.

Referenced by DeExcite(), and InitialiseModel().

◆ minA

G4int G4PreCompoundModel::minA = 5
private

Definition at line 108 of file G4PreCompoundModel.hh.

Referenced by DeExcite(), and InitialiseModel().

◆ minZ

G4int G4PreCompoundModel::minZ = 3
private

Definition at line 107 of file G4PreCompoundModel.hh.

Referenced by DeExcite(), and InitialiseModel().

◆ modelID

G4int G4PreCompoundModel::modelID = -1
private

Definition at line 109 of file G4PreCompoundModel.hh.

Referenced by ApplyYourself(), and G4PreCompoundModel().

◆ neutron

const G4ParticleDefinition* G4PreCompoundModel::neutron
private

Definition at line 98 of file G4PreCompoundModel.hh.

Referenced by ApplyYourself(), and G4PreCompoundModel().

◆ proton

const G4ParticleDefinition* G4PreCompoundModel::proton
private

Definition at line 97 of file G4PreCompoundModel.hh.

Referenced by ApplyYourself(), and G4PreCompoundModel().

◆ recoilEnergyThreshold

G4double G4HadronicInteraction::recoilEnergyThreshold
privateinherited

◆ registry

G4HadronicInteractionRegistry* G4HadronicInteraction::registry
privateinherited

◆ theBlockedList

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

◆ theBlockedListElements

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

◆ theEmission

G4PreCompoundEmission* G4PreCompoundModel::theEmission = nullptr
private

Definition at line 93 of file G4PreCompoundModel.hh.

Referenced by DeExcite(), InitialiseModel(), and ~G4PreCompoundModel().

◆ theExcitationHandler

G4ExcitationHandler* G4VPreCompoundModel::theExcitationHandler
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().

◆ theResult

G4HadFinalState G4PreCompoundModel::theResult
private

Definition at line 111 of file G4PreCompoundModel.hh.

Referenced by ApplyYourself(), and ~G4PreCompoundModel().

◆ theTransition

G4VPreCompoundTransitions* G4PreCompoundModel::theTransition = nullptr
private

Definition at line 94 of file G4PreCompoundModel.hh.

Referenced by DeExcite(), InitialiseModel(), and ~G4PreCompoundModel().

◆ useSCO

G4bool G4PreCompoundModel::useSCO = false
private

Definition at line 103 of file G4PreCompoundModel.hh.

Referenced by DeExcite(), and InitialiseModel().

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