CMS 3D CMS Logo

FTFPCMS_BERT_EMM.cc
Go to the documentation of this file.
1 #include "FTFPCMS_BERT_EMM.h"
5 
6 #include "G4DecayPhysics.hh"
7 #include "G4EmExtraPhysics.hh"
8 #include "G4IonPhysics.hh"
9 #include "G4StoppingPhysics.hh"
10 #include "G4HadronElasticPhysics.hh"
11 #include "G4HadronicParameters.hh"
12 
14  int ver = p.getUntrackedParameter<int>("Verbosity", 0);
15  bool emPhys = p.getUntrackedParameter<bool>("EMPhysics", true);
16  bool hadPhys = p.getUntrackedParameter<bool>("HadPhysics", true);
17  double minFTFP = p.getParameter<double>("EminFTFP") * CLHEP::GeV;
18  double maxBERT = p.getParameter<double>("EmaxBERT") * CLHEP::GeV;
19  double maxBERTpi = p.getParameter<double>("EmaxBERTpi") * CLHEP::GeV;
20  edm::LogVerbatim("PhysicsList") << "CMS Physics List FTFP_BERT_EMM: "
21  << "\n Flags for EM Physics: " << emPhys << "; Hadronic Physics: " << hadPhys
22  << "\n Transition energy Bertini/FTFP from " << minFTFP / CLHEP::GeV << " to "
23  << maxBERT / CLHEP::GeV << "; for pions to " << maxBERTpi / CLHEP::GeV << " GeV";
24 
25  if (emPhys) {
26  // EM Physics
27  RegisterPhysics(new CMSEmStandardPhysics(ver, p));
28 
29  // Synchroton Radiation & GN Physics
30  G4EmExtraPhysics* gn = new G4EmExtraPhysics(ver);
31  RegisterPhysics(gn);
32  bool mu = p.getParameter<bool>("G4MuonPairProductionByMuon");
33  gn->MuonToMuMu(mu);
34  edm::LogVerbatim("PhysicsList") << " Muon pair production by muons: " << mu;
35  }
36 
37  // Decays
38  this->RegisterPhysics(new G4DecayPhysics(ver));
39 
40  if (hadPhys) {
41  bool ngen = p.getParameter<bool>("G4NeutronGeneralProcess");
42  bool bc = p.getParameter<bool>("G4BCHadronicProcess");
43  bool hn = p.getParameter<bool>("G4LightHyperNucleiTracking");
44  auto param = G4HadronicParameters::Instance();
45  param->SetEnableNeutronGeneralProcess(ngen);
46  param->SetEnableBCParticles(bc);
47  param->SetEnableHyperNuclei(hn);
48  edm::LogVerbatim("PhysicsList") << " Eneble neutron general process: " << ngen
49  << "\n Enable b- and c- hadron physics: " << bc
50  << "\n Enable light hyper-nuclei physics: " << hn;
51 
52  // Hadron Elastic scattering
53  RegisterPhysics(new G4HadronElasticPhysics(ver));
54 
55  // Hadron Physics
56  RegisterPhysics(new CMSHadronPhysicsFTFP_BERT(minFTFP, maxBERT, maxBERTpi, minFTFP, maxBERT));
57 
58  // Stopping Physics
59  RegisterPhysics(new G4StoppingPhysics(ver));
60 
61  // Ion Physics
62  RegisterPhysics(new G4IonPhysics(ver));
63  }
64 }
Log< level::Info, true > LogVerbatim
FTFPCMS_BERT_EMM(const edm::ParameterSet &p)