CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HiggsValidation.cc
Go to the documentation of this file.
1 /*class HiggsValidation
2  *
3  * Class to fill dqm monitor elements from existing EDM file
4  *
5  * $Date: 2012/08/12 16:13:29 $
6  * $Revision: 1.1 $
7  */
8 
10 
12 
13 #include "CLHEP/Units/defs.h"
14 #include "CLHEP/Units/PhysicalConstants.h"
15 
17 
20 
21 using namespace edm;
22 
24  _wmanager(iPSet),
25  hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
26  particle_id(iPSet.getParameter<int>("pdg_id")),
27  particle_name(iPSet.getParameter<std::string>("particleName"))
28 {
29  dbe = 0;
31 
32  monitoredDecays = new MonitoredDecays(iPSet);
33 
34 }
35 
37 
39 {
40  if(dbe){
42  TString dir="Generator/";
43  dir+=particle_name;
44  dbe->setCurrentFolder(dir.Data());
45 
46  // Number of analyzed events
47  nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
48 
49  //decay type
50 
51  std::string channel = particle_name+"_DecayChannels";
52  HiggsDecayChannels = dbe->book1D(channel.c_str(),(particle_name+" decay channels").c_str(),monitoredDecays->size(),0,monitoredDecays->size());
53 
54  for(size_t i = 0; i < monitoredDecays->size(); ++i){
56  }
57  }
58 
59  //Kinematics
60  Higgs_pt = dbe->book1D((particle_name+"_pt"),(particle_name+" p_{t}"),50,0,250);
61  Higgs_eta = dbe->book1D((particle_name+"_eta"),(particle_name+" #eta"),50,-5,5);
62  Higgs_mass = dbe->book1D((particle_name+"_m"),(particle_name+" M"),500,0,500);
63 
64  int idx=0;
65  for(unsigned int i=0;i<monitoredDecays->NDecayParticles();i++){
66  HiggsDecayProd_pt.push_back(dbe->book1D((monitoredDecays->ConvertIndex(idx)+"_pt"),(monitoredDecays->ConvertIndex(idx)+" p_{t}"),50,0,250));
67  HiggsDecayProd_eta.push_back(dbe->book1D((monitoredDecays->ConvertIndex(idx)+"_eta"),(monitoredDecays->ConvertIndex(idx)+" #eta"),50,-5,5));
68  idx++;
69  }
70 
71  return;
72 }
73 
75  return;
76 }
77 
78 void HiggsValidation::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup)
79 {
81  iSetup.getData( fPDGTable );
82  return;
83 }
84 void HiggsValidation::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
86 {
87  double weight = _wmanager.weight(iEvent);
88  nEvt->Fill(0.5,weight);
89 
90  //Gathering the HepMCProduct information
92  iEvent.getByLabel(hepmcCollection_, evt);
93 
94  //Get EVENT
95  HepMC::GenEvent *myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
96 
97  // loop over all particles
98  bool filled = false;
99  for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin();
100  iter!= myGenEvent->particles_end() && !filled; ++iter) {
101  if(particle_id == fabs((*iter)->pdg_id())){
102  std::vector<HepMC::GenParticle*> decayprod;
103  int channel = findHiggsDecayChannel(*iter,decayprod);
104  HiggsDecayChannels->Fill(channel,weight);
105  Higgs_pt->Fill((*iter)->momentum().perp(),weight);
106  Higgs_eta->Fill((*iter)->momentum().eta(),weight);
107  Higgs_mass->Fill((*iter)->momentum().m(),weight);
108  for(unsigned int i=0;i<decayprod.size();i++){
109  int idx=monitoredDecays->isDecayParticle(decayprod.at(i)->pdg_id());
110  if(0<=idx && idx<=(int)HiggsDecayProd_pt.size()){
111  HiggsDecayProd_pt.at(idx)->Fill(decayprod.at(i)->momentum().perp(),weight);
112  HiggsDecayProd_eta.at(idx)->Fill(decayprod.at(i)->momentum().eta(),weight);
113  }
114  }
115  filled = true;
116  }
117  }
118 
119  delete myGenEvent;
120 
121 }//analyze
122 
123 int HiggsValidation::findHiggsDecayChannel(const HepMC::GenParticle* genParticle,std::vector<HepMC::GenParticle*> &decayprod){
124  if(genParticle->status() == 1) return monitoredDecays->stable();
125  std::vector<int> children;
126  if ( genParticle->end_vertex() ) {
127  HepMC::GenVertex::particle_iterator des;
128  for(des = genParticle->end_vertex()->particles_begin(HepMC::descendants);
129  des!= genParticle->end_vertex()->particles_end(HepMC::descendants);++des ) {
130 
131  if((*des)->pdg_id() == genParticle->pdg_id()) continue;
132 
133  HepMC::GenVertex::particle_iterator mother = (*des)->production_vertex()->particles_begin(HepMC::parents);
134  if((*mother)->pdg_id() == genParticle->pdg_id()){
135  children.push_back((*des)->pdg_id());
136  decayprod.push_back((*des));
137  }
138  }
139  }
140 
141  if(children.size() == 2 && children.at(0) != 0 && children.at(1) != 0){
142  return monitoredDecays->position(children.at(0),children.at(1));
143  }
144  return monitoredDecays->undetermined();
145 }
146 
147 //define this as a plug-in
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
PDT table.
virtual void analyze(const edm::Event &, const edm::EventSetup &)
int i
Definition: DBlmapReader.cc:9
TPRegexp parents
Definition: eve_filter.cc:24
MonitorElement * Higgs_eta
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:717
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void endRun(const edm::Run &, const edm::EventSetup &)
WeightManager _wmanager
virtual void beginRun(const edm::Run &, const edm::EventSetup &)
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
MonitorElement * nEvt
std::vector< MonitorElement * > HiggsDecayProd_pt
size_t position(int pid1, int pid2)
MonitoredDecays * monitoredDecays
virtual ~HiggsValidation()
void getData(T &iHolder) const
Definition: EventSetup.h:67
void Fill(long long x)
HiggsValidation(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:243
int findHiggsDecayChannel(const HepMC::GenParticle *, std::vector< HepMC::GenParticle * > &decayprod)
std::vector< MonitorElement * > HiggsDecayProd_eta
MonitorElement * Higgs_pt
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
MonitorElement * HiggsDecayChannels
edm::InputTag hepmcCollection_
DQMStore * dbe
ME&#39;s &quot;container&quot;.
virtual void beginJob()
std::string ConvertIndex(int index)
std::string channel(size_t i)
virtual void endJob()
dbl *** dir
Definition: mlp_gen.cc:35
std::string particle_name
double weight(const edm::Event &)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:429
Definition: Run.h:33
MonitorElement * Higgs_mass