Geant4-11
Public Member Functions | Static Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | Static Protected Attributes
G4EmPenelopePhysics Class Reference

#include <G4EmPenelopePhysics.hh>

Inheritance diagram for G4EmPenelopePhysics:
G4VPhysicsConstructor

Public Member Functions

void ConstructParticle () override
 
void ConstructProcess () override
 
 G4EmPenelopePhysics (G4int ver=1, const G4String &name="")
 
G4int GetInstanceID () const
 
const G4StringGetPhysicsName () const
 
G4int GetPhysicsType () const
 
G4int GetVerboseLevel () const
 
void SetPhysicsName (const G4String &="")
 
void SetPhysicsType (G4int)
 
void SetVerboseLevel (G4int value)
 
virtual void TerminateWorker ()
 
 ~G4EmPenelopePhysics () override
 

Static Public Member Functions

static const G4VPCManagerGetSubInstanceManager ()
 

Protected Types

using PhysicsBuilder_V = G4VPCData::PhysicsBuilders_V
 

Protected Member Functions

void AddBuilder (G4PhysicsBuilderInterface *bld)
 
PhysicsBuilder_V GetBuilders () const
 
G4ParticleTable::G4PTblDicIteratorGetParticleIterator () const
 
G4bool RegisterProcess (G4VProcess *process, G4ParticleDefinition *particle)
 

Protected Attributes

G4int g4vpcInstanceID = 0
 
G4String namePhysics = ""
 
G4ParticleTabletheParticleTable = nullptr
 
G4int typePhysics = 0
 
G4int verboseLevel = 0
 

Static Protected Attributes

static G4RUN_DLL G4VPCManager subInstanceManager
 

Detailed Description

Definition at line 35 of file G4EmPenelopePhysics.hh.

Member Typedef Documentation

◆ PhysicsBuilder_V

Definition at line 149 of file G4VPhysicsConstructor.hh.

Constructor & Destructor Documentation

◆ G4EmPenelopePhysics()

G4EmPenelopePhysics::G4EmPenelopePhysics ( G4int  ver = 1,
const G4String name = "" 
)
explicit

Definition at line 107 of file G4EmPenelopePhysics.cc.

108 : G4VPhysicsConstructor("G4EmPenelope")
109{
110 SetVerboseLevel(ver);
112 param->SetDefaults();
113 param->SetVerbose(ver);
114 param->SetMinEnergy(100*CLHEP::eV);
116 param->SetNumberOfBinsPerDecade(20);
118 param->SetStepFunction(0.2, 10*CLHEP::um);
119 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
120 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
121 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
122 param->SetUseMottCorrection(true);
124 param->SetMscSkin(3);
125 param->SetMscRangeFactor(0.08);
126 param->SetMuHadLateralDisplacement(true);
127 param->SetFluo(true);
129 param->SetPIXEElectronCrossSectionModel("Penelope");
131}
@ bElectromagnetic
@ fUseSafetyPlus
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
void SetPIXEElectronCrossSectionModel(const G4String &)
void SetMuHadLateralDisplacement(G4bool val)
void ActivateAngularGeneratorForIonisation(G4bool val)
void SetStepFunction(G4double v1, G4double v2)
void SetFluo(G4bool val)
void SetStepFunctionMuHad(G4double v1, G4double v2)
void SetVerbose(G4int val)
void SetMscSkin(G4double val)
void SetMaxNIELEnergy(G4double val)
void SetStepFunctionIons(G4double v1, G4double v2)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetUseMottCorrection(G4bool val)
void SetMscRangeFactor(G4double val)
G4VPhysicsConstructor(const G4String &="")
void SetVerboseLevel(G4int value)
static constexpr double um
Definition: SystemOfUnits.h:94
static constexpr double MeV
static constexpr double eV

References G4EmParameters::ActivateAngularGeneratorForIonisation(), bElectromagnetic, CLHEP::eV, fUseSafetyPlus, G4EmParameters::Instance(), CLHEP::MeV, G4EmParameters::SetDefaults(), G4EmParameters::SetFluo(), G4EmParameters::SetLowestElectronEnergy(), G4EmParameters::SetMaxNIELEnergy(), G4EmParameters::SetMinEnergy(), G4EmParameters::SetMscRangeFactor(), G4EmParameters::SetMscSkin(), G4EmParameters::SetMscStepLimitType(), G4EmParameters::SetMuHadLateralDisplacement(), G4EmParameters::SetNumberOfBinsPerDecade(), G4VPhysicsConstructor::SetPhysicsType(), G4EmParameters::SetPIXEElectronCrossSectionModel(), G4EmParameters::SetStepFunction(), G4EmParameters::SetStepFunctionIons(), G4EmParameters::SetStepFunctionLightIons(), G4EmParameters::SetStepFunctionMuHad(), G4EmParameters::SetUseMottCorrection(), G4EmParameters::SetVerbose(), G4VPhysicsConstructor::SetVerboseLevel(), and CLHEP::um.

◆ ~G4EmPenelopePhysics()

G4EmPenelopePhysics::~G4EmPenelopePhysics ( )
override

Definition at line 135 of file G4EmPenelopePhysics.cc.

136{}

Member Function Documentation

◆ AddBuilder()

void G4VPhysicsConstructor::AddBuilder ( G4PhysicsBuilderInterface bld)
protectedinherited

Definition at line 99 of file G4VPhysicsConstructor.cc.

100{
101 (subInstanceManager.offset[g4vpcInstanceID])._builders->push_back(bld);
102}
static G4RUN_DLL G4VPCManager subInstanceManager
G4RUN_DLL G4ThreadLocalStatic T * offset

References G4VPhysicsConstructor::g4vpcInstanceID, G4VUPLSplitter< T >::offset, and G4VPhysicsConstructor::subInstanceManager.

Referenced by G4HadronPhysicsFTFP_BERT::Kaon(), G4HadronPhysicsFTF_BIC::Kaon(), G4HadronPhysicsINCLXX::Kaon(), G4HadronPhysicsFTFP_BERT::Neutron(), G4HadronPhysicsQGSP_BERT::Neutron(), G4HadronPhysicsQGSP_BIC::Neutron(), G4HadronPhysicsFTF_BIC::Neutron(), G4HadronPhysicsFTFP_BERT_HP::Neutron(), G4HadronPhysicsINCLXX::Neutron(), G4HadronPhysicsQGS_BIC::Neutron(), G4HadronPhysicsQGSP_BERT_HP::Neutron(), G4HadronPhysicsQGSP_BIC_HP::Neutron(), G4HadronPhysicsShielding::Neutron(), G4HadronPhysicsFTFP_BERT::Pion(), G4HadronPhysicsQGSP_BERT::Pion(), G4HadronPhysicsQGSP_BIC::Pion(), G4HadronPhysicsFTF_BIC::Pion(), G4HadronPhysicsINCLXX::Pion(), G4HadronPhysicsQGS_BIC::Pion(), G4HadronPhysicsFTFP_BERT::Proton(), G4HadronPhysicsQGSP_BERT::Proton(), G4HadronPhysicsQGSP_BIC::Proton(), G4HadronPhysicsFTF_BIC::Proton(), G4HadronPhysicsINCLXX::Proton(), G4HadronPhysicsNuBeam::Proton(), G4HadronPhysicsQGS_BIC::Proton(), and G4HadronPhysicsQGSP_BIC_AllHP::Proton().

◆ ConstructParticle()

void G4EmPenelopePhysics::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 140 of file G4EmPenelopePhysics.cc.

141{
142 // minimal set of particles for EM physics
144}
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:347

References G4EmBuilder::ConstructMinimalEmSet().

◆ ConstructProcess()

void G4EmPenelopePhysics::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 148 of file G4EmPenelopePhysics.cc.

149{
150 if(verboseLevel > 1) {
151 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
152 }
156
157 // processes used by several particles
158 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
159
160 // high energy limit for e+- scattering models
161 G4double highEnergyLimit = param->MscEnergyLimit();
162
163 // nuclear stopping
164 G4double nielEnergyLimit = param->MaxNIELEnergy();
165 G4NuclearStopping* pnuc = nullptr;
166 if(nielEnergyLimit > 0.0) {
167 pnuc = new G4NuclearStopping();
168 pnuc->SetMaxKinEnergy(nielEnergyLimit);
169 }
170
171 //Applicability range for Penelope models
172 //for higher energies, the Standard models are used
173 G4double PenelopeHighEnergyLimit = 1.0*CLHEP::GeV;
174
175 // Add gamma EM Processes
177
178 //Photo-electric effect
179 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
180 thePhotoElectricEffect->SetEmModel(new G4PEEffectFluoModel());
182 thePEPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
183 thePhotoElectricEffect->AddEmModel(0, thePEPenelopeModel);
184
185 //Compton scattering
186 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
187 G4PenelopeComptonModel* theComptonPenelopeModel = new G4PenelopeComptonModel();
188 theComptonScattering->SetEmModel(new G4KleinNishinaModel());
189 theComptonPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
190 theComptonScattering->AddEmModel(0, theComptonPenelopeModel);
191
192 //Gamma conversion
193 G4GammaConversion* theGammaConversion = new G4GammaConversion();
195 theGCPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
196 theGammaConversion->AddEmModel(0, theGCPenelopeModel);
197
198 //Rayleigh scattering
199 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
200 G4PenelopeRayleighModel* theRayleighPenelopeModel = new G4PenelopeRayleighModel();
201 theRayleigh->SetEmModel(theRayleighPenelopeModel);
202
203 ph->RegisterProcess(thePhotoElectricEffect, particle);
204 ph->RegisterProcess(theComptonScattering, particle);
205 ph->RegisterProcess(theGammaConversion, particle);
206 ph->RegisterProcess(theRayleigh, particle);
207
208 // e-
209 particle = G4Electron::Electron();
210
211 // multiple scattering
215 msc1->SetHighEnergyLimit(highEnergyLimit);
216 msc2->SetLowEnergyLimit(highEnergyLimit);
217 msc->SetEmModel(msc1);
218 msc->SetEmModel(msc2);
219
222 ss->SetEmModel(ssm);
223 ss->SetMinKinEnergy(highEnergyLimit);
224 ssm->SetLowEnergyLimit(highEnergyLimit);
225 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
226
227 //Ionisation
228 G4eIonisation* eioni = new G4eIonisation();
230 theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
231 eioni->AddEmModel(0, theIoniPenelope, new G4UniversalFluctuation());
232
233 //Bremsstrahlung
236 theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
237 brem->SetEmModel(theBremPenelope);
238
240
241 // register processes
242 ph->RegisterProcess(msc, particle);
243 ph->RegisterProcess(eioni, particle);
244 ph->RegisterProcess(brem, particle);
245 ph->RegisterProcess(ee, particle);
246 ph->RegisterProcess(ss, particle);
247
248 // e+
249 particle = G4Positron::Positron();
250
251 // multiple scattering
252 msc = new G4eMultipleScattering;
253 msc1 = new G4GoudsmitSaundersonMscModel();
254 msc2 = new G4WentzelVIModel();
255 msc1->SetHighEnergyLimit(highEnergyLimit);
256 msc2->SetLowEnergyLimit(highEnergyLimit);
257 msc->SetEmModel(msc1);
258 msc->SetEmModel(msc2);
259
260 ssm = new G4eCoulombScatteringModel();
261 ss = new G4CoulombScattering();
262 ss->SetEmModel(ssm);
263 ss->SetMinKinEnergy(highEnergyLimit);
264 ssm->SetLowEnergyLimit(highEnergyLimit);
265 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
266
267 //Ionisation
268 eioni = new G4eIonisation();
269 theIoniPenelope = new G4PenelopeIonisationModel();
270 theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
271 eioni->AddEmModel(0,theIoniPenelope, new G4UniversalFluctuation());
272
273 //Bremsstrahlung
274 brem = new G4eBremsstrahlung();
275 theBremPenelope = new G4PenelopeBremsstrahlungModel();
276 theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
277 brem->SetEmModel(theBremPenelope);
278
279 //Annihilation
282 theAnnPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
283 anni->AddEmModel(0, theAnnPenelope);
284
285 // register processes
286 ph->RegisterProcess(msc, particle);
287 ph->RegisterProcess(eioni, particle);
288 ph->RegisterProcess(brem, particle);
289 ph->RegisterProcess(ee, particle);
290 ph->RegisterProcess(anni, particle);
291 ph->RegisterProcess(ss, particle);
292
293 // generic ion
294 particle = G4GenericIon::GenericIon();
295 G4ionIonisation* ionIoni = new G4ionIonisation();
297 ph->RegisterProcess(hmsc, particle);
298 ph->RegisterProcess(ionIoni, particle);
299 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
300
301 // muons, hadrons, ions
303
304 // extra configuration
306}
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4Electron * Electron()
Definition: G4Electron.cc:93
static void ConstructCharged(G4hMultipleScattering *hmsc, G4NuclearStopping *nucStopping, G4bool isWVI=true)
Definition: G4EmBuilder.cc:219
static void PrepareEMPhysics()
Definition: G4EmBuilder.cc:375
G4double MaxNIELEnergy() const
G4double MscEnergyLimit() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Positron * Positron()
Definition: G4Positron.cc:93
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:788
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:774
void SetMinKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetEmModel(G4VMscModel *, G4int idx=0)
const G4String & GetPhysicsName() const
static constexpr double GeV

References G4VEmProcess::AddEmModel(), G4VEnergyLossProcess::AddEmModel(), G4EmBuilder::ConstructCharged(), G4Electron::Electron(), G4cout, G4endl, G4Gamma::Gamma(), G4GenericIon::GenericIon(), G4PhysicsListHelper::GetPhysicsListHelper(), G4VPhysicsConstructor::GetPhysicsName(), CLHEP::GeV, G4EmParameters::Instance(), G4EmParameters::MaxNIELEnergy(), G4EmParameters::MscEnergyLimit(), G4Positron::Positron(), G4EmBuilder::PrepareEMPhysics(), G4PhysicsListHelper::RegisterProcess(), G4VEmModel::SetActivationLowEnergyLimit(), G4VEmProcess::SetEmModel(), G4VEnergyLossProcess::SetEmModel(), G4VMultipleScattering::SetEmModel(), G4VEmModel::SetHighEnergyLimit(), G4VEmModel::SetLowEnergyLimit(), G4VEmProcess::SetMaxKinEnergy(), G4VEmProcess::SetMinKinEnergy(), and G4VPhysicsConstructor::verboseLevel.

◆ GetBuilders()

G4VPhysicsConstructor::PhysicsBuilder_V G4VPhysicsConstructor::GetBuilders ( ) const
protectedinherited

Definition at line 86 of file G4VPhysicsConstructor.cc.

87{
88 const auto& tls = *((subInstanceManager.offset[g4vpcInstanceID])._builders);
89 PhysicsBuilder_V copy(tls.size());
90 G4int i = 0;
91 for(const auto& el : tls)
92 {
93 copy[i++] = el;
94 }
95 return copy;
96}
int G4int
Definition: G4Types.hh:85
G4VPCData::PhysicsBuilders_V PhysicsBuilder_V
void copy(G4double dst[], const G4double src[], size_t size=G4FieldTrack::ncompSVEC)
Definition: G4FieldUtils.cc:98

References field_utils::copy(), G4VPhysicsConstructor::g4vpcInstanceID, G4VUPLSplitter< T >::offset, and G4VPhysicsConstructor::subInstanceManager.

◆ GetInstanceID()

G4int G4VPhysicsConstructor::GetInstanceID ( ) const
inlineinherited

◆ GetParticleIterator()

G4ParticleTable::G4PTblDicIterator * G4VPhysicsConstructor::GetParticleIterator ( ) const
protectedinherited

◆ GetPhysicsName()

const G4String & G4VPhysicsConstructor::GetPhysicsName ( ) const
inlineinherited

Definition at line 191 of file G4VPhysicsConstructor.hh.

192{
193 return namePhysics;
194}

References G4VPhysicsConstructor::namePhysics.

Referenced by G4EmDNAPhysics_option1::ConstructProcess(), G4EmDNAPhysics_option2::ConstructProcess(), G4EmDNAPhysics_option3::ConstructProcess(), G4EmDNAPhysics_option4::ConstructProcess(), G4EmDNAPhysics_option5::ConstructProcess(), G4EmDNAPhysics_option6::ConstructProcess(), G4EmDNAPhysics_option7::ConstructProcess(), G4EmDNAPhysics_option8::ConstructProcess(), G4EmDNAPhysics_stationary_option2::ConstructProcess(), G4EmDNAPhysics_stationary_option4::ConstructProcess(), G4EmDNAPhysics_stationary_option6::ConstructProcess(), G4EmDNAPhysics::ConstructProcess(), G4EmDNAPhysics_stationary::ConstructProcess(), G4EmLivermorePhysics::ConstructProcess(), G4EmLowEPPhysics::ConstructProcess(), ConstructProcess(), G4EmStandardPhysics::ConstructProcess(), G4EmStandardPhysics_option1::ConstructProcess(), G4EmStandardPhysics_option2::ConstructProcess(), G4EmStandardPhysics_option3::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4EmStandardPhysicsGS::ConstructProcess(), G4EmStandardPhysicsSS::ConstructProcess(), G4EmStandardPhysicsWVI::ConstructProcess(), G4ThermalNeutrons::ConstructProcess(), G4HadronPhysicsFTFP_BERT::DumpBanner(), G4HadronPhysicsQGSP_BERT::DumpBanner(), export_G4VPhysicsConstructor(), G4HadronDElasticPhysics::G4HadronDElasticPhysics(), G4HadronElasticPhysics::G4HadronElasticPhysics(), G4HadronElasticPhysicsHP::G4HadronElasticPhysicsHP(), G4HadronElasticPhysicsLEND::G4HadronElasticPhysicsLEND(), G4HadronElasticPhysicsPHP::G4HadronElasticPhysicsPHP(), G4HadronElasticPhysicsXS::G4HadronElasticPhysicsXS(), G4HadronHElasticPhysics::G4HadronHElasticPhysics(), G4IonElasticPhysics::G4IonElasticPhysics(), G4VModularPhysicsList::RegisterPhysics(), and G4VModularPhysicsList::ReplacePhysics().

◆ GetPhysicsType()

G4int G4VPhysicsConstructor::GetPhysicsType ( ) const
inlineinherited

◆ GetSubInstanceManager()

const G4VPCManager & G4VPhysicsConstructor::GetSubInstanceManager ( )
inlinestaticinherited

◆ GetVerboseLevel()

G4int G4VPhysicsConstructor::GetVerboseLevel ( ) const
inlineinherited

◆ RegisterProcess()

G4bool G4VPhysicsConstructor::RegisterProcess ( G4VProcess process,
G4ParticleDefinition particle 
)
inlineprotectedinherited

◆ SetPhysicsName()

void G4VPhysicsConstructor::SetPhysicsName ( const G4String name = "")
inlineinherited

Definition at line 186 of file G4VPhysicsConstructor.hh.

187{
189}
const char * name(G4int ptype)

References G4InuclParticleNames::name(), and G4VPhysicsConstructor::namePhysics.

Referenced by export_G4VPhysicsConstructor().

◆ SetPhysicsType()

void G4VPhysicsConstructor::SetPhysicsType ( G4int  val)
inlineinherited

Definition at line 196 of file G4VPhysicsConstructor.hh.

197{
198 if(val > 0) { typePhysics = val; }
199}

References G4VPhysicsConstructor::typePhysics.

Referenced by G4DecayPhysics::G4DecayPhysics(), G4EmDNAPhysics::G4EmDNAPhysics(), G4EmDNAPhysics_option1::G4EmDNAPhysics_option1(), G4EmDNAPhysics_option2::G4EmDNAPhysics_option2(), G4EmDNAPhysics_option3::G4EmDNAPhysics_option3(), G4EmDNAPhysics_option4::G4EmDNAPhysics_option4(), G4EmDNAPhysics_option5::G4EmDNAPhysics_option5(), G4EmDNAPhysics_option6::G4EmDNAPhysics_option6(), G4EmDNAPhysics_option7::G4EmDNAPhysics_option7(), G4EmDNAPhysics_option8::G4EmDNAPhysics_option8(), G4EmDNAPhysics_stationary_option2::G4EmDNAPhysics_stationary_option2(), G4EmDNAPhysics_stationary_option4::G4EmDNAPhysics_stationary_option4(), G4EmDNAPhysics_stationary_option6::G4EmDNAPhysics_stationary_option6(), G4EmExtraPhysics::G4EmExtraPhysics(), G4EmLivermorePhysics::G4EmLivermorePhysics(), G4EmLowEPPhysics::G4EmLowEPPhysics(), G4EmPenelopePhysics(), G4EmStandardPhysics::G4EmStandardPhysics(), G4EmStandardPhysics_option1::G4EmStandardPhysics_option1(), G4EmStandardPhysics_option2::G4EmStandardPhysics_option2(), G4EmStandardPhysics_option3::G4EmStandardPhysics_option3(), G4EmStandardPhysics_option4::G4EmStandardPhysics_option4(), G4EmStandardPhysicsGS::G4EmStandardPhysicsGS(), G4EmStandardPhysicsSS::G4EmStandardPhysicsSS(), G4EmStandardPhysicsWVI::G4EmStandardPhysicsWVI(), G4HadronElasticPhysics::G4HadronElasticPhysics(), G4HadronInelasticQBBC::G4HadronInelasticQBBC(), G4HadronPhysicsFTFP_BERT::G4HadronPhysicsFTFP_BERT(), G4HadronPhysicsQGSP_BERT::G4HadronPhysicsQGSP_BERT(), G4HadronPhysicsQGSP_BIC::G4HadronPhysicsQGSP_BIC(), G4IonINCLXXPhysics::G4IonINCLXXPhysics(), G4IonPhysics::G4IonPhysics(), G4IonPhysicsPHP::G4IonPhysicsPHP(), G4IonQMDPhysics::G4IonQMDPhysics(), G4NeutronTrackingCut::G4NeutronTrackingCut(), G4StepLimiterPhysics::G4StepLimiterPhysics(), G4StoppingPhysics::G4StoppingPhysics(), and G4StoppingPhysicsFritiofWithBinaryCascade::G4StoppingPhysicsFritiofWithBinaryCascade().

◆ SetVerboseLevel()

void G4VPhysicsConstructor::SetVerboseLevel ( G4int  value)
inlineinherited

◆ TerminateWorker()

void G4VPhysicsConstructor::TerminateWorker ( )
virtualinherited

Definition at line 105 of file G4VPhysicsConstructor.cc.

106{
107 if(subInstanceManager.offset[g4vpcInstanceID]._builders != nullptr)
108 {
109 std::for_each(subInstanceManager.offset[g4vpcInstanceID]._builders->begin(),
110 subInstanceManager.offset[g4vpcInstanceID]._builders->end(),
111 [](PhysicsBuilder_V::value_type bld) { delete bld; });
112 subInstanceManager.offset[g4vpcInstanceID]._builders->clear();
113 }
114}

References G4VPhysicsConstructor::g4vpcInstanceID, G4VUPLSplitter< T >::offset, and G4VPhysicsConstructor::subInstanceManager.

Referenced by G4VPhysicsConstructor::~G4VPhysicsConstructor().

Field Documentation

◆ g4vpcInstanceID

G4int G4VPhysicsConstructor::g4vpcInstanceID = 0
protectedinherited

◆ namePhysics

G4String G4VPhysicsConstructor::namePhysics = ""
protectedinherited

◆ subInstanceManager

G4VPCManager G4VPhysicsConstructor::subInstanceManager
staticprotectedinherited

◆ theParticleTable

G4ParticleTable* G4VPhysicsConstructor::theParticleTable = nullptr
protectedinherited

◆ typePhysics

G4int G4VPhysicsConstructor::typePhysics = 0
protectedinherited

◆ verboseLevel

G4int G4VPhysicsConstructor::verboseLevel = 0
protectedinherited

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