CMS 3D CMS Logo

FTFPCMS_BERT_EMH.cc
Go to the documentation of this file.
1 #include "FTFPCMS_BERT_EMH.h"
5 
6 #include "G4DecayPhysics.hh"
7 #include "G4IonPhysics.hh"
8 #include "G4StoppingPhysics.hh"
9 #include "G4HadronElasticPhysics.hh"
10 
12  int ver = p.getUntrackedParameter<int>("Verbosity", 0);
13  bool emPhys = p.getUntrackedParameter<bool>("EMPhysics", true);
14  bool hadPhys = p.getUntrackedParameter<bool>("HadPhysics", true);
15  double minFTFP = p.getParameter<double>("EminFTFP") * CLHEP::GeV;
16  double maxBERT = p.getParameter<double>("EmaxBERT") * CLHEP::GeV;
17  double maxBERTpi = p.getParameter<double>("EmaxBERTpi") * CLHEP::GeV;
18  edm::LogVerbatim("PhysicsList") << "CMS Physics List FTFP_BERT_EMH: "
19  << "\n Flags for EM Physics: " << emPhys << "; Hadronic Physics: " << hadPhys
20  << "\n Transition energy Bertini/FTFP from " << minFTFP / CLHEP::GeV << " to "
21  << maxBERT / CLHEP::GeV << "; for pions to " << maxBERTpi / CLHEP::GeV << " GeV";
22 
23  if (emPhys) {
24  // EM Physics
25  RegisterPhysics(new CMSEmStandardPhysicsEMH(ver, p));
26  }
27 
28  // Decays
29  this->RegisterPhysics(new G4DecayPhysics(ver));
30 
31  if (hadPhys) {
32  // Hadron Elastic scattering
33  RegisterPhysics(new G4HadronElasticPhysics(ver));
34 
35  // Hadron Physics
36  RegisterPhysics(new CMSHadronPhysicsFTFP_BERT(minFTFP, maxBERT, maxBERTpi, minFTFP, maxBERT));
37 
38  // Stopping Physics
39  RegisterPhysics(new G4StoppingPhysics(ver));
40 
41  // Ion Physics
42  RegisterPhysics(new G4IonPhysics(ver));
43  }
44 }
Log< level::Info, true > LogVerbatim
FTFPCMS_BERT_EMH(const edm::ParameterSet &p)