Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:01:51

0001 /*!
0002   \file SiStripPedestals_PayloadInspector
0003   \Payload Inspector Plugin for SiStrip Noises
0004   \author M. Musich
0005   \version $Revision: 1.0 $
0006   \date $Date: 2017/09/22 11:02:00 $
0007 */
0008 
0009 #include "CondCore/Utilities/interface/PayloadInspectorModule.h"
0010 #include "CondCore/Utilities/interface/PayloadInspector.h"
0011 #include "CondCore/CondDB/interface/Time.h"
0012 
0013 // the data format of the condition to be inspected
0014 #include "CondFormats/SiStripObjects/interface/SiStripPedestals.h"
0015 #include "DataFormats/DetId/interface/DetId.h"
0016 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
0017 #include "CondFormats/SiStripObjects/interface/SiStripDetSummary.h"
0018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0019 
0020 // needed for the tracker map
0021 #include "CommonTools/TrackerMap/interface/TrackerMap.h"
0022 
0023 // auxilliary functions
0024 #include "CondCore/SiStripPlugins/interface/SiStripPayloadInspectorHelper.h"
0025 #include "CalibTracker/StandaloneTrackerTopology/interface/StandaloneTrackerTopology.h"
0026 #include "CalibTracker/SiStripCommon/interface/SiStripDetInfoFileReader.h"
0027 #include "FWCore/ParameterSet/interface/FileInPath.h"
0028 #include "SiStripCondObjectRepresent.h"
0029 
0030 #include <memory>
0031 #include <sstream>
0032 #include <iostream>
0033 #include <iomanip>
0034 #include <boost/tokenizer.hpp>
0035 
0036 // include ROOT
0037 #include "TH2F.h"
0038 #include "TLegend.h"
0039 #include "TCanvas.h"
0040 #include "TLine.h"
0041 #include "TStyle.h"
0042 #include "TLatex.h"
0043 #include "TPave.h"
0044 #include "TGaxis.h"
0045 #include "TPaveStats.h"
0046 
0047 namespace {
0048   using namespace cond::payloadInspector;
0049 
0050   class SiStripPedestalContainer : public SiStripCondObjectRepresent::SiStripDataContainer<SiStripPedestals, float> {
0051   public:
0052     SiStripPedestalContainer(const std::shared_ptr<SiStripPedestals>& payload,
0053                              const SiStripPI::MetaData& metadata,
0054                              const std::string& tagName)
0055         : SiStripCondObjectRepresent::SiStripDataContainer<SiStripPedestals, float>(payload, metadata, tagName) {
0056       payloadType_ = "SiStripPedestals";
0057       setGranularity(SiStripCondObjectRepresent::PERSTRIP);
0058     }
0059 
0060     void storeAllValues() override {
0061       std::vector<uint32_t> detid;
0062       payload_->getDetIds(detid);
0063 
0064       for (const auto& d : detid) {
0065         SiStripPedestals::Range range = payload_->getRange(d);
0066         for (int it = 0; it < (range.second - range.first) * 8 / 10; ++it) {
0067           // to be used to fill the histogram
0068           SiStripCondData_.fillByPushBack(d, payload_->getPed(it, range));
0069         }
0070       }
0071     }
0072   };
0073 
0074   class SiStripPedestalCompareByPartition : public PlotImage<SiStripPedestals, MULTI_IOV, 2> {
0075   public:
0076     SiStripPedestalCompareByPartition()
0077         : PlotImage<SiStripPedestals, MULTI_IOV, 2>("SiStrip Compare Pedestals By Partition") {}
0078 
0079     bool fill() override {
0080       // trick to deal with the multi-ioved tag and two tag case at the same time
0081       auto theIOVs = PlotBase::getTag<0>().iovs;
0082       auto tagname1 = PlotBase::getTag<0>().name;
0083       auto tag2iovs = PlotBase::getTag<1>().iovs;
0084       auto tagname2 = PlotBase::getTag<1>().name;
0085       SiStripPI::MetaData firstiov = theIOVs.front();
0086       SiStripPI::MetaData lastiov = tag2iovs.front();
0087 
0088       std::shared_ptr<SiStripPedestals> last_payload = fetchPayload(std::get<1>(lastiov));
0089       std::shared_ptr<SiStripPedestals> first_payload = fetchPayload(std::get<1>(firstiov));
0090 
0091       SiStripPedestalContainer* l_objContainer = new SiStripPedestalContainer(last_payload, lastiov, tagname1);
0092       SiStripPedestalContainer* f_objContainer = new SiStripPedestalContainer(first_payload, firstiov, tagname2);
0093 
0094       l_objContainer->compare(f_objContainer);
0095 
0096       //l_objContainer->printAll();
0097 
0098       TCanvas canvas("Partition summary", "partition summary", 1400, 1000);
0099       l_objContainer->fillByPartition(canvas, 300, 0.1, 300.);
0100 
0101       std::string fileName(m_imageFileName);
0102       canvas.SaveAs(fileName.c_str());
0103 
0104       return true;
0105     }  // fill
0106   };
0107 
0108   class SiStripPedestalDiffByPartition : public PlotImage<SiStripPedestals, MULTI_IOV, 2> {
0109   public:
0110     SiStripPedestalDiffByPartition()
0111         : PlotImage<SiStripPedestals, MULTI_IOV, 2>("SiStrip Diff Pedestals By Partition") {}
0112 
0113     bool fill() override {
0114       // trick to deal with the multi-ioved tag and two tag case at the same time
0115       auto theIOVs = PlotBase::getTag<0>().iovs;
0116       auto tagname1 = PlotBase::getTag<0>().name;
0117       auto tag2iovs = PlotBase::getTag<1>().iovs;
0118       auto tagname2 = PlotBase::getTag<1>().name;
0119       SiStripPI::MetaData firstiov = theIOVs.front();
0120       SiStripPI::MetaData lastiov = tag2iovs.front();
0121 
0122       std::shared_ptr<SiStripPedestals> last_payload = fetchPayload(std::get<1>(lastiov));
0123       std::shared_ptr<SiStripPedestals> first_payload = fetchPayload(std::get<1>(firstiov));
0124 
0125       SiStripPedestalContainer* l_objContainer = new SiStripPedestalContainer(last_payload, lastiov, tagname1);
0126       SiStripPedestalContainer* f_objContainer = new SiStripPedestalContainer(first_payload, firstiov, tagname2);
0127 
0128       l_objContainer->subtract(f_objContainer);
0129 
0130       //l_objContainer->printAll();
0131 
0132       TCanvas canvas("Partition summary", "partition summary", 1400, 1000);
0133       l_objContainer->fillByPartition(canvas, 100, -30., 30.);
0134 
0135       std::string fileName(m_imageFileName);
0136       canvas.SaveAs(fileName.c_str());
0137 
0138       return true;
0139     }  // fill
0140   };
0141 
0142   class SiStripPedestalCorrelationByPartition : public PlotImage<SiStripPedestals> {
0143   public:
0144     SiStripPedestalCorrelationByPartition()
0145         : PlotImage<SiStripPedestals>("SiStrip Pedestals Correlation By Partition") {
0146       setSingleIov(false);
0147     }
0148 
0149     bool fill(const std::vector<SiStripPI::MetaData>& iovs) override {
0150       std::vector<SiStripPI::MetaData> sorted_iovs = iovs;
0151 
0152       // make absolute sure the IOVs are sortd by since
0153       std::sort(begin(sorted_iovs), end(sorted_iovs), [](auto const& t1, auto const& t2) {
0154         return std::get<0>(t1) < std::get<0>(t2);
0155       });
0156 
0157       auto firstiov = sorted_iovs.front();
0158       auto lastiov = sorted_iovs.back();
0159 
0160       std::shared_ptr<SiStripPedestals> last_payload = fetchPayload(std::get<1>(lastiov));
0161       std::shared_ptr<SiStripPedestals> first_payload = fetchPayload(std::get<1>(firstiov));
0162 
0163       SiStripPedestalContainer* l_objContainer = new SiStripPedestalContainer(last_payload, lastiov, "");
0164       SiStripPedestalContainer* f_objContainer = new SiStripPedestalContainer(first_payload, firstiov, "");
0165 
0166       l_objContainer->compare(f_objContainer);
0167 
0168       TCanvas canvas("Partition summary", "partition summary", 1200, 1200);
0169       l_objContainer->fillCorrelationByPartition(canvas, 100, 0., 300.);
0170 
0171       std::string fileName(m_imageFileName);
0172       canvas.SaveAs(fileName.c_str());
0173 
0174       return true;
0175     }  // fill
0176   };
0177 
0178   /************************************************
0179     test class
0180   *************************************************/
0181 
0182   class SiStripPedestalsTest : public Histogram1D<SiStripPedestals, SINGLE_IOV> {
0183   public:
0184     SiStripPedestalsTest()
0185         : Histogram1D<SiStripPedestals, SINGLE_IOV>("SiStrip Pedestals test", "SiStrip Pedestals test", 10, 0.0, 10.0),
0186           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0187               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath())} {}
0188 
0189     bool fill() override {
0190       auto tag = PlotBase::getTag<0>();
0191       for (auto const& iov : tag.iovs) {
0192         std::shared_ptr<SiStripPedestals> payload = Base::fetchPayload(std::get<1>(iov));
0193         if (payload.get()) {
0194           fillWithValue(1.);
0195 
0196           std::stringstream ss;
0197           ss << "Summary of strips pedestals:" << std::endl;
0198 
0199           //payload->printDebug(ss);
0200           payload->printSummary(ss, &m_trackerTopo);
0201 
0202           std::vector<uint32_t> detid;
0203           payload->getDetIds(detid);
0204 
0205           // for (const auto & d : detid) {
0206           //   int nstrip=0;
0207           //   SiStripPedestals::Range range=payload->getRange(d);
0208           //   for( int it=0; it < (range.second-range.first)*8/10; ++it ){
0209           //     auto ped = payload->getPed(it,range);
0210           //     nstrip++;
0211           //     ss << "DetId="<< d << " Strip=" << nstrip <<": "<< ped << std::endl;
0212           //   } // end of loop on strips
0213           // } // end of loop on detIds
0214 
0215           std::cout << ss.str() << std::endl;
0216 
0217         }  // payload
0218       }    // iovs
0219       return true;
0220     }  // fill
0221   private:
0222     TrackerTopology m_trackerTopo;
0223   };
0224 
0225   /************************************************
0226     SiStrip Pedestals Profile of 1 IOV for one selected DetId
0227   *************************************************/
0228 
0229   class SiStripPedestalPerDetId : public PlotImage<SiStripPedestals, SINGLE_IOV> {
0230   public:
0231     SiStripPedestalPerDetId() : PlotImage<SiStripPedestals, SINGLE_IOV>("SiStrip Pedestal values per DetId") {
0232       PlotBase::addInputParam("DetIds");
0233     }
0234 
0235     bool fill() override {
0236       auto tag = PlotBase::getTag<0>();
0237       auto iov = tag.iovs.front();
0238       auto tagname = tag.name;
0239       std::shared_ptr<SiStripPedestals> payload = fetchPayload(std::get<1>(iov));
0240 
0241       std::vector<uint32_t> the_detids = {};
0242 
0243       auto paramValues = PlotBase::inputParamValues();
0244       auto ip = paramValues.find("DetIds");
0245       if (ip != paramValues.end()) {
0246         auto input = ip->second;
0247         typedef boost::tokenizer<boost::char_separator<char>> tokenizer;
0248         boost::char_separator<char> sep{","};
0249         tokenizer tok{input, sep};
0250         for (const auto& t : tok) {
0251           the_detids.push_back(atoi(t.c_str()));
0252         }
0253       } else {
0254         edm::LogWarning("SiStripNoisePerDetId")
0255             << "\n WARNING!!!! \n The needed parameter DetIds has not been passed. Will use all Strip DetIds! \n\n";
0256         the_detids.push_back(0xFFFFFFFF);
0257       }
0258 
0259       size_t ndets = the_detids.size();
0260       std::vector<std::shared_ptr<TH1F>> hpedestal;
0261       std::vector<std::shared_ptr<TLegend>> legends;
0262       std::vector<unsigned int> v_nAPVs;
0263       std::vector<std::vector<std::shared_ptr<TLine>>> lines;
0264       hpedestal.reserve(ndets);
0265       legends.reserve(ndets);
0266 
0267       auto sides = getClosestFactors(the_detids.size());
0268       edm::LogPrint("SiStripPedestalPerDetId") << "Aspect ratio: " << sides.first << ":" << sides.second << std::endl;
0269 
0270       if (payload.get()) {
0271         //=========================
0272         TCanvas canvas("ByDetId", "ByDetId", sides.second * 800, sides.first * 600);
0273         canvas.Divide(sides.second, sides.first);
0274         const auto detInfo =
0275             SiStripDetInfoFileReader::read(edm::FileInPath(SiStripDetInfoFileReader::kDefaultFile).fullPath());
0276         for (const auto& the_detid : the_detids) {
0277           edm::LogPrint("SiStripNoisePerDetId") << "DetId:" << the_detid << std::endl;
0278 
0279           unsigned int nAPVs = detInfo.getNumberOfApvsAndStripLength(the_detid).first;
0280           if (nAPVs == 0)
0281             nAPVs = 6;
0282           v_nAPVs.push_back(nAPVs);
0283 
0284           auto histo = std::make_shared<TH1F>(
0285               Form("Pedestal profile %s", std::to_string(the_detid).c_str()),
0286               Form("SiStrip Pedestal profile for DetId: %s;Strip number;SiStrip Pedestal [ADC counts]",
0287                    std::to_string(the_detid).c_str()),
0288               sistrip::STRIPS_PER_APV * nAPVs,
0289               -0.5,
0290               (sistrip::STRIPS_PER_APV * nAPVs) - 0.5);
0291 
0292           histo->SetStats(false);
0293           histo->SetTitle("");
0294 
0295           if (the_detid != 0xFFFFFFFF) {
0296             fillHisto(payload, histo, the_detid);
0297           } else {
0298             auto allDetIds = detInfo.getAllDetIds();
0299             for (const auto& id : allDetIds) {
0300               fillHisto(payload, histo, id);
0301             }
0302           }
0303 
0304           SiStripPI::makeNicePlotStyle(histo.get());
0305           histo->GetYaxis()->SetTitleOffset(1.0);
0306           hpedestal.push_back(histo);
0307         }  // loop on the detids
0308 
0309         for (size_t index = 0; index < ndets; index++) {
0310           canvas.cd(index + 1);
0311           canvas.cd(index + 1)->SetBottomMargin(0.11);
0312           canvas.cd(index + 1)->SetTopMargin(0.06);
0313           canvas.cd(index + 1)->SetLeftMargin(0.10);
0314           canvas.cd(index + 1)->SetRightMargin(0.02);
0315           hpedestal.at(index)->Draw();
0316           hpedestal.at(index)->GetYaxis()->SetRangeUser(0, hpedestal.at(index)->GetMaximum() * 1.2);
0317           canvas.cd(index)->Update();
0318 
0319           std::vector<int> boundaries;
0320           for (size_t b = 0; b < v_nAPVs.at(index); b++) {
0321             boundaries.push_back(b * sistrip::STRIPS_PER_APV);
0322           }
0323 
0324           std::vector<std::shared_ptr<TLine>> linesVec;
0325           for (const auto& bound : boundaries) {
0326             auto line = std::make_shared<TLine>(hpedestal.at(index)->GetBinLowEdge(bound),
0327                                                 canvas.cd(index + 1)->GetUymin(),
0328                                                 hpedestal.at(index)->GetBinLowEdge(bound),
0329                                                 canvas.cd(index + 1)->GetUymax());
0330             line->SetLineWidth(1);
0331             line->SetLineStyle(9);
0332             line->SetLineColor(2);
0333             linesVec.push_back(line);
0334           }
0335           lines.push_back(linesVec);
0336 
0337           for (const auto& line : lines.at(index)) {
0338             line->Draw("same");
0339           }
0340 
0341           canvas.cd(index + 1);
0342 
0343           auto ltx = TLatex();
0344           ltx.SetTextFont(62);
0345           ltx.SetTextSize(0.05);
0346           ltx.SetTextAlign(11);
0347           ltx.DrawLatexNDC(gPad->GetLeftMargin(),
0348                            1 - gPad->GetTopMargin() + 0.01,
0349                            Form("SiStrip Pedestals profile for DetId %s", std::to_string(the_detids[index]).c_str()));
0350 
0351           legends.push_back(std::make_shared<TLegend>(0.45, 0.83, 0.95, 0.93));
0352           legends.at(index)->SetHeader(tagname.c_str(), "C");  // option "C" allows to center the header
0353           legends.at(index)->AddEntry(
0354               hpedestal.at(index).get(), ("IOV: " + std::to_string(std::get<0>(iov))).c_str(), "PL");
0355           legends.at(index)->SetTextSize(0.045);
0356           legends.at(index)->Draw("same");
0357         }
0358 
0359         std::string fileName(m_imageFileName);
0360         canvas.SaveAs(fileName.c_str());
0361       }  // payload
0362       return true;
0363     }  // fill
0364 
0365   private:
0366     int nextPerfectSquare(int N) { return std::floor(sqrt(N)) + 1; }
0367 
0368     std::pair<int, int> getClosestFactors(int input) {
0369       if ((input % 2 != 0) && input > 1) {
0370         input += 1;
0371       }
0372 
0373       int testNum = (int)sqrt(input);
0374       while (input % testNum != 0) {
0375         testNum--;
0376       }
0377       return std::make_pair(testNum, input / testNum);
0378     }
0379 
0380     void fillHisto(const std::shared_ptr<SiStripPedestals> payload, std::shared_ptr<TH1F>& histo, uint32_t the_detid) {
0381       int nstrip = 0;
0382       SiStripPedestals::Range range = payload->getRange(the_detid);
0383       for (int it = 0; it < (range.second - range.first) * 8 / 10; ++it) {
0384         auto pedestal = payload->getPed(it, range);
0385         nstrip++;
0386         histo->AddBinContent(nstrip, pedestal);
0387       }  // end of loop on strips
0388     }
0389   };
0390 
0391   /************************************************
0392     1d histogram of SiStripPedestals of 1 IOV 
0393   *************************************************/
0394 
0395   // inherit from one of the predefined plot class: Histogram1D
0396   class SiStripPedestalsValue : public Histogram1D<SiStripPedestals, SINGLE_IOV> {
0397   public:
0398     SiStripPedestalsValue()
0399         : Histogram1D<SiStripPedestals, SINGLE_IOV>(
0400               "SiStrip Pedestals values", "SiStrip Pedestals values", 300, 0.0, 300.0) {}
0401 
0402     bool fill() override {
0403       auto tag = PlotBase::getTag<0>();
0404       for (auto const& iov : tag.iovs) {
0405         std::shared_ptr<SiStripPedestals> payload = Base::fetchPayload(std::get<1>(iov));
0406         if (payload.get()) {
0407           std::vector<uint32_t> detid;
0408           payload->getDetIds(detid);
0409 
0410           for (const auto& d : detid) {
0411             SiStripPedestals::Range range = payload->getRange(d);
0412             for (int it = 0; it < (range.second - range.first) * 8 / 10; ++it) {
0413               auto ped = payload->getPed(it, range);
0414               //to be used to fill the histogram
0415               fillWithValue(ped);
0416             }  // loop over APVs
0417           }    // loop over detIds
0418         }      // payload
0419       }        // iovs
0420       return true;
0421     }  // fill
0422   };
0423 
0424   /************************************************
0425     1d histogram of SiStripPedestals of 1 IOV per Detid
0426   *************************************************/
0427 
0428   // inherit from one of the predefined plot class: Histogram1D
0429   class SiStripPedestalsValuePerDetId : public Histogram1D<SiStripPedestals, SINGLE_IOV> {
0430   public:
0431     SiStripPedestalsValuePerDetId()
0432         : Histogram1D<SiStripPedestals, SINGLE_IOV>(
0433               "SiStrip Pedestal values per DetId", "SiStrip Pedestal values per DetId", 100, 0.0, 10.0) {
0434       PlotBase::addInputParam("DetId");
0435     }
0436 
0437     bool fill() override {
0438       auto tag = PlotBase::getTag<0>();
0439       for (auto const& iov : tag.iovs) {
0440         std::shared_ptr<SiStripPedestals> payload = Base::fetchPayload(std::get<1>(iov));
0441         unsigned int the_detid(0xFFFFFFFF);
0442         auto paramValues = PlotBase::inputParamValues();
0443         auto ip = paramValues.find("DetId");
0444         if (ip != paramValues.end()) {
0445           the_detid = std::stoul(ip->second);
0446         }
0447 
0448         if (payload.get()) {
0449           SiStripPedestals::Range range = payload->getRange(the_detid);
0450           for (int it = 0; it < (range.second - range.first) * 8 / 9; ++it) {
0451             auto noise = payload->getPed(it, range);
0452             //to be used to fill the histogram
0453             fillWithValue(noise);
0454           }  // loop over APVs
0455         }    // payload
0456       }      // iovs
0457       return true;
0458     }  // fill
0459   };
0460 
0461   /************************************************
0462     templated 1d histogram of SiStripPedestals of 1 IOV
0463   *************************************************/
0464 
0465   // inherit from one of the predefined plot class: PlotImage
0466   template <SiStripPI::OpMode op_mode_>
0467   class SiStripPedestalDistribution : public PlotImage<SiStripPedestals, SINGLE_IOV> {
0468   public:
0469     SiStripPedestalDistribution() : PlotImage<SiStripPedestals, SINGLE_IOV>("SiStrip Pedestal values") {}
0470 
0471     bool fill() override {
0472       auto tag = PlotBase::getTag<0>();
0473       auto iov = tag.iovs.front();
0474       TGaxis::SetMaxDigits(3);
0475       gStyle->SetOptStat("emr");
0476 
0477       std::shared_ptr<SiStripPedestals> payload = fetchPayload(std::get<1>(iov));
0478 
0479       auto mon1D = std::unique_ptr<SiStripPI::Monitor1D>(new SiStripPI::Monitor1D(
0480           op_mode_,
0481           "Pedestal",
0482           Form("#LT Pedestal #GT per %s for IOV [%s];#LTStrip Pedestal per %s#GT [ADC counts];n. %ss",
0483                opType(op_mode_).c_str(),
0484                std::to_string(std::get<0>(iov)).c_str(),
0485                opType(op_mode_).c_str(),
0486                opType(op_mode_).c_str()),
0487           300,
0488           0.,
0489           300.0));
0490 
0491       unsigned int prev_det = 0, prev_apv = 0;
0492       SiStripPI::Entry epedestal;
0493 
0494       std::vector<uint32_t> detids;
0495       payload->getDetIds(detids);
0496 
0497       // loop on payload
0498       for (const auto& d : detids) {
0499         SiStripPedestals::Range range = payload->getRange(d);
0500 
0501         unsigned int istrip = 0;
0502 
0503         for (int it = 0; it < (range.second - range.first) * 8 / 10; ++it) {
0504           auto pedestal = payload->getPed(it, range);
0505           bool flush = false;
0506           switch (op_mode_) {
0507             case (SiStripPI::APV_BASED):
0508               flush = (prev_det != 0 && prev_apv != istrip / sistrip::STRIPS_PER_APV);
0509               break;
0510             case (SiStripPI::MODULE_BASED):
0511               flush = (prev_det != 0 && prev_det != d);
0512               break;
0513             case (SiStripPI::STRIP_BASED):
0514               flush = (istrip != 0);
0515               break;
0516           }
0517 
0518           if (flush) {
0519             mon1D->Fill(prev_apv, prev_det, epedestal.mean());
0520             epedestal.reset();
0521           }
0522 
0523           epedestal.add(std::min<float>(pedestal, 300.));
0524           prev_apv = istrip / sistrip::STRIPS_PER_APV;
0525           istrip++;
0526         }
0527         prev_det = d;
0528       }
0529 
0530       //=========================
0531       TCanvas canvas("Partion summary", "partition summary", 1200, 1000);
0532       canvas.cd();
0533       canvas.SetBottomMargin(0.11);
0534       canvas.SetTopMargin(0.07);
0535       canvas.SetLeftMargin(0.13);
0536       canvas.SetRightMargin(0.05);
0537       canvas.Modified();
0538 
0539       auto hist = mon1D->getHist();
0540       SiStripPI::makeNicePlotStyle(&hist);
0541       hist.SetStats(kTRUE);
0542       hist.SetFillColorAlpha(kRed, 0.35);
0543       hist.Draw();
0544 
0545       canvas.Update();
0546 
0547       TPaveStats* st = (TPaveStats*)hist.GetListOfFunctions()->FindObject("stats");
0548       st->SetLineColor(kRed);
0549       st->SetTextColor(kRed);
0550       st->SetX1NDC(.75);
0551       st->SetX2NDC(.95);
0552       st->SetY1NDC(.83);
0553       st->SetY2NDC(.93);
0554 
0555       TLegend legend = TLegend(0.13, 0.83, 0.43, 0.93);
0556       legend.SetHeader(Form("SiStrip Pedestal values per %s", opType(op_mode_).c_str()),
0557                        "C");  // option "C" allows to center the header
0558       legend.AddEntry(&hist, ("IOV: " + std::to_string(std::get<0>(iov))).c_str(), "F");
0559       legend.SetTextSize(0.025);
0560       legend.Draw("same");
0561 
0562       std::string fileName(m_imageFileName);
0563       canvas.SaveAs(fileName.c_str());
0564 
0565       return true;
0566     }
0567 
0568     std::string opType(SiStripPI::OpMode mode) {
0569       std::string types[3] = {"Strip", "APV", "Module"};
0570       return types[mode];
0571     }
0572   };
0573 
0574   typedef SiStripPedestalDistribution<SiStripPI::STRIP_BASED> SiStripPedestalValuePerStrip;
0575   typedef SiStripPedestalDistribution<SiStripPI::APV_BASED> SiStripPedestalValuePerAPV;
0576   typedef SiStripPedestalDistribution<SiStripPI::MODULE_BASED> SiStripPedestalValuePerModule;
0577 
0578   /************************************************
0579   template 1d histogram comparison of SiStripPedestals of 1 IOV
0580   *************************************************/
0581 
0582   // inherit from one of the predefined plot class: PlotImage
0583 
0584   template <SiStripPI::OpMode op_mode_, int ntags, IOVMultiplicity nIOVs>
0585   class SiStripPedestalDistributionComparisonBase : public PlotImage<SiStripPedestals, nIOVs, ntags> {
0586   public:
0587     SiStripPedestalDistributionComparisonBase()
0588         : PlotImage<SiStripPedestals, nIOVs, ntags>("SiStrip Pedestal values comparison") {}
0589 
0590     bool fill() override {
0591       TGaxis::SetExponentOffset(-0.1, 0.01, "y");  // X and Y offset for Y axis
0592 
0593       // trick to deal with the multi-ioved tag and two tag case at the same time
0594       auto theIOVs = PlotBase::getTag<0>().iovs;
0595       auto tagname1 = PlotBase::getTag<0>().name;
0596       std::string tagname2 = "";
0597       auto firstiov = theIOVs.front();
0598       SiStripPI::MetaData lastiov;
0599 
0600       // we don't support (yet) comparison with more than 2 tags
0601       assert(this->m_plotAnnotations.ntags < 3);
0602 
0603       if (this->m_plotAnnotations.ntags == 2) {
0604         auto tag2iovs = PlotBase::getTag<1>().iovs;
0605         tagname2 = PlotBase::getTag<1>().name;
0606         lastiov = tag2iovs.front();
0607       } else {
0608         lastiov = theIOVs.back();
0609       }
0610 
0611       std::shared_ptr<SiStripPedestals> f_payload = this->fetchPayload(std::get<1>(firstiov));
0612       std::shared_ptr<SiStripPedestals> l_payload = this->fetchPayload(std::get<1>(lastiov));
0613 
0614       auto f_mon = std::unique_ptr<SiStripPI::Monitor1D>(new SiStripPI::Monitor1D(
0615           op_mode_,
0616           "f_Pedestal",
0617           Form(";#LTStrip Pedestal per %s#GT [ADC counts];n. %ss", opType(op_mode_).c_str(), opType(op_mode_).c_str()),
0618           300,
0619           0.,
0620           300.));
0621 
0622       auto l_mon = std::unique_ptr<SiStripPI::Monitor1D>(new SiStripPI::Monitor1D(
0623           op_mode_,
0624           "l_Pedestal",
0625           Form(";#LTStrip Pedestal per %s#GT [ADC counts];n. %ss", opType(op_mode_).c_str(), opType(op_mode_).c_str()),
0626           300,
0627           0.,
0628           300.));
0629 
0630       unsigned int prev_det = 0, prev_apv = 0;
0631       SiStripPI::Entry epedestal;
0632 
0633       std::vector<uint32_t> f_detid;
0634       f_payload->getDetIds(f_detid);
0635 
0636       // loop on first payload
0637       for (const auto& d : f_detid) {
0638         SiStripPedestals::Range range = f_payload->getRange(d);
0639 
0640         unsigned int istrip = 0;
0641         for (int it = 0; it < (range.second - range.first) * 8 / 10; ++it) {
0642           float pedestal = f_payload->getPed(it, range);
0643           //to be used to fill the histogram
0644 
0645           bool flush = false;
0646           switch (op_mode_) {
0647             case (SiStripPI::APV_BASED):
0648               flush = (prev_det != 0 && prev_apv != istrip / sistrip::STRIPS_PER_APV);
0649               break;
0650             case (SiStripPI::MODULE_BASED):
0651               flush = (prev_det != 0 && prev_det != d);
0652               break;
0653             case (SiStripPI::STRIP_BASED):
0654               flush = (istrip != 0);
0655               break;
0656           }
0657 
0658           if (flush) {
0659             f_mon->Fill(prev_apv, prev_det, epedestal.mean());
0660             epedestal.reset();
0661           }
0662           epedestal.add(std::min<float>(pedestal, 300.));
0663           prev_apv = istrip / sistrip::STRIPS_PER_APV;
0664           istrip++;
0665         }
0666         prev_det = d;
0667       }
0668 
0669       prev_det = 0, prev_apv = 0;
0670       epedestal.reset();
0671 
0672       std::vector<uint32_t> l_detid;
0673       l_payload->getDetIds(l_detid);
0674 
0675       // loop on first payload
0676       for (const auto& d : l_detid) {
0677         SiStripPedestals::Range range = l_payload->getRange(d);
0678 
0679         unsigned int istrip = 0;
0680         for (int it = 0; it < (range.second - range.first) * 8 / 10; ++it) {
0681           float pedestal = l_payload->getPed(it, range);
0682 
0683           bool flush = false;
0684           switch (op_mode_) {
0685             case (SiStripPI::APV_BASED):
0686               flush = (prev_det != 0 && prev_apv != istrip / sistrip::STRIPS_PER_APV);
0687               break;
0688             case (SiStripPI::MODULE_BASED):
0689               flush = (prev_det != 0 && prev_det != d);
0690               break;
0691             case (SiStripPI::STRIP_BASED):
0692               flush = (istrip != 0);
0693               break;
0694           }
0695 
0696           if (flush) {
0697             l_mon->Fill(prev_apv, prev_det, epedestal.mean());
0698             epedestal.reset();
0699           }
0700 
0701           epedestal.add(std::min<float>(pedestal, 300.));
0702           prev_apv = istrip / sistrip::STRIPS_PER_APV;
0703           istrip++;
0704         }
0705         prev_det = d;
0706       }
0707 
0708       auto h_first = f_mon->getHist();
0709       h_first.SetStats(kFALSE);
0710       auto h_last = l_mon->getHist();
0711       h_last.SetStats(kFALSE);
0712 
0713       SiStripPI::makeNicePlotStyle(&h_first);
0714       SiStripPI::makeNicePlotStyle(&h_last);
0715 
0716       h_first.GetYaxis()->CenterTitle(true);
0717       h_last.GetYaxis()->CenterTitle(true);
0718 
0719       h_first.GetXaxis()->CenterTitle(true);
0720       h_last.GetXaxis()->CenterTitle(true);
0721 
0722       h_first.SetLineWidth(2);
0723       h_last.SetLineWidth(2);
0724 
0725       h_first.SetLineColor(kBlack);
0726       h_last.SetLineColor(kBlue);
0727 
0728       //=========================
0729       TCanvas canvas("Partion summary", "partition summary", 1200, 1000);
0730       canvas.cd();
0731       canvas.SetTopMargin(0.06);
0732       canvas.SetBottomMargin(0.10);
0733       canvas.SetLeftMargin(0.13);
0734       canvas.SetRightMargin(0.05);
0735       canvas.Modified();
0736 
0737       float theMax = (h_first.GetMaximum() > h_last.GetMaximum()) ? h_first.GetMaximum() : h_last.GetMaximum();
0738 
0739       h_first.SetMaximum(theMax * 1.20);
0740       h_last.SetMaximum(theMax * 1.20);
0741 
0742       h_first.Draw();
0743       h_last.SetFillColorAlpha(kBlue, 0.15);
0744       h_last.Draw("same");
0745 
0746       TLegend legend = TLegend(0.13, 0.83, 0.95, 0.94);
0747       if (this->m_plotAnnotations.ntags == 2) {
0748         legend.SetHeader("#bf{Two Tags Comparison}", "C");  // option "C" allows to center the header
0749         legend.AddEntry(&h_first, (tagname1 + " : " + std::to_string(std::get<0>(firstiov))).c_str(), "F");
0750         legend.AddEntry(&h_last, (tagname2 + " : " + std::to_string(std::get<0>(lastiov))).c_str(), "F");
0751       } else {
0752         legend.SetHeader(("tag: #bf{" + tagname1 + "}").c_str(), "C");  // option "C" allows to center the header
0753         legend.AddEntry(&h_first, ("IOV since: " + std::to_string(std::get<0>(firstiov))).c_str(), "F");
0754         legend.AddEntry(&h_last, ("IOV since: " + std::to_string(std::get<0>(lastiov))).c_str(), "F");
0755       }
0756       legend.SetTextSize(0.025);
0757       legend.Draw("same");
0758 
0759       auto ltx = TLatex();
0760       ltx.SetTextFont(62);
0761       ltx.SetTextSize(0.05);
0762       ltx.SetTextAlign(11);
0763       ltx.DrawLatexNDC(gPad->GetLeftMargin(),
0764                        1 - gPad->GetTopMargin() + 0.01,
0765                        Form("#LTSiStrip Pedestals#GT Comparison per %s", opType(op_mode_).c_str()));
0766 
0767       std::string fileName(this->m_imageFileName);
0768       canvas.SaveAs(fileName.c_str());
0769 
0770       return true;
0771     }
0772 
0773     std::string opType(SiStripPI::OpMode mode) {
0774       std::string types[3] = {"Strip", "APV", "Module"};
0775       return types[mode];
0776     }
0777   };
0778 
0779   template <SiStripPI::OpMode op_mode_>
0780   using SiStripPedestalDistributionComparisonSingleTag =
0781       SiStripPedestalDistributionComparisonBase<op_mode_, 1, MULTI_IOV>;
0782 
0783   template <SiStripPI::OpMode op_mode_>
0784   using SiStripPedestalDistributionComparisonTwoTags =
0785       SiStripPedestalDistributionComparisonBase<op_mode_, 2, SINGLE_IOV>;
0786 
0787   typedef SiStripPedestalDistributionComparisonSingleTag<SiStripPI::STRIP_BASED>
0788       SiStripPedestalValueComparisonPerStripSingleTag;
0789   typedef SiStripPedestalDistributionComparisonSingleTag<SiStripPI::APV_BASED>
0790       SiStripPedestalValueComparisonPerAPVSingleTag;
0791   typedef SiStripPedestalDistributionComparisonSingleTag<SiStripPI::MODULE_BASED>
0792       SiStripPedestalValueComparisonPerModuleSingleTag;
0793 
0794   typedef SiStripPedestalDistributionComparisonTwoTags<SiStripPI::STRIP_BASED>
0795       SiStripPedestalValueComparisonPerStripTwoTags;
0796   typedef SiStripPedestalDistributionComparisonTwoTags<SiStripPI::APV_BASED> SiStripPedestalValueComparisonPerAPVTwoTags;
0797   typedef SiStripPedestalDistributionComparisonTwoTags<SiStripPI::MODULE_BASED>
0798       SiStripPedestalValueComparisonPerModuleTwoTags;
0799 
0800   /************************************************
0801     1d histogram of fraction of Zero SiStripPedestals of 1 IOV 
0802   *************************************************/
0803 
0804   // inherit from one of the predefined plot class PlotImage
0805   class SiStripZeroPedestalsFraction_TrackerMap : public PlotImage<SiStripPedestals, SINGLE_IOV> {
0806   public:
0807     SiStripZeroPedestalsFraction_TrackerMap()
0808         : PlotImage<SiStripPedestals, SINGLE_IOV>("Tracker Map of Zero SiStripPedestals fraction per module") {}
0809 
0810     bool fill() override {
0811       auto tag = PlotBase::getTag<0>();
0812       auto iov = tag.iovs.front();
0813       std::shared_ptr<SiStripPedestals> payload = fetchPayload(std::get<1>(iov));
0814 
0815       const auto detInfo =
0816           SiStripDetInfoFileReader::read(edm::FileInPath(SiStripDetInfoFileReader::kDefaultFile).fullPath());
0817 
0818       std::string titleMap =
0819           "Tracker Map of Zero SiStrip Pedestals fraction per module (payload : " + std::get<1>(iov) + ")";
0820 
0821       std::unique_ptr<TrackerMap> tmap = std::make_unique<TrackerMap>("SiStripPedestals");
0822       tmap->setTitle(titleMap);
0823       tmap->setPalette(1);
0824 
0825       std::vector<uint32_t> detid;
0826       payload->getDetIds(detid);
0827 
0828       std::map<uint32_t, int> zeropeds_per_detid;
0829 
0830       for (const auto& d : detid) {
0831         SiStripPedestals::Range range = payload->getRange(d);
0832         for (int it = 0; it < (range.second - range.first) * 8 / 10; ++it) {
0833           auto ped = payload->getPed(it, range);
0834           if (ped == 0.) {
0835             zeropeds_per_detid[d] += 1;
0836           }
0837         }  // end of loop on strips
0838         float fraction =
0839             zeropeds_per_detid[d] / (sistrip::STRIPS_PER_APV * detInfo.getNumberOfApvsAndStripLength(d).first);
0840         if (fraction > 0.) {
0841           tmap->fill(d, fraction);
0842           std::cout << "detid: " << d << " (n. APVs=" << detInfo.getNumberOfApvsAndStripLength(d).first << ") has "
0843                     << std::setw(4) << zeropeds_per_detid[d]
0844                     << " zero-pedestals strips (i.e. a fraction:" << std::setprecision(5) << fraction << ")"
0845                     << std::endl;
0846         }
0847       }  // end of loop on detids
0848 
0849       std::string fileName(m_imageFileName);
0850       tmap->save(true, 0., 0., fileName);
0851 
0852       return true;
0853     }
0854   };
0855 
0856   /************************************************
0857     Tracker Map of SiStrip Pedestals
0858   *************************************************/
0859 
0860   template <SiStripPI::estimator est>
0861   class SiStripPedestalsTrackerMap : public PlotImage<SiStripPedestals, SINGLE_IOV> {
0862   public:
0863     SiStripPedestalsTrackerMap()
0864         : PlotImage<SiStripPedestals, SINGLE_IOV>("Tracker Map of SiStripPedestals " + estimatorType(est) +
0865                                                   " per module") {}
0866 
0867     bool fill() override {
0868       auto tag = PlotBase::getTag<0>();
0869       auto iov = tag.iovs.front();
0870       std::shared_ptr<SiStripPedestals> payload = fetchPayload(std::get<1>(iov));
0871 
0872       std::string titleMap =
0873           "Tracker Map of SiStrip Pedestals " + estimatorType(est) + " per module (payload : " + std::get<1>(iov) + ")";
0874 
0875       std::unique_ptr<TrackerMap> tmap = std::make_unique<TrackerMap>("SiStripPedestals");
0876       tmap->setTitle(titleMap);
0877       tmap->setPalette(1);
0878 
0879       std::vector<uint32_t> detid;
0880       payload->getDetIds(detid);
0881 
0882       std::map<unsigned int, float> info_per_detid;
0883 
0884       for (const auto& d : detid) {
0885         int nstrips = 0;
0886         double mean(0.), rms(0.), min(10000.), max(0.);
0887         SiStripPedestals::Range range = payload->getRange(d);
0888         for (int it = 0; it < (range.second - range.first) * 8 / 10; ++it) {
0889           nstrips++;
0890           auto ped = payload->getPed(it, range);
0891           mean += ped;
0892           rms += (ped * ped);
0893           if (ped < min)
0894             min = ped;
0895           if (ped > max)
0896             max = ped;
0897         }  // end of loop on strips
0898 
0899         mean /= nstrips;
0900         if ((rms / nstrips - mean * mean) > 0.) {
0901           rms = sqrt(rms / nstrips - mean * mean);
0902         } else {
0903           rms = 0.;
0904         }
0905 
0906         switch (est) {
0907           case SiStripPI::min:
0908             info_per_detid[d] = min;
0909             break;
0910           case SiStripPI::max:
0911             info_per_detid[d] = max;
0912             break;
0913           case SiStripPI::mean:
0914             info_per_detid[d] = mean;
0915             break;
0916           case SiStripPI::rms:
0917             info_per_detid[d] = rms;
0918             break;
0919           default:
0920             edm::LogWarning("LogicError") << "Unknown estimator: " << est;
0921             break;
0922         }
0923       }  // end of loop on detids
0924 
0925       for (const auto& d : detid) {
0926         tmap->fill(d, info_per_detid[d]);
0927       }
0928 
0929       std::string fileName(m_imageFileName);
0930       tmap->save(true, 0., 0., fileName);
0931 
0932       return true;
0933     }
0934   };
0935 
0936   typedef SiStripPedestalsTrackerMap<SiStripPI::min> SiStripPedestalsMin_TrackerMap;
0937   typedef SiStripPedestalsTrackerMap<SiStripPI::max> SiStripPedestalsMax_TrackerMap;
0938   typedef SiStripPedestalsTrackerMap<SiStripPI::mean> SiStripPedestalsMean_TrackerMap;
0939   typedef SiStripPedestalsTrackerMap<SiStripPI::rms> SiStripPedestalsRMS_TrackerMap;
0940 
0941   /************************************************
0942     Tracker Map of SiStrip Pedestals Summaries
0943   *************************************************/
0944 
0945   template <SiStripPI::estimator est>
0946   class SiStripPedestalsByRegion : public PlotImage<SiStripPedestals, SINGLE_IOV> {
0947   public:
0948     SiStripPedestalsByRegion()
0949         : PlotImage<SiStripPedestals, SINGLE_IOV>("SiStrip Pedestals " + estimatorType(est) + " by Region"),
0950           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0951               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath())} {}
0952 
0953     bool fill() override {
0954       auto tag = PlotBase::getTag<0>();
0955       auto iov = tag.iovs.front();
0956       std::shared_ptr<SiStripPedestals> payload = fetchPayload(std::get<1>(iov));
0957 
0958       SiStripDetSummary summaryPedestals{&m_trackerTopo};
0959       std::vector<uint32_t> detid;
0960       payload->getDetIds(detid);
0961 
0962       for (const auto& d : detid) {
0963         int nstrips = 0;
0964         double mean(0.), rms(0.), min(10000.), max(0.);
0965         SiStripPedestals::Range range = payload->getRange(d);
0966         for (int it = 0; it < (range.second - range.first) * 8 / 10; ++it) {
0967           nstrips++;
0968           auto ped = payload->getPed(it, range);
0969           mean += ped;
0970           rms += (ped * ped);
0971           if (ped < min)
0972             min = ped;
0973           if (ped > max)
0974             max = ped;
0975         }  // end of loop on strips
0976 
0977         mean /= nstrips;
0978         if ((rms / nstrips - mean * mean) > 0.) {
0979           rms = sqrt(rms / nstrips - mean * mean);
0980         } else {
0981           rms = 0.;
0982         }
0983 
0984         switch (est) {
0985           case SiStripPI::min:
0986             summaryPedestals.add(d, min);
0987             break;
0988           case SiStripPI::max:
0989             summaryPedestals.add(d, max);
0990             break;
0991           case SiStripPI::mean:
0992             summaryPedestals.add(d, mean);
0993             break;
0994           case SiStripPI::rms:
0995             summaryPedestals.add(d, rms);
0996             break;
0997           default:
0998             edm::LogWarning("LogicError") << "Unknown estimator: " << est;
0999             break;
1000         }
1001       }  // loop on the detIds
1002 
1003       std::map<unsigned int, SiStripDetSummary::Values> map = summaryPedestals.getCounts();
1004       //=========================
1005 
1006       TCanvas canvas("Partion summary", "partition summary", 1200, 1000);
1007       canvas.cd();
1008       auto h1 = std::unique_ptr<TH1F>(new TH1F(
1009           "byRegion",
1010           Form("Average by partition of %s SiStrip Pedestals per module;;average SiStrip Pedestals %s [ADC counts]",
1011                estimatorType(est).c_str(),
1012                estimatorType(est).c_str()),
1013           map.size(),
1014           0.,
1015           map.size()));
1016       h1->SetStats(false);
1017       canvas.SetBottomMargin(0.18);
1018       canvas.SetLeftMargin(0.17);
1019       canvas.SetRightMargin(0.05);
1020       canvas.Modified();
1021 
1022       std::vector<int> boundaries;
1023       unsigned int iBin = 0;
1024 
1025       std::string detector;
1026       std::string currentDetector;
1027 
1028       for (const auto& element : map) {
1029         iBin++;
1030         int count = element.second.count;
1031         double mean = (element.second.mean) / count;
1032         double rms = (element.second.rms) / count - mean * mean;
1033 
1034         if (rms <= 0)
1035           rms = 0;
1036         else
1037           rms = sqrt(rms);
1038 
1039         if (currentDetector.empty())
1040           currentDetector = "TIB";
1041 
1042         switch ((element.first) / 1000) {
1043           case 1:
1044             detector = "TIB";
1045             break;
1046           case 2:
1047             detector = "TOB";
1048             break;
1049           case 3:
1050             detector = "TEC";
1051             break;
1052           case 4:
1053             detector = "TID";
1054             break;
1055         }
1056 
1057         h1->SetBinContent(iBin, mean);
1058         h1->GetXaxis()->SetBinLabel(iBin, SiStripPI::regionType(element.first).second);
1059         h1->GetXaxis()->LabelsOption("v");
1060 
1061         if (detector != currentDetector) {
1062           boundaries.push_back(iBin);
1063           currentDetector = detector;
1064         }
1065       }
1066 
1067       h1->SetMarkerStyle(20);
1068       h1->SetMarkerSize(1);
1069       h1->SetMaximum(h1->GetMaximum() * 1.1);
1070       h1->Draw("HIST");
1071       h1->Draw("Psame");
1072 
1073       canvas.Update();
1074 
1075       TLine l[boundaries.size()];
1076       unsigned int i = 0;
1077       for (const auto& line : boundaries) {
1078         l[i] = TLine(h1->GetBinLowEdge(line), canvas.GetUymin(), h1->GetBinLowEdge(line), canvas.GetUymax());
1079         l[i].SetLineWidth(1);
1080         l[i].SetLineStyle(9);
1081         l[i].SetLineColor(2);
1082         l[i].Draw("same");
1083         i++;
1084       }
1085 
1086       TLegend legend = TLegend(0.52, 0.82, 0.95, 0.9);
1087       legend.SetHeader((std::get<1>(iov)).c_str(), "C");  // option "C" allows to center the header
1088       legend.AddEntry(h1.get(), ("IOV: " + std::to_string(std::get<0>(iov))).c_str(), "PL");
1089       legend.SetTextSize(0.025);
1090       legend.Draw("same");
1091 
1092       std::string fileName(m_imageFileName);
1093       canvas.SaveAs(fileName.c_str());
1094 
1095       return true;
1096     }
1097 
1098   private:
1099     TrackerTopology m_trackerTopo;
1100   };
1101 
1102   typedef SiStripPedestalsByRegion<SiStripPI::mean> SiStripPedestalsMeanByRegion;
1103   typedef SiStripPedestalsByRegion<SiStripPI::min> SiStripPedestalsMinByRegion;
1104   typedef SiStripPedestalsByRegion<SiStripPI::max> SiStripPedestalsMaxByRegion;
1105   typedef SiStripPedestalsByRegion<SiStripPI::rms> SiStripPedestalsRMSByRegion;
1106 
1107 }  // namespace
1108 
1109 PAYLOAD_INSPECTOR_MODULE(SiStripPedestals) {
1110   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalCompareByPartition);
1111   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalDiffByPartition);
1112   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalCorrelationByPartition);
1113   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsTest);
1114   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalPerDetId);
1115   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsValue);
1116   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsValuePerDetId);
1117   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalValuePerStrip);
1118   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalValuePerAPV);
1119   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalValuePerModule);
1120   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalValueComparisonPerStripSingleTag);
1121   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalValueComparisonPerAPVSingleTag);
1122   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalValueComparisonPerModuleSingleTag);
1123   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalValueComparisonPerStripTwoTags);
1124   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalValueComparisonPerAPVTwoTags);
1125   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalValueComparisonPerModuleTwoTags);
1126   PAYLOAD_INSPECTOR_CLASS(SiStripZeroPedestalsFraction_TrackerMap);
1127   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsMin_TrackerMap);
1128   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsMax_TrackerMap);
1129   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsMean_TrackerMap);
1130   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsRMS_TrackerMap);
1131   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsMeanByRegion);
1132   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsMinByRegion);
1133   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsMaxByRegion);
1134   PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsRMSByRegion);
1135 }