CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/DQM/Physics/src/QcdHighPtDQM.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2012/01/11 13:53:29 $
00005  *  $Revision: 1.5 $
00006  *  \author S. Bolognesi, Eric - CERN
00007  */
00008 
00009 #include "DQM/Physics/src/QcdHighPtDQM.h"
00010 
00011 #include "FWCore/Framework/interface/Event.h"
00012 #include "FWCore/Framework/interface/EventSetup.h"
00013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00014 #include "DataFormats/Common/interface/Handle.h"
00015 
00016 #include "DQMServices/Core/interface/DQMStore.h"
00017 #include "DQMServices/Core/interface/MonitorElement.h"
00018 #include "FWCore/ServiceRegistry/interface/Service.h"
00019 
00020 #include "DataFormats/METReco/interface/CaloMET.h"
00021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00022 #include <vector>
00023 
00024 #include <string>
00025 #include <cmath>
00026 
00027 using namespace std;
00028 using namespace edm;
00029 using namespace reco;
00030 using namespace math;
00031 
00032 
00033 //Get Jets and MET (no MET plots yet pending converging w/JetMET group)
00034 
00035 QcdHighPtDQM::QcdHighPtDQM(const ParameterSet& iConfig): 
00036   jetLabel_(iConfig.getUntrackedParameter<edm::InputTag>("jetTag")),
00037   metLabel1_(iConfig.getUntrackedParameter<edm::InputTag>("metTag1")),  
00038   metLabel2_(iConfig.getUntrackedParameter<edm::InputTag>("metTag2")),  
00039   metLabel3_(iConfig.getUntrackedParameter<edm::InputTag>("metTag3")),  
00040   metLabel4_(iConfig.getUntrackedParameter<edm::InputTag>("metTag4"))
00041 {
00042 
00043   theDbe = Service<DQMStore>().operator->();
00044   
00045 }
00046 
00047 QcdHighPtDQM::~QcdHighPtDQM() { 
00048   
00049 
00050 }
00051 
00052 
00053 void QcdHighPtDQM::beginJob() {
00054  
00055 
00056   //Book MEs
00057  
00058   theDbe->setCurrentFolder("Physics/QcdHighPt");  
00059 
00060   MEcontainer_["dijet_mass"] = theDbe->book1D("dijet_mass", "dijet resonance invariant mass, barrel region", 100, 0, 1000);
00061   MEcontainer_["njets"] = theDbe->book1D("njets", "jet multiplicity", 10, 0, 10);
00062   MEcontainer_["etaphi"] = theDbe->book2D("etaphi", "eta/phi distribution", 83, -42, 42, 72, -M_PI, M_PI);
00063   MEcontainer_["njets30"] = theDbe->book1D("njets30", "jet multiplicity, pt > 30 GeV", 10, 0, 10);
00064 
00065   //book histograms for inclusive jet quantities
00066   MEcontainer_["inclusive_jet_pt"] = theDbe->book1D("inclusive_jet_pt", "inclusive jet Pt spectrum", 100, 0, 1000);
00067   MEcontainer_["inclusive_jet_pt_barrel"] = theDbe->book1D("inclusive_jet_pt_barrel", "inclusive jet Pt, eta < 1.3", 100, 0, 1000);
00068   MEcontainer_["inclusive_jet_pt_forward"] = theDbe->book1D("inclusive_jet_pt_forward", "inclusive jet Pt, 3.0 < eta < 5.0", 100, 0, 1000);
00069   MEcontainer_["inclusive_jet_pt_endcap"] = theDbe->book1D("inclusive_jet_pt_endcap", "inclusive jet Pt, 1.3 < eta < 3.0", 100, 0, 1000);
00070 
00071   //book histograms for leading jet quantities
00072   MEcontainer_["leading_jet_pt"] = theDbe->book1D("leading_jet_pt", "leading jet Pt", 100, 0, 1000);
00073   MEcontainer_["leading_jet_pt_barrel"] = theDbe->book1D("leading_jet_pt_barrel", "leading jet Pt, eta < 1.3", 100, 0, 1000);
00074   MEcontainer_["leading_jet_pt_forward"] = theDbe->book1D("leading_jet_pt_forward", "leading jet Pt, 3.0 < eta < 5.0", 100, 0, 1000);
00075   MEcontainer_["leading_jet_pt_endcap"] = theDbe->book1D("leading_jet_pt_endcap", "leading jet Pt, 1.3 < eta < 3.0", 100, 0, 1000);
00076 
00077   //book histograms for met over sum et and met over leading jet pt for various
00078   //flavors of MET
00079   MEcontainer_["movers_met"] = theDbe->book1D("movers_met", "MET over Sum ET for basic MET collection", 50, 0, 1);
00080   MEcontainer_["moverl_met"] = theDbe->book1D("moverl_met", "MET over leading jet Pt for basic MET collection", 50, 0, 2);
00081   
00082   MEcontainer_["movers_metho"] = theDbe->book1D("movers_metho", "MET over Sum ET for MET HO collection", 50, 0, 1);
00083   MEcontainer_["moverl_metho"] = theDbe->book1D("moverl_metho", "MET over leading jet Pt for MET HO collection", 50, 0, 2);
00084   
00085   MEcontainer_["movers_metnohf"] = theDbe->book1D("movers_metnohf", "MET over Sum ET for MET no HF collection", 50, 0, 1);
00086   MEcontainer_["moverl_metnohf"] = theDbe->book1D("moverl_metnohf", "MET over leading jet Pt for MET no HF collection", 50, 0, 2);
00087 
00088   MEcontainer_["movers_metnohfho"] = theDbe->book1D("movers_metnohfho", "MET over Sum ET for MET no HF HO collection", 50, 0, 1);
00089   MEcontainer_["moverl_metnohfho"] = theDbe->book1D("moverl_metnohfho", "MET over leading jet Pt for MET no HF HO collection", 50, 0, 2);
00090 
00091   //book histograms for EMF fraction for all jets and first 3 jets
00092   MEcontainer_["inclusive_jet_EMF"] = theDbe->book1D("inclusive_jet_EMF", "inclusive jet EMF", 50, -1, 1);
00093   MEcontainer_["leading_jet_EMF"] = theDbe->book1D("leading_jet_EMF", "leading jet EMF", 50, -1, 1);
00094   MEcontainer_["second_jet_EMF"] = theDbe->book1D("second_jet_EMF", "second jet EMF", 50, -1, 1);
00095   MEcontainer_["third_jet_EMF"] = theDbe->book1D("third_jet_EMF", "third jet EMF", 50, -1, 1);
00096   
00097 
00098 
00099 
00100 }
00101 
00102 void QcdHighPtDQM::endJob(void) {
00103 
00104 }
00105 
00106 // method to calculate MET over Sum ET from a particular MET collection
00107 float QcdHighPtDQM::movers(const CaloMETCollection &metcollection) {
00108 
00109   float metovers = 0;
00110   CaloMETCollection::const_iterator met_iter;
00111   for(met_iter = metcollection.begin(); met_iter!= metcollection.end(); ++met_iter)
00112     {
00113       float mex = met_iter->momentum().x();
00114       float mey = met_iter->momentum().y();
00115       float met = sqrt(mex*mex+mey*mey);
00116       float sumet = met_iter->sumEt();
00117       metovers = (met / sumet);
00118     }
00119   return metovers;
00120 }
00121 
00122 
00123 // method to calculate MET over Leading jet PT for a particular MET collection
00124 float QcdHighPtDQM::moverl(const CaloMETCollection &metcollection, float &ljpt) {
00125  float metoverl = 0;
00126  CaloMETCollection::const_iterator met_iter;
00127  for(met_iter = metcollection.begin(); met_iter!= metcollection.end(); ++met_iter)
00128    {
00129      float mex = met_iter->momentum().x();
00130      float mey = met_iter->momentum().y();
00131      float met = sqrt(mex*mex+mey*mey);
00132      metoverl = (met / ljpt);
00133    }
00134  return metoverl;
00135 }
00136 
00137 
00138 void QcdHighPtDQM::analyze(const Event& iEvent, const EventSetup& iSetup) {
00139   
00140   //Get Jets
00141   edm::Handle<CaloJetCollection> jetHandle;
00142   iEvent.getByLabel(jetLabel_,jetHandle);
00143   const CaloJetCollection & jets = *jetHandle;
00144   CaloJetCollection::const_iterator jet_iter;
00145 
00146   //Get MET collections
00147   edm::Handle<CaloMETCollection> metHandle;
00148   iEvent.getByLabel(metLabel1_, metHandle);
00149   const CaloMETCollection &met = *metHandle;
00150 
00151   edm::Handle<CaloMETCollection> metHOHandle;
00152   iEvent.getByLabel(metLabel2_, metHOHandle);
00153   const CaloMETCollection &metHO = *metHOHandle;
00154   
00155   edm::Handle<CaloMETCollection> metNoHFHandle;
00156   iEvent.getByLabel(metLabel3_, metNoHFHandle);
00157   const CaloMETCollection &metNoHF = *metNoHFHandle;
00158 
00159   edm::Handle<CaloMETCollection> metNoHFHOHandle;
00160   iEvent.getByLabel(metLabel4_, metNoHFHOHandle);
00161   const CaloMETCollection &metNoHFHO = *metNoHFHOHandle;
00162 
00163   //initialize leading jet value and jet multiplicity counter
00164   int njets = 0;
00165   int njets30 = 0;
00166   float leading_jetpt = 0;
00167   float leading_jeteta = 0;
00168 
00169   //initialize variables for picking out leading 2 barrel jets
00170   reco::CaloJet leadingbarreljet; 
00171   reco::CaloJet secondbarreljet;
00172   int nbarreljets = 0;
00173 
00174   //get bins in eta.  
00175   //Bins correspond to calotower regions.
00176 
00177   const float etabins[83] = {-5.191, -4.889, -4.716, -4.538, -4.363, -4.191, -4.013, -3.839, -3.664, -3.489, -3.314, -3.139, -2.964, -2.853, -2.650, -2.500, -2.322, -2.172, -2.043, -1.930, -1.830, -1.740, -1.653, -1.566, -1.479, -1.392, -1.305, -1.218, -1.131, -1.044, -.957, -.879, -.783, -.696, -.609, -.522, -.435, -.348, -.261, -.174, -.087, 0, .087, .174, .261, .348, .435, .522, .609, .696, .783, .879, .957, 1.044, 1.131, 1.218, 1.305, 1.392, 1.479, 1.566, 1.653, 1.740, 1.830, 1.930, 2.043, 2.172, 2.322, 2.500, 2.650, 2.853, 2.964, 3.139, 3.314, 3.489, 3.664, 3.839, 4.013, 4.191, 4.363, 4.538, 4.889, 5.191};
00178 
00179   for(jet_iter = jets.begin(); jet_iter!= jets.end(); ++jet_iter){
00180     njets++;
00181   
00182     //get Jet stats
00183     float jet_pt = jet_iter->pt();
00184     float jet_eta = jet_iter->eta();
00185     float jet_phi = jet_iter->phi();
00186 
00187     //fill jet Pt and jet EMF
00188     MEcontainer_["inclusive_jet_pt"]->Fill(jet_pt);
00189     MEcontainer_["inclusive_jet_EMF"]->Fill(jet_iter->emEnergyFraction());
00190 
00191     // pick out up to the first 2 leading barrel jets
00192     // for use in calculating dijet mass in barrel region
00193     // also fill jet Pt histogram for barrel
00194 
00195     if(jet_eta <= 1.3){
00196       MEcontainer_["inclusive_jet_pt_barrel"]->Fill(jet_pt);
00197       if(nbarreljets == 0){
00198         leadingbarreljet = jets[(njets-1)];
00199         nbarreljets++;
00200       }
00201       else if(nbarreljets == 1){
00202         secondbarreljet = jets[(njets-1)];
00203         nbarreljets++;
00204       }
00205 
00206     }
00207 
00208     // fill jet Pt for endcap and forward regions
00209     else if(jet_eta <= 3.0 && jet_eta > 1.3){MEcontainer_["inclusive_jet_pt_endcap"]->Fill(jet_pt);}
00210 
00211     else if(jet_eta <= 5.0 && jet_eta > 3.0){MEcontainer_["inclusive_jet_pt_forward"]->Fill(jet_pt);}
00212 
00213 
00214     // count jet multiplicity for jets with Pt > 30
00215     if((jet_pt) > 30) njets30++;
00216 
00217     // check leading jet quantities
00218     if(jet_pt > leading_jetpt){
00219       leading_jetpt = jet_pt;
00220       leading_jeteta = jet_eta;
00221     }
00222 
00223     //fill eta-phi plot
00224     for (int eit = 0; eit < 81; eit++){
00225       for(int pit = 0; pit < 72; pit++){
00226         float low_eta = etabins[eit];
00227         float high_eta = etabins[eit+1];
00228         float low_phi = (-M_PI) + pit*(M_PI/36);
00229         float high_phi = low_phi + (M_PI/36);
00230         if(jet_eta > low_eta && jet_eta < high_eta && jet_phi > low_phi && jet_phi < high_phi){
00231           MEcontainer_["etaphi"]->Fill((eit - 41), jet_phi);}
00232 
00233 
00234       }
00235     }
00236   }
00237   
00238   // after iterating over all jets, fill leading jet quantity histograms 
00239   // and jet multiplicity histograms
00240 
00241   MEcontainer_["leading_jet_pt"]->Fill(leading_jetpt);
00242   
00243   if(leading_jeteta <= 1.3){MEcontainer_["leading_jet_pt_barrel"]->Fill(leading_jetpt);}
00244 
00245   else if(leading_jeteta <= 3.0 && leading_jeteta > 1.3){MEcontainer_["leading_jet_pt_endcap"]->Fill(leading_jetpt);}
00246 
00247   else if(leading_jeteta <= 5.0 && leading_jeteta > 3.0){MEcontainer_["leading_jet_pt_forward"]->Fill(leading_jetpt);}
00248 
00249   MEcontainer_["njets"]->Fill(njets);
00250 
00251   MEcontainer_["njets30"]->Fill(njets30);
00252 
00253   // fill MET over Sum ET and Leading jet PT for all MET flavors
00254   MEcontainer_["movers_met"]->Fill(movers(met));
00255   MEcontainer_["moverl_met"]->Fill(movers(met), leading_jetpt);
00256   MEcontainer_["movers_metho"]->Fill(movers(metHO));
00257   MEcontainer_["moverl_metho"]->Fill(movers(metHO), leading_jetpt);
00258   MEcontainer_["movers_metnohf"]->Fill(movers(metNoHF));
00259   MEcontainer_["moverl_metnohf"]->Fill(movers(metNoHF), leading_jetpt);
00260   MEcontainer_["movers_metnohfho"]->Fill(movers(metNoHFHO));
00261   MEcontainer_["moverl_metnohfho"]->Fill(movers(metNoHFHO), leading_jetpt);
00262   
00263 
00264   // fetch first 3 jet EMF
00265 
00266   if(jets.size() >= 1){
00267     MEcontainer_["leading_jet_EMF"]->Fill(jets[0].emEnergyFraction());
00268     if(jets.size() >=2){
00269       MEcontainer_["second_jet_EMF"]->Fill(jets[1].emEnergyFraction());    
00270       if(jets.size() >= 3){
00271         MEcontainer_["third_jet_EMF"]->Fill(jets[2].emEnergyFraction());
00272       }
00273     }
00274   }
00275 
00276   // if 2 nontrivial barrel jets, reconstruct dijet mass
00277 
00278   if(nbarreljets == 2){
00279     if(leadingbarreljet.energy() > 0 && secondbarreljet.energy() > 0){
00280       math::XYZTLorentzVector DiJet = leadingbarreljet.p4() + secondbarreljet.p4();
00281       float dijet_mass = DiJet.mass();
00282         MEcontainer_["dijet_mass"]->Fill(dijet_mass);
00283     }
00284   }
00285       
00286     
00287 }
00288 
00289 
00290 
00291 
00292 
00293 
00294 
00295 
00296 
00297 
00298 
00299 
00300