CMS 3D CMS Logo

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