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;
73  void endJob() override;
74 
75  // ----------member data ---------------------------
77 
78  bool useReco_;
79  std::vector<double> pTthresholds_;
80  float maxSVdist_;
81 
82  // plot configurations
83 
92 
93  // control plots
94 
95  TH1F* hSVProb_;
96  TH1F* hSVDist_;
97  TH1F* hSVDistSig_;
98  TH1F* hSVDist3D_;
100 
101  TH1F* hCosPhi_;
102  TH1F* hCosPhi3D_;
103  TH1F* hCosPhiInv_;
105 
106  TH1F* hInvMass_;
108 
109  TH1F* hCutFlow_;
110 
111  // 2D maps
112 
121 
123 
124  edm::EDGetTokenT<reco::TrackCollection> tracksToken_; //used to select what tracks to read from configuration file
125  edm::EDGetTokenT<reco::VertexCollection> vertexToken_; //used to select what vertices to read from configuration file
126 
127  // either on or the other!
128  edm::EDGetTokenT<reco::MuonCollection> muonsToken_; //used to select what tracks to read from configuration file
130  alcaRecoToken_; //used to select what muon tracks to read from configuration file
131 };
132 
133 //
134 // constants, enums and typedefs
135 //
136 
137 static constexpr float cmToum = 10e4;
138 static constexpr float mumass2 = 0.105658367 * 0.105658367; //mu mass squared (GeV^2/c^4)
139 
140 //
141 // static data member definitions
142 //
143 
144 //
145 // constructors and destructor
146 //
148  : useReco_(iConfig.getParameter<bool>("useReco")),
149  pTthresholds_(iConfig.getParameter<std::vector<double>>("pTThresholds")),
150  maxSVdist_(iConfig.getParameter<double>("maxSVdist")),
151  CosPhiConfiguration_(iConfig.getParameter<edm::ParameterSet>("CosPhiConfig")),
152  CosPhi3DConfiguration_(iConfig.getParameter<edm::ParameterSet>("CosPhi3DConfig")),
153  VtxProbConfiguration_(iConfig.getParameter<edm::ParameterSet>("VtxProbConfig")),
154  VtxDistConfiguration_(iConfig.getParameter<edm::ParameterSet>("VtxDistConfig")),
155  VtxDist3DConfiguration_(iConfig.getParameter<edm::ParameterSet>("VtxDist3DConfig")),
156  VtxDistSigConfiguration_(iConfig.getParameter<edm::ParameterSet>("VtxDistSigConfig")),
157  VtxDist3DSigConfiguration_(iConfig.getParameter<edm::ParameterSet>("VtxDist3DSigConfig")),
158  DiMuMassConfiguration_(iConfig.getParameter<edm::ParameterSet>("DiMuMassConfig")),
159  ttbESToken_(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))),
160  tracksToken_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("tracks"))),
161  vertexToken_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"))) {
162  if (useReco_) {
163  muonsToken_ = mayConsume<reco::MuonCollection>(iConfig.getParameter<edm::InputTag>("muons"));
164  } else {
165  alcaRecoToken_ = mayConsume<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("muonTracks"));
166  }
167 
168  usesResource(TFileService::kSharedResource);
169 
170  // sort the vector of thresholds
171  std::sort(pTthresholds_.begin(), pTthresholds_.end(), [](const double& lhs, const double& rhs) { return lhs > rhs; });
172 
173  edm::LogInfo("DiMuonVertexValidation") << __FUNCTION__;
174  for (const auto& thr : pTthresholds_) {
175  edm::LogInfo("DiMuonVertexValidation") << " Threshold: " << thr << " ";
176  }
177  edm::LogInfo("DiMuonVertexValidation") << "Max SV distance: " << maxSVdist_ << " ";
178 }
179 
181 
182 //
183 // member functions
184 //
185 
186 // ------------ method called for each event ------------
188  using namespace edm;
189 
191 
192  // the di-muon tracks
193  std::vector<const reco::Track*> myTracks;
194 
195  // if we have to start from scratch from RECO data-tier
196  if (useReco_) {
197  // select the good muons
198  std::vector<const reco::Muon*> myGoodMuonVector;
199  for (const auto& muon : iEvent.get(muonsToken_)) {
200  const reco::TrackRef t = muon.innerTrack();
201  if (!t.isNull()) {
202  if (t->quality(reco::TrackBase::highPurity)) {
203  if (t->chi2() / t->ndof() <= 2.5 && t->numberOfValidHits() >= 5 &&
204  t->hitPattern().numberOfValidPixelHits() >= 2 && t->quality(reco::TrackBase::highPurity))
205  myGoodMuonVector.emplace_back(&muon);
206  }
207  }
208  }
209 
210  LogDebug("DiMuonVertexValidation") << "myGoodMuonVector size: " << myGoodMuonVector.size() << std::endl;
211  std::sort(myGoodMuonVector.begin(), myGoodMuonVector.end(), [](const reco::Muon*& lhs, const reco::Muon*& rhs) {
212  return lhs->pt() > rhs->pt();
213  });
214 
215  // just check the ordering
216  for (const auto& muon : myGoodMuonVector) {
217  LogDebug("DiMuonVertexValidation") << "pT: " << muon->pt() << " ";
218  }
219  LogDebug("DiMuonVertexValidation") << std::endl;
220 
221  // reject if there's no Z
222  if (myGoodMuonVector.size() < 2)
223  return;
224 
226 
227  if ((myGoodMuonVector[0]->pt()) < pTthresholds_[0] || (myGoodMuonVector[1]->pt() < pTthresholds_[1]))
228  return;
229 
232 
233  if (myGoodMuonVector[0]->charge() * myGoodMuonVector[1]->charge() > 0)
234  return;
235 
236  const auto& m1 = myGoodMuonVector[1]->p4();
237  const auto& m0 = myGoodMuonVector[0]->p4();
238  const auto& mother = m0 + m1;
239 
240  float invMass = mother.M();
241  hInvMass_->Fill(invMass);
242 
243  // just copy the top two muons
244  std::vector<const reco::Muon*> theZMuonVector;
245  theZMuonVector.reserve(2);
246  theZMuonVector.emplace_back(myGoodMuonVector[1]);
247  theZMuonVector.emplace_back(myGoodMuonVector[0]);
248 
249  // do the matching of Z muons with inner tracks
250  unsigned int i = 0;
251  for (const auto& muon : theZMuonVector) {
252  i++;
253  float minD = 1000.;
254  const reco::Track* theMatch = nullptr;
255  for (const auto& track : iEvent.get(tracksToken_)) {
256  float D = ::deltaR(muon->eta(), muon->phi(), track.eta(), track.phi());
257  if (D < minD) {
258  minD = D;
259  theMatch = &track;
260  }
261  }
262  LogDebug("DiMuonVertexValidation") << "pushing new track: " << i << std::endl;
263  myTracks.emplace_back(theMatch);
264  }
265  } else {
266  // we start directly with the pre-selected ALCARECO tracks
267  for (const auto& muon : iEvent.get(alcaRecoToken_)) {
268  myTracks.emplace_back(&muon);
269  }
270  }
271 
272  LogDebug("DiMuonVertexValidation") << "selected tracks: " << myTracks.size() << std::endl;
273 
274  const TransientTrackBuilder* theB = &iSetup.getData(ttbESToken_);
275  TransientVertex aTransientVertex;
276  std::vector<reco::TransientTrack> tks;
277 
278  if (myTracks.size() != 2)
279  return;
280 
281  const auto& t1 = myTracks[1]->momentum();
282  const auto& t0 = myTracks[0]->momentum();
283  const auto& ditrack = t1 + t0;
284 
285  const auto& tplus = myTracks[0]->charge() > 0 ? myTracks[0] : myTracks[1];
286  const auto& tminus = myTracks[0]->charge() < 0 ? myTracks[0] : myTracks[1];
287 
288  TLorentzVector p4_tplus(tplus->px(), tplus->py(), tplus->pz(), sqrt((tplus->p() * tplus->p()) + mumass2));
289  TLorentzVector p4_tminus(tminus->px(), tminus->py(), tminus->pz(), sqrt((tminus->p() * tminus->p()) + mumass2));
290 
291  const auto& Zp4 = p4_tplus + p4_tminus;
292  float track_invMass = Zp4.M();
293  hTrackInvMass_->Fill(track_invMass);
294 
295  // creat the pair of TLorentVectors used to make the plos
296  std::pair<TLorentzVector, TLorentzVector> tktk_p4 = std::make_pair(p4_tplus, p4_tminus);
297 
298  // fill the z->mm mass plots
299  ZMassPlots.fillPlots(track_invMass, tktk_p4);
300 
301  math::XYZPoint ZpT(ditrack.x(), ditrack.y(), 0);
302  math::XYZPoint Zp(ditrack.x(), ditrack.y(), ditrack.z());
303 
304  for (const auto& track : myTracks) {
305  reco::TransientTrack trajectory = theB->build(track);
306  tks.push_back(trajectory);
307  }
308 
309  KalmanVertexFitter kalman(true);
310  aTransientVertex = kalman.vertex(tks);
311 
312  double SVProb = TMath::Prob(aTransientVertex.totalChiSquared(), (int)aTransientVertex.degreesOfFreedom());
313 
314  LogDebug("DiMuonVertexValidation") << " vertex prob: " << SVProb << std::endl;
315 
316  hSVProb_->Fill(SVProb);
317 
318  if (!aTransientVertex.isValid())
319  return;
320 
322 
323  // fill the VtxProb plots
324  VtxProbPlots.fillPlots(SVProb, tktk_p4);
325 
326  // get collection of reconstructed vertices from event
327  edm::Handle<reco::VertexCollection> vertexHandle = iEvent.getHandle(vertexToken_);
328 
329  math::XYZPoint MainVertex(0, 0, 0);
330  reco::Vertex TheMainVertex = vertexHandle.product()->front();
331 
332  if (vertexHandle.isValid()) {
333  const reco::VertexCollection* vertices = vertexHandle.product();
334  if ((*vertices).at(0).isValid()) {
335  auto theMainVtx = (*vertices).at(0);
336  MainVertex.SetXYZ(theMainVtx.position().x(), theMainVtx.position().y(), theMainVtx.position().z());
337  }
338  }
339 
340  const math::XYZPoint myVertex(
341  aTransientVertex.position().x(), aTransientVertex.position().y(), aTransientVertex.position().z());
342  const math::XYZPoint deltaVtx(
343  MainVertex.x() - myVertex.x(), MainVertex.y() - myVertex.y(), MainVertex.z() - myVertex.z());
344 
345  if (TheMainVertex.isValid()) {
346  // Z Vertex distance in the xy plane
347 
348  VertexDistanceXY vertTool;
349  double distance = vertTool.distance(aTransientVertex, TheMainVertex).value();
350  double dist_err = vertTool.distance(aTransientVertex, TheMainVertex).error();
351 
352  hSVDist_->Fill(distance * cmToum);
353  hSVDistSig_->Fill(distance / dist_err);
354 
355  // fill the VtxDist plots
357 
358  // fill the VtxDisSig plots
359  VtxDistSigPlots.fillPlots(distance / dist_err, tktk_p4);
360 
361  // Z Vertex distance in 3D
362 
363  VertexDistance3D vertTool3D;
364  double distance3D = vertTool3D.distance(aTransientVertex, TheMainVertex).value();
365  double dist3D_err = vertTool3D.distance(aTransientVertex, TheMainVertex).error();
366 
367  hSVDist3D_->Fill(distance3D * cmToum);
368  hSVDist3DSig_->Fill(distance3D / dist3D_err);
369 
370  // fill the VtxDist3D plots
371  VtxDist3DPlots.fillPlots(distance3D * cmToum, tktk_p4);
372 
373  // fill the VtxDisSig plots
374  VtxDist3DSigPlots.fillPlots(distance3D / dist3D_err, tktk_p4);
375 
376  LogDebug("DiMuonVertexValidation") << "distance: " << distance << "+/-" << dist_err << std::endl;
377  // cut on the PV - SV distance
378  if (distance * cmToum < maxSVdist_) {
380 
381  double cosphi = (ZpT.x() * deltaVtx.x() + ZpT.y() * deltaVtx.y()) /
382  (sqrt(ZpT.x() * ZpT.x() + ZpT.y() * ZpT.y()) *
383  sqrt(deltaVtx.x() * deltaVtx.x() + deltaVtx.y() * deltaVtx.y()));
384 
385  double cosphi3D = (Zp.x() * deltaVtx.x() + Zp.y() * deltaVtx.y() + Zp.z() * deltaVtx.z()) /
386  (sqrt(Zp.x() * Zp.x() + Zp.y() * Zp.y() + Zp.z() * Zp.z()) *
387  sqrt(deltaVtx.x() * deltaVtx.x() + deltaVtx.y() * deltaVtx.y() + deltaVtx.z() * deltaVtx.z()));
388 
389  LogDebug("DiMuonVertexValidation") << "cos(phi) = " << cosphi << std::endl;
390 
391  hCosPhi_->Fill(cosphi);
392  hCosPhi3D_->Fill(cosphi3D);
393 
394  // fill the cosphi plots
395  CosPhiPlots.fillPlots(cosphi, tktk_p4);
396 
397  // fill the VtxDisSig plots
398  CosPhi3DPlots.fillPlots(cosphi3D, tktk_p4);
399  }
400  }
401 }
402 
403 // ------------ method called once each job just before starting event loop ------------
406 
407  // clang-format off
408  TH1F::SetDefaultSumw2(kTRUE);
409  hSVProb_ = fs->make<TH1F>("VtxProb", ";ZV vertex probability;N(#mu#mu pairs)", 100, 0., 1.);
410 
411  hSVDist_ = fs->make<TH1F>("VtxDist", ";PV-ZV xy distance [#mum];N(#mu#mu pairs)", 100, 0., 300.);
412  hSVDistSig_ = fs->make<TH1F>("VtxDistSig", ";PV-ZV xy distance signficance;N(#mu#mu pairs)", 100, 0., 5.);
413 
414  hSVDist3D_ = fs->make<TH1F>("VtxDist3D", ";PV-ZV 3D distance [#mum];N(#mu#mu pairs)", 100, 0., 300.);
415  hSVDist3DSig_ = fs->make<TH1F>("VtxDist3DSig", ";PV-ZV 3D distance signficance;N(#mu#mu pairs)", 100, 0., 5.);
416 
417  hInvMass_ = fs->make<TH1F>("InvMass", ";M(#mu#mu) [GeV];N(#mu#mu pairs)", 70., 50., 120.);
418  hTrackInvMass_ = fs->make<TH1F>("TkTkInvMass", ";M(tk,tk) [GeV];N(tk tk pairs)", 70., 50., 120.);
419 
420  hCosPhi_ = fs->make<TH1F>("CosPhi", ";cos(#phi_{xy});N(#mu#mu pairs)", 50, -1., 1.);
421  hCosPhi3D_ = fs->make<TH1F>("CosPhi3D", ";cos(#phi_{3D});N(#mu#mu pairs)", 50, -1., 1.);
422 
423  hCosPhiInv_ = fs->make<TH1F>("CosPhiInv", ";inverted cos(#phi_{xy});N(#mu#mu pairs)", 50, -1., 1.);
424  hCosPhiInv3D_ = fs->make<TH1F>("CosPhiInv3D", ";inverted cos(#phi_{3D});N(#mu#mu pairs)", 50, -1., 1.);
425  // clang-format on
426 
427  // 2D Maps
428 
429  TFileDirectory dirCosPhi = fs->mkdir("CosPhiPlots");
431 
432  TFileDirectory dirCosPhi3D = fs->mkdir("CosPhi3DPlots");
434 
435  TFileDirectory dirVtxProb = fs->mkdir("VtxProbPlots");
437 
438  TFileDirectory dirVtxDist = fs->mkdir("VtxDistPlots");
440 
441  TFileDirectory dirVtxDist3D = fs->mkdir("VtxDist3DPlots");
443 
444  TFileDirectory dirVtxDistSig = fs->mkdir("VtxDistSigPlots");
446 
447  TFileDirectory dirVtxDist3DSig = fs->mkdir("VtxDist3DSigPlots");
449 
450  TFileDirectory dirInvariantMass = fs->mkdir("InvariantMassPlots");
452 
453  // cut flow
454 
455  hCutFlow_ = fs->make<TH1F>("hCutFlow", "cut flow;cut step;events left", 6, -0.5, 5.5);
456  std::string names[6] = {"Total", "Mult.", ">pT", "<eta", "hasVtx", "VtxDist"};
457  for (unsigned int i = 0; i < 6; i++) {
458  hCutFlow_->GetXaxis()->SetBinLabel(i + 1, names[i].c_str());
459  }
460 
461  myCounts.zeroAll();
462 }
463 
464 // ------------ method called once each job just after ending the event loop ------------
467 
468  hCutFlow_->SetBinContent(1, myCounts.eventsTotal);
469  hCutFlow_->SetBinContent(2, myCounts.eventsAfterMult);
470  hCutFlow_->SetBinContent(3, myCounts.eventsAfterPt);
471  hCutFlow_->SetBinContent(4, myCounts.eventsAfterEta);
472  hCutFlow_->SetBinContent(5, myCounts.eventsAfterVtx);
473  hCutFlow_->SetBinContent(6, myCounts.eventsAfterDist);
474 
475  TH1F::SetDefaultSumw2(kTRUE);
476  const unsigned int nBinsX = hCosPhi_->GetNbinsX();
477  for (unsigned int i = 1; i <= nBinsX; i++) {
478  //float binContent = hCosPhi_->GetBinContent(i);
479  float invertedBinContent = hCosPhi_->GetBinContent(nBinsX + 1 - i);
480  float invertedBinError = hCosPhi_->GetBinError(nBinsX + 1 - i);
481  hCosPhiInv_->SetBinContent(i, invertedBinContent);
482  hCosPhiInv_->SetBinError(i, invertedBinError);
483 
484  //float binContent3D = hCosPhi3D_->GetBinContent(i);
485  float invertedBinContent3D = hCosPhi3D_->GetBinContent(nBinsX + 1 - i);
486  float invertedBinError3D = hCosPhi3D_->GetBinError(nBinsX + 1 - i);
487  hCosPhiInv3D_->SetBinContent(i, invertedBinContent3D);
488  hCosPhiInv3D_->SetBinError(i, invertedBinError3D);
489  }
490 }
491 
492 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
495  desc.ifValue(edm::ParameterDescription<bool>("useReco", true, true),
496  true >> edm::ParameterDescription<edm::InputTag>("muons", edm::InputTag("muons"), true) or
498  "muonTracks", edm::InputTag("ALCARECOTkAlDiMuon"), true))
499  ->setComment("If useReco is true need to specify the muon tracks, otherwise take the ALCARECO Inner tracks");
500  //desc.add<bool>("useReco",true);
501  //desc.add<edm::InputTag>("muons", edm::InputTag("muons"));
502  //desc.add<edm::InputTag>("muonTracks", edm::InputTag("ALCARECOTkAlDiMuon"));
503  desc.add<edm::InputTag>("tracks", edm::InputTag("generalTracks"));
504  desc.add<edm::InputTag>("vertices", edm::InputTag("offlinePrimaryVertices"));
505  desc.add<std::vector<double>>("pTThresholds", {30., 10.});
506  desc.add<double>("maxSVdist", 50.);
507 
508  {
509  edm::ParameterSetDescription psDiMuMass;
510  psDiMuMass.add<std::string>("name", "DiMuMass");
511  psDiMuMass.add<std::string>("title", "M(#mu#mu)");
512  psDiMuMass.add<std::string>("yUnits", "[GeV]");
513  psDiMuMass.add<int>("NxBins", 24);
514  psDiMuMass.add<int>("NyBins", 50);
515  psDiMuMass.add<double>("ymin", 70.);
516  psDiMuMass.add<double>("ymax", 120.);
517  desc.add<edm::ParameterSetDescription>("DiMuMassConfig", psDiMuMass);
518  }
519  {
521  psCosPhi.add<std::string>("name", "CosPhi");
522  psCosPhi.add<std::string>("title", "cos(#phi_{xy})");
523  psCosPhi.add<std::string>("yUnits", "");
524  psCosPhi.add<int>("NxBins", 50);
525  psCosPhi.add<int>("NyBins", 50);
526  psCosPhi.add<double>("ymin", -1.);
527  psCosPhi.add<double>("ymax", 1.);
528  desc.add<edm::ParameterSetDescription>("CosPhiConfig", psCosPhi);
529  }
530  {
531  edm::ParameterSetDescription psCosPhi3D;
532  psCosPhi3D.add<std::string>("name", "CosPhi3D");
533  psCosPhi3D.add<std::string>("title", "cos(#phi_{3D})");
534  psCosPhi3D.add<std::string>("yUnits", "");
535  psCosPhi3D.add<int>("NxBins", 50);
536  psCosPhi3D.add<int>("NyBins", 50);
537  psCosPhi3D.add<double>("ymin", -1.);
538  psCosPhi3D.add<double>("ymax", 1.);
539  desc.add<edm::ParameterSetDescription>("CosPhi3DConfig", psCosPhi3D);
540  }
541  {
543  psVtxProb.add<std::string>("name", "VtxProb");
544  psVtxProb.add<std::string>("title", "Prob(#chi^{2}_{SV})");
545  psVtxProb.add<std::string>("yUnits", "");
546  psVtxProb.add<int>("NxBins", 50);
547  psVtxProb.add<int>("NyBins", 50);
548  psVtxProb.add<double>("ymin", 0);
549  psVtxProb.add<double>("ymax", 1.);
550  desc.add<edm::ParameterSetDescription>("VtxProbConfig", psVtxProb);
551  }
552  {
554  psVtxDist.add<std::string>("name", "VtxDist");
555  psVtxDist.add<std::string>("title", "d_{xy}(PV,SV)");
556  psVtxDist.add<std::string>("yUnits", "[#mum]");
557  psVtxDist.add<int>("NxBins", 50);
558  psVtxDist.add<int>("NyBins", 100);
559  psVtxDist.add<double>("ymin", 0);
560  psVtxDist.add<double>("ymax", 300.);
561  desc.add<edm::ParameterSetDescription>("VtxDistConfig", psVtxDist);
562  }
563  {
564  edm::ParameterSetDescription psVtxDist3D;
565  psVtxDist3D.add<std::string>("name", "VtxDist3D");
566  psVtxDist3D.add<std::string>("title", "d_{3D}(PV,SV)");
567  psVtxDist3D.add<std::string>("yUnits", "[#mum]");
568  psVtxDist3D.add<int>("NxBins", 50);
569  psVtxDist3D.add<int>("NyBins", 250);
570  psVtxDist3D.add<double>("ymin", 0);
571  psVtxDist3D.add<double>("ymax", 500.);
572  desc.add<edm::ParameterSetDescription>("VtxDist3DConfig", psVtxDist3D);
573  }
574  {
575  edm::ParameterSetDescription psVtxDistSig;
576  psVtxDistSig.add<std::string>("name", "VtxDistSig");
577  psVtxDistSig.add<std::string>("title", "d_{xy}(PV,SV)/#sigma_{dxy}(PV,SV)");
578  psVtxDistSig.add<std::string>("yUnits", "");
579  psVtxDistSig.add<int>("NxBins", 50);
580  psVtxDistSig.add<int>("NyBins", 100);
581  psVtxDistSig.add<double>("ymin", 0);
582  psVtxDistSig.add<double>("ymax", 5.);
583  desc.add<edm::ParameterSetDescription>("VtxDistSigConfig", psVtxDistSig);
584  }
585  {
586  edm::ParameterSetDescription psVtxDist3DSig;
587  psVtxDist3DSig.add<std::string>("name", "VtxDist3DSig");
588  psVtxDist3DSig.add<std::string>("title", "d_{3D}(PV,SV)/#sigma_{d3D}(PV,SV)");
589  psVtxDist3DSig.add<std::string>("yUnits", "");
590  psVtxDist3DSig.add<int>("NxBins", 50);
591  psVtxDist3DSig.add<int>("NyBins", 100);
592  psVtxDist3DSig.add<double>("ymin", 0);
593  psVtxDist3DSig.add<double>("ymax", 5.);
594  desc.add<edm::ParameterSetDescription>("VtxDist3DSigConfig", psVtxDist3DSig);
595  }
596 
597  descriptions.addWithDefaultLabel(desc);
598 }
599 
600 //define this as a plug-in
DiMuonVertexValidation(const edm::ParameterSet &)
static const std::string kSharedResource
Definition: TFileService.h:76
~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_
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
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
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)
bool getData(T &iHolder) const
Definition: EventSetup.h:122
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
edm::ParameterSet CosPhiConfiguration_
DiLeptonHelp::PlotsVsKinematics CosPhi3DPlots
fixed size matrix
HLT enums.
edm::ParameterSet CosPhi3DConfiguration_
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)
DiLeptonHelp::PlotsVsKinematics ZMassPlots
edm::EDGetTokenT< reco::TrackCollection > alcaRecoToken_
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:72
#define LogDebug(id)