CMS 3D CMS Logo

SplitVertexResolution.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Alignment/OfflineValidation
4 // Class: SplitVertexResolution
5 //
9 //
10 // Original Author: Marco Musich
11 // Created: Mon, 13 Jun 2016 15:07:11 GMT
12 //
13 //
14 
15 // system include files
16 #include <memory>
17 #include <algorithm> // std::sort
18 #include <vector> // std::vector
19 #include <chrono>
20 #include <iostream>
21 #include <random>
22 
23 // ROOT include files
24 #include "TTree.h"
25 #include "TProfile.h"
26 #include "TF1.h"
27 #include "TMath.h"
28 
29 // user include files
32 
40 
42 
47 
51 
57 
60 
61 //
62 // useful code
63 //
64 
65 namespace statmode {
66  using fitParams = std::pair<Measurement1D, Measurement1D>;
67 }
68 
69 //
70 // class declaration
71 //
72 
73 class SplitVertexResolution : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
74 public:
76  ~SplitVertexResolution() override;
77 
78  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
79  static bool mysorter(reco::Track i, reco::Track j) { return (i.pt() > j.pt()); }
80 
81 private:
82  void beginJob() override;
83  void beginRun(edm::Run const& iEvent, edm::EventSetup const&) override;
84  virtual void beginEvent() final;
85  void analyze(const edm::Event&, const edm::EventSetup&) override;
86  void endJob() override;
87  void endRun(edm::Run const&, edm::EventSetup const&) override{};
88 
89  template <std::size_t SIZE>
90  bool checkBinOrdering(std::array<float, SIZE>& bins);
91  std::vector<TH1F*> bookResidualsHistogram(TFileDirectory dir,
92  unsigned int theNOfBins,
93  TString resType,
94  TString varType);
95 
96  void fillTrendPlotByIndex(TH1F* trendPlot, std::vector<TH1F*>& h, PVValHelper::estimator fitPar_);
97  statmode::fitParams fitResiduals(TH1* hist, bool singleTime = false);
99 
100  std::pair<long long, long long> getRunTime(const edm::EventSetup& iSetup) const;
101 
102  // counters
103  int ievt;
104  int itrks;
105 
106  // switch to keep the ntuple
108 
109  // storing integrated luminosity
110  double intLumi_;
111  bool debug_;
112 
115 
118 
119  edm::InputTag triggerResultsTag_ = edm::InputTag("TriggerResults", "", "HLT"); //InputTag tag("TriggerResults");
121 
122  // ES Tokens
125 
126  double minVtxNdf_;
127  double minVtxWgt_;
128 
129  // checks on the run to be processed
130 
132  std::vector<unsigned int> runControlNumbers_;
133 
134  static constexpr double cmToUm = 10000.;
135 
137 
140 
141  std::map<unsigned int, std::pair<long long, long long> > runNumbersTimesLog_;
144 
145  TH1F* h_diffX;
146  TH1F* h_diffY;
147  TH1F* h_diffZ;
148 
152 
153  TH1F* h_errX;
154  TH1F* h_errY;
155  TH1F* h_errZ;
156 
157  TH1F* h_pullX;
158  TH1F* h_pullY;
159  TH1F* h_pullZ;
160 
161  TH1F* h_ntrks;
162  TH1F* h_sumPt;
163  TH1F* h_avgSumPt;
164 
165  TH1F* h_sumPt1;
166  TH1F* h_sumPt2;
167 
168  TH1F* h_wTrks1;
169  TH1F* h_wTrks2;
170 
171  TH1F* h_minWTrks1;
172  TH1F* h_minWTrks2;
173 
176 
177  TH1F* h_runNumber;
178 
180  TH1I* h_nVertices;
183 
184  // trigger results
185 
188 
189  // resolutions
190 
191  std::vector<TH1F*> h_resolX_sumPt_;
192  std::vector<TH1F*> h_resolY_sumPt_;
193  std::vector<TH1F*> h_resolZ_sumPt_;
194 
195  std::vector<TH1F*> h_resolX_Ntracks_;
196  std::vector<TH1F*> h_resolY_Ntracks_;
197  std::vector<TH1F*> h_resolZ_Ntracks_;
198 
199  std::vector<TH1F*> h_resolX_Nvtx_;
200  std::vector<TH1F*> h_resolY_Nvtx_;
201  std::vector<TH1F*> h_resolZ_Nvtx_;
202 
206 
210 
214 
215  // pulls
216  std::vector<TH1F*> h_pullX_sumPt_;
217  std::vector<TH1F*> h_pullY_sumPt_;
218  std::vector<TH1F*> h_pullZ_sumPt_;
219 
220  std::vector<TH1F*> h_pullX_Ntracks_;
221  std::vector<TH1F*> h_pullY_Ntracks_;
222  std::vector<TH1F*> h_pullZ_Ntracks_;
223 
224  std::vector<TH1F*> h_pullX_Nvtx_;
225  std::vector<TH1F*> h_pullY_Nvtx_;
226  std::vector<TH1F*> h_pullZ_Nvtx_;
227 
231 
235 
239 
240  std::mt19937 engine_;
241 
243  TTree* tree_;
244 
245  // ----------member data ---------------------------
246  static const int nPtBins_ = 30;
247  std::array<float, nPtBins_ + 1> mypT_bins_ = PVValHelper::makeLogBins<float, nPtBins_>(1., 1e3);
248 
249  static const int nTrackBins_ = 60;
250  std::array<float, nTrackBins_ + 1> myNTrack_bins_;
251 
252  static const int nVtxBins_ = 40;
253  std::array<float, nVtxBins_ + 1> myNVtx_bins_;
254 
255  std::map<std::string, std::pair<int, int> > triggerMap_;
256 };
257 
259  : storeNtuple_(iConfig.getParameter<bool>("storeNtuple")),
260  intLumi_(iConfig.getUntrackedParameter<double>("intLumi", 0.)),
261  debug_(iConfig.getUntrackedParameter<bool>("Debug", false)),
262  pvsTag_(iConfig.getParameter<edm::InputTag>("vtxCollection")),
263  pvsToken_(consumes<reco::VertexCollection>(pvsTag_)),
264  tracksTag_(iConfig.getParameter<edm::InputTag>("trackCollection")),
265  tracksToken_(consumes<reco::TrackCollection>(tracksTag_)),
266  triggerResultsToken_(consumes<edm::TriggerResults>(triggerResultsTag_)),
267  transientTrackBuilderToken_(
268  esConsumes<TransientTrackBuilder, TransientTrackRecord>(edm::ESInputTag("", "TransientTrackBuilder"))),
269  runInfoToken_(esConsumes<RunInfo, RunInfoRcd, edm::Transition::BeginRun>()),
270  minVtxNdf_(iConfig.getUntrackedParameter<double>("minVertexNdf")),
271  minVtxWgt_(iConfig.getUntrackedParameter<double>("minVertexMeanWeight")),
272  runControl_(iConfig.getUntrackedParameter<bool>("runControl", false)) {
273  usesResource(TFileService::kSharedResource);
274 
275  std::vector<unsigned int> defaultRuns;
276  defaultRuns.push_back(0);
277  runControlNumbers_ = iConfig.getUntrackedParameter<std::vector<unsigned int> >("runControlNumber", defaultRuns);
278 
279  std::vector<float> vect = PVValHelper::generateBins(nTrackBins_ + 1, -0.5, 120.);
280  std::copy(vect.begin(), vect.begin() + nTrackBins_ + 1, myNTrack_bins_.begin());
281 
282  vect.clear();
283  vect = PVValHelper::generateBins(nVtxBins_ + 1, 1., 40.);
284  std::copy(vect.begin(), vect.begin() + nVtxBins_ + 1, myNVtx_bins_.begin());
285 }
286 
288 
289 // ------------ method called for each event ------------
291  using namespace edm;
292 
293  // deterministic seed from the event number
294  // should not bias the result as the event number is already
295  // assigned randomly-enough
296  engine_.seed(iEvent.id().event() + (iEvent.id().luminosityBlock() << 10) + (iEvent.id().run() << 20));
297 
298  // Fill general info
299  h_runNumber->Fill(iEvent.id().run());
300 
301  bool passesRunControl = false;
302 
303  if (runControl_) {
304  for (const auto& runControlNumber : runControlNumbers_) {
305  if (iEvent.eventAuxiliary().run() == runControlNumber) {
306  if (debug_) {
307  edm::LogInfo("SplitVertexResolution")
308  << " run number: " << iEvent.eventAuxiliary().run() << " keeping run:" << runControlNumber;
309  }
310  passesRunControl = true;
311  break;
312  }
313  }
314  if (!passesRunControl)
315  return;
316  }
317 
318  ievt++;
321 
322  const edm::TriggerNames& triggerNames_ = iEvent.triggerNames(*hltresults);
323  int ntrigs = hltresults->size();
324  //const std::vector<std::string>& triggernames = triggerNames_.triggerNames();
325 
326  beginEvent();
327 
328  // Fill general info
329  event_.runNumber = iEvent.id().run();
330  event_.luminosityBlockNumber = iEvent.id().luminosityBlock();
331  event_.eventNumber = iEvent.id().event();
332 
334 
336  iEvent.getByToken(pvsToken_, vertices);
337  const reco::VertexCollection pvtx = *(vertices.product());
338 
339  event_.nVtx = pvtx.size();
340  int nOfflineVtx = pvtx.size();
341  h_nOfflineVertices->Fill(nOfflineVtx);
342 
344  iEvent.getByToken(tracksToken_, tracks);
345  itrks += tracks.product()->size();
346 
347  for (int itrig = 0; itrig != ntrigs; ++itrig) {
348  const std::string& trigName = triggerNames_.triggerName(itrig);
349  bool accept = hltresults->accept(itrig);
350  if (accept == 1) {
351  triggerMap_[trigName].first += 1;
352  triggerMap_[trigName].second += tracks.product()->size();
353  // triggerInfo.push_back(pair <string, int> (trigName, accept));
354  }
355  }
356 
357  int counter = 0;
358  int noFakecounter = 0;
359  int goodcounter = 0;
360 
361  for (auto pvIt = pvtx.cbegin(); pvIt != pvtx.cend(); ++pvIt) {
362  reco::Vertex iPV = *pvIt;
363  counter++;
364  if (iPV.isFake())
365  continue;
366  noFakecounter++;
367 
368  // vertex selection as in bs code
369  if (iPV.ndof() < minVtxNdf_ || (iPV.ndof() + 3.) / iPV.tracksSize() < 2 * minVtxWgt_)
370  continue;
371 
372  goodcounter++;
374  reco::TrackCollection groupOne, groupTwo;
375  for (auto trki = iPV.tracks_begin(); trki != iPV.tracks_end(); ++trki) {
376  if (trki->isNonnull()) {
377  reco::TrackRef trk_now(tracks, (*trki).key());
378  allTracks.push_back(*trk_now);
379  }
380  }
381 
382  if (goodcounter > 1)
383  continue;
384 
385  // order with decreasing pt
386  std::sort(allTracks.begin(), allTracks.end(), mysorter);
387 
388  int ntrks = allTracks.size();
389  h_ntrks->Fill(ntrks);
390 
391  // discard lowest pt track
392  uint even_ntrks;
393  ntrks % 2 == 0 ? even_ntrks = ntrks : even_ntrks = ntrks - 1;
394 
395  // split into two sets equally populated
396  for (uint tracksIt = 0; tracksIt < even_ntrks; tracksIt = tracksIt + 2) {
397  reco::Track firstTrk = allTracks.at(tracksIt);
398  reco::Track secondTrk = allTracks.at(tracksIt + 1);
399  auto dis = std::uniform_int_distribution<>(0, 1); // [0, 1]
400 
401  if (dis(engine_) > 0.5) {
402  groupOne.push_back(firstTrk);
403  groupTwo.push_back(secondTrk);
404  } else {
405  groupOne.push_back(secondTrk);
406  groupTwo.push_back(firstTrk);
407  }
408  }
409 
410  if (!(groupOne.size() >= 2 && groupTwo.size() >= 2))
411  continue;
412 
413  h_OrigVertexErrX->Fill(iPV.xError() * cmToUm);
414  h_OrigVertexErrY->Fill(iPV.yError() * cmToUm);
415  h_OrigVertexErrZ->Fill(iPV.zError() * cmToUm);
416 
417  float sumPt = 0, sumPt1 = 0, sumPt2 = 0, avgSumPt = 0;
418 
419  // refit the two sets of tracks
420  std::vector<reco::TransientTrack> groupOne_ttks;
421  groupOne_ttks.clear();
422  for (auto itrk = groupOne.cbegin(); itrk != groupOne.cend(); itrk++) {
423  reco::TransientTrack tmpTransientTrack = theB.build(*itrk);
424  groupOne_ttks.push_back(tmpTransientTrack);
425  sumPt1 += itrk->pt();
426  sumPt += itrk->pt();
427  }
428 
429  AdaptiveVertexFitter pvFitter;
430  TransientVertex pvOne = pvFitter.vertex(groupOne_ttks);
431  if (!pvOne.isValid())
432  continue;
433 
434  reco::Vertex onePV = pvOne;
435 
436  std::vector<reco::TransientTrack> groupTwo_ttks;
437  groupTwo_ttks.clear();
438  for (auto itrk = groupTwo.cbegin(); itrk != groupTwo.cend(); itrk++) {
439  reco::TransientTrack tmpTransientTrack = theB.build(*itrk);
440  groupTwo_ttks.push_back(tmpTransientTrack);
441  sumPt2 += itrk->pt();
442  sumPt += itrk->pt();
443  }
444 
445  // average sumPt
446  avgSumPt = (sumPt1 + sumPt2) / 2.;
447  h_avgSumPt->Fill(avgSumPt);
448 
449  TransientVertex pvTwo = pvFitter.vertex(groupTwo_ttks);
450  if (!pvTwo.isValid())
451  continue;
452 
453  reco::Vertex twoPV = pvTwo;
454 
455  float theminW1 = 1.;
456  float theminW2 = 1.;
457  for (auto otrk = pvOne.originalTracks().cbegin(); otrk != pvOne.originalTracks().cend(); ++otrk) {
458  h_wTrks1->Fill(pvOne.trackWeight(*otrk));
459  if (pvOne.trackWeight(*otrk) < theminW1) {
460  theminW1 = pvOne.trackWeight(*otrk);
461  }
462  }
463  for (auto otrk = pvTwo.originalTracks().cbegin(); otrk != pvTwo.originalTracks().end(); ++otrk) {
464  h_wTrks2->Fill(pvTwo.trackWeight(*otrk));
465  if (pvTwo.trackWeight(*otrk) < theminW2) {
466  theminW2 = pvTwo.trackWeight(*otrk);
467  }
468  }
469 
470  h_sumPt->Fill(sumPt);
471 
472  int half_trks = twoPV.nTracks();
473 
474  const double invSqrt2 = 1. / std::sqrt(2.);
475 
476  double deltaX = (twoPV.x() - onePV.x());
477  double deltaY = (twoPV.y() - onePV.y());
478  double deltaZ = (twoPV.z() - onePV.z());
479 
480  double resX = deltaX * invSqrt2;
481  double resY = deltaY * invSqrt2;
482  double resZ = deltaZ * invSqrt2;
483 
484  h_diffX->Fill(resX * cmToUm);
485  h_diffY->Fill(resY * cmToUm);
486  h_diffZ->Fill(resZ * cmToUm);
487 
488  double errX = sqrt(pow(twoPV.xError(), 2) + pow(onePV.xError(), 2));
489  double errY = sqrt(pow(twoPV.yError(), 2) + pow(onePV.yError(), 2));
490  double errZ = sqrt(pow(twoPV.zError(), 2) + pow(onePV.zError(), 2));
491 
492  h_errX->Fill(errX * cmToUm);
493  h_errY->Fill(errY * cmToUm);
494  h_errZ->Fill(errZ * cmToUm);
495 
496  h_pullX->Fill(deltaX / errX);
497  h_pullY->Fill(deltaY / errY);
498  h_pullZ->Fill(deltaZ / errZ);
499 
500  // filling the pT-binned distributions
501 
502  for (int ipTBin = 0; ipTBin < nPtBins_; ipTBin++) {
503  float pTF = mypT_bins_[ipTBin];
504  float pTL = mypT_bins_[ipTBin + 1];
505 
506  if (avgSumPt >= pTF && avgSumPt < pTL) {
507  PVValHelper::fillByIndex(h_resolX_sumPt_, ipTBin, resX * cmToUm, "1");
508  PVValHelper::fillByIndex(h_resolY_sumPt_, ipTBin, resY * cmToUm, "2");
509  PVValHelper::fillByIndex(h_resolZ_sumPt_, ipTBin, resZ * cmToUm, "3");
510 
511  PVValHelper::fillByIndex(h_pullX_sumPt_, ipTBin, deltaX / errX, "4");
512  PVValHelper::fillByIndex(h_pullY_sumPt_, ipTBin, deltaY / errY, "5");
513  PVValHelper::fillByIndex(h_pullZ_sumPt_, ipTBin, deltaZ / errZ, "6");
514  }
515  }
516 
517  // filling the track multeplicity binned distributions
518 
519  for (int inTrackBin = 0; inTrackBin < nTrackBins_; inTrackBin++) {
520  float nTrackF = myNTrack_bins_[inTrackBin];
521  float nTrackL = myNTrack_bins_[inTrackBin + 1];
522 
523  if (ntrks >= nTrackF && ntrks < nTrackL) {
524  PVValHelper::fillByIndex(h_resolX_Ntracks_, inTrackBin, resX * cmToUm, "7");
525  PVValHelper::fillByIndex(h_resolY_Ntracks_, inTrackBin, resY * cmToUm, "8");
526  PVValHelper::fillByIndex(h_resolZ_Ntracks_, inTrackBin, resZ * cmToUm, "9");
527 
528  PVValHelper::fillByIndex(h_pullX_Ntracks_, inTrackBin, deltaX / errX, "10");
529  PVValHelper::fillByIndex(h_pullY_Ntracks_, inTrackBin, deltaY / errY, "11");
530  PVValHelper::fillByIndex(h_pullZ_Ntracks_, inTrackBin, deltaZ / errZ, "12");
531  }
532  }
533 
534  // filling the vertex multeplicity binned distributions
535 
536  for (int inVtxBin = 0; inVtxBin < nVtxBins_; inVtxBin++) {
537  /*
538  float nVtxF = myNVtx_bins_[inVtxBin];
539  float nVtxL = myNVtx_bins_[inVtxBin+1];
540  if(nOfflineVtx >= nVtxF && nOfflineVtx < nVtxL){
541  */
542 
543  if (nOfflineVtx == inVtxBin) {
544  PVValHelper::fillByIndex(h_resolX_Nvtx_, inVtxBin, deltaX * cmToUm, "7");
545  PVValHelper::fillByIndex(h_resolY_Nvtx_, inVtxBin, deltaY * cmToUm, "8");
546  PVValHelper::fillByIndex(h_resolZ_Nvtx_, inVtxBin, deltaZ * cmToUm, "9");
547 
548  PVValHelper::fillByIndex(h_pullX_Nvtx_, inVtxBin, deltaX / errX, "10");
549  PVValHelper::fillByIndex(h_pullY_Nvtx_, inVtxBin, deltaY / errY, "11");
550  PVValHelper::fillByIndex(h_pullZ_Nvtx_, inVtxBin, deltaZ / errZ, "12");
551  }
552  }
553 
554  h_sumPt1->Fill(sumPt1);
555  h_sumPt2->Fill(sumPt2);
556 
557  h_minWTrks1->Fill(theminW1);
558  h_minWTrks2->Fill(theminW2);
559 
560  h_PVCL_subVtx1->Fill(TMath::Prob(pvOne.totalChiSquared(), (int)(pvOne.degreesOfFreedom())));
561  h_PVCL_subVtx2->Fill(TMath::Prob(pvTwo.totalChiSquared(), (int)(pvTwo.degreesOfFreedom())));
562 
563  // fill ntuples
564  pvCand thePV;
565  thePV.ipos = counter;
566  thePV.nTrks = ntrks;
567 
568  thePV.x_origVtx = iPV.x();
569  thePV.y_origVtx = iPV.y();
570  thePV.z_origVtx = iPV.z();
571 
572  thePV.xErr_origVtx = iPV.xError();
573  thePV.yErr_origVtx = iPV.yError();
574  thePV.zErr_origVtx = iPV.zError();
575 
576  thePV.n_subVtx1 = half_trks;
577  thePV.x_subVtx1 = onePV.x();
578  thePV.y_subVtx1 = onePV.y();
579  thePV.z_subVtx1 = onePV.z();
580 
581  thePV.xErr_subVtx1 = onePV.xError();
582  thePV.yErr_subVtx1 = onePV.yError();
583  thePV.zErr_subVtx1 = onePV.zError();
584  thePV.sumPt_subVtx1 = sumPt1;
585 
586  thePV.n_subVtx2 = half_trks;
587  thePV.x_subVtx2 = twoPV.x();
588  thePV.y_subVtx2 = twoPV.y();
589  thePV.z_subVtx2 = twoPV.z();
590 
591  thePV.xErr_subVtx2 = twoPV.xError();
592  thePV.yErr_subVtx2 = twoPV.yError();
593  thePV.zErr_subVtx2 = twoPV.zError();
594  thePV.sumPt_subVtx2 = sumPt2;
595 
596  thePV.CL_subVtx1 = TMath::Prob(pvOne.totalChiSquared(), (int)(pvOne.degreesOfFreedom()));
597  thePV.CL_subVtx2 = TMath::Prob(pvTwo.totalChiSquared(), (int)(pvTwo.degreesOfFreedom()));
598 
599  thePV.minW_subVtx1 = theminW1;
600  thePV.minW_subVtx2 = theminW2;
601 
602  event_.pvs.push_back(thePV);
603 
604  } // loop on the vertices
605 
606  // fill the histogram of vertices per event
607  h_nVertices->Fill(counter);
608  h_nNonFakeVertices->Fill(noFakecounter);
609  h_nFinalVertices->Fill(goodcounter);
610 
611  if (storeNtuple_) {
612  tree_->Fill();
613  }
614 }
615 
617  event_.pvs.clear();
618  event_.nVtx = -1;
619 }
620 
622  unsigned int RunNumber_ = run.run();
623 
624  if (!runNumbersTimesLog_.count(RunNumber_)) {
625  auto times = getRunTime(iSetup);
626 
627  if (debug_) {
628  const time_t start_time = times.first / 1000000;
629  edm::LogInfo("SplitVertexResolution")
630  << RunNumber_ << " has start time: " << times.first << " - " << times.second << std::endl;
631  edm::LogInfo("SplitVertexResolution")
632  << "human readable time: " << std::asctime(std::gmtime(&start_time)) << std::endl;
633  }
634  runNumbersTimesLog_[RunNumber_] = times;
635  }
636 }
637 
638 // ------------ method called once each job just before starting event loop ------------
640  ievt = 0;
641  itrks = 0;
642 
643  // luminosity histo
644  TFileDirectory EventFeatures = outfile_->mkdir("EventFeatures");
646  EventFeatures.make<TH1F>("h_lumiFromConfig", "luminosity from config;;luminosity of present run", 1, -0.5, 0.5);
647  h_lumiFromConfig->SetBinContent(1, intLumi_);
648 
649  h_runFromConfig = EventFeatures.make<TH1I>("h_runFromConfig",
650  "run number from config;;run number (from configuration)",
651  runControlNumbers_.size(),
652  0.,
653  runControlNumbers_.size());
654  for (const auto r : runControlNumbers_) {
655  h_runFromConfig->SetBinContent(r + 1, r);
656  }
657 
658  // resolutions
659 
661  edm::LogError("SplitVertexResolution") << " Warning - the vector of pT bins is not ordered " << std::endl;
662  }
663 
665  edm::LogError("SplitVertexResolution") << " Warning -the vector of n. tracks bins is not ordered " << std::endl;
666  }
667 
669  edm::LogError("SplitVertexResolution") << " Warning -the vector of n. vertices bins is not ordered " << std::endl;
670  }
671 
672  TFileDirectory xResolSumPt = outfile_->mkdir("xResolSumPt");
673  h_resolX_sumPt_ = bookResidualsHistogram(xResolSumPt, nPtBins_, "resolX", "sumPt");
674 
675  TFileDirectory yResolSumPt = outfile_->mkdir("yResolSumPt");
676  h_resolY_sumPt_ = bookResidualsHistogram(yResolSumPt, nPtBins_, "resolY", "sumPt");
677 
678  TFileDirectory zResolSumPt = outfile_->mkdir("zResolSumPt");
679  h_resolZ_sumPt_ = bookResidualsHistogram(zResolSumPt, nPtBins_, "resolZ", "sumPt");
680 
681  TFileDirectory xResolNtracks_ = outfile_->mkdir("xResolNtracks");
682  h_resolX_Ntracks_ = bookResidualsHistogram(xResolNtracks_, nTrackBins_, "resolX", "Ntracks");
683 
684  TFileDirectory yResolNtracks_ = outfile_->mkdir("yResolNtracks");
685  h_resolY_Ntracks_ = bookResidualsHistogram(yResolNtracks_, nTrackBins_, "resolY", "Ntracks");
686 
687  TFileDirectory zResolNtracks_ = outfile_->mkdir("zResolNtracks");
688  h_resolZ_Ntracks_ = bookResidualsHistogram(zResolNtracks_, nTrackBins_, "resolZ", "Ntracks");
689 
690  TFileDirectory xResolNvtx_ = outfile_->mkdir("xResolNvtx");
691  h_resolX_Nvtx_ = bookResidualsHistogram(xResolNvtx_, nVtxBins_, "resolX", "Nvtx");
692 
693  TFileDirectory yResolNvtx_ = outfile_->mkdir("yResolNvtx");
694  h_resolY_Nvtx_ = bookResidualsHistogram(yResolNvtx_, nVtxBins_, "resolY", "Nvtx");
695 
696  TFileDirectory zResolNvtx_ = outfile_->mkdir("zResolNvtx");
697  h_resolZ_Nvtx_ = bookResidualsHistogram(zResolNvtx_, nVtxBins_, "resolZ", "Nvtx");
698 
699  // pulls
700 
701  TFileDirectory xPullSumPt = outfile_->mkdir("xPullSumPt");
702  h_pullX_sumPt_ = bookResidualsHistogram(xPullSumPt, nPtBins_, "pullX", "sumPt");
703 
704  TFileDirectory yPullSumPt = outfile_->mkdir("yPullSumPt");
705  h_pullY_sumPt_ = bookResidualsHistogram(yPullSumPt, nPtBins_, "pullY", "sumPt");
706 
707  TFileDirectory zPullSumPt = outfile_->mkdir("zPullSumPt");
708  h_pullZ_sumPt_ = bookResidualsHistogram(zPullSumPt, nPtBins_, "pullZ", "sumPt");
709 
710  TFileDirectory xPullNtracks_ = outfile_->mkdir("xPullNtracks");
711  h_pullX_Ntracks_ = bookResidualsHistogram(xPullNtracks_, nTrackBins_, "pullX", "Ntracks");
712 
713  TFileDirectory yPullNtracks_ = outfile_->mkdir("yPullNtracks");
714  h_pullY_Ntracks_ = bookResidualsHistogram(yPullNtracks_, nTrackBins_, "pullY", "Ntracks");
715 
716  TFileDirectory zPullNtracks_ = outfile_->mkdir("zPullNtracks");
717  h_pullZ_Ntracks_ = bookResidualsHistogram(zPullNtracks_, nTrackBins_, "pullZ", "Ntracks");
718 
719  TFileDirectory xPullNvtx_ = outfile_->mkdir("xPullNvtx");
720  h_pullX_Nvtx_ = bookResidualsHistogram(xPullNvtx_, nVtxBins_, "pullX", "Nvtx");
721 
722  TFileDirectory yPullNvtx_ = outfile_->mkdir("yPullNvtx");
723  h_pullY_Nvtx_ = bookResidualsHistogram(yPullNvtx_, nVtxBins_, "pullY", "Nvtx");
724 
725  TFileDirectory zPullNvtx_ = outfile_->mkdir("zPullNvtx");
726  h_pullZ_Nvtx_ = bookResidualsHistogram(zPullNvtx_, nVtxBins_, "pullZ", "Nvtx");
727 
728  // control plots
729  h_runNumber = outfile_->make<TH1F>("h_runNumber", "run number;run number;n_{events}", 100000, 250000., 350000.);
730 
731  h_nOfflineVertices = outfile_->make<TH1I>("h_nOfflineVertices", "n. of vertices;n. vertices; events", 100, 0, 100);
732  h_nVertices = outfile_->make<TH1I>("h_nVertices", "n. of vertices;n. vertices; events", 100, 0, 100);
734  outfile_->make<TH1I>("h_nRealVertices", "n. of non-fake vertices;n. vertices; events", 100, 0, 100);
736  outfile_->make<TH1I>("h_nSelectedVertices", "n. of selected vertices vertices;n. vertices; events", 100, 0, 100);
737 
738  h_diffX = outfile_->make<TH1F>(
739  "h_diffX", "x-coordinate vertex resolution;vertex resolution (x) [#mum];vertices", 100, -300, 300.);
740  h_diffY = outfile_->make<TH1F>(
741  "h_diffY", "y-coordinate vertex resolution;vertex resolution (y) [#mum];vertices", 100, -300, 300.);
742  h_diffZ = outfile_->make<TH1F>(
743  "h_diffZ", "z-coordinate vertex resolution;vertex resolution (z) [#mum];vertices", 100, -500, 500.);
744 
745  h_OrigVertexErrX = outfile_->make<TH1F>(
746  "h_OrigVertexErrX", "x-coordinate vertex error;vertex error (x) [#mum];vertices", 300, 0., 300.);
747  h_OrigVertexErrY = outfile_->make<TH1F>(
748  "h_OrigVertexErrY", "y-coordinate vertex error;vertex error (y) [#mum];vertices", 300, 0., 300.);
749  h_OrigVertexErrZ = outfile_->make<TH1F>(
750  "h_OrigVertexErrZ", "z-coordinate vertex error;vertex error (z) [#mum];vertices", 500, 0., 500.);
751 
752  h_errX = outfile_->make<TH1F>(
753  "h_errX", "x-coordinate vertex resolution error;vertex resoltuion error (x) [#mum];vertices", 300, 0., 300.);
754  h_errY = outfile_->make<TH1F>(
755  "h_errY", "y-coordinate vertex resolution error;vertex resolution error (y) [#mum];vertices", 300, 0., 300.);
756  h_errZ = outfile_->make<TH1F>(
757  "h_errZ", "z-coordinate vertex resolution error;vertex resolution error (z) [#mum];vertices", 500, 0., 500.);
758 
759  h_pullX = outfile_->make<TH1F>("h_pullX", "x-coordinate vertex pull;vertex pull (x);vertices", 500, -10, 10.);
760  h_pullY = outfile_->make<TH1F>("h_pullY", "y-coordinate vertex pull;vertex pull (y);vertices", 500, -10, 10.);
761  h_pullZ = outfile_->make<TH1F>("h_pullZ", "z-coordinate vertex pull;vertex pull (z);vertices", 500, -10, 10.);
762 
763  h_ntrks = outfile_->make<TH1F>("h_ntrks",
764  "number of tracks in vertex;vertex multeplicity;vertices",
765  myNTrack_bins_.size() - 1,
766  myNTrack_bins_.data());
767 
768  h_sumPt = outfile_->make<TH1F>(
769  "h_sumPt", "#Sigma p_{T};#sum p_{T} [GeV];vertices", mypT_bins_.size() - 1, mypT_bins_.data());
770 
771  h_avgSumPt = outfile_->make<TH1F>(
772  "h_avgSumPt", "#LT #Sigma p_{T} #GT;#LT #sum p_{T} #GT [GeV];vertices", mypT_bins_.size() - 1, mypT_bins_.data());
773 
774  h_sumPt1 = outfile_->make<TH1F>("h_sumPt1",
775  "#Sigma p_{T} sub-vertex 1;#sum p_{T} sub-vertex 1 [GeV];subvertices",
776  mypT_bins_.size() - 1,
777  mypT_bins_.data());
778  h_sumPt2 = outfile_->make<TH1F>("h_sumPt2",
779  "#Sigma p_{T} sub-vertex 2;#sum p_{T} sub-vertex 2 [GeV];subvertices",
780  mypT_bins_.size() - 1,
781  mypT_bins_.data());
782 
783  h_wTrks1 = outfile_->make<TH1F>("h_wTrks1", "weight of track for sub-vertex 1;track weight;subvertices", 500, 0., 1.);
784  h_wTrks2 = outfile_->make<TH1F>("h_wTrks2", "weithg of track for sub-vertex 2;track weight;subvertices", 500, 0., 1.);
785 
786  h_minWTrks1 = outfile_->make<TH1F>(
787  "h_minWTrks1", "minimum weight of track for sub-vertex 1;minimum track weight;subvertices", 500, 0., 1.);
788  h_minWTrks2 = outfile_->make<TH1F>(
789  "h_minWTrks2", "minimum weithg of track for sub-vertex 2;minimum track weight;subvertices", 500, 0., 1.);
790 
792  outfile_->make<TH1F>("h_PVCL_subVtx1",
793  "#chi^{2} probability for sub-vertex 1;Prob(#chi^{2},ndof) for sub-vertex 1;subvertices",
794  100,
795  0.,
796  1);
798  outfile_->make<TH1F>("h_PVCL_subVtx2",
799  "#chi^{2} probability for sub-vertex 2;Prob(#chi^{2},ndof) for sub-vertex 2;subvertices",
800  100,
801  0.,
802  1);
803 
804  // resolutions
805 
806  p_resolX_vsSumPt = outfile_->make<TH1F>("p_resolX_vsSumPt",
807  "x-resolution vs #Sigma p_{T};#sum p_{T}; x vertex resolution [#mum]",
808  mypT_bins_.size() - 1,
809  mypT_bins_.data());
810  p_resolY_vsSumPt = outfile_->make<TH1F>("p_resolY_vsSumPt",
811  "y-resolution vs #Sigma p_{T};#sum p_{T}; y vertex resolution [#mum]",
812  mypT_bins_.size() - 1,
813  mypT_bins_.data());
814  p_resolZ_vsSumPt = outfile_->make<TH1F>("p_resolZ_vsSumPt",
815  "z-resolution vs #Sigma p_{T};#sum p_{T}; z vertex resolution [#mum]",
816  mypT_bins_.size() - 1,
817  mypT_bins_.data());
818 
819  p_resolX_vsNtracks = outfile_->make<TH1F>("p_resolX_vsNtracks",
820  "x-resolution vs n_{tracks};n_{tracks}; x vertex resolution [#mum]",
821  myNTrack_bins_.size() - 1,
822  myNTrack_bins_.data());
823  p_resolY_vsNtracks = outfile_->make<TH1F>("p_resolY_vsNtracks",
824  "y-resolution vs n_{tracks};n_{tracks}; y vertex resolution [#mum]",
825  myNTrack_bins_.size() - 1,
826  myNTrack_bins_.data());
827  p_resolZ_vsNtracks = outfile_->make<TH1F>("p_resolZ_vsNtracks",
828  "z-resolution vs n_{tracks};n_{tracks}; z vertex resolution [#mum]",
829  myNTrack_bins_.size() - 1,
830  myNTrack_bins_.data());
831 
832  p_resolX_vsNvtx = outfile_->make<TH1F>("p_resolX_vsNvtx",
833  "x-resolution vs n_{vertices};n_{vertices}; x vertex resolution [#mum]",
834  myNVtx_bins_.size() - 1,
835  myNVtx_bins_.data());
836  p_resolY_vsNvtx = outfile_->make<TH1F>("p_resolY_vsNvtx",
837  "y-resolution vs n_{vertices};n_{vertices}; y vertex resolution [#mum]",
838  myNVtx_bins_.size() - 1,
839  myNVtx_bins_.data());
840  p_resolZ_vsNvtx = outfile_->make<TH1F>("p_resolZ_vsNvtx",
841  "z-resolution vs n_{vertices};n_{vertices}; z vertex resolution [#mum]",
842  myNVtx_bins_.size() - 1,
843  myNVtx_bins_.data());
844 
845  // pulls
846 
847  p_pullX_vsSumPt = outfile_->make<TH1F>(
848  "p_pullX_vsSumPt", "x-pull vs #Sigma p_{T};#sum p_{T}; x vertex pull", mypT_bins_.size() - 1, mypT_bins_.data());
849  p_pullY_vsSumPt = outfile_->make<TH1F>(
850  "p_pullY_vsSumPt", "y-pull vs #Sigma p_{T};#sum p_{T}; y vertex pull", mypT_bins_.size() - 1, mypT_bins_.data());
851  p_pullZ_vsSumPt = outfile_->make<TH1F>(
852  "p_pullZ_vsSumPt", "z-pull vs #Sigma p_{T};#sum p_{T}; z vertex pull", mypT_bins_.size() - 1, mypT_bins_.data());
853 
854  p_pullX_vsNtracks = outfile_->make<TH1F>("p_pullX_vsNtracks",
855  "x-pull vs n_{tracks};n_{tracks}; x vertex pull",
856  myNTrack_bins_.size() - 1,
857  myNTrack_bins_.data());
858  p_pullY_vsNtracks = outfile_->make<TH1F>("p_pullY_vsNtracks",
859  "y-pull vs n_{tracks};n_{tracks}; y vertex pull",
860  myNTrack_bins_.size() - 1,
861  myNTrack_bins_.data());
862  p_pullZ_vsNtracks = outfile_->make<TH1F>("p_pullZ_vsNtracks",
863  "z-pull vs n_{tracks};n_{tracks}; z vertex pull",
864  myNTrack_bins_.size() - 1,
865  myNTrack_bins_.data());
866 
867  p_pullX_vsNvtx = outfile_->make<TH1F>("p_pullX_vsNvtx",
868  "x-pull vs n_{vertices};n_{vertices}; x vertex pull",
869  myNVtx_bins_.size() - 1,
870  myNVtx_bins_.data());
871  p_pullY_vsNvtx = outfile_->make<TH1F>("p_pullY_vsNvtx",
872  "y-pull vs n_{vertices};n_{vertices}; y vertex pull",
873  myNVtx_bins_.size() - 1,
874  myNVtx_bins_.data());
875  p_pullZ_vsNvtx = outfile_->make<TH1F>("p_pullZ_vsNvtx",
876  "z-pull vs n_{vertices};n_{vertices}; z vertex pull",
877  myNVtx_bins_.size() - 1,
878  myNVtx_bins_.data());
879 
880  tree_ = outfile_->make<TTree>("pvTree", "pvTree");
881  tree_->Branch("event", &event_, 64000, 2);
882 }
883 
884 //*************************************************************
885 // Generic booker function
886 //*************************************************************
888  unsigned int theNOfBins,
889  TString resType,
890  TString varType) {
891  TH1F::SetDefaultSumw2(kTRUE);
892 
893  double up = 500.;
894  double down = -500.;
895 
896  if (resType.Contains("pull")) {
897  up *= 0.01;
898  down *= 0.01;
899  }
900 
901  std::vector<TH1F*> h;
902  h.reserve(theNOfBins);
903 
904  const char* auxResType = (resType.ReplaceAll("_", "")).Data();
905 
906  for (unsigned int i = 0; i < theNOfBins; i++) {
907  TH1F* htemp = dir.make<TH1F>(Form("histo_%s_%s_plot%i", resType.Data(), varType.Data(), i),
908  Form("%s vs %s - bin %i;%s;vertices", auxResType, varType.Data(), i, auxResType),
909  250,
910  down,
911  up);
912  h.push_back(htemp);
913  }
914 
915  return h;
916 }
917 
918 // ------------ method called once each job just after ending the event loop ------------
920  edm::LogVerbatim("SplitVertexResolution") << "*******************************" << std::endl;
921  edm::LogVerbatim("SplitVertexResolution") << "Events run in total: " << ievt << std::endl;
922  edm::LogVerbatim("SplitVertexResolution") << "n. tracks: " << itrks << std::endl;
923  edm::LogVerbatim("SplitVertexResolution") << "*******************************" << std::endl;
924 
925  int nFiringTriggers = triggerMap_.size();
926  edm::LogVerbatim("SplitVertexResolution") << "firing triggers: " << nFiringTriggers << std::endl;
927  edm::LogVerbatim("SplitVertexResolution") << "*******************************" << std::endl;
928 
929  tksByTrigger_ = outfile_->make<TH1D>(
930  "tksByTrigger", "tracks by HLT path;;% of # traks", nFiringTriggers, -0.5, nFiringTriggers - 0.5);
931  evtsByTrigger_ = outfile_->make<TH1D>(
932  "evtsByTrigger", "events by HLT path;;% of # events", nFiringTriggers, -0.5, nFiringTriggers - 0.5);
933 
934  int i = 0;
935  for (std::map<std::string, std::pair<int, int> >::iterator it = triggerMap_.begin(); it != triggerMap_.end(); ++it) {
936  i++;
937 
938  double trkpercent = ((it->second).second) * 100. / double(itrks);
939  double evtpercent = ((it->second).first) * 100. / double(ievt);
940 
941  edm::LogVerbatim("SplitVertexResolution")
942  << "HLT path: " << std::setw(60) << std::left << it->first << " | events firing: " << std::right << std::setw(8)
943  << (it->second).first << " (" << std::setw(8) << std::fixed << std::setprecision(4) << evtpercent << "%)"
944  << " | tracks collected: " << std::setw(8) << (it->second).second << " (" << std::setw(8) << std::fixed
945  << std::setprecision(4) << trkpercent << "%)";
946 
947  tksByTrigger_->SetBinContent(i, trkpercent);
948  tksByTrigger_->GetXaxis()->SetBinLabel(i, (it->first).c_str());
949 
950  evtsByTrigger_->SetBinContent(i, evtpercent);
951  evtsByTrigger_->GetXaxis()->SetBinLabel(i, (it->first).c_str());
952  }
953 
954  TFileDirectory RunFeatures = outfile_->mkdir("RunFeatures");
955  h_runStartTimes = RunFeatures.make<TH1I>(
956  "runStartTimes", "run start times", runNumbersTimesLog_.size(), 0, runNumbersTimesLog_.size());
957  h_runEndTimes =
958  RunFeatures.make<TH1I>("runEndTimes", "run end times", runNumbersTimesLog_.size(), 0, runNumbersTimesLog_.size());
959 
960  unsigned int count = 1;
961  for (const auto& run : runNumbersTimesLog_) {
962  // strip down the microseconds
963  h_runStartTimes->SetBinContent(count, run.second.first / 10e6);
964  h_runStartTimes->GetXaxis()->SetBinLabel(count, (std::to_string(run.first)).c_str());
965 
966  h_runEndTimes->SetBinContent(count, run.second.second / 10e6);
967  h_runEndTimes->GetXaxis()->SetBinLabel(count, (std::to_string(run.first)).c_str());
968 
969  count++;
970  }
971 
972  // resolutions
973 
977 
981 
985 
986  // pulls
987 
991 
995 
999 }
1000 
1001 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
1003  //The following says we do not know what parameters are allowed so do no validation
1004  // Please change this to state exactly what you do use, even if it is no parameters
1006  desc.setUnknown();
1007  descriptions.addDefault(desc);
1008 }
1009 
1010 //*************************************************************
1011 std::pair<long long, long long> SplitVertexResolution::getRunTime(const edm::EventSetup& iSetup) const
1012 //*************************************************************
1013 {
1014  const auto& runInfo = iSetup.getData(runInfoToken_);
1015  if (debug_) {
1016  edm::LogInfo("SplitVertexResolution") << runInfo.m_start_time_str << " " << runInfo.m_stop_time_str << std::endl;
1017  }
1018  return std::make_pair(runInfo.m_start_time_ll, runInfo.m_stop_time_ll);
1019 }
1020 
1021 //*************************************************************
1022 void SplitVertexResolution::fillTrendPlotByIndex(TH1F* trendPlot, std::vector<TH1F*>& h, PVValHelper::estimator fitPar_)
1023 //*************************************************************
1024 {
1025  for (auto iterator = h.begin(); iterator != h.end(); iterator++) {
1026  unsigned int bin = std::distance(h.begin(), iterator) + 1;
1027  statmode::fitParams myFit = fitResiduals((*iterator));
1028 
1029  switch (fitPar_) {
1030  case PVValHelper::MEAN: {
1031  float mean_ = myFit.first.value();
1032  float meanErr_ = myFit.first.error();
1033  trendPlot->SetBinContent(bin, mean_);
1034  trendPlot->SetBinError(bin, meanErr_);
1035  break;
1036  }
1037  case PVValHelper::WIDTH: {
1038  float width_ = myFit.second.value();
1039  float widthErr_ = myFit.second.error();
1040  trendPlot->SetBinContent(bin, width_);
1041  trendPlot->SetBinError(bin, widthErr_);
1042  break;
1043  }
1044  case PVValHelper::MEDIAN: {
1045  float median_ = PVValHelper::getMedian((*iterator)).value();
1046  float medianErr_ = PVValHelper::getMedian((*iterator)).error();
1047  trendPlot->SetBinContent(bin, median_);
1048  trendPlot->SetBinError(bin, medianErr_);
1049  break;
1050  }
1051  case PVValHelper::MAD: {
1052  float mad_ = PVValHelper::getMAD((*iterator)).value();
1053  float madErr_ = PVValHelper::getMAD((*iterator)).error();
1054  trendPlot->SetBinContent(bin, mad_);
1055  trendPlot->SetBinError(bin, madErr_);
1056  break;
1057  }
1058  default:
1059  edm::LogWarning("SplitVertexResolution")
1060  << "fillTrendPlotByIndex() " << fitPar_ << " unknown estimator!" << std::endl;
1061  break;
1062  }
1063  }
1064 }
1065 
1066 //*************************************************************
1068 //*************************************************************
1069 {
1070  if (hist->GetEntries() < 10) {
1071  LogDebug("SplitVertexResolution") << "hist name: " << hist->GetName() << " has less than 10 entries" << std::endl;
1072  return std::make_pair(Measurement1D(0., 0.), Measurement1D(0., 0.));
1073  }
1074 
1075  float maxHist = hist->GetXaxis()->GetXmax();
1076  float minHist = hist->GetXaxis()->GetXmin();
1077  float mean = hist->GetMean();
1078  float sigma = hist->GetRMS();
1079 
1080  if (TMath::IsNaN(mean) || TMath::IsNaN(sigma)) {
1081  mean = 0;
1082  //sigma= - hist->GetXaxis()->GetBinLowEdge(1) + hist->GetXaxis()->GetBinLowEdge(hist->GetNbinsX()+1);
1083  sigma = -minHist + maxHist;
1084  edm::LogWarning("SplitVertexResolution")
1085  << "FitPVResiduals::fitResiduals(): histogram" << hist->GetName() << " mean or sigma are NaN!!" << std::endl;
1086  }
1087 
1088  TF1 func("tmp", "gaus", mean - 2. * sigma, mean + 2. * sigma);
1089  if (0 == hist->Fit(&func, "QNR")) { // N: do not blow up file by storing fit!
1090  mean = func.GetParameter(1);
1091  sigma = func.GetParameter(2);
1092 
1093  if (!singleTime) {
1094  // second fit: three sigma of first fit around mean of first fit
1095  func.SetRange(std::max(mean - 3 * sigma, minHist), std::min(mean + 3 * sigma, maxHist));
1096  // I: integral gives more correct results if binning is too wide
1097  // L: Likelihood can treat empty bins correctly (if hist not weighted...)
1098  if (0 == hist->Fit(&func, "Q0LR")) {
1099  if (hist->GetFunction(func.GetName())) { // Take care that it is later on drawn:
1100  hist->GetFunction(func.GetName())->ResetBit(TF1::kNotDraw);
1101  }
1102  }
1103  }
1104  }
1105 
1106  float res_mean = func.GetParameter(1);
1107  float res_width = func.GetParameter(2);
1108 
1109  float res_mean_err = func.GetParError(1);
1110  float res_width_err = func.GetParError(2);
1111 
1112  Measurement1D resultM(res_mean, res_mean_err);
1113  Measurement1D resultW(res_width, res_width_err);
1114 
1115  statmode::fitParams result = std::make_pair(resultM, resultW);
1116  return result;
1117 }
1118 
1119 //*************************************************************
1121 //*************************************************************
1122 {
1123  float mean = hist->GetMean();
1124  float sigma = hist->GetRMS();
1125 
1126  TF1 func("tmp", "gaus", mean - 1.5 * sigma, mean + 1.5 * sigma);
1127  if (0 == hist->Fit(&func, "QNR")) { // N: do not blow up file by storing fit!
1128  mean = func.GetParameter(1);
1129  sigma = func.GetParameter(2);
1130  // second fit: three sigma of first fit around mean of first fit
1131  func.SetRange(mean - 2 * sigma, mean + 2 * sigma);
1132  // I: integral gives more correct results if binning is too wide
1133  // L: Likelihood can treat empty bins correctly (if hist not weighted...)
1134  if (0 == hist->Fit(&func, "Q0LR")) {
1135  if (hist->GetFunction(func.GetName())) { // Take care that it is later on drawn:
1136  hist->GetFunction(func.GetName())->ResetBit(TF1::kNotDraw);
1137  }
1138  }
1139  }
1140 
1141  float res_mean = func.GetParameter(1);
1142  float res_width = func.GetParameter(2);
1143 
1144  float res_mean_err = func.GetParError(1);
1145  float res_width_err = func.GetParError(2);
1146 
1147  Measurement1D resultM(res_mean, res_mean_err);
1148  Measurement1D resultW(res_width, res_width_err);
1149 
1150  statmode::fitParams result = std::make_pair(resultM, resultW);
1151  return result;
1152 }
1153 
1154 //*************************************************************
1155 template <std::size_t SIZE>
1156 bool SplitVertexResolution::checkBinOrdering(std::array<float, SIZE>& bins)
1157 //*************************************************************
1158 {
1159  int i = 1;
1160 
1161  if (std::is_sorted(bins.begin(), bins.end())) {
1162  return true;
1163  } else {
1164  for (const auto& bin : bins) {
1165  edm::LogInfo("SplitVertexResolution") << "bin: " << i << " : " << bin << std::endl;
1166  i++;
1167  }
1168  edm::LogInfo("SplitVertexResolution") << "--------------------------------" << std::endl;
1169  return false;
1170  }
1171 }
1172 
1173 //define this as a plug-in
AdaptiveVertexFitter
Definition: AdaptiveVertexFitter.h:29
SplitVertexResolution::fitResiduals_v0
statmode::fitParams fitResiduals_v0(TH1 *hist)
Definition: SplitVertexResolution.cc:1120
SplitVertexResolution::outfile_
edm::Service< TFileService > outfile_
Definition: SplitVertexResolution.cc:136
pvCand::zErr_subVtx2
float zErr_subVtx2
Definition: pvTree.h:39
fftjetvertexadder_cfi.errZ
errZ
Definition: fftjetvertexadder_cfi.py:39
SplitVertexResolution::h_runFromConfig
TH1I * h_runFromConfig
Definition: SplitVertexResolution.cc:139
counter
Definition: counter.py:1
SplitVertexResolution::p_pullX_vsNvtx
TH1F * p_pullX_vsNvtx
Definition: SplitVertexResolution.cc:236
alignBH_cfg.fixed
fixed
Definition: alignBH_cfg.py:54
SplitVertexResolution::triggerResultsTag_
edm::InputTag triggerResultsTag_
Definition: SplitVertexResolution.cc:119
PDWG_EXOHSCP_cff.tracks
tracks
Definition: PDWG_EXOHSCP_cff.py:28
Measurement1D
Definition: Measurement1D.h:11
electrons_cff.bool
bool
Definition: electrons_cff.py:372
EDAnalyzer.h
mps_fire.i
i
Definition: mps_fire.py:355
SplitVertexResolution::checkBinOrdering
bool checkBinOrdering(std::array< float, SIZE > &bins)
Definition: SplitVertexResolution.cc:1156
SplitVertexResolution::h_pullX_sumPt_
std::vector< TH1F * > h_pullX_sumPt_
Definition: SplitVertexResolution.cc:216
SplitVertexResolution::h_PVCL_subVtx2
TH1F * h_PVCL_subVtx2
Definition: SplitVertexResolution.cc:175
SplitVertexResolution::h_diffX
TH1F * h_diffX
Definition: SplitVertexResolution.cc:145
SiStripPI::mean
Definition: SiStripPayloadInspectorHelper.h:169
funct::false
false
Definition: Factorize.h:34
filterCSVwithJSON.copy
copy
Definition: filterCSVwithJSON.py:36
TriggerResults.h
ESInputTag
HLTBitAnalyser_cfi.hltresults
hltresults
Definition: HLTBitAnalyser_cfi.py:13
fftjetvertexadder_cfi.errY
errY
Definition: fftjetvertexadder_cfi.py:38
TransientVertex::isValid
bool isValid() const
Definition: TransientVertex.h:195
SplitVertexResolution::p_resolY_vsNvtx
TH1F * p_resolY_vsNvtx
Definition: SplitVertexResolution.cc:212
PVValHelper::WIDTH
Definition: PVValidationHelpers.h:44
edm::Run
Definition: Run.h:45
min
T min(T a, T b)
Definition: MathUtil.h:58
reco::Vertex::z
double z() const
z coordinate
Definition: Vertex.h:120
SplitVertexResolution::p_resolY_vsNtracks
TH1F * p_resolY_vsNtracks
Definition: SplitVertexResolution.cc:208
edm::EDGetTokenT< reco::VertexCollection >
SplitVertexResolution::tracksToken_
edm::EDGetTokenT< reco::TrackCollection > tracksToken_
Definition: SplitVertexResolution.cc:117
SplitVertexResolution::h_runNumber
TH1F * h_runNumber
Definition: SplitVertexResolution.cc:177
TFileDirectory::make
T * make(const Args &... args) const
make new ROOT object
Definition: TFileDirectory.h:53
edm
HLT enums.
Definition: AlignableModifier.h:19
Measurement1D::value
double value() const
Definition: Measurement1D.h:25
pvCand::y_subVtx2
float y_subVtx2
Definition: pvTree.h:34
SplitVertexResolution::nVtxBins_
static const int nVtxBins_
Definition: SplitVertexResolution.cc:252
pvEvent::pvs
std::vector< pvCand > pvs
Definition: pvTree.h:62
SplitVertexResolution::minVtxNdf_
double minVtxNdf_
Definition: SplitVertexResolution.cc:126
SplitVertexResolution::h_ntrks
TH1F * h_ntrks
Definition: SplitVertexResolution.cc:161
h
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
Definition: L1TUtmAlgorithmRcd.h:4
pvCand::CL_subVtx2
float CL_subVtx2
Definition: pvTree.h:43
reco::VertexCollection
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::LogInfo
Definition: MessageLogger.h:254
SplitVertexResolution::bookResidualsHistogram
std::vector< TH1F * > bookResidualsHistogram(TFileDirectory dir, unsigned int theNOfBins, TString resType, TString varType)
Definition: SplitVertexResolution.cc:887
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
pvCand::x_subVtx1
float x_subVtx1
Definition: pvTree.h:23
pvCand::minW_subVtx2
float minW_subVtx2
Definition: pvTree.h:46
SplitVertexResolution::triggerMap_
std::map< std::string, std::pair< int, int > > triggerMap_
Definition: SplitVertexResolution.cc:255
PVValHelper::fillByIndex
void fillByIndex(std::vector< TH1F * > &h, unsigned int index, double x, std::string tag="")
Definition: PVValidationHelpers.cc:42
pvCand::n_subVtx2
int n_subVtx2
Definition: pvTree.h:32
SplitVertexResolution::h_PVCL_subVtx1
TH1F * h_PVCL_subVtx1
Definition: SplitVertexResolution.cc:174
muonTagProbeFilters_cff.allTracks
allTracks
Definition: muonTagProbeFilters_cff.py:22
HLT_2018_cff.distance
distance
Definition: HLT_2018_cff.py:6417
SplitVertexResolution::p_pullZ_vsNvtx
TH1F * p_pullZ_vsNvtx
Definition: SplitVertexResolution.cc:238
edm::second
U second(std::pair< T, U > const &p)
Definition: ParameterSet.cc:215
SplitVertexResolution::h_OrigVertexErrZ
TH1F * h_OrigVertexErrZ
Definition: SplitVertexResolution.cc:151
pvEvent::nVtx
int nVtx
Definition: pvTree.h:60
SplitVertexResolution::h_resolX_Nvtx_
std::vector< TH1F * > h_resolX_Nvtx_
Definition: SplitVertexResolution.cc:199
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
SplitVertexResolution::h_sumPt1
TH1F * h_sumPt1
Definition: SplitVertexResolution.cc:165
SplitVertexResolution::mypT_bins_
std::array< float, nPtBins_+1 > mypT_bins_
Definition: SplitVertexResolution.cc:247
TFileDirectory
Definition: TFileDirectory.h:24
VertexDistance3D.h
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::one::EDAnalyzer
Definition: EDAnalyzer.h:30
tools.TF1
TF1
Definition: tools.py:23
SplitVertexResolution::p_resolX_vsNvtx
TH1F * p_resolX_vsNvtx
Definition: SplitVertexResolution.cc:211
watchdog.const
const
Definition: watchdog.py:83
AdaptiveVertexFitter::vertex
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const override
Definition: AdaptiveVertexFitter.cc:158
edm::Handle< edm::TriggerResults >
pvCand::z_origVtx
float z_origVtx
Definition: pvTree.h:16
parallelization.uint
uint
Definition: parallelization.py:124
SplitVertexResolution::h_resolZ_Ntracks_
std::vector< TH1F * > h_resolZ_Ntracks_
Definition: SplitVertexResolution.cc:197
SplitVertexResolution::h_errY
TH1F * h_errY
Definition: SplitVertexResolution.cc:154
pvCand::zErr_origVtx
float zErr_origVtx
Definition: pvTree.h:20
dqmdumpme.first
first
Definition: dqmdumpme.py:55
SplitVertexResolution::h_pullZ_sumPt_
std::vector< TH1F * > h_pullZ_sumPt_
Definition: SplitVertexResolution.cc:218
SplitVertexResolution::p_pullZ_vsNtracks
TH1F * p_pullZ_vsNtracks
Definition: SplitVertexResolution.cc:234
SplitVertexResolution::p_resolZ_vsSumPt
TH1F * p_resolZ_vsSumPt
Definition: SplitVertexResolution.cc:205
pvCand::sumPt_subVtx2
float sumPt_subVtx2
Definition: pvTree.h:40
edm::Ref< TrackCollection >
TtFullHadEvtBuilder_cfi.sumPt
sumPt
Definition: TtFullHadEvtBuilder_cfi.py:38
SplitVertexResolution::h_resolZ_Nvtx_
std::vector< TH1F * > h_resolZ_Nvtx_
Definition: SplitVertexResolution.cc:201
pvEvent
Definition: pvTree.h:54
SplitVertexResolution::runNumbersTimesLog_
std::map< unsigned int, std::pair< long long, long long > > runNumbersTimesLog_
Definition: SplitVertexResolution.cc:141
SplitVertexResolution::h_pullY_Nvtx_
std::vector< TH1F * > h_pullY_Nvtx_
Definition: SplitVertexResolution.cc:225
SplitVertexResolution::ievt
int ievt
Definition: SplitVertexResolution.cc:103
SplitVertexResolution::p_pullY_vsNvtx
TH1F * p_pullY_vsNvtx
Definition: SplitVertexResolution.cc:237
SplitVertexResolution::h_errX
TH1F * h_errX
Definition: SplitVertexResolution.cc:153
SplitVertexResolution::h_runStartTimes
TH1I * h_runStartTimes
Definition: SplitVertexResolution.cc:142
SplitVertexResolution::p_pullX_vsNtracks
TH1F * p_pullX_vsNtracks
Definition: SplitVertexResolution.cc:232
accept
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
MakerMacros.h
pvCand::n_subVtx1
int n_subVtx1
Definition: pvTree.h:22
pvCand::y_origVtx
float y_origVtx
Definition: pvTree.h:15
down
Definition: BitonicSort.h:7
RunInfo
Definition: RunInfo.h:18
SplitVertexResolution::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: SplitVertexResolution.cc:1002
h
Track.h
pvEvent::eventNumber
int eventNumber
Definition: pvTree.h:58
SplitVertexResolution::h_wTrks1
TH1F * h_wTrks1
Definition: SplitVertexResolution.cc:168
pvCand::sumPt_subVtx1
float sumPt_subVtx1
Definition: pvTree.h:30
reco::Vertex::isFake
bool isFake() const
Definition: Vertex.h:76
PVValHelper::estimator
estimator
Definition: PVValidationHelpers.h:44
TrackFwd.h
pvCand::zErr_subVtx1
float zErr_subVtx1
Definition: pvTree.h:29
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
SplitVertexResolution::itrks
int itrks
Definition: SplitVertexResolution.cc:104
reco::Vertex::xError
double xError() const
error on x
Definition: Vertex.h:124
SplitVertexResolution::h_lumiFromConfig
TH1F * h_lumiFromConfig
Definition: SplitVertexResolution.cc:138
compare.hist
hist
Definition: compare.py:376
SplitVertexResolution::h_minWTrks1
TH1F * h_minWTrks1
Definition: SplitVertexResolution.cc:171
SplitVertexResolution::h_resolZ_sumPt_
std::vector< TH1F * > h_resolZ_sumPt_
Definition: SplitVertexResolution.cc:193
SplitVertexResolution::tree_
TTree * tree_
Definition: SplitVertexResolution.cc:243
Service.h
SplitVertexResolution::p_resolX_vsNtracks
TH1F * p_resolX_vsNtracks
Definition: SplitVertexResolution.cc:207
pvEvent::luminosityBlockNumber
int luminosityBlockNumber
Definition: pvTree.h:57
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
SplitVertexResolution::myNVtx_bins_
std::array< float, nVtxBins_+1 > myNVtx_bins_
Definition: SplitVertexResolution.cc:253
SplitVertexResolution::transientTrackBuilderToken_
edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > transientTrackBuilderToken_
Definition: SplitVertexResolution.cc:123
reco::Vertex::tracks_end
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:73
SplitVertexResolution::runControl_
bool runControl_
Definition: SplitVertexResolution.cc:131
TransientTrackRecord
Definition: TransientTrackRecord.h:11
reco::Track
Definition: Track.h:27
pvCand::minW_subVtx1
float minW_subVtx1
Definition: pvTree.h:45
Event
pvCand
Definition: pvTree.h:9
pvTree.h
SplitVertexResolution::fillTrendPlotByIndex
void fillTrendPlotByIndex(TH1F *trendPlot, std::vector< TH1F * > &h, PVValHelper::estimator fitPar_)
Definition: SplitVertexResolution.cc:1022
SplitVertexResolution::mysorter
static bool mysorter(reco::Track i, reco::Track j)
Definition: SplitVertexResolution.cc:79
SplitVertexResolution::h_nOfflineVertices
TH1I * h_nOfflineVertices
Definition: SplitVertexResolution.cc:179
h
SplitVertexResolution::h_pullY_Ntracks_
std::vector< TH1F * > h_pullY_Ntracks_
Definition: SplitVertexResolution.cc:221
pvCand::z_subVtx1
float z_subVtx1
Definition: pvTree.h:25
SplitVertexResolution
Definition: SplitVertexResolution.cc:73
SplitVertexResolution::debug_
bool debug_
Definition: SplitVertexResolution.cc:111
dumpRecoGeometry_cfg.varType
varType
Definition: dumpRecoGeometry_cfg.py:8
pvCand::x_origVtx
float x_origVtx
Definition: pvTree.h:14
Measurement1D::error
double error() const
Definition: Measurement1D.h:27
GlobalTrackingGeometryRecord.h
RunInfoRcd
Definition: RunSummaryRcd.h:26
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
reco::Vertex::zError
double zError() const
error on z
Definition: Vertex.h:128
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
Vertex.h
EgHLTOffTrigSelection_cfi.trigName
trigName
Definition: EgHLTOffTrigSelection_cfi.py:8
edm::LogWarning
Definition: MessageLogger.h:141
TFileService.h
SplitVertexResolution::p_pullY_vsNtracks
TH1F * p_pullY_vsNtracks
Definition: SplitVertexResolution.cc:233
SplitVertexResolution::fitResiduals
statmode::fitParams fitResiduals(TH1 *hist, bool singleTime=false)
Definition: SplitVertexResolution.cc:1067
TFileService::mkdir
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
TransientTrackBuilder.h
SplitVertexResolution::h_resolY_Ntracks_
std::vector< TH1F * > h_resolY_Ntracks_
Definition: SplitVertexResolution.cc:196
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:670
PVValHelper::getMedian
Measurement1D getMedian(TH1F *histo)
Definition: PVValidationHelpers.cc:178
edm::ParameterSet
Definition: ParameterSet.h:36
SplitVertexResolution::runInfoToken_
edm::ESGetToken< RunInfo, RunInfoRcd > runInfoToken_
Definition: SplitVertexResolution.cc:124
edm::Transition
Transition
Definition: Transition.h:12
SplitVertexResolution::runControlNumbers_
std::vector< unsigned int > runControlNumbers_
Definition: SplitVertexResolution.cc:132
edm::LogError
Definition: MessageLogger.h:183
SplitVertexResolution::h_diffZ
TH1F * h_diffZ
Definition: SplitVertexResolution.cc:147
SplitVertexResolution::h_pullX_Ntracks_
std::vector< TH1F * > h_pullX_Ntracks_
Definition: SplitVertexResolution.cc:220
SplitVertexResolution::h_resolY_Nvtx_
std::vector< TH1F * > h_resolY_Nvtx_
Definition: SplitVertexResolution.cc:200
SplitVertexResolution::p_resolZ_vsNtracks
TH1F * p_resolZ_vsNtracks
Definition: SplitVertexResolution.cc:209
SplitVertexResolution::h_minWTrks2
TH1F * h_minWTrks2
Definition: SplitVertexResolution.cc:172
Event.h
SplitVertexResolution::nTrackBins_
static const int nTrackBins_
Definition: SplitVertexResolution.cc:249
SplitVertexResolution::p_resolZ_vsNvtx
TH1F * p_resolZ_vsNvtx
Definition: SplitVertexResolution.cc:213
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
reco::Vertex::tracks_begin
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:71
PVValidationHelpers.h
pvCand::z_subVtx2
float z_subVtx2
Definition: pvTree.h:35
KineDebug3::count
void count()
Definition: KinematicConstrainedVertexUpdatorT.h:21
SplitVertexResolution::h_pullZ_Ntracks_
std::vector< TH1F * > h_pullZ_Ntracks_
Definition: SplitVertexResolution.cc:222
SplitVertexResolution::intLumi_
double intLumi_
Definition: SplitVertexResolution.cc:110
reco::Vertex::x
double x() const
x coordinate
Definition: Vertex.h:116
SplitVertexResolution::h_sumPt2
TH1F * h_sumPt2
Definition: SplitVertexResolution.cc:166
SplitVertexResolution::nPtBins_
static const int nPtBins_
Definition: SplitVertexResolution.cc:246
edm::Service< TFileService >
PVValHelper::MAD
Definition: PVValidationHelpers.h:44
createfilelist.int
int
Definition: createfilelist.py:10
SplitVertexResolution::h_OrigVertexErrY
TH1F * h_OrigVertexErrY
Definition: SplitVertexResolution.cc:150
TriggerNames.h
reco::Vertex::tracksSize
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:69
iEvent
int iEvent
Definition: GenABIO.cc:224
SplitVertexResolution::triggerResultsToken_
edm::EDGetTokenT< edm::TriggerResults > triggerResultsToken_
Definition: SplitVertexResolution.cc:120
SplitVertexResolution::myNTrack_bins_
std::array< float, nTrackBins_+1 > myNTrack_bins_
Definition: SplitVertexResolution.cc:250
SplitVertexResolution::p_pullY_vsSumPt
TH1F * p_pullY_vsSumPt
Definition: SplitVertexResolution.cc:229
edm::LogVerbatim
Definition: MessageLogger.h:297
TransientTrackBuilder
Definition: TransientTrackBuilder.h:16
edm::ParameterSetDescription::setUnknown
void setUnknown()
Definition: ParameterSetDescription.cc:39
SplitVertexResolution::h_sumPt
TH1F * h_sumPt
Definition: SplitVertexResolution.cc:162
SplitVertexResolution::h_runEndTimes
TH1I * h_runEndTimes
Definition: SplitVertexResolution.cc:143
counter
static std::atomic< unsigned int > counter
Definition: SharedResourceNames.cc:15
SplitVertexResolution::h_OrigVertexErrX
TH1F * h_OrigVertexErrX
Definition: SplitVertexResolution.cc:149
TransientVertex
Definition: TransientVertex.h:18
edm::EventSetup
Definition: EventSetup.h:57
PVValHelper::generateBins
std::vector< float > generateBins(int n, float start, float range)
Definition: PVValidationHelpers.cc:161
TransientTrackRecord.h
TrackCollections2monitor_cff.func
func
Definition: TrackCollections2monitor_cff.py:359
SplitVertexResolution::tracksTag_
edm::InputTag tracksTag_
Definition: SplitVertexResolution.cc:116
PVValHelper::MEAN
Definition: PVValidationHelpers.h:44
SplitVertexResolution::h_pullY_sumPt_
std::vector< TH1F * > h_pullY_sumPt_
Definition: SplitVertexResolution.cc:217
edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord >
pvCand::x_subVtx2
float x_subVtx2
Definition: pvTree.h:33
SplitVertexResolution::getRunTime
std::pair< long long, long long > getRunTime(const edm::EventSetup &iSetup) const
Definition: SplitVertexResolution.cc:1011
InputTag.h
edm::TriggerNames::triggerName
std::string const & triggerName(unsigned int index) const
Definition: TriggerNames.cc:22
edm::EventSetup::getData
bool getData(T &iHolder) const
Definition: EventSetup.h:113
alignCSCRings.r
r
Definition: alignCSCRings.py:93
newFWLiteAna.bin
bin
Definition: newFWLiteAna.py:161
SplitVertexResolution::beginJob
void beginJob() override
Definition: SplitVertexResolution.cc:639
SplitVertexResolution::event_
pvEvent event_
Definition: SplitVertexResolution.cc:242
SiPixelPhase1Clusters_cfi.e3
e3
Definition: SiPixelPhase1Clusters_cfi.py:9
SplitVertexResolution::p_resolY_vsSumPt
TH1F * p_resolY_vsSumPt
Definition: SplitVertexResolution.cc:204
SplitVertexResolution::h_resolX_Ntracks_
std::vector< TH1F * > h_resolX_Ntracks_
Definition: SplitVertexResolution.cc:195
VertexFwd.h
pvCand::xErr_origVtx
float xErr_origVtx
Definition: pvTree.h:18
SplitVertexResolution::h_resolX_sumPt_
std::vector< TH1F * > h_resolX_sumPt_
Definition: SplitVertexResolution.cc:191
pvCand::yErr_subVtx1
float yErr_subVtx1
Definition: pvTree.h:28
RunInfo.h
statmode
Definition: SplitVertexResolution.cc:65
ESInputTag.h
writedatasetfile.run
run
Definition: writedatasetfile.py:27
reco::Vertex::y
double y() const
y coordinate
Definition: Vertex.h:118
SplitVertexResolution::endJob
void endJob() override
Definition: SplitVertexResolution.cc:919
reco::TransientTrack
Definition: TransientTrack.h:19
pvCand::nTrks
int nTrks
Definition: pvTree.h:11
SplitVertexResolution::h_resolY_sumPt_
std::vector< TH1F * > h_resolY_sumPt_
Definition: SplitVertexResolution.cc:192
Frameworkfwd.h
TransientVertex::trackWeight
float trackWeight(const reco::TransientTrack &track) const
Definition: TransientVertex.cc:241
SplitVertexResolution::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: SplitVertexResolution.cc:290
SplitVertexResolution::p_pullX_vsSumPt
TH1F * p_pullX_vsSumPt
Definition: SplitVertexResolution.cc:228
pvCand::xErr_subVtx1
float xErr_subVtx1
Definition: pvTree.h:27
PVValHelper::MEDIAN
Definition: PVValidationHelpers.h:44
SplitVertexResolution::h_pullX
TH1F * h_pullX
Definition: SplitVertexResolution.cc:157
TFileService::kSharedResource
static const std::string kSharedResource
Definition: TFileService.h:76
TransientVertex::originalTracks
std::vector< reco::TransientTrack > const & originalTracks() const
Definition: TransientVertex.h:200
PVValHelper::getMAD
Measurement1D getMAD(TH1F *histo)
Definition: PVValidationHelpers.cc:194
SplitVertexResolution::h_avgSumPt
TH1F * h_avgSumPt
Definition: SplitVertexResolution.cc:163
edm::TriggerNames
Definition: TriggerNames.h:55
SplitVertexResolution::h_nFinalVertices
TH1I * h_nFinalVertices
Definition: SplitVertexResolution.cc:182
pvCand::ipos
int ipos
Definition: pvTree.h:12
EventSetup.h
SplitVertexResolution::minVtxWgt_
double minVtxWgt_
Definition: SplitVertexResolution.cc:127
SplitVertexResolution::p_resolX_vsSumPt
TH1F * p_resolX_vsSumPt
Definition: SplitVertexResolution.cc:203
SplitVertexResolution::h_diffY
TH1F * h_diffY
Definition: SplitVertexResolution.cc:146
SplitVertexResolution::h_nNonFakeVertices
TH1I * h_nNonFakeVertices
Definition: SplitVertexResolution.cc:181
TransientTrackBuilder::build
reco::TransientTrack build(const reco::Track *p) const
Definition: TransientTrackBuilder.cc:20
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:30
SplitVertexResolution::beginRun
void beginRun(edm::Run const &iEvent, edm::EventSetup const &) override
Definition: SplitVertexResolution.cc:621
reco::Vertex::nTracks
unsigned int nTracks(float minWeight=0.5) const
Returns the number of tracks in the vertex with weight above minWeight.
Definition: Vertex.cc:146
mps_fire.result
result
Definition: mps_fire.py:303
reco::Vertex::yError
double yError() const
error on y
Definition: Vertex.h:126
SplitVertexResolution::h_pullX_Nvtx_
std::vector< TH1F * > h_pullX_Nvtx_
Definition: SplitVertexResolution.cc:224
SplitVertexResolution::SplitVertexResolution
SplitVertexResolution(const edm::ParameterSet &)
Definition: SplitVertexResolution.cc:258
genParticles_cff.map
map
Definition: genParticles_cff.py:11
trigObjTnPSource_cfi.bins
bins
Definition: trigObjTnPSource_cfi.py:20
EventSetup
pvCand::yErr_subVtx2
float yErr_subVtx2
Definition: pvTree.h:38
ParameterSet.h
SplitVertexResolution::pvsTag_
edm::InputTag pvsTag_
Definition: SplitVertexResolution.cc:113
pvCand::CL_subVtx1
float CL_subVtx1
Definition: pvTree.h:42
SplitVertexResolution::h_nVertices
TH1I * h_nVertices
Definition: SplitVertexResolution.cc:180
pvCand::xErr_subVtx2
float xErr_subVtx2
Definition: pvTree.h:37
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
SplitVertexResolution::h_pullZ
TH1F * h_pullZ
Definition: SplitVertexResolution.cc:159
SplitVertexResolution::beginEvent
virtual void beginEvent() final
Definition: SplitVertexResolution.cc:616
HLTObjectsMonitor_cfi.TriggerResults
TriggerResults
Definition: HLTObjectsMonitor_cfi.py:9
SplitVertexResolution::h_pullY
TH1F * h_pullY
Definition: SplitVertexResolution.cc:158
edm::Event
Definition: Event.h:73
AdaptiveVertexFitter.h
SplitVertexResolution::engine_
std::mt19937 engine_
Definition: SplitVertexResolution.cc:240
edm::ConfigurationDescriptions::addDefault
void addDefault(ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:99
pvCand::y_subVtx1
float y_subVtx1
Definition: pvTree.h:24
reco::Vertex::ndof
double ndof() const
Definition: Vertex.h:110
SplitVertexResolution::p_pullZ_vsSumPt
TH1F * p_pullZ_vsSumPt
Definition: SplitVertexResolution.cc:230
fftjetvertexadder_cfi.errX
errX
Definition: fftjetvertexadder_cfi.py:37
edm::InputTag
Definition: InputTag.h:15
pvCand::yErr_origVtx
float yErr_origVtx
Definition: pvTree.h:19
SplitVertexResolution::h_wTrks2
TH1F * h_wTrks2
Definition: SplitVertexResolution.cc:169
reco::TrackCollection
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
SplitVertexResolution::h_pullZ_Nvtx_
std::vector< TH1F * > h_pullZ_Nvtx_
Definition: SplitVertexResolution.cc:226
up
Definition: BitonicSort.h:7
reco::Vertex
Definition: Vertex.h:35
statmode::fitParams
std::pair< Measurement1D, Measurement1D > fitParams
Definition: SplitVertexResolution.cc:66
pvEvent::runNumber
int runNumber
Definition: pvTree.h:56
SplitVertexResolution::h_errZ
TH1F * h_errZ
Definition: SplitVertexResolution.cc:155
TFileService::make
T * make(const Args &... args) const
make new ROOT object
Definition: TFileService.h:64
SplitVertexResolution::storeNtuple_
bool storeNtuple_
Definition: SplitVertexResolution.cc:107
SplitVertexResolution::tksByTrigger_
TH1D * tksByTrigger_
Definition: SplitVertexResolution.cc:186
SplitVertexResolution::evtsByTrigger_
TH1D * evtsByTrigger_
Definition: SplitVertexResolution.cc:187
SplitVertexResolution::pvsToken_
edm::EDGetTokenT< reco::VertexCollection > pvsToken_
Definition: SplitVertexResolution.cc:114
DeadROC_duringRun.dir
dir
Definition: DeadROC_duringRun.py:23
SplitVertexResolution::endRun
void endRun(edm::Run const &, edm::EventSetup const &) override
Definition: SplitVertexResolution.cc:87
pwdgSkimBPark_cfi.vertices
vertices
Definition: pwdgSkimBPark_cfi.py:7
SplitVertexResolution::cmToUm
static constexpr double cmToUm
Definition: SplitVertexResolution.cc:134
SplitVertexResolution::~SplitVertexResolution
~SplitVertexResolution() override
Definition: SplitVertexResolution.cc:287
Run