SiPixelQualityBadFractionComparisonBase

SiPixelQualityBadFractionMap

SiPixelQualityBadRocsSummary

SiPixelQualityBadRocsTimeHistory

SiPixelQualityDebugger

SiPixelQualityMap

SiPixelQualityMapComparisonBase

SiPixelQualityTest

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732
/*!
  \file SiPixelQuality_PayloadInspector
  \Payload Inspector Plugin for SiPixelQuality
  \author M. Musich
  \version $Revision: 1.0 $
  \date $Date: 2018/10/18 14:48:00 $
*/

#include "FWCore/MessageLogger/interface/MessageLogger.h"

#include "CondCore/Utilities/interface/PayloadInspectorModule.h"
#include "CondCore/Utilities/interface/PayloadInspector.h"
#include "CondCore/CondDB/interface/Time.h"
#include "CondCore/SiPixelPlugins/interface/SiPixelPayloadInspectorHelper.h"
#include "FWCore/ParameterSet/interface/FileInPath.h"
#include "CalibTracker/StandaloneTrackerTopology/interface/StandaloneTrackerTopology.h"
#include "CalibTracker/SiPixelESProducers/interface/SiPixelDetInfoFileReader.h"

// the data format of the condition to be inspected
#include "CondFormats/SiPixelObjects/interface/SiPixelQuality.h"
#include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
#include "DataFormats/DetId/interface/DetId.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "DQM/TrackerRemapper/interface/Phase1PixelROCMaps.h"
#include "DQM/TrackerRemapper/interface/Phase1PixelSummaryMap.h"

#include <memory>
#include <sstream>
#include <iostream>

// include ROOT
#include "TH2F.h"
#include "TLegend.h"
#include "TCanvas.h"
#include "TLine.h"
#include "TGraph.h"
#include "TStyle.h"
#include "TLatex.h"
#include "TPave.h"
#include "TPaveStats.h"

namespace {

  using namespace cond::payloadInspector;

  /************************************************
    test class
  *************************************************/

  class SiPixelQualityTest : public Histogram1D<SiPixelQuality, SINGLE_IOV> {
  public:
    SiPixelQualityTest()
        : Histogram1D<SiPixelQuality, SINGLE_IOV>("SiPixelQuality test", "SiPixelQuality test", 10, 0.0, 10.0) {}

    bool fill() override {
      auto tag = PlotBase::getTag<0>();
      for (auto const& iov : tag.iovs) {
        std::shared_ptr<SiPixelQuality> payload = Base::fetchPayload(std::get<1>(iov));
        if (payload.get()) {
          fillWithValue(1.);

          auto theDisabledModules = payload->getBadComponentList();
          for (const auto& mod : theDisabledModules) {
            int BadRocCount(0);
            for (unsigned short n = 0; n < 16; n++) {
              unsigned short mask = 1 << n;  // 1 << n = 2^{n} using bitwise shift
              if (mod.BadRocs & mask)
                BadRocCount++;
            }
            COUT << "detId:" << mod.DetID << " error type:" << mod.errorType << " BadRocs:" << BadRocCount << std::endl;
          }
        }  // payload
      }  // iovs
      return true;
    }  // fill
  };

  /************************************************
    Debugging class, to not be displayed
  *************************************************/

  class SiPixelQualityDebugger : public Histogram1D<SiPixelQuality, SINGLE_IOV> {
  public:
    SiPixelQualityDebugger()
        : Histogram1D<SiPixelQuality, SINGLE_IOV>("SiPixelQuality test", "SiPixelQuality test", 10, 0.0, 10.0) {}

    bool fill() override {
      auto tag = PlotBase::getTag<0>();
      for (auto const& iov : tag.iovs) {
        std::shared_ptr<SiPixelQuality> payload = Base::fetchPayload(std::get<1>(iov));
        if (payload.get()) {
          fillWithValue(1.);

          SiPixelPI::topolInfo t_info_fromXML;
          t_info_fromXML.init();

          auto theDisabledModules = payload->getBadComponentList();
          for (const auto& mod : theDisabledModules) {
            DetId detid(mod.DetID);
            auto PhInfo = SiPixelPI::PhaseInfo(SiPixelPI::phase1size);
            const char* path_toTopologyXML = PhInfo.pathToTopoXML();
            auto tTopo =
                StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath(path_toTopologyXML).fullPath());
            t_info_fromXML.fillGeometryInfo(detid, tTopo, PhInfo.phase());
            std::stringstream ss;
            t_info_fromXML.printAll(ss);
            std::bitset<16> bad_rocs(mod.BadRocs);

            if (t_info_fromXML.subDetId() == 1 && t_info_fromXML.layer() == 1) {
              std::cout << ss.str() << " s_module: " << SiPixelPI::signed_module(mod.DetID, tTopo, true)
                        << " s_ladder: " << SiPixelPI::signed_ladder(mod.DetID, tTopo, true)
                        << " error type:" << mod.errorType << " BadRocs: " << bad_rocs.to_string('O', 'X') << std::endl;
            }
          }
        }  // payload
      }  // iovs
      return true;
    }  // fill
  };

  /************************************************
    summary class
  *************************************************/

  class SiPixelQualityBadRocsSummary : public PlotImage<SiPixelQuality, MULTI_IOV, 1> {
  public:
    SiPixelQualityBadRocsSummary() : PlotImage<SiPixelQuality, MULTI_IOV, 1>("SiPixel Quality Summary") {}

    bool fill() override {
      auto tag = PlotBase::getTag<0>();
      for (const auto& iov : tag.iovs) {
        std::shared_ptr<SiPixelQuality> payload = fetchPayload(std::get<1>(iov));
        auto unpacked = SiPixelPI::unpack(std::get<0>(iov));

        COUT << "======================= " << unpacked.first << " : " << unpacked.second << std::endl;
        auto theDisabledModules = payload->getBadComponentList();
        for (const auto& mod : theDisabledModules) {
          COUT << "detId: " << mod.DetID << " |error type: " << mod.errorType << " |BadRocs: " << mod.BadRocs
               << std::endl;
        }
      }

      //=========================
      TCanvas canv("Partition summary", "partition summary", 1200, 1000);
      canv.SetBottomMargin(0.11);
      canv.SetLeftMargin(0.13);
      canv.SetRightMargin(0.05);
      canv.cd();
      SiPixelPI::displayNotSupported(canv, 0);

      std::string fileName(m_imageFileName);
      canv.SaveAs(fileName.c_str());

      return true;
    }
  };

  /************************************************
    time history class
  *************************************************/

  class SiPixelQualityBadRocsTimeHistory : public TimeHistoryPlot<SiPixelQuality, std::pair<double, double> > {
  public:
    SiPixelQualityBadRocsTimeHistory()
        : TimeHistoryPlot<SiPixelQuality, std::pair<double, double> >("bad ROCs count vs time", "bad ROCs count") {}

    std::pair<double, double> getFromPayload(SiPixelQuality& payload) override {
      return std::make_pair(extractBadRocCount(payload), 0.);
    }

    unsigned int extractBadRocCount(SiPixelQuality& payload) {
      unsigned int BadRocCount(0);
      auto theDisabledModules = payload.getBadComponentList();
      for (const auto& mod : theDisabledModules) {
        for (unsigned short n = 0; n < 16; n++) {
          unsigned short mask = 1 << n;  // 1 << n = 2^{n} using bitwise shift
          if (mod.BadRocs & mask)
            BadRocCount++;
        }
      }
      return BadRocCount;
    }
  };

  /************************************************
   occupancy style map whole Pixel
  *************************************************/
  template <SiPixelPI::DetType myType>
  class SiPixelQualityMap : public PlotImage<SiPixelQuality, SINGLE_IOV> {
  public:
    SiPixelQualityMap()
        : PlotImage<SiPixelQuality, SINGLE_IOV>("SiPixelQuality Pixel Map"),
          m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
              edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {}

    bool fill() override {
      auto tag = PlotBase::getTag<0>();
      auto iov = tag.iovs.front();
      auto tagname = tag.name;
      std::shared_ptr<SiPixelQuality> payload = fetchPayload(std::get<1>(iov));

      Phase1PixelROCMaps theMap("");

      auto theDisabledModules = payload->getBadComponentList();
      if (this->isPhase0(theDisabledModules)) {
        edm::LogError("SiPixelQuality_PayloadInspector")
            << "SiPixelQuality maps are not supported for non-Phase1 Pixel geometries !";
        TCanvas canvas("Canv", "Canv", 1200, 1000);
        SiPixelPI::displayNotSupported(canvas, 0);
        std::string fileName(m_imageFileName);
        canvas.SaveAs(fileName.c_str());
        return false;
      }

      for (const auto& mod : theDisabledModules) {
        int subid = DetId(mod.DetID).subdetId();

        if ((subid == PixelSubdetector::PixelBarrel && myType == SiPixelPI::t_barrel) ||
            (subid == PixelSubdetector::PixelEndcap && myType == SiPixelPI::t_forward) ||
            (myType == SiPixelPI::t_all)) {
          std::bitset<16> bad_rocs(mod.BadRocs);
          if (payload->IsModuleBad(mod.DetID)) {
            theMap.fillWholeModule(mod.DetID, 1.);
          } else {
            theMap.fillSelectedRocs(mod.DetID, bad_rocs, 1.);
          }
        }
      }

      gStyle->SetOptStat(0);
      //=========================
      TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
      canvas.cd();

      auto unpacked = SiPixelPI::unpack(std::get<0>(iov));

      std::string IOVstring = (unpacked.first == 0)
                                  ? std::to_string(unpacked.second)
                                  : (std::to_string(unpacked.first) + "," + std::to_string(unpacked.second));

      const auto headerText = fmt::sprintf("#color[4]{%s},  IOV: #color[4]{%s}", tagname, IOVstring);

      switch (myType) {
        case SiPixelPI::t_barrel:
          theMap.drawBarrelMaps(canvas, headerText);
          break;
        case SiPixelPI::t_forward:
          theMap.drawForwardMaps(canvas, headerText);
          break;
        case SiPixelPI::t_all:
          theMap.drawMaps(canvas, headerText);
          break;
        default:
          throw cms::Exception("SiPixelQualityMap") << "\nERROR: unrecognized Pixel Detector part " << std::endl;
      }

      std::string fileName(m_imageFileName);
      canvas.SaveAs(fileName.c_str());
#ifdef MMDEBUG
      canvas.SaveAs("outAll.root");
#endif

      return true;
    }

  private:
    static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
    TrackerTopology m_trackerTopo;

    //_________________________________________________
    bool isPhase0(std::vector<SiPixelQuality::disabledModuleType> mods) {
      SiPixelDetInfoFileReader reader =
          SiPixelDetInfoFileReader(edm::FileInPath(SiPixelDetInfoFileReader::kPh0DefaultFile).fullPath());
      const auto& p0detIds = reader.getAllDetIds();

      std::vector<uint32_t> ownDetIds;
      std::transform(mods.begin(),
                     mods.end(),
                     std::back_inserter(ownDetIds),
                     [](SiPixelQuality::disabledModuleType d) -> uint32_t { return d.DetID; });

      for (const auto& det : ownDetIds) {
        // if found at least one phase-0 detId early return
        if (std::find(p0detIds.begin(), p0detIds.end(), det) != p0detIds.end()) {
          return true;
        }
      }
      return false;
    }
  };

  using SiPixelBPixQualityMap = SiPixelQualityMap<SiPixelPI::t_barrel>;
  using SiPixelFPixQualityMap = SiPixelQualityMap<SiPixelPI::t_forward>;
  using SiPixelFullQualityMap = SiPixelQualityMap<SiPixelPI::t_all>;

  /************************************************
   occupancy style map whole Pixel, difference of payloads
  *************************************************/
  template <SiPixelPI::DetType myType, IOVMultiplicity nIOVs, int ntags>
  class SiPixelQualityMapComparisonBase : public PlotImage<SiPixelQuality, nIOVs, ntags> {
  public:
    SiPixelQualityMapComparisonBase()
        : PlotImage<SiPixelQuality, nIOVs, ntags>(
              Form("SiPixelQuality %s Pixel Map", SiPixelPI::DetNames[myType].c_str())),
          m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
              edm::FileInPath("Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml").fullPath())} {}

    bool fill() override {
      // trick to deal with the multi-ioved tag and two tag case at the same time
      auto theIOVs = PlotBase::getTag<0>().iovs;
      auto f_tagname = PlotBase::getTag<0>().name;
      std::string l_tagname = "";
      auto firstiov = theIOVs.front();
      std::tuple<cond::Time_t, cond::Hash> lastiov;

      // we don't support (yet) comparison with more than 2 tags
      assert(this->m_plotAnnotations.ntags < 3);

      if (this->m_plotAnnotations.ntags == 2) {
        auto tag2iovs = PlotBase::getTag<1>().iovs;
        l_tagname = PlotBase::getTag<1>().name;
        lastiov = tag2iovs.front();
      } else {
        lastiov = theIOVs.back();
      }

      std::shared_ptr<SiPixelQuality> last_payload = this->fetchPayload(std::get<1>(lastiov));
      std::shared_ptr<SiPixelQuality> first_payload = this->fetchPayload(std::get<1>(firstiov));

      if (this->isPhase0(first_payload) || this->isPhase0(last_payload)) {
        edm::LogError("SiPixelQuality_PayloadInspector")
            << "SiPixelQuality comparison maps are not supported for non-Phase1 Pixel geometries !";
        TCanvas canvas("Canv", "Canv", 1200, 1000);
        SiPixelPI::displayNotSupported(canvas, 0);
        std::string fileName(this->m_imageFileName);
        canvas.SaveAs(fileName.c_str());
        return false;
      }

      Phase1PixelROCMaps theMap("", "#Delta payload A - payload B");

      gStyle->SetOptStat(0);
      //=========================
      TCanvas canvas("Summary", "Summary", 1200, k_height[myType]);
      canvas.cd();

      auto f_unpacked = SiPixelPI::unpack(std::get<0>(firstiov));
      auto l_unpacked = SiPixelPI::unpack(std::get<0>(lastiov));

      std::string f_IOVstring = (f_unpacked.first == 0)
                                    ? std::to_string(f_unpacked.second)
                                    : (std::to_string(f_unpacked.first) + "," + std::to_string(f_unpacked.second));

      std::string l_IOVstring = (l_unpacked.first == 0)
                                    ? std::to_string(l_unpacked.second)
                                    : (std::to_string(l_unpacked.first) + "," + std::to_string(l_unpacked.second));

      std::string headerText;

      if (this->m_plotAnnotations.ntags == 2) {
        headerText = fmt::sprintf(
            "#Delta #color[2]{A: %s, %s} - #color[4]{B: %s, %s}", f_tagname, f_IOVstring, l_tagname, l_IOVstring);
      } else {
        headerText =
            fmt::sprintf("%s, #Delta IOV #color[2]{A: %s} - #color[4]{B: %s} ", f_tagname, f_IOVstring, l_IOVstring);
      }

      switch (myType) {
        case SiPixelPI::t_barrel:
          theMap.drawBarrelMaps(canvas, headerText);
          break;
        case SiPixelPI::t_forward:
          theMap.drawForwardMaps(canvas, headerText);
          break;
        case SiPixelPI::t_all:
          theMap.drawMaps(canvas, headerText);
          break;
        default:
          throw cms::Exception("SiPixelQualityMapComparison")
              << "\nERROR: unrecognized Pixel Detector part " << std::endl;
      }

      // first loop on the first payload (newest)
      fillTheMapFromPayload(theMap, first_payload, false);

      // then loop on the second payload (oldest)
      fillTheMapFromPayload(theMap, last_payload, true);  // true will subtract

      std::string fileName(this->m_imageFileName);
      canvas.SaveAs(fileName.c_str());
#ifdef MMDEBUG
      canvas.SaveAs("outAll.root");
#endif

      return true;
    }

  private:
    static constexpr std::array<int, 3> k_height = {{1200, 600, 1600}};
    TrackerTopology m_trackerTopo;

    //_________________________________________________
    bool isPhase0(const std::shared_ptr<SiPixelQuality>& payload) {
      const auto mods = payload->getBadComponentList();
      SiPixelDetInfoFileReader reader =
          SiPixelDetInfoFileReader(edm::FileInPath(SiPixelDetInfoFileReader::kPh0DefaultFile).fullPath());
      const auto& p0detIds = reader.getAllDetIds();

      std::vector<uint32_t> ownDetIds;
      std::transform(mods.begin(),
                     mods.end(),
                     std::back_inserter(ownDetIds),
                     [](SiPixelQuality::disabledModuleType d) -> uint32_t { return d.DetID; });

      for (const auto& det : ownDetIds) {
        // if found at least one phase-0 detId early return
        if (std::find(p0detIds.begin(), p0detIds.end(), det) != p0detIds.end()) {
          return true;
        }
      }
      return false;
    }

    //____________________________________________________________________________________________
    void fillTheMapFromPayload(Phase1PixelROCMaps& theMap,
                               const std::shared_ptr<SiPixelQuality>& payload,
                               bool subtract) {
      const auto theDisabledModules = payload->getBadComponentList();
      for (const auto& mod : theDisabledModules) {
        int subid = DetId(mod.DetID).subdetId();
        if ((subid == PixelSubdetector::PixelBarrel && myType == SiPixelPI::t_barrel) ||
            (subid == PixelSubdetector::PixelEndcap && myType == SiPixelPI::t_forward) ||
            (myType == SiPixelPI::t_all)) {
          std::bitset<16> bad_rocs(mod.BadRocs);
          if (payload->IsModuleBad(mod.DetID)) {
            theMap.fillWholeModule(mod.DetID, (subtract ? -1. : 1.));
          } else {
            theMap.fillSelectedRocs(mod.DetID, bad_rocs, (subtract ? -1. : 1.));
          }
        }
      }
    }
  };

  using SiPixelBPixQualityMapCompareSingleTag = SiPixelQualityMapComparisonBase<SiPixelPI::t_barrel, MULTI_IOV, 1>;
  using SiPixelFPixQualityMapCompareSingleTag = SiPixelQualityMapComparisonBase<SiPixelPI::t_forward, MULTI_IOV, 1>;
  using SiPixelFullQualityMapCompareSingleTag = SiPixelQualityMapComparisonBase<SiPixelPI::t_all, MULTI_IOV, 1>;
  using SiPixelBPixQualityMapCompareTwoTags = SiPixelQualityMapComparisonBase<SiPixelPI::t_barrel, SINGLE_IOV, 2>;
  using SiPixelFPixQualityMapCompareTwoTags = SiPixelQualityMapComparisonBase<SiPixelPI::t_forward, SINGLE_IOV, 2>;
  using SiPixelFullQualityMapCompareTwoTags = SiPixelQualityMapComparisonBase<SiPixelPI::t_all, SINGLE_IOV, 2>;

  /************************************************
   Full Pixel Tracker Map class
  *************************************************/
  class SiPixelQualityBadFractionMap : public PlotImage<SiPixelQuality, SINGLE_IOV> {
  public:
    SiPixelQualityBadFractionMap() : PlotImage<SiPixelQuality, SINGLE_IOV>("SiPixelQuality Map") {
      label_ = "SiPixelQualityFullPixelMap";
      payloadString = "Quality";
    }

    bool fill() override {
      gStyle->SetPalette(1);
      auto tag = PlotBase::getTag<0>();
      auto iov = tag.iovs.front();
      auto unpacked = SiPixelPI::unpack(std::get<0>(iov));
      std::shared_ptr<SiPixelQuality> payload = this->fetchPayload(std::get<1>(iov));

      if (payload.get()) {
        Phase1PixelSummaryMap fullMap(
            "", fmt::sprintf("%s", payloadString), fmt::sprintf("bad %s fraction [%%]", payloadString));
        fullMap.createTrackerBaseMap();

        auto theDisabledModules = payload->getBadComponentList();
        if (this->isPhase0(theDisabledModules)) {
          edm::LogError("SiPixelQuality_PayloadInspector")
              << "SiPixelQuality maps are not supported for non-Phase1 Pixel geometries !";
          TCanvas canvas("Canv", "Canv", 1200, 1000);
          SiPixelPI::displayNotSupported(canvas, 0);
          std::string fileName(m_imageFileName);
          canvas.SaveAs(fileName.c_str());
          return false;
        }

        for (const auto& mod : theDisabledModules) {
          std::bitset<16> bad_rocs(mod.BadRocs);
          fullMap.fillTrackerMap(mod.DetID, (bad_rocs.count() / 16.) * 100);  // 16 ROCs in a Pixel Phase-1 module!
        }

        TCanvas canvas("Canv", "Canv", 3000, 2000);
        fullMap.printTrackerMap(canvas);

        auto ltx = TLatex();
        ltx.SetTextFont(62);
        ltx.SetTextSize(0.025);
        ltx.SetTextAlign(11);
        ltx.DrawLatexNDC(gPad->GetLeftMargin() + 0.01,
                         gPad->GetBottomMargin() + 0.01,
                         ("#color[4]{" + tag.name + "}, IOV: #color[4]{" + std::to_string(unpacked.first) + "," +
                          std::to_string(unpacked.second) + " }")
                             .c_str());

        std::string fileName(this->m_imageFileName);
        canvas.SaveAs(fileName.c_str());
      }
      return true;
    }

  protected:
    std::string payloadString;
    std::string label_;

  private:
    //_________________________________________________
    bool isPhase0(std::vector<SiPixelQuality::disabledModuleType> mods) {
      SiPixelDetInfoFileReader reader =
          SiPixelDetInfoFileReader(edm::FileInPath(SiPixelDetInfoFileReader::kPh0DefaultFile).fullPath());
      const auto& p0detIds = reader.getAllDetIds();

      std::vector<uint32_t> ownDetIds;
      std::transform(mods.begin(),
                     mods.end(),
                     std::back_inserter(ownDetIds),
                     [](SiPixelQuality::disabledModuleType d) -> uint32_t { return d.DetID; });

      for (const auto& det : ownDetIds) {
        // if found at least one phase-0 detId early return
        if (std::find(p0detIds.begin(), p0detIds.end(), det) != p0detIds.end()) {
          return true;
        }
      }
      return false;
    }
  };

  /************************************************
   Full Pixel Tracker Map class difference of two IOVs
  *************************************************/
  template <IOVMultiplicity nIOVs, int ntags>
  class SiPixelQualityBadFractionComparisonBase : public PlotImage<SiPixelQuality, nIOVs, ntags> {
  public:
    SiPixelQualityBadFractionComparisonBase() : PlotImage<SiPixelQuality, nIOVs, ntags>("SiPixelQuality Map") {
      label_ = "SiPixelQualityFullPixelMap";
      payloadString = "Quality";
    }

    bool fill() override {
      // Handle multi-IOV tag and two-tag case
      auto iovs = PlotBase::getTag<0>().iovs;
      std::string firstTagName = PlotBase::getTag<0>().name;
      std::string lastTagName = "";

      auto firstIOV = iovs.front();
      std::tuple<cond::Time_t, cond::Hash> lastIOV;

      // Ensure no more than two tags are supported
      assert(this->m_plotAnnotations.ntags < 3);

      if (this->m_plotAnnotations.ntags == 2) {
        auto secondTagIOVs = PlotBase::getTag<1>().iovs;
        lastTagName = PlotBase::getTag<1>().name;
        lastIOV = secondTagIOVs.front();
      } else {
        lastIOV = iovs.back();
      }

      // Fetch payloads for the first and last IOVs
      auto firstPayload = this->fetchPayload(std::get<1>(firstIOV));
      auto lastPayload = this->fetchPayload(std::get<1>(lastIOV));

      // Early exit if payloads are invalid
      if (!firstPayload || !lastPayload) {
        return false;
      }

      // Create summary map and base tracker map
      Phase1PixelSummaryMap summaryMap(
          "", fmt::sprintf("%s", payloadString), fmt::sprintf("bad %s fraction difference [%%]", payloadString));
      summaryMap.createTrackerBaseMap();

      // Get disabled modules for both IOVs
      auto firstDisabledModules = firstPayload->getBadComponentList();
      auto lastDisabledModules = lastPayload->getBadComponentList();

      // Check if geometry is supported (Phase 1)
      if (this->isPhase0(firstDisabledModules) || this->isPhase0(lastDisabledModules)) {
        edm::LogError("SiPixelQuality_PayloadInspector")
            << "SiPixelQuality maps are not supported for non-Phase1 Pixel geometries!";
        TCanvas canvas("Canv", "Canv", 1200, 1000);
        SiPixelPI::displayNotSupported(canvas, 0);
        std::string fileName(this->m_imageFileName);
        canvas.SaveAs(fileName.c_str());
        return false;
      }

      // Fill tracker map with bad ROC fractions for both IOVs
      const int NumROCsPerModule = 16;
      for (const auto& mod : firstDisabledModules) {
        std::bitset<NumROCsPerModule> badRocs(mod.BadRocs);
        summaryMap.fillTrackerMap(mod.DetID, (badRocs.count() / double(NumROCsPerModule)) * 100);
      }

      for (const auto& mod : lastDisabledModules) {
        std::bitset<NumROCsPerModule> badRocs(mod.BadRocs);
        summaryMap.fillTrackerMap(mod.DetID, -(badRocs.count() / double(NumROCsPerModule)) * 100);
      }

      // Apply custom color palette
      applyCustomPalette(summaryMap);

      // Draw canvas
      TCanvas canvas("Canv", "Canv", 3000, 2000);
      summaryMap.printTrackerMap(canvas);

      // Add header text with IOV information
      auto firstUnpacked = SiPixelPI::unpack(std::get<0>(firstIOV));
      auto lastUnpacked = SiPixelPI::unpack(std::get<0>(lastIOV));

      std::string firstIOVString = formatIOVString(firstUnpacked);
      std::string lastIOVString = formatIOVString(lastUnpacked);

      std::string headerText =
          formatHeaderText(this->m_plotAnnotations.ntags, firstTagName, firstIOVString, lastTagName, lastIOVString);

      TLatex ltx;
      ltx.SetTextFont(62);
      ltx.SetTextSize(0.025);
      ltx.SetTextAlign(11);
      ltx.DrawLatexNDC(gPad->GetLeftMargin() + 0.01, gPad->GetBottomMargin() + 0.01, headerText.c_str());

      // Save canvas
      std::string fileName(this->m_imageFileName);
      canvas.SaveAs(fileName.c_str());

      return true;
    }

  protected:
    std::string payloadString;
    std::string label_;

  private:
    //_________________________________________________
    bool isPhase0(std::vector<SiPixelQuality::disabledModuleType> mods) {
      SiPixelDetInfoFileReader reader =
          SiPixelDetInfoFileReader(edm::FileInPath(SiPixelDetInfoFileReader::kPh0DefaultFile).fullPath());
      const auto& p0detIds = reader.getAllDetIds();

      std::vector<uint32_t> ownDetIds;
      std::transform(mods.begin(),
                     mods.end(),
                     std::back_inserter(ownDetIds),
                     [](SiPixelQuality::disabledModuleType d) -> uint32_t { return d.DetID; });

      for (const auto& det : ownDetIds) {
        // if found at least one phase-0 detId early return
        if (std::find(p0detIds.begin(), p0detIds.end(), det) != p0detIds.end()) {
          return true;
        }
      }
      return false;
    }

    // Helper function to apply custom palette
    //_________________________________________________
    void applyCustomPalette(Phase1PixelSummaryMap& map) {
      double maxVal = map.getZAxisRange().second;
      double minVal = map.getZAxisRange().first;
      double whiteVal = 0.;
      double whitePos = (maxVal != minVal) ? ((whiteVal - minVal) / (maxVal - minVal)) : 0.5;

      const int numColors = 3;
      double red[numColors] = {0., 1., 1.};
      double green[numColors] = {0., 1., 0.};
      double blue[numColors] = {1., 1., 0.};
      double stops[numColors] = {0., whitePos, 1.};

      int numBins = 256;
      TColor::CreateGradientColorTable(numColors, stops, red, green, blue, numBins);
    }

    // Helper function to format IOV string
    //_________________________________________________
    std::string formatIOVString(const std::pair<int, int>& unpacked) {
      return (unpacked.first == 0) ? std::to_string(unpacked.second)
                                   : (std::to_string(unpacked.first) + "," + std::to_string(unpacked.second));
    }

    // Helper function to format header text
    //_________________________________________________
    std::string formatHeaderText(int nTAGs,
                                 const std::string& firstTagName,
                                 const std::string& firstIOVString,
                                 const std::string& lastTagName,
                                 const std::string& lastIOVString) {
      if (nTAGs == 2) {
        return fmt::sprintf("#Delta #color[2]{A: %s, %s} - #color[4]{B: %s, %s}",
                            firstTagName,
                            firstIOVString,
                            lastTagName,
                            lastIOVString);
      } else {
        return fmt::sprintf(
            "%s, #Delta IOV #color[2]{A: %s} - #color[4]{B: %s}", firstTagName, firstIOVString, lastIOVString);
      }
    }
  };

  using SiPixelQualityBadFracCompareSingleTag = SiPixelQualityBadFractionComparisonBase<MULTI_IOV, 1>;
  using SiPixelQualityBadFracCompareTwoTags = SiPixelQualityBadFractionComparisonBase<SINGLE_IOV, 2>;

}  // namespace

// Register the classes as boost python plugin
PAYLOAD_INSPECTOR_MODULE(SiPixelQuality) {
  PAYLOAD_INSPECTOR_CLASS(SiPixelQualityTest);
  PAYLOAD_INSPECTOR_CLASS(SiPixelQualityBadRocsSummary);
  PAYLOAD_INSPECTOR_CLASS(SiPixelQualityBadRocsTimeHistory);
  //PAYLOAD_INSPECTOR_CLASS(SiPixelQualityDebugger);
  PAYLOAD_INSPECTOR_CLASS(SiPixelBPixQualityMap);
  PAYLOAD_INSPECTOR_CLASS(SiPixelFPixQualityMap);
  PAYLOAD_INSPECTOR_CLASS(SiPixelFullQualityMap);
  PAYLOAD_INSPECTOR_CLASS(SiPixelQualityBadFractionMap);
  PAYLOAD_INSPECTOR_CLASS(SiPixelQualityBadFracCompareSingleTag);
  PAYLOAD_INSPECTOR_CLASS(SiPixelQualityBadFracCompareTwoTags);
  PAYLOAD_INSPECTOR_CLASS(SiPixelBPixQualityMapCompareSingleTag);
  PAYLOAD_INSPECTOR_CLASS(SiPixelFPixQualityMapCompareSingleTag);
  PAYLOAD_INSPECTOR_CLASS(SiPixelFullQualityMapCompareSingleTag);
  PAYLOAD_INSPECTOR_CLASS(SiPixelBPixQualityMapCompareTwoTags);
  PAYLOAD_INSPECTOR_CLASS(SiPixelFPixQualityMapCompareTwoTags);
  PAYLOAD_INSPECTOR_CLASS(SiPixelFullQualityMapCompareTwoTags);
}