Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-07-18 00:48:14

0001 #include <string>
0002 
0003 #include "FWCore/Framework/interface/MakerMacros.h"
0004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0005 #include "FWCore/ServiceRegistry/interface/Service.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 
0008 #include "DQMServices/Core/interface/DQMEDHarvester.h"
0009 #include "DQMServices/Core/interface/DQMStore.h"
0010 
0011 #include "DataFormats/ForwardDetId/interface/ETLDetId.h"
0012 
0013 class Primary4DVertexHarvester : public DQMEDHarvester {
0014 public:
0015   explicit Primary4DVertexHarvester(const edm::ParameterSet& iConfig);
0016   ~Primary4DVertexHarvester() override;
0017 
0018   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0019 
0020 protected:
0021   void dqmEndJob(DQMStore::IBooker&, DQMStore::IGetter&) override;
0022 
0023 private:
0024   void computeEfficiency1D(MonitorElement* num, MonitorElement* den, MonitorElement* result);
0025 
0026   void incrementME(MonitorElement* base, MonitorElement* toBeAdded);
0027 
0028   const std::string folder_;
0029 
0030   // --- Histograms
0031   MonitorElement* meTPPtSelEff_;
0032   MonitorElement* meTPEtaSelEff_;
0033   MonitorElement* meTPPtMatchEff_;
0034   MonitorElement* meTPEtaMatchEff_;
0035 
0036   MonitorElement* meBarrelTruePi_;
0037   MonitorElement* meBarrelTrueK_;
0038   MonitorElement* meBarrelTrueP_;
0039 
0040   MonitorElement* meEndcapTruePi_;
0041   MonitorElement* meEndcapTrueK_;
0042   MonitorElement* meEndcapTrueP_;
0043 
0044   MonitorElement* meBarrelAsPi_;
0045   MonitorElement* meBarrelAsK_;
0046   MonitorElement* meBarrelAsP_;
0047   MonitorElement* meBarrelNoPID_;
0048 
0049   MonitorElement* meEndcapAsPi_;
0050   MonitorElement* meEndcapAsK_;
0051   MonitorElement* meEndcapAsP_;
0052   MonitorElement* meEndcapNoPID_;
0053 
0054   MonitorElement* meBarrelPIDPiAsPiEff_;
0055   MonitorElement* meBarrelPIDPiAsKEff_;
0056   MonitorElement* meBarrelPIDPiAsPEff_;
0057   MonitorElement* meBarrelPIDPiNoPIDEff_;
0058 
0059   MonitorElement* meBarrelPIDKAsPiEff_;
0060   MonitorElement* meBarrelPIDKAsKEff_;
0061   MonitorElement* meBarrelPIDKAsPEff_;
0062   MonitorElement* meBarrelPIDKNoPIDEff_;
0063 
0064   MonitorElement* meBarrelPIDPAsPiEff_;
0065   MonitorElement* meBarrelPIDPAsKEff_;
0066   MonitorElement* meBarrelPIDPAsPEff_;
0067   MonitorElement* meBarrelPIDPNoPIDEff_;
0068 
0069   MonitorElement* meEndcapPIDPiAsPiEff_;
0070   MonitorElement* meEndcapPIDPiAsKEff_;
0071   MonitorElement* meEndcapPIDPiAsPEff_;
0072   MonitorElement* meEndcapPIDPiNoPIDEff_;
0073 
0074   MonitorElement* meEndcapPIDKAsPiEff_;
0075   MonitorElement* meEndcapPIDKAsKEff_;
0076   MonitorElement* meEndcapPIDKAsPEff_;
0077   MonitorElement* meEndcapPIDKNoPIDEff_;
0078 
0079   MonitorElement* meEndcapPIDPAsPiEff_;
0080   MonitorElement* meEndcapPIDPAsKEff_;
0081   MonitorElement* meEndcapPIDPAsPEff_;
0082   MonitorElement* meEndcapPIDPNoPIDEff_;
0083 
0084   MonitorElement* meBarrelPiPurity_;
0085   MonitorElement* meBarrelKPurity_;
0086   MonitorElement* meBarrelPPurity_;
0087 
0088   MonitorElement* meEndcapPiPurity_;
0089   MonitorElement* meEndcapKPurity_;
0090   MonitorElement* meEndcapPPurity_;
0091 };
0092 
0093 // ------------ constructor and destructor --------------
0094 Primary4DVertexHarvester::Primary4DVertexHarvester(const edm::ParameterSet& iConfig)
0095     : folder_(iConfig.getParameter<std::string>("folder")) {}
0096 
0097 Primary4DVertexHarvester::~Primary4DVertexHarvester() {}
0098 
0099 // auxiliary method to compute efficiency from the ratio of two 1D MonitorElement
0100 void Primary4DVertexHarvester::computeEfficiency1D(MonitorElement* num, MonitorElement* den, MonitorElement* result) {
0101   for (int ibin = 1; ibin <= den->getNbinsX(); ibin++) {
0102     double eff = num->getBinContent(ibin) / den->getBinContent(ibin);
0103     double bin_err = sqrt((num->getBinContent(ibin) * (den->getBinContent(ibin) - num->getBinContent(ibin))) /
0104                           pow(den->getBinContent(ibin), 3));
0105     if (den->getBinContent(ibin) == 0) {
0106       eff = 0;
0107       bin_err = 0;
0108     }
0109     result->setBinContent(ibin, eff);
0110     result->setBinError(ibin, bin_err);
0111   }
0112 }
0113 
0114 // auxiliary method to add 1D MonitorElement toBeAdded to a base ME
0115 void Primary4DVertexHarvester::incrementME(MonitorElement* base, MonitorElement* toBeAdded) {
0116   for (int ibin = 1; ibin <= base->getNbinsX(); ibin++) {
0117     double newC = base->getBinContent(ibin) + toBeAdded->getBinContent(ibin);
0118     double newE = std::sqrt(newC);
0119     base->setBinContent(ibin, newC);
0120     base->setBinError(ibin, newE);
0121   }
0122 }
0123 
0124 // ------------ endjob tasks ----------------------------
0125 void Primary4DVertexHarvester::dqmEndJob(DQMStore::IBooker& ibook, DQMStore::IGetter& igetter) {
0126   // --- Get the monitoring histograms
0127   MonitorElement* meTrackEffPtTot = igetter.get(folder_ + "EffPtTot");
0128   MonitorElement* meTrackMatchedTPEffPtTot = igetter.get(folder_ + "MatchedTPEffPtTot");
0129   MonitorElement* meTrackMatchedTPEffPtMtd = igetter.get(folder_ + "MatchedTPEffPtMtd");
0130   MonitorElement* meTrackEffEtaTot = igetter.get(folder_ + "EffEtaTot");
0131   MonitorElement* meTrackMatchedTPEffEtaTot = igetter.get(folder_ + "MatchedTPEffEtaTot");
0132   MonitorElement* meTrackMatchedTPEffEtaMtd = igetter.get(folder_ + "MatchedTPEffEtaMtd");
0133   MonitorElement* meRecoVtxVsLineDensity = igetter.get(folder_ + "RecoVtxVsLineDensity");
0134   MonitorElement* meRecVerNumber = igetter.get(folder_ + "RecVerNumber");
0135 
0136   if (!meTrackEffPtTot || !meTrackMatchedTPEffPtTot || !meTrackMatchedTPEffPtMtd || !meTrackEffEtaTot ||
0137       !meTrackMatchedTPEffEtaTot || !meTrackMatchedTPEffEtaMtd || !meRecoVtxVsLineDensity || !meRecVerNumber) {
0138     edm::LogError("Primary4DVertexHarvester") << "Monitoring histograms not found!" << std::endl;
0139     return;
0140   }
0141 
0142   // Normalize line density plot
0143   double nEvt = meRecVerNumber->getEntries();
0144   if (nEvt > 0.) {
0145     nEvt = 1. / nEvt;
0146     double nEntries = meRecoVtxVsLineDensity->getEntries();
0147     for (int ibin = 1; ibin <= meRecoVtxVsLineDensity->getNbinsX(); ibin++) {
0148       double cont = meRecoVtxVsLineDensity->getBinContent(ibin) * nEvt;
0149       double bin_err = meRecoVtxVsLineDensity->getBinError(ibin) * nEvt;
0150       meRecoVtxVsLineDensity->setBinContent(ibin, cont);
0151       meRecoVtxVsLineDensity->setBinError(ibin, bin_err);
0152     }
0153     meRecoVtxVsLineDensity->setEntries(nEntries);
0154   }
0155 
0156   // --- Book  histograms
0157   ibook.cd(folder_);
0158   meTPPtSelEff_ = ibook.book1D("TPPtSelEff",
0159                                "Track associated to LV selected efficiency TP VS Pt;Pt [GeV];Efficiency",
0160                                meTrackEffPtTot->getNbinsX(),
0161                                meTrackEffPtTot->getTH1()->GetXaxis()->GetXmin(),
0162                                meTrackEffPtTot->getTH1()->GetXaxis()->GetXmax());
0163   meTPPtSelEff_->getTH1()->SetMinimum(0.);
0164   computeEfficiency1D(meTrackMatchedTPEffPtTot, meTrackEffPtTot, meTPPtSelEff_);
0165 
0166   meTPEtaSelEff_ = ibook.book1D("TPEtaSelEff",
0167                                 "Track associated to LV selected efficiency TP VS Eta;Eta;Efficiency",
0168                                 meTrackEffEtaTot->getNbinsX(),
0169                                 meTrackEffEtaTot->getTH1()->GetXaxis()->GetXmin(),
0170                                 meTrackEffEtaTot->getTH1()->GetXaxis()->GetXmax());
0171   meTPEtaSelEff_->getTH1()->SetMinimum(0.);
0172   computeEfficiency1D(meTrackMatchedTPEffEtaTot, meTrackEffEtaTot, meTPEtaSelEff_);
0173 
0174   meTPPtMatchEff_ = ibook.book1D("TPPtMatchEff",
0175                                  "Track associated to LV matched to TP efficiency VS Pt;Pt [GeV];Efficiency",
0176                                  meTrackMatchedTPEffPtTot->getNbinsX(),
0177                                  meTrackMatchedTPEffPtTot->getTH1()->GetXaxis()->GetXmin(),
0178                                  meTrackMatchedTPEffPtTot->getTH1()->GetXaxis()->GetXmax());
0179   meTPPtMatchEff_->getTH1()->SetMinimum(0.);
0180   computeEfficiency1D(meTrackMatchedTPEffPtMtd, meTrackMatchedTPEffPtTot, meTPPtMatchEff_);
0181 
0182   meTPEtaMatchEff_ = ibook.book1D("TPEtaMatchEff",
0183                                   "Track associated to LV matched to TP efficiency VS Eta;Eta;Efficiency",
0184                                   meTrackMatchedTPEffEtaTot->getNbinsX(),
0185                                   meTrackMatchedTPEffEtaTot->getTH1()->GetXaxis()->GetXmin(),
0186                                   meTrackMatchedTPEffEtaTot->getTH1()->GetXaxis()->GetXmax());
0187   meTPEtaMatchEff_->getTH1()->SetMinimum(0.);
0188   computeEfficiency1D(meTrackMatchedTPEffEtaMtd, meTrackMatchedTPEffEtaTot, meTPEtaMatchEff_);
0189 
0190   MonitorElement* meBarrelPIDp = igetter.get(folder_ + "BarrelPIDp");
0191   MonitorElement* meEndcapPIDp = igetter.get(folder_ + "EndcapPIDp");
0192 
0193   MonitorElement* meBarrelTruePiNoPID = igetter.get(folder_ + "BarrelTruePiNoPID");
0194   MonitorElement* meBarrelTrueKNoPID = igetter.get(folder_ + "BarrelTrueKNoPID");
0195   MonitorElement* meBarrelTruePNoPID = igetter.get(folder_ + "BarrelTruePNoPID");
0196   MonitorElement* meEndcapTruePiNoPID = igetter.get(folder_ + "EndcapTruePiNoPID");
0197   MonitorElement* meEndcapTrueKNoPID = igetter.get(folder_ + "EndcapTrueKNoPID");
0198   MonitorElement* meEndcapTruePNoPID = igetter.get(folder_ + "EndcapTruePNoPID");
0199 
0200   MonitorElement* meBarrelTruePiAsPi = igetter.get(folder_ + "BarrelTruePiAsPi");
0201   MonitorElement* meBarrelTrueKAsPi = igetter.get(folder_ + "BarrelTrueKAsPi");
0202   MonitorElement* meBarrelTruePAsPi = igetter.get(folder_ + "BarrelTruePAsPi");
0203   MonitorElement* meEndcapTruePiAsPi = igetter.get(folder_ + "EndcapTruePiAsPi");
0204   MonitorElement* meEndcapTrueKAsPi = igetter.get(folder_ + "EndcapTrueKAsPi");
0205   MonitorElement* meEndcapTruePAsPi = igetter.get(folder_ + "EndcapTruePAsPi");
0206 
0207   MonitorElement* meBarrelTruePiAsK = igetter.get(folder_ + "BarrelTruePiAsK");
0208   MonitorElement* meBarrelTrueKAsK = igetter.get(folder_ + "BarrelTrueKAsK");
0209   MonitorElement* meBarrelTruePAsK = igetter.get(folder_ + "BarrelTruePAsK");
0210   MonitorElement* meEndcapTruePiAsK = igetter.get(folder_ + "EndcapTruePiAsK");
0211   MonitorElement* meEndcapTrueKAsK = igetter.get(folder_ + "EndcapTrueKAsK");
0212   MonitorElement* meEndcapTruePAsK = igetter.get(folder_ + "EndcapTruePAsK");
0213 
0214   MonitorElement* meBarrelTruePiAsP = igetter.get(folder_ + "BarrelTruePiAsP");
0215   MonitorElement* meBarrelTrueKAsP = igetter.get(folder_ + "BarrelTrueKAsP");
0216   MonitorElement* meBarrelTruePAsP = igetter.get(folder_ + "BarrelTruePAsP");
0217   MonitorElement* meEndcapTruePiAsP = igetter.get(folder_ + "EndcapTruePiAsP");
0218   MonitorElement* meEndcapTrueKAsP = igetter.get(folder_ + "EndcapTrueKAsP");
0219   MonitorElement* meEndcapTruePAsP = igetter.get(folder_ + "EndcapTruePAsP");
0220 
0221   if (!meBarrelPIDp || !meEndcapPIDp || !meBarrelTruePiNoPID || !meBarrelTrueKNoPID || !meBarrelTruePNoPID ||
0222       !meEndcapTruePiNoPID || !meEndcapTrueKNoPID || !meEndcapTruePNoPID || !meBarrelTruePiAsPi || !meBarrelTrueKAsPi ||
0223       !meBarrelTruePAsPi || !meEndcapTruePiAsPi || !meEndcapTrueKAsPi || !meEndcapTruePAsPi || !meBarrelTruePiAsK ||
0224       !meBarrelTrueKAsK || !meBarrelTruePAsK || !meEndcapTruePiAsK || !meEndcapTrueKAsK || !meEndcapTruePAsK ||
0225       !meBarrelTruePiAsP || !meBarrelTrueKAsP || !meBarrelTruePAsP || !meEndcapTruePiAsP || !meEndcapTrueKAsP ||
0226       !meEndcapTruePAsP) {
0227     edm::LogWarning("Primary4DVertexHarvester") << "PID Monitoring histograms not found!" << std::endl;
0228     return;
0229   }
0230 
0231   meBarrelTruePi_ = ibook.book1D("BarrelTruePi",
0232                                  "Barrel True Pi P;P [GeV]",
0233                                  meBarrelPIDp->getNbinsX(),
0234                                  meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0235                                  meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0236   incrementME(meBarrelTruePi_, meBarrelTruePiAsPi);
0237   incrementME(meBarrelTruePi_, meBarrelTruePiAsK);
0238   incrementME(meBarrelTruePi_, meBarrelTruePiAsP);
0239   incrementME(meBarrelTruePi_, meBarrelTruePiNoPID);
0240 
0241   meEndcapTruePi_ = ibook.book1D("EndcapTruePi",
0242                                  "Endcap True Pi P;P [GeV]",
0243                                  meBarrelPIDp->getNbinsX(),
0244                                  meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0245                                  meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0246   incrementME(meEndcapTruePi_, meEndcapTruePiAsPi);
0247   incrementME(meEndcapTruePi_, meEndcapTruePiAsK);
0248   incrementME(meEndcapTruePi_, meEndcapTruePiAsP);
0249   incrementME(meEndcapTruePi_, meEndcapTruePiNoPID);
0250 
0251   meBarrelTrueK_ = ibook.book1D("BarrelTrueK",
0252                                 "Barrel True K P;P [GeV]",
0253                                 meBarrelPIDp->getNbinsX(),
0254                                 meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0255                                 meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0256   incrementME(meBarrelTrueK_, meBarrelTrueKAsPi);
0257   incrementME(meBarrelTrueK_, meBarrelTrueKAsK);
0258   incrementME(meBarrelTrueK_, meBarrelTrueKAsP);
0259   incrementME(meBarrelTrueK_, meBarrelTrueKNoPID);
0260 
0261   meEndcapTrueK_ = ibook.book1D("EndcapTrueK",
0262                                 "Endcap True K P;P [GeV]",
0263                                 meBarrelPIDp->getNbinsX(),
0264                                 meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0265                                 meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0266   incrementME(meEndcapTrueK_, meEndcapTrueKAsPi);
0267   incrementME(meEndcapTrueK_, meEndcapTrueKAsK);
0268   incrementME(meEndcapTrueK_, meEndcapTrueKAsP);
0269   incrementME(meEndcapTrueK_, meEndcapTrueKNoPID);
0270 
0271   meBarrelTrueP_ = ibook.book1D("BarrelTrueP",
0272                                 "Barrel True P P;P [GeV]",
0273                                 meBarrelPIDp->getNbinsX(),
0274                                 meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0275                                 meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0276   incrementME(meBarrelTrueP_, meBarrelTruePAsPi);
0277   incrementME(meBarrelTrueP_, meBarrelTruePAsK);
0278   incrementME(meBarrelTrueP_, meBarrelTruePAsP);
0279   incrementME(meBarrelTrueP_, meBarrelTruePNoPID);
0280 
0281   meEndcapTrueP_ = ibook.book1D("EndcapTrueP",
0282                                 "Endcap True P P;P [GeV]",
0283                                 meBarrelPIDp->getNbinsX(),
0284                                 meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0285                                 meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0286   incrementME(meEndcapTrueP_, meEndcapTruePAsPi);
0287   incrementME(meEndcapTrueP_, meEndcapTruePAsK);
0288   incrementME(meEndcapTrueP_, meEndcapTruePAsP);
0289   incrementME(meEndcapTrueP_, meEndcapTruePNoPID);
0290 
0291   meBarrelPIDPiAsPiEff_ = ibook.book1D("BarrelPIDPiAsPiEff",
0292                                        "Barrel True pi as pi id. fraction VS P;P [GeV]",
0293                                        meBarrelPIDp->getNbinsX(),
0294                                        meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0295                                        meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0296   meBarrelPIDPiAsPiEff_->getTH1()->SetMinimum(0.);
0297   computeEfficiency1D(meBarrelTruePiAsPi, meBarrelTruePi_, meBarrelPIDPiAsPiEff_);
0298 
0299   meBarrelPIDPiAsKEff_ = ibook.book1D("BarrelPIDPiAsKEff",
0300                                       "Barrel True pi as k id. fraction VS P;P [GeV]",
0301                                       meBarrelPIDp->getNbinsX(),
0302                                       meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0303                                       meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0304   meBarrelPIDPiAsKEff_->getTH1()->SetMinimum(0.);
0305   computeEfficiency1D(meBarrelTruePiAsK, meBarrelTruePi_, meBarrelPIDPiAsKEff_);
0306 
0307   meBarrelPIDPiAsPEff_ = ibook.book1D("BarrelPIDPiAsPEff",
0308                                       "Barrel True pi as p id. fraction VS P;P [GeV]",
0309                                       meBarrelPIDp->getNbinsX(),
0310                                       meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0311                                       meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0312   meBarrelPIDPiAsPEff_->getTH1()->SetMinimum(0.);
0313   computeEfficiency1D(meBarrelTruePiAsP, meBarrelTruePi_, meBarrelPIDPiAsPEff_);
0314 
0315   meBarrelPIDPiNoPIDEff_ = ibook.book1D("BarrelPIDPiNoPIDEff",
0316                                         "Barrel True pi no PID id. fraction VS P;P [GeV]",
0317                                         meBarrelPIDp->getNbinsX(),
0318                                         meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0319                                         meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0320   meBarrelPIDPiNoPIDEff_->getTH1()->SetMinimum(0.);
0321   computeEfficiency1D(meBarrelTruePiNoPID, meBarrelTruePi_, meBarrelPIDPiNoPIDEff_);
0322 
0323   meBarrelPIDKAsPiEff_ = ibook.book1D("BarrelPIDKAsPiEff",
0324                                       "Barrel True k as pi id. fraction VS P;P [GeV]",
0325                                       meBarrelPIDp->getNbinsX(),
0326                                       meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0327                                       meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0328   meBarrelPIDKAsPiEff_->getTH1()->SetMinimum(0.);
0329   computeEfficiency1D(meBarrelTrueKAsPi, meBarrelTrueK_, meBarrelPIDKAsPiEff_);
0330 
0331   meBarrelPIDKAsKEff_ = ibook.book1D("BarrelPIDKAsKEff",
0332                                      "Barrel True k as k id. fraction VS P;P [GeV]",
0333                                      meBarrelPIDp->getNbinsX(),
0334                                      meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0335                                      meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0336   meBarrelPIDKAsKEff_->getTH1()->SetMinimum(0.);
0337   computeEfficiency1D(meBarrelTrueKAsK, meBarrelTrueK_, meBarrelPIDKAsKEff_);
0338 
0339   meBarrelPIDKAsPEff_ = ibook.book1D("BarrelPIDKAsPEff",
0340                                      "Barrel True k as p id. fraction VS P;P [GeV]",
0341                                      meBarrelPIDp->getNbinsX(),
0342                                      meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0343                                      meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0344   meBarrelPIDKAsPEff_->getTH1()->SetMinimum(0.);
0345   computeEfficiency1D(meBarrelTrueKAsP, meBarrelTrueK_, meBarrelPIDKAsPEff_);
0346 
0347   meBarrelPIDKNoPIDEff_ = ibook.book1D("BarrelPIDKNoPIDEff",
0348                                        "Barrel True k no PID id. fraction VS P;P [GeV]",
0349                                        meBarrelPIDp->getNbinsX(),
0350                                        meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0351                                        meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0352   meBarrelPIDKNoPIDEff_->getTH1()->SetMinimum(0.);
0353   computeEfficiency1D(meBarrelTrueKNoPID, meBarrelTrueK_, meBarrelPIDKNoPIDEff_);
0354 
0355   meBarrelPIDPAsPiEff_ = ibook.book1D("BarrelPIDPAsPiEff",
0356                                       "Barrel True p as pi id. fraction VS P;P [GeV]",
0357                                       meBarrelPIDp->getNbinsX(),
0358                                       meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0359                                       meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0360   meBarrelPIDPAsPiEff_->getTH1()->SetMinimum(0.);
0361   computeEfficiency1D(meBarrelTruePAsPi, meBarrelTrueP_, meBarrelPIDPAsPiEff_);
0362 
0363   meBarrelPIDPAsKEff_ = ibook.book1D("BarrelPIDPAsKEff",
0364                                      "Barrel True p as k id. fraction VS P;P [GeV]",
0365                                      meBarrelPIDp->getNbinsX(),
0366                                      meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0367                                      meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0368   meBarrelPIDPAsKEff_->getTH1()->SetMinimum(0.);
0369   computeEfficiency1D(meBarrelTruePAsK, meBarrelTrueP_, meBarrelPIDPAsKEff_);
0370 
0371   meBarrelPIDPAsPEff_ = ibook.book1D("BarrelPIDPAsPEff",
0372                                      "Barrel True p as p id. fraction VS P;P [GeV]",
0373                                      meBarrelPIDp->getNbinsX(),
0374                                      meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0375                                      meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0376   meBarrelPIDPAsPEff_->getTH1()->SetMinimum(0.);
0377   computeEfficiency1D(meBarrelTruePAsP, meBarrelTrueP_, meBarrelPIDPAsPEff_);
0378 
0379   meBarrelPIDPNoPIDEff_ = ibook.book1D("BarrelPIDPNoPIDEff",
0380                                        "Barrel True p no PID id. fraction VS P;P [GeV]",
0381                                        meBarrelPIDp->getNbinsX(),
0382                                        meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0383                                        meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0384   meBarrelPIDPNoPIDEff_->getTH1()->SetMinimum(0.);
0385   computeEfficiency1D(meBarrelTruePNoPID, meBarrelTrueP_, meBarrelPIDPNoPIDEff_);
0386 
0387   meEndcapPIDPiAsPiEff_ = ibook.book1D("EndcapPIDPiAsPiEff",
0388                                        "Endcap True pi as pi id. fraction VS P;P [GeV]",
0389                                        meBarrelPIDp->getNbinsX(),
0390                                        meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0391                                        meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0392   meEndcapPIDPiAsPiEff_->getTH1()->SetMinimum(0.);
0393   computeEfficiency1D(meEndcapTruePiAsPi, meEndcapTruePi_, meEndcapPIDPiAsPiEff_);
0394 
0395   meEndcapPIDPiAsKEff_ = ibook.book1D("EndcapPIDPiAsKEff",
0396                                       "Endcap True pi as k id. fraction VS P;P [GeV]",
0397                                       meBarrelPIDp->getNbinsX(),
0398                                       meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0399                                       meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0400   meEndcapPIDPiAsKEff_->getTH1()->SetMinimum(0.);
0401   computeEfficiency1D(meEndcapTruePiAsK, meEndcapTruePi_, meEndcapPIDPiAsKEff_);
0402 
0403   meEndcapPIDPiAsPEff_ = ibook.book1D("EndcapPIDPiAsPEff",
0404                                       "Endcap True pi as p id. fraction VS P;P [GeV]",
0405                                       meBarrelPIDp->getNbinsX(),
0406                                       meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0407                                       meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0408   meEndcapPIDPiAsPEff_->getTH1()->SetMinimum(0.);
0409   computeEfficiency1D(meEndcapTruePiAsP, meEndcapTruePi_, meEndcapPIDPiAsPEff_);
0410 
0411   meEndcapPIDPiNoPIDEff_ = ibook.book1D("EndcapPIDPiNoPIDEff",
0412                                         "Endcap True pi no PID id. fraction VS P;P [GeV]",
0413                                         meBarrelPIDp->getNbinsX(),
0414                                         meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0415                                         meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0416   meEndcapPIDPiNoPIDEff_->getTH1()->SetMinimum(0.);
0417   computeEfficiency1D(meEndcapTruePiNoPID, meEndcapTruePi_, meEndcapPIDPiNoPIDEff_);
0418 
0419   meEndcapPIDKAsPiEff_ = ibook.book1D("EndcapPIDKAsPiEff",
0420                                       "Endcap True k as pi id. fraction VS P;P [GeV]",
0421                                       meBarrelPIDp->getNbinsX(),
0422                                       meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0423                                       meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0424   meEndcapPIDKAsPiEff_->getTH1()->SetMinimum(0.);
0425   computeEfficiency1D(meEndcapTrueKAsPi, meEndcapTrueK_, meEndcapPIDKAsPiEff_);
0426 
0427   meEndcapPIDKAsKEff_ = ibook.book1D("EndcapPIDKAsKEff",
0428                                      "Endcap True k as k id. fraction VS P;P [GeV]",
0429                                      meBarrelPIDp->getNbinsX(),
0430                                      meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0431                                      meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0432   meEndcapPIDKAsKEff_->getTH1()->SetMinimum(0.);
0433   computeEfficiency1D(meEndcapTrueKAsK, meEndcapTrueK_, meEndcapPIDKAsKEff_);
0434 
0435   meEndcapPIDKAsPEff_ = ibook.book1D("EndcapPIDKAsPEff",
0436                                      "Endcap True k as p id. fraction VS P;P [GeV]",
0437                                      meBarrelPIDp->getNbinsX(),
0438                                      meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0439                                      meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0440   meEndcapPIDKAsPEff_->getTH1()->SetMinimum(0.);
0441   computeEfficiency1D(meEndcapTrueKAsP, meEndcapTrueK_, meEndcapPIDKAsPEff_);
0442 
0443   meEndcapPIDKNoPIDEff_ = ibook.book1D("EndcapPIDKNoPIDEff",
0444                                        "Endcap True k no PID id. fraction VS P;P [GeV]",
0445                                        meBarrelPIDp->getNbinsX(),
0446                                        meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0447                                        meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0448   meEndcapPIDKNoPIDEff_->getTH1()->SetMinimum(0.);
0449   computeEfficiency1D(meEndcapTrueKNoPID, meEndcapTrueK_, meEndcapPIDKNoPIDEff_);
0450 
0451   meEndcapPIDPAsPiEff_ = ibook.book1D("EndcapPIDPAsPiEff",
0452                                       "Endcap True p as pi id. fraction VS P;P [GeV]",
0453                                       meBarrelPIDp->getNbinsX(),
0454                                       meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0455                                       meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0456   meEndcapPIDPAsPiEff_->getTH1()->SetMinimum(0.);
0457   computeEfficiency1D(meEndcapTruePAsPi, meEndcapTrueP_, meEndcapPIDPAsPiEff_);
0458 
0459   meEndcapPIDPAsKEff_ = ibook.book1D("EndcapPIDPAsKEff",
0460                                      "Endcap True p as k id. fraction VS P;P [GeV]",
0461                                      meBarrelPIDp->getNbinsX(),
0462                                      meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0463                                      meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0464   meEndcapPIDPAsKEff_->getTH1()->SetMinimum(0.);
0465   computeEfficiency1D(meEndcapTruePAsK, meEndcapTrueP_, meEndcapPIDPAsKEff_);
0466 
0467   meEndcapPIDPAsPEff_ = ibook.book1D("EndcapPIDPAsPEff",
0468                                      "Endcap True p as p id. fraction VS P;P [GeV]",
0469                                      meBarrelPIDp->getNbinsX(),
0470                                      meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0471                                      meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0472   meEndcapPIDPAsPEff_->getTH1()->SetMinimum(0.);
0473   computeEfficiency1D(meEndcapTruePAsP, meEndcapTrueP_, meEndcapPIDPAsPEff_);
0474 
0475   meEndcapPIDPNoPIDEff_ = ibook.book1D("EndcapPIDPNoPIDEff",
0476                                        "Endcap True p no PID id. fraction VS P;P [GeV]",
0477                                        meBarrelPIDp->getNbinsX(),
0478                                        meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0479                                        meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0480   meEndcapPIDPNoPIDEff_->getTH1()->SetMinimum(0.);
0481   computeEfficiency1D(meEndcapTruePNoPID, meEndcapTrueP_, meEndcapPIDPNoPIDEff_);
0482 
0483   meBarrelAsPi_ = ibook.book1D("BarrelAsPi",
0484                                "Barrel Identified Pi P;P [GeV]",
0485                                meBarrelPIDp->getNbinsX(),
0486                                meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0487                                meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0488   incrementME(meBarrelAsPi_, meBarrelTruePiAsPi);
0489   incrementME(meBarrelAsPi_, meBarrelTrueKAsPi);
0490   incrementME(meBarrelAsPi_, meBarrelTruePAsPi);
0491 
0492   meEndcapAsPi_ = ibook.book1D("EndcapAsPi",
0493                                "Endcap Identified Pi P;P [GeV]",
0494                                meBarrelPIDp->getNbinsX(),
0495                                meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0496                                meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0497   incrementME(meEndcapAsPi_, meEndcapTruePiAsPi);
0498   incrementME(meEndcapAsPi_, meEndcapTrueKAsPi);
0499   incrementME(meEndcapAsPi_, meEndcapTruePAsPi);
0500 
0501   meBarrelAsK_ = ibook.book1D("BarrelAsK",
0502                               "Barrel Identified K P;P [GeV]",
0503                               meBarrelPIDp->getNbinsX(),
0504                               meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0505                               meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0506   incrementME(meBarrelAsK_, meBarrelTruePiAsK);
0507   incrementME(meBarrelAsK_, meBarrelTrueKAsK);
0508   incrementME(meBarrelAsK_, meBarrelTruePAsK);
0509 
0510   meEndcapAsK_ = ibook.book1D("EndcapAsK",
0511                               "Endcap Identified K P;P [GeV]",
0512                               meBarrelPIDp->getNbinsX(),
0513                               meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0514                               meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0515   incrementME(meEndcapAsK_, meEndcapTruePiAsK);
0516   incrementME(meEndcapAsK_, meEndcapTrueKAsK);
0517   incrementME(meEndcapAsK_, meEndcapTruePAsK);
0518 
0519   meBarrelAsP_ = ibook.book1D("BarrelAsP",
0520                               "Barrel Identified P P;P [GeV]",
0521                               meBarrelPIDp->getNbinsX(),
0522                               meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0523                               meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0524   incrementME(meBarrelAsP_, meBarrelTruePiAsP);
0525   incrementME(meBarrelAsP_, meBarrelTrueKAsP);
0526   incrementME(meBarrelAsP_, meBarrelTruePAsP);
0527 
0528   meEndcapAsP_ = ibook.book1D("EndcapAsP",
0529                               "Endcap Identified P P;P [GeV]",
0530                               meBarrelPIDp->getNbinsX(),
0531                               meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0532                               meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0533   incrementME(meEndcapAsP_, meEndcapTruePiAsP);
0534   incrementME(meEndcapAsP_, meEndcapTrueKAsP);
0535   incrementME(meEndcapAsP_, meEndcapTruePAsP);
0536 
0537   meBarrelNoPID_ = ibook.book1D("BarrelNoPID",
0538                                 "Barrel NoPID P;P [GeV]",
0539                                 meBarrelPIDp->getNbinsX(),
0540                                 meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0541                                 meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0542   incrementME(meBarrelNoPID_, meBarrelTruePiNoPID);
0543   incrementME(meBarrelNoPID_, meBarrelTrueKNoPID);
0544   incrementME(meBarrelNoPID_, meBarrelTruePNoPID);
0545 
0546   meEndcapNoPID_ = ibook.book1D("EndcapNoPID",
0547                                 "Endcap NoPID P;P [GeV]",
0548                                 meBarrelPIDp->getNbinsX(),
0549                                 meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0550                                 meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0551   incrementME(meEndcapNoPID_, meEndcapTruePiNoPID);
0552   incrementME(meEndcapNoPID_, meEndcapTrueKNoPID);
0553   incrementME(meEndcapNoPID_, meEndcapTruePNoPID);
0554 
0555   meBarrelPiPurity_ = ibook.book1D("BarrelPiPurity",
0556                                    "Barrel pi id. fraction true pi VS P;P [GeV]",
0557                                    meBarrelPIDp->getNbinsX(),
0558                                    meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0559                                    meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0560   meBarrelPiPurity_->getTH1()->SetMinimum(0.);
0561   computeEfficiency1D(meBarrelTruePiAsPi, meBarrelAsPi_, meBarrelPiPurity_);
0562 
0563   meBarrelKPurity_ = ibook.book1D("BarrelKPurity",
0564                                   "Barrel k id. fraction true k VS P;P [GeV]",
0565                                   meBarrelPIDp->getNbinsX(),
0566                                   meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0567                                   meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0568   meBarrelKPurity_->getTH1()->SetMinimum(0.);
0569   computeEfficiency1D(meBarrelTrueKAsK, meBarrelAsK_, meBarrelKPurity_);
0570 
0571   meBarrelPPurity_ = ibook.book1D("BarrelPPurity",
0572                                   "Barrel p id. fraction true p VS P;P [GeV]",
0573                                   meBarrelPIDp->getNbinsX(),
0574                                   meBarrelPIDp->getTH1()->GetXaxis()->GetXmin(),
0575                                   meBarrelPIDp->getTH1()->GetXaxis()->GetXmax());
0576   meBarrelPPurity_->getTH1()->SetMinimum(0.);
0577   computeEfficiency1D(meBarrelTruePAsP, meBarrelAsP_, meBarrelPPurity_);
0578 
0579   meEndcapPiPurity_ = ibook.book1D("EndcapPiPurity",
0580                                    "Endcap pi id. fraction true pi VS P;P [GeV]",
0581                                    meEndcapPIDp->getNbinsX(),
0582                                    meEndcapPIDp->getTH1()->GetXaxis()->GetXmin(),
0583                                    meEndcapPIDp->getTH1()->GetXaxis()->GetXmax());
0584   meEndcapPiPurity_->getTH1()->SetMinimum(0.);
0585   computeEfficiency1D(meEndcapTruePiAsPi, meEndcapAsPi_, meEndcapPiPurity_);
0586 
0587   meEndcapKPurity_ = ibook.book1D("EndcapKPurity",
0588                                   "Endcap k id. fraction true k VS P;P [GeV]",
0589                                   meEndcapPIDp->getNbinsX(),
0590                                   meEndcapPIDp->getTH1()->GetXaxis()->GetXmin(),
0591                                   meEndcapPIDp->getTH1()->GetXaxis()->GetXmax());
0592   meEndcapKPurity_->getTH1()->SetMinimum(0.);
0593   computeEfficiency1D(meEndcapTrueKAsK, meEndcapAsK_, meEndcapKPurity_);
0594 
0595   meEndcapPPurity_ = ibook.book1D("EndcapPPurity",
0596                                   "Endcap p id. fraction true p VS P;P [GeV]",
0597                                   meEndcapPIDp->getNbinsX(),
0598                                   meEndcapPIDp->getTH1()->GetXaxis()->GetXmin(),
0599                                   meEndcapPIDp->getTH1()->GetXaxis()->GetXmax());
0600   meEndcapPPurity_->getTH1()->SetMinimum(0.);
0601   computeEfficiency1D(meEndcapTruePAsP, meEndcapAsP_, meEndcapPPurity_);
0602 }
0603 
0604 // ------------ method fills 'descriptions' with the allowed parameters for the module  ----------
0605 void Primary4DVertexHarvester::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0606   edm::ParameterSetDescription desc;
0607 
0608   desc.add<std::string>("folder", "MTD/Vertices/");
0609 
0610   descriptions.add("Primary4DVertexPostProcessor", desc);
0611 }
0612 
0613 DEFINE_FWK_MODULE(Primary4DVertexHarvester);