CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/PhysicsTools/PatExamples/plugins/PatVertexAnalyzer.cc

Go to the documentation of this file.
00001 #include <iostream>
00002 #include <cmath>
00003 #include <vector>
00004 #include <string>
00005 
00006 #include <TH1.h>
00007 #include <TProfile.h>
00008 
00009 #include "FWCore/Framework/interface/Frameworkfwd.h"
00010 #include "FWCore/Framework/interface/Event.h"
00011 #include "FWCore/Framework/interface/EventSetup.h"
00012 #include "FWCore/Framework/interface/EDAnalyzer.h"
00013 #include "FWCore/Utilities/interface/InputTag.h"
00014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00015 #include "FWCore/ServiceRegistry/interface/Service.h"
00016 
00017 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00018 
00019 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00020 #include "DataFormats/TrackReco/interface/Track.h"
00021 #include "DataFormats/VertexReco/interface/Vertex.h"
00022 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00023 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00024 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
00025 
00026 class PatVertexAnalyzer : public edm::EDAnalyzer  {
00027     public: 
00029         PatVertexAnalyzer(const edm::ParameterSet &params);
00030         ~PatVertexAnalyzer();
00031 
00032         // virtual methods called from base class EDAnalyzer
00033         virtual void beginJob();
00034         virtual void analyze(const edm::Event &event, const edm::EventSetup &es);
00035 
00036     private:
00037         // configuration parameters
00038         edm::InputTag src_;
00039         edm::InputTag genParticles_;
00040 
00041         TH1 *nVertices_, *nTracks_;
00042         TH1 *x_, *y_, *z_;
00043         TH1 *xErr_, *yErr_, *zErr_;
00044         TH1 *xDelta_, *yDelta_, *zDelta_;
00045         TH1 *xPull_, *yPull_, *zPull_;
00046 };
00047 
00048 PatVertexAnalyzer::PatVertexAnalyzer(const edm::ParameterSet &params) :
00049         src_(params.getParameter<edm::InputTag>("src")),
00050         genParticles_(params.getParameter<edm::InputTag>("mc"))
00051 {
00052 }
00053 
00054 PatVertexAnalyzer::~PatVertexAnalyzer()
00055 {
00056 }
00057 
00058 void PatVertexAnalyzer::beginJob()
00059 {
00060         // retrieve handle to auxiliary service
00061         //  used for storing histograms into ROOT file
00062         edm::Service<TFileService> fs;
00063 
00064         nVertices_ = fs->make<TH1F>("nVertices", "number of reconstructed primary vertices", 50, 0, 50);
00065         nTracks_ = fs->make<TH1F>("nTracks", "number of tracks at primary vertex", 100, 0, 300);
00066         x_ = fs->make<TH1F>("pvX", "primary vertex x", 100, -0.1, 0.1);
00067         y_ = fs->make<TH1F>("pvY", "primary vertex y", 100, -0.1, 0.1);
00068         z_ = fs->make<TH1F>("pvZ", "primary vertex z", 100, -30, 30);
00069         xErr_ = fs->make<TH1F>("pvErrorX", "primary vertex x error", 100, 0, 0.005);
00070         yErr_ = fs->make<TH1F>("pvErrorY", "primary vertex y error", 100, 0, 0.005);
00071         zErr_ = fs->make<TH1F>("pvErrorZ", "primary vertex z error", 100, 0, 0.01);
00072         xDelta_ = fs->make<TH1F>("pvDeltaX", "x shift wrt simulated vertex", 100, -0.01, 0.01);
00073         yDelta_ = fs->make<TH1F>("pvDeltaY", "y shift wrt simulated vertex", 100, -0.01, 0.01);
00074         zDelta_ = fs->make<TH1F>("pvDeltaZ", "z shift wrt simulated vertex", 100, -0.02, 0.02);
00075         xPull_ = fs->make<TH1F>("pvPullX", "primary vertex x pull", 100, -5, 5);
00076         yPull_ = fs->make<TH1F>("pvPullY", "primary vertex y pull", 100, -5, 5);
00077         zPull_ = fs->make<TH1F>("pvPullZ", "primary vertex z pull", 100, -5, 5);
00078 }
00079 
00080 void PatVertexAnalyzer::analyze(const edm::Event &event, const edm::EventSetup &es)
00081 {  
00082         // handle to the primary vertex collection
00083         edm::Handle<reco::VertexCollection> pvHandle;
00084         event.getByLabel(src_, pvHandle);
00085 
00086         // handle to the generator particles (i.e. the MC truth)
00087         edm::Handle<reco::GenParticleCollection> genParticlesHandle;
00088         event.getByLabel(genParticles_, genParticlesHandle);
00089 
00090         // extract the position of the simulated vertex
00091         math::XYZPoint simPV = (*genParticlesHandle)[2].vertex();
00092 
00093         // the number of reconstructed primary vertices
00094         nVertices_->Fill(pvHandle->size());
00095 
00096         // if we have at least one, use the first (highest pt^2 sum)
00097         if (!pvHandle->empty()) {
00098                 const reco::Vertex &pv = (*pvHandle)[0];
00099 
00100                 nTracks_->Fill(pv.tracksSize());
00101 
00102                 x_->Fill(pv.x());
00103                 y_->Fill(pv.y());
00104                 z_->Fill(pv.z());
00105 
00106                 xErr_->Fill(pv.xError());
00107                 yErr_->Fill(pv.yError());
00108                 zErr_->Fill(pv.zError());
00109 
00110                 xDelta_->Fill(pv.x() - simPV.X());
00111                 yDelta_->Fill(pv.y() - simPV.Y());
00112                 zDelta_->Fill(pv.z() - simPV.Z());
00113 
00114                 xPull_->Fill((pv.x() - simPV.X()) / pv.xError());
00115                 yPull_->Fill((pv.y() - simPV.Y()) / pv.yError());
00116                 zPull_->Fill((pv.z() - simPV.Z()) / pv.zError());
00117 
00118                 // we could access the tracks using the
00119                 // pv.tracks_begin() ... pv.tracks_end() iterators
00120         }
00121 }
00122         
00123 #include "FWCore/Framework/interface/MakerMacros.h"
00124 
00125 DEFINE_FWK_MODULE(PatVertexAnalyzer);