CMS 3D CMS Logo

PatVertexAnalyzer.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <cmath>
3 #include <vector>
4 #include <string>
5 
6 #include <TH1.h>
7 #include <TProfile.h>
8 
16 
18 
25 
27  public:
29  PatVertexAnalyzer(const edm::ParameterSet &params);
31 
32  // virtual methods called from base class EDAnalyzer
33  virtual void beginJob() override;
34  virtual void analyze(const edm::Event &event, const edm::EventSetup &es) override;
35 
36  private:
37  // configuration parameters
40 
42  TH1 *x_, *y_, *z_;
43  TH1 *xErr_, *yErr_, *zErr_;
45  TH1 *xPull_, *yPull_, *zPull_;
46 };
47 
49  srcToken_(consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("src"))),
50  genParticlesToken_(consumes<reco::GenParticleCollection>(params.getParameter<edm::InputTag>("mc")))
51 {
52 }
53 
55 {
56 }
57 
59 {
60  // retrieve handle to auxiliary service
61  // used for storing histograms into ROOT file
63 
64  nVertices_ = fs->make<TH1F>("nVertices", "number of reconstructed primary vertices", 50, 0, 50);
65  nTracks_ = fs->make<TH1F>("nTracks", "number of tracks at primary vertex", 100, 0, 300);
66  x_ = fs->make<TH1F>("pvX", "primary vertex x", 100, -0.1, 0.1);
67  y_ = fs->make<TH1F>("pvY", "primary vertex y", 100, -0.1, 0.1);
68  z_ = fs->make<TH1F>("pvZ", "primary vertex z", 100, -30, 30);
69  xErr_ = fs->make<TH1F>("pvErrorX", "primary vertex x error", 100, 0, 0.005);
70  yErr_ = fs->make<TH1F>("pvErrorY", "primary vertex y error", 100, 0, 0.005);
71  zErr_ = fs->make<TH1F>("pvErrorZ", "primary vertex z error", 100, 0, 0.01);
72  xDelta_ = fs->make<TH1F>("pvDeltaX", "x shift wrt simulated vertex", 100, -0.01, 0.01);
73  yDelta_ = fs->make<TH1F>("pvDeltaY", "y shift wrt simulated vertex", 100, -0.01, 0.01);
74  zDelta_ = fs->make<TH1F>("pvDeltaZ", "z shift wrt simulated vertex", 100, -0.02, 0.02);
75  xPull_ = fs->make<TH1F>("pvPullX", "primary vertex x pull", 100, -5, 5);
76  yPull_ = fs->make<TH1F>("pvPullY", "primary vertex y pull", 100, -5, 5);
77  zPull_ = fs->make<TH1F>("pvPullZ", "primary vertex z pull", 100, -5, 5);
78 }
79 
81 {
82  // handle to the primary vertex collection
84  event.getByToken(srcToken_, pvHandle);
85 
86  // handle to the generator particles (i.e. the MC truth)
88  event.getByToken(genParticlesToken_, genParticlesHandle);
89 
90  // extract the position of the simulated vertex
91  math::XYZPoint simPV = (*genParticlesHandle)[2].vertex();
92 
93  // the number of reconstructed primary vertices
94  nVertices_->Fill(pvHandle->size());
95 
96  // if we have at least one, use the first (highest pt^2 sum)
97  if (!pvHandle->empty()) {
98  const reco::Vertex &pv = (*pvHandle)[0];
99 
100  nTracks_->Fill(pv.tracksSize());
101 
102  x_->Fill(pv.x());
103  y_->Fill(pv.y());
104  z_->Fill(pv.z());
105 
106  xErr_->Fill(pv.xError());
107  yErr_->Fill(pv.yError());
108  zErr_->Fill(pv.zError());
109 
110  xDelta_->Fill(pv.x() - simPV.X());
111  yDelta_->Fill(pv.y() - simPV.Y());
112  zDelta_->Fill(pv.z() - simPV.Z());
113 
114  xPull_->Fill((pv.x() - simPV.X()) / pv.xError());
115  yPull_->Fill((pv.y() - simPV.Y()) / pv.yError());
116  zPull_->Fill((pv.z() - simPV.Z()) / pv.zError());
117 
118  // we could access the tracks using the
119  // pv.tracks_begin() ... pv.tracks_end() iterators
120  }
121 }
122 
124 
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
edm::EDGetTokenT< reco::GenParticleCollection > genParticlesToken_
double zError() const
error on z
Definition: Vertex.h:123
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
double y() const
y coordinate
Definition: Vertex.h:113
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
edm::EDGetTokenT< reco::VertexCollection > srcToken_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
virtual void beginJob() override
virtual void analyze(const edm::Event &event, const edm::EventSetup &es) override
def pv(vc)
Definition: MetAnalyzer.py:6
double z() const
z coordinate
Definition: Vertex.h:115
double x() const
x coordinate
Definition: Vertex.h:111
double xError() const
error on x
Definition: Vertex.h:119
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
fixed size matrix
HLT enums.
PatVertexAnalyzer(const edm::ParameterSet &params)
constructor and destructor
double yError() const
error on y
Definition: Vertex.h:121
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:71
Definition: event.py:1