Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LossTableManager Class Reference

#include <G4LossTableManager.hh>

Public Member Functions

 ~G4LossTableManager ()
 
void Clear ()
 
void PreparePhysicsTable (const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
 
void PreparePhysicsTable (const G4ParticleDefinition *aParticle, G4VEmProcess *p)
 
void PreparePhysicsTable (const G4ParticleDefinition *aParticle, G4VMultipleScattering *p)
 
void BuildPhysicsTable (const G4ParticleDefinition *aParticle)
 
void BuildPhysicsTable (const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
 
G4double GetDEDX (const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetSubDEDX (const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetRange (const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetCSDARange (const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetRangeFromRestricteDEDX (const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetEnergy (const G4ParticleDefinition *aParticle, G4double range, const G4MaterialCutsCouple *couple)
 
G4double GetDEDXDispersion (const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double &length)
 
void Register (G4VEnergyLossProcess *p)
 
void DeRegister (G4VEnergyLossProcess *p)
 
void Register (G4VMultipleScattering *p)
 
void DeRegister (G4VMultipleScattering *p)
 
void Register (G4VEmProcess *p)
 
void DeRegister (G4VEmProcess *p)
 
void Register (G4VEmModel *p)
 
void DeRegister (G4VEmModel *p)
 
void Register (G4VEmFluctuationModel *p)
 
void DeRegister (G4VEmFluctuationModel *p)
 
void RegisterIon (const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
 
void RegisterExtraParticle (const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
 
void SetLossFluctuations (G4bool val)
 
void SetSubCutoff (G4bool val, const G4Region *r=0)
 
void SetIntegral (G4bool val)
 
void SetRandomStep (G4bool val)
 
void SetMinSubRange (G4double val)
 
void SetMinEnergy (G4double val)
 
void SetMaxEnergy (G4double val)
 
void SetMaxEnergyForCSDARange (G4double val)
 
void SetMaxEnergyForMuons (G4double val)
 
void SetDEDXBinning (G4int val)
 
void SetDEDXBinningForCSDARange (G4int val)
 
void SetLambdaBinning (G4int val)
 
G4int GetNumberOfBinsPerDecade () const
 
void SetStepFunction (G4double v1, G4double v2)
 
void SetBuildCSDARange (G4bool val)
 
void SetLPMFlag (G4bool val)
 
void SetSplineFlag (G4bool val)
 
void SetLinearLossLimit (G4double val)
 
void SetBremsstrahlungTh (G4double val)
 
void SetFactorForAngleLimit (G4double val)
 
void SetVerbose (G4int val)
 
G4EnergyLossMessengerGetMessenger ()
 
G4bool BuildCSDARange () const
 
G4bool LPMFlag () const
 
G4bool SplineFlag () const
 
G4double BremsstrahlungTh () const
 
G4double FactorForAngleLimit () const
 
G4double MinKinEnergy () const
 
G4double MaxKinEnergy () const
 
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector ()
 
const std::vector< G4VEmProcess * > & GetEmProcessVector ()
 
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector ()
 
G4VEnergyLossProcessGetEnergyLossProcess (const G4ParticleDefinition *)
 
G4EmCorrectionsEmCorrections ()
 
G4EmSaturationEmSaturation ()
 
G4EmConfiguratorEmConfigurator ()
 
G4ElectronIonPairElectronIonPair ()
 
G4VAtomDeexcitationAtomDeexcitation ()
 
G4LossTableBuilderGetTableBuilder ()
 
void SetAtomDeexcitation (G4VAtomDeexcitation *)
 

Static Public Member Functions

static G4LossTableManagerInstance ()
 

Detailed Description

Definition at line 99 of file G4LossTableManager.hh.

Constructor & Destructor Documentation

◆ ~G4LossTableManager()

G4LossTableManager::~G4LossTableManager ( )

Definition at line 118 of file G4LossTableManager.cc.

119{
120 for (G4int i=0; i<n_loss; ++i) {
121 if( loss_vector[i] ) { delete loss_vector[i]; }
122 }
123 size_t msc = msc_vector.size();
124 for (size_t j=0; j<msc; ++j) {
125 if( msc_vector[j] ) { delete msc_vector[j]; }
126 }
127 size_t emp = emp_vector.size();
128 for (size_t k=0; k<emp; ++k) {
129 if( emp_vector[k] ) { delete emp_vector[k]; }
130 }
131 size_t mod = mod_vector.size();
132 for (size_t a=0; a<mod; ++a) {
133 if( mod_vector[a] ) { delete mod_vector[a]; }
134 }
135 size_t fmod = fmod_vector.size();
136 for (size_t b=0; b<fmod; ++b) {
137 if( fmod_vector[b] ) { delete fmod_vector[b]; }
138 }
139 Clear();
140 delete theMessenger;
141 delete tableBuilder;
142 delete emCorrections;
143 delete emSaturation;
144 delete emConfigurator;
145 delete emElectronIonPair;
146 delete atomDeexcitation;
147}
int G4int
Definition: G4Types.hh:66

Member Function Documentation

◆ AtomDeexcitation()

◆ BremsstrahlungTh()

G4double G4LossTableManager::BremsstrahlungTh ( ) const

Definition at line 1043 of file G4LossTableManager.cc.

1044{
1045 return bremsTh;
1046}

Referenced by G4eBremsstrahlung::InitialiseEnergyLossProcess(), and G4eBremsstrahlung::PrintInfo().

◆ BuildCSDARange()

G4bool G4LossTableManager::BuildCSDARange ( ) const

Definition at line 762 of file G4LossTableManager.cc.

763{
764 return buildCSDARange;
765}

Referenced by G4VEnergyLossProcess::PreparePhysicsTable().

◆ BuildPhysicsTable() [1/2]

void G4LossTableManager::BuildPhysicsTable ( const G4ParticleDefinition aParticle)

Definition at line 458 of file G4LossTableManager.cc.

459{
460 if(0 == run && startInitialisation) {
461 emConfigurator->Clear();
462 }
463}

Referenced by G4VEnergyLossProcess::BuildPhysicsTable(), and G4VMultipleScattering::BuildPhysicsTable().

◆ BuildPhysicsTable() [2/2]

void G4LossTableManager::BuildPhysicsTable ( const G4ParticleDefinition aParticle,
G4VEnergyLossProcess p 
)

Definition at line 467 of file G4LossTableManager.cc.

470{
471 if(1 < verbose) {
472 G4cout << "### G4LossTableManager::BuildDEDXTable() is requested for "
473 << aParticle->GetParticleName()
474 << " and process " << p->GetProcessName() << G4endl;
475 }
476 // clear configurator
477 if(0 == run && startInitialisation) {
478 emConfigurator->Clear();
479 firstParticle = aParticle;
480 }
481 if(startInitialisation && atomDeexcitation) {
482 atomDeexcitation->InitialiseAtomicDeexcitation();
483 }
484 startInitialisation = false;
485
486 // initialisation before any table is built
487 if ( aParticle == firstParticle ) {
488 all_tables_are_built = true;
489
490 if(1 < verbose) {
491 G4cout << "### G4LossTableManager start initilisation for first particle "
492 << firstParticle->GetParticleName()
493 << G4endl;
494 }
495 for (G4int i=0; i<n_loss; ++i) {
496 G4VEnergyLossProcess* el = loss_vector[i];
497
498 if(el) {
499 const G4ProcessManager* pm = el->GetProcessManager();
500 isActive[i] = false;
501 if(pm) { isActive[i] = pm->GetProcessActivation(el); }
502 if(0 == run) { base_part_vector[i] = el->BaseParticle(); }
503 tables_are_built[i] = false;
504 all_tables_are_built= false;
505 if(!isActive[i]) { el->SetIonisation(false); }
506
507 if(1 < verbose) {
508 G4cout << i <<". "<< el->GetProcessName();
509 if(el->Particle()) {
510 G4cout << " for " << el->Particle()->GetParticleName();
511 }
512 G4cout << " active= " << isActive[i]
513 << " table= " << tables_are_built[i]
514 << " isIonisation= " << el->IsIonisationProcess();
515 if(base_part_vector[i]) {
516 G4cout << " base particle " << base_part_vector[i]->GetParticleName();
517 }
518 G4cout << G4endl;
519 }
520 } else {
521 tables_are_built[i] = true;
522 part_vector[i] = 0;
523 }
524 }
525 ++run;
526 currentParticle = 0;
527 }
528
529 // Set run time parameters
530 SetParameters(aParticle, p);
531
532 if (all_tables_are_built) { return; }
533
534 // Build tables for given particle
535 all_tables_are_built = true;
536
537 for(G4int i=0; i<n_loss; ++i) {
538 if(p == loss_vector[i] && !tables_are_built[i] && !base_part_vector[i]) {
539 const G4ParticleDefinition* curr_part = part_vector[i];
540 if(1 < verbose) {
541 G4cout << "### BuildPhysicsTable for " << p->GetProcessName()
542 << " and " << curr_part->GetParticleName()
543 << " start BuildTable " << G4endl;
544 }
545 G4VEnergyLossProcess* curr_proc = BuildTables(curr_part);
546 if(curr_proc) { CopyTables(curr_part, curr_proc); }
547 }
548 if ( !tables_are_built[i] ) { all_tables_are_built = false; }
549 }
550
551 if(1 < verbose) {
552 G4cout << "### G4LossTableManager::BuildDEDXTable end: "
553 << "all_tables_are_built= " << all_tables_are_built
554 << G4endl;
555
556 if(all_tables_are_built) {
557 G4cout << "### All dEdx and Range tables are built #####" << G4endl;
558 }
559 }
560}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
const G4String & GetParticleName() const
G4bool GetProcessActivation(G4VProcess *aProcess) const
const G4ParticleDefinition * BaseParticle() const
const G4ParticleDefinition * Particle() const
G4bool IsIonisationProcess() const
void SetIonisation(G4bool val)
virtual const G4ProcessManager * GetProcessManager()
Definition: G4VProcess.hh:485
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ Clear()

void G4LossTableManager::Clear ( )

Definition at line 196 of file G4LossTableManager.cc.

197{
198 all_tables_are_built = false;
199 currentLoss = 0;
200 currentParticle = 0;
201 if(n_loss)
202 {
203 dedx_vector.clear();
204 range_vector.clear();
205 inv_range_vector.clear();
206 loss_map.clear();
207 loss_vector.clear();
208 part_vector.clear();
209 base_part_vector.clear();
210 tables_are_built.clear();
211 isActive.clear();
212 n_loss = 0;
213 }
214}

Referenced by ~G4LossTableManager().

◆ DeRegister() [1/5]

void G4LossTableManager::DeRegister ( G4VEmFluctuationModel p)

Definition at line 345 of file G4LossTableManager.cc.

346{
347 size_t n = fmod_vector.size();
348 for (size_t i=0; i<n; ++i) {
349 if(fmod_vector[i] == p) { fmod_vector[i] = 0; }
350 }
351}

◆ DeRegister() [2/5]

void G4LossTableManager::DeRegister ( G4VEmModel p)

Definition at line 324 of file G4LossTableManager.cc.

325{
326 size_t n = mod_vector.size();
327 for (size_t i=0; i<n; ++i) {
328 if(mod_vector[i] == p) { mod_vector[i] = 0; }
329 }
330}

◆ DeRegister() [3/5]

void G4LossTableManager::DeRegister ( G4VEmProcess p)

Definition at line 302 of file G4LossTableManager.cc.

303{
304 if(!p) { return; }
305 size_t emp = emp_vector.size();
306 for (size_t i=0; i<emp; ++i) {
307 if(emp_vector[i] == p) { emp_vector[i] = 0; }
308 }
309}

◆ DeRegister() [4/5]

void G4LossTableManager::DeRegister ( G4VEnergyLossProcess p)

Definition at line 249 of file G4LossTableManager.cc.

250{
251 if(!p) { return; }
252 for (G4int i=0; i<n_loss; ++i) {
253 if(loss_vector[i] == p) { loss_vector[i] = 0; }
254 }
255}

Referenced by G4VEmFluctuationModel::~G4VEmFluctuationModel(), G4VEmModel::~G4VEmModel(), and G4VMultipleScattering::~G4VMultipleScattering().

◆ DeRegister() [5/5]

void G4LossTableManager::DeRegister ( G4VMultipleScattering p)

Definition at line 275 of file G4LossTableManager.cc.

276{
277 if(!p) { return; }
278 size_t msc = msc_vector.size();
279 for (size_t i=0; i<msc; ++i) {
280 if(msc_vector[i] == p) { msc_vector[i] = 0; }
281 }
282}

◆ ElectronIonPair()

G4ElectronIonPair * G4LossTableManager::ElectronIonPair ( )

Definition at line 1099 of file G4LossTableManager.cc.

1100{
1101 return emElectronIonPair;
1102}

◆ EmConfigurator()

G4EmConfigurator * G4LossTableManager::EmConfigurator ( )

Definition at line 1092 of file G4LossTableManager.cc.

1093{
1094 return emConfigurator;
1095}

◆ EmCorrections()

◆ EmSaturation()

G4EmSaturation * G4LossTableManager::EmSaturation ( )

Definition at line 1085 of file G4LossTableManager.cc.

1086{
1087 return emSaturation;
1088}

Referenced by G4OpticalPhysics::ConstructProcess().

◆ FactorForAngleLimit()

G4double G4LossTableManager::FactorForAngleLimit ( ) const

Definition at line 1057 of file G4LossTableManager.cc.

1058{
1059 return factorForAngleLimit;
1060}

Referenced by G4WentzelOKandVIxSection::Initialise(), G4WentzelVIRelXSection::Initialise(), and G4CoulombScattering::InitialiseProcess().

◆ GetCSDARange()

G4double G4LossTableManager::GetCSDARange ( const G4ParticleDefinition aParticle,
G4double  kineticEnergy,
const G4MaterialCutsCouple couple 
)

Definition at line 1171 of file G4LossTableManager.cc.

1174{
1175 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1176 G4double x = DBL_MAX;
1177 if(currentLoss) { x = currentLoss->GetCSDARange(kineticEnergy, couple); }
1178 return x;
1179}
double G4double
Definition: G4Types.hh:64
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
G4double GetCSDARange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
#define DBL_MAX
Definition: templates.hh:83

Referenced by G4EmCalculator::GetCSDARange().

◆ GetDEDX()

G4double G4LossTableManager::GetDEDX ( const G4ParticleDefinition aParticle,
G4double  kineticEnergy,
const G4MaterialCutsCouple couple 
)

Definition at line 1145 of file G4LossTableManager.cc.

1148{
1149 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1150 G4double x = 0.0;
1151 if(currentLoss) { x = currentLoss->GetDEDX(kineticEnergy, couple); }
1152 else { x = G4EnergyLossTables::GetDEDX(currentParticle,
1153 kineticEnergy,couple,false); }
1154 return x;
1155}
static G4double GetDEDX(const G4ParticleDefinition *aParticle, G4double KineticEnergy, const G4Material *aMaterial)
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)

Referenced by G4EnergyLossTables::GetDEDX(), G4EmCalculator::GetDEDX(), G4EnergyLossTables::GetPreciseDEDX(), and G4EnergyLossTables::GetPreciseRangeFromEnergy().

◆ GetDEDXDispersion()

G4double G4LossTableManager::GetDEDXDispersion ( const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double length 
)

Definition at line 1226 of file G4LossTableManager.cc.

1229{
1230 const G4ParticleDefinition* aParticle = dp->GetParticleDefinition();
1231 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1232 G4double x = 0.0;
1233 if(currentLoss) { currentLoss->GetDEDXDispersion(couple, dp, length); }
1234 return x;
1235}
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)

◆ GetEmProcessVector()

◆ GetEnergy()

G4double G4LossTableManager::GetEnergy ( const G4ParticleDefinition aParticle,
G4double  range,
const G4MaterialCutsCouple couple 
)

Definition at line 1212 of file G4LossTableManager.cc.

1215{
1216 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1217 G4double x;
1218 if(currentLoss) { x = currentLoss->GetKineticEnergy(range, couple); }
1219 else { x = G4EnergyLossTables::GetPreciseEnergyFromRange(currentParticle,range,
1220 couple,false); }
1221 return x;
1222}
static G4double GetPreciseEnergyFromRange(const G4ParticleDefinition *aParticle, G4double range, const G4Material *aMaterial)
G4double GetKineticEnergy(G4double &range, const G4MaterialCutsCouple *)

Referenced by G4EmCalculator::GetKinEnergy(), and G4EnergyLossTables::GetPreciseEnergyFromRange().

◆ GetEnergyLossProcess()

G4VEnergyLossProcess * G4LossTableManager::GetEnergyLossProcess ( const G4ParticleDefinition aParticle)

Definition at line 1128 of file G4LossTableManager.cc.

1129{
1130 if(aParticle != currentParticle) {
1131 currentParticle = aParticle;
1132 std::map<PD,G4VEnergyLossProcess*,std::less<PD> >::const_iterator pos;
1133 if ((pos = loss_map.find(aParticle)) != loss_map.end()) {
1134 currentLoss = (*pos).second;
1135 } else {
1136 currentLoss = 0;
1137 //ParticleHaveNoLoss(aParticle);
1138 }
1139 }
1140 return currentLoss;
1141}

Referenced by GetCSDARange(), GetDEDX(), GetDEDXDispersion(), GetEnergy(), GetRange(), GetRangeFromRestricteDEDX(), GetSubDEDX(), and G4VMultipleScattering::StartTracking().

◆ GetEnergyLossProcessVector()

◆ GetMessenger()

G4EnergyLossMessenger * G4LossTableManager::GetMessenger ( )

Definition at line 742 of file G4LossTableManager.cc.

743{
744 return theMessenger;
745}

◆ GetMultipleScatteringVector()

◆ GetNumberOfBinsPerDecade()

G4int G4LossTableManager::GetNumberOfBinsPerDecade ( ) const

Definition at line 908 of file G4LossTableManager.cc.

909{
910 return nbinsPerDecade;
911}

Referenced by G4LossTableBuilder::BuildTableForModel(), and G4VEmModel::InitialiseElementSelectors().

◆ GetRange()

G4double G4LossTableManager::GetRange ( const G4ParticleDefinition aParticle,
G4double  kineticEnergy,
const G4MaterialCutsCouple couple 
)

Definition at line 1198 of file G4LossTableManager.cc.

1201{
1202 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1203 G4double x;
1204 if(currentLoss) { x = currentLoss->GetRange(kineticEnergy, couple); }
1205 else { x = G4EnergyLossTables::GetRange(currentParticle,kineticEnergy,
1206 couple,false); }
1207 return x;
1208}
static G4double GetRange(const G4ParticleDefinition *aParticle, G4double KineticEnergy, const G4Material *aMaterial)
G4double GetRange(G4double &kineticEnergy, const G4MaterialCutsCouple *)

Referenced by G4ITStepProcessor::ApplyProductionCut(), G4EnergyLossTables::GetRange(), G4EmCalculator::GetRange(), G4UserSpecialCuts::PostStepGetPhysicalInteractionLength(), and G4EmSaturation::VisibleEnergyDeposition().

◆ GetRangeFromRestricteDEDX()

G4double G4LossTableManager::GetRangeFromRestricteDEDX ( const G4ParticleDefinition aParticle,
G4double  kineticEnergy,
const G4MaterialCutsCouple couple 
)

Definition at line 1184 of file G4LossTableManager.cc.

1187{
1188 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1189 G4double x;
1190 if(currentLoss) { x = currentLoss->GetRangeForLoss(kineticEnergy, couple); }
1191 else { x = G4EnergyLossTables::GetRange(currentParticle,kineticEnergy,
1192 couple,false); }
1193 return x;
1194}
G4double GetRangeForLoss(G4double &kineticEnergy, const G4MaterialCutsCouple *)

Referenced by G4EmCalculator::GetRangeFromRestricteDEDX().

◆ GetSubDEDX()

G4double G4LossTableManager::GetSubDEDX ( const G4ParticleDefinition aParticle,
G4double  kineticEnergy,
const G4MaterialCutsCouple couple 
)

Definition at line 1159 of file G4LossTableManager.cc.

1162{
1163 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1164 G4double x = 0.0;
1165 if(currentLoss) { x = currentLoss->GetDEDXForSubsec(kineticEnergy, couple); }
1166 return x;
1167}
G4double GetDEDXForSubsec(G4double &kineticEnergy, const G4MaterialCutsCouple *)

◆ GetTableBuilder()

G4LossTableBuilder * G4LossTableManager::GetTableBuilder ( )

◆ Instance()

G4LossTableManager * G4LossTableManager::Instance ( )
static

Definition at line 107 of file G4LossTableManager.cc.

108{
109 if(0 == theInstance) {
110 static G4LossTableManager manager;
111 theInstance = &manager;
112 }
113 return theInstance;
114}

Referenced by G4ITStepProcessor::ApplyProductionCut(), G4VEnergyLossProcess::BuildDEDXTable(), G4VEnergyLossProcess::BuildLambdaTable(), G4VEmProcess::BuildPhysicsTable(), G4VEnergyLossProcess::BuildPhysicsTable(), G4RadioactiveDecay::BuildPhysicsTable(), G4LossTableBuilder::BuildTableForModel(), G4EmCalculator::ComputeElectronicDEDX(), G4EmDNAPhysics::ConstructProcess(), G4EmDNAPhysicsChemistry::ConstructProcess(), G4EmLivermorePhysics::ConstructProcess(), G4EmLivermorePolarizedPhysics::ConstructProcess(), G4EmLowEPPhysics::ConstructProcess(), G4EmPenelopePhysics::ConstructProcess(), G4EmStandardPhysics::ConstructProcess(), G4EmStandardPhysics_option1::ConstructProcess(), G4EmStandardPhysics_option2::ConstructProcess(), G4EmStandardPhysics_option3::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4OpticalPhysics::ConstructProcess(), G4NuclearDecayChannel::DecayIt(), G4BetheBlochModel::G4BetheBlochModel(), G4BraggModel::G4BraggModel(), G4EmCalculator::G4EmCalculator(), G4EmLivermorePhysics::G4EmLivermorePhysics(), G4EmLivermorePolarizedPhysics::G4EmLivermorePolarizedPhysics(), G4EmLowEPPhysics::G4EmLowEPPhysics(), G4EmPenelopePhysics::G4EmPenelopePhysics(), G4EmProcessOptions::G4EmProcessOptions(), G4EmStandardPhysics::G4EmStandardPhysics(), G4EmStandardPhysics_option1::G4EmStandardPhysics_option1(), G4EmStandardPhysics_option2::G4EmStandardPhysics_option2(), G4EmStandardPhysics_option3::G4EmStandardPhysics_option3(), G4EmStandardPhysics_option4::G4EmStandardPhysics_option4(), G4GoudsmitSaundersonMscModel::G4GoudsmitSaundersonMscModel(), G4ionIonisation::G4ionIonisation(), G4IonParametrisedLossModel::G4IonParametrisedLossModel(), G4MuBetheBlochModel::G4MuBetheBlochModel(), G4QAtomicPhysics::G4QAtomicPhysics(), G4UAtomicDeexcitation::G4UAtomicDeexcitation(), G4UrbanMscModel90::G4UrbanMscModel90(), G4UrbanMscModel92::G4UrbanMscModel92(), G4UrbanMscModel93::G4UrbanMscModel93(), G4UrbanMscModel95::G4UrbanMscModel95(), G4UrbanMscModel96::G4UrbanMscModel96(), G4UserSpecialCuts::G4UserSpecialCuts(), G4VEmFluctuationModel::G4VEmFluctuationModel(), G4VEmModel::G4VEmModel(), G4VEmProcess::G4VEmProcess(), G4VEnergyLossProcess::G4VEnergyLossProcess(), G4VMscModel::G4VMscModel(), G4VMultipleScattering::G4VMultipleScattering(), G4WentzelVIModel::G4WentzelVIModel(), G4WentzelVIRelModel::G4WentzelVIRelModel(), G4EnergyLossTables::GetDEDX(), G4EnergyLossTables::GetPreciseDEDX(), G4EnergyLossTables::GetPreciseEnergyFromRange(), G4EnergyLossTables::GetPreciseRangeFromEnergy(), G4EnergyLossTables::GetRange(), G4DNARuddIonisationExtendedModel::Initialise(), G4DNARuddIonisationModel::Initialise(), G4LivermoreComptonModel::Initialise(), G4LivermoreComptonModifiedModel::Initialise(), G4LivermorePhotoElectricModel::Initialise(), G4LivermorePolarizedPhotoElectricModel::Initialise(), G4LowEPComptonModel::Initialise(), G4MuElecInelasticModel::Initialise(), G4PenelopeComptonModel::Initialise(), G4PenelopeIonisationModel::Initialise(), G4PenelopePhotoElectricModel::Initialise(), G4BraggIonModel::Initialise(), G4KleinNishinaModel::Initialise(), G4PEEffectFluoModel::Initialise(), G4DNABornIonisationModel::Initialise(), G4WentzelOKandVIxSection::Initialise(), G4WentzelVIRelXSection::Initialise(), G4VEmModel::InitialiseElementSelectors(), G4hhIonisation::InitialiseEnergyLossProcess(), G4mplIonisation::InitialiseEnergyLossProcess(), G4eBremsstrahlung::InitialiseEnergyLossProcess(), G4CoulombScattering::InitialiseProcess(), G4VEmProcess::LambdaPhysicsVector(), G4VEnergyLossProcess::LambdaPhysicsVector(), G4Cerenkov::PostStepGetPhysicalInteractionLength(), G4VEmProcess::PreparePhysicsTable(), G4VEnergyLossProcess::PreparePhysicsTable(), G4eBremsstrahlung::PrintInfo(), G4VEnergyLossProcess::PrintInfoDefinition(), G4VEmProcess::RetrievePhysicsTable(), G4VEmFluctuationModel::~G4VEmFluctuationModel(), G4VEmModel::~G4VEmModel(), G4VEmProcess::~G4VEmProcess(), and G4VEnergyLossProcess::~G4VEnergyLossProcess().

◆ LPMFlag()

G4bool G4LossTableManager::LPMFlag ( ) const

Definition at line 1014 of file G4LossTableManager.cc.

1015{
1016 return flagLPM;
1017}

Referenced by G4eBremsstrahlung::InitialiseEnergyLossProcess(), and G4eBremsstrahlung::PrintInfo().

◆ MaxKinEnergy()

G4double G4LossTableManager::MaxKinEnergy ( ) const

Definition at line 1071 of file G4LossTableManager.cc.

1072{
1073 return maxKinEnergy;
1074}

Referenced by G4VMscModel::GetParticleChangeForMSC(), and G4VMultipleScattering::PreparePhysicsTable().

◆ MinKinEnergy()

G4double G4LossTableManager::MinKinEnergy ( ) const

Definition at line 1064 of file G4LossTableManager.cc.

1065{
1066 return minKinEnergy;
1067}

Referenced by G4VMscModel::GetParticleChangeForMSC().

◆ PreparePhysicsTable() [1/3]

void G4LossTableManager::PreparePhysicsTable ( const G4ParticleDefinition aParticle,
G4VEmProcess p 
)

Definition at line 419 of file G4LossTableManager.cc.

421{
422 if (1 < verbose) {
423 G4cout << "G4LossTableManager::PreparePhysicsTable for "
424 << particle->GetParticleName()
425 << " and " << p->GetProcessName() << G4endl;
426 }
427 if(!startInitialisation) { tableBuilder->SetInitialisationFlag(false); }
428
429 // start initialisation for the first run
430 if( 0 == run ) {
431 emConfigurator->PrepareModels(particle, p);
432 }
433 startInitialisation = true;
434}
void PrepareModels(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
void SetInitialisationFlag(G4bool flag)

◆ PreparePhysicsTable() [2/3]

void G4LossTableManager::PreparePhysicsTable ( const G4ParticleDefinition aParticle,
G4VEnergyLossProcess p 
)

Definition at line 390 of file G4LossTableManager.cc.

392{
393 if (1 < verbose) {
394 G4cout << "G4LossTableManager::PreparePhysicsTable for "
395 << particle->GetParticleName()
396 << " and " << p->GetProcessName() << " run= " << run
397 << " loss_vector " << loss_vector.size() << G4endl;
398 }
399 if(!startInitialisation) { tableBuilder->SetInitialisationFlag(false); }
400
401 // start initialisation for the first run
402 startInitialisation = true;
403
404 if( 0 == run ) {
405 emConfigurator->PrepareModels(particle, p);
406
407 // initialise particles for given process
408 for (G4int j=0; j<n_loss; ++j) {
409 if (p == loss_vector[j]) {
410 if (!part_vector[j]) { part_vector[j] = particle; }
411 }
412 }
413 }
414}

Referenced by G4VEmProcess::PreparePhysicsTable(), G4VEnergyLossProcess::PreparePhysicsTable(), and G4VMultipleScattering::PreparePhysicsTable().

◆ PreparePhysicsTable() [3/3]

void G4LossTableManager::PreparePhysicsTable ( const G4ParticleDefinition aParticle,
G4VMultipleScattering p 
)

Definition at line 439 of file G4LossTableManager.cc.

441{
442 if (1 < verbose) {
443 G4cout << "G4LossTableManager::PreparePhysicsTable for "
444 << particle->GetParticleName()
445 << " and " << p->GetProcessName() << G4endl;
446 }
447 if(!startInitialisation) { tableBuilder->SetInitialisationFlag(false); }
448
449 // start initialisation for the first run
450 if( 0 == run ) {
451 emConfigurator->PrepareModels(particle, p);
452 }
453}

◆ Register() [1/5]

void G4LossTableManager::Register ( G4VEmFluctuationModel p)

Definition at line 334 of file G4LossTableManager.cc.

335{
336 fmod_vector.push_back(p);
337 if(verbose > 1) {
338 G4cout << "G4LossTableManager::Register G4VEmFluctuationModel : "
339 << p->GetName() << G4endl;
340 }
341}

◆ Register() [2/5]

void G4LossTableManager::Register ( G4VEmModel p)

Definition at line 313 of file G4LossTableManager.cc.

314{
315 mod_vector.push_back(p);
316 if(verbose > 1) {
317 G4cout << "G4LossTableManager::Register G4VEmModel : "
318 << p->GetName() << G4endl;
319 }
320}
const G4String & GetName() const
Definition: G4VEmModel.hh:655

◆ Register() [3/5]

void G4LossTableManager::Register ( G4VEmProcess p)

Definition at line 286 of file G4LossTableManager.cc.

287{
288 if(!p) { return; }
289 G4int n = emp_vector.size();
290 for (G4int i=0; i<n; ++i) {
291 if(emp_vector[i] == p) { return; }
292 }
293 if(verbose > 1) {
294 G4cout << "G4LossTableManager::Register G4VEmProcess : "
295 << p->GetProcessName() << " idx= " << emp_vector.size() << G4endl;
296 }
297 emp_vector.push_back(p);
298}

◆ Register() [4/5]

void G4LossTableManager::Register ( G4VEnergyLossProcess p)

Definition at line 218 of file G4LossTableManager.cc.

219{
220 if(!p) { return; }
221 for (G4int i=0; i<n_loss; ++i) {
222 if(loss_vector[i] == p) { return; }
223 }
224 if(verbose > 1) {
225 G4cout << "G4LossTableManager::Register G4VEnergyLossProcess : "
226 << p->GetProcessName() << " idx= " << n_loss << G4endl;
227 }
228 ++n_loss;
229 loss_vector.push_back(p);
230 part_vector.push_back(0);
231 base_part_vector.push_back(0);
232 dedx_vector.push_back(0);
233 range_vector.push_back(0);
234 inv_range_vector.push_back(0);
235 tables_are_built.push_back(false);
236 isActive.push_back(true);
237 all_tables_are_built = false;
238 if(!lossFluctuationFlag) { p->SetLossFluctuations(false); }
239 if(subCutoffFlag) { p->ActivateSubCutoff(true); }
240 if(rndmStepFlag) { p->SetRandomStep(true); }
241 if(stepFunctionActive) { p->SetStepFunction(maxRangeVariation, maxFinalStep); }
242 if(integralActive) { p->SetIntegral(integral); }
243 if(minEnergyActive) { p->SetMinKinEnergy(minKinEnergy); }
244 if(maxEnergyActive) { p->SetMaxKinEnergy(maxKinEnergy); }
245}
void SetRandomStep(G4bool val)
void SetMaxKinEnergy(G4double e)
void ActivateSubCutoff(G4bool val, const G4Region *region=0)
void SetStepFunction(G4double v1, G4double v2)
void SetLossFluctuations(G4bool val)
void SetIntegral(G4bool val)
void SetMinKinEnergy(G4double e)

Referenced by G4VEmFluctuationModel::G4VEmFluctuationModel(), G4VEmModel::G4VEmModel(), and G4VMultipleScattering::G4VMultipleScattering().

◆ Register() [5/5]

void G4LossTableManager::Register ( G4VMultipleScattering p)

Definition at line 259 of file G4LossTableManager.cc.

260{
261 if(!p) { return; }
262 G4int n = msc_vector.size();
263 for (G4int i=0; i<n; ++i) {
264 if(msc_vector[i] == p) { return; }
265 }
266 if(verbose > 1) {
267 G4cout << "G4LossTableManager::Register G4VMultipleScattering : "
268 << p->GetProcessName() << " idx= " << msc_vector.size() << G4endl;
269 }
270 msc_vector.push_back(p);
271}

◆ RegisterExtraParticle()

void G4LossTableManager::RegisterExtraParticle ( const G4ParticleDefinition aParticle,
G4VEnergyLossProcess p 
)

Definition at line 363 of file G4LossTableManager.cc.

366{
367 if(!p || !part) { return; }
368 for (G4int i=0; i<n_loss; ++i) {
369 if(loss_vector[i] == p) { return; }
370 }
371 if(verbose > 1) {
372 G4cout << "G4LossTableManager::RegisterExtraParticle "
373 << part->GetParticleName() << " G4VEnergyLossProcess : "
374 << p->GetProcessName() << " idx= " << n_loss << G4endl;
375 }
376 ++n_loss;
377 loss_vector.push_back(p);
378 part_vector.push_back(part);
379 base_part_vector.push_back(p->BaseParticle());
380 dedx_vector.push_back(0);
381 range_vector.push_back(0);
382 inv_range_vector.push_back(0);
383 tables_are_built.push_back(false);
384 all_tables_are_built = false;
385}

Referenced by G4VEnergyLossProcess::PreparePhysicsTable().

◆ RegisterIon()

void G4LossTableManager::RegisterIon ( const G4ParticleDefinition aParticle,
G4VEnergyLossProcess p 
)

Definition at line 355 of file G4LossTableManager.cc.

357{
358 loss_map[ion] = p;
359}

Referenced by G4VEnergyLossProcess::PreparePhysicsTable().

◆ SetAtomDeexcitation()

◆ SetBremsstrahlungTh()

void G4LossTableManager::SetBremsstrahlungTh ( G4double  val)

Definition at line 1036 of file G4LossTableManager.cc.

1037{
1038 bremsTh = val;
1039}

Referenced by G4EmProcessOptions::SetBremsstrahlungTh().

◆ SetBuildCSDARange()

void G4LossTableManager::SetBuildCSDARange ( G4bool  val)

Definition at line 960 of file G4LossTableManager.cc.

961{
962 buildCSDARange = val;
963}

Referenced by G4EmProcessOptions::SetBuildCSDARange().

◆ SetDEDXBinning()

void G4LossTableManager::SetDEDXBinning ( G4int  val)

Definition at line 871 of file G4LossTableManager.cc.

872{
873 for(G4int i=0; i<n_loss; ++i) {
874 if(loss_vector[i]) { loss_vector[i]->SetDEDXBinning(val); }
875 }
876}

Referenced by G4EmProcessOptions::SetDEDXBinning().

◆ SetDEDXBinningForCSDARange()

void G4LossTableManager::SetDEDXBinningForCSDARange ( G4int  val)

Definition at line 880 of file G4LossTableManager.cc.

881{
882 for(G4int i=0; i<n_loss; ++i) {
883 if(loss_vector[i]) { loss_vector[i]->SetDEDXBinningForCSDARange(val); }
884 }
885}

Referenced by G4EmProcessOptions::SetDEDXBinningForCSDARange().

◆ SetFactorForAngleLimit()

void G4LossTableManager::SetFactorForAngleLimit ( G4double  val)

Definition at line 1050 of file G4LossTableManager.cc.

1051{
1052 if(val > 0.0) { factorForAngleLimit = val; }
1053}

Referenced by G4EmProcessOptions::SetFactorForAngleLimit().

◆ SetIntegral()

void G4LossTableManager::SetIntegral ( G4bool  val)

Definition at line 789 of file G4LossTableManager.cc.

790{
791 integral = val;
792 integralActive = true;
793 for(G4int i=0; i<n_loss; ++i) {
794 if(loss_vector[i]) { loss_vector[i]->SetIntegral(val); }
795 }
796 size_t emp = emp_vector.size();
797 for (size_t k=0; k<emp; ++k) {
798 if(emp_vector[k]) { emp_vector[k]->SetIntegral(val); }
799 }
800}

Referenced by G4EmProcessOptions::SetIntegral().

◆ SetLambdaBinning()

void G4LossTableManager::SetLambdaBinning ( G4int  val)

Definition at line 889 of file G4LossTableManager.cc.

890{
891 G4int n = val/G4int(std::log10(maxKinEnergy/minKinEnergy) + 0.5);
892 if(n < 5) {
893 G4cout << "G4LossTableManager::SetLambdaBinning WARNING "
894 << "too small number of bins " << val << " ignored"
895 << G4endl;
896 return;
897 }
898 nbinsLambda = val;
899 nbinsPerDecade = n;
900 size_t emp = emp_vector.size();
901 for (size_t k=0; k<emp; ++k) {
902 if(emp_vector[k]) { emp_vector[k]->SetLambdaBinning(val); }
903 }
904}

Referenced by G4EmProcessOptions::SetLambdaBinning().

◆ SetLinearLossLimit()

void G4LossTableManager::SetLinearLossLimit ( G4double  val)

Definition at line 951 of file G4LossTableManager.cc.

952{
953 for(G4int i=0; i<n_loss; ++i) {
954 if(loss_vector[i]) { loss_vector[i]->SetLinearLossLimit(val); }
955 }
956}

Referenced by G4EmProcessOptions::SetLinearLossLimit().

◆ SetLossFluctuations()

void G4LossTableManager::SetLossFluctuations ( G4bool  val)

Definition at line 769 of file G4LossTableManager.cc.

770{
771 lossFluctuationFlag = val;
772 for(G4int i=0; i<n_loss; ++i) {
773 if(loss_vector[i]) { loss_vector[i]->SetLossFluctuations(val); }
774 }
775}

Referenced by G4EmProcessOptions::SetLossFluctuations().

◆ SetLPMFlag()

void G4LossTableManager::SetLPMFlag ( G4bool  val)

Definition at line 1007 of file G4LossTableManager.cc.

1008{
1009 flagLPM = val;
1010}

Referenced by G4EmProcessOptions::SetLPMFlag().

◆ SetMaxEnergy()

void G4LossTableManager::SetMaxEnergy ( G4double  val)

Definition at line 839 of file G4LossTableManager.cc.

840{
841 maxEnergyActive = true;
842 maxKinEnergy = val;
843 for(G4int i=0; i<n_loss; ++i) {
844 if(loss_vector[i]) { loss_vector[i]->SetMaxKinEnergy(val); }
845 }
846 size_t emp = emp_vector.size();
847 for (size_t k=0; k<emp; ++k) {
848 if(emp_vector[k]) { emp_vector[k]->SetMaxKinEnergy(val); }
849 }
850}

Referenced by G4EmProcessOptions::SetMaxEnergy().

◆ SetMaxEnergyForCSDARange()

void G4LossTableManager::SetMaxEnergyForCSDARange ( G4double  val)

Definition at line 854 of file G4LossTableManager.cc.

855{
856 for(G4int i=0; i<n_loss; ++i) {
857 if(loss_vector[i]) { loss_vector[i]->SetMaxKinEnergyForCSDARange(val); }
858 }
859}

Referenced by G4EmProcessOptions::SetMaxEnergyForCSDARange().

◆ SetMaxEnergyForMuons()

void G4LossTableManager::SetMaxEnergyForMuons ( G4double  val)

Definition at line 863 of file G4LossTableManager.cc.

864{
865 maxEnergyForMuonsActive = true;
866 maxKinEnergyForMuons = val;
867}

Referenced by G4EmProcessOptions::SetMaxEnergyForMuons().

◆ SetMinEnergy()

void G4LossTableManager::SetMinEnergy ( G4double  val)

Definition at line 824 of file G4LossTableManager.cc.

825{
826 minEnergyActive = true;
827 minKinEnergy = val;
828 for(G4int i=0; i<n_loss; ++i) {
829 if(loss_vector[i]) { loss_vector[i]->SetMinKinEnergy(val); }
830 }
831 size_t emp = emp_vector.size();
832 for (size_t k=0; k<emp; ++k) {
833 if(emp_vector[k]) { emp_vector[k]->SetMinKinEnergy(val); }
834 }
835}

Referenced by G4EmProcessOptions::SetMinEnergy().

◆ SetMinSubRange()

void G4LossTableManager::SetMinSubRange ( G4double  val)

Definition at line 804 of file G4LossTableManager.cc.

805{
806 minSubRange = val;
807 for(G4int i=0; i<n_loss; ++i) {
808 if(loss_vector[i]) { loss_vector[i]->SetMinSubRange(val); }
809 }
810}

Referenced by G4EmProcessOptions::SetMinSubRange().

◆ SetRandomStep()

void G4LossTableManager::SetRandomStep ( G4bool  val)

Definition at line 814 of file G4LossTableManager.cc.

815{
816 rndmStepFlag = val;
817 for(G4int i=0; i<n_loss; ++i) {
818 if(loss_vector[i]) { loss_vector[i]->SetRandomStep(val); }
819 }
820}

Referenced by G4EmProcessOptions::SetRandomStep().

◆ SetSplineFlag()

void G4LossTableManager::SetSplineFlag ( G4bool  val)

Definition at line 1021 of file G4LossTableManager.cc.

1022{
1023 splineFlag = val;
1024 tableBuilder->SetSplineFlag(splineFlag);
1025}
void SetSplineFlag(G4bool flag)

Referenced by G4EmProcessOptions::SetSplineFlag().

◆ SetStepFunction()

void G4LossTableManager::SetStepFunction ( G4double  v1,
G4double  v2 
)

Definition at line 939 of file G4LossTableManager.cc.

940{
941 stepFunctionActive = true;
942 maxRangeVariation = v1;
943 maxFinalStep = v2;
944 for(G4int i=0; i<n_loss; ++i) {
945 if(loss_vector[i]) { loss_vector[i]->SetStepFunction(v1, v2); }
946 }
947}

Referenced by G4EmProcessOptions::SetStepFunction().

◆ SetSubCutoff()

void G4LossTableManager::SetSubCutoff ( G4bool  val,
const G4Region r = 0 
)

Definition at line 779 of file G4LossTableManager.cc.

780{
781 subCutoffFlag = val;
782 for(G4int i=0; i<n_loss; ++i) {
783 if(loss_vector[i]) { loss_vector[i]->ActivateSubCutoff(val, r); }
784 }
785}

Referenced by G4EmProcessOptions::SetSubCutoff().

◆ SetVerbose()

void G4LossTableManager::SetVerbose ( G4int  val)

Definition at line 915 of file G4LossTableManager.cc.

916{
917 verbose = val;
918 for(G4int i=0; i<n_loss; ++i) {
919 if(loss_vector[i]) { loss_vector[i]->SetVerboseLevel(val); }
920 }
921 size_t msc = msc_vector.size();
922 for (size_t j=0; j<msc; ++j) {
923 if(msc_vector[j]) { msc_vector[j]->SetVerboseLevel(val); }
924 }
925 size_t emp = emp_vector.size();
926 for (size_t k=0; k<emp; ++k) {
927 if(emp_vector[k]) { emp_vector[k]->SetVerboseLevel(val); }
928 }
929 emConfigurator->SetVerbose(val);
930 //tableBuilder->SetVerbose(val);
931 //emCorrections->SetVerbose(val);
932 emSaturation->SetVerbose(val);
933 emElectronIonPair->SetVerbose(val);
934 if(atomDeexcitation) { atomDeexcitation->SetVerboseLevel(val); }
935}
void SetVerbose(G4int value)
void SetVerbose(G4int)

Referenced by G4EmProcessOptions::SetVerbose().

◆ SplineFlag()

G4bool G4LossTableManager::SplineFlag ( ) const

Definition at line 1029 of file G4LossTableManager.cc.

1030{
1031 return splineFlag;
1032}

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