Geant4 11.2.2
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 98 of file G4LossTableManager.cc.

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

◆ G4LossTableManager()

G4LossTableManager::G4LossTableManager ( G4LossTableManager & )
delete

Member Function Documentation

◆ AtomDeexcitation()

◆ BuildPhysicsTable() [1/2]

void G4LossTableManager::BuildPhysicsTable ( const G4ParticleDefinition * aParticle)

Definition at line 516 of file G4LossTableManager.cc.

517{
518 if(-1 == run && startInitialisation) {
519 if(emConfigurator) { emConfigurator->Clear(); }
520 }
521}

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

◆ BuildPhysicsTable() [2/2]

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

Definition at line 603 of file G4LossTableManager.cc.

606{
607 if(1 < verbose) {
608 G4cout << "### G4LossTableManager::BuildPhysicsTable() for "
609 << aParticle->GetParticleName()
610 << " and process " << p->GetProcessName() << G4endl;
611 }
612 // clear configurator
613 if(-1 == run && startInitialisation) {
614 if(emConfigurator) { emConfigurator->Clear(); }
615 firstParticle = aParticle;
616 }
617 if(startInitialisation) {
618 ++run;
619 if(1 < verbose) {
620 G4cout << "===== G4LossTableManager::BuildPhysicsTable() for run "
621 << run << " ===== " << atomDeexcitation << G4endl;
622 }
623 currentParticle = nullptr;
624 all_tables_are_built= true;
625 }
626
627 // initialisation before any table is built
628 if ( startInitialisation && aParticle == firstParticle ) {
629
630 startInitialisation = false;
631 if(1 < verbose) {
632 G4cout << "### G4LossTableManager start initialisation for first particle "
633 << firstParticle->GetParticleName()
634 << G4endl;
635 }
636
637 if(nielCalculator) { nielCalculator->Initialise(); }
638
639 for (G4int i=0; i<n_loss; ++i) {
640 G4VEnergyLossProcess* el = loss_vector[i];
641
642 if(el) {
643 isActive[i] = true;
644 base_part_vector[i] = el->BaseParticle();
645 tables_are_built[i] = false;
646 all_tables_are_built= false;
647 if(1 < verbose) {
648 G4cout << i <<". "<< el->GetProcessName();
649 if(el->Particle()) {
650 G4cout << " for " << el->Particle()->GetParticleName();
651 }
652 G4cout << " active= " << isActive[i]
653 << " table= " << tables_are_built[i]
654 << " isIonisation= " << el->IsIonisationProcess();
655 if(base_part_vector[i]) {
656 G4cout << " base particle "
657 << base_part_vector[i]->GetParticleName();
658 }
659 G4cout << G4endl;
660 }
661 } else {
662 tables_are_built[i] = true;
663 part_vector[i] = nullptr;
664 isActive[i] = false;
665 }
666 }
667 }
668
669 if (all_tables_are_built) {
670 theParameters->SetIsPrintedFlag(true);
671 return;
672 }
673
674 // Build tables for given particle
675 all_tables_are_built = true;
676
677 for(G4int i=0; i<n_loss; ++i) {
678 if(p == loss_vector[i] && !tables_are_built[i] && !base_part_vector[i]) {
679 const G4ParticleDefinition* curr_part = part_vector[i];
680 if(1 < verbose) {
681 G4cout << "### Build Table for " << p->GetProcessName()
682 << " and " << curr_part->GetParticleName()
683 << " " << tables_are_built[i] << " " << base_part_vector[i]
684 << G4endl;
685 }
686 G4VEnergyLossProcess* curr_proc = BuildTables(curr_part);
687 if(curr_proc) {
688 CopyTables(curr_part, curr_proc);
689 if(p == curr_proc && 0 == run && p->IsIonisationProcess()) {
690 loss_map[aParticle] = p;
691 //G4cout << "G4LossTableManager::BuildPhysicsTable: "
692 // << aParticle->GetParticleName()
693 // << " added to map " << p << G4endl;
694 }
695 }
696 }
697 if ( !tables_are_built[i] ) { all_tables_are_built = false; }
698 }
699 if(1 < verbose) {
700 G4cout << "### G4LossTableManager::BuildPhysicsTable end: "
701 << "all_tables_are_built= " << all_tables_are_built << " "
702 << aParticle->GetParticleName() << " proc: " << p << G4endl;
703 }
704 if(all_tables_are_built) {
705 if(1 < verbose) {
706 G4cout << "%%%%% All dEdx and Range tables are built for master run= "
707 << run << " %%%%%" << G4endl;
708 }
709 }
710}
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void SetIsPrintedFlag(G4bool val)
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 366 of file G4LossTableManager.cc.

367{
368 std::size_t n = fmod_vector.size();
369 for (std::size_t i=0; i<n; ++i) {
370 if(fmod_vector[i] == p) { fmod_vector[i] = nullptr; }
371 }
372}

◆ DeRegister() [2/6]

void G4LossTableManager::DeRegister ( G4VEmModel * p)

Definition at line 341 of file G4LossTableManager.cc.

342{
343 //G4cout << "G4LossTableManager::DeRegister G4VEmModel : " << p << G4endl;
344 std::size_t n = mod_vector.size();
345 for (std::size_t i=0; i<n; ++i) {
346 if(mod_vector[i] == p) {
347 mod_vector[i] = nullptr;
348 break;
349 }
350 }
351}

◆ DeRegister() [3/6]

void G4LossTableManager::DeRegister ( G4VEmProcess * p)

Definition at line 286 of file G4LossTableManager.cc.

287{
288 if(!p) { return; }
289 std::size_t emp = emp_vector.size();
290 for (std::size_t i=0; i<emp; ++i) {
291 if(emp_vector[i] == p) {
292 emp_vector[i] = nullptr;
293 break;
294 }
295 }
296}

◆ DeRegister() [4/6]

◆ DeRegister() [5/6]

void G4LossTableManager::DeRegister ( G4VMultipleScattering * p)

Definition at line 256 of file G4LossTableManager.cc.

257{
258 if(!p) { return; }
259 std::size_t msc = msc_vector.size();
260 for (std::size_t i=0; i<msc; ++i) {
261 if(msc_vector[i] == p) {
262 msc_vector[i] = nullptr;
263 break;
264 }
265 }
266}

◆ DeRegister() [6/6]

void G4LossTableManager::DeRegister ( G4VProcess * p)

Definition at line 316 of file G4LossTableManager.cc.

317{
318 if(!p) { return; }
319 std::size_t emp = p_vector.size();
320 for (std::size_t i=0; i<emp; ++i) {
321 if(p_vector[i] == p) {
322 p_vector[i] = nullptr;
323 break;
324 }
325 }
326}

◆ DumpHtml()

void G4LossTableManager::DumpHtml ( )

Definition at line 1041 of file G4LossTableManager.cc.

1042{
1043 // Automatic generation of html documentation page for physics lists
1044 // List processes and models for the most important
1045 // particles in descending order of importance
1046 // NB. for model names with length > 18 characters the .rst file needs
1047 // to be edited by hand. Or modify G4EmModelManager::DumpModelList
1048
1049 char* dirName = std::getenv("G4PhysListDocDir");
1050 char* physList = std::getenv("G4PhysListName");
1051 if (dirName && physList) {
1052 G4String physListName = G4String(physList);
1053 G4String pathName = G4String(dirName) + "/" + physListName + ".rst";
1054
1055 std::ofstream outFile;
1056 outFile.open(pathName);
1057
1058 outFile << physListName << G4endl;
1059 outFile << std::string(physListName.length(), '=') << G4endl;
1060
1061 std::vector<G4ParticleDefinition*> particles {
1068 };
1069
1070 std::vector<G4VEmProcess*> emproc_vector = GetEmProcessVector();
1071 std::vector<G4VEnergyLossProcess*> enloss_vector =
1073 std::vector<G4VMultipleScattering*> mscat_vector =
1075
1076 for (auto theParticle : particles) {
1077 outFile << G4endl << "**" << theParticle->GetParticleName()
1078 << "**" << G4endl << G4endl << " .. code-block:: none" << G4endl;
1079
1080 G4ProcessManager* pm = theParticle->GetProcessManager();
1081 G4ProcessVector* pv = pm->GetProcessList();
1082 G4int plen = pm->GetProcessListLength();
1083
1084 for (auto emproc : emproc_vector) {
1085 for (G4int i = 0; i < plen; ++i) {
1086 G4VProcess* proc = (*pv)[i];
1087 if (proc == emproc) {
1088 outFile << G4endl;
1089 proc->ProcessDescription(outFile);
1090 break;
1091 }
1092 }
1093 }
1094
1095 for (auto mscproc : mscat_vector) {
1096 for (G4int i = 0; i < plen; ++i) {
1097 G4VProcess* proc = (*pv)[i];
1098 if (proc == mscproc) {
1099 outFile << G4endl;
1100 proc->ProcessDescription(outFile);
1101 break;
1102 }
1103 }
1104 }
1105
1106 for (auto enlossproc : enloss_vector) {
1107 for (G4int i = 0; i < plen; ++i) {
1108 G4VProcess* proc = (*pv)[i];
1109 if (proc == enlossproc) {
1110 outFile << G4endl;
1111 proc->ProcessDescription(outFile);
1112 break;
1113 }
1114 }
1115 }
1116 }
1117 outFile.close();
1118 }
1119}
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 974 of file G4LossTableManager.cc.

975{
976 if(!emElectronIonPair) {
977 emElectronIonPair = new G4ElectronIonPair(verbose);
978 }
979 return emElectronIonPair;
980}

◆ EmConfigurator()

G4EmConfigurator * G4LossTableManager::EmConfigurator ( )

Definition at line 964 of file G4LossTableManager.cc.

965{
966 if(!emConfigurator) {
967 emConfigurator = new G4EmConfigurator(verbose);
968 }
969 return emConfigurator;
970}

Referenced by G4TransportationWithMsc::PreparePhysicsTable().

◆ EmCorrections()

◆ EmSaturation()

G4EmSaturation * G4LossTableManager::EmSaturation ( )

Definition at line 957 of file G4LossTableManager.cc.

958{
959 return theParameters->GetEmSaturation();
960}
G4EmSaturation * GetEmSaturation()

Referenced by G4OpticalPhysics::ConstructProcess().

◆ GetCSDARange()

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

Definition at line 318 of file G4LossTableManager.hh.

321{
322 if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
323 return currentLoss ? currentLoss->GetCSDARange(kineticEnergy, couple) : DBL_MAX;
324}
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
G4double GetCSDARange(G4double kineticEnergy, const G4MaterialCutsCouple *)
#define DBL_MAX
Definition templates.hh:62

Referenced by G4EmCalculator::GetCSDARange().

◆ GetDEDX()

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

Definition at line 307 of file G4LossTableManager.hh.

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

Referenced by G4EmCalculator::GetDEDX(), 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 363 of file G4LossTableManager.hh.

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

◆ GetElectronGeneralProcess()

G4VEmProcess * G4LossTableManager::GetElectronGeneralProcess ( )
inline

Definition at line 431 of file G4LossTableManager.hh.

432{
433 return eGeneral;
434}

◆ GetEmProcessVector()

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

Definition at line 942 of file G4LossTableManager.cc.

943{
944 return emp_vector;
945}

Referenced by DumpHtml().

◆ GetEnergy()

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

Definition at line 352 of file G4LossTableManager.hh.

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

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

◆ GetEnergyLossProcess()

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

Definition at line 403 of file G4LossTableManager.cc.

404{
405 if(aParticle != currentParticle) {
406 currentParticle = aParticle;
407 std::map<PD,G4VEnergyLossProcess*,std::less<PD> >::const_iterator pos;
408 if ((pos = loss_map.find(aParticle)) != loss_map.end()) {
409 currentLoss = (*pos).second;
410 } else {
411 currentLoss = nullptr;
412 if(0.0 != aParticle->GetPDGCharge() &&
413 (pos = loss_map.find(theGenericIon)) != loss_map.end()) {
414 currentLoss = (*pos).second;
415 }
416 }
417 }
418 return currentLoss;
419}

Referenced by GetCSDARange(), GetDEDX(), GetDEDXDispersion(), GetEnergy(), GetRange(), GetRangeFromRestricteDEDX(), G4EmCalculator::PrintDEDXTable(), G4EmCalculator::PrintInverseRangeTable(), G4EmCalculator::PrintRangeTable(), G4TransportationWithMsc::StartTracking(), and G4VMultipleScattering::StartTracking().

◆ GetEnergyLossProcessVector()

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

Definition at line 935 of file G4LossTableManager.cc.

936{
937 return loss_vector;
938}

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

◆ GetGammaGeneralProcess()

G4VEmProcess * G4LossTableManager::GetGammaGeneralProcess ( )
inline

Definition at line 417 of file G4LossTableManager.hh.

418{
419 return gGeneral;
420}

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

◆ GetMultipleScatteringVector()

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

Definition at line 950 of file G4LossTableManager.cc.

951{
952 return msc_vector;
953}

Referenced by DumpHtml().

◆ GetPositronGeneralProcess()

G4VEmProcess * G4LossTableManager::GetPositronGeneralProcess ( )
inline

Definition at line 445 of file G4LossTableManager.hh.

446{
447 return pGeneral;
448}

◆ GetRange()

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

Definition at line 341 of file G4LossTableManager.hh.

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

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

◆ GetRangeFromRestricteDEDX()

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

Definition at line 329 of file G4LossTableManager.hh.

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

Referenced by G4EmCalculator::GetRangeFromRestricteDEDX().

◆ GetTableBuilder()

◆ Instance()

G4LossTableManager * G4LossTableManager::Instance ( )
static

Definition at line 87 of file G4LossTableManager.cc.

88{
89 if(nullptr == instance) {
91 instance = inst.Instance();
92 }
93 return instance;
94}

Referenced by G4ITStepProcessor::ApplyProductionCut(), G4BertiniElectroNuclearBuilder::Build(), 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(), 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(), G4VMscModel::GetParticleChangeForMSC(), G4EnergyLossTables::GetPreciseDEDX(), G4EnergyLossTables::GetPreciseEnergyFromRange(), G4EnergyLossTables::GetPreciseRangeFromEnergy(), G4EnergyLossTables::GetRange(), G4DNABornIonisationModel1::Initialise(), G4DNABornIonisationModel2::Initialise(), G4DNACPA100IonisationModel::Initialise(), G4DNAEmfietzoglouIonisationModel::Initialise(), G4DNARelativisticIonisationModel::Initialise(), G4DNARPWBAIonisationModel::Initialise(), G4DNARuddIonisationExtendedModel::Initialise(), G4DNARuddIonisationModel::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()

◆ LocalPhysicsTables()

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

Definition at line 525 of file G4LossTableManager.cc.

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

Referenced by G4VEnergyLossProcess::BuildPhysicsTable().

◆ NIELCalculator()

G4NIELCalculator * G4LossTableManager::NIELCalculator ( )

Definition at line 994 of file G4LossTableManager.cc.

995{
996 if(!nielCalculator) {
997 nielCalculator = new G4NIELCalculator(nullptr, verbose);
998 }
999 return nielCalculator;
1000}

◆ operator=()

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

◆ PreparePhysicsTable() [1/3]

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

Definition at line 462 of file G4LossTableManager.cc.

464{
465 if (1 < verbose) {
466 G4cout << "G4LossTableManager::PreparePhysicsTable for "
467 << particle->GetParticleName()
468 << " and " << p->GetProcessName() << G4endl;
469 }
470
471 if(!startInitialisation) {
473 if (1 < verbose) {
474 G4cout << "====== G4LossTableManager::PreparePhysicsTable start ====="
475 << G4endl;
476 }
477 }
478
479 // start initialisation for the first run
480 if( -1 == run ) {
481 if(emConfigurator) { emConfigurator->PrepareModels(particle, p); }
482 }
483 startInitialisation = true;
484}
void PrepareModels(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)

◆ PreparePhysicsTable() [2/3]

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

Definition at line 424 of file G4LossTableManager.cc.

426{
427 if (1 < verbose) {
428 G4cout << "G4LossTableManager::PreparePhysicsTable for "
429 << particle->GetParticleName()
430 << " and " << p->GetProcessName() << " run= " << run
431 << " loss_vector " << loss_vector.size() << G4endl;
432 }
433
434 if(!startInitialisation) {
436 if (1 < verbose) {
437 G4cout << "====== G4LossTableManager::PreparePhysicsTable start ====="
438 << G4endl;
439 }
440 }
441
442 // start initialisation for the first run
443 if( -1 == run ) {
444 if(emConfigurator) { emConfigurator->PrepareModels(particle, p); }
445
446 // initialise particles for given process
447 for (G4int j=0; j<n_loss; ++j) {
448 if (p == loss_vector[j] && !part_vector[j]) {
449 part_vector[j] = particle;
450 if(particle->GetParticleName() == "GenericIon") {
451 theGenericIon = particle;
452 }
453 }
454 }
455 }
456 startInitialisation = true;
457}

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

◆ PreparePhysicsTable() [3/3]

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

Definition at line 489 of file G4LossTableManager.cc.

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

◆ Register() [1/6]

void G4LossTableManager::Register ( G4VEmFluctuationModel * p)

Definition at line 355 of file G4LossTableManager.cc.

356{
357 fmod_vector.push_back(p);
358 if(verbose > 1) {
359 G4cout << "G4LossTableManager::Register G4VEmFluctuationModel : "
360 << p->GetName() << " " << fmod_vector.size() << G4endl;
361 }
362}
const G4String & GetName() const

◆ Register() [2/6]

void G4LossTableManager::Register ( G4VEmModel * p)

Definition at line 330 of file G4LossTableManager.cc.

331{
332 mod_vector.push_back(p);
333 if(verbose > 1) {
334 G4cout << "G4LossTableManager::Register G4VEmModel : "
335 << p->GetName() << " " << p << " " << mod_vector.size() << G4endl;
336 }
337}
const G4String & GetName() const

◆ Register() [3/6]

void G4LossTableManager::Register ( G4VEmProcess * p)

Definition at line 270 of file G4LossTableManager.cc.

271{
272 if(!p) { return; }
273 std::size_t n = emp_vector.size();
274 for (std::size_t i=0; i<n; ++i) {
275 if(emp_vector[i] == p) { return; }
276 }
277 if(verbose > 1) {
278 G4cout << "G4LossTableManager::Register G4VEmProcess : "
279 << p->GetProcessName() << " idx= " << emp_vector.size() << G4endl;
280 }
281 emp_vector.push_back(p);
282}

◆ Register() [4/6]

void G4LossTableManager::Register ( G4VEnergyLossProcess * p)

Definition at line 183 of file G4LossTableManager.cc.

184{
185 if(!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}

Referenced by G4AnnihiToMuPair::G4AnnihiToMuPair(), G4GammaConversionToMuons::G4GammaConversionToMuons(), G4SynchrotronRadiation::G4SynchrotronRadiation(), G4VEmFluctuationModel::G4VEmFluctuationModel(), G4VEmModel::G4VEmModel(), G4VEmProcess::G4VEmProcess(), G4VEnergyLossProcess::G4VEnergyLossProcess(), G4VMultipleScattering::G4VMultipleScattering(), and G4VTransitionRadiation::G4VTransitionRadiation().

◆ Register() [5/6]

void G4LossTableManager::Register ( G4VMultipleScattering * p)

Definition at line 240 of file G4LossTableManager.cc.

241{
242 if(!p) { return; }
243 std::size_t n = msc_vector.size();
244 for (std::size_t i=0; i<n; ++i) {
245 if(msc_vector[i] == p) { return; }
246 }
247 if(verbose > 1) {
248 G4cout << "G4LossTableManager::Register G4VMultipleScattering : "
249 << p->GetProcessName() << " idx= " << msc_vector.size() << G4endl;
250 }
251 msc_vector.push_back(p);
252}

◆ Register() [6/6]

void G4LossTableManager::Register ( G4VProcess * p)

Definition at line 300 of file G4LossTableManager.cc.

301{
302 if(!p) { return; }
303 std::size_t n = p_vector.size();
304 for (std::size_t i=0; i<n; ++i) {
305 if(p_vector[i] == p) { return; }
306 }
307 if(verbose > 1) {
308 G4cout << "G4LossTableManager::Register G4VProcess : "
309 << p->GetProcessName() << " idx= " << p_vector.size() << G4endl;
310 }
311 p_vector.push_back(p);
312}

◆ RegisterExtraParticle()

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

Definition at line 376 of file G4LossTableManager.cc.

379{
380 if(!p || !part) { return; }
381 for (G4int i=0; i<n_loss; ++i) {
382 if(loss_vector[i] == p) { return; }
383 }
384 if(verbose > 1) {
385 G4cout << "G4LossTableManager::RegisterExtraParticle "
386 << part->GetParticleName() << " G4VEnergyLossProcess : "
387 << p->GetProcessName() << " idx= " << n_loss << G4endl;
388 }
389 ++n_loss;
390 loss_vector.push_back(p);
391 part_vector.push_back(part);
392 base_part_vector.push_back(p->BaseParticle());
393 dedx_vector.push_back(nullptr);
394 range_vector.push_back(nullptr);
395 inv_range_vector.push_back(nullptr);
396 tables_are_built.push_back(false);
397 all_tables_are_built = false;
398}

Referenced by G4VEnergyLossProcess::PreparePhysicsTable().

◆ ResetParameters()

void G4LossTableManager::ResetParameters ( )

Definition at line 207 of file G4LossTableManager.cc.

208{
209 verbose = theParameters->Verbose();
210 if(!isMaster) {
211 verbose = theParameters->WorkerVerbose();
212 } else {
213 if(verbose > 0) { theParameters->Dump(); }
214 }
215 tableBuilder->SetInitialisationFlag(false);
216 emCorrections->SetVerbose(verbose);
217 if(nullptr != emConfigurator) { emConfigurator->SetVerbose(verbose); };
218 if(nullptr != emElectronIonPair) { emElectronIonPair->SetVerbose(verbose); };
219 if(nullptr != atomDeexcitation) {
220 atomDeexcitation->SetVerboseLevel(verbose);
221 atomDeexcitation->InitialiseAtomicDeexcitation();
222 }
223}
void SetVerbose(G4int value)
void SetVerbose(G4int verb)
G4int Verbose() const
G4int WorkerVerbose() const
void SetInitialisationFlag(G4bool flag)

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

◆ SetAtomDeexcitation()

void G4LossTableManager::SetAtomDeexcitation ( G4VAtomDeexcitation * p)

Definition at line 1004 of file G4LossTableManager.cc.

1005{
1006 if(atomDeexcitation != p) {
1007 delete atomDeexcitation;
1008 atomDeexcitation = p;
1009 }
1010}

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

◆ SetElectronGeneralProcess()

void G4LossTableManager::SetElectronGeneralProcess ( G4VEmProcess * ptr)
inline

Definition at line 424 of file G4LossTableManager.hh.

425{
426 eGeneral = ptr;
427}

◆ SetGammaGeneralProcess()

◆ SetNIELCalculator()

void G4LossTableManager::SetNIELCalculator ( G4NIELCalculator * ptr)

Definition at line 984 of file G4LossTableManager.cc.

985{
986 if(ptr && ptr != nielCalculator) {
987 delete nielCalculator;
988 nielCalculator = ptr;
989 }
990}

Referenced by G4NIELCalculator::G4NIELCalculator().

◆ SetPositronGeneralProcess()

void G4LossTableManager::SetPositronGeneralProcess ( G4VEmProcess * ptr)
inline

Definition at line 438 of file G4LossTableManager.hh.

439{
440 pGeneral = ptr;
441}

◆ SetSubCutProducer()

void G4LossTableManager::SetSubCutProducer ( G4VSubCutProducer * p)

Definition at line 1014 of file G4LossTableManager.cc.

1015{
1016 if(subcutProducer != p) {
1017 delete subcutProducer;
1018 subcutProducer = p;
1019 }
1020}

◆ SetVerbose()

void G4LossTableManager::SetVerbose ( G4int val)

Definition at line 927 of file G4LossTableManager.cc.

928{
929 verbose = val;
930}

◆ SubCutProducer()

G4VSubCutProducer * G4LossTableManager::SubCutProducer ( )
inline

Definition at line 396 of file G4LossTableManager.hh.

397{
398 return subcutProducer;
399}

Referenced by G4VEnergyLossProcess::PreparePhysicsTable().

Friends And Related Symbol Documentation

◆ G4ThreadLocalSingleton< G4LossTableManager >

Definition at line 445 of file G4LossTableManager.hh.


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