CMS 3D CMS Logo

List of all members | Public Member Functions
FTFPCMS_BERT_HP_EML Class Reference

#include <FTFPCMS_BERT_HP_EML.h>

Inheritance diagram for FTFPCMS_BERT_HP_EML:
PhysicsList

Public Member Functions

 FTFPCMS_BERT_HP_EML (G4LogicalVolumeToDDLogicalPartMap &map, const HepPDT::ParticleDataTable *table_, sim::ChordFinderSetter *chordFinderSetter_, const edm::ParameterSet &p)
 
- Public Member Functions inherited from PhysicsList
 PhysicsList (G4LogicalVolumeToDDLogicalPartMap &map, const HepPDT::ParticleDataTable *table_, sim::ChordFinderSetter *chordFinderSetter_, const edm::ParameterSet &p)
 
void SetCuts () override
 
 ~PhysicsList () override
 

Detailed Description

Definition at line 7 of file FTFPCMS_BERT_HP_EML.h.

Constructor & Destructor Documentation

FTFPCMS_BERT_HP_EML::FTFPCMS_BERT_HP_EML ( G4LogicalVolumeToDDLogicalPartMap map,
const HepPDT::ParticleDataTable table_,
sim::ChordFinderSetter chordFinderSetter_,
const edm::ParameterSet p 
)

Definition at line 19 of file FTFPCMS_BERT_HP_EML.cc.

References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), and muons2muons_cfi::photon.

23  : PhysicsList(map, table_, chordFinderSetter_, p) {
24 
25  G4DataQuestionaire it(photon);
26 
27  int ver = p.getUntrackedParameter<int>("Verbosity",0);
28  bool emPhys = p.getUntrackedParameter<bool>("EMPhysics",true);
29  bool hadPhys = p.getUntrackedParameter<bool>("HadPhysics",true);
30  bool tracking= p.getParameter<bool>("TrackingCut");
31  bool thermal = p.getUntrackedParameter<bool>("ThermalNeutrons");
32  double timeLimit = p.getParameter<double>("MaxTrackTime")*ns;
33  edm::LogInfo("PhysicsList") << "You are using the simulation engine: "
34  << "FTFP_BERT_HP_EML \n Flags for EM Physics "
35  << emPhys << ", for Hadronic Physics "
36  << hadPhys << " and tracking cut " << tracking
37  << " t(ns)= " << timeLimit/ns
38  << " ThermalNeutrons: " << thermal;
39 
40  if (emPhys) {
41  // EM Physics
42  RegisterPhysics( new CMSEmStandardPhysicsXS(ver));
43 
44  // Synchroton Radiation & GN Physics
45  G4EmExtraPhysics* gn = new G4EmExtraPhysics(ver);
46  RegisterPhysics(gn);
47  }
48 
49  // Decays
50  this->RegisterPhysics( new G4DecayPhysics(ver) );
51 
52  if (hadPhys) {
53  G4HadronicProcessStore::Instance()->SetVerbose(ver);
54 
55  // Hadron Elastic scattering
56  RegisterPhysics( new G4HadronElasticPhysicsHP(ver));
57 
58  // Hadron Physics
59  RegisterPhysics( new G4HadronPhysicsFTFP_BERT_HP(ver));
60 
61  // Stopping Physics
62  RegisterPhysics( new G4StoppingPhysics(ver));
63 
64  // Ion Physics
65  RegisterPhysics( new G4IonPhysics(ver));
66 
67  // Neutron tracking cut
68  if (tracking) {
69  G4NeutronTrackingCut* ncut= new G4NeutronTrackingCut(ver);
70  ncut->SetTimeLimit(timeLimit);
71  RegisterPhysics(ncut);
72  }
73  if(thermal) {
74  RegisterPhysics(new CMSThermalNeutrons(ver));
75  }
76  }
77 
78  // Monopoles
79  //RegisterPhysics( new CMSMonopolePhysics(table_,chordFinderSetter_,p));
80 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
PhysicsList(G4LogicalVolumeToDDLogicalPartMap &map, const HepPDT::ParticleDataTable *table_, sim::ChordFinderSetter *chordFinderSetter_, const edm::ParameterSet &p)
Definition: PhysicsList.cc:3
Table table_