Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-23 22:48:18

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