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  monitoredDecays = new MonitoredDecays(iPSet);
28  hepmcCollectionToken_=consumes<HepMCProduct>(hepmcCollection_);
29 }
30 
32 
34  c.getData( fPDGTable );
35 }
36 
38 
40  TString dir="Generator/";
41  dir+=particle_name;
42  i.setCurrentFolder(dir.Data());
43 
44  // Number of analyzed events
45  nEvt = i.book1D("nEvt", "n analyzed Events", 1, 0., 1.);
46 
47  //decay type
48 
49  std::string channel = particle_name+"_DecayChannels";
50  HiggsDecayChannels = i.book1D(channel.c_str(),(particle_name+" decay channels").c_str(),monitoredDecays->size(),0,monitoredDecays->size());
51 
52  for(size_t j = 0; j < monitoredDecays->size(); ++j){
54  }
55 
56 
57  //Kinematics
58  Higgs_pt = i.book1D((particle_name+"_pt"),(particle_name+" p_{t}"),50,0,250);
59  Higgs_eta = i.book1D((particle_name+"_eta"),(particle_name+" #eta"),50,-5,5);
60  Higgs_mass = i.book1D((particle_name+"_m"),(particle_name+" M"),500,0,500);
61 
62  int idx=0;
63  for(unsigned int j=0;j<monitoredDecays->NDecayParticles();j++){
64  HiggsDecayProd_pt.push_back(i.book1D((monitoredDecays->ConvertIndex(idx)+"_pt"),(monitoredDecays->ConvertIndex(idx)+" p_{t}"),50,0,250));
65  HiggsDecayProd_eta.push_back(i.book1D((monitoredDecays->ConvertIndex(idx)+"_eta"),(monitoredDecays->ConvertIndex(idx)+" #eta"),50,-5,5));
66  idx++;
67  }
68 
69  return;
70 }
71 
73 {
74  double weight = wmanager_.weight(iEvent);
75  nEvt->Fill(0.5,weight);
76 
77  //Gathering the HepMCProduct information
79  iEvent.getByToken(hepmcCollectionToken_, evt);
80 
81  //Get EVENT
82  HepMC::GenEvent *myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
83 
84  // loop over all particles
85  bool filled = false;
86  for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin();
87  iter!= myGenEvent->particles_end() && !filled; ++iter) {
88  if(particle_id == fabs((*iter)->pdg_id())){
89  std::vector<HepMC::GenParticle*> decayprod;
90  int channel = findHiggsDecayChannel(*iter,decayprod);
91  HiggsDecayChannels->Fill(channel,weight);
92  Higgs_pt->Fill((*iter)->momentum().perp(),weight);
93  Higgs_eta->Fill((*iter)->momentum().eta(),weight);
94  Higgs_mass->Fill((*iter)->momentum().m(),weight);
95  for(unsigned int i=0;i<decayprod.size();i++){
96  int idx=monitoredDecays->isDecayParticle(decayprod.at(i)->pdg_id());
97  if(0<=idx && idx<=(int)HiggsDecayProd_pt.size()){
98  HiggsDecayProd_pt.at(idx)->Fill(decayprod.at(i)->momentum().perp(),weight);
99  HiggsDecayProd_eta.at(idx)->Fill(decayprod.at(i)->momentum().eta(),weight);
100  }
101  }
102  filled = true;
103  }
104  }
105 
106  delete myGenEvent;
107 
108 }//analyze
109 
110 int HiggsValidation::findHiggsDecayChannel(const HepMC::GenParticle* genParticle,std::vector<HepMC::GenParticle*> &decayprod){
111  if(genParticle->status() == 1) return monitoredDecays->stable();
112  std::vector<int> children;
113  if ( genParticle->end_vertex() ) {
114  HepMC::GenVertex::particle_iterator des;
115  for(des = genParticle->end_vertex()->particles_begin(HepMC::descendants);
116  des!= genParticle->end_vertex()->particles_end(HepMC::descendants);++des ) {
117 
118  if((*des)->pdg_id() == genParticle->pdg_id()) continue;
119 
120  HepMC::GenVertex::particle_iterator mother = (*des)->production_vertex()->particles_begin(HepMC::parents);
121  if((*mother)->pdg_id() == genParticle->pdg_id()){
122  children.push_back((*des)->pdg_id());
123  decayprod.push_back((*des));
124  }
125  }
126  }
127 
128  if(children.size() == 2 && children.at(0) != 0 && children.at(1) != 0){
129  return monitoredDecays->position(children.at(0),children.at(1));
130  }
131  return monitoredDecays->undetermined();
132 }
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
PDT table.
int i
Definition: DBlmapReader.cc:9
TPRegexp parents
Definition: eve_filter.cc:24
MonitorElement * Higgs_eta
virtual void dqmBeginRun(const edm::Run &r, const edm::EventSetup &c) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
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:230
int findHiggsDecayChannel(const HepMC::GenParticle *, std::vector< HepMC::GenParticle * > &decayprod)
virtual void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
std::vector< MonitorElement * > HiggsDecayProd_eta
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:113
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
int j
Definition: DBlmapReader.cc:9
MonitorElement * Higgs_pt
MonitorElement * HiggsDecayChannels
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
virtual void analyze(edm::Event const &, edm::EventSetup const &) override
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
edm::InputTag hepmcCollection_
std::string ConvertIndex(int index)
WeightManager wmanager_
std::string channel(size_t i)
dbl *** dir
Definition: mlp_gen.cc:35
int weight
Definition: histoStyle.py:50
std::string particle_name
double weight(const edm::Event &)
Definition: Run.h:41
MonitorElement * Higgs_mass