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

#include <G4EmDNAPhysics_option2.hh>

Inheritance diagram for G4EmDNAPhysics_option2:
G4VPhysicsConstructor

Public Member Functions

virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
 G4EmDNAPhysics_option2 (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 ()
 
virtual ~G4EmDNAPhysics_option2 ()
 

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
 

Private Attributes

G4int verbose
 

Detailed Description

Definition at line 36 of file G4EmDNAPhysics_option2.hh.

Member Typedef Documentation

◆ PhysicsBuilder_V

Definition at line 149 of file G4VPhysicsConstructor.hh.

Constructor & Destructor Documentation

◆ G4EmDNAPhysics_option2()

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

Definition at line 89 of file G4EmDNAPhysics_option2.cc.

90 : G4VPhysicsConstructor("G4EmDNAPhysics_option2"), verbose(ver)
91{
93
94 param->SetDefaults();
95 param->SetFluo(true);
96 param->SetAuger(true);
97 param->SetDeexcitationIgnoreCut(true);
98 param->ActivateDNA();
99
101}
@ bElectromagnetic
static G4EmParameters * Instance()
void SetDeexcitationIgnoreCut(G4bool val)
void SetFluo(G4bool val)
void SetAuger(G4bool val)
G4VPhysicsConstructor(const G4String &="")

References G4EmParameters::ActivateDNA(), bElectromagnetic, G4EmParameters::Instance(), G4EmParameters::SetAuger(), G4EmParameters::SetDeexcitationIgnoreCut(), G4EmParameters::SetDefaults(), G4EmParameters::SetFluo(), and G4VPhysicsConstructor::SetPhysicsType().

◆ ~G4EmDNAPhysics_option2()

G4EmDNAPhysics_option2::~G4EmDNAPhysics_option2 ( )
virtual

Definition at line 105 of file G4EmDNAPhysics_option2.cc.

106{}

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 G4EmDNAPhysics_option2::ConstructParticle ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 110 of file G4EmDNAPhysics_option2.cc.

111{
112// bosons
114
115// leptons
118
119// baryons
121
123
124 G4DNAGenericIonsManager * genericIonsManager;
125 genericIonsManager=G4DNAGenericIonsManager::Instance();
126 genericIonsManager->GetIon("alpha++");
127 genericIonsManager->GetIon("alpha+");
128 genericIonsManager->GetIon("helium");
129 genericIonsManager->GetIon("hydrogen");
130
131}
static G4DNAGenericIonsManager * Instance(void)
G4ParticleDefinition * GetIon(const G4String &name)
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:87
static G4Positron * Positron()
Definition: G4Positron.cc:93
static G4Proton * Proton()
Definition: G4Proton.cc:92

References G4Electron::Electron(), G4Gamma::Gamma(), G4GenericIon::GenericIonDefinition(), G4DNAGenericIonsManager::GetIon(), G4DNAGenericIonsManager::Instance(), G4Positron::Positron(), and G4Proton::Proton().

◆ ConstructProcess()

void G4EmDNAPhysics_option2::ConstructProcess ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 135 of file G4EmDNAPhysics_option2.cc.

136{
137 if(verbose > 1) {
138 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
139 }
141
142 auto myParticleIterator=GetParticleIterator();
143 myParticleIterator->reset();
144 while( (*myParticleIterator)() )
145 {
146 G4ParticleDefinition* particle = myParticleIterator->value();
147 G4String particleName = particle->GetParticleName();
148
149 if (particleName == "e-") {
150
151 G4DNAElectronSolvation* solvation =
152 new G4DNAElectronSolvation("e-_G4DNAElectronSolvation");
153
155 therm->SetHighEnergyLimit(7.4*eV); // limit of the Champion's model
156 //therm->SetHighEnergyLimit(10*eV); // limit of the ELSEPA model
157 solvation->SetEmModel(therm);
158 ph->RegisterProcess(solvation, particle);
159
160 // *** Elastic scattering (two alternative models available) ***
161
162 G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
163 theDNAElasticProcess->SetEmModel(new G4DNAChampionElasticModel());
164
165 // or alternative model
166 //theDNAElasticProcess->SetEmModel(new G4DNAELSEPAElasticModel());
167 //theDNAElasticProcess->SetEmModel(new G4DNAScreenedRutherfordElasticModel());
168
169 ph->RegisterProcess(theDNAElasticProcess, particle);
170
171 // *** Excitation ***
172 ph->RegisterProcess(new G4DNAExcitation("e-_G4DNAExcitation"), particle);
173
174 // *** Ionisation ***
175 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("e-_G4DNAIonisation");
177 mod->SelectFasterComputation(true);
178 theDNAIonisationProcess->SetEmModel(mod);
179 ph->RegisterProcess(theDNAIonisationProcess, particle);
180
181 // *** Vibrational excitation ***
182 ph->RegisterProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"), particle);
183
184 // *** Attachment ***
185 ph->RegisterProcess(new G4DNAAttachment("e-_G4DNAAttachment"), particle);
186
187 } else if ( particleName == "proton" ) {
188
189 ph->RegisterProcess(new G4DNAElastic("proton_G4DNAElastic"), particle);
190
191 ph->RegisterProcess(new G4DNAExcitation("proton_G4DNAExcitation"), particle);
192
193 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("proton_G4DNAIonisation");
194
196 mod1->SetLowEnergyLimit(0*eV);
197 mod1->SetHighEnergyLimit(500*keV);
198
200 mod2->SetLowEnergyLimit(500*keV);
201 mod2->SetHighEnergyLimit(100*MeV);
202 mod2->SelectFasterComputation(true);
203
204 theDNAIonisationProcess->SetEmModel(mod1);
205 theDNAIonisationProcess->SetEmModel(mod2);
206
207 ph->RegisterProcess(theDNAIonisationProcess, particle);
208
209 ph->RegisterProcess(new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"), particle);
210
211 } else if ( particleName == "hydrogen" ) {
212
213 ph->RegisterProcess(new G4DNAElastic("hydrogen_G4DNAElastic"), particle);
214
215 ph->RegisterProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"), particle);
216
217 //ph->RegisterProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"), particle);
218 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("hydrogen_G4DNAIonisation");
219 theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationExtendedModel());
220 ph->RegisterProcess(theDNAIonisationProcess, particle);
221
222 ph->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), particle);
223
224 } else if ( particleName == "alpha" ) {
225
226 ph->RegisterProcess(new G4DNAElastic("alpha_G4DNAElastic"), particle);
227
228 ph->RegisterProcess(new G4DNAExcitation("alpha_G4DNAExcitation"), particle);
229
230 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("alpha_G4DNAIonisation");
231 theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationExtendedModel());
232 ph->RegisterProcess(theDNAIonisationProcess, particle);
233
234 ph->RegisterProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"), particle);
235
236 } else if ( particleName == "alpha+" ) {
237
238 ph->RegisterProcess(new G4DNAElastic("alpha+_G4DNAElastic"), particle);
239
240 ph->RegisterProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"), particle);
241
242 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("alpha+_G4DNAIonisation");
243 theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationExtendedModel());
244 ph->RegisterProcess(theDNAIonisationProcess, particle);
245
246 ph->RegisterProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"), particle);
247 ph->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), particle);
248
249 } else if ( particleName == "helium" ) {
250
251 ph->RegisterProcess(new G4DNAElastic("helium_G4DNAElastic"), particle);
252
253 ph->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), particle);
254
255 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("helium_G4DNAIonisation");
256 theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationExtendedModel());
257 ph->RegisterProcess(theDNAIonisationProcess, particle);
258
259 ph->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), particle);
260
261 } else if ( particleName == "GenericIon" ) {
262 ph->RegisterProcess(new G4DNAIonisation("GenericIon_G4DNAIonisation"), particle);
263 }
264
265 // Warning : the following particles and processes are needed by EM Physics builders
266 // They are taken from the default Livermore Physics list
267 // These particles are currently not handled by Geant4-DNA
268
269 // e+
270
271 else if (particleName == "e+") {
272
275 G4eIonisation* eIoni = new G4eIonisation();
276 eIoni->SetStepFunction(0.2, 100*um);
277
278 ph->RegisterProcess(msc, particle);
279 ph->RegisterProcess(eIoni, particle);
280 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
281 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
282
283 } else if (particleName == "gamma") {
284
285 // photoelectric effect - Livermore model only
286 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
287 thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
288 ph->RegisterProcess(thePhotoElectricEffect, particle);
289
290 // Compton scattering - Livermore model only
291 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
292 theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
293 ph->RegisterProcess(theComptonScattering, particle);
294
295 // gamma conversion - Livermore model below 80 GeV
296 G4GammaConversion* theGammaConversion = new G4GammaConversion();
297 theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
298 ph->RegisterProcess(theGammaConversion, particle);
299
300 // default Rayleigh scattering is Livermore
301 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
302 ph->RegisterProcess(theRayleigh, particle);
303 }
304
305 // Warning : end of particles and processes are needed by EM Physics builders
306
307 }
308
309 // Deexcitation
310 //
313}
#define G4DNABornIonisationModel
@ fUseDistanceToBoundary
static constexpr double keV
Definition: G4SIunits.hh:202
static constexpr double eV
Definition: G4SIunits.hh:201
static constexpr double um
Definition: G4SIunits.hh:93
static constexpr double MeV
Definition: G4SIunits.hh:200
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4VEmModel * GetMacroDefinedModel()
One step thermalization model can be chosen via macro using /process/dna/e-SolvationSubType Ritchie19...
void SetAtomDeexcitation(G4VAtomDeexcitation *)
static G4LossTableManager * Instance()
const G4String & GetParticleName() const
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:774
void SetEmModel(G4VEmModel *, G4int index=0)
void SetStepFunction(G4double v1, G4double v2)
void SetStepLimitType(G4MscStepLimitType val)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
const G4String & GetPhysicsName() const

References eV, fUseDistanceToBoundary, G4cout, G4DNABornIonisationModel, G4endl, G4DNASolvationModelFactory::GetMacroDefinedModel(), G4VPhysicsConstructor::GetParticleIterator(), G4ParticleDefinition::GetParticleName(), G4PhysicsListHelper::GetPhysicsListHelper(), G4VPhysicsConstructor::GetPhysicsName(), G4LossTableManager::Instance(), keV, MeV, G4PhysicsListHelper::RegisterProcess(), G4LossTableManager::SetAtomDeexcitation(), G4VEmProcess::SetEmModel(), G4VEmModel::SetHighEnergyLimit(), G4VEmModel::SetLowEnergyLimit(), G4VEnergyLossProcess::SetStepFunction(), G4VMultipleScattering::SetStepLimitType(), um, and verbose.

◆ 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(), 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(), 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_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::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

◆ verbose

G4int G4EmDNAPhysics_option2::verbose
private

Definition at line 48 of file G4EmDNAPhysics_option2.hh.

Referenced by ConstructProcess().

◆ verboseLevel

G4int G4VPhysicsConstructor::verboseLevel = 0
protectedinherited

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