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

#include <G4EmLowEPPhysics.hh>

Inheritance diagram for G4EmLowEPPhysics:
G4VPhysicsConstructor

Public Member Functions

void ConstructParticle () override
 
void ConstructProcess () override
 
 G4EmLowEPPhysics (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 ()
 
 ~G4EmLowEPPhysics () 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 48 of file G4EmLowEPPhysics.hh.

Member Typedef Documentation

◆ PhysicsBuilder_V

Definition at line 149 of file G4VPhysicsConstructor.hh.

Constructor & Destructor Documentation

◆ G4EmLowEPPhysics()

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

Definition at line 121 of file G4EmLowEPPhysics.cc.

122 : G4VPhysicsConstructor("G4EmLowEPPhysics")
123{
124 SetVerboseLevel(ver);
126 param->SetDefaults();
127 param->SetVerbose(ver);
128 param->SetMinEnergy(100*CLHEP::eV);
130 param->SetNumberOfBinsPerDecade(20);
132 param->SetStepFunction(0.2, 100*CLHEP::um);
133 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
134 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
135 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
136 param->SetUseMottCorrection(true);
137 param->SetMscRangeFactor(0.04);
138 param->SetMuHadLateralDisplacement(true);
139 param->SetFluo(true);
140 param->SetAuger(true);
141 param->SetUseICRU90Data(true);
143}
@ bElectromagnetic
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
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 SetStepFunctionIons(G4double v1, G4double v2)
void SetAuger(G4bool val)
void SetUseICRU90Data(G4bool 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 eV

References G4EmParameters::ActivateAngularGeneratorForIonisation(), bElectromagnetic, CLHEP::eV, G4EmParameters::Instance(), G4EmParameters::SetAuger(), G4EmParameters::SetDefaults(), G4EmParameters::SetFluo(), G4EmParameters::SetLowestElectronEnergy(), G4EmParameters::SetMinEnergy(), G4EmParameters::SetMscRangeFactor(), G4EmParameters::SetMuHadLateralDisplacement(), G4EmParameters::SetNumberOfBinsPerDecade(), G4VPhysicsConstructor::SetPhysicsType(), G4EmParameters::SetStepFunction(), G4EmParameters::SetStepFunctionIons(), G4EmParameters::SetStepFunctionLightIons(), G4EmParameters::SetStepFunctionMuHad(), G4EmParameters::SetUseICRU90Data(), G4EmParameters::SetUseMottCorrection(), G4EmParameters::SetVerbose(), G4VPhysicsConstructor::SetVerboseLevel(), and CLHEP::um.

◆ ~G4EmLowEPPhysics()

G4EmLowEPPhysics::~G4EmLowEPPhysics ( )
override

Definition at line 147 of file G4EmLowEPPhysics.cc.

148{}

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 G4EmLowEPPhysics::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 152 of file G4EmLowEPPhysics.cc.

153{
154 // minimal set of particles for EM physics
156}
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:347

References G4EmBuilder::ConstructMinimalEmSet().

◆ ConstructProcess()

void G4EmLowEPPhysics::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 160 of file G4EmLowEPPhysics.cc.

161{
162 if(verboseLevel > 1) {
163 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
164 }
166
169
170 // common processes
171 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
172
173 // nuclear stopping
174 G4double nielEnergyLimit = param->MaxNIELEnergy();
175 G4NuclearStopping* pnuc = nullptr;
176 if(nielEnergyLimit > 0.0) {
177 pnuc = new G4NuclearStopping();
178 pnuc->SetMaxKinEnergy(nielEnergyLimit);
179 }
180
181 // high energy limit for e+- scattering models and bremsstrahlung
182 G4double highEnergyLimit = param->MscEnergyLimit();
183
184 // Add gamma EM Processes
186
187 // Photoelectric absorption
189 G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
191 pe->SetEmModel(theLivermorePEModel);
192
193 // Compton scattering - Polarised Monash model
195 G4VEmModel* theLowEPCSModel = new G4LowEPPolarizedComptonModel();
196 cs->SetEmModel(theLowEPCSModel);
197
198 // gamma conversion - 5D model below 80 GeV with Livermore x-sections
199 G4GammaConversion* theGammaConversion = new G4GammaConversion();
200 G4VEmModel* conv = new G4BetheHeitler5DModel();
201 theGammaConversion->SetEmModel(conv);
202
203 // default Rayleigh scattering is Livermore
204 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
205 G4VEmModel* theLivermorePRSModel = new G4LivermorePolarizedRayleighModel();
206 theRayleigh->SetEmModel(theLivermorePRSModel);
207
208 ph->RegisterProcess(pe, particle);
209 ph->RegisterProcess(cs, particle);
210 ph->RegisterProcess(theGammaConversion, particle);
211 ph->RegisterProcess(theRayleigh, particle);
212
213 // e-
214 particle = G4Electron::Electron();
215
216 // multiple scattering
220 msc1->SetHighEnergyLimit(highEnergyLimit);
221 msc2->SetLowEnergyLimit(highEnergyLimit);
222 msc->SetEmModel(msc1);
223 msc->SetEmModel(msc2);
224
227 ss->SetEmModel(ssm);
228 ss->SetMinKinEnergy(highEnergyLimit);
229 ssm->SetLowEnergyLimit(highEnergyLimit);
230 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
231
232 // Ionisation - Livermore should be used only for low energies
233 G4eIonisation* eioni = new G4eIonisation();
235 theIoniLivermore->SetHighEnergyLimit(0.1*CLHEP::MeV);
236 eioni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
237
238 // Bremsstrahlung
244 brem->SetEmModel(br1);
245 brem->SetEmModel(br2);
247
249
250 // register processes
251 ph->RegisterProcess(msc, particle);
252 ph->RegisterProcess(ss, particle);
253 ph->RegisterProcess(eioni, particle);
254 ph->RegisterProcess(brem, particle);
255 ph->RegisterProcess(ee, particle);
256
257 // e+
258 particle = G4Positron::Positron();
259
260 // multiple scattering
261 msc = new G4eMultipleScattering();
262 msc1 = new G4GoudsmitSaundersonMscModel();
263 msc2 = new G4WentzelVIModel();
264 msc1->SetHighEnergyLimit(highEnergyLimit);
265 msc2->SetLowEnergyLimit(highEnergyLimit);
266 msc->SetEmModel(msc1);
267 msc->SetEmModel(msc2);
268
269 ssm = new G4eCoulombScatteringModel();
270 ss = new G4CoulombScattering();
271 ss->SetEmModel(ssm);
272 ss->SetMinKinEnergy(highEnergyLimit);
273 ssm->SetLowEnergyLimit(highEnergyLimit);
274 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
275
276 // Standard ionisation
277 eioni = new G4eIonisation();
279 pen->SetHighEnergyLimit(0.1*MeV);
280 eioni->AddEmModel(0, pen, new G4UniversalFluctuation());
281
282 // Bremsstrahlung
283 brem = new G4eBremsstrahlung();
284 br1 = new G4SeltzerBergerModel();
285 br2 = new G4eBremsstrahlungRelModel();
288 brem->SetEmModel(br1);
289 brem->SetEmModel(br2);
291
292 // register processes
293 ph->RegisterProcess(msc, particle);
294 ph->RegisterProcess(ss, particle);
295 ph->RegisterProcess(eioni, particle);
296 ph->RegisterProcess(brem, particle);
297 ph->RegisterProcess(ee, particle);
298 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
299
300 // generic ion
301 particle = G4GenericIon::GenericIon();
302
303 G4ionIonisation* ionIoni = new G4ionIonisation();
305 ionIoni->SetEmModel(mod);
306
307 ph->RegisterProcess(hmsc, particle);
308 ph->RegisterProcess(ionIoni, particle);
309 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
310
311 // muons, hadrons ions
313
314 // extra configuration
316}
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
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 SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:628
void SetMinKinEnergy(G4double e)
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
static constexpr double MeV

References G4VEnergyLossProcess::AddEmModel(), G4EmBuilder::ConstructCharged(), G4Electron::Electron(), G4cout, G4endl, G4Gamma::Gamma(), G4GenericIon::GenericIon(), G4PhysicsListHelper::GetPhysicsListHelper(), G4VPhysicsConstructor::GetPhysicsName(), CLHEP::GeV, GeV, G4EmParameters::Instance(), G4EmParameters::MaxNIELEnergy(), CLHEP::MeV, MeV, G4EmParameters::MscEnergyLimit(), G4Positron::Positron(), G4EmBuilder::PrepareEMPhysics(), G4PhysicsListHelper::RegisterProcess(), G4VEmModel::SetActivationLowEnergyLimit(), G4VEmModel::SetAngularDistribution(), 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(), ConstructProcess(), G4EmPenelopePhysics::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(), G4EmPenelopePhysics::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: