Geant4  10.02.p03
G4EmModelActivator Class Reference

#include <G4EmModelActivator.hh>

Collaboration diagram for G4EmModelActivator:

Public Member Functions

 G4EmModelActivator ()
 
 ~G4EmModelActivator ()
 
void ConstructParticle ()
 
void ConstructProcess ()
 

Private Member Functions

void ConstructDNAParticles ()
 
void ActivatePAI ()
 
void ActivateMicroElec ()
 
void ActivateDNA ()
 
G4bool HasMsc (G4ProcessManager *) const
 
G4EmModelActivatoroperator= (const G4EmModelActivator &right)
 
 G4EmModelActivator (const G4EmModelActivator &)
 

Private Attributes

G4EmParameterstheParameters
 

Detailed Description

Definition at line 57 of file G4EmModelActivator.hh.

Constructor & Destructor Documentation

◆ G4EmModelActivator() [1/2]

G4EmModelActivator::G4EmModelActivator ( )

Definition at line 107 of file G4EmModelActivator.cc.

108 {
110 }
G4EmParameters * theParameters
static G4EmParameters * Instance()
Here is the call graph for this function:

◆ ~G4EmModelActivator()

G4EmModelActivator::~G4EmModelActivator ( )

Definition at line 114 of file G4EmModelActivator.cc.

115 {
116 }

◆ G4EmModelActivator() [2/2]

G4EmModelActivator::G4EmModelActivator ( const G4EmModelActivator )
private

Member Function Documentation

◆ ActivateDNA()

void G4EmModelActivator::ActivateDNA ( )
private

Definition at line 436 of file G4EmModelActivator.cc.

437 {
438  const std::vector<G4String>& regnamesDNA = theParameters->RegionsDNA();
439  G4int nreg = regnamesDNA.size();
440  if(0 == nreg)
441  {
442  return;
443  }
444  const std::vector<G4String> typesDNA = theParameters->TypesDNA();
445  G4int verbose = theParameters->Verbose() - 1;
446  if(verbose > 0)
447  {
448  G4cout << "### G4EmModelActivator::ActivateMicroElec for " << nreg
449  << " regions" << G4endl;
450  }
452 
453  // list of particles
455  const G4ParticleDefinition* prot = G4Proton::Proton();
457 
458  G4DNAGenericIonsManager * genericIonsManager =
460  const G4ParticleDefinition* alpha2 = G4Alpha::Alpha();
461  const G4ParticleDefinition* alpha1 = genericIonsManager->GetIon("alpha+");
462  const G4ParticleDefinition* alpha0 = genericIonsManager->GetIon("helium");
463  const G4ParticleDefinition* h0 = genericIonsManager->GetIon("hydrogen");
464 
465  G4ProcessManager* eman = elec->GetProcessManager();
466  G4ProcessManager* pman = prot->GetProcessManager();
467  G4ProcessManager* iman = gion->GetProcessManager();
468  G4ProcessManager* a2man = alpha2->GetProcessManager();
469  G4ProcessManager* a1man = alpha1->GetProcessManager();
470  G4ProcessManager* a0man = alpha0->GetProcessManager();
471  G4ProcessManager* h0man = h0->GetProcessManager();
472 
473  G4bool emsc = HasMsc(eman);
474  //G4bool pmsc = HasMsc(pman);
475  //G4bool a2msc = HasMsc(a2man);
476  //G4bool a1msc = HasMsc(a1man);
477 
478  // alpha+ standard processes
480  G4ParticleDefinition* alpha11 = const_cast<G4ParticleDefinition*>(alpha1);
481  ph->RegisterProcess(new G4hMultipleScattering(), alpha11);
482  ph->RegisterProcess(new G4ionIonisation(), alpha11);
483 
484  // processes are defined with dummy models for the world
485  // elastic scatetring
486  G4DNAElastic* theDNAeElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
487  theDNAeElasticProcess->AddEmModel(0, new G4DummyModel());
488  eman->AddDiscreteProcess(theDNAeElasticProcess);
489 
490  G4DNAElastic* theDNApElasticProcess = new G4DNAElastic("proton_G4DNAElastic");
491  theDNApElasticProcess->AddEmModel(0, new G4DummyModel());
492  pman->AddDiscreteProcess(theDNApElasticProcess);
493 
494  G4DNAElastic* theDNAa2ElasticProcess = new G4DNAElastic("alpha_G4DNAElastic");
495  theDNAa2ElasticProcess->AddEmModel(0, new G4DummyModel());
496  a2man->AddDiscreteProcess(theDNAa2ElasticProcess);
497 
498  G4DNAElastic* theDNAa1ElasticProcess =
499  new G4DNAElastic("alpha+_G4DNAElastic");
500  theDNAa1ElasticProcess->AddEmModel(0, new G4DummyModel());
501  a1man->AddDiscreteProcess(theDNAa1ElasticProcess);
502 
503  G4DNAElastic* theDNAa0ElasticProcess =
504  new G4DNAElastic("helium_G4DNAElastic");
505  theDNAa0ElasticProcess->AddEmModel(0, new G4DummyModel());
506  a0man->AddDiscreteProcess(theDNAa0ElasticProcess);
507 
508  G4DNAElastic* theDNAh0ElasticProcess =
509  new G4DNAElastic("hydrogen_G4DNAElastic");
510  theDNAh0ElasticProcess->AddEmModel(0, new G4DummyModel());
511  h0man->AddDiscreteProcess(theDNAh0ElasticProcess);
512 
513  // excitation
514  G4DNAExcitation* theDNAeExcProcess =
515  new G4DNAExcitation("e-_G4DNAExcitation");
516  theDNAeExcProcess->AddEmModel(0, new G4DummyModel());
517  eman->AddDiscreteProcess(theDNAeExcProcess);
518 
519  G4DNAExcitation* theDNApExcProcess =
520  new G4DNAExcitation("proton_G4DNAExcitation");
521  theDNApExcProcess->AddEmModel(0, new G4DummyModel());
522  pman->AddDiscreteProcess(theDNApExcProcess);
523 
524  G4DNAExcitation* theDNAa2ExcProcess =
525  new G4DNAExcitation("alpha_G4DNAExcitation");
526  theDNAa2ExcProcess->AddEmModel(0, new G4DummyModel());
527  a2man->AddDiscreteProcess(theDNAa2ExcProcess);
528 
529  G4DNAExcitation* theDNAa1ExcProcess =
530  new G4DNAExcitation("alpha+_G4DNAExcitation");
531  theDNAa1ExcProcess->AddEmModel(0, new G4DummyModel());
532  a1man->AddDiscreteProcess(theDNAa1ExcProcess);
533 
534  G4DNAExcitation* theDNAa0ExcProcess =
535  new G4DNAExcitation("helium_G4DNAExcitation");
536  theDNAa0ExcProcess->AddEmModel(0, new G4DummyModel());
537  a0man->AddDiscreteProcess(theDNAa0ExcProcess);
538 
539  G4DNAExcitation* theDNAh0ExcProcess =
540  new G4DNAExcitation("hydrogen_G4DNAExcitation");
541  theDNAh0ExcProcess->AddEmModel(0, new G4DummyModel());
542  h0man->AddDiscreteProcess(theDNAh0ExcProcess);
543 
544  // vibration excitation
545  G4DNAVibExcitation* theDNAeVibExcProcess =
546  new G4DNAVibExcitation("e-_G4DNAVibExcitation");
547  theDNAeVibExcProcess->AddEmModel(0, new G4DummyModel());
548  eman->AddDiscreteProcess(theDNAeVibExcProcess);
549 
550  // ionisation
551  G4DNAIonisation* theDNAeIoniProcess =
552  new G4DNAIonisation("e-_G4DNAIonisation");
553  theDNAeIoniProcess->AddEmModel(0, new G4DummyModel());
554  eman->AddDiscreteProcess(theDNAeIoniProcess);
555 
556  G4DNAIonisation* theDNApIoniProcess =
557  new G4DNAIonisation("proton_G4DNAIonisation");
558  theDNApIoniProcess->AddEmModel(0, new G4DummyModel());
559  pman->AddDiscreteProcess(theDNApIoniProcess);
560 
561  G4DNAIonisation* theDNAa2IoniProcess =
562  new G4DNAIonisation("alpha_G4DNAIonisation");
563  theDNAa2IoniProcess->AddEmModel(0, new G4DummyModel());
564  a2man->AddDiscreteProcess(theDNAa2IoniProcess);
565 
566  G4DNAIonisation* theDNAa1IoniProcess =
567  new G4DNAIonisation("alpha+_G4DNAIonisation");
568  theDNAa1IoniProcess->AddEmModel(0, new G4DummyModel());
569  a1man->AddDiscreteProcess(theDNAa1IoniProcess);
570 
571  G4DNAIonisation* theDNAa0IoniProcess =
572  new G4DNAIonisation("helium_G4DNAIonisation");
573  theDNAa0IoniProcess->AddEmModel(0, new G4DummyModel());
574  a0man->AddDiscreteProcess(theDNAa0IoniProcess);
575 
576  G4DNAIonisation* theDNAh0IoniProcess =
577  new G4DNAIonisation("hydrogen_G4DNAIonisation");
578  theDNAh0IoniProcess->AddEmModel(0, new G4DummyModel());
579  h0man->AddDiscreteProcess(theDNAh0IoniProcess);
580 
581  G4DNAIonisation* theDNAiIoniProcess =
582  new G4DNAIonisation("GenericIon_G4DNAIonisation");
583  theDNAiIoniProcess->AddEmModel(0, new G4DummyModel());
584  iman->AddDiscreteProcess(theDNAiIoniProcess);
585 
586  // attachment
587  G4DNAAttachment* theDNAAttachProcess =
588  new G4DNAAttachment("e-_G4DNAAttachment");
589  theDNAAttachProcess->AddEmModel(0, new G4DummyModel());
590  eman->AddDiscreteProcess(theDNAAttachProcess);
591 
592  // charge exchange
593  G4DNAChargeDecrease* theDNApChargeDecreaseProcess =
594  new G4DNAChargeDecrease("proton_G4DNAChargeDecrease");
595  theDNApChargeDecreaseProcess->AddEmModel(0, new G4DummyModel());
596  pman->AddDiscreteProcess(theDNApChargeDecreaseProcess);
597 
598  G4DNAChargeDecrease* theDNAa2ChargeDecreaseProcess =
599  new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease");
600  theDNAa2ChargeDecreaseProcess->AddEmModel(0, new G4DummyModel());
601  a2man->AddDiscreteProcess(theDNAa2ChargeDecreaseProcess);
602 
603  G4DNAChargeDecrease* theDNAa1ChargeDecreaseProcess =
604  new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease");
605  theDNAa1ChargeDecreaseProcess->AddEmModel(0, new G4DummyModel());
606  a1man->AddDiscreteProcess(theDNAa1ChargeDecreaseProcess);
607 
608  G4DNAChargeIncrease* theDNAa1ChargeIncreaseProcess =
609  new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease");
610  theDNAa1ChargeIncreaseProcess->AddEmModel(0, new G4DummyModel());
611  a1man->AddDiscreteProcess(theDNAa1ChargeIncreaseProcess);
612 
613  G4DNAChargeIncrease* theDNAa0ChargeIncreaseProcess =
614  new G4DNAChargeIncrease("helium_G4DNAChargeIncrease");
615  theDNAa0ChargeIncreaseProcess->AddEmModel(0, new G4DummyModel());
616  a0man->AddDiscreteProcess(theDNAa0ChargeIncreaseProcess);
617 
618  G4DNAChargeIncrease* theDNAh0ChargeIncreaseProcess =
619  new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease");
620  theDNAh0ChargeIncreaseProcess->AddEmModel(0, new G4DummyModel());
621  h0man->AddDiscreteProcess(theDNAh0ChargeIncreaseProcess);
622 
623  // start configuration of models
624  G4EmConfigurator* em_config = man->EmConfigurator();
625  G4VEmModel* mod;
626 
627  // limits for DNA model applicability
628  G4double elowest = 11 * eV;
629  G4double elimel = 1 * MeV;
630  G4double elimin = 1 * MeV;
631  G4double elimvb = 100 * eV;
632  G4double elimat = 13 * eV;
633  G4double pmin = 10 * keV;
634  G4double pminch = 100 * eV;
635  G4double mgmin = 1 * keV;
636  G4double gmmax = 500 * keV;
637  G4double pmax = 100 * MeV;
638  G4double ionmin = 10 * keV;
639  G4double ionmax = 400 * MeV;
640  G4double hmax = 100 * MeV;
641  G4double gionmax = 1 * TeV;
642 
643  // low-energy capture
644  G4LowECapture* ecap = new G4LowECapture(elowest);
645  eman->AddDiscreteProcess(ecap);
646  G4LowECapture* pcap = new G4LowECapture(pmin);
647  pman->AddDiscreteProcess(pcap);
648  G4LowECapture* icap = new G4LowECapture(ionmin);
649  iman->AddDiscreteProcess(icap);
650  G4LowECapture* a2cap = new G4LowECapture(ionmin);
651  a2man->AddDiscreteProcess(a2cap);
652  G4LowECapture* a1cap = new G4LowECapture(ionmin);
653  a1man->AddDiscreteProcess(a1cap);
654  G4LowECapture* a0cap = new G4LowECapture(ionmin);
655  a0man->AddDiscreteProcess(a0cap);
656  G4LowECapture* h0cap = new G4LowECapture(ionmin);
657  h0man->AddDiscreteProcess(h0cap);
658 
659  // loop over regions
660  for(G4int i = 0; i < nreg; ++i)
661  {
662 
663  G4String reg = regnamesDNA[i];
664  if(0 < verbose)
665  {
666  G4cout << "### DNA models are activated for G4Region " << reg << G4endl
667  << " Energy limits for e- elastic: " << elowest/eV << " eV - "
668  << elimel/MeV << " MeV" << G4endl
669  << " Energy limits for e- inelastic: " << elowest/eV << " eV - "
670  << elimin/MeV << " MeV" << G4endl
671  << " Energy limits for hadrons/ions: " << pmin/MeV << " MeV - "
672  << pmax/MeV << " MeV" << G4endl;
673  }
674  // e-
675  if(emsc)
676  {
677  G4UrbanMscModel* msc = new G4UrbanMscModel();
678  msc->SetActivationLowEnergyLimit(elimel);
679  em_config->SetExtraEmModel("e-", "msc", msc, reg);
680  }
681  else
682  {
683  mod = new G4DummyModel();
684  em_config->SetExtraEmModel("e-", "CoulombScat", mod, reg, 0.0, elimel);
685  }
686 
688  em_config->SetExtraEmModel("e-",
689  "e-_G4DNAElastic",
690  mod,
691  reg,
692  elowest,
693  elimel);
694 
695  mod = new G4MollerBhabhaModel();
696  mod->SetActivationLowEnergyLimit(elimin);
697  em_config->SetExtraEmModel("e-",
698  "eIoni",
699  mod,
700  reg,
701  0.0,
702  10 * TeV,
703  new G4UniversalFluctuation());
704 
705  mod = new G4DNABornIonisationModel();
706  em_config->SetExtraEmModel("e-",
707  "e-_G4DNAIonisation",
708  mod,
709  reg,
710  elowest,
711  elimin);
712 
713  mod = new G4DNABornExcitationModel();
714  em_config->SetExtraEmModel("e-",
715  "e-_G4DNAExcitation",
716  mod,
717  reg,
718  elowest,
719  elimin);
720 
721  mod = new G4DNASancheExcitationModel();
722  em_config->SetExtraEmModel("e-",
723  "e-_G4DNAVibExcitation",
724  mod,
725  reg,
726  elowest,
727  elimvb);
728 
729  mod = new G4DNAMeltonAttachmentModel();
730  em_config->SetExtraEmModel("e-",
731  "e-_G4DNAAttachment",
732  mod,
733  reg,
734  elowest,
735  elimat);
736 
737  // proton
738  mod = new G4BraggModel();
740  em_config->SetExtraEmModel("proton",
741  "hIoni",
742  mod,
743  reg,
744  0.0,
745  2 * MeV,
746  new G4IonFluctuations());
747 
748  mod = new G4BetheBlochModel();
749  mod->SetActivationLowEnergyLimit(pmax);
750  em_config->SetExtraEmModel("proton",
751  "hIoni",
752  mod,
753  reg,
754  2 * MeV,
755  10 * TeV,
756  new G4UniversalFluctuation());
757 
758  mod = new G4DNARuddIonisationModel();
759  em_config->SetExtraEmModel("proton",
760  "proton_G4DNAIonisation",
761  mod,
762  reg,
763  0.0,
764  gmmax);
765 
766  mod = new G4DNABornIonisationModel();
767  em_config->SetExtraEmModel("proton",
768  "proton_G4DNAIonisation",
769  mod,
770  reg,
771  gmmax,
772  pmax);
773 
775  em_config->SetExtraEmModel("proton",
776  "proton_G4DNAExcitation",
777  mod,
778  reg,
779  elowest,
780  gmmax);
781 
782  mod = new G4DNABornExcitationModel();
783  em_config->SetExtraEmModel("proton",
784  "proton_G4DNAExcitation",
785  mod,
786  reg,
787  gmmax,
788  pmax);
789 
791  em_config->SetExtraEmModel("proton",
792  "proton_G4DNAChargeDecrease",
793  mod,
794  reg,
795  pminch,
796  pmax);
797 
798  mod = new G4DNAIonElasticModel();
799  em_config->SetExtraEmModel("proton",
800  "proton_G4DNAIonElasticModel",
801  mod,
802  reg,
803  0.0,
804  1 * MeV);
805 
806  // ions
807  mod = new G4BraggIonModel();
809  em_config->SetExtraEmModel("GenericIon",
810  "ionIoni",
811  mod,
812  reg,
813  0.0,
814  2 * MeV,
815  new G4IonFluctuations());
816 
817  mod = new G4BetheBlochModel();
818  mod->SetActivationLowEnergyLimit(gionmax);
819  em_config->SetExtraEmModel("GenericIon",
820  "ionIoni",
821  mod,
822  reg,
823  2 * MeV,
824  10 * TeV,
825  new G4IonFluctuations());
826 
828  em_config->SetExtraEmModel("GenericIon",
829  "GenericIon_G4DNAIonisation",
830  mod,
831  reg,
832  0.0,
833  gionmax);
834 
835  // alpha++
836  mod = new G4BraggIonModel();
838  em_config->SetExtraEmModel("alpha",
839  "ionIoni",
840  mod,
841  reg,
842  0.0,
843  2 * MeV,
844  new G4IonFluctuations());
845 
846  mod = new G4BetheBlochModel();
847  mod->SetActivationLowEnergyLimit(pmax);
848  em_config->SetExtraEmModel("alpha",
849  "ionIoni",
850  mod,
851  reg,
852  2 * MeV,
853  10 * TeV,
854  new G4IonFluctuations());
855 
857  em_config->SetExtraEmModel("alpha",
858  "alpha_G4DNAIonisation",
859  mod,
860  reg,
861  0.0,
862  pmax);
863 
865  em_config->SetExtraEmModel("alpha",
866  "alpha_G4DNAExcitation",
867  mod,
868  reg,
869  mgmin,
870  pmax);
871 
873  em_config->SetExtraEmModel("alpha",
874  "alpha_G4DNAChargeDecrease",
875  mod,
876  reg,
877  mgmin,
878  ionmax);
879 
880  mod = new G4DNAIonElasticModel();
881  em_config->SetExtraEmModel("alpha",
882  "alpha_G4DNAIonElasticModel",
883  mod,
884  reg,
885  0.0,
886  1 * MeV);
887 
888  // alpha+
889  mod = new G4BraggIonModel();
891  em_config->SetExtraEmModel("alpha+",
892  "ionIoni",
893  mod,
894  reg,
895  0.0,
896  2 * MeV,
897  new G4IonFluctuations());
898 
899  mod = new G4BetheBlochModel();
900  mod->SetActivationLowEnergyLimit(pmax);
901  em_config->SetExtraEmModel("alpha+",
902  "ionIoni",
903  mod,
904  reg,
905  2 * MeV,
906  10 * TeV,
907  new G4IonFluctuations());
908 
909  mod = new G4DNARuddIonisationModel();
910  em_config->SetExtraEmModel("alpha+",
911  "alpha+_G4DNAIonisation",
912  mod,
913  reg,
914  0.0,
915  pmax);
916 
918  em_config->SetExtraEmModel("alpha+",
919  "alpha+_G4DNAExcitation",
920  mod,
921  reg,
922  mgmin,
923  pmax);
924 
926  em_config->SetExtraEmModel("alpha+",
927  "alpha+_G4DNAChargeDecrease",
928  mod,
929  reg,
930  mgmin,
931  ionmax);
932 
934  em_config->SetExtraEmModel("alpha+",
935  "alpha+_G4DNAChargeIncrease",
936  mod,
937  reg,
938  mgmin,
939  ionmax);
940 
941  mod = new G4DNAIonElasticModel();
942  em_config->SetExtraEmModel("alpha+",
943  "alpha+_G4DNAIonElasticModel",
944  mod,
945  reg,
946  0.0,
947  1 * MeV);
948 
949  // helium
950  mod = new G4DNARuddIonisationModel();
951  em_config->SetExtraEmModel("helium",
952  "helium_G4DNAIonisation",
953  mod,
954  reg,
955  0.0,
956  pmax);
957 
959  em_config->SetExtraEmModel("helium",
960  "helium_G4DNAExcitation",
961  mod,
962  reg,
963  mgmin,
964  pmax);
965 
967  em_config->SetExtraEmModel("helium",
968  "helium_G4DNAChargeIncrease",
969  mod,
970  reg,
971  mgmin,
972  ionmax);
973 
974  mod = new G4DNAIonElasticModel();
975  em_config->SetExtraEmModel("helium",
976  "helium_G4DNAIonElasticModel",
977  mod,
978  reg,
979  0.0,
980  1 * MeV);
981 
982  // hydrogen
983  mod = new G4DNARuddIonisationModel();
984  em_config->SetExtraEmModel("hydrogen",
985  "hydrogen_G4DNAIonisation",
986  mod,
987  reg,
988  0.0,
989  hmax);
990 
992  em_config->SetExtraEmModel("hydrogen",
993  "hydrogen_G4DNAExcitation",
994  mod,
995  reg,
996  elowest,
997  gmmax);
998 
1000  em_config->SetExtraEmModel("hydrogen",
1001  "hydrogen_G4DNAChargeIncrease",
1002  mod,
1003  reg,
1004  pminch,
1005  pmax);
1006 
1007  mod = new G4DNAIonElasticModel();
1008  em_config->SetExtraEmModel("hydrogen",
1009  "hydrogen_G4DNAIonElasticModel",
1010  mod,
1011  reg,
1012  0.0,
1013  1 * MeV);
1014  }
1015 }
G4EmConfigurator * EmConfigurator()
void SetActivationHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:739
static const double MeV
Definition: G4SIunits.hh:211
#define G4DNABornIonisationModel
static G4LossTableManager * Instance()
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4ProcessManager * GetProcessManager() const
int G4int
Definition: G4Types.hh:78
G4EmParameters * theParameters
G4DNABornExcitationModel1 G4DNABornExcitationModel
const std::vector< G4String > & RegionsDNA() const
G4int Verbose() const
static const G4double reg
G4GLOB_DLL std::ostream G4cout
const std::vector< G4String > & TypesDNA() const
bool G4bool
Definition: G4Types.hh:79
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4DNAGenericIonsManager * Instance(void)
G4bool HasMsc(G4ProcessManager *) const
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
static const double eV
Definition: G4SIunits.hh:212
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:93
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
static const double TeV
Definition: G4SIunits.hh:215
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
void SetExtraEmModel(const G4String &particleName, const G4String &processName, G4VEmModel *, const G4String &regionName="", G4double emin=0.0, G4double emax=DBL_MAX, G4VEmFluctuationModel *fm=0)
G4ParticleDefinition * GetIon(const G4String &name)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ActivateMicroElec()

void G4EmModelActivator::ActivateMicroElec ( )
private

Definition at line 249 of file G4EmModelActivator.cc.

250 {
251  const std::vector<G4String> regnamesME = theParameters->RegionsMicroElec();
252  G4int nreg = regnamesME.size();
253  if(0 == nreg)
254  {
255  return;
256  }
257  G4int verbose = theParameters->Verbose() - 1;
258  if(verbose > 0)
259  {
260  G4cout << "### G4EmModelActivator::ActivateMicroElec for " << nreg
261  << " regions" << G4endl;
262  }
264 
266  const G4ParticleDefinition* prot = G4Proton::Proton();
268  G4ProcessManager* eman = elec->GetProcessManager();
269  G4ProcessManager* pman = prot->GetProcessManager();
270  G4ProcessManager* iman = gion->GetProcessManager();
271 
272  G4bool emsc = HasMsc(eman);
273 
274  // MicroElec elastic is not active in the world
275  G4MicroElecElastic* eElasProc =
276  new G4MicroElecElastic("e-G4MicroElecElastic");
277  eman->AddDiscreteProcess(eElasProc);
278 
279  G4MicroElecInelastic* eInelProc =
280  new G4MicroElecInelastic("e-G4MicroElecInelastic");
281  eman->AddDiscreteProcess(eInelProc);
282 
283  G4MicroElecInelastic* pInelProc =
284  new G4MicroElecInelastic("p_G4MicroElecInelastic");
285  pman->AddDiscreteProcess(pInelProc);
286 
287  G4MicroElecInelastic* iInelProc =
288  new G4MicroElecInelastic("ion_G4MicroElecInelastic");
289  iman->AddDiscreteProcess(iInelProc);
290 
291  // start configuration of models
292  G4EmConfigurator* em_config = man->EmConfigurator();
293  G4VEmModel* mod;
294 
295  // limits for MicroElec applicability
296  G4double elowest = 16.7 * eV;
297  G4double elimel = 100 * MeV;
298  G4double elimin = 9 * MeV;
299  G4double pmin = 50 * keV;
300  G4double pmax = 99.9 * MeV;
301 
302  G4LowECapture* ecap = new G4LowECapture(elowest);
303  eman->AddDiscreteProcess(ecap);
304 
305  for(G4int i = 0; i < nreg; ++i)
306  {
307 
308  G4String reg = regnamesME[i];
309  G4cout << "### MicroElec models are activated for G4Region " << reg
310  << G4endl
311  << " Energy limits for e- elastic: " << elowest/eV << " eV - "
312  << elimel/MeV << " MeV"
313  << G4endl
314  << " Energy limits for e- inelastic: " << elowest/eV << " eV - "
315  << elimin/MeV << " MeV"
316  << G4endl
317  << " Energy limits for hadrons/ions: " << pmin/MeV << " MeV - "
318  << pmax/MeV << " MeV"
319  << G4endl;
320 
321  // e-
322  if(emsc)
323  {
324  G4UrbanMscModel* msc = new G4UrbanMscModel();
325  msc->SetActivationLowEnergyLimit(elimel);
326  em_config->SetExtraEmModel("e-", "msc", msc, reg);
327  }
328  else
329  {
330  mod = new G4DummyModel();
331  em_config->SetExtraEmModel("e-", "CoulombScat", mod, reg, 0.0, elimel);
332  }
333 
334  mod = new G4MicroElecElasticModel();
335  em_config->SetExtraEmModel("e-",
336  "e-G4MicroElecElastic",
337  mod,
338  reg,
339  elowest,
340  elimin);
341 
342  mod = new G4MollerBhabhaModel();
343  mod->SetActivationLowEnergyLimit(elimin);
344  em_config->SetExtraEmModel("e-",
345  "eIoni",
346  mod,
347  reg,
348  0.0,
349  10 * TeV,
350  new G4UniversalFluctuation());
351 
352  mod = new G4MicroElecInelasticModel();
353  em_config->SetExtraEmModel("e-",
354  "e-G4MicroElecInelastic",
355  mod,
356  reg,
357  elowest,
358  elimin);
359 
360  // proton
361  mod = new G4BraggModel();
362  mod->SetActivationHighEnergyLimit(pmin);
363  em_config->SetExtraEmModel("proton",
364  "hIoni",
365  mod,
366  reg,
367  0.0,
368  2 * MeV,
369  new G4IonFluctuations());
370 
371  mod = new G4BetheBlochModel();
372  mod->SetActivationLowEnergyLimit(pmax);
373  em_config->SetExtraEmModel("proton",
374  "hIoni",
375  mod,
376  reg,
377  2 * MeV,
378  10 * TeV,
379  new G4UniversalFluctuation());
380 
381  mod = new G4MicroElecInelasticModel();
382  em_config->SetExtraEmModel("proton",
383  "p_G4MicroElecInelastic",
384  mod,
385  reg,
386  pmin,
387  pmax);
388 
389  // ions
390  mod = new G4BraggIonModel();
391  mod->SetActivationHighEnergyLimit(pmin);
392  em_config->SetExtraEmModel("GenericIon",
393  "ionIoni",
394  mod,
395  reg,
396  0.0,
397  2 * MeV,
398  new G4IonFluctuations());
399 
400  mod = new G4BetheBlochModel();
401  mod->SetActivationLowEnergyLimit(pmax);
402  em_config->SetExtraEmModel("GenericIon",
403  "ionIoni",
404  mod,
405  reg,
406  2 * MeV,
407  10 * TeV,
408  new G4IonFluctuations());
409 
410  mod = new G4MicroElecInelasticModel();
411  em_config->SetExtraEmModel("GenericIon",
412  "ion_G4MicroElecInelastic",
413  mod,
414  reg,
415  pmin,
416  pmax);
417  }
418 }
G4EmConfigurator * EmConfigurator()
void SetActivationHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:739
static const double MeV
Definition: G4SIunits.hh:211
static G4LossTableManager * Instance()
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4ProcessManager * GetProcessManager() const
int G4int
Definition: G4Types.hh:78
G4EmParameters * theParameters
G4int Verbose() const
static const G4double reg
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
const std::vector< G4String > & RegionsMicroElec() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4bool HasMsc(G4ProcessManager *) const
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
static const double eV
Definition: G4SIunits.hh:212
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:93
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
static const double TeV
Definition: G4SIunits.hh:215
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
void SetExtraEmModel(const G4String &particleName, const G4String &processName, G4VEmModel *, const G4String &regionName="", G4double emin=0.0, G4double emax=DBL_MAX, G4VEmFluctuationModel *fm=0)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ActivatePAI()

void G4EmModelActivator::ActivatePAI ( )
private

Definition at line 152 of file G4EmModelActivator.cc.

153 {
154  const std::vector<G4String> regnamesPAI = theParameters->RegionsPAI();
155  G4int nreg = regnamesPAI.size();
156  if(0 == nreg) return;
157  G4int verbose = theParameters->Verbose() - 1;
158  if(verbose > 0)
159  {
160  G4cout << "### G4EmModelActivator::ActivatePAI for " << nreg << " regions"
161  << G4endl;
162  }
163  const std::vector<G4String> particlesPAI = theParameters->ParticlesPAI();
164  const std::vector<G4String> typesPAI = theParameters->TypesPAI();
165 
166  const std::vector<G4VEnergyLossProcess*>& v = G4LossTableManager::Instance()
168  std::vector<G4VEnergyLossProcess*>::const_iterator itr;
169  G4RegionStore* regionStore = G4RegionStore::GetInstance();
170 
176 
177  for(G4int i = 0; i < nreg; ++i)
178  {
179  G4bool res = true;
180  const G4ParticleDefinition* p = 0;
181  if(particlesPAI[i] != "all")
182  {
183  p = G4ParticleTable::GetParticleTable()->FindParticle(particlesPAI[i]);
184  if(!p)
185  {
186  G4cout << "### WARNING: ActivatePAI::FindParticle fails to find "
187  << particlesPAI[i] << G4endl;
188  res = false;
189  }
190  }
191 
192  if(res)
193  {
194  const G4Region* r = regionStore->GetRegion(regnamesPAI[i], false);
195  if(!r)
196  {
197  G4cout << "### WARNING: ActivatePAI::GetRegion fails to find "
198  << regnamesPAI[i] << G4endl;
199  }
200  else
201  {
202 
203  G4String name = "hIoni";
204  if(p == elec || p == posi)
205  { name = "eIoni";}
206  else if (p == mupl || p == mumi)
207  { name = "muIoni";}
208  else if (p == gion)
209  { name = "ionIoni";}
210 
211  for(itr = v.begin(); itr != v.end(); itr++)
212  {
213  G4VEnergyLossProcess* proc = *itr;
214  if(proc->IsIonisationProcess())
215  {
216  if(p == 0 || (p != 0 && name == proc->GetProcessName()))
217  {
218  G4VEmModel* em = 0;
219  G4VEmFluctuationModel* fm = 0;
220  if(typesPAI[i] == "PAIphoton")
221  {
222  G4PAIPhotModel* mod = new G4PAIPhotModel(p,"PAIPhotModel");
223  em = mod;
224  fm = mod;
225  }
226  else
227  {
228  G4PAIModel* mod = new G4PAIModel(p,"PAIModel");
229  em = mod;
230  fm = mod;
231  }
232  proc->AddEmModel(0, em, fm, r);
233  if(verbose > 0)
234  {
235  G4cout << "### G4EmModelActivator: add <" << typesPAI[i]
236  << "> model for " << particlesPAI[i]
237  << " in the " << regnamesPAI[i] << G4endl;
238  }
239  }
240  }
241  }
242  }
243  }
244  }
245 }
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4LossTableManager * Instance()
G4bool IsIonisationProcess() const
G4String name
Definition: TRTMaterials.hh:40
const std::vector< G4String > & RegionsPAI() const
int G4int
Definition: G4Types.hh:78
G4EmParameters * theParameters
G4int Verbose() const
static G4RegionStore * GetInstance()
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
const std::vector< G4String > & TypesPAI() const
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:93
static G4ParticleTable * GetParticleTable()
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
const std::vector< G4String > & ParticlesPAI() const
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructDNAParticles()

void G4EmModelActivator::ConstructDNAParticles ( )
private

Definition at line 422 of file G4EmModelActivator.cc.

423 {
424  // this is only addition on top of any Physics Lists
425  G4Alpha::Alpha();
426 
427  G4DNAGenericIonsManager * genericIonsManager =
429  genericIonsManager->GetIon("alpha+");
430  genericIonsManager->GetIon("helium");
431  genericIonsManager->GetIon("hydrogen");
432 }
static G4DNAGenericIonsManager * Instance(void)
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
G4ParticleDefinition * GetIon(const G4String &name)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructParticle()

void G4EmModelActivator::ConstructParticle ( void  )

Definition at line 120 of file G4EmModelActivator.cc.

121 {
122  const std::vector<G4String> regnamesDNA = theParameters->RegionsDNA();
123  if(regnamesDNA.size() > 0)
124  {
126  }
127 }
G4EmParameters * theParameters
const std::vector< G4String > & RegionsDNA() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructProcess()

void G4EmModelActivator::ConstructProcess ( void  )

Definition at line 131 of file G4EmModelActivator.cc.

132 {
133  const std::vector<G4String> regnamesPAI = theParameters->RegionsPAI();
134  if(regnamesPAI.size() > 0)
135  {
136  ActivatePAI();
137  }
138  const std::vector<G4String> regnamesME = theParameters->RegionsMicroElec();
139  if(regnamesME.size() > 0)
140  {
142  }
143  const std::vector<G4String> regnamesDNA = theParameters->RegionsDNA();
144  if(regnamesDNA.size() > 0)
145  {
146  ActivateDNA();
147  }
148 }
const std::vector< G4String > & RegionsPAI() const
G4EmParameters * theParameters
const std::vector< G4String > & RegionsDNA() const
const std::vector< G4String > & RegionsMicroElec() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ HasMsc()

G4bool G4EmModelActivator::HasMsc ( G4ProcessManager pm) const
private

Definition at line 1019 of file G4EmModelActivator.cc.

1020 {
1021  G4bool res = false;
1022  G4ProcessVector* pv = pm->GetProcessList();
1023  G4int nproc = pm->GetProcessListLength();
1024  for(G4int i = 0; i < nproc; ++i)
1025  {
1026  if(((*pv)[i])->GetProcessSubType() == fMultipleScattering)
1027  {
1028  res = true;
1029  break;
1030  }
1031  }
1032  return res;
1033 }
G4ProcessVector * GetProcessList() const
int G4int
Definition: G4Types.hh:78
bool G4bool
Definition: G4Types.hh:79
G4int GetProcessListLength() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

G4EmModelActivator& G4EmModelActivator::operator= ( const G4EmModelActivator right)
private

Member Data Documentation

◆ theParameters

G4EmParameters* G4EmModelActivator::theParameters
private

Definition at line 83 of file G4EmModelActivator.hh.


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