CMS 3D CMS Logo

DiMuonVertexValidation.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Alignment/OfflineValidation
4 // Class: DiMuonVertexValidation
5 //
11 //
12 // Original Author: Marco Musich
13 // Created: Wed, 21 Apr 2021 09:06:25 GMT
14 //
15 //
16 
17 // system include files
18 #include <memory>
19 #include <utility>
20 
21 // user include files
28 
29 // muons
33 
34 // utils
37 #include "TLorentzVector.h"
38 
39 // tracks
46 
47 // vertices
52 
53 // ROOT
54 #include "TH1F.h"
55 #include "TH2F.h"
56 
57 //#define LogDebug(X) std::cout << X <<
58 
59 //
60 // class declaration
61 //
62 
63 class DiMuonVertexValidation : public edm::one::EDAnalyzer<edm::one::SharedResources> {
64 public:
66  ~DiMuonVertexValidation() override;
67 
68  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
69 
70 private:
71  void beginJob() override;
72  void analyze(const edm::Event&, const edm::EventSetup&) override;
74  void endJob() override;
75 
76  // ----------member data ---------------------------
78 
79  const bool useReco_;
80  const bool useClosestVertex_;
81  std::vector<double> pTthresholds_;
82  float maxSVdist_;
83 
84  // plot configurations
85 
94 
95  // control plots
96 
97  TH1F* hSVProb_;
98  TH1F* hSVDist_;
99  TH1F* hSVDistSig_;
100  TH1F* hSVDist3D_;
102 
103  TH1F* hCosPhi_;
104  TH1F* hCosPhi3D_;
105  TH1F* hCosPhiInv_;
107 
108  TH1F* hInvMass_;
110 
111  TH1F* hCutFlow_;
112 
113  // 2D maps
114 
123 
125 
126  edm::EDGetTokenT<reco::TrackCollection> tracksToken_; //used to select what tracks to read from configuration file
127  edm::EDGetTokenT<reco::VertexCollection> vertexToken_; //used to select what vertices to read from configuration file
128 
129  // either on or the other!
130  edm::EDGetTokenT<reco::MuonCollection> muonsToken_; //used to select what tracks to read from configuration file
132  alcaRecoToken_; //used to select what muon tracks to read from configuration file
133 };
134 
135 //
136 // constants, enums and typedefs
137 //
138 
139 static constexpr float cmToum = 10e4;
140 static constexpr float mumass2 = 0.105658367 * 0.105658367; //mu mass squared (GeV^2/c^4)
141 
142 //
143 // static data member definitions
144 //
145 
146 //
147 // constructors and destructor
148 //
150  : useReco_(iConfig.getParameter<bool>("useReco")),
151  useClosestVertex_(iConfig.getParameter<bool>("useClosestVertex")),
152  pTthresholds_(iConfig.getParameter<std::vector<double>>("pTThresholds")),
153  maxSVdist_(iConfig.getParameter<double>("maxSVdist")),
154  CosPhiConfiguration_(iConfig.getParameter<edm::ParameterSet>("CosPhiConfig")),
155  CosPhi3DConfiguration_(iConfig.getParameter<edm::ParameterSet>("CosPhi3DConfig")),
156  VtxProbConfiguration_(iConfig.getParameter<edm::ParameterSet>("VtxProbConfig")),
157  VtxDistConfiguration_(iConfig.getParameter<edm::ParameterSet>("VtxDistConfig")),
158  VtxDist3DConfiguration_(iConfig.getParameter<edm::ParameterSet>("VtxDist3DConfig")),
159  VtxDistSigConfiguration_(iConfig.getParameter<edm::ParameterSet>("VtxDistSigConfig")),
160  VtxDist3DSigConfiguration_(iConfig.getParameter<edm::ParameterSet>("VtxDist3DSigConfig")),
161  DiMuMassConfiguration_(iConfig.getParameter<edm::ParameterSet>("DiMuMassConfig")),
162  ttbESToken_(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))),
163  tracksToken_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("tracks"))),
164  vertexToken_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"))) {
165  if (useReco_) {
166  muonsToken_ = mayConsume<reco::MuonCollection>(iConfig.getParameter<edm::InputTag>("muons"));
167  } else {
168  alcaRecoToken_ = mayConsume<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("muonTracks"));
169  }
170 
171  usesResource(TFileService::kSharedResource);
172 
173  // sort the vector of thresholds
174  std::sort(pTthresholds_.begin(), pTthresholds_.end(), [](const double& lhs, const double& rhs) { return lhs > rhs; });
175 
176  edm::LogInfo("DiMuonVertexValidation") << __FUNCTION__;
177  for (const auto& thr : pTthresholds_) {
178  edm::LogInfo("DiMuonVertexValidation") << " Threshold: " << thr << " ";
179  }
180  edm::LogInfo("DiMuonVertexValidation") << "Max SV distance: " << maxSVdist_ << " ";
181 }
182 
184 
185 //
186 // member functions
187 //
188 
189 // ------------ method called for each event ------------
191  using namespace edm;
192 
194 
195  // the di-muon tracks
196  std::vector<const reco::Track*> myTracks;
197 
198  // if we have to start from scratch from RECO data-tier
199  if (useReco_) {
200  // select the good muons
201  std::vector<const reco::Muon*> myGoodMuonVector;
202  for (const auto& muon : iEvent.get(muonsToken_)) {
203  const reco::TrackRef t = muon.innerTrack();
204  if (!t.isNull()) {
205  if (t->quality(reco::TrackBase::highPurity)) {
206  if (t->chi2() / t->ndof() <= 2.5 && t->numberOfValidHits() >= 5 &&
207  t->hitPattern().numberOfValidPixelHits() >= 2 && t->quality(reco::TrackBase::highPurity))
208  myGoodMuonVector.emplace_back(&muon);
209  }
210  }
211  }
212 
213  LogDebug("DiMuonVertexValidation") << "myGoodMuonVector size: " << myGoodMuonVector.size() << std::endl;
214  std::sort(myGoodMuonVector.begin(), myGoodMuonVector.end(), [](const reco::Muon*& lhs, const reco::Muon*& rhs) {
215  return lhs->pt() > rhs->pt();
216  });
217 
218  // just check the ordering
219  for (const auto& muon : myGoodMuonVector) {
220  LogDebug("DiMuonVertexValidation") << "pT: " << muon->pt() << " ";
221  }
222  LogDebug("DiMuonVertexValidation") << std::endl;
223 
224  // reject if there's no Z
225  if (myGoodMuonVector.size() < 2)
226  return;
227 
229 
230  if ((myGoodMuonVector[0]->pt()) < pTthresholds_[0] || (myGoodMuonVector[1]->pt() < pTthresholds_[1]))
231  return;
232 
235 
236  if (myGoodMuonVector[0]->charge() * myGoodMuonVector[1]->charge() > 0)
237  return;
238 
239  const auto& m1 = myGoodMuonVector[1]->p4();
240  const auto& m0 = myGoodMuonVector[0]->p4();
241  const auto& mother = m0 + m1;
242 
243  float invMass = mother.M();
244  hInvMass_->Fill(invMass);
245 
246  // just copy the top two muons
247  std::vector<const reco::Muon*> theZMuonVector;
248  theZMuonVector.reserve(2);
249  theZMuonVector.emplace_back(myGoodMuonVector[1]);
250  theZMuonVector.emplace_back(myGoodMuonVector[0]);
251 
252  // do the matching of Z muons with inner tracks
253  unsigned int i = 0;
254  for (const auto& muon : theZMuonVector) {
255  i++;
256  float minD = 1000.;
257  const reco::Track* theMatch = nullptr;
258  for (const auto& track : iEvent.get(tracksToken_)) {
259  float D = ::deltaR(muon->eta(), muon->phi(), track.eta(), track.phi());
260  if (D < minD) {
261  minD = D;
262  theMatch = &track;
263  }
264  }
265  LogDebug("DiMuonVertexValidation") << "pushing new track: " << i << std::endl;
266  myTracks.emplace_back(theMatch);
267  }
268  } else {
269  // we start directly with the pre-selected ALCARECO tracks
270  for (const auto& muon : iEvent.get(alcaRecoToken_)) {
271  myTracks.emplace_back(&muon);
272  }
273  }
274 
275 #ifdef EDM_ML_DEBUG
276  for (const auto& track : myTracks) {
277  edm::LogVerbatim("DiMuonVertexValidation") << __PRETTY_FUNCTION__ << " pT: " << track->pt() << " GeV"
278  << " , pT error: " << track->ptError() << " GeV"
279  << " , eta: " << track->eta() << " , phi: " << track->phi() << std::endl;
280  }
281 #endif
282 
283  LogDebug("DiMuonVertexValidation") << "selected tracks: " << myTracks.size() << std::endl;
284 
285  const TransientTrackBuilder* theB = &iSetup.getData(ttbESToken_);
286  TransientVertex aTransientVertex;
287  std::vector<reco::TransientTrack> tks;
288 
289  if (myTracks.size() != 2)
290  return;
291 
292  const auto& t1 = myTracks[1]->momentum();
293  const auto& t0 = myTracks[0]->momentum();
294  const auto& ditrack = t1 + t0;
295 
296  const auto& tplus = myTracks[0]->charge() > 0 ? myTracks[0] : myTracks[1];
297  const auto& tminus = myTracks[0]->charge() < 0 ? myTracks[0] : myTracks[1];
298 
299  TLorentzVector p4_tplus(tplus->px(), tplus->py(), tplus->pz(), sqrt((tplus->p() * tplus->p()) + mumass2));
300  TLorentzVector p4_tminus(tminus->px(), tminus->py(), tminus->pz(), sqrt((tminus->p() * tminus->p()) + mumass2));
301 
302 #ifdef EDM_ML_DEBUG
303  // Define a lambda function to convert TLorentzVector to a string
304  auto tLorentzVectorToString = [](const TLorentzVector& vector) {
305  return std::to_string(vector.Px()) + " " + std::to_string(vector.Py()) + " " + std::to_string(vector.Pz()) + " " +
306  std::to_string(vector.E());
307  };
308 
309  edm::LogVerbatim("DiMuonVertexValidation") << "mu+" << tLorentzVectorToString(p4_tplus) << std::endl;
310  edm::LogVerbatim("DiMuonVertexValidation") << "mu-" << tLorentzVectorToString(p4_tminus) << std::endl;
311 #endif
312 
313  const auto& Zp4 = p4_tplus + p4_tminus;
314  float track_invMass = Zp4.M();
315  hTrackInvMass_->Fill(track_invMass);
316 
317  // creat the pair of TLorentVectors used to make the plos
318  std::pair<TLorentzVector, TLorentzVector> tktk_p4 = std::make_pair(p4_tplus, p4_tminus);
319 
320  // fill the z->mm mass plots
321  ZMassPlots.fillPlots(track_invMass, tktk_p4);
322 
323  math::XYZPoint ZpT(ditrack.x(), ditrack.y(), 0);
324  math::XYZPoint Zp(ditrack.x(), ditrack.y(), ditrack.z());
325 
326  for (const auto& track : myTracks) {
327  reco::TransientTrack trajectory = theB->build(track);
328  tks.push_back(trajectory);
329  }
330 
331  KalmanVertexFitter kalman(true);
332  aTransientVertex = kalman.vertex(tks);
333 
334  double SVProb = TMath::Prob(aTransientVertex.totalChiSquared(), (int)aTransientVertex.degreesOfFreedom());
335 
336  LogDebug("DiMuonVertexValidation") << " vertex prob: " << SVProb << std::endl;
337 
338  hSVProb_->Fill(SVProb);
339 
340  if (!aTransientVertex.isValid())
341  return;
342 
344 
345  // fill the VtxProb plots
346  VtxProbPlots.fillPlots(SVProb, tktk_p4);
347 
348  math::XYZPoint MainVertex(0, 0, 0);
349  const reco::Vertex* theClosestVertex = nullptr;
350  // get collection of reconstructed vertices from event
351  edm::Handle<reco::VertexCollection> vertexHandle = iEvent.getHandle(vertexToken_);
352  if (vertexHandle.isValid()) {
353  const reco::VertexCollection* vertices = vertexHandle.product();
354  theClosestVertex = this->findClosestVertex(aTransientVertex, vertices);
355  } else {
356  edm::LogWarning("DiMuonVertexMonitor") << "invalid vertex collection encountered Skipping event!";
357  return;
358  }
359 
360  reco::Vertex TheMainVertex;
361  if (!useClosestVertex_ || theClosestVertex == nullptr) {
362  // if the closest vertex is not available, or explicitly not chosen
363  TheMainVertex = vertexHandle.product()->front();
364  } else {
365  TheMainVertex = *theClosestVertex;
366  }
367 
368  MainVertex.SetXYZ(TheMainVertex.position().x(), TheMainVertex.position().y(), TheMainVertex.position().z());
369  const math::XYZPoint myVertex(
370  aTransientVertex.position().x(), aTransientVertex.position().y(), aTransientVertex.position().z());
371  const math::XYZPoint deltaVtx(
372  MainVertex.x() - myVertex.x(), MainVertex.y() - myVertex.y(), MainVertex.z() - myVertex.z());
373 
374 #ifdef EDM_ML_DEBUG
375  edm::LogVerbatim("DiMuonVertexValidation")
376  << "mm vertex position:" << aTransientVertex.position().x() << "," << aTransientVertex.position().y() << ","
377  << aTransientVertex.position().z();
378 
379  edm::LogVerbatim("DiMuonVertexValidation") << "main vertex position:" << TheMainVertex.position().x() << ","
380  << TheMainVertex.position().y() << "," << TheMainVertex.position().z();
381 #endif
382 
383  if (TheMainVertex.isValid()) {
384  // Z Vertex distance in the xy plane
385 
386  VertexDistanceXY vertTool;
387  double distance = vertTool.distance(aTransientVertex, TheMainVertex).value();
388  double dist_err = vertTool.distance(aTransientVertex, TheMainVertex).error();
389 
390  hSVDist_->Fill(distance * cmToum);
391  hSVDistSig_->Fill(distance / dist_err);
392 
393  // fill the VtxDist plots
395 
396  // fill the VtxDisSig plots
397  VtxDistSigPlots.fillPlots(distance / dist_err, tktk_p4);
398 
399  // Z Vertex distance in 3D
400 
401  VertexDistance3D vertTool3D;
402  double distance3D = vertTool3D.distance(aTransientVertex, TheMainVertex).value();
403  double dist3D_err = vertTool3D.distance(aTransientVertex, TheMainVertex).error();
404 
405  hSVDist3D_->Fill(distance3D * cmToum);
406  hSVDist3DSig_->Fill(distance3D / dist3D_err);
407 
408  // fill the VtxDist3D plots
409  VtxDist3DPlots.fillPlots(distance3D * cmToum, tktk_p4);
410 
411  // fill the VtxDisSig plots
412  VtxDist3DSigPlots.fillPlots(distance3D / dist3D_err, tktk_p4);
413 
414  LogDebug("DiMuonVertexValidation") << "distance: " << distance << "+/-" << dist_err << std::endl;
415  // cut on the PV - SV distance
416  if (distance * cmToum < maxSVdist_) {
418 
419  double cosphi = (ZpT.x() * deltaVtx.x() + ZpT.y() * deltaVtx.y()) /
420  (sqrt(ZpT.x() * ZpT.x() + ZpT.y() * ZpT.y()) *
421  sqrt(deltaVtx.x() * deltaVtx.x() + deltaVtx.y() * deltaVtx.y()));
422 
423  double cosphi3D = (Zp.x() * deltaVtx.x() + Zp.y() * deltaVtx.y() + Zp.z() * deltaVtx.z()) /
424  (sqrt(Zp.x() * Zp.x() + Zp.y() * Zp.y() + Zp.z() * Zp.z()) *
425  sqrt(deltaVtx.x() * deltaVtx.x() + deltaVtx.y() * deltaVtx.y() + deltaVtx.z() * deltaVtx.z()));
426 
427  LogDebug("DiMuonVertexValidation") << "cos(phi) = " << cosphi << std::endl;
428 
429  hCosPhi_->Fill(cosphi);
430  hCosPhi3D_->Fill(cosphi3D);
431 
432 #ifdef EDM_ML_DEBUG
433  edm::LogVerbatim("DiMuonVertexValidation")
434  << "distance " << distance * cmToum << " cosphi3D:" << cosphi3D << std::endl;
435 #endif
436 
437  // fill the cosphi plots
438  CosPhiPlots.fillPlots(cosphi, tktk_p4);
439 
440  // fill the VtxDisSig plots
441  CosPhi3DPlots.fillPlots(cosphi3D, tktk_p4);
442  }
443  }
444 }
445 
446 // ------------ method called once each job just before starting event loop ------------
449 
450  // clang-format off
451  TH1F::SetDefaultSumw2(kTRUE);
452  hSVProb_ = fs->make<TH1F>("VtxProb", ";ZV vertex probability;N(#mu#mu pairs)", 100, 0., 1.);
453 
454  hSVDist_ = fs->make<TH1F>("VtxDist", ";PV-ZV xy distance [#mum];N(#mu#mu pairs)", 100, 0., 300.);
455  hSVDistSig_ = fs->make<TH1F>("VtxDistSig", ";PV-ZV xy distance signficance;N(#mu#mu pairs)", 100, 0., 5.);
456 
457  hSVDist3D_ = fs->make<TH1F>("VtxDist3D", ";PV-ZV 3D distance [#mum];N(#mu#mu pairs)", 100, 0., 300.);
458  hSVDist3DSig_ = fs->make<TH1F>("VtxDist3DSig", ";PV-ZV 3D distance signficance;N(#mu#mu pairs)", 100, 0., 5.);
459 
460  hInvMass_ = fs->make<TH1F>("InvMass", ";M(#mu#mu) [GeV];N(#mu#mu pairs)", 70., 50., 120.);
461  hTrackInvMass_ = fs->make<TH1F>("TkTkInvMass", ";M(tk,tk) [GeV];N(tk tk pairs)", 70., 50., 120.);
462 
463  hCosPhi_ = fs->make<TH1F>("CosPhi", ";cos(#phi_{xy});N(#mu#mu pairs)", 50, -1., 1.);
464  hCosPhi3D_ = fs->make<TH1F>("CosPhi3D", ";cos(#phi_{3D});N(#mu#mu pairs)", 50, -1., 1.);
465 
466  hCosPhiInv_ = fs->make<TH1F>("CosPhiInv", ";inverted cos(#phi_{xy});N(#mu#mu pairs)", 50, -1., 1.);
467  hCosPhiInv3D_ = fs->make<TH1F>("CosPhiInv3D", ";inverted cos(#phi_{3D});N(#mu#mu pairs)", 50, -1., 1.);
468  // clang-format on
469 
470  // 2D Maps
471 
472  TFileDirectory dirCosPhi = fs->mkdir("CosPhiPlots");
474 
475  TFileDirectory dirCosPhi3D = fs->mkdir("CosPhi3DPlots");
477 
478  TFileDirectory dirVtxProb = fs->mkdir("VtxProbPlots");
480 
481  TFileDirectory dirVtxDist = fs->mkdir("VtxDistPlots");
483 
484  TFileDirectory dirVtxDist3D = fs->mkdir("VtxDist3DPlots");
486 
487  TFileDirectory dirVtxDistSig = fs->mkdir("VtxDistSigPlots");
489 
490  TFileDirectory dirVtxDist3DSig = fs->mkdir("VtxDist3DSigPlots");
492 
493  TFileDirectory dirInvariantMass = fs->mkdir("InvariantMassPlots");
495 
496  // cut flow
497 
498  hCutFlow_ = fs->make<TH1F>("hCutFlow", "cut flow;cut step;events left", 6, -0.5, 5.5);
499  std::string names[6] = {"Total", "Mult.", ">pT", "<eta", "hasVtx", "VtxDist"};
500  for (unsigned int i = 0; i < 6; i++) {
501  hCutFlow_->GetXaxis()->SetBinLabel(i + 1, names[i].c_str());
502  }
503 
504  myCounts.zeroAll();
505 }
506 
507 // ------------ method called once each job just after ending the event loop ------------
510 
511  hCutFlow_->SetBinContent(1, myCounts.eventsTotal);
512  hCutFlow_->SetBinContent(2, myCounts.eventsAfterMult);
513  hCutFlow_->SetBinContent(3, myCounts.eventsAfterPt);
514  hCutFlow_->SetBinContent(4, myCounts.eventsAfterEta);
515  hCutFlow_->SetBinContent(5, myCounts.eventsAfterVtx);
516  hCutFlow_->SetBinContent(6, myCounts.eventsAfterDist);
517 
518  TH1F::SetDefaultSumw2(kTRUE);
519  const unsigned int nBinsX = hCosPhi_->GetNbinsX();
520  for (unsigned int i = 1; i <= nBinsX; i++) {
521  //float binContent = hCosPhi_->GetBinContent(i);
522  float invertedBinContent = hCosPhi_->GetBinContent(nBinsX + 1 - i);
523  float invertedBinError = hCosPhi_->GetBinError(nBinsX + 1 - i);
524  hCosPhiInv_->SetBinContent(i, invertedBinContent);
525  hCosPhiInv_->SetBinError(i, invertedBinError);
526 
527  //float binContent3D = hCosPhi3D_->GetBinContent(i);
528  float invertedBinContent3D = hCosPhi3D_->GetBinContent(nBinsX + 1 - i);
529  float invertedBinError3D = hCosPhi3D_->GetBinError(nBinsX + 1 - i);
530  hCosPhiInv3D_->SetBinContent(i, invertedBinContent3D);
531  hCosPhiInv3D_->SetBinError(i, invertedBinError3D);
532  }
533 }
534 
535 // compute the closest vertex to di-lepton ------------------------------------
537  const reco::VertexCollection* vertices) const {
538  reco::Vertex* defaultVtx = nullptr;
539 
540  if (!aTransVtx.isValid())
541  return defaultVtx;
542 
543  // find the closest vertex to the secondary vertex in 3D
544  VertexDistance3D vertTool3D;
545  float minD = 9999.;
546  int closestVtxIndex = 0;
547  int counter = 0;
548  for (const auto& vtx : *vertices) {
549  double dist3D = vertTool3D.distance(aTransVtx, vtx).value();
550  if (dist3D < minD) {
551  minD = dist3D;
552  closestVtxIndex = counter;
553  }
554  counter++;
555  }
556 
557  if ((*vertices).at(closestVtxIndex).isValid()) {
558  return &(vertices->at(closestVtxIndex));
559  } else {
560  return defaultVtx;
561  }
562 }
563 
564 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
567  desc.ifValue(edm::ParameterDescription<bool>("useReco", true, true),
568  true >> edm::ParameterDescription<edm::InputTag>("muons", edm::InputTag("muons"), true) or
570  "muonTracks", edm::InputTag("ALCARECOTkAlDiMuon"), true))
571  ->setComment("If useReco is true need to specify the muon tracks, otherwise take the ALCARECO Inner tracks");
572  //desc.add<bool>("useReco",true);
573  //desc.add<edm::InputTag>("muons", edm::InputTag("muons"));
574  //desc.add<edm::InputTag>("muonTracks", edm::InputTag("ALCARECOTkAlDiMuon"));
575  desc.add<edm::InputTag>("tracks", edm::InputTag("generalTracks"));
576  desc.add<edm::InputTag>("vertices", edm::InputTag("offlinePrimaryVertices"));
577  desc.add<std::vector<double>>("pTThresholds", {30., 10.});
578  desc.add<bool>("useClosestVertex", true);
579  desc.add<double>("maxSVdist", 50.);
580 
581  {
582  edm::ParameterSetDescription psDiMuMass;
583  psDiMuMass.add<std::string>("name", "DiMuMass");
584  psDiMuMass.add<std::string>("title", "M(#mu#mu)");
585  psDiMuMass.add<std::string>("yUnits", "[GeV]");
586  psDiMuMass.add<int>("NxBins", 24);
587  psDiMuMass.add<int>("NyBins", 50);
588  psDiMuMass.add<double>("ymin", 70.);
589  psDiMuMass.add<double>("ymax", 120.);
590  desc.add<edm::ParameterSetDescription>("DiMuMassConfig", psDiMuMass);
591  }
592  {
594  psCosPhi.add<std::string>("name", "CosPhi");
595  psCosPhi.add<std::string>("title", "cos(#phi_{xy})");
596  psCosPhi.add<std::string>("yUnits", "");
597  psCosPhi.add<int>("NxBins", 50);
598  psCosPhi.add<int>("NyBins", 50);
599  psCosPhi.add<double>("ymin", -1.);
600  psCosPhi.add<double>("ymax", 1.);
601  desc.add<edm::ParameterSetDescription>("CosPhiConfig", psCosPhi);
602  }
603  {
604  edm::ParameterSetDescription psCosPhi3D;
605  psCosPhi3D.add<std::string>("name", "CosPhi3D");
606  psCosPhi3D.add<std::string>("title", "cos(#phi_{3D})");
607  psCosPhi3D.add<std::string>("yUnits", "");
608  psCosPhi3D.add<int>("NxBins", 50);
609  psCosPhi3D.add<int>("NyBins", 50);
610  psCosPhi3D.add<double>("ymin", -1.);
611  psCosPhi3D.add<double>("ymax", 1.);
612  desc.add<edm::ParameterSetDescription>("CosPhi3DConfig", psCosPhi3D);
613  }
614  {
616  psVtxProb.add<std::string>("name", "VtxProb");
617  psVtxProb.add<std::string>("title", "Prob(#chi^{2}_{SV})");
618  psVtxProb.add<std::string>("yUnits", "");
619  psVtxProb.add<int>("NxBins", 50);
620  psVtxProb.add<int>("NyBins", 50);
621  psVtxProb.add<double>("ymin", 0);
622  psVtxProb.add<double>("ymax", 1.);
623  desc.add<edm::ParameterSetDescription>("VtxProbConfig", psVtxProb);
624  }
625  {
627  psVtxDist.add<std::string>("name", "VtxDist");
628  psVtxDist.add<std::string>("title", "d_{xy}(PV,SV)");
629  psVtxDist.add<std::string>("yUnits", "[#mum]");
630  psVtxDist.add<int>("NxBins", 50);
631  psVtxDist.add<int>("NyBins", 100);
632  psVtxDist.add<double>("ymin", 0);
633  psVtxDist.add<double>("ymax", 300.);
634  desc.add<edm::ParameterSetDescription>("VtxDistConfig", psVtxDist);
635  }
636  {
637  edm::ParameterSetDescription psVtxDist3D;
638  psVtxDist3D.add<std::string>("name", "VtxDist3D");
639  psVtxDist3D.add<std::string>("title", "d_{3D}(PV,SV)");
640  psVtxDist3D.add<std::string>("yUnits", "[#mum]");
641  psVtxDist3D.add<int>("NxBins", 50);
642  psVtxDist3D.add<int>("NyBins", 250);
643  psVtxDist3D.add<double>("ymin", 0);
644  psVtxDist3D.add<double>("ymax", 500.);
645  desc.add<edm::ParameterSetDescription>("VtxDist3DConfig", psVtxDist3D);
646  }
647  {
648  edm::ParameterSetDescription psVtxDistSig;
649  psVtxDistSig.add<std::string>("name", "VtxDistSig");
650  psVtxDistSig.add<std::string>("title", "d_{xy}(PV,SV)/#sigma_{dxy}(PV,SV)");
651  psVtxDistSig.add<std::string>("yUnits", "");
652  psVtxDistSig.add<int>("NxBins", 50);
653  psVtxDistSig.add<int>("NyBins", 100);
654  psVtxDistSig.add<double>("ymin", 0);
655  psVtxDistSig.add<double>("ymax", 5.);
656  desc.add<edm::ParameterSetDescription>("VtxDistSigConfig", psVtxDistSig);
657  }
658  {
659  edm::ParameterSetDescription psVtxDist3DSig;
660  psVtxDist3DSig.add<std::string>("name", "VtxDist3DSig");
661  psVtxDist3DSig.add<std::string>("title", "d_{3D}(PV,SV)/#sigma_{d3D}(PV,SV)");
662  psVtxDist3DSig.add<std::string>("yUnits", "");
663  psVtxDist3DSig.add<int>("NxBins", 50);
664  psVtxDist3DSig.add<int>("NyBins", 100);
665  psVtxDist3DSig.add<double>("ymin", 0);
666  psVtxDist3DSig.add<double>("ymax", 5.);
667  desc.add<edm::ParameterSetDescription>("VtxDist3DSigConfig", psVtxDist3DSig);
668  }
669 
670  descriptions.addWithDefaultLabel(desc);
671 }
672 
673 //define this as a plug-in
DiMuonVertexValidation(const edm::ParameterSet &)
static const std::string kSharedResource
Definition: TFileService.h:76
Log< level::Info, true > LogVerbatim
~DiMuonVertexValidation() override
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
float totalChiSquared() const
edm::ParameterSet VtxProbConfiguration_
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
double pt() const final
transverse momentum
GlobalPoint position() const
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
DiLeptonHelp::Counts myCounts
T z() const
Definition: PV3DBase.h:61
void analyze(const edm::Event &, const edm::EventSetup &) override
static constexpr float mumass2
const Point & position() const
position
Definition: Vertex.h:127
T const * product() const
Definition: Handle.h:70
DiLeptonHelp::PlotsVsKinematics CosPhiPlots
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
edm::ParameterSet DiMuMassConfiguration_
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::ParameterSet VtxDist3DConfiguration_
DiLeptonHelp::PlotsVsKinematics VtxDist3DPlots
edm::ParameterSet VtxDistConfiguration_
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
float degreesOfFreedom() const
const std::string names[nVars_]
DiLeptonHelp::PlotsVsKinematics VtxDist3DSigPlots
static std::string to_string(const XMLCh *ch)
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const override
reco::TransientTrack build(const reco::Track *p) const
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
int iEvent
Definition: GenABIO.cc:224
edm::ParameterSet VtxDist3DSigConfiguration_
bool isValid() const
T sqrt(T t)
Definition: SSEVec.h:19
DiLeptonHelp::PlotsVsKinematics VtxDistSigPlots
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
static constexpr float cmToum
edm::EDGetTokenT< reco::MuonCollection > muonsToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > ttbESToken_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
Log< level::Info, false > LogInfo
edm::EDGetTokenT< reco::TrackCollection > tracksToken_
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:141
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
edm::ParameterSet VtxDistSigConfiguration_
DiLeptonHelp::PlotsVsKinematics VtxProbPlots
bool isValid() const
Definition: HandleBase.h:70
double value() const
Definition: Measurement1D.h:25
const reco::Vertex * findClosestVertex(const TransientVertex aTransVtx, const reco::VertexCollection *vertices) const
edm::ParameterSet CosPhiConfiguration_
DiLeptonHelp::PlotsVsKinematics CosPhi3DPlots
fixed size matrix
HLT enums.
edm::ParameterSet CosPhi3DConfiguration_
static std::atomic< unsigned int > counter
double error() const
Definition: Measurement1D.h:27
std::vector< double > pTthresholds_
DiLeptonHelp::PlotsVsKinematics VtxDistPlots
void fillPlots(const float val, const std::pair< TLorentzVector, TLorentzVector > &momenta)
void bookFromPSet(const TFileDirectory &fs, const edm::ParameterSet &hpar)
Log< level::Warning, false > LogWarning
DiLeptonHelp::PlotsVsKinematics ZMassPlots
edm::EDGetTokenT< reco::TrackCollection > alcaRecoToken_
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:72
#define LogDebug(id)