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