Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VUserPhysicsList Class Referenceabstract

#include <G4VUserPhysicsList.hh>

+ Inheritance diagram for G4VUserPhysicsList:

Public Member Functions

 G4VUserPhysicsList ()
 
virtual ~G4VUserPhysicsList ()
 
 G4VUserPhysicsList (const G4VUserPhysicsList &)
 
G4VUserPhysicsListoperator= (const G4VUserPhysicsList &)
 
virtual void ConstructParticle ()=0
 
void Construct ()
 
virtual void ConstructProcess ()=0
 
virtual void SetCuts ()
 
void SetDefaultCutValue (G4double newCutValue)
 
G4double GetDefaultCutValue () const
 
void BuildPhysicsTable ()
 
void PreparePhysicsTable (G4ParticleDefinition *)
 
void BuildPhysicsTable (G4ParticleDefinition *)
 
G4bool StorePhysicsTable (const G4String &directory=".")
 
G4bool IsPhysicsTableRetrieved () const
 
G4bool IsStoredInAscii () const
 
const G4StringGetPhysicsTableDirectory () const
 
void SetPhysicsTableRetrieved (const G4String &directory="")
 
void SetStoredInAscii ()
 
void ResetPhysicsTableRetrieved ()
 
void ResetStoredInAscii ()
 
void DumpList () const
 
void DumpCutValuesTable (G4int flag=1)
 
void DumpCutValuesTableIfRequested ()
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
void UseCoupledTransportation (G4bool vl=true)
 
void SetCutsWithDefault ()
 
void SetCutValue (G4double aCut, const G4String &pname)
 
G4double GetCutValue (const G4String &pname) const
 
void SetCutValue (G4double aCut, const G4String &pname, const G4String &rname)
 
void SetParticleCuts (G4double cut, G4ParticleDefinition *particle, G4Region *region=nullptr)
 
void SetParticleCuts (G4double cut, const G4String &particleName, G4Region *region=nullptr)
 
void SetCutsForRegion (G4double aCut, const G4String &rname)
 
void SetApplyCuts (G4bool value, const G4String &name)
 
G4bool GetApplyCuts (const G4String &name) const
 
void RemoveProcessManager ()
 
void RemoveTrackingManager ()
 
void AddProcessManager (G4ParticleDefinition *newParticle, G4ProcessManager *newManager=nullptr)
 
void CheckParticleList ()
 
void DisableCheckParticleList ()
 
G4int GetInstanceID () const
 
virtual void InitializeWorker ()
 
virtual void TerminateWorker ()
 

Static Public Member Functions

static const G4VUPLManagerGetSubInstanceManager ()
 

Protected Member Functions

void AddTransportation ()
 
G4bool RegisterProcess (G4VProcess *process, G4ParticleDefinition *particle)
 
void BuildIntegralPhysicsTable (G4VProcess *, G4ParticleDefinition *)
 
virtual void RetrievePhysicsTable (G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
void InitializeProcessManager ()
 
G4ParticleTable::G4PTblDicIteratorGetParticleIterator () const
 

Protected Attributes

G4ParticleTabletheParticleTable = nullptr
 
G4int verboseLevel = 1
 
G4double defaultCutValue = 1.0
 
G4bool isSetDefaultCutValue = false
 
G4ProductionCutsTablefCutsTable = nullptr
 
G4bool fRetrievePhysicsTable = false
 
G4bool fStoredInAscii = true
 
G4bool fIsCheckedForRetrievePhysicsTable = false
 
G4bool fIsRestoredCutValues = false
 
G4String directoryPhysicsTable = "."
 
G4bool fDisableCheckParticleList = false
 
G4int g4vuplInstanceID = 0
 

Static Protected Attributes

static G4RUN_DLL G4VUPLManager subInstanceManager
 

Detailed Description

Definition at line 103 of file G4VUserPhysicsList.hh.

Constructor & Destructor Documentation

◆ G4VUserPhysicsList() [1/2]

G4VUserPhysicsList::G4VUserPhysicsList ( )

Definition at line 82 of file G4VUserPhysicsList.cc.

83{
85 // default cut value (1.0mm)
86 defaultCutValue = 1.0 * mm;
87
88 // pointer to the particle table
90
91 // pointer to the cuts table
93
94 // set energy range for SetCut calcuration
95 fCutsTable->SetEnergyRange(0.99 * keV, 100 * TeV);
96
97 // UI Messenger
99
100 // PhysicsListHelper
101 G4MT_thePLHelper->SetVerboseLevel(verboseLevel); // AND
102
103 fIsPhysicsTableBuilt = false;
105}
#define G4MT_theMessenger
#define G4MT_thePLHelper
#define fIsPhysicsTableBuilt
#define fDisplayThreshold
static G4ParticleTable * GetParticleTable()
void SetEnergyRange(G4double lowedge, G4double highedge)
static G4ProductionCutsTable * GetProductionCutsTable()
G4int CreateSubInstance()
G4ProductionCutsTable * fCutsTable
G4ParticleTable * theParticleTable
static G4RUN_DLL G4VUPLManager subInstanceManager

◆ ~G4VUserPhysicsList()

G4VUserPhysicsList::~G4VUserPhysicsList ( )
virtual

Definition at line 126 of file G4VUserPhysicsList.cc.

127{
128 delete G4MT_theMessenger;
129 G4MT_theMessenger = nullptr;
130
133
134 // invoke DeleteAllParticle
136}

◆ G4VUserPhysicsList() [2/2]

G4VUserPhysicsList::G4VUserPhysicsList ( const G4VUserPhysicsList & right)

Definition at line 139 of file G4VUserPhysicsList.cc.

140 : verboseLevel(right.verboseLevel),
149{
151 // pointer to the particle table
154 // pointer to the cuts table
156
157 // UI Messenger
159
160 // PhysicsListHelper
162 G4MT_thePLHelper->SetVerboseLevel(verboseLevel);
163
165 right.GetSubInstanceManager().offset[right.GetInstanceID()]._fIsPhysicsTableBuilt;
167 right.GetSubInstanceManager().offset[right.GetInstanceID()]._fDisplayThreshold;
168}
#define theParticleIterator
G4PTblDicIterator * GetIterator() const
static G4PhysicsListHelper * GetPhysicsListHelper()
G4RUN_DLL G4ThreadLocalStatic T * offset
static const G4VUPLManager & GetSubInstanceManager()

Member Function Documentation

◆ AddProcessManager()

void G4VUserPhysicsList::AddProcessManager ( G4ParticleDefinition * newParticle,
G4ProcessManager * newManager = nullptr )

Definition at line 192 of file G4VUserPhysicsList.cc.

193{
194 if (newParticle == nullptr) return;
195 G4Exception("G4VUserPhysicsList::AddProcessManager", "Run0252", JustWarning,
196 "This method is obsolete");
197}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)

Referenced by G4UserPhysicsListMessenger::SetNewValue().

◆ AddTransportation()

void G4VUserPhysicsList::AddTransportation ( )
protected

Definition at line 920 of file G4VUserPhysicsList.cc.

921{
922 G4MT_thePLHelper->AddTransportation();
923}

Referenced by LBE::AddTransportation(), and G4VModularPhysicsList::ConstructProcess().

◆ BuildIntegralPhysicsTable()

void G4VUserPhysicsList::BuildIntegralPhysicsTable ( G4VProcess * process,
G4ParticleDefinition * particle )
protected

Definition at line 745 of file G4VUserPhysicsList.cc.

747{
748 // TODO Should we change this function?
749 //*******************************************************************
750 // Temporary addition to make the integral schema of electromagnetic
751 // processes work.
752
753 if ((process->GetProcessName() == "Imsc") || (process->GetProcessName() == "IeIoni")
754 || (process->GetProcessName() == "IeBrems") || (process->GetProcessName() == "Iannihil")
755 || (process->GetProcessName() == "IhIoni") || (process->GetProcessName() == "IMuIoni")
756 || (process->GetProcessName() == "IMuBrems") || (process->GetProcessName() == "IMuPairProd"))
757 {
758#ifdef G4VERBOSE
759 if (verboseLevel > 2) {
760 G4cout << "G4VUserPhysicsList::BuildIntegralPhysicsTable "
761 << " BuildPhysicsTable is invoked for " << process->GetProcessName() << "("
762 << particle->GetParticleName() << ")" << G4endl;
763 }
764#endif
765 process->BuildPhysicsTable(*particle);
766 }
767}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
const G4String & GetParticleName() const
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
const G4String & GetProcessName() const

Referenced by RetrievePhysicsTable().

◆ BuildPhysicsTable() [1/2]

void G4VUserPhysicsList::BuildPhysicsTable ( )

Definition at line 501 of file G4VUserPhysicsList.cc.

502{
503 // Prepare Physics table for all particles
504 theParticleIterator->reset();
505 while ((*theParticleIterator)()) {
506 G4ParticleDefinition* particle = theParticleIterator->value();
507 PreparePhysicsTable(particle);
508 }
509
510 // ask processes to prepare physics table
513 // check if retrieve Cut Table successfully
515#ifdef G4VERBOSE
516 if (verboseLevel > 0) {
517 G4cout << "G4VUserPhysicsList::BuildPhysicsTable"
518 << " Retrieve Cut Table failed !!" << G4endl;
519 }
520#endif
521 G4Exception("G4VUserPhysicsList::BuildPhysicsTable", "Run0255", RunMustBeAborted,
522 "Fail to retrieve Production Cut Table");
523 }
524 else {
525#ifdef G4VERBOSE
526 if (verboseLevel > 2) {
527 G4cout << "G4VUserPhysicsList::BuildPhysicsTable"
528 << " Retrieve Cut Table successfully " << G4endl;
529 }
530#endif
531 }
532 }
533 else {
534#ifdef G4VERBOSE
535 if (verboseLevel > 2) {
536 G4cout << "G4VUserPhysicsList::BuildPhysicsTable"
537 << " does not retrieve Cut Table but calculate " << G4endl;
538 }
539#endif
540 }
541
542 // Sets a value to particle
543 // set cut values for gamma at first and for e- and e+
544 G4String particleName;
546 if (GammaP != nullptr) BuildPhysicsTable(GammaP);
548 if (EMinusP != nullptr) BuildPhysicsTable(EMinusP);
550 if (EPlusP != nullptr) BuildPhysicsTable(EPlusP);
552 if (ProtonP != nullptr) BuildPhysicsTable(ProtonP);
553
554 theParticleIterator->reset();
555 while ((*theParticleIterator)()) {
556 G4ParticleDefinition* particle = theParticleIterator->value();
557 if (particle != GammaP && particle != EMinusP && particle != EPlusP && particle != ProtonP) {
558 BuildPhysicsTable(particle);
559 }
560 }
561
562 // Set flag
564}
@ RunMustBeAborted
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
G4bool RetrieveCutsTable(const G4String &directory, G4bool ascii=false)
void PreparePhysicsTable(G4ParticleDefinition *)

Referenced by BuildPhysicsTable(), G4RunManagerKernel::BuildPhysicsTables(), and G4UserPhysicsListMessenger::SetNewValue().

◆ BuildPhysicsTable() [2/2]

void G4VUserPhysicsList::BuildPhysicsTable ( G4ParticleDefinition * particle)

Definition at line 567 of file G4VUserPhysicsList.cc.

568{
569 if (auto* trackingManager = particle->GetTrackingManager()) {
570#ifdef G4VERBOSE
571 if (verboseLevel > 2) {
572 G4cout << "G4VUserPhysicsList::BuildPhysicsTable "
573 << "Calculate Physics Table for " << particle->GetParticleName()
574 << " via custom TrackingManager" << G4endl;
575 }
576#endif
577 trackingManager->BuildPhysicsTable(*particle);
578 return;
579 }
580
581 // Change in order to share physics tables for two kind of process.
582
583 if (particle->GetMasterProcessManager() == nullptr) {
584#ifdef G4VERBOSE
585 if (verboseLevel > 0) {
586 G4cout << "#### G4VUserPhysicsList::BuildPhysicsTable() - BuildPhysicsTable("
587 << particle->GetParticleName() << ") skipped..." << G4endl;
588 }
589#endif
590 return;
591 }
594 // fail to retrieve cut tables
595#ifdef G4VERBOSE
596 if (verboseLevel > 0) {
597 G4cout << "G4VUserPhysicsList::BuildPhysicsTable "
598 << "Physics table can not be retrieved and will be calculated " << G4endl;
599 }
600#endif
601 fRetrievePhysicsTable = false;
602 }
603 else {
604#ifdef G4VERBOSE
605 if (verboseLevel > 2) {
606 G4cout << "G4VUserPhysicsList::BuildPhysicsTable "
607 << " Retrieve Physics Table for " << particle->GetParticleName() << G4endl;
608 }
609#endif
610 // Retrieve PhysicsTable from files for proccesses
612 }
613 }
614
615#ifdef G4VERBOSE
616 if (verboseLevel > 2) {
617 G4cout << "G4VUserPhysicsList::BuildPhysicsTable "
618 << "Calculate Physics Table for " << particle->GetParticleName() << G4endl;
619 }
620#endif
621 // Rebuild the physics tables for every process for this particle type
622 // if particle is not ShortLived
623 if (!particle->IsShortLived()) {
624 G4ProcessManager* pManager = particle->GetProcessManager();
625 if (pManager == nullptr) {
626#ifdef G4VERBOSE
627 if (verboseLevel > 0) {
628 G4cout << "G4VUserPhysicsList::BuildPhysicsTable "
629 << " : No Process Manager for " << particle->GetParticleName() << G4endl;
630 G4cout << particle->GetParticleName() << " should be created in your PhysicsList" << G4endl;
631 }
632#endif
633 G4Exception("G4VUserPhysicsList::BuildPhysicsTable", "Run0271", FatalException,
634 "No process manager");
635 return;
636 }
637
638 // Get processes from master thread;
639 G4ProcessManager* pManagerShadow = particle->GetMasterProcessManager();
640
641 G4ProcessVector* pVector = pManager->GetProcessList();
642 if (pVector == nullptr) {
643#ifdef G4VERBOSE
644 if (verboseLevel > 0) {
645 G4cout << "G4VUserPhysicsList::BuildPhysicsTable "
646 << " : No Process Vector for " << particle->GetParticleName() << G4endl;
647 }
648#endif
649 G4Exception("G4VUserPhysicsList::BuildPhysicsTable", "Run0272", FatalException,
650 "No process Vector");
651 return;
652 }
653#ifdef G4VERBOSE
654 if (verboseLevel > 2) {
655 G4cout << "G4VUserPhysicsList::BuildPhysicsTable %%%%%% " << particle->GetParticleName()
656 << G4endl;
657 G4cout << " ProcessManager : " << pManager << " ProcessManagerShadow : " << pManagerShadow
658 << G4endl;
659 for (G4int iv1 = 0; iv1 < (G4int)pVector->size(); ++iv1) {
660 G4cout << " " << iv1 << " - " << (*pVector)[iv1]->GetProcessName() << G4endl;
661 }
662 G4cout << "--------------------------------------------------------------" << G4endl;
663 G4ProcessVector* pVectorShadow = pManagerShadow->GetProcessList();
664
665 for (G4int iv2 = 0; iv2 < (G4int)pVectorShadow->size(); ++iv2) {
666 G4cout << " " << iv2 << " - " << (*pVectorShadow)[iv2]->GetProcessName() << G4endl;
667 }
668 }
669#endif
670 for (G4int j = 0; j < (G4int)pVector->size(); ++j) {
671 // Andrea July 16th 2013 : migration to new interface...
672 // Infer if we are in a worker thread or master thread
673 // Master thread is the one in which the process manager
674 // and process manager shadow pointers are the same
675 if (pManagerShadow == pManager) {
676 (*pVector)[j]->BuildPhysicsTable(*particle);
677 }
678 else {
679 (*pVector)[j]->BuildWorkerPhysicsTable(*particle);
680 }
681
682 } // End loop on processes vector
683 } // End if short-lived
684}
@ FatalException
int G4int
Definition G4Types.hh:85
G4ProcessManager * GetProcessManager() const
G4VTrackingManager * GetTrackingManager() const
G4ProcessManager * GetMasterProcessManager() const
G4ProcessVector * GetProcessList() const
std::size_t size() const
virtual void RetrievePhysicsTable(G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)

◆ CheckParticleList()

void G4VUserPhysicsList::CheckParticleList ( )

Definition at line 912 of file G4VUserPhysicsList.cc.

913{
915 G4MT_thePLHelper->CheckParticleList();
916 }
917}

Referenced by G4RunManagerKernel::InitializePhysics().

◆ Construct()

void G4VUserPhysicsList::Construct ( )
inline

Definition at line 302 of file G4VUserPhysicsList.hh.

303{
304#ifdef G4VERBOSE
305 if (verboseLevel > 1) G4cout << "G4VUserPhysicsList::Construct()" << G4endl;
306#endif
307
309
311
312#ifdef G4VERBOSE
313 if (verboseLevel > 1) G4cout << "Construct processes " << G4endl;
314#endif
316}
virtual void ConstructProcess()=0
G4bool IsMasterThread()

Referenced by G4RunManagerKernel::InitializePhysics().

◆ ConstructParticle()

virtual void G4VUserPhysicsList::ConstructParticle ( )
pure virtual

◆ ConstructProcess()

virtual void G4VUserPhysicsList::ConstructProcess ( )
pure virtual

Implemented in G4ErrorPhysicsList, G4VModularPhysicsList, and LBE.

Referenced by Construct().

◆ DisableCheckParticleList()

void G4VUserPhysicsList::DisableCheckParticleList ( )
inline

Definition at line 360 of file G4VUserPhysicsList.hh.

361{
363}

◆ DumpCutValuesTable()

void G4VUserPhysicsList::DumpCutValuesTable ( G4int flag = 1)

◆ DumpCutValuesTableIfRequested()

void G4VUserPhysicsList::DumpCutValuesTableIfRequested ( )

◆ DumpList()

void G4VUserPhysicsList::DumpList ( ) const

Definition at line 770 of file G4VUserPhysicsList.cc.

771{
772 theParticleIterator->reset();
773 G4int idx = 0;
774 while ((*theParticleIterator)()) {
775 G4ParticleDefinition* particle = theParticleIterator->value();
776 G4cout << particle->GetParticleName();
777 if ((idx++ % 4) == 3) {
778 G4cout << G4endl;
779 }
780 else {
781 G4cout << ", ";
782 }
783 }
784 G4cout << G4endl;
785}

Referenced by G4UserPhysicsListMessenger::SetNewValue().

◆ GetApplyCuts()

G4bool G4VUserPhysicsList::GetApplyCuts ( const G4String & name) const

Definition at line 906 of file G4VUserPhysicsList.cc.

907{
909}
G4bool GetApplyCutsFlag() const

◆ GetCutValue()

G4double G4VUserPhysicsList::GetCutValue ( const G4String & pname) const

Definition at line 379 of file G4VUserPhysicsList.cc.

380{
381 std::size_t nReg = (G4RegionStore::GetInstance())->size();
382 if (nReg == 0) {
383#ifdef G4VERBOSE
384 if (verboseLevel > 0) {
385 G4cout << "G4VUserPhysicsList::GetCutValue "
386 << " : No Default Region " << G4endl;
387 }
388#endif
389 G4Exception("G4VUserPhysicsList::GetCutValue", "Run0253", FatalException, "No Default Region");
390 return -1. * mm;
391 }
392 G4Region* region = G4RegionStore::GetInstance()->GetRegion("DefaultRegionForTheWorld", false);
393 return region->GetProductionCuts()->GetProductionCut(name);
394}
G4double GetProductionCut(G4int index) const
static G4RegionStore * GetInstance()
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
G4ProductionCuts * GetProductionCuts() const

Referenced by SetCuts(), and G4UserPhysicsListMessenger::SetNewValue().

◆ GetDefaultCutValue()

G4double G4VUserPhysicsList::GetDefaultCutValue ( ) const
inline

Definition at line 318 of file G4VUserPhysicsList.hh.

319{
320 return defaultCutValue;
321}

Referenced by G4UserPhysicsListMessenger::GetCurrentValue().

◆ GetInstanceID()

G4int G4VUserPhysicsList::GetInstanceID ( ) const
inline

Definition at line 365 of file G4VUserPhysicsList.hh.

366{
367 return g4vuplInstanceID;
368}

Referenced by G4VUserPhysicsList(), and operator=().

◆ GetParticleIterator()

G4ParticleTable::G4PTblDicIterator * G4VUserPhysicsList::GetParticleIterator ( ) const
protected

Definition at line 938 of file G4VUserPhysicsList.cc.

939{
940 return (subInstanceManager.offset[g4vuplInstanceID])._theParticleIterator;
941}

Referenced by G4ErrorPhysicsList::ConstructEM(), and G4ErrorPhysicsList::ConstructProcess().

◆ GetPhysicsTableDirectory()

const G4String & G4VUserPhysicsList::GetPhysicsTableDirectory ( ) const
inline

Definition at line 338 of file G4VUserPhysicsList.hh.

339{
341}

Referenced by G4UserPhysicsListMessenger::GetCurrentValue().

◆ GetSubInstanceManager()

const G4VUPLManager & G4VUserPhysicsList::GetSubInstanceManager ( )
inlinestatic

◆ GetVerboseLevel()

G4int G4VUserPhysicsList::GetVerboseLevel ( ) const
inline

Definition at line 323 of file G4VUserPhysicsList.hh.

324{
325 return verboseLevel;
326}

Referenced by G4UserPhysicsListMessenger::GetCurrentValue().

◆ InitializeProcessManager()

void G4VUserPhysicsList::InitializeProcessManager ( )
protected

Definition at line 200 of file G4VUserPhysicsList.cc.

201{
202 // Request lock for particle table access. Some changes are inside
203 // this critical region.
204#ifdef G4MULTITHREADED
205 G4MUTEXLOCK(&G4ParticleTable::particleTableMutex());
206 G4ParticleTable::lockCount()++;
207#endif
209
210 // loop over all particles in G4ParticleTable
211 theParticleIterator->reset();
212 while ((*theParticleIterator)()) {
213 G4ParticleDefinition* particle = theParticleIterator->value();
214 G4ProcessManager* pmanager = particle->GetProcessManager();
215
216 if (pmanager == nullptr) {
217 // create process manager if the particle does not have its own.
218 pmanager = new G4ProcessManager(particle);
219 particle->SetProcessManager(pmanager);
220 if (particle->GetMasterProcessManager() == nullptr)
221 particle->SetMasterProcessManager(pmanager);
222#ifdef G4VERBOSE
223 if (verboseLevel > 2) {
224 G4cout << "G4VUserPhysicsList::InitializeProcessManager: creating "
225 "ProcessManager to "
226 << particle->GetParticleName() << G4endl;
227 }
228#endif
229 }
230 }
231
232 if (gion != nullptr) {
233 G4ProcessManager* gionPM = gion->GetProcessManager();
234 // loop over all particles once again (this time, with all general ions)
235 theParticleIterator->reset(false);
236 while ((*theParticleIterator)()) {
237 G4ParticleDefinition* particle = theParticleIterator->value();
238 if (particle->IsGeneralIon()) {
239 particle->SetProcessManager(gionPM);
240#ifdef G4VERBOSE
241 if (verboseLevel > 2) {
242 G4cout << "G4VUserPhysicsList::InitializeProcessManager: copying "
243 "ProcessManager to "
244 << particle->GetParticleName() << G4endl;
245 }
246#endif
247 }
248 }
249 }
250
251 // release lock for particle table access
252#ifdef G4MULTITHREADED
253 G4MUTEXUNLOCK(&G4ParticleTable::particleTableMutex());
254#endif
255}
#define G4MUTEXLOCK(mutex)
#define G4MUTEXUNLOCK(mutex)
void SetMasterProcessManager(G4ProcessManager *aNewPM)
G4bool IsGeneralIon() const
void SetProcessManager(G4ProcessManager *aProcessManager)
G4ParticleDefinition * GetGenericIon() const

Referenced by Construct().

◆ InitializeWorker()

void G4VUserPhysicsList::InitializeWorker ( )
virtual

Definition at line 108 of file G4VUserPhysicsList.cc.

109{
110 // Remember messengers are per-thread, so this needs to be done by each
111 // worker and due to the presence of "this" cannot be done in
112 // G4VUPLData::initialize()
114}

Referenced by G4WorkerRunManager::SetUserInitialization().

◆ IsPhysicsTableRetrieved()

G4bool G4VUserPhysicsList::IsPhysicsTableRetrieved ( ) const
inline

Definition at line 328 of file G4VUserPhysicsList.hh.

329{
331}

Referenced by G4UserPhysicsListMessenger::GetCurrentValue().

◆ IsStoredInAscii()

G4bool G4VUserPhysicsList::IsStoredInAscii ( ) const
inline

Definition at line 333 of file G4VUserPhysicsList.hh.

334{
335 return fStoredInAscii;
336}

Referenced by G4UserPhysicsListMessenger::GetCurrentValue().

◆ operator=()

G4VUserPhysicsList & G4VUserPhysicsList::operator= ( const G4VUserPhysicsList & right)

◆ PreparePhysicsTable()

void G4VUserPhysicsList::PreparePhysicsTable ( G4ParticleDefinition * particle)

Definition at line 687 of file G4VUserPhysicsList.cc.

688{
689 if (auto* trackingManager = particle->GetTrackingManager()) {
690 trackingManager->PreparePhysicsTable(*particle);
691 return;
692 }
693
694 if ((particle->GetMasterProcessManager()) == nullptr) {
695 return;
696 }
697 // Prepare the physics tables for every process for this particle type
698 // if particle is not ShortLived
699 if (!particle->IsShortLived()) {
700 G4ProcessManager* pManager = particle->GetProcessManager();
701 if (pManager == nullptr) {
702#ifdef G4VERBOSE
703 if (verboseLevel > 0) {
704 G4cout << "G4VUserPhysicsList::PreparePhysicsTable "
705 << ": No Process Manager for " << particle->GetParticleName() << G4endl;
706 G4cout << particle->GetParticleName() << " should be created in your PhysicsList" << G4endl;
707 }
708#endif
709 G4Exception("G4VUserPhysicsList::PreparePhysicsTable", "Run0273", FatalException,
710 "No process manager");
711 return;
712 }
713
714 // Get processes from master thread
715 G4ProcessManager* pManagerShadow = particle->GetMasterProcessManager();
716
717 G4ProcessVector* pVector = pManager->GetProcessList();
718 if (pVector == nullptr) {
719#ifdef G4VERBOSE
720 if (verboseLevel > 0) {
721 G4cout << "G4VUserPhysicsList::PreparePhysicsTable "
722 << ": No Process Vector for " << particle->GetParticleName() << G4endl;
723 }
724#endif
725 G4Exception("G4VUserPhysicsList::PreparePhysicsTable", "Run0274", FatalException,
726 "No process Vector");
727 return;
728 }
729 for (G4int j = 0; j < (G4int)pVector->size(); ++j) {
730 // Andrea July 16th 2013 : migration to new interface...
731 // Infer if we are in a worker thread or master thread
732 // Master thread is the one in which the process manager
733 // and process manager shadow pointers are the same
734 if (pManagerShadow == pManager) {
735 (*pVector)[j]->PreparePhysicsTable(*particle);
736 }
737 else {
738 (*pVector)[j]->PrepareWorkerPhysicsTable(*particle);
739 }
740 } // End loop on processes vector
741 } // End if pn ShortLived
742}

Referenced by BuildPhysicsTable(), and G4UserPhysicsListMessenger::SetNewValue().

◆ RegisterProcess()

G4bool G4VUserPhysicsList::RegisterProcess ( G4VProcess * process,
G4ParticleDefinition * particle )
protected

Definition at line 932 of file G4VUserPhysicsList.cc.

933{
934 return G4MT_thePLHelper->RegisterProcess(process, particle);
935}

◆ RemoveProcessManager()

void G4VUserPhysicsList::RemoveProcessManager ( )

Definition at line 258 of file G4VUserPhysicsList.cc.

259{
260 // Request lock for particle table access. Some changes are inside
261 // this critical region.
262#ifdef G4MULTITHREADED
263 G4MUTEXLOCK(&G4ParticleTable::particleTableMutex());
264 G4ParticleTable::lockCount()++;
265#endif
266
267 // loop over all particles in G4ParticleTable
268 theParticleIterator->reset();
269 while ((*theParticleIterator)()) {
270 G4ParticleDefinition* particle = theParticleIterator->value();
272 if (particle->GetParticleSubType() != "generic"
273 || particle->GetParticleName() == "GenericIon")
274 {
275 G4ProcessManager* pmanager = particle->GetProcessManager();
276
277 delete pmanager;
278#ifdef G4VERBOSE
279 if (verboseLevel > 2) {
280 G4cout << "G4VUserPhysicsList::RemoveProcessManager: ";
281 G4cout << "remove ProcessManager from ";
282 G4cout << particle->GetParticleName() << G4endl;
283 }
284#endif
285 }
286 particle->SetProcessManager(nullptr);
287 }
288 }
289
290 // release lock for particle table access
291#ifdef G4MULTITHREADED
292 G4MUTEXUNLOCK(&G4ParticleTable::particleTableMutex());
293#endif
294}
static G4PART_DLL G4int & slavetotalspace()
G4int GetInstanceID() const
const G4String & GetParticleSubType() const

Referenced by TerminateWorker(), and ~G4VUserPhysicsList().

◆ RemoveTrackingManager()

void G4VUserPhysicsList::RemoveTrackingManager ( )

Definition at line 297 of file G4VUserPhysicsList.cc.

298{
299 // One tracking manager may be registered for multiple particles, make sure
300 // to delete every object only once.
301 std::unordered_set<G4VTrackingManager*> trackingManagers;
302
303 // loop over all particles in G4ParticleTable
304 theParticleIterator->reset();
305 while ((*theParticleIterator)()) {
306 G4ParticleDefinition* particle = theParticleIterator->value();
307 if (auto* trackingManager = particle->GetTrackingManager()) {
308#ifdef G4VERBOSE
309 if (verboseLevel > 2) {
310 G4cout << "G4VUserPhysicsList::RemoveTrackingManager: ";
311 G4cout << "remove TrackingManager from ";
312 G4cout << particle->GetParticleName() << G4endl;
313 }
314#endif
315 trackingManagers.insert(trackingManager);
316 particle->SetTrackingManager(nullptr);
317 }
318 }
319
320 for (G4VTrackingManager* tm : trackingManagers) {
321 delete tm;
322 }
323}
void SetTrackingManager(G4VTrackingManager *aTrackingManager)

Referenced by TerminateWorker(), and ~G4VUserPhysicsList().

◆ ResetPhysicsTableRetrieved()

void G4VUserPhysicsList::ResetPhysicsTableRetrieved ( )
inline

Definition at line 348 of file G4VUserPhysicsList.hh.

349{
350 fRetrievePhysicsTable = false;
351 fIsRestoredCutValues = false;
353}

Referenced by G4UserPhysicsListMessenger::SetNewValue().

◆ ResetStoredInAscii()

void G4VUserPhysicsList::ResetStoredInAscii ( )
inline

Definition at line 355 of file G4VUserPhysicsList.hh.

356{
357 fStoredInAscii = false;
358}

Referenced by G4UserPhysicsListMessenger::SetNewValue().

◆ RetrievePhysicsTable()

void G4VUserPhysicsList::RetrievePhysicsTable ( G4ParticleDefinition * particle,
const G4String & directory,
G4bool ascii = false )
protectedvirtual

Definition at line 859 of file G4VUserPhysicsList.cc.

861{
862 G4bool success[100];
863 // Retrieve physics tables for every process for this particle type
864 G4ProcessVector* pVector = (particle->GetProcessManager())->GetProcessList();
865 for (G4int j = 0; j < (G4int)pVector->size(); ++j) {
866 success[j] = (*pVector)[j]->RetrievePhysicsTable(particle, directory, ascii);
867
868 if (!success[j]) {
869#ifdef G4VERBOSE
870 if (verboseLevel > 2) {
871 G4cout << "G4VUserPhysicsList::RetrievePhysicsTable "
872 << " Fail to retrieve Physics Table for " << (*pVector)[j]->GetProcessName()
873 << G4endl;
874 G4cout << "Calculate Physics Table for " << particle->GetParticleName() << G4endl;
875 }
876#endif
877 (*pVector)[j]->BuildPhysicsTable(*particle);
878 }
879 }
880 for (G4int j = 0; j < (G4int)pVector->size(); ++j) {
881 // temporary addition to make the integral schema
882 if (!success[j]) BuildIntegralPhysicsTable((*pVector)[j], particle);
883 }
884}
bool G4bool
Definition G4Types.hh:86
void BuildIntegralPhysicsTable(G4VProcess *, G4ParticleDefinition *)

Referenced by BuildPhysicsTable().

◆ SetApplyCuts()

void G4VUserPhysicsList::SetApplyCuts ( G4bool value,
const G4String & name )

Definition at line 887 of file G4VUserPhysicsList.cc.

888{
889#ifdef G4VERBOSE
890 if (verboseLevel > 2) {
891 G4cout << "G4VUserPhysicsList::SetApplyCuts for " << name << G4endl;
892 }
893#endif
894 if (name == "all") {
899 }
900 else {
902 }
903}
const char * name(G4int ptype)

Referenced by G4UserPhysicsListMessenger::SetNewValue().

◆ SetCuts()

void G4VUserPhysicsList::SetCuts ( )
virtual

Reimplemented in FTFP_BERT_HP, G4ErrorPhysicsList, LBE, and QGSP_BERT_HP.

Definition at line 326 of file G4VUserPhysicsList.cc.

327{
330 }
331
332#ifdef G4VERBOSE
333 if (verboseLevel > 1) {
334 G4cout << "G4VUserPhysicsList::SetCuts: " << G4endl;
335 G4cout << "Cut for gamma: " << GetCutValue("gamma") / mm << "[mm]" << G4endl;
336 G4cout << "Cut for e-: " << GetCutValue("e-") / mm << "[mm]" << G4endl;
337 G4cout << "Cut for e+: " << GetCutValue("e+") / mm << "[mm]" << G4endl;
338 G4cout << "Cut for proton: " << GetCutValue("proton") / mm << "[mm]" << G4endl;
339 }
340#endif
341
342 // dump Cut values if verboseLevel==3
343 if (verboseLevel > 2) {
345 }
346}
G4double GetCutValue(const G4String &pname) const
void SetDefaultCutValue(G4double newCutValue)
void DumpCutValuesTable(G4int flag=1)

Referenced by G4RunManagerKernel::InitializePhysics(), SetCutsWithDefault(), and G4UserPhysicsListMessenger::SetNewValue().

◆ SetCutsForRegion()

void G4VUserPhysicsList::SetCutsForRegion ( G4double aCut,
const G4String & rname )

Definition at line 428 of file G4VUserPhysicsList.cc.

429{
430 // set cut values for gamma at first and for e- and e+
431 SetCutValue(aCut, "gamma", rname);
432 SetCutValue(aCut, "e-", rname);
433 SetCutValue(aCut, "e+", rname);
434 SetCutValue(aCut, "proton", rname);
435}
void SetCutValue(G4double aCut, const G4String &pname)

Referenced by G4UserPhysicsListMessenger::SetNewValue().

◆ SetCutsWithDefault()

void G4VUserPhysicsList::SetCutsWithDefault ( )

◆ SetCutValue() [1/2]

void G4VUserPhysicsList::SetCutValue ( G4double aCut,
const G4String & pname )

◆ SetCutValue() [2/2]

void G4VUserPhysicsList::SetCutValue ( G4double aCut,
const G4String & pname,
const G4String & rname )

Definition at line 403 of file G4VUserPhysicsList.cc.

404{
406 if (region != nullptr) {
407 // set cut value
408 SetParticleCuts(aCut, pname, region);
409 }
410 else {
411#ifdef G4VERBOSE
412 if (verboseLevel > 0) {
413 G4cout << "G4VUserPhysicsList::SetCutValue "
414 << " : No Region of " << rname << G4endl;
415 }
416#endif
417 }
418}

◆ SetDefaultCutValue()

void G4VUserPhysicsList::SetDefaultCutValue ( G4double newCutValue)

Definition at line 349 of file G4VUserPhysicsList.cc.

350{
351 if (value < 0.0) {
352#ifdef G4VERBOSE
353 if (verboseLevel > 0) {
354 G4cout << "G4VUserPhysicsList::SetDefaultCutValue: negative cut values"
355 << " :" << value / mm << "[mm]" << G4endl;
356 }
357#endif
358 return;
359 }
360
361 defaultCutValue = value;
363
364 // set cut values for gamma at first and for e- and e+
368 SetCutValue(defaultCutValue, "proton");
369
370#ifdef G4VERBOSE
371 if (verboseLevel > 1) {
372 G4cout << "G4VUserPhysicsList::SetDefaultCutValue:"
373 << "default cut value is changed to :" << defaultCutValue / mm << "[mm]" << G4endl;
374 }
375#endif
376}

Referenced by SetCuts(), SetCutsWithDefault(), G4UserPhysicsListMessenger::SetNewValue(), and SetParticleCuts().

◆ SetParticleCuts() [1/2]

void G4VUserPhysicsList::SetParticleCuts ( G4double cut,
const G4String & particleName,
G4Region * region = nullptr )

Definition at line 445 of file G4VUserPhysicsList.cc.

447{
448 if (cut < 0.0) {
449#ifdef G4VERBOSE
450 if (verboseLevel > 0) {
451 G4cout << "G4VUserPhysicsList::SetParticleCuts: negative cut values"
452 << " :" << cut / mm << "[mm]"
453 << " for " << particleName << G4endl;
454 }
455#endif
456 return;
457 }
458
459 G4Region* world_region =
460 G4RegionStore::GetInstance()->GetRegion("DefaultRegionForTheWorld", false);
461 if (region == nullptr) {
462 std::size_t nReg = G4RegionStore::GetInstance()->size();
463 if (nReg == 0) {
464#ifdef G4VERBOSE
465 if (verboseLevel > 0) {
466 G4cout << "G4VUserPhysicsList::SetParticleCuts "
467 << " : No Default Region " << G4endl;
468 }
469#endif
470 G4Exception("G4VUserPhysicsList::SetParticleCuts ", "Run0254", FatalException,
471 "No Default Region");
472 return;
473 }
474 region = world_region;
475 }
476
479 }
480
481 G4ProductionCuts* pcuts = region->GetProductionCuts();
482 if (region != world_region
483 && pcuts == G4ProductionCutsTable::GetProductionCutsTable()->GetDefaultProductionCuts())
484 { // This region had no unique cuts yet but shares the default cuts.
485 // Need to create a new object before setting the value.
486 pcuts = new G4ProductionCuts(
487 *(G4ProductionCutsTable::GetProductionCutsTable()->GetDefaultProductionCuts()));
488 region->SetProductionCuts(pcuts);
489 }
490 pcuts->SetProductionCut(cut, particleName);
491#ifdef G4VERBOSE
492 if (verboseLevel > 2) {
493 G4cout << "G4VUserPhysicsList::SetParticleCuts: "
494 << " :" << cut / mm << "[mm]"
495 << " for " << particleName << G4endl;
496 }
497#endif
498}
void SetProductionCut(G4double cut, G4int index)
void SetProductionCuts(G4ProductionCuts *cut)

◆ SetParticleCuts() [2/2]

void G4VUserPhysicsList::SetParticleCuts ( G4double cut,
G4ParticleDefinition * particle,
G4Region * region = nullptr )

Definition at line 438 of file G4VUserPhysicsList.cc.

440{
441 SetParticleCuts(cut, particle->GetParticleName(), region);
442}

Referenced by SetCutValue(), SetCutValue(), and SetParticleCuts().

◆ SetPhysicsTableRetrieved()

void G4VUserPhysicsList::SetPhysicsTableRetrieved ( const G4String & directory = "")

Definition at line 848 of file G4VUserPhysicsList.cc.

849{
851 if (!directory.empty()) {
852 directoryPhysicsTable = directory;
853 }
855 fIsRestoredCutValues = false;
856}

Referenced by G4UserPhysicsListMessenger::SetNewValue().

◆ SetStoredInAscii()

void G4VUserPhysicsList::SetStoredInAscii ( )
inline

Definition at line 343 of file G4VUserPhysicsList.hh.

344{
345 fStoredInAscii = true;
346}

Referenced by G4UserPhysicsListMessenger::SetNewValue().

◆ SetVerboseLevel()

void G4VUserPhysicsList::SetVerboseLevel ( G4int value)

Definition at line 944 of file G4VUserPhysicsList.cc.

945{
946 verboseLevel = value;
947 // set verboseLevel for G4ProductionCutsTable same as one for
948 // G4VUserPhysicsList:
950
951 G4MT_thePLHelper->SetVerboseLevel(verboseLevel);
952
953#ifdef G4VERBOSE
954 if (verboseLevel > 1) {
955 G4cout << "G4VUserPhysicsList::SetVerboseLevel :"
956 << " Verbose level is set to " << verboseLevel << G4endl;
957 }
958#endif
959}

Referenced by G4UserPhysicsListMessenger::SetNewValue().

◆ StorePhysicsTable()

G4bool G4VUserPhysicsList::StorePhysicsTable ( const G4String & directory = ".")

Definition at line 802 of file G4VUserPhysicsList.cc.

803{
804 G4bool ascii = fStoredInAscii;
805 G4String dir = directory;
806 if (dir.empty())
808 else
810
811 // store CutsTable info
812 if (!fCutsTable->StoreCutsTable(dir, ascii)) {
813 G4Exception("G4VUserPhysicsList::StorePhysicsTable", "Run0281", JustWarning,
814 "Fail to store Cut Table");
815 return false;
816 }
817#ifdef G4VERBOSE
818 if (verboseLevel > 2) {
819 G4cout << "G4VUserPhysicsList::StorePhysicsTable "
820 << " Store material and cut values successfully" << G4endl;
821 }
822#endif
823
824 G4bool success = true;
825
826 // loop over all particles in G4ParticleTable
827 theParticleIterator->reset();
828 while ((*theParticleIterator)()) {
829 G4ParticleDefinition* particle = theParticleIterator->value();
830 // Store physics tables for every process for this particle type
831 G4ProcessVector* pVector = (particle->GetProcessManager())->GetProcessList();
832 for (G4int j = 0; j < (G4int)pVector->size(); ++j) {
833 if (!(*pVector)[j]->StorePhysicsTable(particle, dir, ascii)) {
834 G4String comment = "Fail to store physics table for ";
835 comment += (*pVector)[j]->GetProcessName();
836 comment += "(" + particle->GetParticleName() + ")";
837 G4Exception("G4VUserPhysicsList::StorePhysicsTable", "Run0282", JustWarning, comment);
838 success = false;
839 }
840 }
841 // end loop over processes
842 }
843 // end loop over particles
844 return success;
845}
G4bool StoreCutsTable(const G4String &directory, G4bool ascii=false)
G4bool StorePhysicsTable(const G4String &directory=".")

Referenced by G4UserPhysicsListMessenger::SetNewValue(), and StorePhysicsTable().

◆ TerminateWorker()

void G4VUserPhysicsList::TerminateWorker ( )
virtual

Reimplemented in G4VModularPhysicsList.

Definition at line 117 of file G4VUserPhysicsList.cc.

118{
121 delete G4MT_theMessenger;
122 G4MT_theMessenger = nullptr;
123}

Referenced by G4VModularPhysicsList::TerminateWorker().

◆ UseCoupledTransportation()

void G4VUserPhysicsList::UseCoupledTransportation ( G4bool vl = true)

Definition at line 926 of file G4VUserPhysicsList.cc.

927{
928 G4MT_thePLHelper->UseCoupledTransportation(vl);
929}

Referenced by G4RunManagerKernel::InitializePhysics().

Member Data Documentation

◆ defaultCutValue

◆ directoryPhysicsTable

◆ fCutsTable

G4ProductionCutsTable* G4VUserPhysicsList::fCutsTable = nullptr
protected

◆ fDisableCheckParticleList

G4bool G4VUserPhysicsList::fDisableCheckParticleList = false
protected

◆ fIsCheckedForRetrievePhysicsTable

G4bool G4VUserPhysicsList::fIsCheckedForRetrievePhysicsTable = false
protected

◆ fIsRestoredCutValues

G4bool G4VUserPhysicsList::fIsRestoredCutValues = false
protected

◆ fRetrievePhysicsTable

◆ fStoredInAscii

◆ g4vuplInstanceID

G4int G4VUserPhysicsList::g4vuplInstanceID = 0
protected

◆ isSetDefaultCutValue

G4bool G4VUserPhysicsList::isSetDefaultCutValue = false
protected

◆ subInstanceManager

G4VUPLManager G4VUserPhysicsList::subInstanceManager
staticprotected

◆ theParticleTable

G4ParticleTable* G4VUserPhysicsList::theParticleTable = nullptr
protected

◆ verboseLevel


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