CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
QGSPCMS_FTFP_BERT_EML.cc
Go to the documentation of this file.
1 #include "QGSPCMS_FTFP_BERT_EML.hh"
6 
7 #include "G4DecayPhysics.hh"
8 #include "G4EmExtraPhysics.hh"
9 #include "G4IonPhysics.hh"
10 #include "G4StoppingPhysics.hh"
11 #include "G4HadronElasticPhysics.hh"
12 #include "G4NeutronTrackingCut.hh"
13 #include "G4HadronicProcessStore.hh"
14 
15 #include "G4DataQuestionaire.hh"
16 #include "G4HadronPhysicsQGSP_FTFP_BERT.hh"
17 
18 QGSPCMS_FTFP_BERT_EML::QGSPCMS_FTFP_BERT_EML(G4LogicalVolumeToDDLogicalPartMap& map,
20  sim::ChordFinderSetter *chordFinderSetter_,
21  const edm::ParameterSet & p) : PhysicsList(map, table_, chordFinderSetter_, p) {
22 
23  G4DataQuestionaire it(photon);
24 
25  int ver = p.getUntrackedParameter<int>("Verbosity",0);
26  bool emPhys = p.getUntrackedParameter<bool>("EMPhysics",true);
27  bool hadPhys = p.getUntrackedParameter<bool>("HadPhysics",true);
28  bool tracking= p.getParameter<bool>("TrackingCut");
29  double timeLimit = p.getParameter<double>("MaxTrackTime")*ns;
30  edm::LogInfo("PhysicsList") << "You are using the simulation engine: "
31  << "QGSP_FTFP_BERT_EML \n Flags for EM Physics "
32  << emPhys << ", for Hadronic Physics "
33  << hadPhys << " and tracking cut " << tracking
34  << " t(ns)= " << timeLimit/ns;
35 
36  if (emPhys) {
37  // EM Physics
38  RegisterPhysics(new CMSEmStandardPhysics(ver));
39 
40  // Synchroton Radiation & GN Physics
41  G4EmExtraPhysics* gn = new G4EmExtraPhysics(ver);
42  RegisterPhysics(gn);
43  }
44 
45  // Decays
46  RegisterPhysics(new G4DecayPhysics(ver) );
47 
48  if (hadPhys) {
49  G4HadronicProcessStore::Instance()->SetVerbose(ver);
50 
51  // Hadron Elastic scattering
52  RegisterPhysics(new G4HadronElasticPhysics(ver));
53 
54  // Hadron Physics
55  RegisterPhysics(new G4HadronPhysicsQGSP_FTFP_BERT(ver));
56 
57  // Stopping Physics
58  RegisterPhysics(new G4StoppingPhysics(ver));
59 
60  // Ion Physics
61  RegisterPhysics(new G4IonPhysics(ver));
62 
63  // Neutron tracking cut
64  if (tracking) {
65  G4NeutronTrackingCut* ncut= new G4NeutronTrackingCut(ver);
66  ncut->SetTimeLimit(timeLimit);
67  RegisterPhysics(ncut);
68  }
69  }
70 
71  // Monopoles
72  RegisterPhysics( new CMSMonopolePhysics(table_,chordFinderSetter_,p));
73 }
74 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
HepPDT::ParticleDataTable ParticleDataTable
Table table_