CMS 3D CMS Logo

GeneralPurposeVertexAnalyzer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Alignment/GeneralPurposeVertexAnalyzer
4 // Class: GeneralPurposeVertexAnalyzer
5 //
10 //
11 // Original Author: Marco Musich
12 // Created: Thu, 13 Apr 2023 14:16:43 GMT
13 //
14 //
15 
16 // ROOT includes files
17 #include "TMath.h"
18 #include "TFile.h"
19 #include "TH1I.h"
20 #include "TH1D.h"
21 #include "TH2D.h"
22 #include "TProfile.h"
23 #include "TProfile2D.h"
24 
25 // system include files
26 #include <memory>
27 #include <fmt/format.h>
28 #include <fmt/printf.h>
29 
30 // user include files
48 
49 //
50 // class declaration
51 //
52 
54 
55 class GeneralPurposeVertexAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
56 public:
58  ~GeneralPurposeVertexAnalyzer() override = default;
59 
60  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
61 
62 private:
63  void pvTracksPlots(const reco::Vertex &v);
64  void vertexPlots(const reco::Vertex &v, const reco::BeamSpot &beamSpot, int i);
65  template <typename T, typename... Args>
66  T *book(const Args &...args) const;
67  void beginJob() override;
68  void analyze(const edm::Event &, const edm::EventSetup &) override;
69 
70  // ----------member data ---------------------------
72 
73  const int ndof_;
75 
81 
82  static constexpr int cmToUm = 10000;
83 
84  const double vposx_;
85  const double vposy_;
86  const int tkNoBin_;
87  const double tkNoMin_;
88  const double tkNoMax_;
89 
90  const int dxyBin_;
91  const double dxyMin_;
92  const double dxyMax_;
93 
94  const int dzBin_;
95  const double dzMin_;
96  const double dzMax_;
97 
98  const int phiBin_;
99  const int phiBin2D_;
100  const double phiMin_;
101  const double phiMax_;
102 
103  const int etaBin_;
104  const int etaBin2D_;
105  const double etaMin_;
106  const double etaMax_;
107 
108  // the histos
109  TH1I *nbvtx, *nbgvtx;
110  TH1D *nbtksinvtx[2], *trksWeight[2], *score[2];
111  TH1D *tt[2];
112  TH1D *xrec[2], *yrec[2], *zrec[2], *xDiff[2], *yDiff[2], *xerr[2], *yerr[2], *zerr[2];
113  TH2D *xerrVsTrks[2], *yerrVsTrks[2], *zerrVsTrks[2];
114  TH1D *ntracksVsZ[2];
115  TH1D *vtxchi2[2], *vtxndf[2], *vtxprob[2], *nans[2];
116  TH1D *type[2];
118 
120  TH1D *dxy, *dxy2, *dz, *dxyErr, *dzErr;
121  TH1D *phi_pt1, *eta_pt1;
129 };
130 
131 //
132 // constructors and destructor
133 //
135  : ndof_(iConfig.getParameter<int>("ndof")),
136  errorPrinted_(false),
137  vertexInputTag_(iConfig.getParameter<edm::InputTag>("vertexLabel")),
138  beamSpotInputTag_(iConfig.getParameter<edm::InputTag>("beamSpotLabel")),
139  vertexToken_(consumes<reco::VertexCollection>(vertexInputTag_)),
140  scoreToken_(consumes<VertexScore>(vertexInputTag_)),
141  beamspotToken_(consumes<reco::BeamSpot>(beamSpotInputTag_)),
142  // to be configured for each year...
143  vposx_(iConfig.getParameter<double>("Xpos")),
144  vposy_(iConfig.getParameter<double>("Ypos")),
145  tkNoBin_(iConfig.getParameter<int>("TkSizeBin")),
146  tkNoMin_(iConfig.getParameter<double>("TkSizeMin")),
147  tkNoMax_(iConfig.getParameter<double>("TkSizeMax")),
148  dxyBin_(iConfig.getParameter<int>("DxyBin")),
149  dxyMin_(iConfig.getParameter<double>("DxyMin")),
150  dxyMax_(iConfig.getParameter<double>("DxyMax")),
151  dzBin_(iConfig.getParameter<int>("DzBin")),
152  dzMin_(iConfig.getParameter<double>("DzMin")),
153  dzMax_(iConfig.getParameter<double>("DzMax")),
154  phiBin_(iConfig.getParameter<int>("PhiBin")),
155  phiBin2D_(iConfig.getParameter<int>("PhiBin2D")),
156  phiMin_(iConfig.getParameter<double>("PhiMin")),
157  phiMax_(iConfig.getParameter<double>("PhiMax")),
158  etaBin_(iConfig.getParameter<int>("EtaBin")),
159  etaBin2D_(iConfig.getParameter<int>("EtaBin2D")),
160  etaMin_(iConfig.getParameter<double>("EtaMin")),
161  etaMax_(iConfig.getParameter<double>("EtaMax")),
162  // histograms
163  nbvtx(nullptr),
164  bsX(nullptr),
165  bsY(nullptr),
166  bsZ(nullptr),
167  bsSigmaZ(nullptr),
168  bsDxdz(nullptr),
169  bsDydz(nullptr),
170  bsBeamWidthX(nullptr),
171  bsBeamWidthY(nullptr),
172  bsType(nullptr),
173  sumpt(nullptr),
174  ntracks(nullptr),
175  weight(nullptr),
176  chi2ndf(nullptr),
177  chi2prob(nullptr),
178  dxy(nullptr),
179  dxy2(nullptr),
180  dz(nullptr),
181  dxyErr(nullptr),
182  dzErr(nullptr),
183  phi_pt1(nullptr),
184  eta_pt1(nullptr),
185  phi_pt10(nullptr),
186  eta_pt10(nullptr),
187  dxyVsPhi_pt1(nullptr),
188  dzVsPhi_pt1(nullptr),
189  dxyVsEta_pt1(nullptr),
190  dzVsEta_pt1(nullptr),
191  dxyVsEtaVsPhi_pt1(nullptr),
192  dzVsEtaVsPhi_pt1(nullptr),
193  dxyVsPhi_pt10(nullptr),
194  dzVsPhi_pt10(nullptr),
195  dxyVsEta_pt10(nullptr),
196  dzVsEta_pt10(nullptr),
197  dxyVsEtaVsPhi_pt10(nullptr),
198  dzVsEtaVsPhi_pt10(nullptr) {
199  usesResource(TFileService::kSharedResource);
200 }
201 
202 //
203 // member functions
204 //
205 
206 // ------------ method called for each event ------------
208  using namespace edm;
209 
210  const auto &recVtxs = iEvent.getHandle(vertexToken_);
211  const auto &scores = iEvent.getHandle(scoreToken_);
212  const auto &beamSpotHandle = iEvent.getHandle(beamspotToken_);
213  reco::BeamSpot beamSpot = *beamSpotHandle;
214 
215  // check for absent products and simply "return" in that case
216  if (recVtxs.isValid() == false || beamSpotHandle.isValid() == false) {
217  edm::LogWarning("GeneralPurposeVertexAnalyzer")
218  << " Some products not available in the event: VertexCollection " << vertexInputTag_ << " " << recVtxs.isValid()
219  << " BeamSpot " << beamSpotInputTag_ << " " << beamSpotHandle.isValid() << ". Skipping plots for this event";
220  return;
221  }
222 
223  // check upfront that refs to track are (likely) to be valid
224  bool ok{true};
225  for (const auto &v : *recVtxs) {
226  if (v.tracksSize() > 0) {
227  const auto &ref = v.trackRefAt(0);
228  if (ref.isNull() || !ref.isAvailable()) {
229  if (!errorPrinted_) {
230  edm::LogWarning("GeneralPurposeVertexAnalyzer")
231  << "Skipping vertex collection: " << vertexInputTag_
232  << " since likely the track collection the vertex has refs pointing to is missing (at least the first "
233  "TrackBaseRef is null or not available)";
234  } else {
235  errorPrinted_ = true;
236  }
237  ok = false;
238  }
239  }
240  }
241 
242  if (!ok) {
243  return;
244  }
245 
246  nbvtx->Fill(recVtxs->size());
247  int ng = 0;
248  for (auto const &vx : (*recVtxs)) {
249  if (vx.isValid() && !vx.isFake() && vx.ndof() >= ndof_) {
250  ++ng;
251  }
252  }
253  nbgvtx->Fill(ng);
254 
255  if (scores.isValid() && !(*scores).empty()) {
256  auto pvScore = (*scores).get(0);
257  score[1]->Fill(std::sqrt(pvScore));
258  for (unsigned int i = 1; i < (*scores).size(); ++i) {
259  score[0]->Fill(std::sqrt((*scores).get(i)));
260  }
261  }
262 
263  // fill PV tracks MEs (as now, for alignment)
264  if (!recVtxs->empty()) {
265  vertexPlots(recVtxs->front(), beamSpot, 1);
266  pvTracksPlots(recVtxs->front());
267 
268  for (reco::VertexCollection::const_iterator v = recVtxs->begin() + 1; v != recVtxs->end(); ++v) {
269  vertexPlots(*v, beamSpot, 0);
270  }
271  }
272 
273  // Beamline plots:
274  bsX->Fill(beamSpot.x0());
275  bsY->Fill(beamSpot.y0());
276  bsZ->Fill(beamSpot.z0());
277  bsSigmaZ->Fill(beamSpot.sigmaZ());
278  bsDxdz->Fill(beamSpot.dxdz());
279  bsDydz->Fill(beamSpot.dydz());
280  bsBeamWidthX->Fill(beamSpot.BeamWidthX() * cmToUm);
281  bsBeamWidthY->Fill(beamSpot.BeamWidthY() * cmToUm);
282  bsType->Fill(beamSpot.type());
283 }
284 
286  if (!v.isValid())
287  return;
288  if (v.isFake())
289  return;
290 
291  if (v.tracksSize() == 0) {
292  ntracks->Fill(0);
293  return;
294  }
295 
296  const math::XYZPoint myVertex(v.position().x(), v.position().y(), v.position().z());
297 
298  float sumPT = 0.;
299  for (const auto &t : v.tracks()) {
300  const bool isHighPurity = t->quality(reco::TrackBase::highPurity);
301  if (!isHighPurity) {
302  continue;
303  }
304 
305  const float pt = t->pt();
306  if (pt < 1.f) {
307  continue;
308  }
309 
310  const float pt2 = pt * pt;
311  const float eta = t->eta();
312  const float phi = t->phi();
313  const float w = v.trackWeight(t);
314  const float chi2NDF = t->normalizedChi2();
315  const float chi2Prob = TMath::Prob(t->chi2(), static_cast<int>(t->ndof()));
316  const float Dxy = t->dxy(myVertex) * cmToUm;
317  const float Dz = t->dz(myVertex) * cmToUm;
318  const float DxyErr = t->dxyError() * cmToUm;
319  const float DzErr = t->dzError() * cmToUm;
320 
321  sumPT += pt2;
322 
323  weight->Fill(w);
324  chi2ndf->Fill(chi2NDF);
325  chi2prob->Fill(chi2Prob);
326  dxy->Fill(Dxy);
327  dxy2->Fill(Dxy);
328  dz->Fill(Dz);
329  dxyErr->Fill(DxyErr);
330  dzErr->Fill(DzErr);
331  phi_pt1->Fill(phi);
332  eta_pt1->Fill(eta);
333  dxyVsPhi_pt1->Fill(phi, Dxy);
334  dzVsPhi_pt1->Fill(phi, Dz);
335  dxyVsEta_pt1->Fill(eta, Dxy);
336  dzVsEta_pt1->Fill(eta, Dz);
337  dxyVsEtaVsPhi_pt1->Fill(eta, phi, Dxy);
338  dzVsEtaVsPhi_pt1->Fill(eta, phi, Dz);
339 
340  if (pt >= 10.f) {
341  phi_pt10->Fill(phi);
342  eta_pt10->Fill(eta);
343  dxyVsPhi_pt10->Fill(phi, Dxy);
344  dzVsPhi_pt10->Fill(phi, Dz);
345  dxyVsEta_pt10->Fill(eta, Dxy);
346  dzVsEta_pt10->Fill(eta, Dz);
347  dxyVsEtaVsPhi_pt10->Fill(eta, phi, Dxy);
348  dzVsEtaVsPhi_pt10->Fill(eta, phi, Dz);
349  }
350  }
351 
352  ntracks->Fill(static_cast<float>(v.tracks().size()));
353  sumpt->Fill(sumPT);
354 }
355 
357  if (i < 0 || i > 1)
358  return;
359  if (!v.isValid())
360  type[i]->Fill(2.);
361  else if (v.isFake())
362  type[i]->Fill(1.);
363  else
364  type[i]->Fill(0.);
365 
366  if (v.isValid() && !v.isFake()) {
367  float weight = 0;
368  for (reco::Vertex::trackRef_iterator t = v.tracks_begin(); t != v.tracks_end(); t++)
369  weight += v.trackWeight(*t);
370  trksWeight[i]->Fill(weight);
371  nbtksinvtx[i]->Fill(v.tracksSize());
372  ntracksVsZ[i]->Fill(v.position().z() - beamSpot.z0(), v.tracksSize());
373 
374  vtxchi2[i]->Fill(v.chi2());
375  vtxndf[i]->Fill(v.ndof());
376  vtxprob[i]->Fill(ChiSquaredProbability(v.chi2(), v.ndof()));
377 
378  xrec[i]->Fill(v.position().x());
379  yrec[i]->Fill(v.position().y());
380  zrec[i]->Fill(v.position().z());
381 
382  float xb = beamSpot.x0() + beamSpot.dxdz() * (v.position().z() - beamSpot.z0());
383  float yb = beamSpot.y0() + beamSpot.dydz() * (v.position().z() - beamSpot.z0());
384  xDiff[i]->Fill((v.position().x() - xb) * cmToUm);
385  yDiff[i]->Fill((v.position().y() - yb) * cmToUm);
386 
387  xerr[i]->Fill(v.xError() * cmToUm);
388  yerr[i]->Fill(v.yError() * cmToUm);
389  zerr[i]->Fill(v.zError() * cmToUm);
390  xerrVsTrks[i]->Fill(weight, v.xError() * cmToUm);
391  yerrVsTrks[i]->Fill(weight, v.yError() * cmToUm);
392  zerrVsTrks[i]->Fill(weight, v.zError() * cmToUm);
393 
394  nans[i]->Fill(1., edm::isNotFinite(v.position().x()) * 1.);
395  nans[i]->Fill(2., edm::isNotFinite(v.position().y()) * 1.);
396  nans[i]->Fill(3., edm::isNotFinite(v.position().z()) * 1.);
397 
398  int index = 3;
399  for (int k = 0; k != 3; k++) {
400  for (int j = k; j != 3; j++) {
401  index++;
402  nans[i]->Fill(index * 1., edm::isNotFinite(v.covariance(k, j)) * 1.);
403  // in addition, diagonal element must be positive
404  if (j == k && v.covariance(k, j) < 0) {
405  nans[i]->Fill(index * 1., 1.);
406  }
407  }
408  }
409  }
410 }
411 
412 template <typename T, typename... Args>
413 T *GeneralPurposeVertexAnalyzer::book(const Args &...args) const {
414  T *t = fs_->make<T>(args...);
415  return t;
416 }
417 
418 // ------------ method called once each job just before starting event loop ------------
420  nbvtx = book<TH1I>("vtxNbr", "Reconstructed Vertices in Event", 80, -0.5, 79.5);
421  nbgvtx = book<TH1I>("goodvtxNbr", "Reconstructed Good Vertices in Event", 80, -0.5, 79.5);
422 
423  nbtksinvtx[0] = book<TH1D>("otherVtxTrksNbr", "Reconstructed Tracks in Vertex (other Vtx)", 40, -0.5, 99.5);
424  ntracksVsZ[0] =
425  book<TProfile>("otherVtxTrksVsZ", "Reconstructed Tracks in Vertex (other Vtx) vs Z", 80, -20., 20., 0., 100., "");
426  ntracksVsZ[0]->SetXTitle("z-bs");
427  ntracksVsZ[0]->SetYTitle("#tracks");
428 
429  score[0] = book<TH1D>("otherVtxScore", "sqrt(score) (other Vtx)", 100, 0., 400.);
430  trksWeight[0] = book<TH1D>("otherVtxTrksWeight", "Total weight of Tracks in Vertex (other Vtx)", 40, 0, 100.);
431  vtxchi2[0] = book<TH1D>("otherVtxChi2", "#chi^{2} (other Vtx)", 100, 0., 200.);
432  vtxndf[0] = book<TH1D>("otherVtxNdf", "ndof (other Vtx)", 100, 0., 200.);
433  vtxprob[0] = book<TH1D>("otherVtxProb", "#chi^{2} probability (other Vtx)", 100, 0., 1.);
434  nans[0] = book<TH1D>("otherVtxNans", "Illegal values for x,y,z,xx,xy,xz,yy,yz,zz (other Vtx)", 9, 0.5, 9.5);
435 
436  nbtksinvtx[1] = book<TH1D>("tagVtxTrksNbr", "Reconstructed Tracks in Vertex (tagged Vtx)", 100, -0.5, 99.5);
437  ntracksVsZ[1] =
438  book<TProfile>("tagVtxTrksVsZ", "Reconstructed Tracks in Vertex (tagged Vtx) vs Z", 80, -20., 20., 0., 100., "");
439  ntracksVsZ[1]->SetXTitle("z-bs");
440  ntracksVsZ[1]->SetYTitle("#tracks");
441 
442  score[1] = book<TH1D>("tagVtxScore", "sqrt(score) (tagged Vtx)", 100, 0., 400.);
443  trksWeight[1] = book<TH1D>("tagVtxTrksWeight", "Total weight of Tracks in Vertex (tagged Vtx)", 100, 0, 100.);
444  vtxchi2[1] = book<TH1D>("tagVtxChi2", "#chi^{2} (tagged Vtx)", 100, 0., 200.);
445  vtxndf[1] = book<TH1D>("tagVtxNdf", "ndof (tagged Vtx)", 100, 0., 200.);
446  vtxprob[1] = book<TH1D>("tagVtxProb", "#chi^{2} probability (tagged Vtx)", 100, 0., 1.);
447  nans[1] = book<TH1D>("tagVtxNans", "Illegal values for x,y,z,xx,xy,xz,yy,yz,zz (tagged Vtx)", 9, 0.5, 9.5);
448 
449  xrec[0] = book<TH1D>("otherPosX", "Position x Coordinate (other Vtx)", 100, vposx_ - 0.1, vposx_ + 0.1);
450  yrec[0] = book<TH1D>("otherPosY", "Position y Coordinate (other Vtx)", 100, vposy_ - 0.1, vposy_ + 0.1);
451  zrec[0] = book<TH1D>("otherPosZ", "Position z Coordinate (other Vtx)", 100, -20., 20.);
452  xDiff[0] = book<TH1D>("otherDiffX", "X distance from BeamSpot (other Vtx)", 100, -500, 500);
453  yDiff[0] = book<TH1D>("otherDiffY", "Y distance from BeamSpot (other Vtx)", 100, -500, 500);
454  xerr[0] = book<TH1D>("otherErrX", "Uncertainty x Coordinate (other Vtx)", 100, 0., 100);
455  yerr[0] = book<TH1D>("otherErrY", "Uncertainty y Coordinate (other Vtx)", 100, 0., 100);
456  zerr[0] = book<TH1D>("otherErrZ", "Uncertainty z Coordinate (other Vtx)", 100, 0., 100);
457  xerrVsTrks[0] = book<TH2D>(
458  "otherErrVsWeightX", "Uncertainty x Coordinate vs. track weight (other Vtx)", 100, 0, 100., 100, 0., 100);
459  yerrVsTrks[0] = book<TH2D>(
460  "otherErrVsWeightY", "Uncertainty y Coordinate vs. track weight (other Vtx)", 100, 0, 100., 100, 0., 100);
461  zerrVsTrks[0] = book<TH2D>(
462  "otherErrVsWeightZ", "Uncertainty z Coordinate vs. track weight (other Vtx)", 100, 0, 100., 100, 0., 100);
463 
464  xrec[1] = book<TH1D>("tagPosX", "Position x Coordinate (tagged Vtx)", 100, vposx_ - 0.1, vposx_ + 0.1);
465  yrec[1] = book<TH1D>("tagPosY", "Position y Coordinate (tagged Vtx)", 100, vposy_ - 0.1, vposy_ + 0.1);
466  zrec[1] = book<TH1D>("tagPosZ", "Position z Coordinate (tagged Vtx)", 100, -20., 20.);
467  xDiff[1] = book<TH1D>("tagDiffX", "X distance from BeamSpot (tagged Vtx)", 100, -500, 500);
468  yDiff[1] = book<TH1D>("tagDiffY", "Y distance from BeamSpot (tagged Vtx)", 100, -500, 500);
469  xerr[1] = book<TH1D>("tagErrX", "Uncertainty x Coordinate (tagged Vtx)", 100, 0., 100);
470  yerr[1] = book<TH1D>("tagErrY", "Uncertainty y Coordinate (tagged Vtx)", 100, 0., 100);
471  zerr[1] = book<TH1D>("tagErrZ", "Uncertainty z Coordinate (tagged Vtx)", 100, 0., 100);
472  xerrVsTrks[1] = book<TH2D>(
473  "tagErrVsWeightX", "Uncertainty x Coordinate vs. track weight (tagged Vtx)", 100, 0, 100., 100, 0., 100);
474  yerrVsTrks[1] = book<TH2D>(
475  "tagErrVsWeightY", "Uncertainty y Coordinate vs. track weight (tagged Vtx)", 100, 0, 100., 100, 0., 100);
476  zerrVsTrks[1] = book<TH2D>(
477  "tagErrVsWeightZ", "Uncertainty z Coordinate vs. track weight (tagged Vtx)", 100, 0, 100., 100, 0., 100);
478 
479  type[0] = book<TH1D>("otherType", "Vertex type (other Vtx)", 3, -0.5, 2.5);
480  type[1] = book<TH1D>("tagType", "Vertex type (tagged Vtx)", 3, -0.5, 2.5);
481 
482  for (int i = 0; i < 2; ++i) {
483  type[i]->GetXaxis()->SetBinLabel(1, "Valid, real");
484  type[i]->GetXaxis()->SetBinLabel(2, "Valid, fake");
485  type[i]->GetXaxis()->SetBinLabel(3, "Invalid");
486  }
487 
488  bsX = book<TH1D>("bsX", "BeamSpot x0", 100, vposx_ - 0.1, vposx_ + 0.1);
489  bsY = book<TH1D>("bsY", "BeamSpot y0", 100, vposy_ - 0.1, vposy_ + 0.1);
490  bsZ = book<TH1D>("bsZ", "BeamSpot z0", 100, -2., 2.);
491  bsSigmaZ = book<TH1D>("bsSigmaZ", "BeamSpot sigmaZ", 100, 0., 10.);
492  bsDxdz = book<TH1D>("bsDxdz", "BeamSpot dxdz", 100, -0.0003, 0.0003);
493  bsDydz = book<TH1D>("bsDydz", "BeamSpot dydz", 100, -0.0003, 0.0003);
494  bsBeamWidthX = book<TH1D>("bsBeamWidthX", "BeamSpot BeamWidthX", 100, 0., 100.);
495  bsBeamWidthY = book<TH1D>("bsBeamWidthY", "BeamSpot BeamWidthY", 100, 0., 100.);
496  bsType = book<TH1D>("bsType", "BeamSpot type", 4, -1.5, 2.5);
497  bsType->GetXaxis()->SetBinLabel(1, "Unknown");
498  bsType->GetXaxis()->SetBinLabel(2, "Fake");
499  bsType->GetXaxis()->SetBinLabel(3, "LHC");
500  bsType->GetXaxis()->SetBinLabel(4, "Tracker");
501 
502  // repeated strings in titles
503  std::string s_1 = "PV Tracks (p_{T} > 1 GeV)";
504  std::string s_10 = "PV Tracks (p_{T} > 10 GeV)";
505 
506  ntracks = book<TH1D>("ntracks", fmt::sprintf("number of %s", s_1).c_str(), tkNoBin_, tkNoMin_, tkNoMax_);
507  ntracks->SetXTitle(fmt::sprintf("Number of %s per Event", s_1).c_str());
508  ntracks->SetYTitle("Number of Events");
509 
510  weight = book<TH1D>("weight", fmt::sprintf("weight of %s", s_1).c_str(), 100, 0., 1.);
511  weight->SetXTitle(fmt::sprintf("weight of %s per Event", s_1).c_str());
512  weight->SetYTitle("Number of Events");
513 
514  sumpt = book<TH1D>("sumpt", fmt::sprintf("#Sum p_{T} of %s", s_1).c_str(), 100, -0.5, 249.5);
515  chi2ndf = book<TH1D>("chi2ndf", fmt::sprintf("%s #chi^{2}/ndof", s_1).c_str(), 100, 0., 20.);
516  chi2prob = book<TH1D>("chi2prob", fmt::sprintf("%s #chi^{2} probability", s_1).c_str(), 100, 0., 1.);
517  dxy = book<TH1D>("dxy", fmt::sprintf("%s d_{xy} (#mum)", s_1).c_str(), dxyBin_, dxyMin_, dxyMax_);
518  dxy2 = book<TH1D>("dxyzoom", fmt::sprintf("%s d_{xy} (#mum)", s_1).c_str(), dxyBin_, dxyMin_ / 5., dxyMax_ / 5.);
519  dxyErr = book<TH1D>("dxyErr", fmt::sprintf("%s d_{xy} error (#mum)", s_1).c_str(), 100, 0., 2000.);
520  dz = book<TH1D>("dz", fmt::sprintf("%s d_{z} (#mum)", s_1).c_str(), dzBin_, dzMin_, dzMax_);
521  dzErr = book<TH1D>("dzErr", fmt::sprintf("%s d_{z} error(#mum)", s_1).c_str(), 100, 0., 10000.);
522 
523  phi_pt1 =
524  book<TH1D>("phi_pt1", fmt::sprintf("%s #phi; PV tracks #phi;#tracks", s_1).c_str(), phiBin_, phiMin_, phiMax_);
525  eta_pt1 =
526  book<TH1D>("eta_pt1", fmt::sprintf("%s #eta; PV tracks #eta;#tracks", s_1).c_str(), etaBin_, etaMin_, etaMax_);
527  phi_pt10 =
528  book<TH1D>("phi_pt10", fmt::sprintf("%s #phi; PV tracks #phi;#tracks", s_10).c_str(), phiBin_, phiMin_, phiMax_);
529  eta_pt10 =
530  book<TH1D>("eta_pt10", fmt::sprintf("%s #phi; PV tracks #eta;#tracks", s_10).c_str(), etaBin_, etaMin_, etaMax_);
531 
532  dxyVsPhi_pt1 = book<TProfile>("dxyVsPhi_pt1",
533  fmt::sprintf("%s d_{xy} (#mum) VS track #phi", s_1).c_str(),
534  phiBin_,
535  phiMin_,
536  phiMax_,
537  dxyMin_,
538  dxyMax_);
539  dxyVsPhi_pt1->SetXTitle("PV track (p_{T} > 1 GeV) #phi");
540  dxyVsPhi_pt1->SetYTitle("PV track (p_{T} > 1 GeV) d_{xy} (#mum)");
541 
542  dzVsPhi_pt1 = book<TProfile>("dzVsPhi_pt1",
543  fmt::sprintf("%s d_{z} (#mum) VS track #phi", s_1).c_str(),
544  phiBin_,
545  phiMin_,
546  phiMax_,
547  dzMin_,
548  dzMax_);
549  dzVsPhi_pt1->SetXTitle("PV track (p_{T} > 1 GeV) #phi");
550  dzVsPhi_pt1->SetYTitle("PV track (p_{T} > 1 GeV) d_{z} (#mum)");
551 
552  dxyVsEta_pt1 = book<TProfile>("dxyVsEta_pt1",
553  fmt::sprintf("%s d_{xy} (#mum) VS track #eta", s_1).c_str(),
554  etaBin_,
555  etaMin_,
556  etaMax_,
557  dxyMin_,
558  dxyMax_);
559  dxyVsEta_pt1->SetXTitle("PV track (p_{T} > 1 GeV) #eta");
560  dxyVsEta_pt1->SetYTitle("PV track (p_{T} > 1 GeV) d_{xy} (#mum)");
561 
562  dzVsEta_pt1 = book<TProfile>("dzVsEta_pt1",
563  fmt::sprintf("%s d_{z} (#mum) VS track #eta", s_1).c_str(),
564  etaBin_,
565  etaMin_,
566  etaMax_,
567  dzMin_,
568  dzMax_);
569  dzVsEta_pt1->SetXTitle("PV track (p_{T} > 1 GeV) #eta");
570  dzVsEta_pt1->SetYTitle("PV track (p_{T} > 1 GeV) d_{z} (#mum)");
571 
572  dxyVsEtaVsPhi_pt1 = book<TProfile2D>("dxyVsEtaVsPhi_pt1",
573  fmt::sprintf("%s d_{xy} (#mum) VS track #eta VS track #phi", s_1).c_str(),
574  etaBin2D_,
575  etaMin_,
576  etaMax_,
577  phiBin2D_,
578  phiMin_,
579  phiMax_,
580  dxyMin_,
581  dxyMax_);
582  dxyVsEtaVsPhi_pt1->SetXTitle("PV track (p_{T} > 1 GeV) #eta");
583  dxyVsEtaVsPhi_pt1->SetYTitle("PV track (p_{T} > 1 GeV) #phi");
584  dxyVsEtaVsPhi_pt1->SetZTitle("PV track (p_{T} > 1 GeV) d_{xy} (#mum)");
585 
586  dzVsEtaVsPhi_pt1 = book<TProfile2D>("dzVsEtaVsPhi_pt1",
587  fmt::sprintf("%s d_{z} (#mum) VS track #eta VS track #phi", s_1).c_str(),
588  etaBin2D_,
589  etaMin_,
590  etaMax_,
591  phiBin2D_,
592  phiMin_,
593  phiMax_,
594  dzMin_,
595  dzMax_);
596  dzVsEtaVsPhi_pt1->SetXTitle("PV track (p_{T} > 1 GeV) #eta");
597  dzVsEtaVsPhi_pt1->SetYTitle("PV track (p_{T} > 1 GeV) #phi");
598  dzVsEtaVsPhi_pt1->SetZTitle("PV track (p_{T} > 1 GeV) d_{z} (#mum)");
599 
600  dxyVsPhi_pt10 = book<TProfile>("dxyVsPhi_pt10",
601  fmt::sprintf("%s d_{xy} (#mum) VS track #phi", s_10).c_str(),
602  phiBin_,
603  phiMin_,
604  phiMax_,
605  dxyMin_,
606  dxyMax_);
607  dxyVsPhi_pt10->SetXTitle("PV track (p_{T} > 10 GeV) #phi");
608  dxyVsPhi_pt10->SetYTitle("PV track (p_{T} > 10 GeV) d_{xy} (#mum)");
609 
610  dzVsPhi_pt10 = book<TProfile>("dzVsPhi_pt10",
611  fmt::sprintf("%s d_{z} (#mum) VS track #phi", s_10).c_str(),
612  phiBin_,
613  phiMin_,
614  phiMax_,
615  dzMin_,
616  dzMax_);
617  dzVsPhi_pt10->SetXTitle("PV track (p_{T} > 10 GeV) #phi");
618  dzVsPhi_pt10->SetYTitle("PV track (p_{T} > 10 GeV) d_{z} (#mum)");
619 
620  dxyVsEta_pt10 = book<TProfile>("dxyVsEta_pt10",
621  fmt::sprintf("%s d_{xy} (#mum) VS track #eta", s_10).c_str(),
622  etaBin_,
623  etaMin_,
624  etaMax_,
625  dxyMin_,
626  dxyMax_);
627  dxyVsEta_pt10->SetXTitle("PV track (p_{T} > 10 GeV) #eta");
628  dxyVsEta_pt10->SetYTitle("PV track (p_{T} > 10 GeV) d_{xy} (#mum)");
629 
630  dzVsEta_pt10 = book<TProfile>("dzVsEta_pt10",
631  fmt::sprintf("%s d_{z} (#mum) VS track #eta", s_10).c_str(),
632  etaBin_,
633  etaMin_,
634  etaMax_,
635  dzMin_,
636  dzMax_);
637  dzVsEta_pt10->SetXTitle("PV track (p_{T} > 10 GeV) #eta");
638  dzVsEta_pt10->SetYTitle("PV track (p_{T} > 10 GeV) d_{z} (#mum)");
639 
640  dxyVsEtaVsPhi_pt10 = book<TProfile2D>("dxyVsEtaVsPhi_pt10",
641  fmt::sprintf("%s d_{xy} (#mum) VS track #eta VS track #phi", s_10).c_str(),
642  etaBin2D_,
643  etaMin_,
644  etaMax_,
645  phiBin2D_,
646  phiMin_,
647  phiMax_,
648  dxyMin_,
649  dxyMax_);
650  dxyVsEtaVsPhi_pt10->SetXTitle("PV track (p_{T} > 10 GeV) #eta");
651  dxyVsEtaVsPhi_pt10->SetYTitle("PV track (p_{T} > 10 GeV) #phi");
652  dxyVsEtaVsPhi_pt10->SetZTitle("PV track (p_{T} > 10 GeV) d_{xy} (#mum)");
653 
654  dzVsEtaVsPhi_pt10 = book<TProfile2D>("dzVsEtaVsPhi_pt10",
655  fmt::sprintf("%s d_{z} (#mum) VS track #eta VS track #phi", s_10).c_str(),
656  etaBin2D_,
657  etaMin_,
658  etaMax_,
659  phiBin2D_,
660  phiMin_,
661  phiMax_,
662  dzMin_,
663  dzMax_);
664  dzVsEtaVsPhi_pt10->SetXTitle("PV track (p_{T} > 10 GeV) #eta");
665  dzVsEtaVsPhi_pt10->SetYTitle("PV track (p_{T} > 10 GeV) #phi");
666  dzVsEtaVsPhi_pt10->SetZTitle("PV track (p_{T} > 10 GeV) d_{z} (#mum)");
667 }
668 
669 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
672  desc.add<int>("ndof", 4);
673  desc.add<edm::InputTag>("vertexLabel", edm::InputTag("offlinePrimaryVertices"));
674  desc.add<edm::InputTag>("beamSpotLabel", edm::InputTag("offlineBeamSpot"));
675  desc.add<double>("Xpos", 0.1);
676  desc.add<double>("Ypos", 0.0);
677  desc.add<int>("TkSizeBin", 100);
678  desc.add<double>("TkSizeMin", 499.5);
679  desc.add<double>("TkSizeMax", -0.5);
680  desc.add<int>("DxyBin", 100);
681  desc.add<double>("DxyMin", 5000.);
682  desc.add<double>("DxyMax", -5000.);
683  desc.add<int>("DzBin", 100);
684  desc.add<double>("DzMin", -2000.0);
685  desc.add<double>("DzMax", 2000.0);
686  desc.add<int>("PhiBin", 32);
687  desc.add<int>("PhiBin2D", 12);
688  desc.add<double>("PhiMin", -M_PI);
689  desc.add<double>("PhiMax", M_PI);
690  desc.add<int>("EtaBin", 26);
691  desc.add<int>("EtaBin2D", 8);
692  desc.add<double>("EtaMin", -2.7);
693  desc.add<double>("EtaMax", 2.7);
694  descriptions.addWithDefaultLabel(desc);
695 }
696 
697 //define this as a plug-in
static const std::string kSharedResource
Definition: TFileService.h:76
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
T * book(const Args &...args) const
const edm::EDGetTokenT< reco::VertexCollection > vertexToken_
const edm::EDGetTokenT< reco::BeamSpot > beamspotToken_
T w() const
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
const edm::EDGetTokenT< VertexScore > scoreToken_
~GeneralPurposeVertexAnalyzer() override=default
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
Definition: weight.py:1
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
void vertexPlots(const reco::Vertex &v, const reco::BeamSpot &beamSpot, int i)
int iEvent
Definition: GenABIO.cc:224
Definition: TTTypes.h:54
void pvTracksPlots(const reco::Vertex &v)
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
T sqrt(T t)
Definition: SSEVec.h:19
double f[11][100]
float ChiSquaredProbability(double chiSquared, double nrDOF)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
#define M_PI
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
fixed size matrix
HLT enums.
GeneralPurposeVertexAnalyzer(const edm::ParameterSet &)
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
Log< level::Warning, false > LogWarning
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
long double T
void analyze(const edm::Event &, const edm::EventSetup &) override
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:38