Geant4 11.3.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 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 77 of file G4LossTableManager.hh.

Constructor & Destructor Documentation

◆ ~G4LossTableManager()

G4LossTableManager::~G4LossTableManager ( )

Definition at line 99 of file G4LossTableManager.cc.

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

◆ G4LossTableManager()

G4LossTableManager::G4LossTableManager ( G4LossTableManager & )
delete

Member Function Documentation

◆ AtomDeexcitation()

◆ BuildPhysicsTable() [1/2]

void G4LossTableManager::BuildPhysicsTable ( const G4ParticleDefinition * aParticle)

Definition at line 513 of file G4LossTableManager.cc.

514{
515 if(-1 == run && startInitialisation) {
516 if (nullptr != emConfigurator) { emConfigurator->Clear(); }
517 }
518 if (startInitialisation) { resetParam = true; }
519}

◆ BuildPhysicsTable() [2/2]

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

Definition at line 602 of file G4LossTableManager.cc.

605{
606 if(1 < verbose) {
607 G4cout << "### G4LossTableManager::BuildPhysicsTable() for "
608 << aParticle->GetParticleName()
609 << " and process " << p->GetProcessName() << G4endl;
610 }
611 // clear configurator
612 if(-1 == run && startInitialisation) {
613 if( nullptr != emConfigurator) { emConfigurator->Clear(); }
614 firstParticle = aParticle;
615 }
616 if(startInitialisation) {
617 ++run;
618 resetParam = true;
619 startInitialisation = false;
620 if(1 < verbose) {
621 G4cout << "===== G4LossTableManager::BuildPhysicsTable() for run "
622 << run << " ===== " << atomDeexcitation << G4endl;
623 }
624 currentParticle = nullptr;
625 all_tables_are_built = false;
626
627 for (G4int i=0; i<n_loss; ++i) {
628 G4VEnergyLossProcess* el = loss_vector[i];
629
630 if(nullptr != el) {
631 isActive[i] = true;
632 part_vector[i] = el->Particle();
633 base_part_vector[i] = el->BaseParticle();
634 tables_are_built[i] = false;
635 if(1 < verbose) {
636 G4cout << i <<". "<< el->GetProcessName();
637 if(el->Particle()) {
638 G4cout << " for " << el->Particle()->GetParticleName();
639 }
640 G4cout << " active= " << isActive[i]
641 << " table= " << tables_are_built[i]
642 << " isIonisation= " << el->IsIonisationProcess();
643 if(base_part_vector[i]) {
644 G4cout << " base particle "
645 << base_part_vector[i]->GetParticleName();
646 }
647 G4cout << G4endl;
648 }
649 } else {
650 tables_are_built[i] = true;
651 part_vector[i] = nullptr;
652 isActive[i] = false;
653 }
654 }
655 }
656
657 if (all_tables_are_built) {
658 theParameters->SetIsPrintedFlag(true);
659 return;
660 }
661
662 // Build tables for given particle
663 all_tables_are_built = true;
664
665 for(G4int i=0; i<n_loss; ++i) {
666 if(p == loss_vector[i] && !tables_are_built[i] && nullptr == base_part_vector[i]) {
667 const G4ParticleDefinition* curr_part = part_vector[i];
668 if(1 < verbose) {
669 G4cout << "### Build Table for " << p->GetProcessName()
670 << " and " << curr_part->GetParticleName()
671 << " " << tables_are_built[i] << " " << base_part_vector[i]
672 << G4endl;
673 }
674 G4VEnergyLossProcess* curr_proc = BuildTables(curr_part);
675 if(curr_proc) {
676 CopyTables(curr_part, curr_proc);
677 if(p == curr_proc && 0 == run && p->IsIonisationProcess()) {
678 loss_map[aParticle] = p;
679 //G4cout << "G4LossTableManager::BuildPhysicsTable: "
680 // << aParticle->GetParticleName()
681 // << " added to map " << p << G4endl;
682 }
683 }
684 }
685 if ( !tables_are_built[i] ) { all_tables_are_built = false; }
686 }
687 if(1 < verbose) {
688 G4cout << "### G4LossTableManager::BuildPhysicsTable end: "
689 << "all_tables_are_built= " << all_tables_are_built << " "
690 << aParticle->GetParticleName() << " proc: " << p << G4endl;
691 }
692}
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/6]

void G4LossTableManager::DeRegister ( G4VEmFluctuationModel * p)

Definition at line 379 of file G4LossTableManager.cc.

380{
381 std::size_t n = fmod_vector.size();
382 for (std::size_t i=0; i<n; ++i) {
383 if(fmod_vector[i] == p) { fmod_vector[i] = nullptr; }
384 }
385}

◆ DeRegister() [2/6]

void G4LossTableManager::DeRegister ( G4VEmModel * p)

Definition at line 354 of file G4LossTableManager.cc.

355{
356 //G4cout << "G4LossTableManager::DeRegister G4VEmModel : " << p << G4endl;
357 std::size_t n = mod_vector.size();
358 for (std::size_t i=0; i<n; ++i) {
359 if(mod_vector[i] == p) {
360 mod_vector[i] = nullptr;
361 break;
362 }
363 }
364}

◆ DeRegister() [3/6]

void G4LossTableManager::DeRegister ( G4VEmProcess * p)

Definition at line 299 of file G4LossTableManager.cc.

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

◆ DeRegister() [4/6]

void G4LossTableManager::DeRegister ( G4VEnergyLossProcess * p)

Definition at line 240 of file G4LossTableManager.cc.

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

◆ DeRegister() [5/6]

void G4LossTableManager::DeRegister ( G4VMultipleScattering * p)

Definition at line 269 of file G4LossTableManager.cc.

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

◆ DeRegister() [6/6]

void G4LossTableManager::DeRegister ( G4VProcess * p)

Definition at line 329 of file G4LossTableManager.cc.

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

◆ DumpHtml()

void G4LossTableManager::DumpHtml ( )

Definition at line 1007 of file G4LossTableManager.cc.

1008{
1009 // Automatic generation of html documentation page for physics lists
1010 // List processes and models for the most important
1011 // particles in descending order of importance
1012 // NB. for model names with length > 18 characters the .rst file needs
1013 // to be edited by hand. Or modify G4EmModelManager::DumpModelList
1014
1015 char* dirName = std::getenv("G4PhysListDocDir");
1016 char* physList = std::getenv("G4PhysListName");
1017 if (dirName && physList) {
1018 G4String physListName = G4String(physList);
1019 G4String pathName = G4String(dirName) + "/" + physListName + ".rst";
1020
1021 std::ofstream outFile;
1022 outFile.open(pathName);
1023
1024 outFile << physListName << G4endl;
1025 outFile << std::string(physListName.length(), '=') << G4endl;
1026
1027 std::vector<G4ParticleDefinition*> particles {
1034 };
1035
1036 std::vector<G4VEmProcess*> emproc_vector = GetEmProcessVector();
1037 std::vector<G4VEnergyLossProcess*> enloss_vector =
1039 std::vector<G4VMultipleScattering*> mscat_vector =
1041
1042 for (auto theParticle : particles) {
1043 outFile << G4endl << "**" << theParticle->GetParticleName()
1044 << "**" << G4endl << G4endl << " .. code-block:: none" << G4endl;
1045
1046 G4ProcessManager* pm = theParticle->GetProcessManager();
1047 G4ProcessVector* pv = pm->GetProcessList();
1048 G4int plen = pm->GetProcessListLength();
1049
1050 for (auto emproc : emproc_vector) {
1051 for (G4int i = 0; i < plen; ++i) {
1052 G4VProcess* proc = (*pv)[i];
1053 if (proc == emproc) {
1054 outFile << G4endl;
1055 proc->ProcessDescription(outFile);
1056 break;
1057 }
1058 }
1059 }
1060
1061 for (auto mscproc : mscat_vector) {
1062 for (G4int i = 0; i < plen; ++i) {
1063 G4VProcess* proc = (*pv)[i];
1064 if (proc == mscproc) {
1065 outFile << G4endl;
1066 proc->ProcessDescription(outFile);
1067 break;
1068 }
1069 }
1070 }
1071
1072 for (auto enlossproc : enloss_vector) {
1073 for (G4int i = 0; i < plen; ++i) {
1074 G4VProcess* proc = (*pv)[i];
1075 if (proc == enlossproc) {
1076 outFile << G4endl;
1077 proc->ProcessDescription(outFile);
1078 break;
1079 }
1080 }
1081 }
1082 }
1083 outFile.close();
1084 }
1085}
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 940 of file G4LossTableManager.cc.

941{
942 if(!emElectronIonPair) {
943 emElectronIonPair = new G4ElectronIonPair(verbose);
944 }
945 return emElectronIonPair;
946}

◆ EmConfigurator()

G4EmConfigurator * G4LossTableManager::EmConfigurator ( )

Definition at line 930 of file G4LossTableManager.cc.

931{
932 if(!emConfigurator) {
933 emConfigurator = new G4EmConfigurator(verbose);
934 }
935 return emConfigurator;
936}

◆ EmCorrections()

◆ EmSaturation()

G4EmSaturation * G4LossTableManager::EmSaturation ( )

Definition at line 923 of file G4LossTableManager.cc.

924{
925 return theParameters->GetEmSaturation();
926}

Referenced by G4OpticalPhysics::ConstructProcess().

◆ GetCSDARange()

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

Definition at line 319 of file G4LossTableManager.hh.

322{
323 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
324 return currentLoss ? currentLoss->GetCSDARange(kineticEnergy, couple) : DBL_MAX;
325}
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 308 of file G4LossTableManager.hh.

311{
312 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
313 return currentLoss ? currentLoss->GetDEDX(kineticEnergy, couple) : 0.0;
314}

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

◆ GetDEDXDispersion()

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

Definition at line 364 of file G4LossTableManager.hh.

368{
369 const G4ParticleDefinition* aParticle = dp->GetParticleDefinition();
370 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
371 return currentLoss ? currentLoss->GetDEDXDispersion(couple, dp, length) : 0.0;
372}
const G4ParticleDefinition * GetParticleDefinition() const

◆ GetElectronGeneralProcess()

G4VEmProcess * G4LossTableManager::GetElectronGeneralProcess ( )
inline

Definition at line 432 of file G4LossTableManager.hh.

433{
434 return eGeneral;
435}

◆ GetEmProcessVector()

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

Definition at line 908 of file G4LossTableManager.cc.

909{
910 return emp_vector;
911}

Referenced by DumpHtml().

◆ GetEnergy()

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

Definition at line 353 of file G4LossTableManager.hh.

356{
357 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
358 return currentLoss ? currentLoss->GetKineticEnergy(range, couple) : 0.0;
359}

Referenced by G4EnergyLossTables::GetPreciseEnergyFromRange().

◆ GetEnergyLossProcess()

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

Definition at line 416 of file G4LossTableManager.cc.

417{
418 if(aParticle != currentParticle) {
419 currentParticle = aParticle;
420 std::map<PD,G4VEnergyLossProcess*,std::less<PD> >::const_iterator pos;
421 if ((pos = loss_map.find(aParticle)) != loss_map.end()) {
422 currentLoss = (*pos).second;
423 } else {
424 currentLoss = nullptr;
425 if(0.0 != aParticle->GetPDGCharge() &&
426 (pos = loss_map.find(theGenericIon)) != loss_map.end()) {
427 currentLoss = (*pos).second;
428 }
429 }
430 }
431 return currentLoss;
432}

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

◆ GetEnergyLossProcessVector()

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

Definition at line 901 of file G4LossTableManager.cc.

902{
903 return loss_vector;
904}

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

◆ GetGammaGeneralProcess()

G4VEmProcess * G4LossTableManager::GetGammaGeneralProcess ( )
inline

Definition at line 418 of file G4LossTableManager.hh.

419{
420 return gGeneral;
421}

Referenced by G4BertiniElectroNuclearBuilder::Build(), and G4EmExtraPhysics::ConstructProcess().

◆ GetMultipleScatteringVector()

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

Definition at line 916 of file G4LossTableManager.cc.

917{
918 return msc_vector;
919}

Referenced by DumpHtml().

◆ GetPositronGeneralProcess()

G4VEmProcess * G4LossTableManager::GetPositronGeneralProcess ( )
inline

Definition at line 446 of file G4LossTableManager.hh.

447{
448 return pGeneral;
449}

◆ GetRange()

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

Definition at line 342 of file G4LossTableManager.hh.

345{
346 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
347 return currentLoss ? currentLoss->GetRange(kineticEnergy, couple) : DBL_MAX;
348}

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

◆ GetRangeFromRestricteDEDX()

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

Definition at line 330 of file G4LossTableManager.hh.

334{
335 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
336 return currentLoss ? currentLoss->GetRange(kineticEnergy, couple) : DBL_MAX;
337}

◆ GetTableBuilder()

G4LossTableBuilder * G4LossTableManager::GetTableBuilder ( )
inline

◆ Instance()

G4LossTableManager * G4LossTableManager::Instance ( )
static

Definition at line 88 of file G4LossTableManager.cc.

89{
90 if(nullptr == instance) {
92 instance = inst.Instance();
93 }
94 return instance;
95}
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(), G4EmStandardPhysicsSS::ConstructProcess(), G4OpticalPhysics::ConstructProcess(), G4RadioactiveDecayPhysics::ConstructProcess(), G4ECDecay::DecayIt(), G4ITDecay::DecayIt(), G4AnnihiToMuPair::G4AnnihiToMuPair(), G4AtimaEnergyLossModel::G4AtimaEnergyLossModel(), G4BetheBlochModel::G4BetheBlochModel(), G4BraggModel::G4BraggModel(), G4DNARuddIonisationExtendedModel::G4DNARuddIonisationExtendedModel(), G4DynamicParticleIonisation::G4DynamicParticleIonisation(), G4DynamicParticleMSC::G4DynamicParticleMSC(), G4EmCalculator::G4EmCalculator(), G4GammaConversionToMuons::G4GammaConversionToMuons(), G4GoudsmitSaundersonMscModel::G4GoudsmitSaundersonMscModel(), 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(), G4DNACPA100IonisationModel::Initialise(), G4DNADoubleIonisationModel::Initialise(), G4DNAEmfietzoglouIonisationModel::Initialise(), G4DNAQuadrupleIonisationModel::Initialise(), G4DNARelativisticIonisationModel::Initialise(), G4DNARPWBAIonisationModel::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(), G4EmBuilder::PrepareEMPhysics(), G4GammaGeneralProcess::PreparePhysicsTable(), and G4EmSaturation::VisibleEnergyDeposition().

◆ IsMaster()

G4bool G4LossTableManager::IsMaster ( ) const
inline

Definition at line 376 of file G4LossTableManager.hh.

377{
378 return isMaster;
379}

◆ LocalPhysicsTables()

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

Definition at line 523 of file G4LossTableManager.cc.

526{
527 if (1 < verbose) {
528 G4cout << "### G4LossTableManager::LocalPhysicsTable() for "
529 << aParticle->GetParticleName()
530 << " and process " << p->GetProcessName()
531 << G4endl;
532 }
533
534 if(-1 == run && startInitialisation) {
535 if (nullptr != emConfigurator) { emConfigurator->Clear(); }
536 firstParticle = aParticle;
537 }
538
539 if (startInitialisation) {
540 ++run;
541 if (1 < verbose) {
542 G4cout << "===== G4LossTableManager::LocalPhysicsTable() for run "
543 << run << " =====" << G4endl;
544 }
545 currentParticle = nullptr;
546 startInitialisation = false;
547 resetParam = true;
548 for (G4int i=0; i<n_loss; ++i) {
549 if (nullptr != loss_vector[i]) {
550 tables_are_built[i] = false;
551 } else {
552 tables_are_built[i] = true;
553 part_vector[i] = nullptr;
554 }
555 }
556 }
557
558 all_tables_are_built= true;
559 for (G4int i=0; i<n_loss; ++i) {
560 if(p == loss_vector[i]) {
561 tables_are_built[i] = true;
562 isActive[i] = true;
563 part_vector[i] = p->Particle();
564 base_part_vector[i] = p->BaseParticle();
565 dedx_vector[i] = p->DEDXTable();
566 range_vector[i] = p->RangeTableForLoss();
567 inv_range_vector[i] = p->InverseRangeTable();
568 if(0 == run && p->IsIonisationProcess()) {
569 loss_map[part_vector[i]] = p;
570 }
571
572 if(1 < verbose) {
573 G4cout << i <<". "<< p->GetProcessName();
574 if(part_vector[i]) {
575 G4cout << " for " << part_vector[i]->GetParticleName();
576 }
577 G4cout << " active= " << isActive[i]
578 << " table= " << tables_are_built[i]
579 << " isIonisation= " << p->IsIonisationProcess()
580 << G4endl;
581 }
582 break;
583 } else if(!tables_are_built[i]) {
584 all_tables_are_built = false;
585 }
586 }
587
588 if(1 < verbose) {
589 G4cout << "### G4LossTableManager::LocalPhysicsTable end"
590 << G4endl;
591 }
592 if(all_tables_are_built) {
593 if(1 < verbose) {
594 G4cout << "%%%%% All dEdx and Range tables for worker are ready for run "
595 << run << " %%%%%" << G4endl;
596 }
597 }
598}
G4PhysicsTable * RangeTableForLoss() const
G4PhysicsTable * InverseRangeTable() const
G4PhysicsTable * DEDXTable() const

◆ NIELCalculator()

G4NIELCalculator * G4LossTableManager::NIELCalculator ( )

Definition at line 960 of file G4LossTableManager.cc.

961{
962 if(!nielCalculator) {
963 nielCalculator = new G4NIELCalculator(nullptr, verbose);
964 }
965 return nielCalculator;
966}

◆ operator=()

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

◆ PreparePhysicsTable() [1/3]

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

Definition at line 469 of file G4LossTableManager.cc.

471{
472 if (1 < verbose) {
473 G4cout << "G4LossTableManager::PreparePhysicsTable for "
474 << particle->GetParticleName()
475 << " and " << p->GetProcessName()
476 << " run=" << run << " master=" << isMaster
477 << G4endl;
478 }
479
480 // start initialisation for the first run
481 if( -1 == run ) {
482 if (nullptr != emConfigurator) { emConfigurator->PrepareModels(particle, p); }
483 }
484
486}

◆ PreparePhysicsTable() [2/3]

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

Definition at line 437 of file G4LossTableManager.cc.

439{
440 if (1 < verbose) {
441 G4cout << "G4LossTableManager::PreparePhysicsTable for "
442 << particle->GetParticleName()
443 << " and " << p->GetProcessName() << " run= " << run
444 << " loss_vector " << loss_vector.size()
445 << " run=" << run << " master=" << isMaster
446 << G4endl;
447 }
448
449 // start initialisation for the first run
450 if( -1 == run ) {
451 if (nullptr != emConfigurator) { emConfigurator->PrepareModels(particle, p); }
452
453 // initialise particles for given process
454 for (G4int j=0; j<n_loss; ++j) {
455 if (p == loss_vector[j] && nullptr == part_vector[j]) {
456 part_vector[j] = particle;
457 if (particle->GetParticleName() == "GenericIon") {
458 theGenericIon = particle;
459 }
460 }
461 }
462 }
464}

◆ PreparePhysicsTable() [3/3]

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

Definition at line 491 of file G4LossTableManager.cc.

493{
494 if (1 < verbose) {
495 G4cout << "G4LossTableManager::PreparePhysicsTable for "
496 << particle->GetParticleName()
497 << " and " << p->GetProcessName()
498 << " run=" << run << " master=" << isMaster
499 << G4endl;
500 }
501
502 // start initialisation for the first run
503 if ( -1 == run ) {
504 if (nullptr != emConfigurator) { emConfigurator->PrepareModels(particle, p); }
505 }
506
508}

◆ Register() [1/6]

void G4LossTableManager::Register ( G4VEmFluctuationModel * p)

Definition at line 368 of file G4LossTableManager.cc.

369{
370 fmod_vector.push_back(p);
371 if(verbose > 1) {
372 G4cout << "G4LossTableManager::Register G4VEmFluctuationModel : "
373 << p->GetName() << " " << fmod_vector.size() << G4endl;
374 }
375}
const G4String & GetName() const

◆ Register() [2/6]

void G4LossTableManager::Register ( G4VEmModel * p)

Definition at line 343 of file G4LossTableManager.cc.

344{
345 mod_vector.push_back(p);
346 if(verbose > 1) {
347 G4cout << "G4LossTableManager::Register G4VEmModel : "
348 << p->GetName() << " " << p << " " << mod_vector.size() << G4endl;
349 }
350}
const G4String & GetName() const

◆ Register() [3/6]

void G4LossTableManager::Register ( G4VEmProcess * p)

Definition at line 283 of file G4LossTableManager.cc.

284{
285 if (nullptr == p) { return; }
286 std::size_t n = emp_vector.size();
287 for (std::size_t i=0; i<n; ++i) {
288 if(emp_vector[i] == p) { return; }
289 }
290 if(verbose > 1) {
291 G4cout << "G4LossTableManager::Register G4VEmProcess : "
292 << p->GetProcessName() << " idx= " << emp_vector.size() << G4endl;
293 }
294 emp_vector.push_back(p);
295}

◆ Register() [4/6]

void G4LossTableManager::Register ( G4VEnergyLossProcess * p)

Definition at line 183 of file G4LossTableManager.cc.

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

◆ Register() [5/6]

void G4LossTableManager::Register ( G4VMultipleScattering * p)

Definition at line 253 of file G4LossTableManager.cc.

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

◆ Register() [6/6]

void G4LossTableManager::Register ( G4VProcess * p)

Definition at line 313 of file G4LossTableManager.cc.

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

◆ RegisterExtraParticle()

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

Definition at line 389 of file G4LossTableManager.cc.

392{
393 if (nullptr == p || nullptr == part) { return; }
394 for (G4int i=0; i<n_loss; ++i) {
395 if(loss_vector[i] == p) { return; }
396 }
397 if(verbose > 1) {
398 G4cout << "G4LossTableManager::RegisterExtraParticle "
399 << part->GetParticleName() << " G4VEnergyLossProcess : "
400 << p->GetProcessName() << " idx= " << n_loss << G4endl;
401 }
402 ++n_loss;
403 loss_vector.push_back(p);
404 part_vector.push_back(part);
405 base_part_vector.push_back(p->BaseParticle());
406 dedx_vector.push_back(nullptr);
407 range_vector.push_back(nullptr);
408 inv_range_vector.push_back(nullptr);
409 tables_are_built.push_back(false);
410 all_tables_are_built = false;
411}

◆ ResetParameters()

void G4LossTableManager::ResetParameters ( )

Definition at line 207 of file G4LossTableManager.cc.

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

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

◆ SetAtomDeexcitation()

void G4LossTableManager::SetAtomDeexcitation ( G4VAtomDeexcitation * p)

Definition at line 970 of file G4LossTableManager.cc.

971{
972 if(atomDeexcitation != p) {
973 delete atomDeexcitation;
974 atomDeexcitation = p;
975 }
976}

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

◆ SetElectronGeneralProcess()

void G4LossTableManager::SetElectronGeneralProcess ( G4VEmProcess * ptr)
inline

Definition at line 425 of file G4LossTableManager.hh.

426{
427 eGeneral = ptr;
428}

◆ SetGammaGeneralProcess()

◆ SetNIELCalculator()

void G4LossTableManager::SetNIELCalculator ( G4NIELCalculator * ptr)

Definition at line 950 of file G4LossTableManager.cc.

951{
952 if(nullptr != ptr && ptr != nielCalculator) {
953 delete nielCalculator;
954 nielCalculator = ptr;
955 }
956}

Referenced by G4NIELCalculator::G4NIELCalculator().

◆ SetPositronGeneralProcess()

void G4LossTableManager::SetPositronGeneralProcess ( G4VEmProcess * ptr)
inline

Definition at line 439 of file G4LossTableManager.hh.

440{
441 pGeneral = ptr;
442}

◆ SetSubCutProducer()

void G4LossTableManager::SetSubCutProducer ( G4VSubCutProducer * p)

Definition at line 980 of file G4LossTableManager.cc.

981{
982 if(subcutProducer != p) {
983 delete subcutProducer;
984 subcutProducer = p;
985 }
986}

◆ SetVerbose()

void G4LossTableManager::SetVerbose ( G4int val)

Definition at line 893 of file G4LossTableManager.cc.

894{
895 verbose = val;
896}

◆ SubCutProducer()

G4VSubCutProducer * G4LossTableManager::SubCutProducer ( )
inline

Definition at line 397 of file G4LossTableManager.hh.

398{
399 return subcutProducer;
400}

Friends And Related Symbol Documentation

◆ G4ThreadLocalSingleton< G4LossTableManager >

Definition at line 446 of file G4LossTableManager.hh.


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