Geant4 11.1.1
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, G4bool theMaster)
 
void PreparePhysicsTable (const G4ParticleDefinition *aParticle, G4VEmProcess *p, G4bool theMaster)
 
void PreparePhysicsTable (const G4ParticleDefinition *aParticle, G4VMultipleScattering *p, G4bool theMaster)
 
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 (G4int i=0; i<n_loss; ++i) {
101 delete loss_vector[i];
102 }
103 std::size_t msc = msc_vector.size();
104 for (std::size_t j=0; j<msc; ++j) {
105 delete msc_vector[j];
106 }
107 std::size_t emp = emp_vector.size();
108 for (std::size_t k=0; k<emp; ++k) {
109 delete emp_vector[k];
110 }
111 emp = p_vector.size();
112 for (std::size_t k=0; k<emp; ++k) {
113 delete p_vector[k];
114 }
115 std::size_t mod = mod_vector.size();
116 std::size_t fmod = fmod_vector.size();
117 //G4cout << " Nmod" << mod << " Nfluc= " << fmod << G4endl;
118 for (std::size_t a=0; a<mod; ++a) {
119 //G4cout << "Delete model #" << a << " " << mod_vector[a] << G4endl;
120 if( nullptr != mod_vector[a] ) {
121 for (std::size_t b=0; b<fmod; ++b) {
122 if((G4VEmModel*)(fmod_vector[b]) == mod_vector[a]) {
123 fmod_vector[b] = nullptr;
124 }
125 }
126 delete mod_vector[a];
127 mod_vector[a] = nullptr;
128 }
129 }
130 for (std::size_t b=0; b<fmod; ++b) {
131 delete fmod_vector[b];
132 }
133 Clear();
134 delete tableBuilder;
135 delete emCorrections;
136 delete emConfigurator;
137 delete emElectronIonPair;
138 delete nielCalculator;
139 delete atomDeexcitation;
140 delete subcutProducer;
141}
int G4int
Definition: G4Types.hh:85

◆ G4LossTableManager()

G4LossTableManager::G4LossTableManager ( G4LossTableManager )
delete

Member Function Documentation

◆ AtomDeexcitation()

◆ BuildPhysicsTable() [1/2]

void G4LossTableManager::BuildPhysicsTable ( const G4ParticleDefinition aParticle)

Definition at line 537 of file G4LossTableManager.cc.

538{
539 if(-1 == run && startInitialisation) {
540 if(emConfigurator) { emConfigurator->Clear(); }
541 }
542}

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

◆ BuildPhysicsTable() [2/2]

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

Definition at line 626 of file G4LossTableManager.cc.

629{
630 if(1 < verbose) {
631 G4cout << "### G4LossTableManager::BuildPhysicsTable() for "
632 << aParticle->GetParticleName()
633 << " and process " << p->GetProcessName() << G4endl;
634 }
635 // clear configurator
636 if(-1 == run && startInitialisation) {
637 if(emConfigurator) { emConfigurator->Clear(); }
638 firstParticle = aParticle;
639 }
640 if(startInitialisation) {
641 ++run;
642 if(1 < verbose) {
643 G4cout << "===== G4LossTableManager::BuildPhysicsTable() for run "
644 << run << " ===== " << atomDeexcitation << G4endl;
645 }
646 currentParticle = nullptr;
647 all_tables_are_built= true;
648 }
649
650 // initialisation before any table is built
651 if ( startInitialisation && aParticle == firstParticle ) {
652
653 startInitialisation = false;
654 if(1 < verbose) {
655 G4cout << "### G4LossTableManager start initialisation for first particle "
656 << firstParticle->GetParticleName()
657 << G4endl;
658 }
659
660 if(nielCalculator) { nielCalculator->Initialise(); }
661
662 for (G4int i=0; i<n_loss; ++i) {
663 G4VEnergyLossProcess* el = loss_vector[i];
664
665 if(el) {
666 isActive[i] = true;
667 base_part_vector[i] = el->BaseParticle();
668 tables_are_built[i] = false;
669 all_tables_are_built= false;
670 if(!isActive[i]) {
671 el->SetIonisation(false);
672 tables_are_built[i] = true;
673 }
674
675 if(1 < verbose) {
676 G4cout << i <<". "<< el->GetProcessName();
677 if(el->Particle()) {
678 G4cout << " for " << el->Particle()->GetParticleName();
679 }
680 G4cout << " active= " << isActive[i]
681 << " table= " << tables_are_built[i]
682 << " isIonisation= " << el->IsIonisationProcess();
683 if(base_part_vector[i]) {
684 G4cout << " base particle "
685 << base_part_vector[i]->GetParticleName();
686 }
687 G4cout << G4endl;
688 }
689 } else {
690 tables_are_built[i] = true;
691 part_vector[i] = nullptr;
692 isActive[i] = false;
693 }
694 }
695 }
696
697 if (all_tables_are_built) {
698 theParameters->SetIsPrintedFlag(true);
699 return;
700 }
701
702 // Build tables for given particle
703 all_tables_are_built = true;
704
705 for(G4int i=0; i<n_loss; ++i) {
706 if(p == loss_vector[i] && !tables_are_built[i] && !base_part_vector[i]) {
707 const G4ParticleDefinition* curr_part = part_vector[i];
708 if(1 < verbose) {
709 G4cout << "### Build Table for " << p->GetProcessName()
710 << " and " << curr_part->GetParticleName()
711 << " " << tables_are_built[i] << " " << base_part_vector[i]
712 << G4endl;
713 }
714 G4VEnergyLossProcess* curr_proc = BuildTables(curr_part);
715 if(curr_proc) {
716 CopyTables(curr_part, curr_proc);
717 if(p == curr_proc && 0 == run && p->IsIonisationProcess()) {
718 loss_map[aParticle] = p;
719 //G4cout << "G4LossTableManager::BuildPhysicsTable: "
720 // << aParticle->GetParticleName()
721 // << " added to map " << p << G4endl;
722 }
723 }
724 }
725 if ( !tables_are_built[i] ) { all_tables_are_built = false; }
726 }
727 if(1 < verbose) {
728 G4cout << "### G4LossTableManager::BuildPhysicsTable end: "
729 << "all_tables_are_built= " << all_tables_are_built << " "
730 << aParticle->GetParticleName() << " proc: " << p << G4endl;
731 }
732 if(all_tables_are_built) {
733 if(1 < verbose) {
734 G4cout << "%%%%% All dEdx and Range tables are built for master run= "
735 << run << " %%%%%" << G4endl;
736 }
737 }
738}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void SetIsPrintedFlag(G4bool val)
const G4String & GetParticleName() const
const G4ParticleDefinition * BaseParticle() const
const G4ParticleDefinition * Particle() const
G4bool IsIonisationProcess() const
void SetIonisation(G4bool val)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386

◆ DeRegister() [1/6]

void G4LossTableManager::DeRegister ( G4VEmFluctuationModel p)

Definition at line 380 of file G4LossTableManager.cc.

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

◆ DeRegister() [2/6]

void G4LossTableManager::DeRegister ( G4VEmModel p)

Definition at line 355 of file G4LossTableManager.cc.

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

◆ DeRegister() [3/6]

void G4LossTableManager::DeRegister ( G4VEmProcess p)

Definition at line 300 of file G4LossTableManager.cc.

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

◆ DeRegister() [4/6]

◆ DeRegister() [5/6]

void G4LossTableManager::DeRegister ( G4VMultipleScattering p)

Definition at line 270 of file G4LossTableManager.cc.

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

◆ DeRegister() [6/6]

void G4LossTableManager::DeRegister ( G4VProcess p)

Definition at line 330 of file G4LossTableManager.cc.

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

◆ DumpHtml()

void G4LossTableManager::DumpHtml ( )

Definition at line 1071 of file G4LossTableManager.cc.

1072{
1073 // Automatic generation of html documentation page for physics lists
1074 // List processes and models for the most important
1075 // particles in descending order of importance
1076 // NB. for model names with length > 18 characters the .rst file needs
1077 // to be edited by hand. Or modify G4EmModelManager::DumpModelList
1078
1079 char* dirName = std::getenv("G4PhysListDocDir");
1080 char* physList = std::getenv("G4PhysListName");
1081 if (dirName && physList) {
1082 G4String physListName = G4String(physList);
1083 G4String pathName = G4String(dirName) + "/" + physListName + ".rst";
1084
1085 std::ofstream outFile;
1086 outFile.open(pathName);
1087
1088 outFile << physListName << G4endl;
1089 outFile << std::string(physListName.length(), '=') << G4endl;
1090
1091 std::vector<G4ParticleDefinition*> particles {
1098 };
1099
1100 std::vector<G4VEmProcess*> emproc_vector = GetEmProcessVector();
1101 std::vector<G4VEnergyLossProcess*> enloss_vector =
1103 std::vector<G4VMultipleScattering*> mscat_vector =
1105
1106 for (auto theParticle : particles) {
1107 outFile << G4endl << "**" << theParticle->GetParticleName()
1108 << "**" << G4endl << G4endl << " .. code-block:: none" << G4endl;
1109
1110 G4ProcessManager* pm = theParticle->GetProcessManager();
1111 G4ProcessVector* pv = pm->GetProcessList();
1112 G4int plen = pm->GetProcessListLength();
1113
1114 for (auto emproc : emproc_vector) {
1115 for (G4int i = 0; i < plen; ++i) {
1116 G4VProcess* proc = (*pv)[i];
1117 if (proc == emproc) {
1118 outFile << G4endl;
1119 proc->ProcessDescription(outFile);
1120 break;
1121 }
1122 }
1123 }
1124
1125 for (auto mscproc : mscat_vector) {
1126 for (G4int i = 0; i < plen; ++i) {
1127 G4VProcess* proc = (*pv)[i];
1128 if (proc == mscproc) {
1129 outFile << G4endl;
1130 proc->ProcessDescription(outFile);
1131 break;
1132 }
1133 }
1134 }
1135
1136 for (auto enlossproc : enloss_vector) {
1137 for (G4int i = 0; i < plen; ++i) {
1138 G4VProcess* proc = (*pv)[i];
1139 if (proc == enlossproc) {
1140 outFile << G4endl;
1141 proc->ProcessDescription(outFile);
1142 break;
1143 }
1144 }
1145 }
1146 }
1147 outFile.close();
1148 }
1149}
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
const std::vector< G4VEmProcess * > & GetEmProcessVector()
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector()
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
static G4MuonMinus * MuonMinusDefinition()
Definition: G4MuonMinus.cc:94
static G4MuonPlus * MuonPlusDefinition()
Definition: G4MuonPlus.cc:93
static G4Positron * Positron()
Definition: G4Positron.cc:93
G4int GetProcessListLength() const
G4ProcessVector * GetProcessList() const
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:87
virtual void ProcessDescription(std::ostream &outfile) const
Definition: G4VProcess.cc:181

◆ ElectronIonPair()

G4ElectronIonPair * G4LossTableManager::ElectronIonPair ( )

Definition at line 1004 of file G4LossTableManager.cc.

1005{
1006 if(!emElectronIonPair) {
1007 emElectronIonPair = new G4ElectronIonPair(verbose);
1008 }
1009 return emElectronIonPair;
1010}

◆ EmConfigurator()

G4EmConfigurator * G4LossTableManager::EmConfigurator ( )

Definition at line 994 of file G4LossTableManager.cc.

995{
996 if(!emConfigurator) {
997 emConfigurator = new G4EmConfigurator(verbose);
998 }
999 return emConfigurator;
1000}

Referenced by G4TransportationWithMsc::PreparePhysicsTable().

◆ EmCorrections()

◆ EmSaturation()

G4EmSaturation * G4LossTableManager::EmSaturation ( )

Definition at line 987 of file G4LossTableManager.cc.

988{
989 return theParameters->GetEmSaturation();
990}
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 G4EnergyLossTables::GetDEDX(), G4EmCalculator::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 972 of file G4LossTableManager.cc.

973{
974 return emp_vector;
975}

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 417 of file G4LossTableManager.cc.

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

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

◆ GetEnergyLossProcessVector()

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

Definition at line 965 of file G4LossTableManager.cc.

966{
967 return loss_vector;
968}

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 980 of file G4LossTableManager.cc.

981{
982 return msc_vector;
983}

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(!instance) {
91 instance = inst.Instance();
92 }
93 return instance;
94}

Referenced by G4ITStepProcessor::ApplyProductionCut(), G4BertiniElectroNuclearBuilder::Build(), G4GammaGeneralProcess::BuildPhysicsTable(), G4EmCalculator::ComputeDEDXForCutInRange(), G4EmCalculator::ComputeElectronicDEDX(), LBE::ConstructGeneral(), G4RadioactiveDecayPhysics::ConstructProcess(), G4EmStandardPhysics::ConstructProcess(), G4EmStandardPhysics_option1::ConstructProcess(), G4EmStandardPhysics_option2::ConstructProcess(), G4EmStandardPhysics_option3::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4EmStandardPhysicsSS::ConstructProcess(), G4OpticalPhysics::ConstructProcess(), G4EmExtraPhysics::ConstructProcess(), G4ECDecay::DecayIt(), G4ITDecay::DecayIt(), G4AnnihiToMuPair::G4AnnihiToMuPair(), G4AtimaEnergyLossModel::G4AtimaEnergyLossModel(), G4BetheBlochModel::G4BetheBlochModel(), G4BraggModel::G4BraggModel(), G4EmCalculator::G4EmCalculator(), G4GammaConversionToMuons::G4GammaConversionToMuons(), G4GoudsmitSaundersonMscModel::G4GoudsmitSaundersonMscModel(), G4ionIonisation::G4ionIonisation(), 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(), G4DNARuddIonisationExtendedModel::Initialise(), G4DNARuddIonisationModel::Initialise(), G4LivermoreComptonModel::Initialise(), G4LivermorePhotoElectricModel::Initialise(), G4LivermorePolarizedComptonModel::Initialise(), G4LowEPComptonModel::Initialise(), G4LowEPPolarizedComptonModel::Initialise(), G4MicroElecInelasticModel::Initialise(), G4MicroElecInelasticModel_new::Initialise(), G4PenelopeComptonModel::Initialise(), G4PenelopeIonisationModel::Initialise(), G4PenelopePhotoElectricModel::Initialise(), G4BraggIonModel::Initialise(), G4KleinNishinaModel::Initialise(), G4PEEffectFluoModel::Initialise(), G4DNABornIonisationModel1::Initialise(), G4DNABornIonisationModel2::Initialise(), G4DNACPA100IonisationModel::Initialise(), G4DNAEmfietzoglouIonisationModel::Initialise(), G4DNARelativisticIonisationModel::Initialise(), G4DNARPWBAIonisationModel::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 546 of file G4LossTableManager.cc.

549{
550 if(1 < verbose) {
551 G4cout << "### G4LossTableManager::LocalPhysicsTable() for "
552 << aParticle->GetParticleName()
553 << " and process " << p->GetProcessName()
554 << G4endl;
555 }
556
557 if(-1 == run && startInitialisation) {
558 if(emConfigurator) { emConfigurator->Clear(); }
559 firstParticle = aParticle;
560 }
561
562 if(startInitialisation) {
563 ++run;
564 if(1 < verbose) {
565 G4cout << "===== G4LossTableManager::LocalPhysicsTable() for run "
566 << run << " =====" << G4endl;
567 }
568 currentParticle = nullptr;
569 startInitialisation = false;
570 for (G4int i=0; i<n_loss; ++i) {
571 if(loss_vector[i]) {
572 tables_are_built[i] = false;
573 } else {
574 tables_are_built[i] = true;
575 part_vector[i] = nullptr;
576 }
577 }
578 }
579
580 all_tables_are_built= true;
581 for (G4int i=0; i<n_loss; ++i) {
582 if(p == loss_vector[i]) {
583 tables_are_built[i] = true;
584 isActive[i] = true;
585 part_vector[i] = p->Particle();
586 base_part_vector[i] = p->BaseParticle();
587 dedx_vector[i] = p->DEDXTable();
588 range_vector[i] = p->RangeTableForLoss();
589 inv_range_vector[i] = p->InverseRangeTable();
590 if(0 == run && p->IsIonisationProcess()) {
591 loss_map[part_vector[i]] = p;
592 //G4cout << "G4LossTableManager::LocalPhysicsTable " << part_vector[i]->GetParticleName()
593 // << " added to map " << p << G4endl;
594 }
595
596 if(1 < verbose) {
597 G4cout << i <<". "<< p->GetProcessName();
598 if(part_vector[i]) {
599 G4cout << " for " << part_vector[i]->GetParticleName();
600 }
601 G4cout << " active= " << isActive[i]
602 << " table= " << tables_are_built[i]
603 << " isIonisation= " << p->IsIonisationProcess()
604 << G4endl;
605 }
606 break;
607 } else if(!tables_are_built[i]) {
608 all_tables_are_built = false;
609 }
610 }
611
612 if(1 < verbose) {
613 G4cout << "### G4LossTableManager::LocalPhysicsTable end"
614 << G4endl;
615 }
616 if(all_tables_are_built) {
617 if(1 < verbose) {
618 G4cout << "%%%%% All dEdx and Range tables for worker are ready for run "
619 << run << " %%%%%" << G4endl;
620 }
621 }
622}
G4PhysicsTable * RangeTableForLoss() const
G4PhysicsTable * InverseRangeTable() const
G4PhysicsTable * DEDXTable() const

Referenced by G4VEnergyLossProcess::BuildPhysicsTable().

◆ NIELCalculator()

G4NIELCalculator * G4LossTableManager::NIELCalculator ( )

Definition at line 1024 of file G4LossTableManager.cc.

1025{
1026 if(!nielCalculator) {
1027 nielCalculator = new G4NIELCalculator(nullptr, verbose);
1028 }
1029 return nielCalculator;
1030}

◆ operator=()

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

◆ PreparePhysicsTable() [1/3]

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

Definition at line 479 of file G4LossTableManager.cc.

481{
482 if (1 < verbose) {
483 G4cout << "G4LossTableManager::PreparePhysicsTable for "
484 << particle->GetParticleName()
485 << " and " << p->GetProcessName() << G4endl;
486 }
487 isMaster = theMaster;
488
489 if(!startInitialisation) {
491 if (1 < verbose) {
492 G4cout << "====== G4LossTableManager::PreparePhysicsTable start ====="
493 << G4endl;
494 }
495 }
496
497 // start initialisation for the first run
498 if( -1 == run ) {
499 if(emConfigurator) { emConfigurator->PrepareModels(particle, p); }
500 }
501 startInitialisation = true;
502}
void PrepareModels(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)

◆ PreparePhysicsTable() [2/3]

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

Definition at line 438 of file G4LossTableManager.cc.

441{
442 if (1 < verbose) {
443 G4cout << "G4LossTableManager::PreparePhysicsTable for "
444 << particle->GetParticleName()
445 << " and " << p->GetProcessName() << " run= " << run
446 << " loss_vector " << loss_vector.size() << G4endl;
447 }
448
449 isMaster = theMaster;
450
451 if(!startInitialisation) {
453 if (1 < verbose) {
454 G4cout << "====== G4LossTableManager::PreparePhysicsTable start ====="
455 << G4endl;
456 }
457 }
458
459 // start initialisation for the first run
460 if( -1 == run ) {
461 if(emConfigurator) { emConfigurator->PrepareModels(particle, p); }
462
463 // initialise particles for given process
464 for (G4int j=0; j<n_loss; ++j) {
465 if (p == loss_vector[j] && !part_vector[j]) {
466 part_vector[j] = particle;
467 if(particle->GetParticleName() == "GenericIon") {
468 theGenericIon = particle;
469 }
470 }
471 }
472 }
473 startInitialisation = true;
474}

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

◆ PreparePhysicsTable() [3/3]

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

Definition at line 507 of file G4LossTableManager.cc.

510{
511 if (1 < verbose) {
512 G4cout << "G4LossTableManager::PreparePhysicsTable for "
513 << particle->GetParticleName()
514 << " and " << p->GetProcessName() << G4endl;
515 }
516
517 isMaster = theMaster;
518
519 if(!startInitialisation) {
521 if (1 < verbose) {
522 G4cout << "====== G4LossTableManager::PreparePhysicsTable start ====="
523 << G4endl;
524 }
525 }
526
527 // start initialisation for the first run
528 if( -1 == run ) {
529 if(emConfigurator) { emConfigurator->PrepareModels(particle, p); }
530 }
531 startInitialisation = true;
532}

◆ Register() [1/6]

void G4LossTableManager::Register ( G4VEmFluctuationModel p)

Definition at line 369 of file G4LossTableManager.cc.

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

◆ Register() [2/6]

void G4LossTableManager::Register ( G4VEmModel p)

Definition at line 344 of file G4LossTableManager.cc.

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

◆ Register() [3/6]

void G4LossTableManager::Register ( G4VEmProcess p)

Definition at line 284 of file G4LossTableManager.cc.

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

◆ Register() [4/6]

void G4LossTableManager::Register ( G4VEnergyLossProcess p)

Definition at line 197 of file G4LossTableManager.cc.

198{
199 if(!p) { return; }
200 for (G4int i=0; i<n_loss; ++i) {
201 if(loss_vector[i] == p) { return; }
202 }
203 if(verbose > 1) {
204 G4cout << "G4LossTableManager::Register G4VEnergyLossProcess : "
205 << p->GetProcessName() << " idx= " << n_loss << G4endl;
206 }
207 ++n_loss;
208 loss_vector.push_back(p);
209 part_vector.push_back(nullptr);
210 base_part_vector.push_back(nullptr);
211 dedx_vector.push_back(nullptr);
212 range_vector.push_back(nullptr);
213 inv_range_vector.push_back(nullptr);
214 tables_are_built.push_back(false);
215 isActive.push_back(true);
216 all_tables_are_built = false;
217}

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 254 of file G4LossTableManager.cc.

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

◆ Register() [6/6]

void G4LossTableManager::Register ( G4VProcess p)

Definition at line 314 of file G4LossTableManager.cc.

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

◆ RegisterExtraParticle()

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

Definition at line 390 of file G4LossTableManager.cc.

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

Referenced by G4VEnergyLossProcess::PreparePhysicsTable().

◆ ResetParameters()

void G4LossTableManager::ResetParameters ( )

Definition at line 221 of file G4LossTableManager.cc.

222{
223 verbose = theParameters->Verbose();
224 if(!isMaster) {
225 verbose = theParameters->WorkerVerbose();
226 } else {
227 if(verbose > 0) { theParameters->Dump(); }
228 }
229 tableBuilder->SetInitialisationFlag(false);
230 emCorrections->SetVerbose(verbose);
231 if(nullptr != emConfigurator) { emConfigurator->SetVerbose(verbose); };
232 if(nullptr != emElectronIonPair) { emElectronIonPair->SetVerbose(verbose); };
233 if(nullptr != atomDeexcitation) {
234 atomDeexcitation->SetVerboseLevel(verbose);
235 atomDeexcitation->InitialiseAtomicDeexcitation();
236 }
237}
void SetVerbose(G4int value)
void SetVerbose(G4int verb)
G4int Verbose() const
G4int WorkerVerbose() const
void SetInitialisationFlag(G4bool flag)

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

◆ SetAtomDeexcitation()

void G4LossTableManager::SetAtomDeexcitation ( G4VAtomDeexcitation p)

Definition at line 1034 of file G4LossTableManager.cc.

1035{
1036 if(atomDeexcitation != p) {
1037 delete atomDeexcitation;
1038 atomDeexcitation = p;
1039 }
1040}

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 1014 of file G4LossTableManager.cc.

1015{
1016 if(ptr && ptr != nielCalculator) {
1017 delete nielCalculator;
1018 nielCalculator = ptr;
1019 }
1020}

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 1044 of file G4LossTableManager.cc.

1045{
1046 if(subcutProducer != p) {
1047 delete subcutProducer;
1048 subcutProducer = p;
1049 }
1050}

◆ SetVerbose()

void G4LossTableManager::SetVerbose ( G4int  val)

Definition at line 957 of file G4LossTableManager.cc.

958{
959 verbose = val;
960}

◆ 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 Function Documentation

◆ G4ThreadLocalSingleton< G4LossTableManager >

Definition at line 445 of file G4LossTableManager.hh.


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