CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/Validation/RecoVertex/src/MCvsRecoVerticesAnalyzer.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    MCvsRecoVerticesAnalyzer
00004 // Class:      MCvsRecoVerticesAnalyzer
00005 // 
00013 //
00014 // Original Author:  Andrea Venturi
00015 //         Created:  Thu Dec 16 16:32:56 CEST 2010
00016 // $Id: MCvsRecoVerticesAnalyzer.cc,v 1.4 2011/11/26 00:51:42 venturia Exp $
00017 //
00018 //
00019 
00020 
00021 // system include files
00022 #include <memory>
00023 #include <numeric>
00024 
00025 #include <vector>
00026 
00027 // user include files
00028 #include "FWCore/Framework/interface/Frameworkfwd.h"
00029 #include "FWCore/Framework/interface/EDAnalyzer.h"
00030 
00031 #include "FWCore/Framework/interface/Event.h"
00032 #include "FWCore/Framework/interface/Run.h"
00033 #include "FWCore/Framework/interface/MakerMacros.h"
00034 
00035 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00036 
00037 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00038 
00039 #include "FWCore/ServiceRegistry/interface/Service.h"
00040 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00041 
00042 #include "FWCore/Utilities/interface/InputTag.h"
00043 
00044 #include "SimDataFormats/PileupSummaryInfo/interface/PileupSummaryInfo.h"
00045 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00046 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
00047 
00048 #include "DataFormats/VertexReco/interface/Vertex.h"
00049 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00050 
00051 #include "TH1F.h"
00052 #include "TH2F.h"
00053 #include "TProfile.h"
00054 
00055 //
00056 // class decleration
00057 //
00058 
00059 
00060 class MCvsRecoVerticesAnalyzer : public edm::EDAnalyzer {
00061 public:
00062   explicit MCvsRecoVerticesAnalyzer(const edm::ParameterSet&);
00063   ~MCvsRecoVerticesAnalyzer();
00064   
00065 private:
00066   virtual void beginJob() ;
00067   virtual void beginRun(const edm::Run&, const edm::EventSetup&);
00068   virtual void endRun(const edm::Run&, const edm::EventSetup&);
00069   virtual void analyze(const edm::Event&, const edm::EventSetup&);
00070   virtual void endJob() ;
00071   
00072       // ----------member data ---------------------------
00073 
00074   
00075 
00076   edm::InputTag m_pileupcollection;
00077   edm::InputTag m_mctruthcollection;
00078   edm::InputTag m_pvcollection;
00079   const bool m_useweight;
00080   edm::InputTag m_weight;
00081   const bool m_useVisibleVertices;
00082   const edm::ParameterSet m_histoParameters;
00083 
00084   TH2F* m_hrecovsmcnvtx2d;
00085   TProfile* m_hrecovsmcnvtxprof;
00086   TProfile* m_hrecovsmcnvtxweightedprof;
00087   TH2F* m_hrecovsmclumi2d;
00088   TProfile* m_hrecovsmclumiprof;
00089   TProfile* m_hrecovsmclumiweightedprof;
00090   TH1F* m_hdeltazfirst;
00091   TH1F* m_hdeltazclose;
00092   TH1F* m_hclosestvtx;
00093  
00094   TH2F* m_hdeltazfirstvsnpu;
00095   TH2F* m_hdeltazclosevsnpu;
00096   TH2F* m_hclosestvtxvsnpu;
00097  
00098 
00099 };
00100 
00101 //
00102 // constants, enums and typedefs
00103 //
00104 
00105 //
00106 // static data member definitions
00107 //
00108 
00109 //
00110 // constructors and destructor
00111 //
00112 MCvsRecoVerticesAnalyzer::MCvsRecoVerticesAnalyzer(const edm::ParameterSet& iConfig):
00113   m_pileupcollection(iConfig.getParameter<edm::InputTag>("pileupSummaryCollection")),
00114   m_mctruthcollection(iConfig.getParameter<edm::InputTag>("mcTruthCollection")),
00115   m_pvcollection(iConfig.getParameter<edm::InputTag>("pvCollection")),
00116   m_useweight(iConfig.getParameter<bool>("useWeight")),
00117   m_weight(iConfig.getParameter<edm::InputTag>("weightProduct")),
00118   m_useVisibleVertices(iConfig.getParameter<bool>("useVisibleVertices")),
00119   m_histoParameters(iConfig.getUntrackedParameter<edm::ParameterSet>("histoParameters",edm::ParameterSet()))
00120 {
00121    //now do what ever initialization is needed
00122 
00123   if(m_useVisibleVertices) edm::LogInfo("UseVisibleVertices") << "Only visible vertices will be used to compute Npileup";
00124 
00125   edm::Service<TFileService> tfserv;
00126 
00127   m_hrecovsmcnvtx2d = tfserv->make<TH2F>("recovsmcnvtx2d","Number of reco vertices vs pileup interactions",60,-0.5,59.5,60,-0.5,59.5);
00128   m_hrecovsmcnvtx2d->GetXaxis()->SetTitle("Pileup Interactions");  m_hrecovsmcnvtx2d->GetYaxis()->SetTitle("Reco Vertices");
00129   m_hrecovsmcnvtxprof = tfserv->make<TProfile>("recovsmcnvtxprof","Mean number of reco vs pileup vertices",60,-0.5,59.5);
00130   m_hrecovsmcnvtxprof->GetXaxis()->SetTitle("Pileup Interactions");  m_hrecovsmcnvtxprof->GetYaxis()->SetTitle("Reco Vertices");
00131 
00132   m_hrecovsmclumi2d = tfserv->make<TH2F>("recovsmclumi2d","Number of reco vertices vs ave pileup interactions",200,0.,50.,60,-0.5,59.5);
00133   m_hrecovsmclumi2d->GetXaxis()->SetTitle("Average Pileup Interactions");  m_hrecovsmclumi2d->GetYaxis()->SetTitle("Reco Vertices");
00134   m_hrecovsmclumiprof = tfserv->make<TProfile>("recovsmclumiprof","Mean number of reco vs ave pileup vertices",200,0.,50.);
00135   m_hrecovsmclumiprof->GetXaxis()->SetTitle("Average Pileup Interactions");  m_hrecovsmclumiprof->GetYaxis()->SetTitle("Reco Vertices");
00136 
00137   if(m_useweight) {
00138     m_hrecovsmcnvtxweightedprof = tfserv->make<TProfile>("recovsmcnvtxweightedprof","Mean number of reco vs pileup vertices (1-w) weight",60,-0.5,59.5);
00139     m_hrecovsmcnvtxweightedprof->GetXaxis()->SetTitle("Pileup Interactions");  m_hrecovsmcnvtxweightedprof->GetYaxis()->SetTitle("Reco Vertices (1-w)");
00140 
00141     m_hrecovsmclumiweightedprof = tfserv->make<TProfile>("recovsmclumiweightedprof","Mean number of reco vs ave pileup vertices (1-w) weight",
00142                                                          200,0.,50.);
00143     m_hrecovsmclumiweightedprof->GetXaxis()->SetTitle("Average Pileup Interactions");  
00144     m_hrecovsmclumiweightedprof->GetYaxis()->SetTitle("Reco Vertices (1-w)");
00145   }
00146 
00147   m_hdeltazfirst = tfserv->make<TH1F>("deltazfirst","Reco-MC vertex z position (first vertex)",
00148                                       m_histoParameters.getUntrackedParameter<unsigned int>("zBins",1000),
00149                                       m_histoParameters.getUntrackedParameter<double>("zMin",-1.),
00150                                       m_histoParameters.getUntrackedParameter<double>("zMax",1.));
00151   m_hdeltazfirst->GetXaxis()->SetTitle("#Delta z (cm)");  m_hdeltazfirst->GetYaxis()->SetTitle("Events");
00152 
00153   m_hdeltazclose = tfserv->make<TH1F>("deltazclose","Reco-MC vertex z position (closest vertex)",
00154                                       m_histoParameters.getUntrackedParameter<unsigned int>("zBins",1000),
00155                                       m_histoParameters.getUntrackedParameter<double>("zMin",-1.),
00156                                       m_histoParameters.getUntrackedParameter<double>("zMax",1.));
00157   m_hdeltazclose->GetXaxis()->SetTitle("#Delta z (cm)");  m_hdeltazclose->GetYaxis()->SetTitle("Events");
00158 
00159   m_hclosestvtx = tfserv->make<TH1F>("closestvtx","Closest reco vtx ID",30,-0.5,29.5);
00160   m_hclosestvtx->GetXaxis()->SetTitle("Vtx ID");  m_hclosestvtx->GetYaxis()->SetTitle("Events");
00161 
00162   m_hdeltazfirstvsnpu = tfserv->make<TH2F>("deltazfirstvsnpu","Reco-MC vertex z position (first vertex) vs Npileup",30,-0.5,29.5,
00163                                       m_histoParameters.getUntrackedParameter<unsigned int>("zBins",1000),
00164                                       m_histoParameters.getUntrackedParameter<double>("zMin",-1.),
00165                                       m_histoParameters.getUntrackedParameter<double>("zMax",1.));
00166   m_hdeltazfirstvsnpu->GetXaxis()->SetTitle("pileup Interactions");  m_hdeltazfirstvsnpu->GetYaxis()->SetTitle("#Delta z (cm)");  
00167 
00168   m_hdeltazclosevsnpu = tfserv->make<TH2F>("deltazclosevsnpu","Reco-MC vertex z position (closest vertex) v Npileup",30,-0.5,29.5,
00169                                       m_histoParameters.getUntrackedParameter<unsigned int>("zBins",1000),
00170                                       m_histoParameters.getUntrackedParameter<double>("zMin",-1.),
00171                                       m_histoParameters.getUntrackedParameter<double>("zMax",1.));
00172   m_hdeltazclosevsnpu->GetXaxis()->SetTitle("Pileup Interactions");  m_hdeltazclosevsnpu->GetYaxis()->SetTitle("#Delta z (cm)");
00173 
00174   m_hclosestvtxvsnpu = tfserv->make<TH2F>("closestvtxvsnpu","Closest reco vtx ID vs Npileup",30,-0.5,29.5,30,-0.5,29.5);
00175   m_hclosestvtxvsnpu->GetXaxis()->SetTitle("Pileup Interactions");  m_hclosestvtxvsnpu->GetYaxis()->SetTitle("Vtx ID");
00176 
00177 }
00178 
00179 
00180 MCvsRecoVerticesAnalyzer::~MCvsRecoVerticesAnalyzer()
00181 {
00182 
00183    // do anything here that needs to be done at desctruction time
00184    // (e.g. close files, deallocate resources etc.)
00185 
00186 }
00187 
00188 
00189 //
00190 // member functions
00191 //
00192 
00193 // ------------ method called to for each event  ------------
00194 void
00195 MCvsRecoVerticesAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00196 {
00197   using namespace edm;
00198   
00199   double weight = 1.;
00200   
00201   if(m_useweight) {
00202     Handle<double> weightprod;
00203     iEvent.getByLabel(m_weight,weightprod);
00204     
00205     weight = *weightprod;
00206     
00207   }
00208   
00209   Handle<std::vector<PileupSummaryInfo> > pileupinfos;
00210   iEvent.getByLabel(m_pileupcollection,pileupinfos);
00211 
00212   // look for the intime PileupSummaryInfo
00213 
00214   std::vector<PileupSummaryInfo>::const_iterator pileupinfo;
00215 
00216   for(pileupinfo = pileupinfos->begin(); pileupinfo != pileupinfos->end() ; ++pileupinfo) {
00217 
00218     if(pileupinfo->getBunchCrossing()==0) break;
00219 
00220   } 
00221   
00222   //
00223   
00224   Handle<reco::VertexCollection> pvcoll;
00225   iEvent.getByLabel(m_pvcollection,pvcoll);
00226   
00227 
00228    //
00229 
00230   if(pileupinfo->getBunchCrossing()!=0) {
00231 
00232     edm::LogError("NoInTimePileUpInfo") << "Cannot find the in-time pileup info " << pileupinfo->getBunchCrossing();
00233 
00234   }
00235   else {
00236 
00237     int npileup = pileupinfo->getPU_NumInteractions();
00238 
00239     if(m_useVisibleVertices) npileup = pileupinfo->getPU_zpositions().size();
00240 
00241     m_hrecovsmcnvtx2d->Fill(npileup,pvcoll->size(),weight);
00242     m_hrecovsmcnvtxprof->Fill(npileup,pvcoll->size(),weight);
00243 
00244     m_hrecovsmclumi2d->Fill(pileupinfo->getTrueNumInteractions(),pvcoll->size(),weight);
00245     m_hrecovsmclumiprof->Fill(pileupinfo->getTrueNumInteractions(),pvcoll->size(),weight);
00246 
00247     if(m_useweight) {
00248       m_hrecovsmcnvtxweightedprof->Fill(npileup,pvcoll->size(),1.-weight);
00249       m_hrecovsmclumiweightedprof->Fill(pileupinfo->getTrueNumInteractions(),pvcoll->size(),1.-weight);
00250     }
00251     //
00252     
00253     Handle< HepMCProduct > EvtHandle ;
00254     iEvent.getByLabel(m_mctruthcollection, EvtHandle ) ;
00255     
00256     const HepMC::GenEvent* Evt = EvtHandle->GetEvent();
00257     
00258     // compute the difference between the main interaction vertex z position and the first vertex of the collection
00259     
00260     if(pvcoll->size() !=0) {
00261       if(!(*pvcoll)[0].isFake()) {
00262         // get the first vertex
00263         if(Evt->vertices_begin() != Evt->vertices_end()) {
00264           m_hdeltazfirst->Fill((*pvcoll)[0].z()-(*Evt->vertices_begin())->point3d().z()/10.,weight);
00265           m_hdeltazfirstvsnpu->Fill(npileup,(*pvcoll)[0].z()-(*Evt->vertices_begin())->point3d().z()/10.,weight);
00266         }
00267       }
00268     }
00269     
00270     // compute the difference between the main interaction vertex z position and the closest reco vertex  
00271     
00272     double minabsdist = -1.;
00273     double mindist = -999.;
00274     int closestvtx = -1;
00275     
00276     for(unsigned int ivtx = 0 ; ivtx < pvcoll->size() ; ++ivtx) {
00277       
00278       if(closestvtx < 0 || minabsdist > std::abs((*pvcoll)[ivtx].z()-(*Evt->vertices_begin())->point3d().z()/10.)) {
00279         mindist = (*pvcoll)[ivtx].z()-(*Evt->vertices_begin())->point3d().z()/10.;
00280         closestvtx = ivtx;
00281         minabsdist = std::abs(mindist);
00282       }
00283       
00284     }
00285     if(closestvtx >= 0) {
00286       m_hdeltazclose->Fill(mindist,weight);
00287       m_hdeltazclosevsnpu->Fill(npileup,mindist,weight);
00288       m_hclosestvtx->Fill(closestvtx,weight);
00289       m_hclosestvtxvsnpu->Fill(npileup,closestvtx,weight);
00290     }
00291     
00292   }
00293 }
00294   
00295   void 
00296 MCvsRecoVerticesAnalyzer::beginRun(const edm::Run& iRun, const edm::EventSetup&)
00297 {
00298 }
00299 
00300 void 
00301 MCvsRecoVerticesAnalyzer::endRun(const edm::Run& iRun, const edm::EventSetup&)
00302 {
00303 }
00304 
00305 
00306 
00307 // ------------ method called once each job just before starting event loop  ------------
00308 void 
00309 MCvsRecoVerticesAnalyzer::beginJob()
00310 {
00311 }
00312 
00313 // ------------ method called once each job just after ending the event loop  ------------
00314 void 
00315 MCvsRecoVerticesAnalyzer::endJob() 
00316 {
00317 }
00318 
00319 //define this as a plug-in
00320 DEFINE_FWK_MODULE(MCvsRecoVerticesAnalyzer);