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

#include <G4EmDNAPhysicsActivator.hh>

Inheritance diagram for G4EmDNAPhysicsActivator:
G4VPhysicsConstructor

Public Member Functions

void ConstructParticle () override
 
void ConstructProcess () override
 
 G4EmDNAPhysicsActivator (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 ()
 
 ~G4EmDNAPhysicsActivator () 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
 

Private Member Functions

void AddElectronModels0 (const G4String &region, G4bool emsc, G4double elimel)
 
void AddElectronModels2 (const G4String &region, G4bool emsc, G4double elimel)
 
void AddElectronModels4 (const G4String &region, G4bool emsc, G4double elimel)
 
void AddElectronModels4a (const G4String &region, G4bool emsc, G4double elimel)
 
void AddElectronModels6 (const G4String &region, G4bool emsc, G4double elimel)
 
void AddElectronModels6a (const G4String &region, G4bool emsc, G4double elimel)
 
void AddElectronModels7 (const G4String &region, G4bool emsc, G4double elimel)
 
void AddGenericIonModels0 (const G4String &region, G4double pminbb)
 
void AddHeliumModels0 (const G4String &region, G4bool a1msc, G4bool a2msc, G4double elimel, G4double pminbb, G4double pmax)
 
void AddProtonModels0 (const G4String &region, G4bool pmsc, G4double elimel, G4double pminbb, G4double pmax)
 
void DeactivateNuclearStopping (G4ProcessManager *, G4double elimel)
 
G4bool HasMsc (G4ProcessManager *) const
 
G4bool IsVerbose () const
 

Private Attributes

G4EmParameterstheParameters
 
G4int verbose
 

Detailed Description

Definition at line 39 of file G4EmDNAPhysicsActivator.hh.

Member Typedef Documentation

◆ PhysicsBuilder_V

Definition at line 149 of file G4VPhysicsConstructor.hh.

Constructor & Destructor Documentation

◆ G4EmDNAPhysicsActivator()

G4EmDNAPhysicsActivator::G4EmDNAPhysicsActivator ( G4int  ver = 1)
explicit

Definition at line 101 of file G4EmDNAPhysicsActivator.cc.

102 : G4VPhysicsConstructor("G4EmDNAPhysicsActivator"), verbose(ver)
103{
106}
static G4EmParameters * Instance()
G4VPhysicsConstructor(const G4String &="")

References G4EmParameters::ActivateDNA(), G4EmParameters::Instance(), and theParameters.

◆ ~G4EmDNAPhysicsActivator()

G4EmDNAPhysicsActivator::~G4EmDNAPhysicsActivator ( )
override

Definition at line 110 of file G4EmDNAPhysicsActivator.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().

◆ AddElectronModels0()

void G4EmDNAPhysicsActivator::AddElectronModels0 ( const G4String region,
G4bool  emsc,
G4double  elimel 
)
private

Definition at line 416 of file G4EmDNAPhysicsActivator.cc.

419{
420 G4EmConfigurator* em_config =
422 G4VEmModel* mod;
423
424 static const G4double elowest = 7.4 * eV; // seems to be option dependent, so I moved it here - MJPietrzak
425 static const G4double elimin = 1 * MeV;
426 static const G4double elimvb = 100 * eV;
427 static const G4double elimat = 13 * eV;
428 static const G4double elim1 = 10 * keV;
429
430 // for e- 100 MeV is a limit between different msc models
432
433 if(emsc) {
434 G4UrbanMscModel* msc = new G4UrbanMscModel();
435 msc->SetActivationLowEnergyLimit(elimel);
436 G4double emaxmsc = std::min(100*MeV, emax);
437 em_config->SetExtraEmModel("e-", "msc", msc, reg, 0.0, emaxmsc);
438 } else {
439 mod = new G4eCoulombScatteringModel();
440 mod->SetActivationLowEnergyLimit(elimel);
441 em_config->SetExtraEmModel("e-", "CoulombScat", mod, reg, 0.0, emax);
442 }
443
444 // cuts and solvation
446 em_config->SetExtraEmModel("e-", "e-_G4DNAElectronSolvation",
447 mod, reg, 0., elowest);
448
449 // elastic
450 mod = new G4DNAChampionElasticModel();
451 em_config->SetExtraEmModel("e-", "e-_G4DNAElastic",
452 mod, reg, 0.0, elimel);
453 // ionisation
454 mod = new G4MollerBhabhaModel();
455 mod->SetActivationLowEnergyLimit(elimin);
456 em_config->SetExtraEmModel("e-", "eIoni",
457 mod, reg, 0.0, emax,
459
460 mod = new G4DNABornIonisationModel();
461 em_config->SetExtraEmModel("e-", "e-_G4DNAIonisation",
462 mod, reg, elim1, elimin);
463
465 em_config->SetExtraEmModel("e-", "e-_G4DNAIonisation",
466 mod, reg, 0.0, elim1);
467
468 // exc
470 em_config->SetExtraEmModel("e-", "e-_G4DNAExcitation",
471 mod, reg, 0.0, elim1);
472
473 mod = new G4DNABornExcitationModel();
474 em_config->SetExtraEmModel("e-", "e-_G4DNAExcitation",
475 mod, reg, elim1, elimin);
476
477 mod = new G4DNASancheExcitationModel();
478 em_config->SetExtraEmModel("e-", "e-_G4DNAVibExcitation",
479 mod, reg, 0.0, elimvb);
480
481 mod = new G4DNAMeltonAttachmentModel();
482 em_config->SetExtraEmModel("e-", "e-_G4DNAAttachment",
483 mod, reg, 0.0, elimat);
484
485}
static const G4double emax
G4DNABornExcitationModel1 G4DNABornExcitationModel
#define G4DNABornIonisationModel
G4TDNAOneStepThermalizationModel< DNA::Penetration::Meesungnoen2002 > G4DNAOneStepThermalizationModel
static const G4double reg
static constexpr double keV
Definition: G4SIunits.hh:202
static constexpr double eV
Definition: G4SIunits.hh:201
static constexpr double MeV
Definition: G4SIunits.hh:200
double G4double
Definition: G4Types.hh:83
void SetExtraEmModel(const G4String &particleName, const G4String &processName, G4VEmModel *, const G4String &regionName="", G4double emin=0.0, G4double emax=DBL_MAX, G4VEmFluctuationModel *fm=nullptr)
G4double MaxKinEnergy() const
static G4LossTableManager * Instance()
G4EmConfigurator * EmConfigurator()
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:788
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References emax, G4LossTableManager::EmConfigurator(), eV, G4DNABornIonisationModel, G4LossTableManager::Instance(), keV, G4EmParameters::MaxKinEnergy(), MeV, G4INCL::Math::min(), reg, G4VEmModel::SetActivationLowEnergyLimit(), G4EmConfigurator::SetExtraEmModel(), and theParameters.

Referenced by ConstructProcess().

◆ AddElectronModels2()

void G4EmDNAPhysicsActivator::AddElectronModels2 ( const G4String region,
G4bool  emsc,
G4double  elimel 
)
private

Definition at line 489 of file G4EmDNAPhysicsActivator.cc.

492{
493 G4EmConfigurator *em_config =
495 G4VEmModel *mod;
496
497 static const G4double elowest = 7.4 * eV; // seems to be option dependent, so I moved it here - MJPietrzak
498 static const G4double elimin = 1 * MeV;
499 static const G4double elimvb = 100 * eV;
500 static const G4double elimat = 13 * eV;
501
502
503 // for e- 100 MeV is a limit between different msc models
505
506 if (emsc) {
507 G4UrbanMscModel *msc = new G4UrbanMscModel();
508 msc->SetActivationLowEnergyLimit(elimel);
509 G4double emaxmsc = std::min(100 * MeV, emax);
510 em_config->SetExtraEmModel("e-", "msc", msc, reg, 0.0, emaxmsc);
511 } else {
512 mod = new G4eCoulombScatteringModel();
513 mod->SetActivationLowEnergyLimit(elimel);
514 em_config->SetExtraEmModel("e-", "CoulombScat", mod, reg, 0.0, emax);
515 }
516
517 // cuts and solvation
519 em_config->SetExtraEmModel("e-", "e-_G4DNAElectronSolvation",
520 mod, reg, 0., elowest);
521
522 // elastic
523 mod = new G4DNAChampionElasticModel();
524 em_config->SetExtraEmModel("e-", "e-_G4DNAElastic",
525 mod, reg, 0.0, elimel);
526
527
528 // ionisation
529 mod = new G4MollerBhabhaModel();
530 mod->SetActivationLowEnergyLimit(elimin);
531 em_config->SetExtraEmModel("e-", "eIoni",
532 mod, reg, 0.0, emax,
534
535 mod = new G4DNABornIonisationModel();
536 em_config->SetExtraEmModel("e-", "e-_G4DNAIonisation",
537 mod, reg, elowest, elimin);
538
539 // excitation
540 mod = new G4DNABornExcitationModel();
541 em_config->SetExtraEmModel("e-", "e-_G4DNAExcitation",
542 mod, reg, 0., elimin);
543
544 // vib excitation
545 mod = new G4DNASancheExcitationModel();
546 em_config->SetExtraEmModel("e-", "e-_G4DNAVibExcitation",
547 mod, reg, 0.0, elimvb);
548
549 // attachment
550 mod = new G4DNAMeltonAttachmentModel();
551 em_config->SetExtraEmModel("e-", "e-_G4DNAAttachment",
552 mod, reg, 0.0, elimat);
553
554}

References emax, G4LossTableManager::EmConfigurator(), eV, G4DNABornIonisationModel, G4LossTableManager::Instance(), G4EmParameters::MaxKinEnergy(), MeV, G4INCL::Math::min(), reg, G4VEmModel::SetActivationLowEnergyLimit(), G4EmConfigurator::SetExtraEmModel(), and theParameters.

Referenced by ConstructProcess().

◆ AddElectronModels4()

void G4EmDNAPhysicsActivator::AddElectronModels4 ( const G4String region,
G4bool  emsc,
G4double  elimel 
)
private

Definition at line 559 of file G4EmDNAPhysicsActivator.cc.

562{
563 G4EmConfigurator *em_config =
565 G4VEmModel *mod;
566
567 static const G4double elowest = 10 * eV; // seems to be option dependent, so I moved it here - MJPietrzak
568 static const G4double elimin = 1 * MeV;
569// static const G4double elimvb = 100 * eV;
570// static const G4double elimat = 13 * eV;
571
572
573 // for e- 100 MeV is a limit between different msc models
575
576 if (emsc) {
577 G4UrbanMscModel *msc = new G4UrbanMscModel();
578 msc->SetActivationLowEnergyLimit(elimel);
579 G4double emaxmsc = std::min(100 * MeV, emax);
580 em_config->SetExtraEmModel("e-", "msc", msc, reg, 0.0, emaxmsc);
581 } else {
582 mod = new G4eCoulombScatteringModel();
583 mod->SetActivationLowEnergyLimit(elimel);
584 em_config->SetExtraEmModel("e-", "CoulombScat", mod, reg, 0.0, emax);
585 }
586
587 // cuts and solvation
589 em_config->SetExtraEmModel("e-", "e-_G4DNAElectronSolvation",
590 mod, reg, 0., elowest);
591
592 // elastic
594 em_config->SetExtraEmModel("e-", "e-_G4DNAElastic",
595 mod, reg, 0.0, elimel);
596
597 // ionisation
598 mod = new G4MollerBhabhaModel();
599 mod->SetActivationLowEnergyLimit(elimin);
600 em_config->SetExtraEmModel("e-", "eIoni",
601 mod, reg, 0.0, emax,
603
605 em_config->SetExtraEmModel("e-", "e-_G4DNAIonisation",
606 mod, reg, elowest, elimin);
607
608 // excitation
610 em_config->SetExtraEmModel("e-", "e-_G4DNAExcitation",
611 mod, reg, 0., elimin);
612
613
614// todo - MJPietrzak
615// I don't understand why vib excit. and attachment models are turned off in option 4 (and 6)
616// therefore I have created option 4a (and 6a), which has these models turned on
617//
618// // vib excitation
619// mod = new G4DNASancheExcitationModel();
620// em_config->SetExtraEmModel("e-", "e-_G4DNAVibExcitation",
621// mod, reg, 0.0, elimvb);
622//
623// // attachment
624// mod = new G4DNAMeltonAttachmentModel();
625// em_config->SetExtraEmModel("e-", "e-_G4DNAAttachment",
626// mod, reg, 0.0, elimat);
627
628}

References emax, G4LossTableManager::EmConfigurator(), eV, G4LossTableManager::Instance(), G4EmParameters::MaxKinEnergy(), MeV, G4INCL::Math::min(), reg, G4VEmModel::SetActivationLowEnergyLimit(), G4EmConfigurator::SetExtraEmModel(), and theParameters.

Referenced by ConstructProcess().

◆ AddElectronModels4a()

void G4EmDNAPhysicsActivator::AddElectronModels4a ( const G4String region,
G4bool  emsc,
G4double  elimel 
)
private

Definition at line 632 of file G4EmDNAPhysicsActivator.cc.

635{
636 G4EmConfigurator *em_config =
638 G4VEmModel *mod;
639
640 static const G4double elowest = 10 * eV; // seems to be option dependent, so I moved it here - MJPietrzak
641 static const G4double elimin = 1 * MeV;
642 static const G4double elimvb = 100 * eV;
643 static const G4double elimat = 13 * eV;
644
645
646 // for e- 100 MeV is a limit between different msc models
648
649 if (emsc) {
650 G4UrbanMscModel *msc = new G4UrbanMscModel();
651 msc->SetActivationLowEnergyLimit(elimel);
652 G4double emaxmsc = std::min(100 * MeV, emax);
653 em_config->SetExtraEmModel("e-", "msc", msc, reg, 0.0, emaxmsc);
654 } else {
655 mod = new G4eCoulombScatteringModel();
656 mod->SetActivationLowEnergyLimit(elimel);
657 em_config->SetExtraEmModel("e-", "CoulombScat", mod, reg, 0.0, emax);
658 }
659
660 // cuts and solvation
662 em_config->SetExtraEmModel("e-", "e-_G4DNAElectronSolvation",
663 mod, reg, 0., elowest);
664
665 // elastic
667 em_config->SetExtraEmModel("e-", "e-_G4DNAElastic",
668 mod, reg, 0.0, elimel);
669
670 // ionisation
671 mod = new G4MollerBhabhaModel();
672 mod->SetActivationLowEnergyLimit(elimin);
673 em_config->SetExtraEmModel("e-", "eIoni",
674 mod, reg, 0.0, emax,
676
678 em_config->SetExtraEmModel("e-", "e-_G4DNAIonisation",
679 mod, reg, elowest, elimin);
680
681 // excitation
683 em_config->SetExtraEmModel("e-", "e-_G4DNAExcitation",
684 mod, reg, 0., elimin);
685
686
687 // I don't understand why vib excit. and attachment models are turned off in option 4
688 // therefore I have created option 4a, which has these models turned on
689 // and here it is - MJPietrzak
690
691 // vib excitation
692 mod = new G4DNASancheExcitationModel();
693 em_config->SetExtraEmModel("e-", "e-_G4DNAVibExcitation",
694 mod, reg, 0.0, elimvb);
695
696 // attachment
697 mod = new G4DNAMeltonAttachmentModel();
698 em_config->SetExtraEmModel("e-", "e-_G4DNAAttachment",
699 mod, reg, 0.0, elimat);
700
701}

References emax, G4LossTableManager::EmConfigurator(), eV, G4LossTableManager::Instance(), G4EmParameters::MaxKinEnergy(), MeV, G4INCL::Math::min(), reg, G4VEmModel::SetActivationLowEnergyLimit(), G4EmConfigurator::SetExtraEmModel(), and theParameters.

Referenced by ConstructProcess().

◆ AddElectronModels6()

void G4EmDNAPhysicsActivator::AddElectronModels6 ( const G4String region,
G4bool  emsc,
G4double  elimel 
)
private

Definition at line 705 of file G4EmDNAPhysicsActivator.cc.

708{
709 G4EmConfigurator *em_config =
711 G4VEmModel *mod;
712
713 static const G4double elowest = 11 * eV; // seems to be option dependent, so I moved it here - MJPietrzak
714 static const G4double elimin = 1 * MeV;
715 // static const G4double elimvb = 100 * eV;
716 // static const G4double elimat = 13 * eV;
717
718 // for e- 100 MeV is a limit between different msc models
720
721 if (emsc) {
722 G4UrbanMscModel *msc = new G4UrbanMscModel();
723 msc->SetActivationLowEnergyLimit(elimel);
724 G4double emaxmsc = std::min(100 * MeV, emax);
725 em_config->SetExtraEmModel("e-", "msc", msc, reg, 0.0, emaxmsc);
726 } else {
727 mod = new G4eCoulombScatteringModel();
728 mod->SetActivationLowEnergyLimit(elimel);
729 em_config->SetExtraEmModel("e-", "CoulombScat", mod, reg, 0.0, emax);
730 }
731
732 // cuts and solvation
734 em_config->SetExtraEmModel("e-", "e-_G4DNAElectronSolvation",
735 mod, reg, 0., elowest);
736
737 // elastic
738 mod = new G4DNACPA100ElasticModel();
739 em_config->SetExtraEmModel("e-", "e-_G4DNAElastic",
740 mod, reg, 0.0, elimel);
741
742 // ionisation
743 mod = new G4MollerBhabhaModel();
744 mod->SetActivationLowEnergyLimit(elimin);
745 em_config->SetExtraEmModel("e-", "eIoni",
746 mod, reg, 0.0, emax,
748
749 mod = new G4DNACPA100IonisationModel();
750 em_config->SetExtraEmModel("e-", "e-_G4DNAIonisation",
751 mod, reg, elowest, elimin);
752
753 // excitation
754 mod = new G4DNACPA100ExcitationModel();
755 em_config->SetExtraEmModel("e-", "e-_G4DNAExcitation",
756 mod, reg, 0., elimin);
757
758
759// I don't understand why vib excit. and attachment models are turned off in option 6 (and 4)
760// therefore I have created option 6a (and 4a), which has these models turned on
761// and here it is - MJPietrzak
762
763// // vib excitation
764// mod = new G4DNASancheExcitationModel();
765// em_config->SetExtraEmModel("e-", "e-_G4DNAVibExcitation",
766// mod, reg, 0.0, elimvb);
767//
768// // attachment
769// mod = new G4DNAMeltonAttachmentModel();
770// em_config->SetExtraEmModel("e-", "e-_G4DNAAttachment",
771// mod, reg, 0.0, elimat);
772
773}

References emax, G4LossTableManager::EmConfigurator(), eV, G4LossTableManager::Instance(), G4EmParameters::MaxKinEnergy(), MeV, G4INCL::Math::min(), reg, G4VEmModel::SetActivationLowEnergyLimit(), G4EmConfigurator::SetExtraEmModel(), and theParameters.

Referenced by ConstructProcess().

◆ AddElectronModels6a()

void G4EmDNAPhysicsActivator::AddElectronModels6a ( const G4String region,
G4bool  emsc,
G4double  elimel 
)
private

Definition at line 776 of file G4EmDNAPhysicsActivator.cc.

779{
780 G4EmConfigurator *em_config =
782 G4VEmModel *mod;
783
784 static const G4double elowest = 11 * eV; // seems to be option dependent, so I moved it here - MJPietrzak
785 static const G4double elimin = 1 * MeV;
786 static const G4double elimvb = 100 * eV;
787 static const G4double elimat = 13 * eV;
788
789
790 // for e- 100 MeV is a limit between different msc models
792
793 if (emsc) {
794 G4UrbanMscModel *msc = new G4UrbanMscModel();
795 msc->SetActivationLowEnergyLimit(elimel);
796 G4double emaxmsc = std::min(100 * MeV, emax);
797 em_config->SetExtraEmModel("e-", "msc", msc, reg, 0.0, emaxmsc);
798 } else {
799 mod = new G4eCoulombScatteringModel();
800 mod->SetActivationLowEnergyLimit(elimel);
801 em_config->SetExtraEmModel("e-", "CoulombScat", mod, reg, 0.0, emax);
802 }
803
804 // cuts and solvation
806 em_config->SetExtraEmModel("e-", "e-_G4DNAElectronSolvation",
807 mod, reg, 0., elowest);
808
809 // elastic
810 mod = new G4DNACPA100ElasticModel();
811 em_config->SetExtraEmModel("e-", "e-_G4DNAElastic",
812 mod, reg, 0.0, elimel);
813
814 // ionisation
815 mod = new G4MollerBhabhaModel();
816 mod->SetActivationLowEnergyLimit(elimin);
817 em_config->SetExtraEmModel("e-", "eIoni",
818 mod, reg, 0.0, emax,
820
821 mod = new G4DNACPA100IonisationModel();
822 em_config->SetExtraEmModel("e-", "e-_G4DNAIonisation",
823 mod, reg, elowest, elimin);
824
825 // excitation
826 mod = new G4DNACPA100ExcitationModel();
827 em_config->SetExtraEmModel("e-", "e-_G4DNAExcitation",
828 mod, reg, 0., elimin);
829
830 // I don't understand why vib excit. and attachment models are turned off in option 6 (and 4)
831 // therefore I have created option 6a (and 4a), which has these models turned on
832 // and here it is - MJPietrzak
833
834 // vib excitation
835 mod = new G4DNASancheExcitationModel();
836 em_config->SetExtraEmModel("e-", "e-_G4DNAVibExcitation",
837 mod, reg, 0.0, elimvb);
838
839 // attachment
840 mod = new G4DNAMeltonAttachmentModel();
841 em_config->SetExtraEmModel("e-", "e-_G4DNAAttachment",
842 mod, reg, 0.0, elimat);
843
844}

References emax, G4LossTableManager::EmConfigurator(), eV, G4LossTableManager::Instance(), G4EmParameters::MaxKinEnergy(), MeV, G4INCL::Math::min(), reg, G4VEmModel::SetActivationLowEnergyLimit(), G4EmConfigurator::SetExtraEmModel(), and theParameters.

Referenced by ConstructProcess().

◆ AddElectronModels7()

void G4EmDNAPhysicsActivator::AddElectronModels7 ( const G4String region,
G4bool  emsc,
G4double  elimel 
)
private

Definition at line 848 of file G4EmDNAPhysicsActivator.cc.

851{
852 G4EmConfigurator *em_config =
854 G4VEmModel *mod;
855
856 static const G4double elowest = 10 * eV; // seems to be option dependent, so I moved it here - MJPietrzak
857 static const G4double elimin = 1 * MeV;
858 static const G4double elimvb = 100 * eV;
859 static const G4double elimat = 13 * eV;
860 static const G4double elim1 = 10 * keV;
861
862 // for e- 100 MeV is a limit between different msc models
864
865 if (emsc) {
866 G4UrbanMscModel *msc = new G4UrbanMscModel();
867 msc->SetActivationLowEnergyLimit(elimel);
868 G4double emaxmsc = std::min(100 * MeV, emax);
869 em_config->SetExtraEmModel("e-", "msc", msc, reg, 0.0, emaxmsc);
870 } else {
871 mod = new G4eCoulombScatteringModel();
872 mod->SetActivationLowEnergyLimit(elimel);
873 em_config->SetExtraEmModel("e-", "CoulombScat", mod, reg, 0.0, emax);
874 }
875
876 // cuts and solvation
878 em_config->SetExtraEmModel("e-", "e-_G4DNAElectronSolvation",
879 mod, reg, 0., elowest);
880
881 // elastic
882 mod = new G4DNAUeharaScreenedRutherfordElasticModel(); // G4DNAChampionElasticModel(); in Opt_0 - the main difference between Opt0 and Opt7
883 em_config->SetExtraEmModel("e-", "e-_G4DNAElastic",
884 mod, reg, 0.0, elimel);
885
886
887 // ionisation
888 mod = new G4MollerBhabhaModel();
889 mod->SetActivationLowEnergyLimit(elimin);
890 em_config->SetExtraEmModel("e-", "eIoni",
891 mod, reg, 0.0, emax,
893 // todo - MJPietrzak
894 // I don't understand why the MollerBhabhaModel is here, as there is no sign of this model in regular DNA lists.
895 // However, it seems that it is necessary for all the options.
896
897 mod = new G4DNABornIonisationModel();
898 em_config->SetExtraEmModel("e-", "e-_G4DNAIonisation",
899 mod, reg, elim1, elimin);
900
902 em_config->SetExtraEmModel("e-", "e-_G4DNAIonisation",
903 mod, reg, elowest, elim1); // mod, reg, 0.0, elim1);
904
905
906 // excitation
908 em_config->SetExtraEmModel("e-", "e-_G4DNAExcitation",
909 mod, reg, 8 * eV, elim1);
910
911 mod = new G4DNABornExcitationModel();
912 em_config->SetExtraEmModel("e-", "e-_G4DNAExcitation",
913 mod, reg, elim1, elimin);
914
915 // vib excitation
916 mod = new G4DNASancheExcitationModel();
917 em_config->SetExtraEmModel("e-", "e-_G4DNAVibExcitation",
918 mod, reg, 0.0, elimvb);
919
920 // attachment
921 mod = new G4DNAMeltonAttachmentModel();
922 em_config->SetExtraEmModel("e-", "e-_G4DNAAttachment",
923 mod, reg, 0.0, elimat);
924
925}

References emax, G4LossTableManager::EmConfigurator(), eV, G4DNABornIonisationModel, G4LossTableManager::Instance(), keV, G4EmParameters::MaxKinEnergy(), MeV, G4INCL::Math::min(), reg, G4VEmModel::SetActivationLowEnergyLimit(), G4EmConfigurator::SetExtraEmModel(), and theParameters.

Referenced by ConstructProcess().

◆ AddGenericIonModels0()

void G4EmDNAPhysicsActivator::AddGenericIonModels0 ( const G4String region,
G4double  pminbb 
)
private

Definition at line 1015 of file G4EmDNAPhysicsActivator.cc.

1017{
1018 G4EmConfigurator* em_config =
1020 G4VEmModel* mod;
1021
1023 G4double iemax = std::min(10*MeV, emax);
1024 //G4double iemin = 100*eV;
1025
1026 mod = new G4BraggIonModel();
1027 mod->SetActivationLowEnergyLimit(iemax);
1028 em_config->SetExtraEmModel("GenericIon", "ionIoni",
1029 mod, reg, 0.0, pminbb,
1030 new G4IonFluctuations());
1031
1032 mod = new G4BetheBlochModel();
1033 mod->SetActivationLowEnergyLimit(iemax);
1034 em_config->SetExtraEmModel("GenericIon", "ionIoni",
1035 mod, reg, pminbb, emax,
1036 new G4IonFluctuations());
1037
1039 em_config->SetExtraEmModel("GenericIon", "GenericIon_G4DNAIonisation",
1040 mod, reg, 0.0, iemax);
1041
1042}

References emax, G4LossTableManager::EmConfigurator(), G4LossTableManager::Instance(), G4EmParameters::MaxKinEnergy(), MeV, G4INCL::Math::min(), reg, G4VEmModel::SetActivationLowEnergyLimit(), G4EmConfigurator::SetExtraEmModel(), and theParameters.

Referenced by ConstructProcess().

◆ AddHeliumModels0()

void G4EmDNAPhysicsActivator::AddHeliumModels0 ( const G4String region,
G4bool  a1msc,
G4bool  a2msc,
G4double  elimel,
G4double  pminbb,
G4double  pmax 
)
private

Definition at line 1046 of file G4EmDNAPhysicsActivator.cc.

1050{
1051 G4EmConfigurator* em_config =
1053 G4VEmModel* mod;
1054
1055 static const G4double hemax = 400 * MeV;
1056 static const G4double massRatio = G4Alpha::Alpha()->GetPDGMass()/CLHEP::proton_mass_c2;
1057
1059 G4double pminbba = massRatio*pminbb;
1060 if(IsVerbose()) {
1061 G4cout << "AddHeliumModels0 for <" << reg << "> a1msc: " << a1msc <<" a2msc: " << a2msc
1062 << " elimel= " << elimel << " pminbba= " << pminbba << G4endl;
1063 }
1064 // alpha++
1065 if(elimel < emax) {
1066 if(a2msc) {
1067 G4UrbanMscModel* msc = new G4UrbanMscModel();
1068 msc->SetActivationLowEnergyLimit(elimel);
1069 em_config->SetExtraEmModel("alpha", "msc", msc, reg, 0.0, emax);
1070 } else {
1071 mod = new G4IonCoulombScatteringModel();
1072 mod->SetActivationLowEnergyLimit(elimel);
1073 em_config->SetExtraEmModel("alpha", "CoulombScat", mod, reg, 0.0, emax);
1074 }
1075 }
1076
1077 mod = new G4BraggIonModel();
1078 mod->SetActivationLowEnergyLimit(hemax/massRatio);
1079 em_config->SetExtraEmModel("alpha", "ionIoni",
1080 mod, reg, 0.0, pminbba,
1081 new G4IonFluctuations());
1082
1083 mod = new G4BetheBlochModel();
1084 mod->SetActivationLowEnergyLimit(hemax/massRatio);
1085 em_config->SetExtraEmModel("alpha", "ionIoni",
1086 mod, reg, pminbba, emax,
1087 new G4IonFluctuations());
1088
1089 mod = new G4DNARuddIonisationModel();
1090 em_config->SetExtraEmModel("alpha", "alpha_G4DNAIonisation",
1091 mod, reg, 0.0, hemax);
1092
1094 em_config->SetExtraEmModel("alpha", "alpha_G4DNAExcitation",
1095 mod, reg, 0.0, hemax);
1096
1098 em_config->SetExtraEmModel("alpha", "alpha_G4DNAChargeDecrease",
1099 mod, reg, 0.0, hemax);
1100
1101 mod = new G4DNAIonElasticModel();
1102 em_config->SetExtraEmModel("alpha", "alpha_G4DNAElastic",
1103 mod, reg, 0.0, elimel);
1104
1105 // ---
1106 // alpha+
1107 if(elimel < emax) {
1108 if(a1msc) {
1109 G4UrbanMscModel* msc = new G4UrbanMscModel();
1110 msc->SetActivationLowEnergyLimit(elimel);
1111 em_config->SetExtraEmModel("alpha+", "msc", msc, reg, 0.0, emax);
1112 } else {
1113 mod = new G4IonCoulombScatteringModel();
1114 mod->SetActivationLowEnergyLimit(elimel);
1115 em_config->SetExtraEmModel("alpha+", "CoulombScat", mod, reg, 0.0, emax);
1116 }
1117 }
1118
1119 mod = new G4BraggIonModel();
1120 mod->SetActivationLowEnergyLimit(hemax/massRatio);
1121 em_config->SetExtraEmModel("alpha+", "hIoni",
1122 mod, reg, 0.0, pminbba,
1123 new G4IonFluctuations());
1124
1125 mod = new G4BetheBlochModel();
1126 mod->SetActivationLowEnergyLimit(hemax/massRatio);
1127 em_config->SetExtraEmModel("alpha+", "hIoni",
1128 mod, reg, pminbba, emax,
1129 new G4IonFluctuations());
1130
1131 mod = new G4DNARuddIonisationModel();
1132 em_config->SetExtraEmModel("alpha+", "alpha+_G4DNAIonisation",
1133 mod, reg, 0.0, hemax);
1134
1136 em_config->SetExtraEmModel("alpha+", "alpha+_G4DNAExcitation",
1137 mod, reg, 0.0, hemax);
1138
1140 em_config->SetExtraEmModel("alpha+", "alpha+_G4DNAChargeDecrease",
1141 mod, reg, 0.0, hemax);
1142
1144 em_config->SetExtraEmModel("alpha+", "alpha+_G4DNAChargeIncrease",
1145 mod, reg, 0.0, hemax);
1146
1147 mod = new G4DNAIonElasticModel();
1148 em_config->SetExtraEmModel("alpha+", "alpha+_G4DNAElastic",
1149 mod, reg, 0.0, elimel);
1150
1151 // ---
1152 // helium
1153 mod = new G4DNARuddIonisationModel();
1154 em_config->SetExtraEmModel("helium", "helium_G4DNAIonisation",
1155 mod, reg, 0.0, hemax);
1156
1158 em_config->SetExtraEmModel("helium", "helium_G4DNAExcitation",
1159 mod, reg, 0.0, hemax);
1160
1162 em_config->SetExtraEmModel("helium", "helium_G4DNAChargeIncrease",
1163 mod, reg, 0.0, hemax);
1164
1165 mod = new G4DNAIonElasticModel();
1166 em_config->SetExtraEmModel("helium", "helium_G4DNAElastic",
1167 mod, reg, 0.0, elimel);
1168}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
static constexpr double proton_mass_c2

References G4Alpha::Alpha(), emax, G4LossTableManager::EmConfigurator(), G4cout, G4endl, G4ParticleDefinition::GetPDGMass(), G4LossTableManager::Instance(), IsVerbose(), G4EmParameters::MaxKinEnergy(), MeV, CLHEP::proton_mass_c2, reg, G4VEmModel::SetActivationLowEnergyLimit(), G4EmConfigurator::SetExtraEmModel(), and theParameters.

Referenced by ConstructProcess().

◆ AddProtonModels0()

void G4EmDNAPhysicsActivator::AddProtonModels0 ( const G4String region,
G4bool  pmsc,
G4double  elimel,
G4double  pminbb,
G4double  pmax 
)
private

Definition at line 933 of file G4EmDNAPhysicsActivator.cc.

936{
937 G4EmConfigurator* em_config =
939 G4VEmModel* mod;
940
941 static const G4double gmmax = 500 * keV;
942
944
945 // proton
946
947 // if SS physics list msc process does not exist
948 if(pmsc) {
950 msc->SetActivationLowEnergyLimit(elimel);
951 em_config->SetExtraEmModel("proton", "msc", msc, reg, 0.0, emax);
952 }
953 // single scattering always applied
954 mod = new G4eCoulombScatteringModel();
955 mod->SetActivationLowEnergyLimit(elimel);
956 em_config->SetExtraEmModel("proton", "CoulombScat", mod, reg, 0.0, emax);
957
958 mod = new G4BraggModel();
959 mod->SetActivationLowEnergyLimit(std::min(pminbb, pmax));
960 em_config->SetExtraEmModel("proton", "hIoni",
961 mod, reg, 0.0, pminbb,
963
964 mod = new G4BetheBlochModel();
966 em_config->SetExtraEmModel("proton", "hIoni",
967 mod, reg, pminbb, emax,
969
970 mod = new G4DNARuddIonisationModel();
971 em_config->SetExtraEmModel("proton", "proton_G4DNAIonisation",
972 mod, reg, 0.0, gmmax);
973
974 mod = new G4DNABornIonisationModel();
975 em_config->SetExtraEmModel("proton", "proton_G4DNAIonisation",
976 mod, reg, gmmax, pmax);
977
979 em_config->SetExtraEmModel("proton", "proton_G4DNAExcitation",
980 mod, reg, 0.0, gmmax);
981
982 mod = new G4DNABornExcitationModel();
983 em_config->SetExtraEmModel("proton", "proton_G4DNAExcitation",
984 mod, reg, gmmax, pmax);
985
987 em_config->SetExtraEmModel("proton", "proton_G4DNAChargeDecrease",
988 mod, reg, 0.0, pmax);
989
990 mod = new G4DNAIonElasticModel();
991 em_config->SetExtraEmModel("proton", "proton_G4DNAElastic",
992 mod, reg, 0.0, elimel);
993
994 // hydrogen
995 mod = new G4DNARuddIonisationModel();
996 em_config->SetExtraEmModel("hydrogen", "hydrogen_G4DNAIonisation",
997 mod, reg, 0.0, pmax);
998
1000 em_config->SetExtraEmModel("hydrogen", "hydrogen_G4DNAExcitation",
1001 mod, reg, 0.0, gmmax);
1002
1004 em_config->SetExtraEmModel("hydrogen", "hydrogen_G4DNAChargeIncrease",
1005 mod, reg, 0.0, pmax);
1006
1007 mod = new G4DNAIonElasticModel();
1008 em_config->SetExtraEmModel("hydrogen", "hydrogen_G4DNAElastic",
1009 mod, reg, 0.0, elimel);
1010
1011}

References emax, G4LossTableManager::EmConfigurator(), G4DNABornIonisationModel, G4LossTableManager::Instance(), keV, G4EmParameters::MaxKinEnergy(), G4INCL::Math::min(), reg, G4VEmModel::SetActivationLowEnergyLimit(), G4EmConfigurator::SetExtraEmModel(), and theParameters.

Referenced by ConstructProcess().

◆ ConstructParticle()

void G4EmDNAPhysicsActivator::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 120 of file G4EmDNAPhysicsActivator.cc.

121{
122 // bosons
124
125 // leptons
128
129 // baryons
131
134
135 G4DNAGenericIonsManager * genericIonsManager;
136 genericIonsManager=G4DNAGenericIonsManager::Instance();
137 genericIonsManager->GetIon("alpha+");
138 genericIonsManager->GetIon("helium");
139 genericIonsManager->GetIon("hydrogen");
140 //genericIonsManager->GetIon("carbon");
141 //genericIonsManager->GetIon("nitrogen");
142 //genericIonsManager->GetIon("oxygen");
143 //genericIonsManager->GetIon("iron");
144
145}
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 G4Alpha::Alpha(), G4Electron::Electron(), G4Gamma::Gamma(), G4GenericIon::GenericIonDefinition(), G4DNAGenericIonsManager::GetIon(), G4DNAGenericIonsManager::Instance(), G4Positron::Positron(), and G4Proton::Proton().

◆ ConstructProcess()

void G4EmDNAPhysicsActivator::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 149 of file G4EmDNAPhysicsActivator.cc.

150{
151 const std::vector<G4String>& regnamesDNA = theParameters->RegionsDNA();
152 G4int nreg = regnamesDNA.size();
153 if(0 == nreg)
154 {
155 return;
156 }
157 const std::vector<G4String>& typesDNA = theParameters->TypesDNA();
158
159 if(IsVerbose()) {
160 G4cout << "### G4EmDNAPhysicsActivator::ConstructProcess for " << nreg
161 << " regions; DNA physics type " << typesDNA[0] << G4endl;
162 }
163
164 // list of particles
168
169 G4DNAGenericIonsManager * genericIonsManager =
172 const G4ParticleDefinition* alpha1 = genericIonsManager->GetIon("alpha+");
173 const G4ParticleDefinition* alpha0 = genericIonsManager->GetIon("helium");
174 const G4ParticleDefinition* h0 = genericIonsManager->GetIon("hydrogen");
175
176 G4ProcessManager* eman = elec->GetProcessManager();
177 G4ProcessManager* pman = prot->GetProcessManager();
178 G4ProcessManager* iman = gion->GetProcessManager();
179 G4ProcessManager* a2man = alpha2->GetProcessManager();
180 G4ProcessManager* a1man = alpha1->GetProcessManager();
181 G4ProcessManager* a0man = alpha0->GetProcessManager();
182 G4ProcessManager* h0man = h0->GetProcessManager();
183
184 // alpha+ standard processes
186 G4ParticleDefinition* alpha11 = const_cast<G4ParticleDefinition*>(alpha1);
187 ph->RegisterProcess(new G4hMultipleScattering(), alpha11);
188 ph->RegisterProcess(new G4hIonisation(), alpha11);
189
190 G4bool emsc = HasMsc(eman);
191 G4bool pmsc = HasMsc(pman);
192 G4bool a2msc = HasMsc(a2man);
193 G4bool a1msc = HasMsc(a1man);
194 // G4bool imsc = HasMsc(iman);
195
196 // processes are defined with dummy models for the world
197 // elastic scatetring
198 G4DNAElastic* theDNAeElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
199 theDNAeElasticProcess->SetEmModel(new G4DummyModel());
200 eman->AddDiscreteProcess(theDNAeElasticProcess);
201
202 G4DNAElastic* theDNApElasticProcess =
203 new G4DNAElastic("proton_G4DNAElastic");
204 theDNApElasticProcess->SetEmModel(new G4DummyModel());
205 pman->AddDiscreteProcess(theDNApElasticProcess);
206
207 G4DNAElastic* theDNAa2ElasticProcess =
208 new G4DNAElastic("alpha_G4DNAElastic");
209 theDNAa2ElasticProcess->SetEmModel(new G4DummyModel());
210 a2man->AddDiscreteProcess(theDNAa2ElasticProcess);
211
212 G4DNAElastic* theDNAa1ElasticProcess =
213 new G4DNAElastic("alpha+_G4DNAElastic");
214 theDNAa1ElasticProcess->SetEmModel(new G4DummyModel());
215 a1man->AddDiscreteProcess(theDNAa1ElasticProcess);
216
217 G4DNAElastic* theDNAa0ElasticProcess =
218 new G4DNAElastic("helium_G4DNAElastic");
219 theDNAa0ElasticProcess->SetEmModel(new G4DummyModel());
220 a0man->AddDiscreteProcess(theDNAa0ElasticProcess);
221
222 G4DNAElastic* theDNAh0ElasticProcess =
223 new G4DNAElastic("hydrogen_G4DNAElastic");
224 theDNAh0ElasticProcess->SetEmModel(new G4DummyModel());
225 h0man->AddDiscreteProcess(theDNAh0ElasticProcess);
226
227 // excitation
228 G4DNAExcitation* theDNAeExcProcess =
229 new G4DNAExcitation("e-_G4DNAExcitation");
230 theDNAeExcProcess->SetEmModel(new G4DummyModel());
231 eman->AddDiscreteProcess(theDNAeExcProcess);
232
233 G4DNAExcitation* theDNApExcProcess =
234 new G4DNAExcitation("proton_G4DNAExcitation");
235 theDNApExcProcess->SetEmModel(new G4DummyModel());
236 pman->AddDiscreteProcess(theDNApExcProcess);
237
238 G4DNAExcitation* theDNAa2ExcProcess =
239 new G4DNAExcitation("alpha_G4DNAExcitation");
240 theDNAa2ExcProcess->SetEmModel(new G4DummyModel());
241 a2man->AddDiscreteProcess(theDNAa2ExcProcess);
242
243 G4DNAExcitation* theDNAa1ExcProcess =
244 new G4DNAExcitation("alpha+_G4DNAExcitation");
245 theDNAa1ExcProcess->SetEmModel(new G4DummyModel());
246 a1man->AddDiscreteProcess(theDNAa1ExcProcess);
247
248 G4DNAExcitation* theDNAa0ExcProcess =
249 new G4DNAExcitation("helium_G4DNAExcitation");
250 theDNAa0ExcProcess->SetEmModel(new G4DummyModel());
251 a0man->AddDiscreteProcess(theDNAa0ExcProcess);
252
253 G4DNAExcitation* theDNAh0ExcProcess =
254 new G4DNAExcitation("hydrogen_G4DNAExcitation");
255 theDNAh0ExcProcess->SetEmModel(new G4DummyModel());
256 h0man->AddDiscreteProcess(theDNAh0ExcProcess);
257
258 // vibration excitation
259 G4DNAVibExcitation* theDNAeVibExcProcess =
260 new G4DNAVibExcitation("e-_G4DNAVibExcitation");
261 theDNAeVibExcProcess->SetEmModel(new G4DummyModel());
262 eman->AddDiscreteProcess(theDNAeVibExcProcess);
263
264 // ionisation
265 G4DNAIonisation* theDNAeIoniProcess =
266 new G4DNAIonisation("e-_G4DNAIonisation");
267 theDNAeIoniProcess->SetEmModel(new G4DummyModel());
268 eman->AddDiscreteProcess(theDNAeIoniProcess);
269
270 G4DNAIonisation* theDNApIoniProcess =
271 new G4DNAIonisation("proton_G4DNAIonisation");
272 theDNApIoniProcess->SetEmModel(new G4DummyModel());
273 pman->AddDiscreteProcess(theDNApIoniProcess);
274
275 G4DNAIonisation* theDNAa2IoniProcess =
276 new G4DNAIonisation("alpha_G4DNAIonisation");
277 theDNAa2IoniProcess->SetEmModel(new G4DummyModel());
278 a2man->AddDiscreteProcess(theDNAa2IoniProcess);
279
280 G4DNAIonisation* theDNAa1IoniProcess =
281 new G4DNAIonisation("alpha+_G4DNAIonisation");
282 theDNAa1IoniProcess->SetEmModel(new G4DummyModel());
283 a1man->AddDiscreteProcess(theDNAa1IoniProcess);
284
285 G4DNAIonisation* theDNAa0IoniProcess =
286 new G4DNAIonisation("helium_G4DNAIonisation");
287 theDNAa0IoniProcess->SetEmModel(new G4DummyModel());
288 a0man->AddDiscreteProcess(theDNAa0IoniProcess);
289
290 G4DNAIonisation* theDNAh0IoniProcess =
291 new G4DNAIonisation("hydrogen_G4DNAIonisation");
292 theDNAh0IoniProcess->SetEmModel(new G4DummyModel());
293 h0man->AddDiscreteProcess(theDNAh0IoniProcess);
294
295 G4DNAIonisation* theDNAiIoniProcess =
296 new G4DNAIonisation("GenericIon_G4DNAIonisation");
297 theDNAiIoniProcess->SetEmModel(new G4DummyModel());
298 iman->AddDiscreteProcess(theDNAiIoniProcess);
299
300 // attachment
301 G4DNAAttachment* theDNAAttachProcess =
302 new G4DNAAttachment("e-_G4DNAAttachment");
303 theDNAAttachProcess->SetEmModel(new G4DummyModel());
304 eman->AddDiscreteProcess(theDNAAttachProcess);
305
306 // charge exchange
307 G4DNAChargeDecrease* theDNApChargeDecreaseProcess =
308 new G4DNAChargeDecrease("proton_G4DNAChargeDecrease");
309 theDNApChargeDecreaseProcess->SetEmModel(new G4DummyModel());
310 pman->AddDiscreteProcess(theDNApChargeDecreaseProcess);
311
312 G4DNAChargeDecrease* theDNAa2ChargeDecreaseProcess =
313 new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease");
314 theDNAa2ChargeDecreaseProcess->SetEmModel(new G4DummyModel());
315 a2man->AddDiscreteProcess(theDNAa2ChargeDecreaseProcess);
316
317 G4DNAChargeDecrease* theDNAa1ChargeDecreaseProcess =
318 new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease");
319 theDNAa1ChargeDecreaseProcess->SetEmModel(new G4DummyModel());
320 a1man->AddDiscreteProcess(theDNAa1ChargeDecreaseProcess);
321
322 G4DNAChargeIncrease* theDNAa1ChargeIncreaseProcess =
323 new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease");
324 theDNAa1ChargeIncreaseProcess->SetEmModel(new G4DummyModel());
325 a1man->AddDiscreteProcess(theDNAa1ChargeIncreaseProcess);
326
327 G4DNAChargeIncrease* theDNAa0ChargeIncreaseProcess =
328 new G4DNAChargeIncrease("helium_G4DNAChargeIncrease");
329 theDNAa0ChargeIncreaseProcess->SetEmModel(new G4DummyModel());
330 a0man->AddDiscreteProcess(theDNAa0ChargeIncreaseProcess);
331
332 G4DNAChargeIncrease* theDNAh0ChargeIncreaseProcess =
333 new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease");
334 theDNAh0ChargeIncreaseProcess->SetEmModel(new G4DummyModel());
335 h0man->AddDiscreteProcess(theDNAh0ChargeIncreaseProcess);
336
337 // limits for DNA model applicability
338 // static const G4double elowest= 7.4 * eV; // seems to be option dependent - MJPietrzak
339 static const G4double elimel = 1 * MeV;
340 static const G4double pminbb = 2 * MeV;
341 static const G4double pmin = 0.1 * keV;
342 static const G4double pmax = 100 * MeV;
343 static const G4double hemin = 1 * keV;
344 static const G4double ionmin = 0.5 * MeV;
345
346
347 // G4DNAElectronSolvation used instead of G4LowECapture (always for simplicity) - MJPietrzak
348
349 // When chemistry is activated: G4DNAElectronSolvation turns the electron
350 // to a solvated electron, otherwise it kills the electron at the
351 // corresponding high energy limit of the model
352 auto pSolvatation = new G4DNAElectronSolvation("e-_G4DNAElectronSolvation");
353 pSolvatation->SetEmModel(G4DNASolvationModelFactory::GetMacroDefinedModel());
354 eman->AddDiscreteProcess(pSolvatation);
355
356 G4LowECapture* pcap = new G4LowECapture(pmin);
357 pman->AddDiscreteProcess(pcap);
358 G4LowECapture* icap = new G4LowECapture(ionmin);
359 iman->AddDiscreteProcess(icap);
360 G4LowECapture* a2cap = new G4LowECapture(hemin);
361 a2man->AddDiscreteProcess(a2cap);
362 G4LowECapture* a1cap = new G4LowECapture(hemin);
363 a1man->AddDiscreteProcess(a1cap);
364 G4LowECapture* a0cap = new G4LowECapture(hemin);
365 a0man->AddDiscreteProcess(a0cap);
366 G4LowECapture* h0cap = new G4LowECapture(ionmin);
367 h0man->AddDiscreteProcess(h0cap);
368
369 // loop over regions
370 for(G4int i = 0; i < nreg; ++i)
371 {
372 G4String reg = regnamesDNA[i];
373 if(IsVerbose())
374 {
375 G4cout << "### DNA models type " << typesDNA[i]
376 << " are activated for G4Region " << reg << G4endl;
377 }
378
379 if(typesDNA[i] == "DNA_Opt0") {
380 AddElectronModels0(reg, emsc, elimel);
381 }
382 else if(typesDNA[i] == "DNA_Opt2"){
383 AddElectronModels2(reg, emsc, elimel);
384 }
385 else if(typesDNA[i] == "DNA_Opt4"){
386 AddElectronModels4(reg, emsc, elimel);
387 }
388 else if(typesDNA[i] == "DNA_Opt4a"){
389 AddElectronModels4a(reg, emsc, elimel);
390 }
391 else if(typesDNA[i] == "DNA_Opt6"){
392 AddElectronModels6(reg, emsc, elimel);
393 }
394 else if(typesDNA[i] == "DNA_Opt6a"){
395 AddElectronModels6a(reg, emsc, elimel);
396 }
397 else if(typesDNA[i] == "DNA_Opt7") {
398 AddElectronModels7(reg, emsc, elimel);
399 }
400
401 // models for for other particles seems to be option independent so I moved them here - MJPietrzak
402 AddProtonModels0(reg, pmsc, elimel, pminbb, pmax);
403 AddHeliumModels0(reg, a1msc, a2msc, elimel, pminbb, pmax);
404 AddGenericIonModels0(reg, pminbb);
405 DeactivateNuclearStopping(pman, elimel);
406 DeactivateNuclearStopping(a1man, elimel);
407 DeactivateNuclearStopping(a2man, elimel);
408 }
409
410 auto ltman = G4LossTableManager::Instance();
411 ltman->EmConfigurator()->AddModels();
412}
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4double alpha2
static G4VEmModel * GetMacroDefinedModel()
One step thermalization model can be chosen via macro using /process/dna/e-SolvationSubType Ritchie19...
void AddElectronModels2(const G4String &region, G4bool emsc, G4double elimel)
G4bool HasMsc(G4ProcessManager *) const
void AddHeliumModels0(const G4String &region, G4bool a1msc, G4bool a2msc, G4double elimel, G4double pminbb, G4double pmax)
void AddProtonModels0(const G4String &region, G4bool pmsc, G4double elimel, G4double pminbb, G4double pmax)
void AddElectronModels4a(const G4String &region, G4bool emsc, G4double elimel)
void AddElectronModels4(const G4String &region, G4bool emsc, G4double elimel)
void AddElectronModels7(const G4String &region, G4bool emsc, G4double elimel)
void DeactivateNuclearStopping(G4ProcessManager *, G4double elimel)
void AddGenericIonModels0(const G4String &region, G4double pminbb)
void AddElectronModels6a(const G4String &region, G4bool emsc, G4double elimel)
void AddElectronModels0(const G4String &region, G4bool emsc, G4double elimel)
void AddElectronModels6(const G4String &region, G4bool emsc, G4double elimel)
const std::vector< G4String > & TypesDNA() const
const std::vector< G4String > & RegionsDNA() const
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
G4ProcessManager * GetProcessManager() const
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
void SetEmModel(G4VEmModel *, G4int index=0)

References G4ProcessManager::AddDiscreteProcess(), AddElectronModels0(), AddElectronModels2(), AddElectronModels4(), AddElectronModels4a(), AddElectronModels6(), AddElectronModels6a(), AddElectronModels7(), AddGenericIonModels0(), AddHeliumModels0(), AddProtonModels0(), G4Alpha::Alpha(), alpha2, DeactivateNuclearStopping(), G4Electron::Electron(), G4cout, G4endl, G4GenericIon::GenericIon(), G4DNAGenericIonsManager::GetIon(), G4DNASolvationModelFactory::GetMacroDefinedModel(), G4PhysicsListHelper::GetPhysicsListHelper(), G4ParticleDefinition::GetProcessManager(), HasMsc(), G4LossTableManager::Instance(), G4DNAGenericIonsManager::Instance(), IsVerbose(), keV, MeV, G4Proton::Proton(), reg, G4EmParameters::RegionsDNA(), G4PhysicsListHelper::RegisterProcess(), G4VEmProcess::SetEmModel(), theParameters, and G4EmParameters::TypesDNA().

◆ DeactivateNuclearStopping()

void G4EmDNAPhysicsActivator::DeactivateNuclearStopping ( G4ProcessManager pman,
G4double  elimel 
)
private

Definition at line 1172 of file G4EmDNAPhysicsActivator.cc.

1174{
1175 G4ProcessVector* pv = pman->GetProcessList();
1176 G4int nproc = pman->GetProcessListLength();
1177 for(G4int i = 0; i < nproc; ++i) {
1178 if(((*pv)[i])->GetProcessSubType() == fNuclearStopping) {
1179 G4VEmProcess* proc = static_cast<G4VEmProcess*>((*pv)[i]);
1180 if(proc) {
1182 mod->SetActivationLowEnergyLimit(elimel);
1183 proc->SetEmModel(mod);
1184 }
1185 break;
1186 }
1187 }
1188}
@ fNuclearStopping
G4int GetProcessListLength() const
G4ProcessVector * GetProcessList() const

References fNuclearStopping, G4ProcessManager::GetProcessList(), G4ProcessManager::GetProcessListLength(), G4VEmModel::SetActivationLowEnergyLimit(), and G4VEmProcess::SetEmModel().

Referenced by ConstructProcess().

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

◆ HasMsc()

G4bool G4EmDNAPhysicsActivator::HasMsc ( G4ProcessManager pman) const
private

Definition at line 1192 of file G4EmDNAPhysicsActivator.cc.

1193{
1194 G4bool res = false;
1195 G4ProcessVector* pv = pman->GetProcessList();
1196 G4int nproc = pman->GetProcessListLength();
1197 for(G4int i = 0; i < nproc; ++i)
1198 {
1199 if(((*pv)[i])->GetProcessSubType() == fMultipleScattering)
1200 {
1201 res = true;
1202 break;
1203 }
1204 }
1205 return res;
1206}
@ fMultipleScattering

References fMultipleScattering, G4ProcessManager::GetProcessList(), and G4ProcessManager::GetProcessListLength().

Referenced by ConstructProcess().

◆ IsVerbose()

G4bool G4EmDNAPhysicsActivator::IsVerbose ( ) const
private

Definition at line 113 of file G4EmDNAPhysicsActivator.cc.

114{
115 return (0 < verbose && G4Threading::IsMasterThread());
116}
G4bool IsMasterThread()
Definition: G4Threading.cc:124

References G4Threading::IsMasterThread(), and verbose.

Referenced by AddHeliumModels0(), and ConstructProcess().

◆ 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::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

◆ theParameters

G4EmParameters* G4EmDNAPhysicsActivator::theParameters
private

◆ theParticleTable

G4ParticleTable* G4VPhysicsConstructor::theParticleTable = nullptr
protectedinherited

◆ typePhysics

G4int G4VPhysicsConstructor::typePhysics = 0
protectedinherited

◆ verbose

G4int G4EmDNAPhysicsActivator::verbose
private

Definition at line 74 of file G4EmDNAPhysicsActivator.hh.

Referenced by IsVerbose().

◆ verboseLevel

G4int G4VPhysicsConstructor::verboseLevel = 0
protectedinherited

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