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  */
6 
8 
10 
11 #include "CLHEP/Units/defs.h"
12 #include "CLHEP/Units/PhysicalConstants.h"
13 
15 
18 
19 using namespace edm;
20 
22  wmanager_(iPSet,consumesCollector()),
23  hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
24  particle_id(iPSet.getParameter<int>("pdg_id")),
25  particle_name(iPSet.getParameter<std::string>("particleName"))
26 {
27  dbe = 0;
29 
30  monitoredDecays = new MonitoredDecays(iPSet);
31 
32  hepmcCollectionToken_=consumes<HepMCProduct>(hepmcCollection_);
33 }
34 
36 
38 {
39  if(dbe){
41  TString dir="Generator/";
42  dir+=particle_name;
43  dbe->setCurrentFolder(dir.Data());
44 
45  // Number of analyzed events
46  nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
47 
48  //decay type
49 
50  std::string channel = particle_name+"_DecayChannels";
51  HiggsDecayChannels = dbe->book1D(channel.c_str(),(particle_name+" decay channels").c_str(),monitoredDecays->size(),0,monitoredDecays->size());
52 
53  for(size_t i = 0; i < monitoredDecays->size(); ++i){
55  }
56  }
57 
58  //Kinematics
59  Higgs_pt = dbe->book1D((particle_name+"_pt"),(particle_name+" p_{t}"),50,0,250);
60  Higgs_eta = dbe->book1D((particle_name+"_eta"),(particle_name+" #eta"),50,-5,5);
61  Higgs_mass = dbe->book1D((particle_name+"_m"),(particle_name+" M"),500,0,500);
62 
63  int idx=0;
64  for(unsigned int i=0;i<monitoredDecays->NDecayParticles();i++){
65  HiggsDecayProd_pt.push_back(dbe->book1D((monitoredDecays->ConvertIndex(idx)+"_pt"),(monitoredDecays->ConvertIndex(idx)+" p_{t}"),50,0,250));
66  HiggsDecayProd_eta.push_back(dbe->book1D((monitoredDecays->ConvertIndex(idx)+"_eta"),(monitoredDecays->ConvertIndex(idx)+" #eta"),50,-5,5));
67  idx++;
68  }
69 
70  return;
71 }
72 
74  return;
75 }
76 
77 void HiggsValidation::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup)
78 {
80  iSetup.getData( fPDGTable );
81  return;
82 }
83 void HiggsValidation::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
85 {
86  double weight = wmanager_.weight(iEvent);
87  nEvt->Fill(0.5,weight);
88 
89  //Gathering the HepMCProduct information
91  iEvent.getByToken(hepmcCollectionToken_, evt);
92 
93  //Get EVENT
94  HepMC::GenEvent *myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
95 
96  // loop over all particles
97  bool filled = false;
98  for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin();
99  iter!= myGenEvent->particles_end() && !filled; ++iter) {
100  if(particle_id == fabs((*iter)->pdg_id())){
101  std::vector<HepMC::GenParticle*> decayprod;
102  int channel = findHiggsDecayChannel(*iter,decayprod);
103  HiggsDecayChannels->Fill(channel,weight);
104  Higgs_pt->Fill((*iter)->momentum().perp(),weight);
105  Higgs_eta->Fill((*iter)->momentum().eta(),weight);
106  Higgs_mass->Fill((*iter)->momentum().m(),weight);
107  for(unsigned int i=0;i<decayprod.size();i++){
108  int idx=monitoredDecays->isDecayParticle(decayprod.at(i)->pdg_id());
109  if(0<=idx && idx<=(int)HiggsDecayProd_pt.size()){
110  HiggsDecayProd_pt.at(idx)->Fill(decayprod.at(i)->momentum().perp(),weight);
111  HiggsDecayProd_eta.at(idx)->Fill(decayprod.at(i)->momentum().eta(),weight);
112  }
113  }
114  filled = true;
115  }
116  }
117 
118  delete myGenEvent;
119 
120 }//analyze
121 
122 int HiggsValidation::findHiggsDecayChannel(const HepMC::GenParticle* genParticle,std::vector<HepMC::GenParticle*> &decayprod){
123  if(genParticle->status() == 1) return monitoredDecays->stable();
124  std::vector<int> children;
125  if ( genParticle->end_vertex() ) {
126  HepMC::GenVertex::particle_iterator des;
127  for(des = genParticle->end_vertex()->particles_begin(HepMC::descendants);
128  des!= genParticle->end_vertex()->particles_end(HepMC::descendants);++des ) {
129 
130  if((*des)->pdg_id() == genParticle->pdg_id()) continue;
131 
132  HepMC::GenVertex::particle_iterator mother = (*des)->production_vertex()->particles_begin(HepMC::parents);
133  if((*mother)->pdg_id() == genParticle->pdg_id()){
134  children.push_back((*des)->pdg_id());
135  decayprod.push_back((*des));
136  }
137  }
138  }
139 
140  if(children.size() == 2 && children.at(0) != 0 && children.at(1) != 0){
141  return monitoredDecays->position(children.at(0),children.at(1));
142  }
143  return monitoredDecays->undetermined();
144 }
145 
146 //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:872
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void endRun(const edm::Run &, const edm::EventSetup &)
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
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
MonitorElement * Higgs_pt
MonitorElement * HiggsDecayChannels
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
edm::InputTag hepmcCollection_
DQMStore * dbe
ME&#39;s &quot;container&quot;.
virtual void beginJob()
std::string ConvertIndex(int index)
WeightManager wmanager_
std::string channel(size_t i)
virtual void endJob()
dbl *** dir
Definition: mlp_gen.cc:35
int weight
Definition: histoStyle.py:50
std::string particle_name
double weight(const edm::Event &)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:584
Definition: Run.h:41
MonitorElement * Higgs_mass