CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/RecoJets/JetAnalyzers/src/CMSDAS11DijetAnalyzer.cc

Go to the documentation of this file.
00001 // CMSDAS11DijetAnalyzer.cc
00002 // Description: A basic dijet analyzer for the CMSDAS 2011
00003 // Author: John Paul Chou
00004 // Date: January 12, 2011
00005 
00006 #include "RecoJets/JetAnalyzers/interface/CMSDAS11DijetAnalyzer.h"
00007 
00008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00009 #include "FWCore/Framework/interface/MakerMacros.h"
00010 #include "FWCore/ServiceRegistry/interface/Service.h"
00011 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00012 
00013 #include "DataFormats/VertexReco/interface/Vertex.h"
00014 #include "JetMETCorrections/Objects/interface/JetCorrector.h"
00015 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
00016 #include "SimDataFormats/GeneratorProducts/interface/GenRunInfoProduct.h"
00017 
00018 #include <TH1D.h>
00019 
00020 CMSDAS11DijetAnalyzer::CMSDAS11DijetAnalyzer(edm::ParameterSet const& params) :
00021   edm::EDAnalyzer(),
00022   jetSrc(params.getParameter<edm::InputTag>("jetSrc")),
00023   vertexSrc(params.getParameter<edm::InputTag>("vertexSrc")),
00024   jetCorrections(params.getParameter<std::string>("jetCorrections")),
00025   innerDeltaEta(params.getParameter<double>("innerDeltaEta")),
00026   outerDeltaEta(params.getParameter<double>("outerDeltaEta")),
00027   JESbias(params.getParameter<double>("JESbias"))
00028 {
00029   // setup file service
00030   edm::Service<TFileService> fs;
00031 
00032   const int NBINS=36;
00033   Double_t BOUNDARIES[NBINS] = { 220, 244, 270, 296, 325, 354, 386, 419, 453,
00034                                  489, 526, 565, 606, 649, 693, 740, 788, 838,
00035                                  890, 944, 1000, 1058, 1118, 1181, 1246, 1313, 1383,
00036                                  1455, 1530, 1607, 1687, 1770, 1856, 1945, 2037, 2132 };
00037   
00038   // setup histograms
00039   hVertexZ = fs->make<TH1D>("hVertexZ", "Z position of the Vertex",50,-20,20);
00040   hJetRawPt    = fs->make<TH1D>("hJetRawPt","Raw Jet Pt",50,0,1000);
00041   hJetCorrPt = fs->make<TH1D>("hJetCorrPt","Corrected Jet Pt",50,0,1000);
00042   hJet1Pt      = fs->make<TH1D>("hJet1Pt","Corrected Jet1 Pt",50,0,1000);
00043   hJet2Pt      = fs->make<TH1D>("hJet2Pt","Corrected Jet2 Pt",50,0,1000);
00044 
00045   hJetEta      = fs->make<TH1D>("hJetEta","Corrected Jet Eta",   10,-5,5);
00046   hJet1Eta      = fs->make<TH1D>("hJet1Eta","Corrected Jet1 Eta",10,-5,5);
00047   hJet2Eta      = fs->make<TH1D>("hJet2Eta","Corrected Jet2 Eta",10,-5,5);
00048 
00049   hJetPhi      = fs->make<TH1D>("hJetPhi","Corrected Jet Phi",   10,-3.1415,3.1415);
00050   hJet1Phi      = fs->make<TH1D>("hJet1Phi","Corrected Jet1 Phi",10,-3.1415,3.1415);
00051   hJet2Phi      = fs->make<TH1D>("hJet2Phi","Corrected Jet2 Phi",10,-3.1415,3.1415);
00052 
00053   hJetEMF       = fs->make<TH1D>("hJetEMF","EM Fraction of Jets",50,0,1);
00054   hJet1EMF      = fs->make<TH1D>("hJet1EMF","EM Fraction of Jet1",50,0,1);
00055   hJet2EMF      = fs->make<TH1D>("hJet2EMF","EM Fraction of Jet2",50,0,1);
00056 
00057   hCorDijetMass = fs->make<TH1D>("hCorDijetMass","Corrected Dijet Mass",NBINS-1,BOUNDARIES);
00058   hDijetDeltaPhi= fs->make<TH1D>("hDijetDeltaPhi","Dijet |#Delta #phi|",50,0,3.1415);
00059   hDijetDeltaEta= fs->make<TH1D>("hDijetDeltaEta","Dijet |#Delta #eta|",50,0,6);
00060 
00061   hInnerDijetMass = fs->make<TH1D>("hInnerDijetMass","Corrected Inner Dijet Mass",NBINS-1,BOUNDARIES);
00062   hOuterDijetMass = fs->make<TH1D>("hOuterDijetMass","Corrected Outer Dijet Mass",NBINS-1,BOUNDARIES);
00063 }
00064 
00065 void CMSDAS11DijetAnalyzer::endJob(void){
00066 }
00067 
00068 
00069 void CMSDAS11DijetAnalyzer::analyze( const edm::Event& iEvent, const edm::EventSetup& iSetup)
00070 {
00071 
00073   double mWeight;
00074   edm::Handle<GenEventInfoProduct> hEventInfo;
00075   iEvent.getByLabel("generator", hEventInfo);
00076   if (hEventInfo.isValid())
00077     mWeight = hEventInfo->weight();
00078   else 
00079     mWeight = 1.;
00080 
00082   // Get event ID information
00084   //int nrun=iEvent.id().run();
00085   //int nlumi=iEvent.luminosityBlock();
00086   //int nevent=iEvent.id().event();
00087 
00089   // Get Primary Vertex Information
00091 
00092   // magic to get the vertices from EDM
00093   edm::Handle< std::vector<reco::Vertex> > vertices_h;
00094   iEvent.getByLabel(vertexSrc, vertices_h);
00095   if (!vertices_h.isValid()) {
00096     std::cout<<"Didja hear the one about the empty vertex collection?\n";
00097     return;
00098   }
00099   
00100   // require in the event that there is at least one reconstructed vertex
00101   if(vertices_h->size()<=0) return;
00102 
00103   // pick the first (i.e. highest sum pt) verte
00104   const reco::Vertex* theVertex=&(vertices_h->front());
00105 
00106   // require that the vertex meets certain criteria
00107   if(theVertex->ndof()<5) return;
00108   if(fabs(theVertex->z())>24.0) return;
00109   if(fabs(theVertex->position().rho())>2.0) return;
00110 
00112   // Get Jet Information
00114 
00115   // magic to get the jets from EDM
00116   edm::Handle<reco::CaloJetCollection> jets_h;
00117   iEvent.getByLabel(jetSrc, jets_h);
00118   if (!jets_h.isValid()) {
00119     std::cout<<"Didja hear the one about the empty jet collection?\n";
00120     return;
00121   }
00122 
00123   // magic to get the jet energy corrections
00124   const JetCorrector* corrector = JetCorrector::getJetCorrector(jetCorrections,iSetup);
00125 
00126   // collection of selected jets
00127   std::vector<reco::CaloJet> selectedJets;
00128 
00129   // loop over the jet collection
00130   for(reco::CaloJetCollection::const_iterator j_it = jets_h->begin(); j_it!=jets_h->end(); j_it++) {
00131     reco::CaloJet jet = *j_it;
00132 
00133     // calculate and apply the correction
00134     double scale = corrector->correction(jet.p4());
00135 
00136     // Introduce a purposeful bias to the correction, to show what happens
00137     scale *= JESbias;
00138 
00139     // fill the histograms
00140     hJetRawPt->Fill(jet.pt(),mWeight);
00141     hJetEta->Fill(jet.eta(),mWeight);
00142     hJetPhi->Fill(jet.phi(),mWeight);
00143     hJetEMF->Fill(jet.emEnergyFraction(),mWeight);
00144 
00145     // correct the jet energy
00146     jet.scaleEnergy(scale);
00147     
00148     // now fill the correct jet pt after the correction
00149     hJetCorrPt->Fill(jet.pt(),mWeight);
00150 
00151     // put the selected jets into a collection
00152     selectedJets.push_back(jet);
00153   }
00154 
00155   // require at least two jets to continue
00156   if(selectedJets.size()<2) return;
00157 
00158   //sort by corrected pt (not the same order as raw pt, sometimes)
00159   sort(selectedJets.begin(), selectedJets.end(), compare_JetPt);
00160 
00161   // select high pt, central, non-noise-like jets
00162   if (selectedJets[0].pt()<50.0) return;
00163   if (fabs(selectedJets[0].eta())>2.5) return;
00164   if (selectedJets[0].emEnergyFraction()<0.01) return;
00165   if (selectedJets[1].pt()<50.0) return;
00166   if (fabs(selectedJets[1].eta())>2.5) return;
00167   if (selectedJets[1].emEnergyFraction()<0.01) return;
00168   
00169 
00170   // fill histograms for the jets in our dijets, only
00171   hJet1Pt ->Fill(selectedJets[0].pt(),mWeight);
00172   hJet1Eta->Fill(selectedJets[0].eta(),mWeight);
00173   hJet1Phi->Fill(selectedJets[0].phi(),mWeight);
00174   hJet1EMF->Fill(selectedJets[0].emEnergyFraction(),mWeight);
00175   hJet2Pt ->Fill(selectedJets[1].pt(),mWeight);
00176   hJet2Eta->Fill(selectedJets[1].eta(),mWeight);
00177   hJet2Phi->Fill(selectedJets[1].phi(),mWeight);
00178   hJet2EMF->Fill(selectedJets[1].emEnergyFraction(),mWeight);
00179 
00180   //Get the mass of the two leading jets.  Needs their 4-vectors...
00181   double corMass = (selectedJets[0].p4()+selectedJets[1].p4()).M();
00182   double deltaEta = fabs(selectedJets[0].eta()-selectedJets[1].eta());
00183   if (corMass < 489) return;
00184   if (deltaEta > 1.3) return;
00185   hCorDijetMass->Fill(corMass,mWeight);
00186   hDijetDeltaPhi->Fill(fabs(selectedJets[0].phi()-selectedJets[1].phi()) ,mWeight);
00187   hDijetDeltaEta->Fill(deltaEta ,mWeight);
00188 
00189   //Fill the inner and outer dijet mass spectra, to make the ratio from
00190   if (deltaEta < innerDeltaEta) hInnerDijetMass->Fill(corMass,mWeight);
00191   else if (deltaEta < outerDeltaEta) hOuterDijetMass->Fill(corMass,mWeight);
00192 
00193   hVertexZ->Fill(theVertex->z(),mWeight);
00194   return;
00195 }
00196 
00197 
00198 DEFINE_FWK_MODULE(CMSDAS11DijetAnalyzer);