Geant4  10.02.p03
G4LossTableManager Class Reference

#include <G4LossTableManager.hh>

Collaboration diagram for G4LossTableManager:

Public Member Functions

 ~G4LossTableManager ()
 
void Clear ()
 
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)
 
G4double GetDEDX (const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetSubDEDX (const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetRange (const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetCSDARange (const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetRangeFromRestricteDEDX (const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetEnergy (const G4ParticleDefinition *aParticle, G4double range, const G4MaterialCutsCouple *couple)
 
G4double GetDEDXDispersion (const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double &length)
 
void Register (G4VEnergyLossProcess *p)
 
void DeRegister (G4VEnergyLossProcess *p)
 
void Register (G4VMultipleScattering *p)
 
void DeRegister (G4VMultipleScattering *p)
 
void Register (G4VEmProcess *p)
 
void DeRegister (G4VEmProcess *p)
 
void Register (G4VEmModel *p)
 
void DeRegister (G4VEmModel *p)
 
void Register (G4VEmFluctuationModel *p)
 
void DeRegister (G4VEmFluctuationModel *p)
 
void RegisterExtraParticle (const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
 
void SetLossFluctuations (G4bool val)
 
void SetSubCutoff (G4bool val, const G4Region *r=0)
 
void SetIntegral (G4bool val)
 
void SetRandomStep (G4bool val)
 
void SetMinSubRange (G4double val)
 
void SetMinEnergy (G4double val)
 
void SetMaxEnergy (G4double val)
 
void SetMaxEnergyForCSDARange (G4double val)
 
void SetMaxEnergyForMuons (G4double val)
 
void SetDEDXBinning (G4int val)
 
void SetDEDXBinningForCSDARange (G4int val)
 
void SetLambdaBinning (G4int val)
 
void SetStepFunction (G4double v1, G4double v2)
 
void SetBuildCSDARange (G4bool val)
 
void SetLinearLossLimit (G4double val)
 
void SetVerbose (G4int val)
 
void SetAtomDeexcitation (G4VAtomDeexcitation *)
 
void SetSubCutProducer (G4VSubCutProducer *)
 
G4EnergyLossMessengerGetMessenger ()
 
G4bool BuildCSDARange () const
 
G4bool IsMaster () const
 
G4double MinKinEnergy () const
 
G4double MaxKinEnergy () const
 
G4int GetNumberOfBinsPerDecade () const
 
G4VEnergyLossProcessGetEnergyLossProcess (const G4ParticleDefinition *)
 
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector ()
 
const std::vector< G4VEmProcess * > & GetEmProcessVector ()
 
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector ()
 
G4EmSaturationEmSaturation ()
 
G4EmConfiguratorEmConfigurator ()
 
G4ElectronIonPairElectronIonPair ()
 
G4EmCorrectionsEmCorrections ()
 
G4VAtomDeexcitationAtomDeexcitation ()
 
G4VSubCutProducerSubCutProducer ()
 
G4LossTableBuilderGetTableBuilder ()
 

Static Public Member Functions

static G4LossTableManagerInstance ()
 

Private Types

typedef const G4ParticleDefinitionPD
 

Private Member Functions

 G4LossTableManager ()
 
void ResetParameters ()
 
G4VEnergyLossProcessBuildTables (const G4ParticleDefinition *aParticle)
 
void CopyTables (const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *)
 
void ParticleHaveNoLoss (const G4ParticleDefinition *aParticle)
 
void SetParameters (const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *)
 
void CopyDEDXTables ()
 
void PrintEWarning (G4String, G4double)
 
 G4LossTableManager (G4LossTableManager &)
 
G4LossTableManageroperator= (const G4LossTableManager &right)
 

Private Attributes

std::map< PD, G4VEnergyLossProcess *, std::less< PD > > loss_map
 
std::vector< G4VEnergyLossProcess * > loss_vector
 
std::vector< PDpart_vector
 
std::vector< PDbase_part_vector
 
std::vector< G4booltables_are_built
 
std::vector< G4boolisActive
 
std::vector< G4PhysicsTable * > dedx_vector
 
std::vector< G4PhysicsTable * > range_vector
 
std::vector< G4PhysicsTable * > inv_range_vector
 
std::vector< G4VMultipleScattering * > msc_vector
 
std::vector< G4VEmProcess * > emp_vector
 
std::vector< G4VEmModel * > mod_vector
 
std::vector< G4VEmFluctuationModel * > fmod_vector
 
G4VEnergyLossProcesscurrentLoss
 
PD currentParticle
 
PD theElectron
 
PD theGenericIon
 
PD firstParticle
 
G4int n_loss
 
G4int run
 
G4bool all_tables_are_built
 
G4bool startInitialisation
 
G4bool subCutoffFlag
 
G4bool integral
 
G4bool integralActive
 
G4bool stepFunctionActive
 
G4bool isMaster
 
G4double maxRangeVariation
 
G4double maxFinalStep
 
G4LossTableBuildertableBuilder
 
G4EnergyLossMessengertheMessenger
 
G4EmCorrectionsemCorrections
 
G4EmSaturationemSaturation
 
G4EmConfiguratoremConfigurator
 
G4ElectronIonPairemElectronIonPair
 
G4VAtomDeexcitationatomDeexcitation
 
G4VSubCutProducersubcutProducer
 
G4EmParameterstheParameters
 
G4int verbose
 

Static Private Attributes

static G4ThreadLocal G4LossTableManagerinstance = 0
 

Friends

class G4ThreadLocalSingleton< G4LossTableManager >
 

Detailed Description

Definition at line 102 of file G4LossTableManager.hh.

Member Typedef Documentation

◆ PD

Definition at line 311 of file G4LossTableManager.hh.

Constructor & Destructor Documentation

◆ ~G4LossTableManager()

G4LossTableManager::~G4LossTableManager ( )

Definition at line 123 of file G4LossTableManager.cc.

124 {
125  //G4cout << "### G4LossTableManager::~G4LossTableManager()" << G4endl;
126  for (G4int i=0; i<n_loss; ++i) {
127  //G4cout << "### eloss #" << i << G4endl;
128  if( loss_vector[i] ) {
129  delete loss_vector[i];
130  }
131  }
132  size_t msc = msc_vector.size();
133  for (size_t j=0; j<msc; ++j) {
134  if( msc_vector[j] ) { delete msc_vector[j]; }
135  }
136  size_t emp = emp_vector.size();
137  for (size_t k=0; k<emp; ++k) {
138  if( emp_vector[k] ) { delete emp_vector[k]; }
139  }
140  size_t mod = mod_vector.size();
141  size_t fmod = fmod_vector.size();
142  for (size_t a=0; a<mod; ++a) {
143  if( mod_vector[a] ) {
144  for (size_t b=0; b<fmod; ++b) {
145  if((G4VEmModel*)(fmod_vector[b]) == mod_vector[a]) {
146  fmod_vector[b] = 0;
147  }
148  }
149  delete mod_vector[a];
150  }
151  }
152  for (size_t b=0; b<fmod; ++b) {
153  if( fmod_vector[b] ) { delete fmod_vector[b]; }
154  }
155  Clear();
156  delete theMessenger;
157  delete tableBuilder;
158  delete emCorrections;
159  delete emSaturation;
160  delete emConfigurator;
161  delete emElectronIonPair;
162  delete atomDeexcitation;
163  delete subcutProducer;
164 }
G4VAtomDeexcitation * atomDeexcitation
G4VSubCutProducer * subcutProducer
std::vector< G4VEmModel * > mod_vector
G4ElectronIonPair * emElectronIonPair
G4EmSaturation * emSaturation
int G4int
Definition: G4Types.hh:78
G4EmConfigurator * emConfigurator
std::vector< G4VEmFluctuationModel * > fmod_vector
G4EmCorrections * emCorrections
G4LossTableBuilder * tableBuilder
std::vector< G4VEnergyLossProcess * > loss_vector
std::vector< G4VEmProcess * > emp_vector
std::vector< G4VMultipleScattering * > msc_vector
G4EnergyLossMessenger * theMessenger
Here is the call graph for this function:

◆ G4LossTableManager() [1/2]

G4LossTableManager::G4LossTableManager ( )
private

Definition at line 168 of file G4LossTableManager.cc.

169 {
171  n_loss = 0;
172  run = -1;
173  startInitialisation = false;
174  all_tables_are_built = false;
175  currentLoss = nullptr;
176  currentParticle = nullptr;
177  firstParticle = nullptr;
178  subCutoffFlag = false;
179  maxRangeVariation = 0.2;
181  integral = true;
182  integralActive = false;
183  stepFunctionActive = false;
184  isMaster = true;
188  theGenericIon= nullptr;
191  isMaster = false;
192  }
195  emSaturation = nullptr;
196  emConfigurator = nullptr;
197  emElectronIonPair = nullptr;
198  atomDeexcitation = nullptr;
199  subcutProducer = nullptr;
200 }
G4VAtomDeexcitation * atomDeexcitation
G4VSubCutProducer * subcutProducer
G4ElectronIonPair * emElectronIonPair
G4EmSaturation * emSaturation
G4int Verbose() const
G4EmConfigurator * emConfigurator
G4VEnergyLossProcess * currentLoss
G4bool IsWorkerThread()
Definition: G4Threading.cc:135
G4EmParameters * theParameters
G4EmCorrections * emCorrections
static const double mm
Definition: SystemOfUnits.h:94
G4LossTableBuilder * tableBuilder
static G4EmParameters * Instance()
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4int WorkerVerbose() const
G4EnergyLossMessenger * theMessenger
Here is the call graph for this function:

◆ G4LossTableManager() [2/2]

G4LossTableManager::G4LossTableManager ( G4LossTableManager )
private

Member Function Documentation

◆ AtomDeexcitation()

G4VAtomDeexcitation * G4LossTableManager::AtomDeexcitation ( )
inline

Definition at line 528 of file G4LossTableManager.hh.

529 {
530  return atomDeexcitation;
531 }
G4VAtomDeexcitation * atomDeexcitation
Here is the caller graph for this function:

◆ BuildCSDARange()

G4bool G4LossTableManager::BuildCSDARange ( ) const
inline

Definition at line 486 of file G4LossTableManager.hh.

487 {
488  return theParameters->BuildCSDARange();
489 }
G4EmParameters * theParameters
G4bool BuildCSDARange() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ BuildPhysicsTable() [1/2]

void G4LossTableManager::BuildPhysicsTable ( const G4ParticleDefinition aParticle)

Definition at line 506 of file G4LossTableManager.cc.

507 {
508  if(-1 == run && startInitialisation) {
510  }
511 }
G4EmConfigurator * emConfigurator
Here is the call graph for this function:
Here is the caller graph for this function:

◆ BuildPhysicsTable() [2/2]

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

Definition at line 600 of file G4LossTableManager.cc.

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

◆ BuildTables()

G4VEnergyLossProcess * G4LossTableManager::BuildTables ( const G4ParticleDefinition aParticle)
private

Definition at line 753 of file G4LossTableManager.cc.

755 {
756  if(1 < verbose) {
757  G4cout << "G4LossTableManager::BuildTables() for "
758  << aParticle->GetParticleName() << G4endl;
759  }
760 
761  std::vector<G4PhysicsTable*> t_list;
762  std::vector<G4VEnergyLossProcess*> loss_list;
763  std::vector<G4bool> build_flags;
764  G4VEnergyLossProcess* em = nullptr;
765  G4VEnergyLossProcess* p = nullptr;
766  G4int iem = 0;
767  G4PhysicsTable* dedx = 0;
768  G4int i;
769 
770  G4ProcessVector* pvec =
771  aParticle->GetProcessManager()->GetProcessList();
772  G4int nvec = pvec->size();
773 
774  for (i=0; i<n_loss; ++i) {
775  p = loss_vector[i];
776  if (p) {
777  G4bool yes = (aParticle == part_vector[i]);
778 
779  // possible case of process sharing between particle/anti-particle
780  if(!yes) {
781  G4VProcess* ptr = static_cast<G4VProcess*>(p);
782  for(G4int j=0; j<nvec; ++j) {
783  //G4cout << "j= " << j << " " << (*pvec)[j] << " " << ptr << G4endl;
784  if(ptr == (*pvec)[j]) {
785  yes = true;
786  break;
787  }
788  }
789  }
790  // process belong to this particle
791  if(yes && isActive[i]) {
792  if (p->IsIonisationProcess() || !em) {
793  em = p;
794  iem= i;
795  }
796  // tables may be shared between particle/anti-particle
797  G4bool val = false;
798  if (!tables_are_built[i]) {
799  val = true;
800  dedx = p->BuildDEDXTable(fRestricted);
801  //G4cout << "Build DEDX table for " << p->GetProcessName()
802  // << " idx= " << i << dedx << " " << dedx->length() << G4endl;
803  p->SetDEDXTable(dedx,fRestricted);
804  tables_are_built[i] = true;
805  } else {
806  dedx = p->DEDXTable();
807  }
808  t_list.push_back(dedx);
809  loss_list.push_back(p);
810  build_flags.push_back(val);
811  }
812  }
813  }
814 
815  G4int n_dedx = t_list.size();
816  if (0 == n_dedx || !em) {
817  G4cout << "G4LossTableManager WARNING: no DEDX processes for "
818  << aParticle->GetParticleName() << G4endl;
819  return 0;
820  }
821  G4int nSubRegions = em->NumberOfSubCutoffRegions();
822 
823  if (1 < verbose) {
824  G4cout << "G4LossTableManager::BuildTables() start to build range tables"
825  << " and the sum of " << n_dedx << " processes"
826  << " iem= " << iem << " em= " << em->GetProcessName()
827  << " buildCSDARange= " << theParameters->BuildCSDARange()
828  << " nSubRegions= " << nSubRegions;
829  if(subcutProducer) {
830  G4cout << " SubCutProducer " << subcutProducer->GetName();
831  }
832  G4cout << G4endl;
833  }
834  // do not build tables if producer class is defined
835  if(subcutProducer) { nSubRegions = 0; }
836 
837  dedx = em->DEDXTable();
838  em->SetIonisation(true);
839  em->SetDEDXTable(dedx, fIsIonisation);
840 
841  if (1 < n_dedx) {
842  dedx = 0;
844  tableBuilder->BuildDEDXTable(dedx, t_list);
845  em->SetDEDXTable(dedx, fRestricted);
846  }
847 
848  /*
849  if(2==run && "e-" == aParticle->GetParticleName()) {
850  G4cout << "G4LossTableManager::BuildTables for e- " << dedx << G4endl;
851  G4cout << (*dedx) << G4endl;
852  G4cout << "%%%%% Instance ID= " << (*dedx)[0]->GetInstanceID() << G4endl;
853  G4cout << "%%%%% LastValue= " << (*dedx)[0]->GetLastValue() << G4endl;
854  G4cout << "%%%%% 1.2 " << (*(dedx))[0]->Value(1.2) << G4endl;
855  }
856  */
857  dedx_vector[iem] = dedx;
858 
859  G4PhysicsTable* range = em->RangeTableForLoss();
860  if(!range) range = G4PhysicsTableHelper::PreparePhysicsTable(range);
861  range_vector[iem] = range;
862 
863  G4PhysicsTable* invrange = em->InverseRangeTable();
864  if(!invrange) invrange = G4PhysicsTableHelper::PreparePhysicsTable(invrange);
865  inv_range_vector[iem] = invrange;
866 
867  tableBuilder->BuildRangeTable(dedx, range, true);
868  tableBuilder->BuildInverseRangeTable(range, invrange, true);
869 
870  // if(1<verbose) G4cout << *dedx << G4endl;
871 
872  em->SetRangeTableForLoss(range);
873  em->SetInverseRangeTable(invrange);
874 
875  // if(1<verbose) G4cout << *range << G4endl;
876 
877  std::vector<G4PhysicsTable*> listSub;
878  std::vector<G4PhysicsTable*> listCSDA;
879 
880  for (i=0; i<n_dedx; ++i) {
881  p = loss_list[i];
882  if(p != em) { p->SetIonisation(false); }
883  if(build_flags[i]) {
885  }
886  if (0 < nSubRegions) {
887  dedx = p->BuildDEDXTable(fSubRestricted);
888  p->SetDEDXTable(dedx,fSubRestricted);
889  listSub.push_back(dedx);
890  if(build_flags[i]) {
892  if(p != em) { em->AddCollaborativeProcess(p); }
893  }
894  }
895  if(theParameters->BuildCSDARange()) {
896  dedx = p->BuildDEDXTable(fTotal);
897  p->SetDEDXTable(dedx,fTotal);
898  listCSDA.push_back(dedx);
899  }
900  }
901 
902  if (0 < nSubRegions) {
903  G4PhysicsTable* dedxSub = em->IonisationTableForSubsec();
904  if (1 < listSub.size()) {
905  em->SetDEDXTable(dedxSub, fIsSubIonisation);
906  dedxSub = 0;
908  tableBuilder->BuildDEDXTable(dedxSub, listSub);
909  em->SetDEDXTable(dedxSub, fSubRestricted);
910  }
911  }
913  G4PhysicsTable* dedxCSDA = em->DEDXunRestrictedTable();
914  if (1 < n_dedx) {
915  dedxCSDA = nullptr;
916  dedxCSDA = G4PhysicsTableHelper::PreparePhysicsTable(dedxCSDA);
917  tableBuilder->BuildDEDXTable(dedxCSDA, listCSDA);
918  em->SetDEDXTable(dedxCSDA,fTotal);
919  }
920  G4PhysicsTable* rCSDA = em->CSDARangeTable();
921  if(!rCSDA) { rCSDA = G4PhysicsTableHelper::PreparePhysicsTable(rCSDA); }
922  tableBuilder->BuildRangeTable(dedxCSDA, rCSDA, true);
923  em->SetCSDARangeTable(rCSDA);
924  }
925 
926  if (1 < verbose) {
927  G4cout << "G4LossTableManager::BuildTables: Tables are built for "
928  << aParticle->GetParticleName()
929  << "; ionisation process: " << em->GetProcessName()
930  << " " << em
931  << G4endl;
932  }
933  return em;
934 }
void BuildRangeTable(const G4PhysicsTable *dedxTable, G4PhysicsTable *rangeTable, G4bool isIonisation=false)
void SetDEDXTable(G4PhysicsTable *p, G4EmTableType tType)
G4PhysicsTable * InverseRangeTable() const
G4ProcessVector * GetProcessList() const
void SetIonisation(G4bool val)
G4VSubCutProducer * subcutProducer
G4bool IsIonisationProcess() const
void push_back(G4PhysicsVector *)
void AddCollaborativeProcess(G4VEnergyLossProcess *)
std::vector< G4PhysicsTable * > range_vector
G4ProcessManager * GetProcessManager() const
G4PhysicsTable * IonisationTableForSubsec() const
G4int NumberOfSubCutoffRegions() const
int G4int
Definition: G4Types.hh:78
G4PhysicsTable * BuildLambdaTable(G4EmTableType tType=fRestricted)
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
void SetInverseRangeTable(G4PhysicsTable *p)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
G4PhysicsTable * RangeTableForLoss() const
void BuildDEDXTable(G4PhysicsTable *dedxTable, const std::vector< G4PhysicsTable *> &)
bool G4bool
Definition: G4Types.hh:79
G4int size() const
void BuildInverseRangeTable(const G4PhysicsTable *rangeTable, G4PhysicsTable *invRangeTable, G4bool isIonisation=false)
std::vector< G4PhysicsTable * > inv_range_vector
G4EmParameters * theParameters
void SetLambdaTable(G4PhysicsTable *p)
G4PhysicsTable * DEDXunRestrictedTable() const
G4LossTableBuilder * tableBuilder
G4PhysicsTable * DEDXTable() const
std::vector< G4PhysicsTable * > dedx_vector
std::vector< G4VEnergyLossProcess * > loss_vector
void SetCSDARangeTable(G4PhysicsTable *pRange)
std::vector< PD > part_vector
#define G4endl
Definition: G4ios.hh:61
std::vector< G4bool > isActive
G4PhysicsTable * CSDARangeTable() const
void SetSubLambdaTable(G4PhysicsTable *p)
void SetRangeTableForLoss(G4PhysicsTable *p)
G4PhysicsTable * BuildDEDXTable(G4EmTableType tType=fRestricted)
std::vector< G4bool > tables_are_built
G4bool BuildCSDARange() const
const G4String & GetName() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Clear()

void G4LossTableManager::Clear ( )

Definition at line 204 of file G4LossTableManager.cc.

205 {
206  all_tables_are_built = false;
207  currentLoss = nullptr;
208  currentParticle = nullptr;
209  if(n_loss)
210  {
211  dedx_vector.clear();
212  range_vector.clear();
213  inv_range_vector.clear();
214  loss_map.clear();
215  loss_vector.clear();
216  part_vector.clear();
217  base_part_vector.clear();
218  tables_are_built.clear();
219  isActive.clear();
220  n_loss = 0;
221  }
222 }
std::vector< G4PhysicsTable * > range_vector
std::map< PD, G4VEnergyLossProcess *, std::less< PD > > loss_map
G4VEnergyLossProcess * currentLoss
std::vector< G4PhysicsTable * > inv_range_vector
std::vector< G4PhysicsTable * > dedx_vector
std::vector< G4VEnergyLossProcess * > loss_vector
std::vector< PD > part_vector
std::vector< G4bool > isActive
std::vector< PD > base_part_vector
std::vector< G4bool > tables_are_built
Here is the caller graph for this function:

◆ CopyDEDXTables()

void G4LossTableManager::CopyDEDXTables ( )
private

◆ CopyTables()

void G4LossTableManager::CopyTables ( const G4ParticleDefinition aParticle,
G4VEnergyLossProcess base_proc 
)
private

Definition at line 711 of file G4LossTableManager.cc.

713 {
714  for (G4int j=0; j<n_loss; ++j) {
715 
717 
718  if (!tables_are_built[j] && part == base_part_vector[j]) {
719  tables_are_built[j] = true;
720  proc->SetDEDXTable(base_proc->IonisationTable(),fRestricted);
721  proc->SetDEDXTable(base_proc->DEDXTableForSubsec(),fSubRestricted);
722  proc->SetDEDXTable(base_proc->DEDXunRestrictedTable(),fTotal);
723  proc->SetCSDARangeTable(base_proc->CSDARangeTable());
724  proc->SetRangeTableForLoss(base_proc->RangeTableForLoss());
725  proc->SetInverseRangeTable(base_proc->InverseRangeTable());
726  proc->SetLambdaTable(base_proc->LambdaTable());
727  proc->SetSubLambdaTable(base_proc->SubLambdaTable());
728  proc->SetIonisation(base_proc->IsIonisationProcess());
729  if(proc->IsIonisationProcess()) {
730  range_vector[j] = base_proc->RangeTableForLoss();
731  inv_range_vector[j] = base_proc->InverseRangeTable();
732  loss_map[part_vector[j]] = proc;
733  //G4cout << "G4LossTableManager::CopyTable " << part_vector[j]->GetParticleName()
734  // << " added to map " << proc << G4endl;
735  }
736  if (1 < verbose) {
737  G4cout << "For " << proc->GetProcessName()
738  << " for " << part_vector[j]->GetParticleName()
739  << " base_part= " << part->GetParticleName()
740  << " tables are assigned"
741  << G4endl;
742  }
743  }
744 
745  if (theElectron == part && theElectron == proc->SecondaryParticle() ) {
746  proc->SetSecondaryRangeTable(base_proc->RangeTableForLoss());
747  }
748  }
749 }
void SetDEDXTable(G4PhysicsTable *p, G4EmTableType tType)
G4PhysicsTable * InverseRangeTable() const
void SetIonisation(G4bool val)
G4bool IsIonisationProcess() const
G4PhysicsTable * DEDXTableForSubsec() const
std::vector< G4PhysicsTable * > range_vector
int G4int
Definition: G4Types.hh:78
std::map< PD, G4VEnergyLossProcess *, std::less< PD > > loss_map
G4PhysicsTable * LambdaTable() const
void SetInverseRangeTable(G4PhysicsTable *p)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
G4PhysicsTable * RangeTableForLoss() const
TString part[npart]
std::vector< G4PhysicsTable * > inv_range_vector
void SetLambdaTable(G4PhysicsTable *p)
G4PhysicsTable * DEDXunRestrictedTable() const
const G4ParticleDefinition * SecondaryParticle() const
std::vector< G4VEnergyLossProcess * > loss_vector
void SetCSDARangeTable(G4PhysicsTable *pRange)
std::vector< PD > part_vector
#define G4endl
Definition: G4ios.hh:61
void SetSecondaryRangeTable(G4PhysicsTable *p)
G4PhysicsTable * CSDARangeTable() const
G4PhysicsTable * SubLambdaTable() const
G4PhysicsTable * IonisationTable() const
void SetSubLambdaTable(G4PhysicsTable *p)
std::vector< PD > base_part_vector
void SetRangeTableForLoss(G4PhysicsTable *p)
std::vector< G4bool > tables_are_built
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DeRegister() [1/5]

void G4LossTableManager::DeRegister ( G4VEnergyLossProcess p)

Definition at line 274 of file G4LossTableManager.cc.

275 {
276  if(!p) { return; }
277  for (G4int i=0; i<n_loss; ++i) {
278  if(loss_vector[i] == p) { loss_vector[i] = nullptr; }
279  }
280 }
int G4int
Definition: G4Types.hh:78
std::vector< G4VEnergyLossProcess * > loss_vector
Here is the caller graph for this function:

◆ DeRegister() [2/5]

void G4LossTableManager::DeRegister ( G4VMultipleScattering p)

Definition at line 300 of file G4LossTableManager.cc.

301 {
302  if(!p) { return; }
303  size_t msc = msc_vector.size();
304  for (size_t i=0; i<msc; ++i) {
305  if(msc_vector[i] == p) { msc_vector[i] = nullptr; }
306  }
307 }
std::vector< G4VMultipleScattering * > msc_vector

◆ DeRegister() [3/5]

void G4LossTableManager::DeRegister ( G4VEmProcess p)

Definition at line 327 of file G4LossTableManager.cc.

328 {
329  if(!p) { return; }
330  size_t emp = emp_vector.size();
331  for (size_t i=0; i<emp; ++i) {
332  if(emp_vector[i] == p) { emp_vector[i] = nullptr; }
333  }
334 }
std::vector< G4VEmProcess * > emp_vector

◆ DeRegister() [4/5]

void G4LossTableManager::DeRegister ( G4VEmModel p)

Definition at line 349 of file G4LossTableManager.cc.

350 {
351  size_t n = mod_vector.size();
352  for (size_t i=0; i<n; ++i) {
353  if(mod_vector[i] == p) { mod_vector[i] = nullptr; }
354  }
355 }
std::vector< G4VEmModel * > mod_vector
Char_t n[5]

◆ DeRegister() [5/5]

void G4LossTableManager::DeRegister ( G4VEmFluctuationModel p)

Definition at line 370 of file G4LossTableManager.cc.

371 {
372  size_t n = fmod_vector.size();
373  for (size_t i=0; i<n; ++i) {
374  if(fmod_vector[i] == p) { fmod_vector[i] = nullptr; }
375  }
376 }
Char_t n[5]
std::vector< G4VEmFluctuationModel * > fmod_vector

◆ ElectronIonPair()

G4ElectronIonPair * G4LossTableManager::ElectronIonPair ( )

Definition at line 1136 of file G4LossTableManager.cc.

1137 {
1138  if(!emElectronIonPair) {
1140  }
1141  return emElectronIonPair;
1142 }
G4ElectronIonPair * emElectronIonPair

◆ EmConfigurator()

G4EmConfigurator * G4LossTableManager::EmConfigurator ( )

Definition at line 1128 of file G4LossTableManager.cc.

1129 {
1131  return emConfigurator;
1132 }
G4EmConfigurator * emConfigurator
Here is the caller graph for this function:

◆ EmCorrections()

G4EmCorrections * G4LossTableManager::EmCorrections ( )
inline

Definition at line 521 of file G4LossTableManager.hh.

522 {
523  return emCorrections;
524 }
G4EmCorrections * emCorrections
Here is the caller graph for this function:

◆ EmSaturation()

G4EmSaturation * G4LossTableManager::EmSaturation ( )

Definition at line 1120 of file G4LossTableManager.cc.

1121 {
1123  return emSaturation;
1124 }
G4EmSaturation * emSaturation
Here is the caller graph for this function:

◆ GetCSDARange()

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

Definition at line 417 of file G4LossTableManager.hh.

420 {
421  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
422  G4double x = DBL_MAX;
423  if(currentLoss) { x = currentLoss->GetCSDARange(kineticEnergy, couple); }
424  return x;
425 }
G4double GetCSDARange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4VEnergyLossProcess * currentLoss
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetDEDX()

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

Definition at line 391 of file G4LossTableManager.hh.

394 {
395  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
396  G4double x = 0.0;
397  if(currentLoss) { x = currentLoss->GetDEDX(kineticEnergy, couple); }
398  return x;
399 }
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4VEnergyLossProcess * currentLoss
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetDEDXDispersion()

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

Definition at line 472 of file G4LossTableManager.hh.

476 {
477  const G4ParticleDefinition* aParticle = dp->GetParticleDefinition();
478  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
479  G4double x = 0.0;
480  if(currentLoss) { currentLoss->GetDEDXDispersion(couple, dp, length); }
481  return x;
482 }
G4VEnergyLossProcess * currentLoss
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)
const G4ParticleDefinition * GetParticleDefinition() const
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ GetEmProcessVector()

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

Definition at line 1105 of file G4LossTableManager.cc.

1106 {
1107  return emp_vector;
1108 }
std::vector< G4VEmProcess * > emp_vector
Here is the caller graph for this function:

◆ GetEnergy()

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

Definition at line 459 of file G4LossTableManager.hh.

462 {
463  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
464  G4double x = 0;
465  if(currentLoss) { x = currentLoss->GetKineticEnergy(range, couple); }
466  return x;
467 }
G4double GetKineticEnergy(G4double &range, const G4MaterialCutsCouple *)
G4VEnergyLossProcess * currentLoss
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetEnergyLossProcess()

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

Definition at line 369 of file G4LossTableManager.hh.

370 {
371  //G4cout << "G4LossTableManager::GetEnergyLossProcess: "
372  // << aParticle << " " << currentParticle << " " << currentLoss << G4endl;
373  if(aParticle != currentParticle) {
374  currentParticle = aParticle;
375  std::map<PD,G4VEnergyLossProcess*,std::less<PD> >::const_iterator pos;
376  if ((pos = loss_map.find(aParticle)) != loss_map.end()) {
377  currentLoss = (*pos).second;
378  } else {
379  currentLoss = 0;
380  if ((pos = loss_map.find(theGenericIon)) != loss_map.end()) {
381  currentLoss = (*pos).second;
382  }
383  }
384  }
385  return currentLoss;
386 }
std::map< PD, G4VEnergyLossProcess *, std::less< PD > > loss_map
G4VEnergyLossProcess * currentLoss
static const G4double pos
Here is the caller graph for this function:

◆ GetEnergyLossProcessVector()

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

Definition at line 1098 of file G4LossTableManager.cc.

1099 {
1100  return loss_vector;
1101 }
std::vector< G4VEnergyLossProcess * > loss_vector
Here is the caller graph for this function:

◆ GetMessenger()

G4EnergyLossMessenger * G4LossTableManager::GetMessenger ( )

Definition at line 938 of file G4LossTableManager.cc.

939 {
940  return theMessenger;
941 }
G4EnergyLossMessenger * theMessenger

◆ GetMultipleScatteringVector()

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

Definition at line 1113 of file G4LossTableManager.cc.

1114 {
1115  return msc_vector;
1116 }
std::vector< G4VMultipleScattering * > msc_vector
Here is the caller graph for this function:

◆ GetNumberOfBinsPerDecade()

G4int G4LossTableManager::GetNumberOfBinsPerDecade ( ) const
inline

Definition at line 514 of file G4LossTableManager.hh.

515 {
517 }
G4int NumberOfBinsPerDecade() const
G4EmParameters * theParameters
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetRange()

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

Definition at line 444 of file G4LossTableManager.hh.

447 {
448  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
449  G4double x = DBL_MAX;
450  if(currentLoss) {
451  x = currentLoss->GetRange(kineticEnergy, couple);
452  }
453  return x;
454 }
G4VEnergyLossProcess * currentLoss
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
G4double GetRange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetRangeFromRestricteDEDX()

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

Definition at line 430 of file G4LossTableManager.hh.

434 {
435  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
436  G4double x = DBL_MAX;
437  if(currentLoss) { x = currentLoss->GetRangeForLoss(kineticEnergy, couple); }
438  return x;
439 }
G4double GetRangeForLoss(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4VEnergyLossProcess * currentLoss
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetSubDEDX()

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

Definition at line 404 of file G4LossTableManager.hh.

407 {
408  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
409  G4double x = 0.0;
410  if(currentLoss) { x = currentLoss->GetDEDXForSubsec(kineticEnergy, couple); }
411  return x;
412 }
G4double GetDEDXForSubsec(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4VEnergyLossProcess * currentLoss
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ GetTableBuilder()

G4LossTableBuilder * G4LossTableManager::GetTableBuilder ( )
inline

Definition at line 542 of file G4LossTableManager.hh.

543 {
544  return tableBuilder;
545 }
G4LossTableBuilder * tableBuilder
Here is the caller graph for this function:

◆ Instance()

G4LossTableManager * G4LossTableManager::Instance ( void  )
static

Definition at line 112 of file G4LossTableManager.cc.

113 {
114  if(!instance) {
116  instance = inst.Instance();
117  }
118  return instance;
119 }
static G4ThreadLocal G4LossTableManager * instance
Here is the call graph for this function:

◆ IsMaster()

G4bool G4LossTableManager::IsMaster ( ) const
inline

Definition at line 493 of file G4LossTableManager.hh.

494 {
495  return isMaster;
496 }
Here is the caller graph for this function:

◆ LocalPhysicsTables()

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

Definition at line 515 of file G4LossTableManager.cc.

518 {
519  if(1 < verbose) {
520  G4cout << "### G4LossTableManager::LocalPhysicsTable() for "
521  << aParticle->GetParticleName()
522  << " and process " << p->GetProcessName()
523  << G4endl;
524  }
525 
526  if(-1 == run && startInitialisation) {
528  firstParticle = aParticle;
529  }
530 
531  if(startInitialisation) {
532  ++run;
534  if(1 < verbose) {
535  G4cout << "===== G4LossTableManager::LocalPhysicsTable() for run "
536  << run << " =====" << G4endl;
537  }
539  currentParticle = nullptr;
540  startInitialisation = false;
541  for (G4int i=0; i<n_loss; ++i) {
542  if(loss_vector[i]) {
543  tables_are_built[i] = false;
544  } else {
545  tables_are_built[i] = true;
546  part_vector[i] = nullptr;
547  }
548  }
549  }
550 
551  all_tables_are_built= true;
552  for (G4int i=0; i<n_loss; ++i) {
553  if(p == loss_vector[i]) {
554  tables_are_built[i] = true;
555  isActive[i] = true;
556  part_vector[i] = p->Particle();
557  base_part_vector[i] = p->BaseParticle();
558  dedx_vector[i] = p->DEDXTable();
561  if(0 == run && p->IsIonisationProcess()) {
562  loss_map[part_vector[i]] = p;
563  //G4cout << "G4LossTableManager::LocalPhysicsTable " << part_vector[i]->GetParticleName()
564  // << " added to map " << p << G4endl;
565  }
566 
567  if(1 < verbose) {
568  G4cout << i <<". "<< p->GetProcessName();
569  if(part_vector[i]) {
570  G4cout << " for " << part_vector[i]->GetParticleName();
571  }
572  G4cout << " active= " << isActive[i]
573  << " table= " << tables_are_built[i]
574  << " isIonisation= " << p->IsIonisationProcess()
575  << G4endl;
576  }
577  break;
578  } else if(!tables_are_built[i]) {
579  all_tables_are_built = false;
580  }
581  }
582 
583  // Set run time parameters
584  SetParameters(aParticle, p);
585 
586  if(1 < verbose) {
587  G4cout << "### G4LossTableManager::LocalPhysicsTable end"
588  << G4endl;
589  }
590  if(all_tables_are_built) {
591  if(1 < verbose) {
592  G4cout << "%%%%% All dEdx and Range tables for worker are ready for run "
593  << run << " %%%%%" << G4endl;
594  }
595  }
596 }
G4PhysicsTable * InverseRangeTable() const
G4bool IsIonisationProcess() const
G4EmSaturation * emSaturation
const G4ParticleDefinition * Particle() const
std::vector< G4PhysicsTable * > range_vector
int G4int
Definition: G4Types.hh:78
std::map< PD, G4VEnergyLossProcess *, std::less< PD > > loss_map
G4EmConfigurator * emConfigurator
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void SetParameters(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *)
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
G4PhysicsTable * RangeTableForLoss() const
std::vector< G4PhysicsTable * > inv_range_vector
void DumpG4BirksCoefficients()
G4PhysicsTable * DEDXTable() const
std::vector< G4PhysicsTable * > dedx_vector
std::vector< G4VEnergyLossProcess * > loss_vector
void SetVerbose(G4int val)
std::vector< PD > part_vector
#define G4endl
Definition: G4ios.hh:61
std::vector< G4bool > isActive
std::vector< PD > base_part_vector
const G4ParticleDefinition * BaseParticle() const
std::vector< G4bool > tables_are_built
Here is the call graph for this function:
Here is the caller graph for this function:

◆ MaxKinEnergy()

G4double G4LossTableManager::MaxKinEnergy ( ) const
inline

Definition at line 507 of file G4LossTableManager.hh.

508 {
509  return theParameters->MaxKinEnergy();
510 }
G4EmParameters * theParameters
G4double MaxKinEnergy() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ MinKinEnergy()

G4double G4LossTableManager::MinKinEnergy ( ) const
inline

Definition at line 500 of file G4LossTableManager.hh.

501 {
502  return theParameters->MinKinEnergy();
503 }
G4double MinKinEnergy() const
G4EmParameters * theParameters
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

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

◆ ParticleHaveNoLoss()

void G4LossTableManager::ParticleHaveNoLoss ( const G4ParticleDefinition aParticle)
private

Definition at line 945 of file G4LossTableManager.cc.

947 {
949  ed << "Energy loss process not found for " << aParticle->GetParticleName()
950  << " !";
951  G4Exception("G4LossTableManager::ParticleHaveNoLoss", "em0001",
952  FatalException, ed);
953 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const G4String & GetParticleName() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
Here is the call graph for this function:

◆ PreparePhysicsTable() [1/3]

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

Definition at line 407 of file G4LossTableManager.cc.

410 {
411  if (1 < verbose) {
412  G4cout << "G4LossTableManager::PreparePhysicsTable for "
413  << particle->GetParticleName()
414  << " and " << p->GetProcessName() << " run= " << run
415  << " loss_vector " << loss_vector.size() << G4endl;
416  }
417 
418  isMaster = theMaster;
419 
420  if(!startInitialisation) {
421  ResetParameters();
422  if (1 < verbose) {
423  G4cout << "====== G4LossTableManager::PreparePhysicsTable start ====="
424  << G4endl;
425  }
426  }
427 
428  // start initialisation for the first run
429  if( -1 == run ) {
430  if(emConfigurator) { emConfigurator->PrepareModels(particle, p); }
431 
432  // initialise particles for given process
433  for (G4int j=0; j<n_loss; ++j) {
434  if (p == loss_vector[j] && !part_vector[j]) {
435  part_vector[j] = particle;
436  if(particle->GetParticleName() == "GenericIon") {
437  theGenericIon = particle;
438  }
439  }
440  }
441  }
442  startInitialisation = true;
443 }
int G4int
Definition: G4Types.hh:78
G4EmConfigurator * emConfigurator
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
void PrepareModels(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
std::vector< G4VEnergyLossProcess * > loss_vector
std::vector< PD > part_vector
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:
Here is the caller graph for this function:

◆ PreparePhysicsTable() [2/3]

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

Definition at line 448 of file G4LossTableManager.cc.

450 {
451  if (1 < verbose) {
452  G4cout << "G4LossTableManager::PreparePhysicsTable for "
453  << particle->GetParticleName()
454  << " and " << p->GetProcessName() << G4endl;
455  }
456  isMaster = theMaster;
457 
458  if(!startInitialisation) {
459  ResetParameters();
460  if (1 < verbose) {
461  G4cout << "====== G4LossTableManager::PreparePhysicsTable start ====="
462  << G4endl;
463  }
464  }
465 
466  // start initialisation for the first run
467  if( -1 == run ) {
468  if(emConfigurator) { emConfigurator->PrepareModels(particle, p); }
469  }
470  startInitialisation = true;
471 }
G4EmConfigurator * emConfigurator
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
void PrepareModels(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:

◆ PreparePhysicsTable() [3/3]

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

Definition at line 476 of file G4LossTableManager.cc.

479 {
480  if (1 < verbose) {
481  G4cout << "G4LossTableManager::PreparePhysicsTable for "
482  << particle->GetParticleName()
483  << " and " << p->GetProcessName() << G4endl;
484  }
485 
486  isMaster = theMaster;
487 
488  if(!startInitialisation) {
489  ResetParameters();
490  if (1 < verbose) {
491  G4cout << "====== G4LossTableManager::PreparePhysicsTable start ====="
492  << G4endl;
493  }
494  }
495 
496  // start initialisation for the first run
497  if( -1 == run ) {
498  if(emConfigurator) { emConfigurator->PrepareModels(particle, p); }
499  }
500  startInitialisation = true;
501 }
G4EmConfigurator * emConfigurator
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
void PrepareModels(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:

◆ PrintEWarning()

void G4LossTableManager::PrintEWarning ( G4String  tit,
G4double   
)
private

Definition at line 1166 of file G4LossTableManager.cc.

1167 {
1168  G4String ss = "G4LossTableManager::" + tit;
1170  /*
1171  ed << "Parameter is out of range: " << val
1172  << " it will have no effect!\n" << " ## "
1173  << " nbins= " << nbinsLambda
1174  << " nbinsPerDecade= " << nbinsPerDecade
1175  << " Emin(keV)= " << minKinEnergy/keV
1176  << " Emax(GeV)= " << maxKinEnergy/GeV;
1177  */
1178  G4Exception(ss, "em0044", JustWarning, ed);
1179 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Register() [1/5]

void G4LossTableManager::Register ( G4VEnergyLossProcess p)

Definition at line 226 of file G4LossTableManager.cc.

227 {
228  if(!p) { return; }
229  for (G4int i=0; i<n_loss; ++i) {
230  if(loss_vector[i] == p) { return; }
231  }
232  if(verbose > 1) {
233  G4cout << "G4LossTableManager::Register G4VEnergyLossProcess : "
234  << p->GetProcessName() << " idx= " << n_loss << G4endl;
235  }
236  ++n_loss;
237  loss_vector.push_back(p);
238  part_vector.push_back(nullptr);
239  base_part_vector.push_back(nullptr);
240  dedx_vector.push_back(nullptr);
241  range_vector.push_back(nullptr);
242  inv_range_vector.push_back(nullptr);
243  tables_are_built.push_back(false);
244  isActive.push_back(true);
245  all_tables_are_built = false;
246  if(subCutoffFlag) { p->ActivateSubCutoff(true); }
248  maxFinalStep); }
249  if(integralActive) { p->SetIntegral(integral); }
250 }
void SetIntegral(G4bool val)
void SetStepFunction(G4double v1, G4double v2)
std::vector< G4PhysicsTable * > range_vector
int G4int
Definition: G4Types.hh:78
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
std::vector< G4PhysicsTable * > inv_range_vector
std::vector< G4PhysicsTable * > dedx_vector
std::vector< G4VEnergyLossProcess * > loss_vector
std::vector< PD > part_vector
#define G4endl
Definition: G4ios.hh:61
std::vector< G4bool > isActive
std::vector< PD > base_part_vector
void ActivateSubCutoff(G4bool val, const G4Region *region=0)
std::vector< G4bool > tables_are_built
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Register() [2/5]

void G4LossTableManager::Register ( G4VMultipleScattering p)

Definition at line 284 of file G4LossTableManager.cc.

285 {
286  if(!p) { return; }
287  G4int n = msc_vector.size();
288  for (G4int i=0; i<n; ++i) {
289  if(msc_vector[i] == p) { return; }
290  }
291  if(verbose > 1) {
292  G4cout << "G4LossTableManager::Register G4VMultipleScattering : "
293  << p->GetProcessName() << " idx= " << msc_vector.size() << G4endl;
294  }
295  msc_vector.push_back(p);
296 }
int G4int
Definition: G4Types.hh:78
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
Char_t n[5]
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
std::vector< G4VMultipleScattering * > msc_vector
Here is the call graph for this function:

◆ Register() [3/5]

void G4LossTableManager::Register ( G4VEmProcess p)

Definition at line 311 of file G4LossTableManager.cc.

312 {
313  if(!p) { return; }
314  G4int n = emp_vector.size();
315  for (G4int i=0; i<n; ++i) {
316  if(emp_vector[i] == p) { return; }
317  }
318  if(verbose > 1) {
319  G4cout << "G4LossTableManager::Register G4VEmProcess : "
320  << p->GetProcessName() << " idx= " << emp_vector.size() << G4endl;
321  }
322  emp_vector.push_back(p);
323 }
int G4int
Definition: G4Types.hh:78
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
Char_t n[5]
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
std::vector< G4VEmProcess * > emp_vector
Here is the call graph for this function:

◆ Register() [4/5]

void G4LossTableManager::Register ( G4VEmModel p)

Definition at line 338 of file G4LossTableManager.cc.

339 {
340  mod_vector.push_back(p);
341  if(verbose > 1) {
342  G4cout << "G4LossTableManager::Register G4VEmModel : "
343  << p->GetName() << G4endl;
344  }
345 }
const G4String & GetName() const
Definition: G4VEmModel.hh:795
std::vector< G4VEmModel * > mod_vector
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:

◆ Register() [5/5]

void G4LossTableManager::Register ( G4VEmFluctuationModel p)

Definition at line 359 of file G4LossTableManager.cc.

360 {
361  fmod_vector.push_back(p);
362  if(verbose > 1) {
363  G4cout << "G4LossTableManager::Register G4VEmFluctuationModel : "
364  << p->GetName() << G4endl;
365  }
366 }
const G4String & GetName() const
G4GLOB_DLL std::ostream G4cout
std::vector< G4VEmFluctuationModel * > fmod_vector
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:

◆ RegisterExtraParticle()

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

Definition at line 380 of file G4LossTableManager.cc.

383 {
384  if(!p || !part) { return; }
385  for (G4int i=0; i<n_loss; ++i) {
386  if(loss_vector[i] == p) { return; }
387  }
388  if(verbose > 1) {
389  G4cout << "G4LossTableManager::RegisterExtraParticle "
390  << part->GetParticleName() << " G4VEnergyLossProcess : "
391  << p->GetProcessName() << " idx= " << n_loss << G4endl;
392  }
393  ++n_loss;
394  loss_vector.push_back(p);
395  part_vector.push_back(part);
396  base_part_vector.push_back(p->BaseParticle());
397  dedx_vector.push_back(nullptr);
398  range_vector.push_back(nullptr);
399  inv_range_vector.push_back(nullptr);
400  tables_are_built.push_back(false);
401  all_tables_are_built = false;
402 }
std::vector< G4PhysicsTable * > range_vector
int G4int
Definition: G4Types.hh:78
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
TString part[npart]
std::vector< G4PhysicsTable * > inv_range_vector
std::vector< G4PhysicsTable * > dedx_vector
std::vector< G4VEnergyLossProcess * > loss_vector
std::vector< PD > part_vector
#define G4endl
Definition: G4ios.hh:61
std::vector< PD > base_part_vector
const G4ParticleDefinition * BaseParticle() const
std::vector< G4bool > tables_are_built
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ResetParameters()

void G4LossTableManager::ResetParameters ( )
private

Definition at line 254 of file G4LossTableManager.cc.

255 {
257  if(!isMaster) {
259  }
266  if(atomDeexcitation) {
269  }
270 }
G4VAtomDeexcitation * atomDeexcitation
void SetVerbose(G4int verb)
G4ElectronIonPair * emElectronIonPair
G4EmSaturation * emSaturation
void SetVerbose(G4int)
G4int Verbose() const
G4EmConfigurator * emConfigurator
void SetInitialisationFlag(G4bool flag)
G4bool Spline() const
void SetVerbose(G4int value)
G4EmParameters * theParameters
G4EmCorrections * emCorrections
G4LossTableBuilder * tableBuilder
G4int WorkerVerbose() const
void SetSplineFlag(G4bool flag)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetAtomDeexcitation()

void G4LossTableManager::SetAtomDeexcitation ( G4VAtomDeexcitation p)

Definition at line 1146 of file G4LossTableManager.cc.

1147 {
1148  if(atomDeexcitation != p) {
1149  delete atomDeexcitation;
1150  atomDeexcitation = p;
1151  }
1152 }
G4VAtomDeexcitation * atomDeexcitation
Here is the caller graph for this function:

◆ SetBuildCSDARange()

void G4LossTableManager::SetBuildCSDARange ( G4bool  val)

Definition at line 1078 of file G4LossTableManager.cc.

1079 {
1080 }
Here is the caller graph for this function:

◆ SetDEDXBinning()

void G4LossTableManager::SetDEDXBinning ( G4int  val)

Definition at line 1028 of file G4LossTableManager.cc.

1029 {
1031 }
void SetNumberOfBins(G4int val)
G4EmParameters * theParameters
Here is the call graph for this function:

◆ SetDEDXBinningForCSDARange()

void G4LossTableManager::SetDEDXBinningForCSDARange ( G4int  val)

Definition at line 1035 of file G4LossTableManager.cc.

1036 {
1037 }

◆ SetIntegral()

void G4LossTableManager::SetIntegral ( G4bool  val)

Definition at line 973 of file G4LossTableManager.cc.

974 {
975  integral = val;
976  integralActive = true;
977  for(G4int i=0; i<n_loss; ++i) {
978  if(loss_vector[i]) { loss_vector[i]->SetIntegral(val); }
979  }
980  size_t emp = emp_vector.size();
981  for (size_t k=0; k<emp; ++k) {
982  if(emp_vector[k]) { emp_vector[k]->SetIntegral(val); }
983  }
984 }
int G4int
Definition: G4Types.hh:78
std::vector< G4VEnergyLossProcess * > loss_vector
std::vector< G4VEmProcess * > emp_vector
Here is the caller graph for this function:

◆ SetLambdaBinning()

void G4LossTableManager::SetLambdaBinning ( G4int  val)

Definition at line 1041 of file G4LossTableManager.cc.

1042 {
1044 }
void SetNumberOfBins(G4int val)
G4EmParameters * theParameters
Here is the call graph for this function:

◆ SetLinearLossLimit()

void G4LossTableManager::SetLinearLossLimit ( G4double  val)

Definition at line 1072 of file G4LossTableManager.cc.

1073 {
1074 }
Here is the caller graph for this function:

◆ SetLossFluctuations()

void G4LossTableManager::SetLossFluctuations ( G4bool  val)

Definition at line 957 of file G4LossTableManager.cc.

958 {
959 }
Here is the caller graph for this function:

◆ SetMaxEnergy()

void G4LossTableManager::SetMaxEnergy ( G4double  val)

Definition at line 1007 of file G4LossTableManager.cc.

1008 {
1010 }
void SetMaxEnergy(G4double val)
G4EmParameters * theParameters
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetMaxEnergyForCSDARange()

void G4LossTableManager::SetMaxEnergyForCSDARange ( G4double  val)

Definition at line 1014 of file G4LossTableManager.cc.

1015 {
1017 }
void SetMaxEnergyForCSDARange(G4double val)
G4EmParameters * theParameters
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetMaxEnergyForMuons()

void G4LossTableManager::SetMaxEnergyForMuons ( G4double  val)

Definition at line 1021 of file G4LossTableManager.cc.

1022 {
1024 }
void SetMaxEnergy(G4double val)
G4EmParameters * theParameters
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetMinEnergy()

void G4LossTableManager::SetMinEnergy ( G4double  val)

Definition at line 1000 of file G4LossTableManager.cc.

1001 {
1003 }
G4EmParameters * theParameters
void SetMinEnergy(G4double val)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetMinSubRange()

void G4LossTableManager::SetMinSubRange ( G4double  val)

Definition at line 988 of file G4LossTableManager.cc.

989 {
990 }
Here is the caller graph for this function:

◆ SetParameters()

void G4LossTableManager::SetParameters ( const G4ParticleDefinition aParticle,
G4VEnergyLossProcess p 
)
private

Definition at line 1085 of file G4LossTableManager.cc.

1087 {
1088  if(stepFunctionActive) {
1090  }
1091  if(integralActive) { p->SetIntegral(integral); }
1092 }
void SetIntegral(G4bool val)
void SetStepFunction(G4double v1, G4double v2)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetRandomStep()

void G4LossTableManager::SetRandomStep ( G4bool  val)

Definition at line 994 of file G4LossTableManager.cc.

995 {
996 }
Here is the caller graph for this function:

◆ SetStepFunction()

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

Definition at line 1054 of file G4LossTableManager.cc.

1055 {
1056  if(0.0 < v1 && 0.0 < v2 && v2 < 1.e+50) {
1057  stepFunctionActive = true;
1058  maxRangeVariation = v1;
1059  maxFinalStep = v2;
1060  for(G4int i=0; i<n_loss; ++i) {
1061  if(loss_vector[i]) { loss_vector[i]->SetStepFunction(v1, v2); }
1062  }
1063  } else if(v1 <= 0.0) {
1064  PrintEWarning("SetStepFunction", v1);
1065  } else {
1066  PrintEWarning("SetStepFunction", v2);
1067  }
1068 }
int G4int
Definition: G4Types.hh:78
std::vector< G4VEnergyLossProcess * > loss_vector
void PrintEWarning(G4String, G4double)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetSubCutoff()

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

Definition at line 963 of file G4LossTableManager.cc.

964 {
965  subCutoffFlag = val;
966  for(G4int i=0; i<n_loss; ++i) {
967  if(loss_vector[i]) { loss_vector[i]->ActivateSubCutoff(val, r); }
968  }
969 }
int G4int
Definition: G4Types.hh:78
std::vector< G4VEnergyLossProcess * > loss_vector
Here is the caller graph for this function:

◆ SetSubCutProducer()

void G4LossTableManager::SetSubCutProducer ( G4VSubCutProducer p)

Definition at line 1156 of file G4LossTableManager.cc.

1157 {
1158  if(subcutProducer != p) {
1159  delete subcutProducer;
1160  subcutProducer = p;
1161  }
1162 }
G4VSubCutProducer * subcutProducer

◆ SetVerbose()

void G4LossTableManager::SetVerbose ( G4int  val)

Definition at line 1048 of file G4LossTableManager.cc.

1049 {
1050 }
Here is the caller graph for this function:

◆ SubCutProducer()

G4VSubCutProducer * G4LossTableManager::SubCutProducer ( )
inline

Definition at line 535 of file G4LossTableManager.hh.

536 {
537  return subcutProducer;
538 }
G4VSubCutProducer * subcutProducer
Here is the caller graph for this function:

Friends And Related Function Documentation

◆ G4ThreadLocalSingleton< G4LossTableManager >

Definition at line 105 of file G4LossTableManager.hh.

Member Data Documentation

◆ all_tables_are_built

G4bool G4LossTableManager::all_tables_are_built
private

Definition at line 338 of file G4LossTableManager.hh.

◆ atomDeexcitation

G4VAtomDeexcitation* G4LossTableManager::atomDeexcitation
private

Definition at line 356 of file G4LossTableManager.hh.

◆ base_part_vector

std::vector<PD> G4LossTableManager::base_part_vector
private

Definition at line 317 of file G4LossTableManager.hh.

◆ currentLoss

G4VEnergyLossProcess* G4LossTableManager::currentLoss
private

Definition at line 329 of file G4LossTableManager.hh.

◆ currentParticle

PD G4LossTableManager::currentParticle
private

Definition at line 330 of file G4LossTableManager.hh.

◆ dedx_vector

std::vector<G4PhysicsTable*> G4LossTableManager::dedx_vector
private

Definition at line 320 of file G4LossTableManager.hh.

◆ emConfigurator

G4EmConfigurator* G4LossTableManager::emConfigurator
private

Definition at line 354 of file G4LossTableManager.hh.

◆ emCorrections

G4EmCorrections* G4LossTableManager::emCorrections
private

Definition at line 352 of file G4LossTableManager.hh.

◆ emElectronIonPair

G4ElectronIonPair* G4LossTableManager::emElectronIonPair
private

Definition at line 355 of file G4LossTableManager.hh.

◆ emp_vector

std::vector<G4VEmProcess*> G4LossTableManager::emp_vector
private

Definition at line 324 of file G4LossTableManager.hh.

◆ emSaturation

G4EmSaturation* G4LossTableManager::emSaturation
private

Definition at line 353 of file G4LossTableManager.hh.

◆ firstParticle

PD G4LossTableManager::firstParticle
private

Definition at line 333 of file G4LossTableManager.hh.

◆ fmod_vector

std::vector<G4VEmFluctuationModel*> G4LossTableManager::fmod_vector
private

Definition at line 326 of file G4LossTableManager.hh.

◆ instance

G4ThreadLocal G4LossTableManager * G4LossTableManager::instance = 0
staticprivate

Definition at line 309 of file G4LossTableManager.hh.

◆ integral

G4bool G4LossTableManager::integral
private

Definition at line 342 of file G4LossTableManager.hh.

◆ integralActive

G4bool G4LossTableManager::integralActive
private

Definition at line 343 of file G4LossTableManager.hh.

◆ inv_range_vector

std::vector<G4PhysicsTable*> G4LossTableManager::inv_range_vector
private

Definition at line 322 of file G4LossTableManager.hh.

◆ isActive

std::vector<G4bool> G4LossTableManager::isActive
private

Definition at line 319 of file G4LossTableManager.hh.

◆ isMaster

G4bool G4LossTableManager::isMaster
private

Definition at line 345 of file G4LossTableManager.hh.

◆ loss_map

std::map<PD,G4VEnergyLossProcess*,std::less<PD> > G4LossTableManager::loss_map
private

Definition at line 313 of file G4LossTableManager.hh.

◆ loss_vector

std::vector<G4VEnergyLossProcess*> G4LossTableManager::loss_vector
private

Definition at line 315 of file G4LossTableManager.hh.

◆ maxFinalStep

G4double G4LossTableManager::maxFinalStep
private

Definition at line 348 of file G4LossTableManager.hh.

◆ maxRangeVariation

G4double G4LossTableManager::maxRangeVariation
private

Definition at line 347 of file G4LossTableManager.hh.

◆ mod_vector

std::vector<G4VEmModel*> G4LossTableManager::mod_vector
private

Definition at line 325 of file G4LossTableManager.hh.

◆ msc_vector

std::vector<G4VMultipleScattering*> G4LossTableManager::msc_vector
private

Definition at line 323 of file G4LossTableManager.hh.

◆ n_loss

G4int G4LossTableManager::n_loss
private

Definition at line 335 of file G4LossTableManager.hh.

◆ part_vector

std::vector<PD> G4LossTableManager::part_vector
private

Definition at line 316 of file G4LossTableManager.hh.

◆ range_vector

std::vector<G4PhysicsTable*> G4LossTableManager::range_vector
private

Definition at line 321 of file G4LossTableManager.hh.

◆ run

G4int G4LossTableManager::run
private

Definition at line 336 of file G4LossTableManager.hh.

◆ startInitialisation

G4bool G4LossTableManager::startInitialisation
private

Definition at line 339 of file G4LossTableManager.hh.

◆ stepFunctionActive

G4bool G4LossTableManager::stepFunctionActive
private

Definition at line 344 of file G4LossTableManager.hh.

◆ subCutoffFlag

G4bool G4LossTableManager::subCutoffFlag
private

Definition at line 341 of file G4LossTableManager.hh.

◆ subcutProducer

G4VSubCutProducer* G4LossTableManager::subcutProducer
private

Definition at line 357 of file G4LossTableManager.hh.

◆ tableBuilder

G4LossTableBuilder* G4LossTableManager::tableBuilder
private

Definition at line 350 of file G4LossTableManager.hh.

◆ tables_are_built

std::vector<G4bool> G4LossTableManager::tables_are_built
private

Definition at line 318 of file G4LossTableManager.hh.

◆ theElectron

PD G4LossTableManager::theElectron
private

Definition at line 331 of file G4LossTableManager.hh.

◆ theGenericIon

PD G4LossTableManager::theGenericIon
private

Definition at line 332 of file G4LossTableManager.hh.

◆ theMessenger

G4EnergyLossMessenger* G4LossTableManager::theMessenger
private

Definition at line 351 of file G4LossTableManager.hh.

◆ theParameters

G4EmParameters* G4LossTableManager::theParameters
private

Definition at line 359 of file G4LossTableManager.hh.

◆ verbose

G4int G4LossTableManager::verbose
private

Definition at line 361 of file G4LossTableManager.hh.


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