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

#include <G4EmDNAPhysics_stationary_option4.hh>

Inheritance diagram for G4EmDNAPhysics_stationary_option4:
G4VPhysicsConstructor

Public Member Functions

virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
 G4EmDNAPhysics_stationary_option4 (G4int ver, const G4String &name)
 
 G4EmDNAPhysics_stationary_option4 (G4int ver=1)
 
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_stationary_option4 ()
 

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 34 of file G4EmDNAPhysics_stationary_option4.hh.

Member Typedef Documentation

◆ PhysicsBuilder_V

Definition at line 149 of file G4VPhysicsConstructor.hh.

Constructor & Destructor Documentation

◆ G4EmDNAPhysics_stationary_option4() [1/2]

G4EmDNAPhysics_stationary_option4::G4EmDNAPhysics_stationary_option4 ( G4int  ver = 1)

Definition at line 88 of file G4EmDNAPhysics_stationary_option4.cc.

89 : G4VPhysicsConstructor("G4EmDNAPhysics_stationary_option4"), verbose(ver)
90{
92 param->SetDefaults();
93 param->SetFluo(true);
94 param->SetAuger(true);
95 param->SetDeexcitationIgnoreCut(true);
96 param->ActivateDNA();
97
99}
@ 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_stationary_option4() [2/2]

G4EmDNAPhysics_stationary_option4::G4EmDNAPhysics_stationary_option4 ( G4int  ver,
const G4String name 
)

◆ ~G4EmDNAPhysics_stationary_option4()

G4EmDNAPhysics_stationary_option4::~G4EmDNAPhysics_stationary_option4 ( )
virtual

Definition at line 110 of file G4EmDNAPhysics_stationary_option4.cc.

111{}

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

Implements G4VPhysicsConstructor.

Definition at line 115 of file G4EmDNAPhysics_stationary_option4.cc.

116{
117// bosons
119
120// leptons
123
124// baryons
126
128
129 G4DNAGenericIonsManager * genericIonsManager;
130 genericIonsManager=G4DNAGenericIonsManager::Instance();
131 genericIonsManager->GetIon("alpha++");
132 genericIonsManager->GetIon("alpha+");
133 genericIonsManager->GetIon("helium");
134 genericIonsManager->GetIon("hydrogen");
135
136}
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_stationary_option4::ConstructProcess ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 140 of file G4EmDNAPhysics_stationary_option4.cc.

141{
142 if(verbose > 1) {
143 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
144 }
146
147 auto myParticleIterator=GetParticleIterator();
148 myParticleIterator->reset();
149 while( (*myParticleIterator)() )
150 {
151 G4ParticleDefinition* particle = myParticleIterator->value();
152 G4String particleName = particle->GetParticleName();
153
154 if (particleName == "e-") {
155
156 // *** Elastic scattering (two alternative models available) ***
157
158 G4DNAElastic* theDNAElasticProcess =
159 new G4DNAElastic("e-_G4DNAElastic");
160 theDNAElasticProcess->SetEmModel(
162 ph->RegisterProcess(theDNAElasticProcess, particle);
163
164 // *** Excitation ***
165
166 G4DNAExcitation* theDNAExcitationProcess =
167 new G4DNAExcitation("e-_G4DNAExcitation");
168 theDNAExcitationProcess->SetEmModel(
171 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
172 ph->RegisterProcess(theDNAExcitationProcess, particle);
173
174 // *** Ionisation ***
175
176 G4DNAIonisation* theDNAIonisationProcess =
177 new G4DNAIonisation("e-_G4DNAIonisation");
178 theDNAIonisationProcess->SetEmModel(
181 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
182 ph->RegisterProcess(theDNAIonisationProcess, particle);
183
184 // *** Vibrational excitation ***
185
186 G4DNAVibExcitation* theDNAVibExcitationProcess =
187 new G4DNAVibExcitation("e-_G4DNAVibExcitation");
188 theDNAVibExcitationProcess->SetEmModel(
191 (theDNAVibExcitationProcess->EmModel()))->SelectStationary(true);
192 ph->RegisterProcess(theDNAVibExcitationProcess, particle);
193
194 // *** Attachment ***
195
196 G4DNAAttachment* theDNAAttachmentProcess =
197 new G4DNAAttachment("e-_G4DNAAttachment");
198 theDNAAttachmentProcess->SetEmModel(
201 (theDNAAttachmentProcess->EmModel()))->SelectStationary(true);
202 ph->RegisterProcess(theDNAAttachmentProcess, particle);
203
204 } else if ( particleName == "proton" ) {
205
206 // *** Elastic ***
207
208 G4DNAElastic* theDNAElasticProcess =
209 new G4DNAElastic("proton_G4DNAElastic");
210 theDNAElasticProcess->SetEmModel(
213 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
214 ph->RegisterProcess(theDNAElasticProcess, particle);
215
216 // *** Excitation ***
217
218 G4DNAExcitation* theDNAExcitationProcess =
219 new G4DNAExcitation("proton_G4DNAExcitation");
220
221 theDNAExcitationProcess->SetEmModel(
223 theDNAExcitationProcess->SetEmModel(
225
227 (theDNAExcitationProcess->EmModel()))->SetLowEnergyLimit(10*eV);
229 (theDNAExcitationProcess->EmModel()))->SetHighEnergyLimit(500*keV);
231 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
232
234 (theDNAExcitationProcess->EmModel(1)))->SetLowEnergyLimit(500*keV);
236 (theDNAExcitationProcess->EmModel(1)))->SetHighEnergyLimit(100*MeV);
238 (theDNAExcitationProcess->EmModel(1)))->SelectStationary(true);
239
240 ph->RegisterProcess(theDNAExcitationProcess, particle);
241
242 // *** Ionisation ***
243
244 G4DNAIonisation* theDNAIonisationProcess =
245 new G4DNAIonisation("proton_G4DNAIonisation");
246
247 theDNAIonisationProcess->SetEmModel(
249 theDNAIonisationProcess->SetEmModel(
251
253 (theDNAIonisationProcess->EmModel()))->SetLowEnergyLimit(0*eV);
255 (theDNAIonisationProcess->EmModel()))->SetHighEnergyLimit(500*keV);
257 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
258
260 (theDNAIonisationProcess->EmModel(1)))->SetLowEnergyLimit(500*keV);
262 (theDNAIonisationProcess->EmModel(1)))->SetHighEnergyLimit(100*MeV);
264 (theDNAIonisationProcess->EmModel(1)))->SelectStationary(true);
265 //
267 (theDNAIonisationProcess->EmModel(1)))->SelectFasterComputation(true);
268 //
269
270 ph->RegisterProcess(theDNAIonisationProcess, particle);
271
272 // *** Charge decrease ***
273
274 G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
275 new G4DNAChargeDecrease("proton_G4DNAChargeDecrease");
276 theDNAChargeDecreaseProcess->SetEmModel(
279 (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
280 ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
281
282 } else if ( particleName == "hydrogen" ) {
283
284 // *** Elastic ***
285
286 G4DNAElastic* theDNAElasticProcess =
287 new G4DNAElastic("hydrogen_G4DNAElastic");
288 theDNAElasticProcess->SetEmModel(
291 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
292 ph->RegisterProcess(theDNAElasticProcess, particle);
293
294 // *** Excitation ***
295
296 G4DNAExcitation* theDNAExcitationProcess =
297 new G4DNAExcitation("hydrogen_G4DNAExcitation");
298 theDNAExcitationProcess->SetEmModel(
301 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
302 ph->RegisterProcess(theDNAExcitationProcess, particle);
303
304 // *** Ionisation ***
305
306 G4DNAIonisation* theDNAIonisationProcess =
307 new G4DNAIonisation("hydrogen_G4DNAIonisation");
308 theDNAIonisationProcess->SetEmModel(
311 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
312 ph->RegisterProcess(theDNAIonisationProcess, particle);
313
314 // *** Charge increase ***
315
316 G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
317 new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease");
318 theDNAChargeIncreaseProcess->SetEmModel(
321 (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
322 ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
323
324 } else if ( particleName == "alpha" ) {
325
326 // *** Elastic ***
327
328 G4DNAElastic* theDNAElasticProcess =
329 new G4DNAElastic("alpha_G4DNAElastic");
330 theDNAElasticProcess->SetEmModel(
333 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
334 ph->RegisterProcess(theDNAElasticProcess, particle);
335
336 // *** Excitation ***
337
338 G4DNAExcitation* theDNAExcitationProcess =
339 new G4DNAExcitation("alpha_G4DNAExcitation");
340 theDNAExcitationProcess->SetEmModel(
343 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
344 ph->RegisterProcess(theDNAExcitationProcess, particle);
345
346 // *** Ionisation ***
347
348 G4DNAIonisation* theDNAIonisationProcess =
349 new G4DNAIonisation("alpha_G4DNAIonisation");
350 theDNAIonisationProcess->SetEmModel(
353 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
354 ph->RegisterProcess(theDNAIonisationProcess, particle);
355
356 // *** Charge decrease ***
357
358 G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
359 new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease");
360 theDNAChargeDecreaseProcess->SetEmModel(
363 (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
364 ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
365
366 } else if ( particleName == "alpha+" ) {
367
368 // *** Elastic ***
369
370 G4DNAElastic* theDNAElasticProcess =
371 new G4DNAElastic("alpha+_G4DNAElastic");
372 theDNAElasticProcess->SetEmModel(
375 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
376 ph->RegisterProcess(theDNAElasticProcess, particle);
377
378 // *** Excitation ***
379
380 G4DNAExcitation* theDNAExcitationProcess =
381 new G4DNAExcitation("alpha+_G4DNAExcitation");
382 theDNAExcitationProcess->SetEmModel(
385 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
386 ph->RegisterProcess(theDNAExcitationProcess, particle);
387
388 // *** Ionisation ***
389
390 G4DNAIonisation* theDNAIonisationProcess =
391 new G4DNAIonisation("alpha+_G4DNAIonisation");
392 theDNAIonisationProcess->SetEmModel(
395 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
396 ph->RegisterProcess(theDNAIonisationProcess, particle);
397
398 // *** Charge decrease ***
399
400 G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
401 new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease");
402 theDNAChargeDecreaseProcess->SetEmModel(
405 (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
406 ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
407
408 // *** Charge increase ***
409
410 G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
411 new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease");
412 theDNAChargeIncreaseProcess->SetEmModel(
415 (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
416 ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
417
418 } else if ( particleName == "helium" ) {
419
420 // *** Elastic ***
421
422 G4DNAElastic* theDNAElasticProcess =
423 new G4DNAElastic("helium_G4DNAElastic");
424 theDNAElasticProcess->SetEmModel(
427 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
428 ph->RegisterProcess(theDNAElasticProcess, particle);
429
430 // *** Excitation ***
431
432 G4DNAExcitation* theDNAExcitationProcess =
433 new G4DNAExcitation("helium_G4DNAExcitation");
434 theDNAExcitationProcess->SetEmModel(
437 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
438 ph->RegisterProcess(theDNAExcitationProcess, particle);
439
440 // *** Ionisation ***
441
442 G4DNAIonisation* theDNAIonisationProcess =
443 new G4DNAIonisation("helium_G4DNAIonisation");
444 theDNAIonisationProcess->SetEmModel(
447 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
448 ph->RegisterProcess(theDNAIonisationProcess, particle);
449
450 // *** Charge increase ***
451
452 G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
453 new G4DNAChargeIncrease("helium_G4DNAChargeIncrease");
454 theDNAChargeIncreaseProcess->SetEmModel(
457 (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
458 ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
459
460 } else if ( particleName == "GenericIon" ) {
461
462 // *** Ionisation ***
463
464 G4DNAIonisation* theDNAIonisationProcess =
465 new G4DNAIonisation("GenericIon_G4DNAIonisation");
466 theDNAIonisationProcess->SetEmModel(
469 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
470 ph->RegisterProcess(theDNAIonisationProcess, particle);
471
472 }
473
474 // Warning : the following particles and processes are needed by EM Physics
475 // builders
476 // They are taken from the default Livermore Physics list
477 // These particles are currently not handled by Geant4-DNA
478
479 // e+
480
481 else if (particleName == "e+") {
482
483 // Identical to G4EmStandardPhysics_stationary
484
487 G4eIonisation* eIoni = new G4eIonisation();
488 eIoni->SetStepFunction(0.2, 100*um);
489
490 ph->RegisterProcess(msc, particle);
491 ph->RegisterProcess(eIoni, particle);
492 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
493 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
494
495 } else if (particleName == "gamma") {
496
497 // photoelectric effect - Livermore model only
498 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
499 thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
500 ph->RegisterProcess(thePhotoElectricEffect, particle);
501
502 // Compton scattering - Livermore model only
503 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
504 theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
505 ph->RegisterProcess(theComptonScattering, particle);
506
507 // gamma conversion - Livermore model below 80 GeV
508 G4GammaConversion* theGammaConversion = new G4GammaConversion();
509 theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
510 ph->RegisterProcess(theGammaConversion, particle);
511
512 // default Rayleigh scattering is Livermore
513 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
514 ph->RegisterProcess(theRayleigh, particle);
515 }
516
517 // Warning : end of particles and processes are needed by EM Physics build.
518
519 }
520
521 // Deexcitation
522 //
525}
G4DNABornExcitationModel1 G4DNABornExcitationModel
#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
void SetAtomDeexcitation(G4VAtomDeexcitation *)
static G4LossTableManager * Instance()
const G4String & GetParticleName() const
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
G4VEmModel * EmModel(size_t index=0) const
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 G4VEmProcess::EmModel(), eV, fUseDistanceToBoundary, G4cout, G4DNABornIonisationModel, G4endl, G4VPhysicsConstructor::GetParticleIterator(), G4ParticleDefinition::GetParticleName(), G4PhysicsListHelper::GetPhysicsListHelper(), G4VPhysicsConstructor::GetPhysicsName(), G4LossTableManager::Instance(), keV, MeV, G4PhysicsListHelper::RegisterProcess(), G4LossTableManager::SetAtomDeexcitation(), G4VEmProcess::SetEmModel(), 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(), 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(), 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_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_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_stationary_option4::verbose
private

Definition at line 48 of file G4EmDNAPhysics_stationary_option4.hh.

Referenced by ConstructProcess().

◆ verboseLevel

G4int G4VPhysicsConstructor::verboseLevel = 0
protectedinherited

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