Geant4 11.4.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 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)
void LocalPhysicsTables (const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
void DumpHtml ()
G4double GetDEDX (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 (G4VProcess *p)
void DeRegister (G4VProcess *p)
void Register (G4VEmModel *p)
void DeRegister (G4VEmModel *p)
void Register (G4VEmFluctuationModel *p)
void DeRegister (G4VEmFluctuationModel *p)
void Register (G4VXRayModel *p)
void DeRegister (G4VXRayModel *p)
void RegisterExtraParticle (const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
void SetVerbose (G4int val)
void ResetParameters ()
void SetAtomDeexcitation (G4VAtomDeexcitation *)
void SetSubCutProducer (G4VSubCutProducer *)
void SetNIELCalculator (G4NIELCalculator *)
G4bool IsMaster () const
G4VEnergyLossProcessGetEnergyLossProcess (const G4ParticleDefinition *)
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector ()
const std::vector< G4VEmProcess * > & GetEmProcessVector ()
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector ()
G4EmSaturationEmSaturation ()
G4EmConfiguratorEmConfigurator ()
G4ElectronIonPairElectronIonPair ()
G4NIELCalculatorNIELCalculator ()
G4EmCorrectionsEmCorrections ()
G4VAtomDeexcitationAtomDeexcitation ()
G4VSubCutProducerSubCutProducer ()
G4LossTableBuilderGetTableBuilder ()
void SetGammaGeneralProcess (G4VEmProcess *)
G4VEmProcessGetGammaGeneralProcess ()
void SetElectronGeneralProcess (G4VEmProcess *)
G4VEmProcessGetElectronGeneralProcess ()
void SetPositronGeneralProcess (G4VEmProcess *)
G4VEmProcessGetPositronGeneralProcess ()
 G4LossTableManager (G4LossTableManager &)=delete
G4LossTableManageroperator= (const G4LossTableManager &right)=delete

Static Public Member Functions

static G4LossTableManagerInstance ()

Friends

class G4ThreadLocalSingleton< G4LossTableManager >

Detailed Description

Definition at line 78 of file G4LossTableManager.hh.

Constructor & Destructor Documentation

◆ ~G4LossTableManager()

G4LossTableManager::~G4LossTableManager ( )

Definition at line 101 of file G4LossTableManager.cc.

102{
103 for (auto const & p : loss_vector) { delete p; }
104 for (auto const & p : msc_vector) { delete p; }
105 for (auto const & p : emp_vector) { delete p; }
106 for (auto const & p : p_vector) { delete p; }
107 for (auto const & p : xray_vector) { delete p; }
108
109 std::size_t mod = mod_vector.size();
110 std::size_t fmod = fmod_vector.size();
111 for (std::size_t a=0; a<mod; ++a) {
112 if( nullptr != mod_vector[a] ) {
113 for (std::size_t b=0; b<fmod; ++b) {
114 if((G4VEmModel*)(fmod_vector[b]) == mod_vector[a]) {
115 fmod_vector[b] = nullptr;
116 }
117 }
118 delete mod_vector[a];
119 mod_vector[a] = nullptr;
120 }
121 }
122 for (auto const & p : fmod_vector) { delete p; }
123
124 Clear();
125 delete tableBuilder;
126 delete emCorrections;
127 delete emConfigurator;
128 delete emElectronIonPair;
129 delete nielCalculator;
130 delete atomDeexcitation;
131 delete subcutProducer;
132}

◆ G4LossTableManager()

G4LossTableManager::G4LossTableManager ( G4LossTableManager & )
delete

Member Function Documentation

◆ AtomDeexcitation()

◆ BuildPhysicsTable() [1/2]

void G4LossTableManager::BuildPhysicsTable ( const G4ParticleDefinition * aParticle)

Definition at line 555 of file G4LossTableManager.cc.

556{
557 if(-1 == run && startInitialisation) {
558 if (nullptr != emConfigurator) { emConfigurator->Clear(); }
559 }
560 if (startInitialisation) { resetParam = true; }
561}

◆ BuildPhysicsTable() [2/2]

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

Definition at line 644 of file G4LossTableManager.cc.

647{
648 if(1 < verbose) {
649 G4cout << "### G4LossTableManager::BuildPhysicsTable() for "
650 << aParticle->GetParticleName()
651 << " and process " << p->GetProcessName() << G4endl;
652 }
653 // clear configurator
654 if(-1 == run && startInitialisation) {
655 if( nullptr != emConfigurator) { emConfigurator->Clear(); }
656 firstParticle = aParticle;
657 }
658 if(startInitialisation) {
659 ++run;
660 resetParam = true;
661 startInitialisation = false;
662 if(1 < verbose) {
663 G4cout << "===== G4LossTableManager::BuildPhysicsTable() for run "
664 << run << " ===== " << atomDeexcitation << G4endl;
665 }
666 currentParticle = nullptr;
667 all_tables_are_built = false;
668
669 for (G4int i=0; i<n_loss; ++i) {
670 G4VEnergyLossProcess* el = loss_vector[i];
671
672 if(nullptr != el) {
673 isActive[i] = true;
674 part_vector[i] = el->Particle();
675 base_part_vector[i] = el->BaseParticle();
676 tables_are_built[i] = false;
677 if(1 < verbose) {
678 G4cout << i <<". "<< el->GetProcessName();
679 if(el->Particle()) {
680 G4cout << " for " << el->Particle()->GetParticleName();
681 }
682 G4cout << " active= " << isActive[i]
683 << " table= " << tables_are_built[i]
684 << " isIonisation= " << el->IsIonisationProcess();
685 if(base_part_vector[i]) {
686 G4cout << " base particle "
687 << base_part_vector[i]->GetParticleName();
688 }
689 G4cout << G4endl;
690 }
691 } else {
692 tables_are_built[i] = true;
693 part_vector[i] = nullptr;
694 isActive[i] = false;
695 }
696 }
697 }
698
699 if (all_tables_are_built) {
700 theParameters->SetIsPrintedFlag(true);
701 return;
702 }
703
704 // Build tables for given particle
705 all_tables_are_built = true;
706
707 for(G4int i=0; i<n_loss; ++i) {
708 if(p == loss_vector[i] && !tables_are_built[i] && nullptr == base_part_vector[i]) {
709 const G4ParticleDefinition* curr_part = part_vector[i];
710 if(1 < verbose) {
711 G4cout << "### Build Table for " << p->GetProcessName()
712 << " and " << curr_part->GetParticleName()
713 << " " << tables_are_built[i] << " " << base_part_vector[i]
714 << G4endl;
715 }
716 G4VEnergyLossProcess* curr_proc = BuildTables(curr_part);
717 if(curr_proc) {
718 CopyTables(curr_part, curr_proc);
719 if(p == curr_proc && 0 == run && p->IsIonisationProcess()) {
720 loss_map[aParticle] = p;
721 //G4cout << "G4LossTableManager::BuildPhysicsTable: "
722 // << aParticle->GetParticleName()
723 // << " added to map " << p << G4endl;
724 }
725 }
726 }
727 if ( !tables_are_built[i] ) { all_tables_are_built = false; }
728 }
729 if(1 < verbose) {
730 G4cout << "### G4LossTableManager::BuildPhysicsTable end: "
731 << "all_tables_are_built= " << all_tables_are_built << " "
732 << aParticle->GetParticleName() << " proc: " << p << G4endl;
733 }
734}
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
const G4String & GetParticleName() const
const G4ParticleDefinition * BaseParticle() const
const G4ParticleDefinition * Particle() const
const G4String & GetProcessName() const

◆ DeRegister() [1/7]

void G4LossTableManager::DeRegister ( G4VEmFluctuationModel * p)

Definition at line 387 of file G4LossTableManager.cc.

388{
389 std::size_t n = fmod_vector.size();
390 for (std::size_t i=0; i<n; ++i) {
391 if(fmod_vector[i] == p) {
392 fmod_vector[i] = nullptr;
393 break;
394 }
395 }
396}

◆ DeRegister() [2/7]

void G4LossTableManager::DeRegister ( G4VEmModel * p)

Definition at line 360 of file G4LossTableManager.cc.

361{
362 std::size_t n = mod_vector.size();
363 for (std::size_t i=0; i<n; ++i) {
364 if(mod_vector[i] == p) {
365 mod_vector[i] = nullptr;
366 break;
367 }
368 }
369}

◆ DeRegister() [3/7]

void G4LossTableManager::DeRegister ( G4VEmProcess * p)

Definition at line 302 of file G4LossTableManager.cc.

303{
304 if (nullptr == p) { return; }
305 std::size_t emp = emp_vector.size();
306 for (std::size_t i=0; i<emp; ++i) {
307 if(emp_vector[i] == p) {
308 emp_vector[i] = nullptr;
309 break;
310 }
311 }
312}

◆ DeRegister() [4/7]

void G4LossTableManager::DeRegister ( G4VEnergyLossProcess * p)

Definition at line 243 of file G4LossTableManager.cc.

244{
245 if (nullptr == p) { return; }
246 for (G4int i=0; i<n_loss; ++i) {
247 if(loss_vector[i] == p) {
248 loss_vector[i] = nullptr;
249 break;
250 }
251 }
252}

◆ DeRegister() [5/7]

void G4LossTableManager::DeRegister ( G4VMultipleScattering * p)

Definition at line 272 of file G4LossTableManager.cc.

273{
274 if (nullptr == p) { return; }
275 std::size_t msc = msc_vector.size();
276 for (std::size_t i=0; i<msc; ++i) {
277 if(msc_vector[i] == p) {
278 msc_vector[i] = nullptr;
279 break;
280 }
281 }
282}

◆ DeRegister() [6/7]

void G4LossTableManager::DeRegister ( G4VProcess * p)

Definition at line 332 of file G4LossTableManager.cc.

333{
334 if (nullptr == p) { return; }
335 std::size_t emp = p_vector.size();
336 for (std::size_t i=0; i<emp; ++i) {
337 if(p_vector[i] == p) {
338 p_vector[i] = nullptr;
339 break;
340 }
341 }
342}

◆ DeRegister() [7/7]

void G4LossTableManager::DeRegister ( G4VXRayModel * p)

Definition at line 414 of file G4LossTableManager.cc.

415{
416 std::size_t n = xray_vector.size();
417 for (std::size_t i=0; i<n; ++i) {
418 if (xray_vector[i] == p) {
419 xray_vector[i] = nullptr;
420 break;
421 }
422 }
423}

◆ DumpHtml()

void G4LossTableManager::DumpHtml ( )

Definition at line 1049 of file G4LossTableManager.cc.

1050{
1051 // Automatic generation of html documentation page for physics lists
1052 // List processes and models for the most important
1053 // particles in descending order of importance
1054 // NB. for model names with length > 18 characters the .rst file needs
1055 // to be edited by hand. Or modify G4EmModelManager::DumpModelList
1056
1057 char* dirName = std::getenv("G4PhysListDocDir");
1058 char* physList = std::getenv("G4PhysListName");
1059 if (dirName && physList) {
1060 G4String physListName = G4String(physList);
1061 G4String pathName = G4String(dirName) + "/" + physListName + ".rst";
1062
1063 std::ofstream outFile;
1064 outFile.open(pathName);
1065
1066 outFile << physListName << G4endl;
1067 outFile << std::string(physListName.length(), '=') << G4endl;
1068
1069 std::vector<G4ParticleDefinition*> particles {
1076 };
1077
1078 std::vector<G4VEmProcess*> emproc_vector = GetEmProcessVector();
1079 std::vector<G4VEnergyLossProcess*> enloss_vector =
1081 std::vector<G4VMultipleScattering*> mscat_vector =
1083
1084 for (auto theParticle : particles) {
1085 outFile << G4endl << "**" << theParticle->GetParticleName()
1086 << "**" << G4endl << G4endl << " .. code-block:: none" << G4endl;
1087
1088 G4ProcessManager* pm = theParticle->GetProcessManager();
1089 G4ProcessVector* pv = pm->GetProcessList();
1090 G4int plen = pm->GetProcessListLength();
1091
1092 for (auto emproc : emproc_vector) {
1093 for (G4int i = 0; i < plen; ++i) {
1094 G4VProcess* proc = (*pv)[i];
1095 if (proc == emproc) {
1096 outFile << G4endl;
1097 proc->ProcessDescription(outFile);
1098 break;
1099 }
1100 }
1101 }
1102
1103 for (auto mscproc : mscat_vector) {
1104 for (G4int i = 0; i < plen; ++i) {
1105 G4VProcess* proc = (*pv)[i];
1106 if (proc == mscproc) {
1107 outFile << G4endl;
1108 proc->ProcessDescription(outFile);
1109 break;
1110 }
1111 }
1112 }
1113
1114 for (auto enlossproc : enloss_vector) {
1115 for (G4int i = 0; i < plen; ++i) {
1116 G4VProcess* proc = (*pv)[i];
1117 if (proc == enlossproc) {
1118 outFile << G4endl;
1119 proc->ProcessDescription(outFile);
1120 break;
1121 }
1122 }
1123 }
1124 }
1125 outFile.close();
1126 }
1127}
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
const std::vector< G4VEmProcess * > & GetEmProcessVector()
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector()
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
static G4MuonMinus * MuonMinusDefinition()
static G4MuonPlus * MuonPlusDefinition()
Definition G4MuonPlus.cc:93
static G4Positron * Positron()
Definition G4Positron.cc:90
G4int GetProcessListLength() const
G4ProcessVector * GetProcessList() const
static G4Proton * ProtonDefinition()
Definition G4Proton.cc:85
virtual void ProcessDescription(std::ostream &outfile) const

◆ ElectronIonPair()

G4ElectronIonPair * G4LossTableManager::ElectronIonPair ( )

Definition at line 982 of file G4LossTableManager.cc.

983{
984 if(!emElectronIonPair) {
985 emElectronIonPair = new G4ElectronIonPair(verbose);
986 }
987 return emElectronIonPair;
988}

◆ EmConfigurator()

G4EmConfigurator * G4LossTableManager::EmConfigurator ( )

Definition at line 972 of file G4LossTableManager.cc.

973{
974 if(!emConfigurator) {
975 emConfigurator = new G4EmConfigurator(verbose);
976 }
977 return emConfigurator;
978}

◆ EmCorrections()

◆ EmSaturation()

G4EmSaturation * G4LossTableManager::EmSaturation ( )

Definition at line 965 of file G4LossTableManager.cc.

966{
967 return theParameters->GetEmSaturation();
968}

Referenced by G4OpticalPhysics::ConstructProcess().

◆ GetCSDARange()

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

Definition at line 325 of file G4LossTableManager.hh.

328{
329 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
330 return currentLoss ? currentLoss->GetCSDARange(kineticEnergy, couple) : DBL_MAX;
331}
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
#define DBL_MAX
Definition templates.hh:62

◆ GetDEDX()

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

Definition at line 314 of file G4LossTableManager.hh.

317{
318 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
319 return currentLoss ? currentLoss->GetDEDX(kineticEnergy, couple) : 0.0;
320}

Referenced by G4EnergyLossTables::GetDEDX(), G4EnergyLossTables::GetPreciseDEDX(), G4EnergyLossTables::GetPreciseRangeFromEnergy(), G4Cerenkov::PostStepGetPhysicalInteractionLength(), and G4QuasiCerenkov::PostStepGetPhysicalInteractionLength().

◆ GetDEDXDispersion()

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

Definition at line 370 of file G4LossTableManager.hh.

374{
375 const G4ParticleDefinition* aParticle = dp->GetParticleDefinition();
376 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
377 return currentLoss ? currentLoss->GetDEDXDispersion(couple, dp, length) : 0.0;
378}
const G4ParticleDefinition * GetParticleDefinition() const

◆ GetElectronGeneralProcess()

G4VEmProcess * G4LossTableManager::GetElectronGeneralProcess ( )
inline

Definition at line 438 of file G4LossTableManager.hh.

439{
440 return eGeneral;
441}

◆ GetEmProcessVector()

const std::vector< G4VEmProcess * > & G4LossTableManager::GetEmProcessVector ( )

Definition at line 950 of file G4LossTableManager.cc.

951{
952 return emp_vector;
953}

Referenced by DumpHtml().

◆ GetEnergy()

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

Definition at line 359 of file G4LossTableManager.hh.

362{
363 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
364 return currentLoss ? currentLoss->GetKineticEnergy(range, couple) : 0.0;
365}

Referenced by G4EnergyLossTables::GetPreciseEnergyFromRange().

◆ GetEnergyLossProcess()

G4VEnergyLossProcess * G4LossTableManager::GetEnergyLossProcess ( const G4ParticleDefinition * aParticle)

Definition at line 454 of file G4LossTableManager.cc.

455{
456 if(aParticle != currentParticle) {
457 currentParticle = aParticle;
458 std::map<PD,G4VEnergyLossProcess*,std::less<PD> >::const_iterator pos;
459 if ((pos = loss_map.find(aParticle)) != loss_map.end()) {
460 currentLoss = (*pos).second;
461 } else {
462 currentLoss = nullptr;
463 if(0.0 != aParticle->GetPDGCharge() &&
464 (pos = loss_map.find(theGenericIon)) != loss_map.end()) {
465 currentLoss = (*pos).second;
466 }
467 }
468 }
469 return currentLoss;
470}

Referenced by GetCSDARange(), GetDEDX(), GetDEDXDispersion(), GetEnergy(), GetRange(), and GetRangeFromRestricteDEDX().

◆ GetEnergyLossProcessVector()

const std::vector< G4VEnergyLossProcess * > & G4LossTableManager::GetEnergyLossProcessVector ( )

Definition at line 943 of file G4LossTableManager.cc.

944{
945 return loss_vector;
946}

Referenced by G4EmCalculator::ComputeDEDXForCutInRange(), G4EmCalculator::ComputeElectronicDEDX(), and DumpHtml().

◆ GetGammaGeneralProcess()

G4VEmProcess * G4LossTableManager::GetGammaGeneralProcess ( )
inline

◆ GetMultipleScatteringVector()

const std::vector< G4VMultipleScattering * > & G4LossTableManager::GetMultipleScatteringVector ( )

Definition at line 958 of file G4LossTableManager.cc.

959{
960 return msc_vector;
961}

Referenced by DumpHtml().

◆ GetPositronGeneralProcess()

G4VEmProcess * G4LossTableManager::GetPositronGeneralProcess ( )
inline

Definition at line 452 of file G4LossTableManager.hh.

453{
454 return pGeneral;
455}

◆ GetRange()

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

Definition at line 348 of file G4LossTableManager.hh.

351{
352 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
353 return currentLoss ? currentLoss->GetRange(kineticEnergy, couple) : DBL_MAX;
354}

Referenced by G4ITStepProcessor::ApplyProductionCut(), G4EnergyLossTables::GetRange(), G4Cerenkov::PostStepGetPhysicalInteractionLength(), G4QuasiCerenkov::PostStepGetPhysicalInteractionLength(), and G4EmSaturation::VisibleEnergyDeposition().

◆ GetRangeFromRestricteDEDX()

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

Definition at line 336 of file G4LossTableManager.hh.

340{
341 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
342 return currentLoss ? currentLoss->GetRange(kineticEnergy, couple) : DBL_MAX;
343}

◆ GetTableBuilder()

G4LossTableBuilder * G4LossTableManager::GetTableBuilder ( )
inline

◆ Instance()

G4LossTableManager * G4LossTableManager::Instance ( )
static

Definition at line 90 of file G4LossTableManager.cc.

91{
92 if(nullptr == instance) {
94 instance = inst.Instance();
95 }
96 return instance;
97}
friend class G4ThreadLocalSingleton< G4LossTableManager >

Referenced by G4ITStepProcessor::ApplyProductionCut(), G4BertiniElectroNuclearBuilder::Build(), G4EmTableUtil::BuildMscProcess(), G4GammaGeneralProcess::BuildPhysicsTable(), G4EmCalculator::ComputeDEDXForCutInRange(), G4EmCalculator::ComputeElectronicDEDX(), LBE::ConstructGeneral(), G4EmExtraPhysics::ConstructProcess(), G4EmStandardPhysics::ConstructProcess(), G4EmStandardPhysics_option1::ConstructProcess(), G4EmStandardPhysics_option2::ConstructProcess(), G4EmStandardPhysics_option3::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4EmStandardPhysicsGS::ConstructProcess(), G4EmStandardPhysicsSS::ConstructProcess(), G4HadronPhysicsLEND::ConstructProcess(), G4OpticalPhysics::ConstructProcess(), G4RadioactiveDecayPhysics::ConstructProcess(), G4ECDecay::DecayIt(), G4ITDecay::DecayIt(), G4AnnihiToMuPair::G4AnnihiToMuPair(), G4AtimaEnergyLossModel::G4AtimaEnergyLossModel(), G4BetheBlochModel::G4BetheBlochModel(), G4BraggModel::G4BraggModel(), G4DNARuddIonisationDynamicModel::G4DNARuddIonisationDynamicModel(), G4DNARuddIonisationExtendedModel::G4DNARuddIonisationExtendedModel(), G4DynamicParticleIonisation::G4DynamicParticleIonisation(), G4DynamicParticleMSC::G4DynamicParticleMSC(), G4EmCalculator::G4EmCalculator(), G4GammaConversionToMuons::G4GammaConversionToMuons(), G4IonParametrisedLossModel::G4IonParametrisedLossModel(), G4LindhardSorensenIonModel::G4LindhardSorensenIonModel(), G4MuBetheBlochModel::G4MuBetheBlochModel(), G4NIELCalculator::G4NIELCalculator(), G4SynchrotronRadiation::G4SynchrotronRadiation(), G4TransportationWithMsc::G4TransportationWithMsc(), G4UAtomicDeexcitation::G4UAtomicDeexcitation(), G4UrbanAdjointMscModel::G4UrbanAdjointMscModel(), G4UserSpecialCuts::G4UserSpecialCuts(), G4VEmFluctuationModel::G4VEmFluctuationModel(), G4VEmModel::G4VEmModel(), G4VEmProcess::G4VEmProcess(), G4VEnergyLossProcess::G4VEnergyLossProcess(), G4VMultipleScattering::G4VMultipleScattering(), G4VTransitionRadiation::G4VTransitionRadiation(), G4EnergyLossTables::GetDEDX(), G4EnergyLossTables::GetPreciseDEDX(), G4EnergyLossTables::GetPreciseEnergyFromRange(), G4EnergyLossTables::GetPreciseRangeFromEnergy(), G4EnergyLossTables::GetRange(), G4DNABornIonisationModel1::Initialise(), G4DNABornIonisationModel2::Initialise(), G4DNABornIonisationModel::Initialise(), G4DNACPA100IonisationModel::Initialise(), G4DNADoubleIonisationModel::Initialise(), G4DNAEmfietzoglouIonisationModel::Initialise(), G4DNAQuadrupleIonisationModel::Initialise(), G4DNARelativisticIonisationModel::Initialise(), G4DNARPWBAIonisationModel::Initialise(), G4DNARuddIonisationDynamicModel::Initialise(), G4DNARuddIonisationExtendedModel::Initialise(), G4DNARuddIonisationModel::Initialise(), G4DNATripleIonisationModel::Initialise(), G4KleinNishinaModel::Initialise(), G4LivermoreComptonModel::Initialise(), G4LivermorePhotoElectricModel::Initialise(), G4LivermorePolarizedComptonModel::Initialise(), G4LowEPComptonModel::Initialise(), G4LowEPPolarizedComptonModel::Initialise(), G4MicroElecInelasticModel::Initialise(), G4MicroElecInelasticModel_new::Initialise(), G4PEEffectFluoModel::Initialise(), G4PenelopeComptonModel::Initialise(), G4PenelopeIonisationModel::Initialise(), G4PenelopePhotoElectricModel::Initialise(), G4GammaGeneralProcess::InitialiseProcess(), G4Cerenkov::PostStepGetPhysicalInteractionLength(), G4QuasiCerenkov::PostStepGetPhysicalInteractionLength(), G4EmBuilder::PrepareEMPhysics(), G4GammaGeneralProcess::PreparePhysicsTable(), and G4EmSaturation::VisibleEnergyDeposition().

◆ IsMaster()

G4bool G4LossTableManager::IsMaster ( ) const
inline

Definition at line 382 of file G4LossTableManager.hh.

383{
384 return isMaster;
385}

◆ LocalPhysicsTables()

void G4LossTableManager::LocalPhysicsTables ( const G4ParticleDefinition * aParticle,
G4VEnergyLossProcess * p )

Definition at line 565 of file G4LossTableManager.cc.

568{
569 if (1 < verbose) {
570 G4cout << "### G4LossTableManager::LocalPhysicsTable() for "
571 << aParticle->GetParticleName()
572 << " and process " << p->GetProcessName()
573 << G4endl;
574 }
575
576 if(-1 == run && startInitialisation) {
577 if (nullptr != emConfigurator) { emConfigurator->Clear(); }
578 firstParticle = aParticle;
579 }
580
581 if (startInitialisation) {
582 ++run;
583 if (1 < verbose) {
584 G4cout << "===== G4LossTableManager::LocalPhysicsTable() for run "
585 << run << " =====" << G4endl;
586 }
587 currentParticle = nullptr;
588 startInitialisation = false;
589 resetParam = true;
590 for (G4int i=0; i<n_loss; ++i) {
591 if (nullptr != loss_vector[i]) {
592 tables_are_built[i] = false;
593 } else {
594 tables_are_built[i] = true;
595 part_vector[i] = nullptr;
596 }
597 }
598 }
599
600 all_tables_are_built= true;
601 for (G4int i=0; i<n_loss; ++i) {
602 if(p == loss_vector[i]) {
603 tables_are_built[i] = true;
604 isActive[i] = true;
605 part_vector[i] = p->Particle();
606 base_part_vector[i] = p->BaseParticle();
607 dedx_vector[i] = p->DEDXTable();
608 range_vector[i] = p->RangeTableForLoss();
609 inv_range_vector[i] = p->InverseRangeTable();
610 if(0 == run && p->IsIonisationProcess()) {
611 loss_map[part_vector[i]] = p;
612 }
613
614 if(1 < verbose) {
615 G4cout << i <<". "<< p->GetProcessName();
616 if(part_vector[i]) {
617 G4cout << " for " << part_vector[i]->GetParticleName();
618 }
619 G4cout << " active= " << isActive[i]
620 << " table= " << tables_are_built[i]
621 << " isIonisation= " << p->IsIonisationProcess()
622 << G4endl;
623 }
624 break;
625 } else if(!tables_are_built[i]) {
626 all_tables_are_built = false;
627 }
628 }
629
630 if(1 < verbose) {
631 G4cout << "### G4LossTableManager::LocalPhysicsTable end"
632 << G4endl;
633 }
634 if(all_tables_are_built) {
635 if(1 < verbose) {
636 G4cout << "%%%%% All dEdx and Range tables for worker are ready for run "
637 << run << " %%%%%" << G4endl;
638 }
639 }
640}
G4PhysicsTable * RangeTableForLoss() const
G4PhysicsTable * InverseRangeTable() const
G4PhysicsTable * DEDXTable() const

◆ NIELCalculator()

G4NIELCalculator * G4LossTableManager::NIELCalculator ( )

Definition at line 1002 of file G4LossTableManager.cc.

1003{
1004 if(!nielCalculator) {
1005 nielCalculator = new G4NIELCalculator(nullptr, verbose);
1006 }
1007 return nielCalculator;
1008}

◆ operator=()

G4LossTableManager & G4LossTableManager::operator= ( const G4LossTableManager & right)
delete

◆ PreparePhysicsTable() [1/3]

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

Definition at line 509 of file G4LossTableManager.cc.

511{
512 if (1 < verbose) {
513 G4cout << "G4LossTableManager::PreparePhysicsTable for "
514 << particle->GetParticleName()
515 << " and " << p->GetProcessName()
516 << " run=" << run << " master=" << isMaster
517 << G4endl;
518 }
519
520 // start initialisation for the first run
521 if( -1 == run ) {
522 if (nullptr != emConfigurator) { emConfigurator->PrepareModels(particle, p); }
523 }
524
526}

◆ PreparePhysicsTable() [2/3]

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

Definition at line 475 of file G4LossTableManager.cc.

477{
478 if (1 < verbose) {
479 G4cout << "G4LossTableManager::PreparePhysicsTable for "
480 << particle->GetParticleName()
481 << " and " << p->GetProcessName() << " run= " << run
482 << " loss_vector " << loss_vector.size()
483 << " run=" << run << " master=" << isMaster
484 << G4endl;
485 }
486
487 // start initialisation for the first run
488 if ( -1 == run ) {
489 if (nullptr != emConfigurator) {
490 emConfigurator->PrepareModels(particle, p);
491 }
492
493 // initialise particles for given process
494 for (G4int j=0; j<n_loss; ++j) {
495 if (p == loss_vector[j] && nullptr == part_vector[j]) {
496 part_vector[j] = particle;
497 if (particle->GetParticleName() == "GenericIon") {
498 theGenericIon = particle;
499 }
500 }
501 }
502 }
504}

◆ PreparePhysicsTable() [3/3]

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

Definition at line 531 of file G4LossTableManager.cc.

533{
534 if (1 < verbose) {
535 G4cout << "G4LossTableManager::PreparePhysicsTable for "
536 << particle->GetParticleName()
537 << " and " << p->GetProcessName()
538 << " run=" << run << " master=" << isMaster
539 << G4endl;
540 }
541
542 // start initialisation for the first run
543 if ( -1 == run ) {
544 if (nullptr != emConfigurator) {
545 emConfigurator->PrepareModels(particle, p);
546 }
547 }
548
550}

◆ Register() [1/7]

void G4LossTableManager::Register ( G4VEmFluctuationModel * p)

Definition at line 373 of file G4LossTableManager.cc.

374{
375 if (nullptr == p) { return; }
376 std::size_t n = fmod_vector.size();
377 for (std::size_t i=0; i<n; ++i) { if (fmod_vector[i] == p) { return; } }
378 fmod_vector.push_back(p);
379 if(verbose > 1) {
380 G4cout << "G4LossTableManager::Register G4VEmFluctuationModel : "
381 << p->GetName() << " " << fmod_vector.size() << G4endl;
382 }
383}
const G4String & GetName() const

◆ Register() [2/7]

void G4LossTableManager::Register ( G4VEmModel * p)

Definition at line 346 of file G4LossTableManager.cc.

347{
348 if (nullptr == p) { return; }
349 std::size_t n = mod_vector.size();
350 for (std::size_t i=0; i<n; ++i) { if (mod_vector[i] == p) { return; } }
351 mod_vector.push_back(p);
352 if(verbose > 1) {
353 G4cout << "G4LossTableManager::Register G4VEmModel : "
354 << p->GetName() << " " << p << " " << mod_vector.size() << G4endl;
355 }
356}
const G4String & GetName() const

◆ Register() [3/7]

void G4LossTableManager::Register ( G4VEmProcess * p)

Definition at line 286 of file G4LossTableManager.cc.

287{
288 if (nullptr == p) { return; }
289 std::size_t n = emp_vector.size();
290 for (std::size_t 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/7]

void G4LossTableManager::Register ( G4VEnergyLossProcess * p)

Definition at line 186 of file G4LossTableManager.cc.

187{
188 if (nullptr == p) { return; }
189 for (G4int i=0; i<n_loss; ++i) {
190 if(loss_vector[i] == p) { return; }
191 }
192 if(verbose > 1) {
193 G4cout << "G4LossTableManager::Register G4VEnergyLossProcess : "
194 << p->GetProcessName() << " idx= " << n_loss << G4endl;
195 }
196 ++n_loss;
197 loss_vector.push_back(p);
198 part_vector.push_back(nullptr);
199 base_part_vector.push_back(nullptr);
200 dedx_vector.push_back(nullptr);
201 range_vector.push_back(nullptr);
202 inv_range_vector.push_back(nullptr);
203 tables_are_built.push_back(false);
204 isActive.push_back(true);
205 all_tables_are_built = false;
206}

◆ Register() [5/7]

void G4LossTableManager::Register ( G4VMultipleScattering * p)

Definition at line 256 of file G4LossTableManager.cc.

257{
258 if (nullptr == p) { return; }
259 std::size_t n = msc_vector.size();
260 for (std::size_t i=0; i<n; ++i) {
261 if(msc_vector[i] == p) { return; }
262 }
263 if(verbose > 1) {
264 G4cout << "G4LossTableManager::Register G4VMultipleScattering : "
265 << p->GetProcessName() << " idx= " << msc_vector.size() << G4endl;
266 }
267 msc_vector.push_back(p);
268}

◆ Register() [6/7]

void G4LossTableManager::Register ( G4VProcess * p)

Definition at line 316 of file G4LossTableManager.cc.

317{
318 if (nullptr == p) { return; }
319 std::size_t n = p_vector.size();
320 for (std::size_t i=0; i<n; ++i) {
321 if(p_vector[i] == p) { return; }
322 }
323 if(verbose > 1) {
324 G4cout << "G4LossTableManager::Register G4VProcess : "
325 << p->GetProcessName() << " idx= " << p_vector.size() << G4endl;
326 }
327 p_vector.push_back(p);
328}

◆ Register() [7/7]

void G4LossTableManager::Register ( G4VXRayModel * p)

Definition at line 400 of file G4LossTableManager.cc.

401{
402 if (nullptr == p) { return; }
403 std::size_t n = xray_vector.size();
404 for (std::size_t i=0; i<n; ++i) { if (xray_vector[i] == p) { return; } }
405 xray_vector.push_back(p);
406 if (verbose > 1) {
407 G4cout << "G4LossTableManager::Register G4VXRayModel : "
408 << p->GetName() << " " << xray_vector.size() << G4endl;
409 }
410}
const G4String & GetName() const

◆ RegisterExtraParticle()

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

Definition at line 427 of file G4LossTableManager.cc.

430{
431 if (nullptr == p || nullptr == part) { return; }
432 for (G4int i=0; i<n_loss; ++i) {
433 if(loss_vector[i] == p) { return; }
434 }
435 if(verbose > 1) {
436 G4cout << "G4LossTableManager::RegisterExtraParticle "
437 << part->GetParticleName() << " G4VEnergyLossProcess : "
438 << p->GetProcessName() << " idx= " << n_loss << G4endl;
439 }
440 ++n_loss;
441 loss_vector.push_back(p);
442 part_vector.push_back(part);
443 base_part_vector.push_back(p->BaseParticle());
444 dedx_vector.push_back(nullptr);
445 range_vector.push_back(nullptr);
446 inv_range_vector.push_back(nullptr);
447 tables_are_built.push_back(false);
448 all_tables_are_built = false;
449}

◆ ResetParameters()

void G4LossTableManager::ResetParameters ( )

Definition at line 210 of file G4LossTableManager.cc.

211{
212 // initialisation once per run
213 if (!resetParam) { return; }
214 resetParam = false;
215 startInitialisation = true;
216 verbose = theParameters->Verbose();
217 if(!isMaster) {
218 verbose = theParameters->WorkerVerbose();
219 } else {
220 if(verbose > 0) { theParameters->Dump(); }
221 }
222
223 tableBuilder->InitialiseBaseMaterials();
224 if (nullptr != nielCalculator) { nielCalculator->Initialise(); }
225
226 emCorrections->SetVerbose(verbose);
227 if(nullptr != emConfigurator) { emConfigurator->SetVerbose(verbose); };
228 if(nullptr != emElectronIonPair) { emElectronIonPair->SetVerbose(verbose); };
229 if(nullptr != atomDeexcitation) {
230 atomDeexcitation->SetVerboseLevel(verbose);
231 atomDeexcitation->InitialiseAtomicDeexcitation();
232 }
233 if (1 < verbose) {
234 G4cout << "====== G4LossTableManager::ResetParameters "
235 << " Nloss=" << loss_vector.size()
236 << " run=" << run << " master=" << isMaster
237 << G4endl;
238 }
239}

Referenced by G4RadioactiveDecayPhysics::ConstructProcess(), PreparePhysicsTable(), PreparePhysicsTable(), and PreparePhysicsTable().

◆ SetAtomDeexcitation()

void G4LossTableManager::SetAtomDeexcitation ( G4VAtomDeexcitation * p)

Definition at line 1012 of file G4LossTableManager.cc.

1013{
1014 if(atomDeexcitation != p) {
1015 delete atomDeexcitation;
1016 atomDeexcitation = p;
1017 }
1018}

Referenced by LBE::ConstructGeneral(), G4RadioactiveDecayPhysics::ConstructProcess(), and G4EmBuilder::PrepareEMPhysics().

◆ SetElectronGeneralProcess()

void G4LossTableManager::SetElectronGeneralProcess ( G4VEmProcess * ptr)
inline

Definition at line 431 of file G4LossTableManager.hh.

432{
433 eGeneral = ptr;
434}

◆ SetGammaGeneralProcess()

◆ SetNIELCalculator()

void G4LossTableManager::SetNIELCalculator ( G4NIELCalculator * ptr)

Definition at line 992 of file G4LossTableManager.cc.

993{
994 if(nullptr != ptr && ptr != nielCalculator) {
995 delete nielCalculator;
996 nielCalculator = ptr;
997 }
998}

Referenced by G4NIELCalculator::G4NIELCalculator().

◆ SetPositronGeneralProcess()

void G4LossTableManager::SetPositronGeneralProcess ( G4VEmProcess * ptr)
inline

Definition at line 445 of file G4LossTableManager.hh.

446{
447 pGeneral = ptr;
448}

◆ SetSubCutProducer()

void G4LossTableManager::SetSubCutProducer ( G4VSubCutProducer * p)

Definition at line 1022 of file G4LossTableManager.cc.

1023{
1024 if(subcutProducer != p) {
1025 delete subcutProducer;
1026 subcutProducer = p;
1027 }
1028}

◆ SetVerbose()

void G4LossTableManager::SetVerbose ( G4int val)

Definition at line 935 of file G4LossTableManager.cc.

936{
937 verbose = val;
938}

◆ SubCutProducer()

G4VSubCutProducer * G4LossTableManager::SubCutProducer ( )
inline

Definition at line 403 of file G4LossTableManager.hh.

404{
405 return subcutProducer;
406}

◆ G4ThreadLocalSingleton< G4LossTableManager >

Definition at line 452 of file G4LossTableManager.hh.


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