Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:46

0001 #ifndef DQMOffline_RecoB_IPTagPlotter_cc_h
0002 #define DQMOffline_RecoB_IPTagPlotter_cc_h
0003 
0004 #include <cstddef>
0005 #include <string>
0006 
0007 #include "DataFormats/VertexReco/interface/Vertex.h"
0008 #include "DQMOffline/RecoB/interface/IPTagPlotter.h"
0009 
0010 template <class Container, class Base>
0011 IPTagPlotter<Container, Base>::IPTagPlotter(const std::string& tagName,
0012                                             const EtaPtBin& etaPtBin,
0013                                             const edm::ParameterSet& pSet,
0014                                             unsigned int mc,
0015                                             bool wf,
0016                                             DQMStore::IBooker& ibook_)
0017     : BaseTagInfoPlotter(tagName, etaPtBin),
0018       nBinEffPur_(pSet.getParameter<int>("nBinEffPur")),
0019       startEffPur_(pSet.getParameter<double>("startEffPur")),
0020       endEffPur_(pSet.getParameter<double>("endEffPur")),
0021       mcPlots_(mc),
0022       willFinalize_(wf),
0023       makeQualityPlots_(pSet.getParameter<bool>("QualityPlots")),
0024       lowerIPSBound(pSet.getParameter<double>("LowerIPSBound")),
0025       upperIPSBound(pSet.getParameter<double>("UpperIPSBound")),
0026       lowerIPBound(pSet.getParameter<double>("LowerIPBound")),
0027       upperIPBound(pSet.getParameter<double>("UpperIPBound")),
0028       lowerIPEBound(pSet.getParameter<double>("LowerIPEBound")),
0029       upperIPEBound(pSet.getParameter<double>("UpperIPEBound")),
0030       nBinsIPS(pSet.getParameter<int>("NBinsIPS")),
0031       nBinsIP(pSet.getParameter<int>("NBinsIP")),
0032       nBinsIPE(pSet.getParameter<int>("NBinsIPE")),
0033       minDecayLength(pSet.getParameter<double>("MinDecayLength")),
0034       maxDecayLength(pSet.getParameter<double>("MaxDecayLength")),
0035       minJetDistance(pSet.getParameter<double>("MinJetDistance")),
0036       maxJetDistance(pSet.getParameter<double>("MaxJetDistance")) {
0037   const std::string trackIPDir(theExtensionString.substr(1));
0038 
0039   if (willFinalize_)
0040     return;
0041 
0042   // Number of tracks
0043   // 3D
0044   trkNbr3D = std::make_unique<TrackIPHistograms<int>>("selTrksNbr_3D" + theExtensionString,
0045                                                       "Number of selected tracks for 3D IPS",
0046                                                       31,
0047                                                       -0.5,
0048                                                       30.5,
0049                                                       false,
0050                                                       true,
0051                                                       true,
0052                                                       "b",
0053                                                       trackIPDir,
0054                                                       mc,
0055                                                       makeQualityPlots_,
0056                                                       ibook_);
0057 
0058   // 2D
0059   trkNbr2D = std::make_unique<TrackIPHistograms<int>>("selTrksNbr_2D" + theExtensionString,
0060                                                       "Number of selected tracks for 2D IPS",
0061                                                       31,
0062                                                       -0.5,
0063                                                       30.5,
0064                                                       false,
0065                                                       true,
0066                                                       true,
0067                                                       "b",
0068                                                       trackIPDir,
0069                                                       mc,
0070                                                       makeQualityPlots_,
0071                                                       ibook_);
0072 
0073   // IP significance
0074   // 3D
0075   for (unsigned int i = 1; i <= 4; i++) {
0076     tkcntHistosSig3D.push_back(
0077         std::make_unique<TrackIPHistograms<double>>("ips" + std::to_string(i) + "_3D" + theExtensionString,
0078                                                     "3D IP significance " + std::to_string(i) + ".trk",
0079                                                     nBinsIPS,
0080                                                     lowerIPSBound,
0081                                                     upperIPSBound,
0082                                                     false,
0083                                                     true,
0084                                                     true,
0085                                                     "b",
0086                                                     trackIPDir,
0087                                                     mc,
0088                                                     makeQualityPlots_,
0089                                                     ibook_));
0090   }
0091   tkcntHistosSig3D.push_back(std::make_unique<TrackIPHistograms<double>>("ips_3D" + theExtensionString,
0092                                                                          "3D IP significance",
0093                                                                          nBinsIPS,
0094                                                                          lowerIPSBound,
0095                                                                          upperIPSBound,
0096                                                                          false,
0097                                                                          true,
0098                                                                          true,
0099                                                                          "b",
0100                                                                          trackIPDir,
0101                                                                          mc,
0102                                                                          makeQualityPlots_,
0103                                                                          ibook_));
0104 
0105   //2D
0106   for (unsigned int i = 1; i <= 4; i++) {
0107     tkcntHistosSig2D.push_back(
0108         std::make_unique<TrackIPHistograms<double>>("ips" + std::to_string(i) + "_2D" + theExtensionString,
0109                                                     "2D IP significance " + std::to_string(i) + ".trk",
0110                                                     nBinsIPS,
0111                                                     lowerIPSBound,
0112                                                     upperIPSBound,
0113                                                     false,
0114                                                     true,
0115                                                     true,
0116                                                     "b",
0117                                                     trackIPDir,
0118                                                     mc,
0119                                                     makeQualityPlots_,
0120                                                     ibook_));
0121   }
0122 
0123   tkcntHistosSig2D.push_back(std::make_unique<TrackIPHistograms<double>>("ips_2D" + theExtensionString,
0124                                                                          "2D IP significance",
0125                                                                          nBinsIPS,
0126                                                                          lowerIPSBound,
0127                                                                          upperIPSBound,
0128                                                                          false,
0129                                                                          true,
0130                                                                          true,
0131                                                                          "b",
0132                                                                          trackIPDir,
0133                                                                          mc,
0134                                                                          makeQualityPlots_,
0135                                                                          ibook_));
0136 
0137   // IP value
0138   //3D
0139   for (unsigned int i = 1; i <= 4; i++) {
0140     tkcntHistosVal3D.push_back(
0141         std::make_unique<TrackIPHistograms<double>>("ip" + std::to_string(i) + "_3D" + theExtensionString,
0142                                                     "3D IP value " + std::to_string(i) + ".trk",
0143                                                     nBinsIP,
0144                                                     lowerIPBound,
0145                                                     upperIPBound,
0146                                                     false,
0147                                                     true,
0148                                                     true,
0149                                                     "b",
0150                                                     trackIPDir,
0151                                                     mc,
0152                                                     makeQualityPlots_,
0153                                                     ibook_));
0154   }
0155 
0156   tkcntHistosVal3D.push_back(std::make_unique<TrackIPHistograms<double>>("ip_3D" + theExtensionString,
0157                                                                          "3D IP value",
0158                                                                          nBinsIP,
0159                                                                          lowerIPBound,
0160                                                                          upperIPBound,
0161                                                                          false,
0162                                                                          true,
0163                                                                          true,
0164                                                                          "b",
0165                                                                          trackIPDir,
0166                                                                          mc,
0167                                                                          makeQualityPlots_,
0168                                                                          ibook_));
0169 
0170   //2D
0171   for (unsigned int i = 1; i <= 4; i++) {
0172     tkcntHistosVal2D.push_back(
0173         std::make_unique<TrackIPHistograms<double>>("ip" + std::to_string(i) + "_2D" + theExtensionString,
0174                                                     "2D IP value " + std::to_string(i) + ".trk",
0175                                                     nBinsIP,
0176                                                     lowerIPBound,
0177                                                     upperIPBound,
0178                                                     false,
0179                                                     true,
0180                                                     true,
0181                                                     "b",
0182                                                     trackIPDir,
0183                                                     mc,
0184                                                     makeQualityPlots_,
0185                                                     ibook_));
0186   }
0187 
0188   tkcntHistosVal2D.push_back(std::make_unique<TrackIPHistograms<double>>("ip_2D" + theExtensionString,
0189                                                                          "2D IP value",
0190                                                                          nBinsIP,
0191                                                                          lowerIPBound,
0192                                                                          upperIPBound,
0193                                                                          false,
0194                                                                          true,
0195                                                                          true,
0196                                                                          "b",
0197                                                                          trackIPDir,
0198                                                                          mc,
0199                                                                          makeQualityPlots_,
0200                                                                          ibook_));
0201 
0202   // IP error
0203   // 3D
0204   for (unsigned int i = 1; i <= 4; i++) {
0205     tkcntHistosErr3D.push_back(
0206         std::make_unique<TrackIPHistograms<double>>("ipe" + std::to_string(i) + "_3D" + theExtensionString,
0207                                                     "3D IP error " + std::to_string(i) + ".trk",
0208                                                     nBinsIPE,
0209                                                     lowerIPEBound,
0210                                                     upperIPEBound,
0211                                                     false,
0212                                                     true,
0213                                                     true,
0214                                                     "b",
0215                                                     trackIPDir,
0216                                                     mc,
0217                                                     makeQualityPlots_,
0218                                                     ibook_));
0219   }
0220 
0221   tkcntHistosErr3D.push_back(std::make_unique<TrackIPHistograms<double>>("ipe_3D" + theExtensionString,
0222                                                                          "3D IP error",
0223                                                                          nBinsIPE,
0224                                                                          lowerIPEBound,
0225                                                                          upperIPEBound,
0226                                                                          false,
0227                                                                          true,
0228                                                                          true,
0229                                                                          "b",
0230                                                                          trackIPDir,
0231                                                                          mc,
0232                                                                          makeQualityPlots_,
0233                                                                          ibook_));
0234 
0235   //2D
0236   for (unsigned int i = 1; i <= 4; i++) {
0237     tkcntHistosErr2D.push_back(
0238         std::make_unique<TrackIPHistograms<double>>("ipe" + std::to_string(i) + "_2D" + theExtensionString,
0239                                                     "2D IP error " + std::to_string(i) + ".trk",
0240                                                     nBinsIPE,
0241                                                     lowerIPEBound,
0242                                                     upperIPEBound,
0243                                                     false,
0244                                                     true,
0245                                                     true,
0246                                                     "b",
0247                                                     trackIPDir,
0248                                                     mc,
0249                                                     makeQualityPlots_,
0250                                                     ibook_));
0251   }
0252 
0253   tkcntHistosErr2D.push_back(std::make_unique<TrackIPHistograms<double>>("ipe_2D" + theExtensionString,
0254                                                                          "2D IP error",
0255                                                                          nBinsIPE,
0256                                                                          lowerIPEBound,
0257                                                                          upperIPEBound,
0258                                                                          false,
0259                                                                          true,
0260                                                                          true,
0261                                                                          "b",
0262                                                                          trackIPDir,
0263                                                                          mc,
0264                                                                          makeQualityPlots_,
0265                                                                          ibook_));
0266 
0267   // decay length
0268   //2D
0269   for (unsigned int i = 1; i <= 4; i++) {
0270     tkcntHistosDecayLengthVal2D.push_back(
0271         std::make_unique<TrackIPHistograms<double>>("decLen" + std::to_string(i) + "_2D" + theExtensionString,
0272                                                     "2D Decay Length " + std::to_string(i) + ".trk",
0273                                                     50,
0274                                                     0.0,
0275                                                     5.0,
0276                                                     false,
0277                                                     true,
0278                                                     true,
0279                                                     "b",
0280                                                     trackIPDir,
0281                                                     mc,
0282                                                     makeQualityPlots_,
0283                                                     ibook_));
0284   }
0285 
0286   tkcntHistosDecayLengthVal2D.push_back(std::make_unique<TrackIPHistograms<double>>("decLen_2D" + theExtensionString,
0287                                                                                     "Decay Length 2D",
0288                                                                                     50,
0289                                                                                     0.0,
0290                                                                                     5.0,
0291                                                                                     false,
0292                                                                                     true,
0293                                                                                     true,
0294                                                                                     "b",
0295                                                                                     trackIPDir,
0296                                                                                     mc,
0297                                                                                     makeQualityPlots_,
0298                                                                                     ibook_));
0299 
0300   // 3D
0301   for (unsigned int i = 1; i <= 4; i++) {
0302     tkcntHistosDecayLengthVal3D.push_back(
0303         std::make_unique<TrackIPHistograms<double>>("decLen" + std::to_string(i) + "_3D" + theExtensionString,
0304                                                     "3D Decay Length " + std::to_string(i) + ".trk",
0305                                                     50,
0306                                                     0.0,
0307                                                     5.0,
0308                                                     false,
0309                                                     true,
0310                                                     true,
0311                                                     "b",
0312                                                     trackIPDir,
0313                                                     mc,
0314                                                     makeQualityPlots_,
0315                                                     ibook_));
0316   }
0317 
0318   tkcntHistosDecayLengthVal3D.push_back(std::make_unique<TrackIPHistograms<double>>("decLen_3D" + theExtensionString,
0319                                                                                     "3D Decay Length",
0320                                                                                     50,
0321                                                                                     0.0,
0322                                                                                     5.0,
0323                                                                                     false,
0324                                                                                     true,
0325                                                                                     true,
0326                                                                                     "b",
0327                                                                                     trackIPDir,
0328                                                                                     mc,
0329                                                                                     makeQualityPlots_,
0330                                                                                     ibook_));
0331 
0332   // jet distance
0333   //2D
0334   for (unsigned int i = 1; i <= 4; i++) {
0335     tkcntHistosJetDistVal2D.push_back(
0336         std::make_unique<TrackIPHistograms<double>>("jetDist" + std::to_string(i) + "_2D" + theExtensionString,
0337                                                     "JetDistance 2D " + std::to_string(i) + ".trk",
0338                                                     50,
0339                                                     -0.1,
0340                                                     0.0,
0341                                                     false,
0342                                                     true,
0343                                                     true,
0344                                                     "b",
0345                                                     trackIPDir,
0346                                                     mc,
0347                                                     makeQualityPlots_,
0348                                                     ibook_));
0349   }
0350 
0351   tkcntHistosJetDistVal2D.push_back(std::make_unique<TrackIPHistograms<double>>("jetDist_2D" + theExtensionString,
0352                                                                                 "JetDistance 2D",
0353                                                                                 50,
0354                                                                                 -0.1,
0355                                                                                 0.0,
0356                                                                                 false,
0357                                                                                 true,
0358                                                                                 true,
0359                                                                                 "b",
0360                                                                                 trackIPDir,
0361                                                                                 mc,
0362                                                                                 makeQualityPlots_,
0363                                                                                 ibook_));
0364 
0365   // 3D
0366   for (unsigned int i = 1; i <= 4; i++) {
0367     tkcntHistosJetDistVal3D.push_back(
0368         std::make_unique<TrackIPHistograms<double>>("jetDist" + std::to_string(i) + "_3D" + theExtensionString,
0369                                                     "JetDistance 3D " + std::to_string(i) + ".trk",
0370                                                     50,
0371                                                     -0.1,
0372                                                     0.0,
0373                                                     false,
0374                                                     true,
0375                                                     true,
0376                                                     "b",
0377                                                     trackIPDir,
0378                                                     mc,
0379                                                     makeQualityPlots_,
0380                                                     ibook_));
0381   }
0382 
0383   tkcntHistosJetDistVal3D.push_back(std::make_unique<TrackIPHistograms<double>>("jetDist_3D" + theExtensionString,
0384                                                                                 "JetDistance 3D",
0385                                                                                 50,
0386                                                                                 -0.1,
0387                                                                                 0.0,
0388                                                                                 false,
0389                                                                                 true,
0390                                                                                 true,
0391                                                                                 "b",
0392                                                                                 trackIPDir,
0393                                                                                 mc,
0394                                                                                 makeQualityPlots_,
0395                                                                                 ibook_));
0396 
0397   // track chi-squared
0398   // 2D
0399   for (unsigned int i = 1; i <= 4; i++) {
0400     tkcntHistosTkNChiSqr2D.push_back(
0401         std::make_unique<TrackIPHistograms<double>>("tkNChiSqr" + std::to_string(i) + "_2D" + theExtensionString,
0402                                                     "Normalized Chi Squared 2D " + std::to_string(i) + ".trk",
0403                                                     50,
0404                                                     -0.1,
0405                                                     10.0,
0406                                                     false,
0407                                                     true,
0408                                                     true,
0409                                                     "b",
0410                                                     trackIPDir,
0411                                                     mc,
0412                                                     makeQualityPlots_,
0413                                                     ibook_));
0414   }
0415 
0416   tkcntHistosTkNChiSqr2D.push_back(std::make_unique<TrackIPHistograms<double>>("tkNChiSqr_2D" + theExtensionString,
0417                                                                                "Normalized Chi Squared 2D",
0418                                                                                50,
0419                                                                                -0.1,
0420                                                                                10.0,
0421                                                                                false,
0422                                                                                true,
0423                                                                                true,
0424                                                                                "b",
0425                                                                                trackIPDir,
0426                                                                                mc,
0427                                                                                makeQualityPlots_,
0428                                                                                ibook_));
0429 
0430   // 3D
0431   for (unsigned int i = 1; i <= 4; i++) {
0432     tkcntHistosTkNChiSqr3D.push_back(
0433         std::make_unique<TrackIPHistograms<double>>("tkNChiSqr" + std::to_string(i) + "_3D" + theExtensionString,
0434                                                     "Normalized Chi Squared 3D " + std::to_string(i) + ".trk",
0435                                                     50,
0436                                                     -0.1,
0437                                                     10.0,
0438                                                     false,
0439                                                     true,
0440                                                     true,
0441                                                     "b",
0442                                                     trackIPDir,
0443                                                     mc,
0444                                                     makeQualityPlots_,
0445                                                     ibook_));
0446   }
0447 
0448   tkcntHistosTkNChiSqr3D.push_back(std::make_unique<TrackIPHistograms<double>>("tkNChiSqr_3D" + theExtensionString,
0449                                                                                "Normalized Chi Squared 3D",
0450                                                                                50,
0451                                                                                -0.1,
0452                                                                                10.0,
0453                                                                                false,
0454                                                                                true,
0455                                                                                true,
0456                                                                                "b",
0457                                                                                trackIPDir,
0458                                                                                mc,
0459                                                                                makeQualityPlots_,
0460                                                                                ibook_));
0461 
0462   // track pT
0463   // 2D
0464   for (unsigned int i = 1; i <= 4; i++) {
0465     tkcntHistosTkPt2D.push_back(
0466         std::make_unique<TrackIPHistograms<double>>("tkPt" + std::to_string(i) + "_2D" + theExtensionString,
0467                                                     "Track Pt 2D " + std::to_string(i) + ".trk",
0468                                                     50,
0469                                                     -0.1,
0470                                                     50.1,
0471                                                     false,
0472                                                     true,
0473                                                     true,
0474                                                     "b",
0475                                                     trackIPDir,
0476                                                     mc,
0477                                                     makeQualityPlots_,
0478                                                     ibook_));
0479   }
0480 
0481   tkcntHistosTkPt2D.push_back(std::make_unique<TrackIPHistograms<double>>("tkPt_2D" + theExtensionString,
0482                                                                           "Track Pt 2D",
0483                                                                           50,
0484                                                                           -0.1,
0485                                                                           50.1,
0486                                                                           false,
0487                                                                           true,
0488                                                                           true,
0489                                                                           "b",
0490                                                                           trackIPDir,
0491                                                                           mc,
0492                                                                           makeQualityPlots_,
0493                                                                           ibook_));
0494 
0495   // 3D
0496   for (unsigned int i = 1; i <= 4; i++) {
0497     tkcntHistosTkPt3D.push_back(
0498         std::make_unique<TrackIPHistograms<double>>("tkPt" + std::to_string(i) + "_3D" + theExtensionString,
0499                                                     "Track Pt 3D " + std::to_string(i) + ".trk",
0500                                                     50,
0501                                                     -0.1,
0502                                                     50.1,
0503                                                     false,
0504                                                     true,
0505                                                     true,
0506                                                     "b",
0507                                                     trackIPDir,
0508                                                     mc,
0509                                                     makeQualityPlots_,
0510                                                     ibook_));
0511   }
0512 
0513   tkcntHistosTkPt3D.push_back(std::make_unique<TrackIPHistograms<double>>("tkPt_3D" + theExtensionString,
0514                                                                           "Track Pt 3D",
0515                                                                           50,
0516                                                                           -0.1,
0517                                                                           50.1,
0518                                                                           false,
0519                                                                           true,
0520                                                                           true,
0521                                                                           "b",
0522                                                                           trackIPDir,
0523                                                                           mc,
0524                                                                           makeQualityPlots_,
0525                                                                           ibook_));
0526 
0527   // track nHits
0528   // 2D
0529   for (unsigned int i = 1; i <= 4; i++) {
0530     tkcntHistosTkNHits2D.push_back(
0531         std::make_unique<TrackIPHistograms<int>>("tkNHits" + std::to_string(i) + "_2D" + theExtensionString,
0532                                                  "Track NHits 2D " + std::to_string(i) + ".trk",
0533                                                  31,
0534                                                  -0.5,
0535                                                  30.5,
0536                                                  false,
0537                                                  true,
0538                                                  true,
0539                                                  "b",
0540                                                  trackIPDir,
0541                                                  mc,
0542                                                  makeQualityPlots_,
0543                                                  ibook_));
0544   }
0545 
0546   tkcntHistosTkNHits2D.push_back(std::make_unique<TrackIPHistograms<int>>("tkNHits_2D" + theExtensionString,
0547                                                                           "Track NHits 2D",
0548                                                                           31,
0549                                                                           -0.5,
0550                                                                           30.5,
0551                                                                           false,
0552                                                                           true,
0553                                                                           true,
0554                                                                           "b",
0555                                                                           trackIPDir,
0556                                                                           mc,
0557                                                                           makeQualityPlots_,
0558                                                                           ibook_));
0559 
0560   // 3D
0561   for (unsigned int i = 1; i <= 4; i++) {
0562     tkcntHistosTkNHits3D.push_back(
0563         std::make_unique<TrackIPHistograms<int>>("tkNHits" + std::to_string(i) + "_3D" + theExtensionString,
0564                                                  "Track NHits 3D " + std::to_string(i) + ".trk",
0565                                                  31,
0566                                                  -0.5,
0567                                                  30.5,
0568                                                  false,
0569                                                  true,
0570                                                  true,
0571                                                  "b",
0572                                                  trackIPDir,
0573                                                  mc,
0574                                                  makeQualityPlots_,
0575                                                  ibook_));
0576   }
0577 
0578   tkcntHistosTkNHits3D.push_back(std::make_unique<TrackIPHistograms<int>>("tkNHits_3D" + theExtensionString,
0579                                                                           "Track NHits 3D",
0580                                                                           31,
0581                                                                           -0.5,
0582                                                                           30.5,
0583                                                                           false,
0584                                                                           true,
0585                                                                           true,
0586                                                                           "b",
0587                                                                           trackIPDir,
0588                                                                           mc,
0589                                                                           makeQualityPlots_,
0590                                                                           ibook_));
0591 
0592   //Pixel hits
0593   // 2D
0594   for (unsigned int i = 1; i <= 4; i++) {
0595     tkcntHistosTkNPixelHits2D.push_back(
0596         std::make_unique<TrackIPHistograms<int>>("tkNPixelHits" + std::to_string(i) + "_2D" + theExtensionString,
0597                                                  "Track NPixelHits 2D " + std::to_string(i) + ".trk",
0598                                                  11,
0599                                                  -0.5,
0600                                                  10.5,
0601                                                  false,
0602                                                  true,
0603                                                  true,
0604                                                  "b",
0605                                                  trackIPDir,
0606                                                  mc,
0607                                                  makeQualityPlots_,
0608                                                  ibook_));
0609   }
0610 
0611   tkcntHistosTkNPixelHits2D.push_back(std::make_unique<TrackIPHistograms<int>>("tkNPixelHits_2D" + theExtensionString,
0612                                                                                "Track NPixelHits 2D",
0613                                                                                11,
0614                                                                                -0.5,
0615                                                                                10.5,
0616                                                                                false,
0617                                                                                true,
0618                                                                                true,
0619                                                                                "b",
0620                                                                                trackIPDir,
0621                                                                                mc,
0622                                                                                makeQualityPlots_,
0623                                                                                ibook_));
0624 
0625   // 3D
0626   for (unsigned int i = 1; i <= 4; i++) {
0627     tkcntHistosTkNPixelHits3D.push_back(
0628         std::make_unique<TrackIPHistograms<int>>("tkNPixelHits" + std::to_string(i) + "_3D" + theExtensionString,
0629                                                  "Track NPixelHits 3D " + std::to_string(i) + ".trk",
0630                                                  11,
0631                                                  -0.5,
0632                                                  10.5,
0633                                                  false,
0634                                                  true,
0635                                                  true,
0636                                                  "b",
0637                                                  trackIPDir,
0638                                                  mc,
0639                                                  makeQualityPlots_,
0640                                                  ibook_));
0641   }
0642 
0643   tkcntHistosTkNPixelHits3D.push_back(std::make_unique<TrackIPHistograms<int>>("tkNPixelHits_3D" + theExtensionString,
0644                                                                                "Track NPixelHits 3D",
0645                                                                                11,
0646                                                                                -0.5,
0647                                                                                10.5,
0648                                                                                false,
0649                                                                                true,
0650                                                                                true,
0651                                                                                "b",
0652                                                                                trackIPDir,
0653                                                                                mc,
0654                                                                                makeQualityPlots_,
0655                                                                                ibook_));
0656 
0657   // probability
0658   // 3D
0659   for (unsigned int i = 1; i <= 4; i++) {
0660     tkcntHistosProb3D.push_back(
0661         std::make_unique<TrackIPHistograms<float>>("prob" + std::to_string(i) + "_3D" + theExtensionString,
0662                                                    "3D IP probability " + std::to_string(i) + ".trk",
0663                                                    52,
0664                                                    -1.04,
0665                                                    1.04,
0666                                                    false,
0667                                                    true,
0668                                                    true,
0669                                                    "b",
0670                                                    trackIPDir,
0671                                                    mc,
0672                                                    makeQualityPlots_,
0673                                                    ibook_));
0674   }
0675 
0676   tkcntHistosProb3D.push_back(std::make_unique<TrackIPHistograms<float>>("prob_3D" + theExtensionString,
0677                                                                          "3D IP probability",
0678                                                                          50,
0679                                                                          -1.04,
0680                                                                          1.04,
0681                                                                          false,
0682                                                                          true,
0683                                                                          true,
0684                                                                          "b",
0685                                                                          trackIPDir,
0686                                                                          mc,
0687                                                                          makeQualityPlots_,
0688                                                                          ibook_));
0689 
0690   // 2D
0691   for (unsigned int i = 1; i <= 4; i++) {
0692     tkcntHistosProb2D.push_back(
0693         std::make_unique<TrackIPHistograms<float>>("prob" + std::to_string(i) + "_2D" + theExtensionString,
0694                                                    "2D IP probability " + std::to_string(i) + ".trk",
0695                                                    52,
0696                                                    -1.04,
0697                                                    1.04,
0698                                                    false,
0699                                                    true,
0700                                                    true,
0701                                                    "b",
0702                                                    trackIPDir,
0703                                                    mc,
0704                                                    makeQualityPlots_,
0705                                                    ibook_));
0706   }
0707 
0708   tkcntHistosProb2D.push_back(std::make_unique<TrackIPHistograms<float>>("prob_2D" + theExtensionString,
0709                                                                          "2D IP probability",
0710                                                                          52,
0711                                                                          -1.04,
0712                                                                          1.04,
0713                                                                          false,
0714                                                                          true,
0715                                                                          true,
0716                                                                          "b",
0717                                                                          trackIPDir,
0718                                                                          mc,
0719                                                                          makeQualityPlots_,
0720                                                                          ibook_));
0721 
0722   // probability for tracks with IP value < 0 or IP value >> 0
0723   tkcntHistosTkProbIPneg2D = std::make_unique<TrackIPHistograms<float>>("probIPneg_2D" + theExtensionString,
0724                                                                         "2D negative IP probability",
0725                                                                         52,
0726                                                                         -1.04,
0727                                                                         1.04,
0728                                                                         false,
0729                                                                         true,
0730                                                                         true,
0731                                                                         "b",
0732                                                                         trackIPDir,
0733                                                                         mc,
0734                                                                         makeQualityPlots_,
0735                                                                         ibook_);
0736   tkcntHistosTkProbIPpos2D = std::make_unique<TrackIPHistograms<float>>("probIPpos_2D" + theExtensionString,
0737                                                                         "2D positive IP probability",
0738                                                                         52,
0739                                                                         -1.04,
0740                                                                         1.04,
0741                                                                         false,
0742                                                                         true,
0743                                                                         true,
0744                                                                         "b",
0745                                                                         trackIPDir,
0746                                                                         mc,
0747                                                                         makeQualityPlots_,
0748                                                                         ibook_);
0749   tkcntHistosTkProbIPneg3D = std::make_unique<TrackIPHistograms<float>>("probIPneg_3D" + theExtensionString,
0750                                                                         "3D negative IP probability",
0751                                                                         52,
0752                                                                         -1.04,
0753                                                                         1.04,
0754                                                                         false,
0755                                                                         true,
0756                                                                         true,
0757                                                                         "b",
0758                                                                         trackIPDir,
0759                                                                         mc,
0760                                                                         makeQualityPlots_,
0761                                                                         ibook_);
0762   tkcntHistosTkProbIPpos3D = std::make_unique<TrackIPHistograms<float>>("probIPpos_3D" + theExtensionString,
0763                                                                         "3D positive IP probability",
0764                                                                         52,
0765                                                                         -1.04,
0766                                                                         1.04,
0767                                                                         false,
0768                                                                         true,
0769                                                                         true,
0770                                                                         "b",
0771                                                                         trackIPDir,
0772                                                                         mc,
0773                                                                         makeQualityPlots_,
0774                                                                         ibook_);
0775 
0776   // ghost Tracks and others
0777   ghostTrackDistanceValuHisto = std::make_unique<TrackIPHistograms<double>>("ghostTrackDist" + theExtensionString,
0778                                                                             "GhostTrackDistance",
0779                                                                             50,
0780                                                                             0.0,
0781                                                                             0.1,
0782                                                                             false,
0783                                                                             true,
0784                                                                             true,
0785                                                                             "b",
0786                                                                             trackIPDir,
0787                                                                             mc,
0788                                                                             makeQualityPlots_,
0789                                                                             ibook_);
0790   ghostTrackDistanceSignHisto = std::make_unique<TrackIPHistograms<double>>("ghostTrackDistSign" + theExtensionString,
0791                                                                             "GhostTrackDistance significance",
0792                                                                             50,
0793                                                                             -5.0,
0794                                                                             15.0,
0795                                                                             false,
0796                                                                             true,
0797                                                                             true,
0798                                                                             "b",
0799                                                                             trackIPDir,
0800                                                                             mc,
0801                                                                             makeQualityPlots_,
0802                                                                             ibook_);
0803   ghostTrackWeightHisto = std::make_unique<TrackIPHistograms<double>>("ghostTrackWeight" + theExtensionString,
0804                                                                       "GhostTrack fit participation weight",
0805                                                                       50,
0806                                                                       0.0,
0807                                                                       1.0,
0808                                                                       false,
0809                                                                       false,
0810                                                                       true,
0811                                                                       "b",
0812                                                                       trackIPDir,
0813                                                                       mc,
0814                                                                       makeQualityPlots_,
0815                                                                       ibook_);
0816 
0817   trackQualHisto = std::make_unique<FlavourHistograms<int>>("trackQual" + theExtensionString,
0818                                                             "Track Quality of Tracks Associated to Jets",
0819                                                             4,
0820                                                             -1.5,
0821                                                             2.5,
0822                                                             false,
0823                                                             true,
0824                                                             true,
0825                                                             "b",
0826                                                             trackIPDir,
0827                                                             mc,
0828                                                             ibook_);
0829 
0830   selectedTrackQualHisto =
0831       std::make_unique<FlavourHistograms<int>>("selectedTrackQual" + theExtensionString,
0832                                                "Track Quality of Selected Tracks Associated to Jets",
0833                                                4,
0834                                                -1.5,
0835                                                2.5,
0836                                                false,
0837                                                true,
0838                                                true,
0839                                                "b",
0840                                                trackIPDir,
0841                                                mc,
0842                                                ibook_);
0843 
0844   trackMultVsJetPtHisto =
0845       std::make_unique<FlavourHistograms2D<double, int>>("trackMultVsJetPt" + theExtensionString,
0846                                                          "Track Multiplicity vs Jet Pt for Tracks Associated to Jets",
0847                                                          50,
0848                                                          0.0,
0849                                                          250.0,
0850                                                          21,
0851                                                          -0.5,
0852                                                          30.5,
0853                                                          false,
0854                                                          trackIPDir,
0855                                                          mc,
0856                                                          true,
0857                                                          ibook_);
0858 
0859   selectedTrackMultVsJetPtHisto = std::make_unique<FlavourHistograms2D<double, int>>(
0860       "selectedTrackMultVsJetPt" + theExtensionString,
0861       "Track Multiplicity vs Jet Pt for Selected Tracks Associated to Jets",
0862       50,
0863       0.0,
0864       250.0,
0865       21,
0866       -0.5,
0867       20.5,
0868       false,
0869       trackIPDir,
0870       mc,
0871       true,
0872       ibook_);
0873 }
0874 
0875 template <class Container, class Base>
0876 IPTagPlotter<Container, Base>::~IPTagPlotter() {}
0877 
0878 template <class Container, class Base>
0879 void IPTagPlotter<Container, Base>::analyzeTag(const reco::BaseTagInfo* baseTagInfo,
0880                                                double jec,
0881                                                int jetFlavour,
0882                                                float w /*=1*/) {
0883   //  const reco::TrackIPTagInfo * tagInfo =
0884   //    dynamic_cast<const reco::TrackIPTagInfo *>(baseTagInfo);
0885   const reco::IPTagInfo<Container, Base>* tagInfo = dynamic_cast<const reco::IPTagInfo<Container, Base>*>(baseTagInfo);
0886 
0887   if (!tagInfo) {
0888     throw cms::Exception("Configuration")
0889         << "BTagPerformanceAnalyzer: Extended TagInfo not of type TrackIPTagInfo. " << std::endl;
0890   }
0891 
0892   const GlobalPoint pv(tagInfo->primaryVertex()->position().x(),
0893                        tagInfo->primaryVertex()->position().y(),
0894                        tagInfo->primaryVertex()->position().z());
0895 
0896   const std::vector<reco::btag::TrackIPData>& ip = tagInfo->impactParameterData();
0897 
0898   std::vector<float> prob2d, prob3d;
0899   if (tagInfo->hasProbabilities()) {
0900     prob2d = tagInfo->probabilities(0);
0901     prob3d = tagInfo->probabilities(1);
0902   }
0903 
0904   std::vector<std::size_t> sortedIndices = tagInfo->sortedIndexes(reco::btag::IP2DSig);
0905   std::vector<std::size_t> selectedIndices;
0906   Container sortedTracks = tagInfo->sortedTracks(sortedIndices);
0907   Container selectedTracks;
0908   for (unsigned int n = 0; n != sortedIndices.size(); ++n) {
0909     double decayLength = (ip[sortedIndices[n]].closestToJetAxis - pv).mag();
0910     double jetDistance = ip[sortedIndices[n]].distanceToJetAxis.value();
0911     if (decayLength > minDecayLength && decayLength < maxDecayLength && fabs(jetDistance) >= minJetDistance &&
0912         fabs(jetDistance) < maxJetDistance) {
0913       selectedIndices.push_back(sortedIndices[n]);
0914       selectedTracks.push_back(sortedTracks[n]);
0915     }
0916   }
0917 
0918   trkNbr2D->fill(jetFlavour, selectedIndices.size(), w);
0919 
0920   for (unsigned int n = 0; n != selectedIndices.size(); ++n) {
0921     const reco::Track* track = reco::btag::toTrack(selectedTracks[n]);
0922     const reco::TrackBase::TrackQuality& trackQual = highestTrackQual(track);
0923     tkcntHistosSig2D[4]->fill(jetFlavour, trackQual, ip[selectedIndices[n]].ip2d.significance(), true, w);
0924     tkcntHistosVal2D[4]->fill(jetFlavour, trackQual, ip[selectedIndices[n]].ip2d.value(), true, w);
0925     tkcntHistosErr2D[4]->fill(jetFlavour, trackQual, ip[selectedIndices[n]].ip2d.error(), true, w);
0926     const double& decayLen = (ip[selectedIndices[n]].closestToJetAxis - pv).mag();
0927     tkcntHistosDecayLengthVal2D[4]->fill(jetFlavour, trackQual, decayLen, true, w);
0928     tkcntHistosJetDistVal2D[4]->fill(jetFlavour, trackQual, ip[selectedIndices[n]].distanceToJetAxis.value(), true, w);
0929     tkcntHistosTkNChiSqr2D[4]->fill(jetFlavour, trackQual, track->normalizedChi2(), true, w);
0930     tkcntHistosTkPt2D[4]->fill(jetFlavour, trackQual, track->pt(), true, w);
0931     tkcntHistosTkNHits2D[4]->fill(jetFlavour, trackQual, track->found(), true, w);
0932     tkcntHistosTkNPixelHits2D[4]->fill(jetFlavour, trackQual, track->hitPattern().numberOfValidPixelHits(), true, w);
0933     if (n >= 4)
0934       continue;
0935     tkcntHistosSig2D[n]->fill(jetFlavour, trackQual, ip[selectedIndices[n]].ip2d.significance(), true, w);
0936     tkcntHistosVal2D[n]->fill(jetFlavour, trackQual, ip[selectedIndices[n]].ip2d.value(), true, w);
0937     tkcntHistosErr2D[n]->fill(jetFlavour, trackQual, ip[selectedIndices[n]].ip2d.error(), true, w);
0938     tkcntHistosDecayLengthVal2D[n]->fill(jetFlavour, trackQual, decayLen, true, w);
0939     tkcntHistosJetDistVal2D[n]->fill(jetFlavour, trackQual, ip[selectedIndices[n]].distanceToJetAxis.value(), true, w);
0940     tkcntHistosTkNChiSqr2D[n]->fill(jetFlavour, trackQual, track->normalizedChi2(), true, w);
0941     tkcntHistosTkPt2D[n]->fill(jetFlavour, trackQual, track->pt(), true, w);
0942     tkcntHistosTkNHits2D[n]->fill(jetFlavour, trackQual, track->found(), true, w);
0943     tkcntHistosTkNPixelHits2D[n]->fill(jetFlavour, trackQual, track->hitPattern().numberOfValidPixelHits(), true, w);
0944   }
0945   sortedIndices = tagInfo->sortedIndexes(reco::btag::Prob2D);
0946   selectedIndices.clear();
0947   sortedTracks = tagInfo->sortedTracks(sortedIndices);
0948   selectedTracks.clear();
0949   for (unsigned int n = 0; n != sortedIndices.size(); ++n) {
0950     double decayLength = (ip[sortedIndices[n]].closestToJetAxis - pv).mag();
0951     double jetDistance = ip[sortedIndices[n]].distanceToJetAxis.value();
0952     if (decayLength > minDecayLength && decayLength < maxDecayLength && fabs(jetDistance) >= minJetDistance &&
0953         fabs(jetDistance) < maxJetDistance) {
0954       selectedIndices.push_back(sortedIndices[n]);
0955       selectedTracks.push_back(sortedTracks[n]);
0956     }
0957   }
0958   for (unsigned int n = 0; n != selectedIndices.size(); ++n) {
0959     const reco::Track* track = reco::btag::toTrack(selectedTracks[n]);
0960     const reco::TrackBase::TrackQuality& trackQual = highestTrackQual(track);
0961     tkcntHistosProb2D[4]->fill(jetFlavour, trackQual, prob2d[selectedIndices[n]], true, w);
0962     if (ip[selectedIndices[n]].ip2d.value() < 0)
0963       tkcntHistosTkProbIPneg2D->fill(jetFlavour, trackQual, prob2d[selectedIndices[n]], true, w);
0964     else
0965       tkcntHistosTkProbIPpos2D->fill(jetFlavour, trackQual, prob2d[selectedIndices[n]], true, w);
0966     if (n >= 4)
0967       continue;
0968     tkcntHistosProb2D[n]->fill(jetFlavour, trackQual, prob2d[selectedIndices[n]], true, w);
0969   }
0970   for (unsigned int n = selectedIndices.size(); n < 4; ++n) {
0971     const reco::TrackBase::TrackQuality trackQual = reco::TrackBase::undefQuality;
0972     tkcntHistosSig2D[n]->fill(jetFlavour, trackQual, lowerIPSBound - 1.0, false, w);
0973     tkcntHistosVal2D[n]->fill(jetFlavour, trackQual, lowerIPBound - 1.0, false, w);
0974     tkcntHistosErr2D[n]->fill(jetFlavour, trackQual, lowerIPEBound - 1.0, false, w);
0975   }
0976   sortedIndices = tagInfo->sortedIndexes(reco::btag::IP3DSig);
0977   selectedIndices.clear();
0978   sortedTracks = tagInfo->sortedTracks(sortedIndices);
0979   selectedTracks.clear();
0980   for (unsigned int n = 0; n != sortedIndices.size(); ++n) {
0981     double decayLength = (ip[sortedIndices[n]].closestToJetAxis - pv).mag();
0982     double jetDistance = ip[sortedIndices[n]].distanceToJetAxis.value();
0983     if (decayLength > minDecayLength && decayLength < maxDecayLength && fabs(jetDistance) >= minJetDistance &&
0984         fabs(jetDistance) < maxJetDistance) {
0985       selectedIndices.push_back(sortedIndices[n]);
0986       selectedTracks.push_back(sortedTracks[n]);
0987     }
0988   }
0989 
0990   trkNbr3D->fill(jetFlavour, selectedIndices.size(), w);
0991   int nSelectedTracks = selectedIndices.size();
0992 
0993   for (unsigned int n = 0; n != selectedIndices.size(); ++n) {
0994     const reco::Track* track = reco::btag::toTrack(selectedTracks[n]);
0995     const reco::TrackBase::TrackQuality& trackQual = highestTrackQual(track);
0996     tkcntHistosSig3D[4]->fill(jetFlavour, trackQual, ip[selectedIndices[n]].ip3d.significance(), true, w);
0997     tkcntHistosVal3D[4]->fill(jetFlavour, trackQual, ip[selectedIndices[n]].ip3d.value(), true, w);
0998     tkcntHistosErr3D[4]->fill(jetFlavour, trackQual, ip[selectedIndices[n]].ip3d.error(), true, w);
0999     const double& decayLen = (ip[selectedIndices[n]].closestToJetAxis - pv).mag();
1000     tkcntHistosDecayLengthVal3D[4]->fill(jetFlavour, trackQual, decayLen, true, w);
1001     tkcntHistosJetDistVal3D[4]->fill(jetFlavour, trackQual, ip[selectedIndices[n]].distanceToJetAxis.value(), true, w);
1002     tkcntHistosTkNChiSqr3D[4]->fill(jetFlavour, trackQual, track->normalizedChi2(), true, w);
1003     tkcntHistosTkPt3D[4]->fill(jetFlavour, trackQual, track->pt(), true, w);
1004     tkcntHistosTkNHits3D[4]->fill(jetFlavour, trackQual, track->found(), true, w);
1005     tkcntHistosTkNPixelHits3D[4]->fill(jetFlavour, trackQual, track->hitPattern().numberOfValidPixelHits(), true, w);
1006     //ghostTrack infos
1007     ghostTrackDistanceValuHisto->fill(
1008         jetFlavour, trackQual, ip[selectedIndices[n]].distanceToGhostTrack.value(), true, w);
1009     ghostTrackDistanceSignHisto->fill(
1010         jetFlavour, trackQual, ip[selectedIndices[n]].distanceToGhostTrack.significance(), true, w);
1011     ghostTrackWeightHisto->fill(jetFlavour, trackQual, ip[selectedIndices[n]].ghostTrackWeight, true, w);
1012     selectedTrackQualHisto->fill(jetFlavour, trackQual, w);
1013     if (n >= 4)
1014       continue;
1015     tkcntHistosSig3D[n]->fill(jetFlavour, trackQual, ip[selectedIndices[n]].ip3d.significance(), true, w);
1016     tkcntHistosVal3D[n]->fill(jetFlavour, trackQual, ip[selectedIndices[n]].ip3d.value(), true, w);
1017     tkcntHistosErr3D[n]->fill(jetFlavour, trackQual, ip[selectedIndices[n]].ip3d.error(), true, w);
1018     tkcntHistosDecayLengthVal3D[n]->fill(jetFlavour, trackQual, decayLen, true, w);
1019     tkcntHistosJetDistVal3D[n]->fill(jetFlavour, trackQual, ip[selectedIndices[n]].distanceToJetAxis.value(), true, w);
1020     tkcntHistosTkNChiSqr3D[n]->fill(jetFlavour, trackQual, track->normalizedChi2(), true, w);
1021     tkcntHistosTkPt3D[n]->fill(jetFlavour, trackQual, track->pt(), true, w);
1022     tkcntHistosTkNHits3D[n]->fill(jetFlavour, trackQual, track->found(), true, w);
1023     tkcntHistosTkNPixelHits3D[n]->fill(jetFlavour, trackQual, track->hitPattern().numberOfValidPixelHits(), true, w);
1024   }
1025   sortedIndices = tagInfo->sortedIndexes(reco::btag::Prob3D);
1026   selectedIndices.clear();
1027   sortedTracks = tagInfo->sortedTracks(sortedIndices);
1028   selectedTracks.clear();
1029   for (unsigned int n = 0; n != sortedIndices.size(); ++n) {
1030     double decayLength = (ip[sortedIndices[n]].closestToJetAxis - pv).mag();
1031     double jetDistance = ip[sortedIndices[n]].distanceToJetAxis.value();
1032     if (decayLength > minDecayLength && decayLength < maxDecayLength && fabs(jetDistance) >= minJetDistance &&
1033         fabs(jetDistance) < maxJetDistance) {
1034       selectedIndices.push_back(sortedIndices[n]);
1035       selectedTracks.push_back(sortedTracks[n]);
1036     }
1037   }
1038   for (unsigned int n = 0; n != selectedIndices.size(); ++n) {
1039     const reco::Track* track = reco::btag::toTrack(selectedTracks[n]);
1040     const reco::TrackBase::TrackQuality& trackQual = highestTrackQual(track);
1041     tkcntHistosProb3D[4]->fill(jetFlavour, trackQual, prob3d[selectedIndices[n]], true, w);
1042     if (ip[selectedIndices[n]].ip3d.value() < 0)
1043       tkcntHistosTkProbIPneg3D->fill(jetFlavour, trackQual, prob3d[selectedIndices[n]], true, w);
1044     else
1045       tkcntHistosTkProbIPpos3D->fill(jetFlavour, trackQual, prob3d[selectedIndices[n]], true, w);
1046     if (n >= 4)
1047       continue;
1048     tkcntHistosProb3D[n]->fill(jetFlavour, trackQual, prob3d[selectedIndices[n]], true, w);
1049   }
1050   for (unsigned int n = selectedIndices.size(); n < 4; ++n) {
1051     const reco::TrackBase::TrackQuality trackQual = reco::TrackBase::undefQuality;
1052     tkcntHistosSig3D[n]->fill(jetFlavour, trackQual, lowerIPSBound - 1.0, false, w);
1053     tkcntHistosVal3D[n]->fill(jetFlavour, trackQual, lowerIPBound - 1.0, false, w);
1054     tkcntHistosErr3D[n]->fill(jetFlavour, trackQual, lowerIPEBound - 1.0, false, w);
1055   }
1056   for (unsigned int n = 0; n != sortedTracks.size(); ++n) {
1057     trackQualHisto->fill(jetFlavour, highestTrackQual(reco::btag::toTrack(sortedTracks[n])), w);
1058   }
1059 
1060   //still need to implement weights in FlavourHistograms2D
1061   trackMultVsJetPtHisto->fill(jetFlavour, tagInfo->jet()->pt() * jec, sortedTracks.size());
1062   selectedTrackMultVsJetPtHisto->fill(
1063       jetFlavour, tagInfo->jet()->pt() * jec, nSelectedTracks);  //tagInfo->selectedTracks().size());
1064 }
1065 
1066 template <class Container, class Base>
1067 void IPTagPlotter<Container, Base>::finalize(DQMStore::IBooker& ibook_, DQMStore::IGetter& igetter_) {
1068   //
1069   // final processing:
1070   // produce the misid. vs. eff histograms
1071   //
1072   const std::string trackIPDir(theExtensionString.substr(1));
1073 
1074   tkcntHistosSig3D.clear();
1075   tkcntHistosSig2D.clear();
1076   effPurFromHistos.clear();
1077 
1078   for (unsigned int i = 2; i <= 3; i++) {
1079     tkcntHistosSig3D.push_back(
1080         std::make_unique<TrackIPHistograms<double>>("ips" + std::to_string(i) + "_3D" + theExtensionString,
1081                                                     "3D IP significance " + std::to_string(i) + ".trk",
1082                                                     nBinsIPS,
1083                                                     lowerIPSBound,
1084                                                     upperIPSBound,
1085                                                     "b",
1086                                                     trackIPDir,
1087                                                     mcPlots_,
1088                                                     makeQualityPlots_,
1089                                                     igetter_));
1090     effPurFromHistos.push_back(std::make_unique<EffPurFromHistos>(
1091         *tkcntHistosSig3D.back(), trackIPDir, mcPlots_, ibook_, nBinEffPur_, startEffPur_, endEffPur_));
1092   }
1093 
1094   for (unsigned int i = 2; i <= 3; i++) {
1095     tkcntHistosSig2D.push_back(
1096         std::make_unique<TrackIPHistograms<double>>("ips" + std::to_string(i) + "_2D" + theExtensionString,
1097                                                     "2D IP significance " + std::to_string(i) + ".trk",
1098                                                     nBinsIPS,
1099                                                     lowerIPSBound,
1100                                                     upperIPSBound,
1101                                                     "b",
1102                                                     trackIPDir,
1103                                                     mcPlots_,
1104                                                     makeQualityPlots_,
1105                                                     igetter_));
1106     effPurFromHistos.push_back(std::make_unique<EffPurFromHistos>(
1107         *tkcntHistosSig2D.back(), trackIPDir, mcPlots_, ibook_, nBinEffPur_, startEffPur_, endEffPur_));
1108   }
1109 
1110   for (int n = 0; n != 4; ++n)
1111     effPurFromHistos[n]->compute(ibook_);
1112 }
1113 
1114 template <class Container, class Base>
1115 void IPTagPlotter<Container, Base>::psPlot(const std::string& name) {
1116   const std::string cName("TrackIPPlots" + theExtensionString);
1117   RecoBTag::setTDRStyle()->cd();
1118   TCanvas canvas(cName.c_str(), cName.c_str(), 600, 900);
1119   canvas.UseCurrentStyle();
1120   if (willFinalize_) {
1121     for (int n = 0; n != 2; ++n) {
1122       canvas.Print((name + cName + ".ps").c_str());
1123       canvas.Clear();
1124       canvas.Divide(2, 3);
1125       canvas.cd(1);
1126       effPurFromHistos[0 + n]->discriminatorNoCutEffic().plot();
1127       canvas.cd(2);
1128       effPurFromHistos[0 + n]->discriminatorCutEfficScan().plot();
1129       canvas.cd(3);
1130       effPurFromHistos[0 + n]->plot();
1131       canvas.cd(4);
1132       effPurFromHistos[1 + n]->discriminatorNoCutEffic().plot();
1133       canvas.cd(5);
1134       effPurFromHistos[1 + n]->discriminatorCutEfficScan().plot();
1135       canvas.cd(6);
1136       effPurFromHistos[1 + n]->plot();
1137     }
1138     return;
1139   }
1140 
1141   canvas.Print((name + cName + ".ps[").c_str());
1142   canvas.Clear();
1143   canvas.Divide(2, 3);
1144 
1145   canvas.cd(1);
1146   trkNbr3D->plot();
1147   canvas.cd(2);
1148   tkcntHistosSig3D[4]->plot();
1149   for (int n = 0; n < 4; n++) {
1150     canvas.cd(3 + n);
1151     tkcntHistosSig3D[n]->plot();
1152   }
1153 
1154   canvas.Print((name + cName + ".ps").c_str());
1155   canvas.Clear();
1156   canvas.Divide(2, 3);
1157 
1158   canvas.cd(1);
1159   trkNbr3D->plot();
1160   canvas.cd(2);
1161   tkcntHistosProb3D[4]->plot();
1162   for (int n = 0; n < 4; n++) {
1163     canvas.cd(3 + n);
1164     tkcntHistosProb3D[n]->plot();
1165   }
1166 
1167   canvas.Print((name + cName + ".ps").c_str());
1168   canvas.Clear();
1169   canvas.Divide(2, 3);
1170   canvas.cd(1);
1171   trkNbr2D->plot();
1172   canvas.cd(2);
1173   tkcntHistosSig2D[4]->plot();
1174   for (int n = 0; n != 4; ++n) {
1175     canvas.cd(3 + n);
1176     tkcntHistosSig2D[n]->plot();
1177   }
1178 
1179   canvas.Print((name + cName + ".ps").c_str());
1180   canvas.Clear();
1181   canvas.Divide(2, 3);
1182   canvas.cd(1);
1183   trkNbr2D->plot();
1184   canvas.cd(2);
1185   tkcntHistosProb2D[4]->plot();
1186   for (int n = 0; n != 4; ++n) {
1187     canvas.cd(3 + n);
1188     tkcntHistosProb2D[n]->plot();
1189   }
1190 
1191   canvas.Print((name + cName + ".ps").c_str());
1192   canvas.Clear();
1193   canvas.Divide(1, 3);
1194   canvas.cd(1);
1195   ghostTrackDistanceValuHisto->plot();
1196   canvas.cd(2);
1197   ghostTrackDistanceSignHisto->plot();
1198   canvas.cd(3);
1199   ghostTrackWeightHisto->plot();
1200 
1201   canvas.Print((name + cName + ".ps").c_str());
1202   canvas.Print((name + cName + ".ps]").c_str());
1203 }
1204 
1205 template <class Container, class Base>
1206 void IPTagPlotter<Container, Base>::epsPlot(const std::string& name) {
1207   if (willFinalize_) {
1208     for (int n = 0; n != 4; ++n)
1209       effPurFromHistos[n]->epsPlot(name);
1210     return;
1211   }
1212   trkNbr2D->epsPlot(name);
1213   trkNbr3D->epsPlot(name);
1214   tkcntHistosTkProbIPneg2D->epsPlot(name);
1215   tkcntHistosTkProbIPpos2D->epsPlot(name);
1216   tkcntHistosTkProbIPneg3D->epsPlot(name);
1217   tkcntHistosTkProbIPpos3D->epsPlot(name);
1218   ghostTrackDistanceValuHisto->epsPlot(name);
1219   ghostTrackDistanceSignHisto->epsPlot(name);
1220   ghostTrackWeightHisto->epsPlot(name);
1221   for (int n = 0; n != 5; ++n) {
1222     tkcntHistosSig2D[n]->epsPlot(name);
1223     tkcntHistosSig3D[n]->epsPlot(name);
1224     tkcntHistosVal2D[n]->epsPlot(name);
1225     tkcntHistosVal3D[n]->epsPlot(name);
1226     tkcntHistosErr2D[n]->epsPlot(name);
1227     tkcntHistosErr3D[n]->epsPlot(name);
1228     tkcntHistosProb2D[n]->epsPlot(name);
1229     tkcntHistosProb3D[n]->epsPlot(name);
1230     tkcntHistosDecayLengthVal2D[n]->epsPlot(name);
1231     tkcntHistosDecayLengthVal3D[n]->epsPlot(name);
1232     tkcntHistosJetDistVal2D[n]->epsPlot(name);
1233     tkcntHistosJetDistVal3D[n]->epsPlot(name);
1234     tkcntHistosTkNChiSqr2D[n]->epsPlot(name);
1235     tkcntHistosTkNChiSqr3D[n]->epsPlot(name);
1236     tkcntHistosTkPt2D[n]->epsPlot(name);
1237     tkcntHistosTkPt3D[n]->epsPlot(name);
1238     tkcntHistosTkNHits2D[n]->epsPlot(name);
1239     tkcntHistosTkNHits3D[n]->epsPlot(name);
1240     tkcntHistosTkNPixelHits2D[n]->epsPlot(name);
1241     tkcntHistosTkNPixelHits3D[n]->epsPlot(name);
1242   }
1243 }
1244 
1245 template <class Container, class Base>
1246 reco::TrackBase::TrackQuality IPTagPlotter<Container, Base>::highestTrackQual(const reco::Track* track) const {
1247   for (reco::TrackBase::TrackQuality i = reco::TrackBase::highPurity; i != reco::TrackBase::undefQuality;
1248        i = reco::TrackBase::TrackQuality(i - 1)) {
1249     if (track->quality(i))
1250       return i;
1251   }
1252 
1253   return reco::TrackBase::undefQuality;
1254 }
1255 
1256 #endif  // DQMOffline_RecoB_IPTagPlotter_cc_h