Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-14 22:39:29

0001 /*!
0002   \file SiStripNoises_PayloadInspector
0003   \Payload Inspector Plugin for SiStrip Noises
0004   \author M. Musich
0005   \version $Revision: 1.0 $
0006   \date $Date: 2017/09/21 13:59:56 $
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/SiStripNoises.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 <boost/tokenizer.hpp>
0034 
0035 // include ROOT
0036 #include "TH2F.h"
0037 #include "TF1.h"
0038 #include "TGraphErrors.h"
0039 #include "TLegend.h"
0040 #include "TCanvas.h"
0041 #include "TLine.h"
0042 #include "TStyle.h"
0043 #include "TLatex.h"
0044 #include "TPave.h"
0045 #include "TPaveStats.h"
0046 #include "TGaxis.h"
0047 
0048 namespace {
0049 
0050   using namespace cond::payloadInspector;
0051 
0052   class SiStripNoiseContainer : public SiStripCondObjectRepresent::SiStripDataContainer<SiStripNoises, float> {
0053   public:
0054     SiStripNoiseContainer(const std::shared_ptr<SiStripNoises>& payload,
0055                           const SiStripPI::MetaData& metadata,
0056                           const std::string& tagName)
0057         : SiStripCondObjectRepresent::SiStripDataContainer<SiStripNoises, float>(payload, metadata, tagName) {
0058       payloadType_ = "SiStripNoises";
0059       setGranularity(SiStripCondObjectRepresent::PERSTRIP);
0060     }
0061 
0062     void storeAllValues() override {
0063       std::vector<uint32_t> detid;
0064       payload_->getDetIds(detid);
0065 
0066       for (const auto& d : detid) {
0067         SiStripNoises::Range range = payload_->getRange(d);
0068         for (int it = 0; it < (range.second - range.first) * 8 / 9; ++it) {
0069           // to be used to fill the histogram
0070           SiStripCondData_.fillByPushBack(d, payload_->getNoise(it, range));
0071         }
0072       }
0073     }
0074   };
0075 
0076   template <int ntags, IOVMultiplicity nIOVs>
0077   class SiStripNoiseCompareByPartitionBase : public PlotImage<SiStripNoises, nIOVs, ntags> {
0078   public:
0079     SiStripNoiseCompareByPartitionBase()
0080         : PlotImage<SiStripNoises, nIOVs, ntags>("SiStrip Compare Noises By Partition") {}
0081 
0082     bool fill() override {
0083       // trick to deal with the multi-ioved tag and two tag case at the same time
0084       auto theIOVs = PlotBase::getTag<0>().iovs;
0085       auto tagname1 = PlotBase::getTag<0>().name;
0086       std::string tagname2{tagname1};
0087       auto firstiov = theIOVs.front();
0088       SiStripPI::MetaData lastiov;
0089 
0090       // we don't support (yet) comparison with more than 2 tags
0091       assert(this->m_plotAnnotations.ntags < 3);
0092 
0093       if (this->m_plotAnnotations.ntags == 2) {
0094         auto tag2iovs = PlotBase::getTag<1>().iovs;
0095         tagname2 = PlotBase::getTag<1>().name;
0096         lastiov = tag2iovs.front();
0097       } else {
0098         lastiov = theIOVs.back();
0099       }
0100 
0101       std::shared_ptr<SiStripNoises> last_payload = this->fetchPayload(std::get<1>(lastiov));
0102       std::shared_ptr<SiStripNoises> first_payload = this->fetchPayload(std::get<1>(firstiov));
0103 
0104       SiStripNoiseContainer* l_objContainer = new SiStripNoiseContainer(last_payload, lastiov, tagname1);
0105       SiStripNoiseContainer* f_objContainer = new SiStripNoiseContainer(first_payload, firstiov, tagname2);
0106 
0107       l_objContainer->compare(f_objContainer);
0108 
0109       //l_objContainer->printAll();
0110 
0111       TCanvas canvas("Partition summary", "partition summary", 1400, 1000);
0112       l_objContainer->fillByPartition(canvas, 100, 0.1, 10.);
0113 
0114       std::string fileName(this->m_imageFileName);
0115       canvas.SaveAs(fileName.c_str());
0116 
0117       return true;
0118     }  // fill
0119   };
0120 
0121   using SiStripNoiseCompareByPartitionSingleTag = SiStripNoiseCompareByPartitionBase<1, MULTI_IOV>;
0122   using SiStripNoiseCompareByPartitionTwoTags = SiStripNoiseCompareByPartitionBase<2, SINGLE_IOV>;
0123 
0124   template <int ntags, IOVMultiplicity nIOVs>
0125   class SiStripNoiseDiffByPartitionBase : public PlotImage<SiStripNoises, nIOVs, ntags> {
0126   public:
0127     SiStripNoiseDiffByPartitionBase() : PlotImage<SiStripNoises, nIOVs, ntags>("SiStrip Diff Noises By Partition") {}
0128 
0129     bool fill() override {
0130       // trick to deal with the multi-ioved tag and two tag case at the same time
0131       auto theIOVs = PlotBase::getTag<0>().iovs;
0132       auto tagname1 = PlotBase::getTag<0>().name;
0133       std::string tagname2{tagname1};
0134       auto firstiov = theIOVs.front();
0135       SiStripPI::MetaData lastiov;
0136 
0137       // we don't support (yet) comparison with more than 2 tags
0138       assert(this->m_plotAnnotations.ntags < 3);
0139 
0140       if (this->m_plotAnnotations.ntags == 2) {
0141         auto tag2iovs = PlotBase::getTag<1>().iovs;
0142         tagname2 = PlotBase::getTag<1>().name;
0143         lastiov = tag2iovs.front();
0144       } else {
0145         lastiov = theIOVs.back();
0146       }
0147 
0148       std::shared_ptr<SiStripNoises> last_payload = this->fetchPayload(std::get<1>(lastiov));
0149       std::shared_ptr<SiStripNoises> first_payload = this->fetchPayload(std::get<1>(firstiov));
0150 
0151       SiStripNoiseContainer* l_objContainer = new SiStripNoiseContainer(last_payload, lastiov, tagname1);
0152       SiStripNoiseContainer* f_objContainer = new SiStripNoiseContainer(first_payload, firstiov, tagname2);
0153 
0154       l_objContainer->subtract(f_objContainer);
0155 
0156       //l_objContainer->printAll();
0157 
0158       TCanvas canvas("Partition summary", "partition summary", 1400, 1000);
0159       l_objContainer->fillByPartition(canvas, 100, -1., 1.);
0160 
0161       std::string fileName(this->m_imageFileName);
0162       canvas.SaveAs(fileName.c_str());
0163 
0164       return true;
0165     }  // fill
0166   };
0167 
0168   using SiStripNoiseDiffByPartitionSingleTag = SiStripNoiseDiffByPartitionBase<1, MULTI_IOV>;
0169   using SiStripNoiseDiffByPartitionTwoTags = SiStripNoiseDiffByPartitionBase<2, SINGLE_IOV>;
0170 
0171   class SiStripNoiseCorrelationByPartition : public PlotImage<SiStripNoises> {
0172   public:
0173     SiStripNoiseCorrelationByPartition() : PlotImage<SiStripNoises>("SiStrip Noises Correlation By Partition") {
0174       setSingleIov(false);
0175     }
0176 
0177     bool fill(const std::vector<SiStripPI::MetaData>& iovs) override {
0178       std::vector<SiStripPI::MetaData> sorted_iovs = iovs;
0179 
0180       // make absolute sure the IOVs are sortd by since
0181       std::sort(begin(sorted_iovs), end(sorted_iovs), [](auto const& t1, auto const& t2) {
0182         return std::get<0>(t1) < std::get<0>(t2);
0183       });
0184 
0185       auto firstiov = sorted_iovs.front();
0186       auto lastiov = sorted_iovs.back();
0187 
0188       std::shared_ptr<SiStripNoises> last_payload = fetchPayload(std::get<1>(lastiov));
0189       std::shared_ptr<SiStripNoises> first_payload = fetchPayload(std::get<1>(firstiov));
0190 
0191       SiStripNoiseContainer* l_objContainer = new SiStripNoiseContainer(last_payload, lastiov, "");
0192       SiStripNoiseContainer* f_objContainer = new SiStripNoiseContainer(first_payload, firstiov, "");
0193 
0194       l_objContainer->compare(f_objContainer);
0195 
0196       TCanvas canvas("Partition summary", "partition summary", 1200, 1200);
0197       l_objContainer->fillCorrelationByPartition(canvas, 100, 0.1, 10.);
0198 
0199       std::string fileName(m_imageFileName);
0200       canvas.SaveAs(fileName.c_str());
0201 
0202       return true;
0203     }  // fill
0204   };
0205 
0206   class SiStripNoiseConsistencyCheck : public PlotImage<SiStripNoises> {
0207   public:
0208     SiStripNoiseConsistencyCheck() : PlotImage<SiStripNoises>("SiStrip Noise Consistency Check") {
0209       setSingleIov(false);
0210     }
0211 
0212     bool fill(const std::vector<SiStripPI::MetaData>& iovs) override {
0213       std::vector<SiStripPI::MetaData> sorted_iovs = iovs;
0214 
0215       // make absolute sure the IOVs are sortd by since
0216       std::sort(begin(sorted_iovs), end(sorted_iovs), [](auto const& t1, auto const& t2) {
0217         return std::get<0>(t1) < std::get<0>(t2);
0218       });
0219 
0220       auto firstiov = sorted_iovs.front();
0221       auto lastiov = sorted_iovs.back();
0222 
0223       std::shared_ptr<SiStripNoises> last_payload = fetchPayload(std::get<1>(lastiov));
0224       std::shared_ptr<SiStripNoises> first_payload = fetchPayload(std::get<1>(firstiov));
0225 
0226       SiStripNoiseContainer* f_objContainer = new SiStripNoiseContainer(first_payload, firstiov, "");
0227       SiStripNoiseContainer* l_objContainer = new SiStripNoiseContainer(last_payload, lastiov, "");
0228 
0229       f_objContainer->compare(l_objContainer);
0230 
0231       //l_objContainer->printAll();
0232 
0233       TCanvas canvas("Partition summary", "partition summary", 1200, 1000);
0234       //f_objContainer->fillValuePlot(canvas,SiStripPI::STRIP_BASED,100,0.1,10);
0235       f_objContainer->fillValuePlot(canvas, SiStripPI::APV_BASED, 100, 0.1, 10);
0236       //f_objContainer->fillValuePlot(canvas,SiStripPI::MODULE_BASED,100,0.1,10);
0237 
0238       std::string fileName(m_imageFileName);
0239       canvas.SaveAs(fileName.c_str());
0240 
0241       return true;
0242     }  // fill
0243   };
0244 
0245   /************************************************
0246     test class
0247   *************************************************/
0248 
0249   class SiStripNoisesTest : public Histogram1D<SiStripNoises, SINGLE_IOV> {
0250   public:
0251     SiStripNoisesTest()
0252         : Histogram1D<SiStripNoises, SINGLE_IOV>("SiStrip Noise test", "SiStrip Noise test", 10, 0.0, 10.0),
0253           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0254               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath())} {}
0255 
0256     bool fill() override {
0257       auto tag = PlotBase::getTag<0>();
0258       for (auto const& iov : tag.iovs) {
0259         std::shared_ptr<SiStripNoises> payload = Base::fetchPayload(std::get<1>(iov));
0260         if (payload.get()) {
0261           fillWithValue(1.);
0262 
0263           std::stringstream ss;
0264           ss << "Summary of strips noises:" << std::endl;
0265 
0266           //payload->printDebug(ss);
0267           payload->printSummary(ss, &m_trackerTopo);
0268 
0269           std::vector<uint32_t> detid;
0270           payload->getDetIds(detid);
0271 
0272           // for (const auto & d : detid) {
0273           //   int nstrip=0;
0274           //   SiStripNoises::Range range=payload->getRange(d);
0275           //   for( int it=0; it < (range.second-range.first)*8/9; ++it ){
0276           //     auto noise = payload->getNoise(it,range);
0277           //     nstrip++;
0278           //     ss << "DetId="<< d << " Strip=" << nstrip <<": "<< noise << std::endl;
0279           //   } // end of loop on strips
0280           // } // end of loop on detIds
0281 
0282           std::cout << ss.str() << std::endl;
0283 
0284         }  // payload
0285       }  // iovs
0286       return true;
0287     }  // fill
0288   private:
0289     TrackerTopology m_trackerTopo;
0290   };
0291 
0292   /************************************************
0293     SiStrip Noise Profile of 1 IOV for one selected DetId
0294   *************************************************/
0295 
0296   class SiStripNoisePerDetId : public PlotImage<SiStripNoises, SINGLE_IOV> {
0297   public:
0298     SiStripNoisePerDetId() : PlotImage<SiStripNoises, SINGLE_IOV>("SiStrip Noise values Per DetId") {
0299       PlotBase::addInputParam("DetIds");
0300     }
0301 
0302     bool fill() override {
0303       auto tag = PlotBase::getTag<0>();
0304       auto iov = tag.iovs.front();
0305       auto tagname = tag.name;
0306       std::shared_ptr<SiStripNoises> payload = fetchPayload(std::get<1>(iov));
0307 
0308       std::vector<uint32_t> the_detids = {};
0309 
0310       auto paramValues = PlotBase::inputParamValues();
0311       auto ip = paramValues.find("DetIds");
0312       if (ip != paramValues.end()) {
0313         auto input = ip->second;
0314         typedef boost::tokenizer<boost::char_separator<char>> tokenizer;
0315         boost::char_separator<char> sep{","};
0316         tokenizer tok{input, sep};
0317         for (const auto& t : tok) {
0318           the_detids.push_back(atoi(t.c_str()));
0319         }
0320       } else {
0321         edm::LogWarning("SiStripNoisePerDetId")
0322             << "\n WARNING!!!! \n The needed parameter DetIds has not been passed. Will use all Strip DetIds! \n\n";
0323         the_detids.push_back(0xFFFFFFFF);
0324       }
0325 
0326       size_t ndets = the_detids.size();
0327       std::vector<std::shared_ptr<TH1F>> hnoise;
0328       std::vector<std::shared_ptr<TLegend>> legends;
0329       std::vector<unsigned int> v_nAPVs;
0330       std::vector<std::vector<std::shared_ptr<TLine>>> lines;
0331       hnoise.reserve(ndets);
0332       legends.reserve(ndets);
0333 
0334       // determine how the plot will be paginated
0335       auto sides = getClosestFactors(the_detids.size());
0336       edm::LogPrint("SiStripNoisePerDetId") << "Aspect ratio: " << sides.first << ":" << sides.second << std::endl;
0337 
0338       if (payload.get()) {
0339         //=========================
0340         TCanvas canvas("ByDetId", "ByDetId", sides.second * 800, sides.first * 600);
0341         canvas.Divide(sides.second, sides.first);
0342         const auto detInfo =
0343             SiStripDetInfoFileReader::read(edm::FileInPath(SiStripDetInfoFileReader::kDefaultFile).fullPath());
0344         for (const auto& the_detid : the_detids) {
0345           edm::LogPrint("SiStripNoisePerDetId") << "DetId:" << the_detid << std::endl;
0346 
0347           unsigned int nAPVs = detInfo.getNumberOfApvsAndStripLength(the_detid).first;
0348           if (nAPVs == 0)
0349             nAPVs = 6;
0350           v_nAPVs.push_back(nAPVs);
0351 
0352           auto histo =
0353               std::make_shared<TH1F>(Form("Noise profile_%s", std::to_string(the_detid).c_str()),
0354                                      Form("SiStrip Noise profile for DetId: %s;Strip number;SiStrip Noise [ADC counts]",
0355                                           std::to_string(the_detid).c_str()),
0356                                      sistrip::STRIPS_PER_APV * nAPVs,
0357                                      -0.5,
0358                                      (sistrip::STRIPS_PER_APV * nAPVs) - 0.5);
0359 
0360           histo->SetStats(false);
0361           histo->SetTitle("");
0362 
0363           if (the_detid != 0xFFFFFFFF) {
0364             fillHisto(payload, histo, the_detid);
0365           } else {
0366             auto allDetIds = detInfo.getAllDetIds();
0367             for (const auto& id : allDetIds) {
0368               fillHisto(payload, histo, id);
0369             }
0370           }
0371 
0372           SiStripPI::makeNicePlotStyle(histo.get());
0373           histo->GetYaxis()->SetTitleOffset(1.0);
0374           hnoise.push_back(histo);
0375         }  // loop on the detids
0376 
0377         for (size_t index = 0; index < ndets; index++) {
0378           canvas.cd(index + 1);
0379           canvas.cd(index + 1)->SetBottomMargin(0.11);
0380           canvas.cd(index + 1)->SetTopMargin(0.06);
0381           canvas.cd(index + 1)->SetLeftMargin(0.10);
0382           canvas.cd(index + 1)->SetRightMargin(0.02);
0383           hnoise.at(index)->Draw();
0384           hnoise.at(index)->GetYaxis()->SetRangeUser(0, hnoise.at(index)->GetMaximum() * 1.2);
0385           canvas.cd(index)->Update();
0386 
0387           std::vector<int> boundaries;
0388           for (size_t b = 0; b < v_nAPVs.at(index); b++) {
0389             boundaries.push_back(b * sistrip::STRIPS_PER_APV);
0390           }
0391 
0392           std::vector<std::shared_ptr<TLine>> linesVec;
0393           for (const auto& bound : boundaries) {
0394             auto line = std::make_shared<TLine>(hnoise.at(index)->GetBinLowEdge(bound),
0395                                                 canvas.cd(index + 1)->GetUymin(),
0396                                                 hnoise.at(index)->GetBinLowEdge(bound),
0397                                                 canvas.cd(index + 1)->GetUymax());
0398             line->SetLineWidth(1);
0399             line->SetLineStyle(9);
0400             line->SetLineColor(2);
0401             linesVec.push_back(line);
0402           }
0403           lines.push_back(linesVec);
0404 
0405           for (const auto& line : lines.at(index)) {
0406             line->Draw("same");
0407           }
0408 
0409           canvas.cd(index + 1);
0410 
0411           auto ltx = TLatex();
0412           ltx.SetTextFont(62);
0413           ltx.SetTextSize(0.05);
0414           ltx.SetTextAlign(11);
0415           ltx.DrawLatexNDC(gPad->GetLeftMargin(),
0416                            1 - gPad->GetTopMargin() + 0.01,
0417                            Form("SiStrip Noise profile for DetId %s", std::to_string(the_detids[index]).c_str()));
0418 
0419           legends.push_back(std::make_shared<TLegend>(0.55, 0.83, 0.95, 0.93));
0420           legends.at(index)->SetHeader(tagname.c_str(), "C");  // option "C" allows to center the header
0421           legends.at(index)->AddEntry(
0422               hnoise.at(index).get(), ("IOV: " + std::to_string(std::get<0>(iov))).c_str(), "PL");
0423           legends.at(index)->SetTextSize(0.045);
0424           legends.at(index)->Draw("same");
0425         }
0426 
0427         std::string fileName(m_imageFileName);
0428         canvas.SaveAs(fileName.c_str());
0429       }  // payload
0430       return true;
0431     }  // fill
0432 
0433   private:
0434     int nextPerfectSquare(int N) { return std::floor(sqrt(N)) + 1; }
0435 
0436     std::pair<int, int> getClosestFactors(int input) {
0437       if ((input % 2 != 0) && input > 1) {
0438         input += 1;
0439       }
0440 
0441       int testNum = (int)sqrt(input);
0442       while (input % testNum != 0) {
0443         testNum--;
0444       }
0445       return std::make_pair(testNum, input / testNum);
0446     }
0447 
0448     void fillHisto(const std::shared_ptr<SiStripNoises> payload, std::shared_ptr<TH1F>& histo, uint32_t the_detid) {
0449       int nstrip = 0;
0450       SiStripNoises::Range range = payload->getRange(the_detid);
0451       for (int it = 0; it < (range.second - range.first) * 8 / 9; ++it) {
0452         auto noise = payload->getNoise(it, range);
0453         nstrip++;
0454         histo->AddBinContent(nstrip, noise);
0455       }  // end of loop on strips
0456     }
0457   };
0458 
0459   /************************************************
0460     1d histogram of SiStripNoises of 1 IOV 
0461   *************************************************/
0462 
0463   // inherit from one of the predefined plot class: Histogram1D
0464   class SiStripNoiseValue : public Histogram1D<SiStripNoises, SINGLE_IOV> {
0465   public:
0466     SiStripNoiseValue()
0467         : Histogram1D<SiStripNoises, SINGLE_IOV>("SiStrip Noise values", "SiStrip Noise values", 100, 0.0, 10.0) {}
0468 
0469     bool fill() override {
0470       auto tag = PlotBase::getTag<0>();
0471       for (auto const& iov : tag.iovs) {
0472         std::shared_ptr<SiStripNoises> payload = Base::fetchPayload(std::get<1>(iov));
0473         if (payload.get()) {
0474           std::vector<uint32_t> detid;
0475           payload->getDetIds(detid);
0476 
0477           for (const auto& d : detid) {
0478             SiStripNoises::Range range = payload->getRange(d);
0479             for (int it = 0; it < (range.second - range.first) * 8 / 9; ++it) {
0480               auto noise = payload->getNoise(it, range);
0481               //to be used to fill the histogram
0482               fillWithValue(noise);
0483             }  // loop over APVs
0484           }  // loop over detIds
0485         }  // payload
0486       }  // iovs
0487       return true;
0488     }  // fill
0489   };
0490 
0491   /************************************************
0492     1d histogram of SiStripNoises of 1 IOV per Detid
0493   *************************************************/
0494 
0495   // inherit from one of the predefined plot class: Histogram1D
0496   class SiStripNoiseValuePerDetId : public Histogram1D<SiStripNoises, SINGLE_IOV> {
0497   public:
0498     SiStripNoiseValuePerDetId()
0499         : Histogram1D<SiStripNoises, SINGLE_IOV>(
0500               "SiStrip Noise values per DetId", "SiStrip Noise values per DetId", 100, 0.0, 10.0) {
0501       PlotBase::addInputParam("DetId");
0502     }
0503 
0504     bool fill() override {
0505       auto tag = PlotBase::getTag<0>();
0506       for (auto const& iov : tag.iovs) {
0507         std::shared_ptr<SiStripNoises> payload = Base::fetchPayload(std::get<1>(iov));
0508         unsigned int the_detid(0xFFFFFFFF);
0509         auto paramValues = PlotBase::inputParamValues();
0510         auto ip = paramValues.find("DetId");
0511         if (ip != paramValues.end()) {
0512           the_detid = std::stoul(ip->second);
0513         }
0514 
0515         if (payload.get()) {
0516           SiStripNoises::Range range = payload->getRange(the_detid);
0517           for (int it = 0; it < (range.second - range.first) * 8 / 9; ++it) {
0518             auto noise = payload->getNoise(it, range);
0519             //to be used to fill the histogram
0520             fillWithValue(noise);
0521           }  // loop over APVs
0522         }  // payload
0523       }  // iovs
0524       return true;
0525     }  // fill
0526   };
0527 
0528   /************************************************
0529     templated 1d histogram of SiStripNoises of 1 IOV
0530   *************************************************/
0531 
0532   // inherit from one of the predefined plot class: PlotImage
0533   template <SiStripPI::OpMode op_mode_>
0534   class SiStripNoiseDistribution : public PlotImage<SiStripNoises, SINGLE_IOV> {
0535   public:
0536     SiStripNoiseDistribution() : PlotImage<SiStripNoises, SINGLE_IOV>("SiStrip Noise values") {}
0537 
0538     bool fill() override {
0539       auto tag = PlotBase::getTag<0>();
0540       auto iov = tag.iovs.front();
0541 
0542       TGaxis::SetMaxDigits(3);
0543       gStyle->SetOptStat("emr");
0544 
0545       std::shared_ptr<SiStripNoises> payload = fetchPayload(std::get<1>(iov));
0546 
0547       auto mon1D = std::unique_ptr<SiStripPI::Monitor1D>(new SiStripPI::Monitor1D(
0548           op_mode_,
0549           "Noise",
0550           Form("#LT Strip Noise #GT per %s for IOV [%s];#LTStrip Noise per %s#GT [ADC counts];n. %ss",
0551                opType(op_mode_).c_str(),
0552                std::to_string(std::get<0>(iov)).c_str(),
0553                opType(op_mode_).c_str(),
0554                opType(op_mode_).c_str()),
0555           100,
0556           0.1,
0557           10.));
0558 
0559       unsigned int prev_det = 0, prev_apv = 0;
0560       SiStripPI::Entry enoise;
0561 
0562       std::vector<uint32_t> detids;
0563       payload->getDetIds(detids);
0564 
0565       // loop on payload
0566 
0567       for (const auto& d : detids) {
0568         SiStripNoises::Range range = payload->getRange(d);
0569 
0570         unsigned int istrip = 0;
0571 
0572         for (int it = 0; it < (range.second - range.first) * 8 / 9; ++it) {
0573           auto noise = payload->getNoise(it, range);
0574           bool flush = false;
0575           switch (op_mode_) {
0576             case (SiStripPI::APV_BASED):
0577               flush = (prev_det != 0 && prev_apv != istrip / sistrip::STRIPS_PER_APV);
0578               break;
0579             case (SiStripPI::MODULE_BASED):
0580               flush = (prev_det != 0 && prev_det != d);
0581               break;
0582             case (SiStripPI::STRIP_BASED):
0583               flush = (istrip != 0);
0584               break;
0585           }
0586 
0587           if (flush) {
0588             mon1D->Fill(prev_apv, prev_det, enoise.mean());
0589             enoise.reset();
0590           }
0591 
0592           enoise.add(std::min<float>(noise, 30.5));
0593           prev_apv = istrip / sistrip::STRIPS_PER_APV;
0594           istrip++;
0595         }
0596         prev_det = d;
0597       }
0598 
0599       //=========================
0600       TCanvas canvas("Partition summary", "partition summary", 1200, 1000);
0601       canvas.cd();
0602       canvas.SetBottomMargin(0.11);
0603       canvas.SetTopMargin(0.07);
0604       canvas.SetLeftMargin(0.13);
0605       canvas.SetRightMargin(0.05);
0606       canvas.Modified();
0607 
0608       auto hist = mon1D->getHist();
0609       SiStripPI::makeNicePlotStyle(&hist);
0610       hist.SetStats(kTRUE);
0611       hist.SetFillColorAlpha(kRed, 0.35);
0612       hist.Draw();
0613 
0614       canvas.Update();
0615 
0616       TPaveStats* st = (TPaveStats*)hist.GetListOfFunctions()->FindObject("stats");
0617       st->SetLineColor(kRed);
0618       st->SetTextColor(kRed);
0619       st->SetX1NDC(.75);
0620       st->SetX2NDC(.95);
0621       st->SetY1NDC(.83);
0622       st->SetY2NDC(.93);
0623 
0624       TLegend legend = TLegend(0.13, 0.83, 0.43, 0.93);
0625       legend.SetHeader(Form("SiStrip Noise values per %s", opType(op_mode_).c_str()),
0626                        "C");  // option "C" allows to center the header
0627       legend.AddEntry(&hist, ("IOV: " + std::to_string(std::get<0>(iov))).c_str(), "F");
0628       legend.SetTextSize(0.025);
0629       legend.Draw("same");
0630 
0631       std::string fileName(m_imageFileName);
0632       canvas.SaveAs(fileName.c_str());
0633 
0634       return true;
0635     }
0636 
0637     std::string opType(SiStripPI::OpMode mode) {
0638       std::string types[3] = {"Strip", "APV", "Module"};
0639       return types[mode];
0640     }
0641   };
0642 
0643   typedef SiStripNoiseDistribution<SiStripPI::STRIP_BASED> SiStripNoiseValuePerStrip;
0644   typedef SiStripNoiseDistribution<SiStripPI::APV_BASED> SiStripNoiseValuePerAPV;
0645   typedef SiStripNoiseDistribution<SiStripPI::MODULE_BASED> SiStripNoiseValuePerModule;
0646 
0647   /************************************************
0648   template 1d histogram comparison of SiStripNoises of 1 IOV
0649   *************************************************/
0650 
0651   // inherit from one of the predefined plot class: PlotImage
0652   template <SiStripPI::OpMode op_mode_, int ntags, IOVMultiplicity nIOVs>
0653   class SiStripNoiseDistributionComparisonBase : public PlotImage<SiStripNoises, nIOVs, ntags> {
0654   public:
0655     SiStripNoiseDistributionComparisonBase()
0656         : PlotImage<SiStripNoises, nIOVs, ntags>("SiStrip Noise values comparison") {}
0657 
0658     bool fill() override {
0659       TGaxis::SetExponentOffset(-0.1, 0.01, "y");  // X and Y offset for Y axis
0660 
0661       // trick to deal with the multi-ioved tag and two tag case at the same time
0662       auto theIOVs = PlotBase::getTag<0>().iovs;
0663       auto tagname1 = PlotBase::getTag<0>().name;
0664       std::string tagname2 = "";
0665       auto firstiov = theIOVs.front();
0666       SiStripPI::MetaData lastiov;
0667 
0668       // we don't support (yet) comparison with more than 2 tags
0669       assert(this->m_plotAnnotations.ntags < 3);
0670 
0671       if (this->m_plotAnnotations.ntags == 2) {
0672         auto tag2iovs = PlotBase::getTag<1>().iovs;
0673         tagname2 = PlotBase::getTag<1>().name;
0674         lastiov = tag2iovs.front();
0675       } else {
0676         lastiov = theIOVs.back();
0677       }
0678 
0679       std::shared_ptr<SiStripNoises> f_payload = this->fetchPayload(std::get<1>(firstiov));
0680       std::shared_ptr<SiStripNoises> l_payload = this->fetchPayload(std::get<1>(lastiov));
0681 
0682       auto f_mon = std::unique_ptr<SiStripPI::Monitor1D>(new SiStripPI::Monitor1D(
0683           op_mode_,
0684           "f_Noise",
0685           Form(";#LTStrip Noise per %s#GT [ADC counts];n. %ss", opType(op_mode_).c_str(), opType(op_mode_).c_str()),
0686           100,
0687           0.1,
0688           10.));
0689 
0690       auto l_mon = std::unique_ptr<SiStripPI::Monitor1D>(new SiStripPI::Monitor1D(
0691           op_mode_,
0692           "l_Noise",
0693           Form(";#LTStrip Noise per %s#GT [ADC counts];n. %ss", opType(op_mode_).c_str(), opType(op_mode_).c_str()),
0694           100,
0695           0.1,
0696           10.));
0697 
0698       unsigned int prev_det = 0, prev_apv = 0;
0699       SiStripPI::Entry enoise;
0700 
0701       std::vector<uint32_t> f_detid;
0702       f_payload->getDetIds(f_detid);
0703 
0704       // loop on first payload
0705       for (const auto& d : f_detid) {
0706         SiStripNoises::Range range = f_payload->getRange(d);
0707 
0708         unsigned int istrip = 0;
0709         for (int it = 0; it < (range.second - range.first) * 8 / 9; ++it) {
0710           float noise = f_payload->getNoise(it, range);
0711           //to be used to fill the histogram
0712 
0713           bool flush = false;
0714           switch (op_mode_) {
0715             case (SiStripPI::APV_BASED):
0716               flush = (prev_det != 0 && prev_apv != istrip / sistrip::STRIPS_PER_APV);
0717               break;
0718             case (SiStripPI::MODULE_BASED):
0719               flush = (prev_det != 0 && prev_det != d);
0720               break;
0721             case (SiStripPI::STRIP_BASED):
0722               flush = (istrip != 0);
0723               break;
0724           }
0725 
0726           if (flush) {
0727             f_mon->Fill(prev_apv, prev_det, enoise.mean());
0728             enoise.reset();
0729           }
0730           enoise.add(std::min<float>(noise, 30.5));
0731           prev_apv = istrip / sistrip::STRIPS_PER_APV;
0732           istrip++;
0733         }
0734         prev_det = d;
0735       }
0736 
0737       prev_det = 0, prev_apv = 0;
0738       enoise.reset();
0739 
0740       std::vector<uint32_t> l_detid;
0741       l_payload->getDetIds(l_detid);
0742 
0743       // loop on first payload
0744 
0745       for (const auto& d : l_detid) {
0746         SiStripNoises::Range range = l_payload->getRange(d);
0747 
0748         unsigned int istrip = 0;
0749         for (int it = 0; it < (range.second - range.first) * 8 / 9; ++it) {
0750           float noise = l_payload->getNoise(it, range);
0751 
0752           bool flush = false;
0753           switch (op_mode_) {
0754             case (SiStripPI::APV_BASED):
0755               flush = (prev_det != 0 && prev_apv != istrip / sistrip::STRIPS_PER_APV);
0756               break;
0757             case (SiStripPI::MODULE_BASED):
0758               flush = (prev_det != 0 && prev_det != d);
0759               break;
0760             case (SiStripPI::STRIP_BASED):
0761               flush = (istrip != 0);
0762               break;
0763           }
0764 
0765           if (flush) {
0766             l_mon->Fill(prev_apv, prev_det, enoise.mean());
0767             enoise.reset();
0768           }
0769 
0770           enoise.add(std::min<float>(noise, 30.5));
0771           prev_apv = istrip / sistrip::STRIPS_PER_APV;
0772           istrip++;
0773         }
0774         prev_det = d;
0775       }
0776 
0777       auto h_first = f_mon->getHist();
0778       h_first.SetStats(kFALSE);
0779       auto h_last = l_mon->getHist();
0780       h_last.SetStats(kFALSE);
0781 
0782       SiStripPI::makeNicePlotStyle(&h_first);
0783       SiStripPI::makeNicePlotStyle(&h_last);
0784 
0785       h_first.GetYaxis()->CenterTitle(true);
0786       h_last.GetYaxis()->CenterTitle(true);
0787 
0788       h_first.GetXaxis()->CenterTitle(true);
0789       h_last.GetXaxis()->CenterTitle(true);
0790 
0791       h_first.SetLineWidth(2);
0792       h_last.SetLineWidth(2);
0793 
0794       h_first.SetLineColor(kBlack);
0795       h_last.SetLineColor(kBlue);
0796 
0797       //=========================
0798       TCanvas canvas("Partition summary", "partition summary", 1200, 1000);
0799       canvas.cd();
0800       canvas.SetTopMargin(0.06);
0801       canvas.SetBottomMargin(0.10);
0802       canvas.SetLeftMargin(0.13);
0803       canvas.SetRightMargin(0.05);
0804       canvas.Modified();
0805 
0806       float theMax = (h_first.GetMaximum() > h_last.GetMaximum()) ? h_first.GetMaximum() : h_last.GetMaximum();
0807 
0808       h_first.SetMaximum(theMax * 1.20);
0809       h_last.SetMaximum(theMax * 1.20);
0810 
0811       h_first.Draw();
0812       h_last.SetFillColorAlpha(kBlue, 0.15);
0813       h_last.Draw("same");
0814 
0815       TLegend legend = TLegend(0.13, 0.83, 0.95, 0.94);
0816       if (this->m_plotAnnotations.ntags == 2) {
0817         legend.SetHeader("#bf{Two Tags Comparison}", "C");  // option "C" allows to center the header
0818         legend.AddEntry(&h_first, (tagname1 + " : " + std::to_string(std::get<0>(firstiov))).c_str(), "F");
0819         legend.AddEntry(&h_last, (tagname2 + " : " + std::to_string(std::get<0>(lastiov))).c_str(), "F");
0820       } else {
0821         legend.SetHeader(("tag: #bf{" + tagname1 + "}").c_str(), "C");  // option "C" allows to center the header
0822         legend.AddEntry(&h_first, ("IOV since: " + std::to_string(std::get<0>(firstiov))).c_str(), "F");
0823         legend.AddEntry(&h_last, ("IOV since: " + std::to_string(std::get<0>(lastiov))).c_str(), "F");
0824       }
0825       legend.SetTextSize(0.025);
0826       legend.Draw("same");
0827 
0828       auto ltx = TLatex();
0829       ltx.SetTextFont(62);
0830       ltx.SetTextSize(0.05);
0831       ltx.SetTextAlign(11);
0832       ltx.DrawLatexNDC(gPad->GetLeftMargin(),
0833                        1 - gPad->GetTopMargin() + 0.01,
0834                        Form("#LTSiStrip Noise#GT Comparison per %s", opType(op_mode_).c_str()));
0835 
0836       std::string fileName(this->m_imageFileName);
0837       canvas.SaveAs(fileName.c_str());
0838 
0839       return true;
0840     }
0841 
0842     std::string opType(SiStripPI::OpMode mode) {
0843       std::string types[3] = {"Strip", "APV", "Module"};
0844       return types[mode];
0845     }
0846   };
0847 
0848   template <SiStripPI::OpMode op_mode_>
0849   using SiStripNoiseDistributionComparisonSingleTag = SiStripNoiseDistributionComparisonBase<op_mode_, 1, MULTI_IOV>;
0850 
0851   template <SiStripPI::OpMode op_mode_>
0852   using SiStripNoiseDistributionComparisonTwoTags = SiStripNoiseDistributionComparisonBase<op_mode_, 2, SINGLE_IOV>;
0853 
0854   typedef SiStripNoiseDistributionComparisonSingleTag<SiStripPI::STRIP_BASED>
0855       SiStripNoiseValueComparisonPerStripSingleTag;
0856   typedef SiStripNoiseDistributionComparisonSingleTag<SiStripPI::APV_BASED> SiStripNoiseValueComparisonPerAPVSingleTag;
0857   typedef SiStripNoiseDistributionComparisonSingleTag<SiStripPI::MODULE_BASED>
0858       SiStripNoiseValueComparisonPerModuleSingleTag;
0859 
0860   typedef SiStripNoiseDistributionComparisonTwoTags<SiStripPI::STRIP_BASED> SiStripNoiseValueComparisonPerStripTwoTags;
0861   typedef SiStripNoiseDistributionComparisonTwoTags<SiStripPI::APV_BASED> SiStripNoiseValueComparisonPerAPVTwoTags;
0862   typedef SiStripNoiseDistributionComparisonTwoTags<SiStripPI::MODULE_BASED> SiStripNoiseValueComparisonPerModuleTwoTags;
0863 
0864   /************************************************
0865     1d histogram comparison of SiStripNoises of 1 IOV 
0866   *************************************************/
0867 
0868   // inherit from one of the predefined plot class: PlotImage
0869 
0870   template <int ntags, IOVMultiplicity nIOVs>
0871   class SiStripNoiseValueComparisonBase : public PlotImage<SiStripNoises, nIOVs, ntags> {
0872   public:
0873     SiStripNoiseValueComparisonBase() : PlotImage<SiStripNoises, nIOVs, ntags>("SiStrip Noise values comparison") {}
0874 
0875     bool fill() override {
0876       TGaxis::SetExponentOffset(-0.1, 0.01, "y");  // X and Y offset for Y axis
0877 
0878       // trick to deal with the multi-ioved tag and two tag case at the same time
0879       auto theIOVs = PlotBase::getTag<0>().iovs;
0880       auto tagname1 = PlotBase::getTag<0>().name;
0881       std::string tagname2 = "";
0882       auto firstiov = theIOVs.front();
0883       SiStripPI::MetaData lastiov;
0884 
0885       // we don't support (yet) comparison with more than 2 tags
0886       assert(this->m_plotAnnotations.ntags < 3);
0887 
0888       if (this->m_plotAnnotations.ntags == 2) {
0889         auto tag2iovs = PlotBase::getTag<1>().iovs;
0890         tagname2 = PlotBase::getTag<1>().name;
0891         lastiov = tag2iovs.front();
0892       } else {
0893         lastiov = theIOVs.back();
0894       }
0895 
0896       std::shared_ptr<SiStripNoises> f_payload = this->fetchPayload(std::get<1>(firstiov));
0897       std::shared_ptr<SiStripNoises> l_payload = this->fetchPayload(std::get<1>(lastiov));
0898 
0899       auto h_first = std::make_unique<TH1F>("f_Noise", ";Strip Noise [ADC counts];n. strips", 100, 0.1, 10.);
0900       h_first->SetStats(false);
0901 
0902       auto h_last = std::make_unique<TH1F>("l_Noise", ";Strip Noise [ADC counts];n. strips", 100, 0.1, 10.);
0903       h_last->SetStats(false);
0904 
0905       std::vector<uint32_t> f_detid;
0906       f_payload->getDetIds(f_detid);
0907 
0908       // loop on first payload
0909       for (const auto& d : f_detid) {
0910         SiStripNoises::Range range = f_payload->getRange(d);
0911         for (int it = 0; it < (range.second - range.first) * 8 / 9; ++it) {
0912           float noise = f_payload->getNoise(it, range);
0913           //to be used to fill the histogram
0914           h_first->Fill(noise);
0915         }  // loop over strips
0916       }
0917 
0918       std::vector<uint32_t> l_detid;
0919       l_payload->getDetIds(l_detid);
0920 
0921       // loop on first payload
0922       for (const auto& d : l_detid) {
0923         SiStripNoises::Range range = l_payload->getRange(d);
0924         for (int it = 0; it < (range.second - range.first) * 8 / 9; ++it) {
0925           float noise = l_payload->getNoise(it, range);
0926           //to be used to fill the histogram
0927           h_last->Fill(noise);
0928         }  // loop over strips
0929       }
0930 
0931       SiStripPI::makeNicePlotStyle(h_first.get());
0932       SiStripPI::makeNicePlotStyle(h_last.get());
0933 
0934       h_first->GetYaxis()->CenterTitle(true);
0935       h_last->GetYaxis()->CenterTitle(true);
0936 
0937       h_first->GetXaxis()->CenterTitle(true);
0938       h_last->GetXaxis()->CenterTitle(true);
0939 
0940       h_first->SetLineWidth(2);
0941       h_last->SetLineWidth(2);
0942 
0943       h_first->SetLineColor(kBlack);
0944       h_last->SetLineColor(kBlue);
0945 
0946       //=========================
0947       TCanvas canvas("Partition summary", "partition summary", 1200, 1000);
0948       canvas.cd();
0949       canvas.SetTopMargin(0.06);
0950       canvas.SetBottomMargin(0.10);
0951       canvas.SetLeftMargin(0.13);
0952       canvas.SetRightMargin(0.05);
0953       canvas.Modified();
0954 
0955       float theMax = (h_first->GetMaximum() > h_last->GetMaximum()) ? h_first->GetMaximum() : h_last->GetMaximum();
0956 
0957       h_first->SetMaximum(theMax * 1.20);
0958       h_last->SetMaximum(theMax * 1.20);
0959 
0960       h_first->Draw();
0961       h_last->SetFillColorAlpha(kBlue, 0.15);
0962       h_last->Draw("same");
0963 
0964       TLegend legend = TLegend(0.13, 0.83, 0.95, 0.94);
0965       if (this->m_plotAnnotations.ntags == 2) {
0966         legend.SetHeader("#bf{Two Tags Comparison}", "C");  // option "C" allows to center the header
0967         legend.AddEntry(h_first.get(), (tagname1 + " : " + std::to_string(std::get<0>(firstiov))).c_str(), "F");
0968         legend.AddEntry(h_last.get(), (tagname2 + " : " + std::to_string(std::get<0>(lastiov))).c_str(), "F");
0969       } else {
0970         legend.SetHeader(("tag: #bf{" + tagname1 + "}").c_str(), "C");  // option "C" allows to center the header
0971         legend.AddEntry(h_first.get(), ("IOV since: " + std::to_string(std::get<0>(firstiov))).c_str(), "F");
0972         legend.AddEntry(h_last.get(), ("IOV since: " + std::to_string(std::get<0>(lastiov))).c_str(), "F");
0973       }
0974       legend.SetTextSize(0.025);
0975       legend.Draw("same");
0976 
0977       auto ltx = TLatex();
0978       ltx.SetTextFont(62);
0979       ltx.SetTextSize(0.05);
0980       ltx.SetTextAlign(11);
0981       ltx.DrawLatexNDC(gPad->GetLeftMargin(), 1 - gPad->GetTopMargin() + 0.01, "SiStrip Noise Values Comparison");
0982 
0983       std::string fileName(this->m_imageFileName);
0984       canvas.SaveAs(fileName.c_str());
0985 
0986       return true;
0987     }
0988   };
0989 
0990   using SiStripNoiseValueComparisonSingleTag = SiStripNoiseValueComparisonBase<1, MULTI_IOV>;
0991   using SiStripNoiseValueComparisonTwoTags = SiStripNoiseValueComparisonBase<2, SINGLE_IOV>;
0992 
0993   /************************************************
0994     SiStrip Noise Tracker Map 
0995   *************************************************/
0996 
0997   template <SiStripPI::estimator est>
0998   class SiStripNoiseTrackerMap : public PlotImage<SiStripNoises, SINGLE_IOV> {
0999   public:
1000     SiStripNoiseTrackerMap()
1001         : PlotImage<SiStripNoises, SINGLE_IOV>("Tracker Map of SiStripNoise " + estimatorType(est) + " per module") {}
1002 
1003     bool fill() override {
1004       auto tag = PlotBase::getTag<0>();
1005       auto iov = tag.iovs.front();
1006       std::shared_ptr<SiStripNoises> payload = fetchPayload(std::get<1>(iov));
1007 
1008       std::string titleMap =
1009           "Tracker Map of Noise " + estimatorType(est) + " per module (payload : " + std::get<1>(iov) + ")";
1010 
1011       std::unique_ptr<TrackerMap> tmap = std::make_unique<TrackerMap>("SiStripNoises");
1012       tmap->setTitle(titleMap);
1013       tmap->setPalette(1);
1014 
1015       // storage of info
1016       std::map<unsigned int, float> info_per_detid;
1017 
1018       SiStripNoises::RegistryIterator rit = payload->getRegistryVectorBegin(), erit = payload->getRegistryVectorEnd();
1019       uint16_t Nstrips;
1020       std::vector<float> vstripnoise;
1021       double mean, rms, min, max;
1022       for (; rit != erit; ++rit) {
1023         Nstrips =
1024             (rit->iend - rit->ibegin) * 8 / 9;  //number of strips = number of chars * char size / strip noise size
1025         vstripnoise.resize(Nstrips);
1026         payload->allNoises(
1027             vstripnoise,
1028             make_pair(payload->getDataVectorBegin() + rit->ibegin, payload->getDataVectorBegin() + rit->iend));
1029         mean = 0;
1030         rms = 0;
1031         min = 10000;
1032         max = 0;
1033 
1034         DetId detId(rit->detid);
1035 
1036         for (size_t i = 0; i < Nstrips; ++i) {
1037           mean += vstripnoise[i];
1038           rms += vstripnoise[i] * vstripnoise[i];
1039           if (vstripnoise[i] < min)
1040             min = vstripnoise[i];
1041           if (vstripnoise[i] > max)
1042             max = vstripnoise[i];
1043         }
1044 
1045         mean /= Nstrips;
1046         if ((rms / Nstrips - mean * mean) > 0.) {
1047           rms = sqrt(rms / Nstrips - mean * mean);
1048         } else {
1049           rms = 0.;
1050         }
1051 
1052         switch (est) {
1053           case SiStripPI::min:
1054             info_per_detid[rit->detid] = min;
1055             break;
1056           case SiStripPI::max:
1057             info_per_detid[rit->detid] = max;
1058             break;
1059           case SiStripPI::mean:
1060             info_per_detid[rit->detid] = mean;
1061             break;
1062           case SiStripPI::rms:
1063             info_per_detid[rit->detid] = rms;
1064             break;
1065           default:
1066             edm::LogWarning("LogicError") << "Unknown estimator: " << est;
1067             break;
1068         }
1069       }
1070 
1071       // loop on the map
1072       for (const auto& item : info_per_detid) {
1073         tmap->fill(item.first, item.second);
1074       }
1075 
1076       auto range = SiStripPI::getTheRange(info_per_detid, 2);
1077 
1078       //=========================
1079 
1080       std::string fileName(m_imageFileName);
1081       if (est == SiStripPI::rms && (range.first < 0.)) {
1082         tmap->save(true, 0., range.second, fileName);
1083       } else {
1084         tmap->save(true, range.first, range.second, fileName);
1085       }
1086 
1087       return true;
1088     }
1089   };
1090 
1091   typedef SiStripNoiseTrackerMap<SiStripPI::min> SiStripNoiseMin_TrackerMap;
1092   typedef SiStripNoiseTrackerMap<SiStripPI::max> SiStripNoiseMax_TrackerMap;
1093   typedef SiStripNoiseTrackerMap<SiStripPI::mean> SiStripNoiseMean_TrackerMap;
1094   typedef SiStripNoiseTrackerMap<SiStripPI::rms> SiStripNoiseRMS_TrackerMap;
1095 
1096   /************************************************
1097     SiStrip Noise Tracker Map  (ratio with previous gain per detid)
1098   *************************************************/
1099 
1100   template <SiStripPI::estimator est, int ntags, IOVMultiplicity nIOVs>
1101   class SiStripNoiseRatioWithPreviousIOVTrackerMapBase : public PlotImage<SiStripNoises, nIOVs, ntags> {
1102   public:
1103     SiStripNoiseRatioWithPreviousIOVTrackerMapBase()
1104         : PlotImage<SiStripNoises, nIOVs, ntags>("Tracker Map of ratio of SiStripNoises " + estimatorType(est) +
1105                                                  "with previous IOV") {
1106       PlotBase::addInputParam("nsigma");
1107     }
1108 
1109     std::map<unsigned int, float> computeEstimator(std::shared_ptr<SiStripNoises> payload) {
1110       std::map<unsigned int, float> info_per_detid;
1111       SiStripNoises::RegistryIterator rit = payload->getRegistryVectorBegin(), erit = payload->getRegistryVectorEnd();
1112       uint16_t Nstrips;
1113       std::vector<float> vstripnoise;
1114       double mean, rms, min, max;
1115       for (; rit != erit; ++rit) {
1116         Nstrips =
1117             (rit->iend - rit->ibegin) * 8 / 9;  //number of strips = number of chars * char size / strip noise size
1118         vstripnoise.resize(Nstrips);
1119         payload->allNoises(
1120             vstripnoise,
1121             make_pair(payload->getDataVectorBegin() + rit->ibegin, payload->getDataVectorBegin() + rit->iend));
1122         mean = 0;
1123         rms = 0;
1124         min = 10000;
1125         max = 0;
1126 
1127         DetId detId(rit->detid);
1128 
1129         for (size_t i = 0; i < Nstrips; ++i) {
1130           mean += vstripnoise[i];
1131           rms += vstripnoise[i] * vstripnoise[i];
1132           if (vstripnoise[i] < min)
1133             min = vstripnoise[i];
1134           if (vstripnoise[i] > max)
1135             max = vstripnoise[i];
1136         }
1137 
1138         mean /= Nstrips;
1139         if ((rms / Nstrips - mean * mean) > 0.) {
1140           rms = sqrt(rms / Nstrips - mean * mean);
1141         } else {
1142           rms = 0.;
1143         }
1144         switch (est) {
1145           case SiStripPI::min:
1146             info_per_detid[rit->detid] = min;
1147             break;
1148           case SiStripPI::max:
1149             info_per_detid[rit->detid] = max;
1150             break;
1151           case SiStripPI::mean:
1152             info_per_detid[rit->detid] = mean;
1153             break;
1154           case SiStripPI::rms:
1155             info_per_detid[rit->detid] = rms;
1156             break;
1157           default:
1158             edm::LogWarning("LogicError") << "Unknown estimator: " << est;
1159             break;
1160         }
1161       }
1162       return info_per_detid;
1163     }
1164 
1165     bool fill() override {
1166       unsigned int nsigma(1);
1167 
1168       auto paramValues = PlotBase::inputParamValues();
1169       auto ip = paramValues.find("nsigma");
1170       if (ip != paramValues.end()) {
1171         nsigma = std::stoul(ip->second);
1172       }
1173 
1174       // trick to deal with the multi-ioved tag and two tag case at the same time
1175       auto theIOVs = PlotBase::getTag<0>().iovs;
1176       auto tagname1 = PlotBase::getTag<0>().name;
1177       std::string tagname2 = "";
1178       auto firstiov = theIOVs.front();
1179       SiStripPI::MetaData lastiov;
1180 
1181       // we don't support (yet) comparison with more than 2 tags
1182       assert(this->m_plotAnnotations.ntags < 3);
1183 
1184       if (this->m_plotAnnotations.ntags == 2) {
1185         auto tag2iovs = PlotBase::getTag<1>().iovs;
1186         tagname2 = PlotBase::getTag<1>().name;
1187         lastiov = tag2iovs.front();
1188       } else {
1189         lastiov = theIOVs.back();
1190       }
1191 
1192       std::shared_ptr<SiStripNoises> last_payload = this->fetchPayload(std::get<1>(lastiov));
1193       std::shared_ptr<SiStripNoises> first_payload = this->fetchPayload(std::get<1>(firstiov));
1194 
1195       std::string titleMap = "SiStripNoise " + estimatorType(est) + " ratio per module average (IOV: ";
1196 
1197       titleMap += std::to_string(std::get<0>(firstiov));
1198       titleMap += "/ IOV:";
1199       titleMap += std::to_string(std::get<0>(lastiov));
1200       titleMap += ")";
1201       titleMap += +" " + std::to_string(nsigma) + " std. dev. saturation";
1202 
1203       std::unique_ptr<TrackerMap> tmap = std::make_unique<TrackerMap>("SiStripNoises");
1204       tmap->setTitle(titleMap);
1205       tmap->setPalette(1);
1206 
1207       std::map<unsigned int, float> map_first, map_second;
1208 
1209       map_first = computeEstimator(first_payload);
1210       map_second = computeEstimator(last_payload);
1211       std::map<unsigned int, float> cachedRatio;
1212 
1213       for (auto entry : map_first) {
1214         auto it2 = map_second.find(entry.first);
1215         if (it2 == map_second.end() || it2->second == 0)
1216           continue;
1217         tmap->fill(entry.first, entry.second / it2->second);
1218         cachedRatio[entry.first] = (entry.second / it2->second);
1219       }
1220 
1221       auto range = SiStripPI::getTheRange(cachedRatio, nsigma);
1222       std::string fileName(this->m_imageFileName);
1223       if (est == SiStripPI::rms && (range.first < 0.)) {
1224         tmap->save(true, 0., range.second, fileName);
1225       } else {
1226         tmap->save(true, range.first, range.second, fileName);
1227       }
1228 
1229       return true;
1230     }
1231   };
1232 
1233   template <SiStripPI::estimator est>
1234   using SiStripNoiseRatioWithPreviousIOVTrackerMapSingleTag =
1235       SiStripNoiseRatioWithPreviousIOVTrackerMapBase<est, 1, MULTI_IOV>;
1236 
1237   template <SiStripPI::estimator est>
1238   using SiStripNoiseRatioWithPreviousIOVTrackerMapTwoTags =
1239       SiStripNoiseRatioWithPreviousIOVTrackerMapBase<est, 2, SINGLE_IOV>;
1240 
1241   typedef SiStripNoiseRatioWithPreviousIOVTrackerMapSingleTag<SiStripPI::min>
1242       SiStripNoiseMin_RatioWithPreviousIOVTrackerMapSingleTag;
1243   typedef SiStripNoiseRatioWithPreviousIOVTrackerMapSingleTag<SiStripPI::max>
1244       SiStripNoiseMax_RatioWithPreviousIOVTrackerMapSingleTag;
1245   typedef SiStripNoiseRatioWithPreviousIOVTrackerMapSingleTag<SiStripPI::mean>
1246       SiStripNoiseMean_RatioWithPreviousIOVTrackerMapSingleTag;
1247   typedef SiStripNoiseRatioWithPreviousIOVTrackerMapSingleTag<SiStripPI::rms>
1248       SiStripNoiseRms_RatioWithPreviousIOVTrackerMapSingleTag;
1249 
1250   typedef SiStripNoiseRatioWithPreviousIOVTrackerMapTwoTags<SiStripPI::min>
1251       SiStripNoiseMin_RatioWithPreviousIOVTrackerMapTwoTags;
1252   typedef SiStripNoiseRatioWithPreviousIOVTrackerMapTwoTags<SiStripPI::max>
1253       SiStripNoiseMax_RatioWithPreviousIOVTrackerMapTwoTags;
1254   typedef SiStripNoiseRatioWithPreviousIOVTrackerMapTwoTags<SiStripPI::mean>
1255       SiStripNoiseMean_RatioWithPreviousIOVTrackerMapTwoTags;
1256   typedef SiStripNoiseRatioWithPreviousIOVTrackerMapTwoTags<SiStripPI::rms>
1257       SiStripNoiseRms_RatioWithPreviousIOVTrackerMapTwoTags;
1258 
1259   /************************************************
1260   SiStrip Noise Tracker Summaries 
1261   *************************************************/
1262 
1263   template <SiStripPI::estimator est>
1264   class SiStripNoiseByRegion : public PlotImage<SiStripNoises, SINGLE_IOV> {
1265   public:
1266     SiStripNoiseByRegion()
1267         : PlotImage<SiStripNoises, SINGLE_IOV>("SiStrip Noise " + estimatorType(est) + " by Region"),
1268           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
1269               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath())} {}
1270 
1271     bool fill() override {
1272       auto tag = PlotBase::getTag<0>();
1273       auto iov = tag.iovs.front();
1274       std::shared_ptr<SiStripNoises> payload = fetchPayload(std::get<1>(iov));
1275 
1276       SiStripDetSummary summaryNoise{&m_trackerTopo};
1277 
1278       SiStripPI::fillNoiseDetSummary(summaryNoise, payload, est);
1279 
1280       std::map<unsigned int, SiStripDetSummary::Values> map = summaryNoise.getCounts();
1281       //=========================
1282 
1283       TCanvas canvas("Partition summary", "partition summary", 1200, 1000);
1284       canvas.cd();
1285       auto h1 = std::unique_ptr<TH1F>(
1286           new TH1F("byRegion",
1287                    Form("Average by partition of %s SiStrip Noise per module;;average SiStrip Noise %s [ADC counts]",
1288                         estimatorType(est).c_str(),
1289                         estimatorType(est).c_str()),
1290                    map.size(),
1291                    0.,
1292                    map.size()));
1293       h1->SetStats(false);
1294       canvas.SetBottomMargin(0.18);
1295       canvas.SetLeftMargin(0.17);
1296       canvas.SetRightMargin(0.05);
1297       canvas.Modified();
1298 
1299       std::vector<int> boundaries;
1300       unsigned int iBin = 0;
1301 
1302       std::string detector;
1303       std::string currentDetector;
1304 
1305       for (const auto& element : map) {
1306         iBin++;
1307         int count = element.second.count;
1308         double mean = (element.second.mean) / count;
1309         double rms = (element.second.rms) / count - mean * mean;
1310 
1311         if (rms <= 0)
1312           rms = 0;
1313         else
1314           rms = sqrt(rms);
1315 
1316         if (currentDetector.empty())
1317           currentDetector = "TIB";
1318 
1319         switch ((element.first) / 1000) {
1320           case 1:
1321             detector = "TIB";
1322             break;
1323           case 2:
1324             detector = "TOB";
1325             break;
1326           case 3:
1327             detector = "TEC";
1328             break;
1329           case 4:
1330             detector = "TID";
1331             break;
1332         }
1333 
1334         h1->SetBinContent(iBin, mean);
1335         h1->GetXaxis()->SetBinLabel(iBin, SiStripPI::regionType(element.first).second);
1336         h1->GetXaxis()->LabelsOption("v");
1337 
1338         if (detector != currentDetector) {
1339           boundaries.push_back(iBin);
1340           currentDetector = detector;
1341         }
1342       }
1343 
1344       h1->SetMarkerStyle(20);
1345       h1->SetMarkerSize(1);
1346       h1->SetMaximum(h1->GetMaximum() * 1.1);
1347       h1->Draw("HIST");
1348       h1->Draw("Psame");
1349 
1350       canvas.Update();
1351 
1352       TLine l[boundaries.size()];
1353       unsigned int i = 0;
1354 
1355       for (const auto& line : boundaries) {
1356         l[i] = TLine(h1->GetBinLowEdge(line), canvas.GetUymin(), h1->GetBinLowEdge(line), canvas.GetUymax());
1357         l[i].SetLineWidth(1);
1358         l[i].SetLineStyle(9);
1359         l[i].SetLineColor(2);
1360         l[i].Draw("same");
1361         i++;
1362       }
1363 
1364       TLegend legend = TLegend(0.52, 0.82, 0.95, 0.9);
1365       legend.SetHeader((std::get<1>(iov)).c_str(), "C");  // option "C" allows to center the header
1366       legend.AddEntry(h1.get(), ("IOV: " + std::to_string(std::get<0>(iov))).c_str(), "PL");
1367       legend.SetTextSize(0.025);
1368       legend.Draw("same");
1369 
1370       std::string fileName(m_imageFileName);
1371       canvas.SaveAs(fileName.c_str());
1372 
1373       return true;
1374     }
1375 
1376   private:
1377     TrackerTopology m_trackerTopo;
1378   };
1379 
1380   typedef SiStripNoiseByRegion<SiStripPI::mean> SiStripNoiseMeanByRegion;
1381   typedef SiStripNoiseByRegion<SiStripPI::min> SiStripNoiseMinByRegion;
1382   typedef SiStripNoiseByRegion<SiStripPI::max> SiStripNoiseMaxByRegion;
1383   typedef SiStripNoiseByRegion<SiStripPI::rms> SiStripNoiseRMSByRegion;
1384 
1385   /************************************************
1386   SiStrip Noise Comparator
1387   *************************************************/
1388 
1389   template <SiStripPI::estimator est, int ntags, IOVMultiplicity nIOVs>
1390   class SiStripNoiseComparatorByRegionBase : public PlotImage<SiStripNoises, nIOVs, ntags> {
1391   public:
1392     SiStripNoiseComparatorByRegionBase()
1393         : PlotImage<SiStripNoises, nIOVs, ntags>("SiStrip Noise " + estimatorType(est) + " comparator by Region"),
1394           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
1395               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath())} {}
1396 
1397     bool fill() override {
1398       // trick to deal with the multi-ioved tag and two tag case at the same time
1399       auto theIOVs = PlotBase::getTag<0>().iovs;
1400       auto tagname1 = PlotBase::getTag<0>().name;
1401       std::string tagname2 = "";
1402       auto firstiov = theIOVs.front();
1403       SiStripPI::MetaData lastiov;
1404 
1405       // we don't support (yet) comparison with more than 2 tags
1406       assert(this->m_plotAnnotations.ntags < 3);
1407 
1408       if (this->m_plotAnnotations.ntags == 2) {
1409         auto tag2iovs = PlotBase::getTag<1>().iovs;
1410         tagname2 = PlotBase::getTag<1>().name;
1411         lastiov = tag2iovs.front();
1412       } else {
1413         lastiov = theIOVs.back();
1414       }
1415 
1416       std::shared_ptr<SiStripNoises> f_payload = this->fetchPayload(std::get<1>(firstiov));
1417       std::shared_ptr<SiStripNoises> l_payload = this->fetchPayload(std::get<1>(lastiov));
1418 
1419       SiStripDetSummary f_summaryNoise{&m_trackerTopo};
1420       SiStripDetSummary l_summaryNoise{&m_trackerTopo};
1421 
1422       SiStripPI::fillNoiseDetSummary(f_summaryNoise, f_payload, est);
1423       SiStripPI::fillNoiseDetSummary(l_summaryNoise, l_payload, est);
1424 
1425       std::map<unsigned int, SiStripDetSummary::Values> f_map = f_summaryNoise.getCounts();
1426       std::map<unsigned int, SiStripDetSummary::Values> l_map = l_summaryNoise.getCounts();
1427 
1428       //=========================
1429       TCanvas canvas("Partition summary", "partition summary", 1200, 1000);
1430       canvas.cd();
1431 
1432       auto hfirst = std::unique_ptr<TH1F>(
1433           new TH1F("f_byRegion",
1434                    Form("Average by partition of %s SiStrip Noise per module;;average SiStrip Noise %s [ADC counts]",
1435                         estimatorType(est).c_str(),
1436                         estimatorType(est).c_str()),
1437                    f_map.size(),
1438                    0.,
1439                    f_map.size()));
1440       hfirst->SetStats(false);
1441 
1442       auto hlast = std::unique_ptr<TH1F>(
1443           new TH1F("l_byRegion",
1444                    Form("Average by partition of %s SiStrip Noise per module;;average SiStrip Noise %s [ADC counts]",
1445                         estimatorType(est).c_str(),
1446                         estimatorType(est).c_str()),
1447                    l_map.size(),
1448                    0.,
1449                    l_map.size()));
1450       hlast->SetStats(false);
1451 
1452       canvas.SetBottomMargin(0.18);
1453       canvas.SetLeftMargin(0.17);
1454       canvas.SetRightMargin(0.05);
1455       canvas.Modified();
1456 
1457       std::vector<int> boundaries;
1458       unsigned int iBin = 0;
1459 
1460       std::string detector;
1461       std::string currentDetector;
1462 
1463       for (const auto& element : f_map) {
1464         iBin++;
1465         int count = element.second.count;
1466         double mean = (element.second.mean) / count;
1467         double rms = (element.second.rms) / count - mean * mean;
1468 
1469         if (rms <= 0)
1470           rms = 0;
1471         else
1472           rms = sqrt(rms);
1473 
1474         if (currentDetector.empty())
1475           currentDetector = "TIB";
1476 
1477         switch ((element.first) / 1000) {
1478           case 1:
1479             detector = "TIB";
1480             break;
1481           case 2:
1482             detector = "TOB";
1483             break;
1484           case 3:
1485             detector = "TEC";
1486             break;
1487           case 4:
1488             detector = "TID";
1489             break;
1490         }
1491 
1492         hfirst->SetBinContent(iBin, mean);
1493         //hfirst->SetBinError(iBin,rms);
1494         hfirst->GetXaxis()->SetBinLabel(iBin, SiStripPI::regionType(element.first).second);
1495         hfirst->GetXaxis()->LabelsOption("v");
1496 
1497         if (detector != currentDetector) {
1498           boundaries.push_back(iBin);
1499           currentDetector = detector;
1500         }
1501       }
1502 
1503       // second payload
1504       // reset the counter
1505       iBin = 0;
1506 
1507       for (const auto& element : l_map) {
1508         iBin++;
1509         int count = element.second.count;
1510         double mean = (element.second.mean) / count;
1511         double rms = (element.second.rms) / count - mean * mean;
1512 
1513         if (rms <= 0)
1514           rms = 0;
1515         else
1516           rms = sqrt(rms);
1517 
1518         hlast->SetBinContent(iBin, mean);
1519         //hlast->SetBinError(iBin,rms);
1520         hlast->GetXaxis()->SetBinLabel(iBin, SiStripPI::regionType(element.first).second);
1521         hlast->GetXaxis()->LabelsOption("v");
1522       }
1523 
1524       float theMax = (hfirst->GetMaximum() > hlast->GetMaximum()) ? hfirst->GetMaximum() : hlast->GetMaximum();
1525       float theMin = (hfirst->GetMinimum() < hlast->GetMinimum()) ? hfirst->GetMinimum() : hlast->GetMinimum();
1526 
1527       hfirst->SetMarkerStyle(20);
1528       hfirst->SetMarkerSize(1);
1529       hfirst->GetYaxis()->SetTitleOffset(1.3);
1530       hfirst->GetYaxis()->SetRangeUser(theMin * 0.9, theMax * 1.1);
1531       hfirst->Draw("HIST");
1532       hfirst->Draw("Psame");
1533 
1534       hlast->SetMarkerStyle(21);
1535       hlast->SetMarkerSize(1);
1536       hlast->SetMarkerColor(kBlue);
1537       hlast->SetLineColor(kBlue);
1538       hlast->GetYaxis()->SetTitleOffset(1.3);
1539       hlast->GetYaxis()->SetRangeUser(theMin * 0.9, theMax * 1.1);
1540       hlast->Draw("HISTsame");
1541       hlast->Draw("Psame");
1542 
1543       canvas.Update();
1544 
1545       TLine l[boundaries.size()];
1546       unsigned int i = 0;
1547 
1548       for (const auto& line : boundaries) {
1549         l[i] = TLine(hfirst->GetBinLowEdge(line), canvas.GetUymin(), hfirst->GetBinLowEdge(line), canvas.GetUymax());
1550         l[i].SetLineWidth(1);
1551         l[i].SetLineStyle(9);
1552         l[i].SetLineColor(2);
1553         l[i].Draw("same");
1554         i++;
1555       }
1556 
1557       TLegend legend = TLegend(0.52, 0.82, 0.95, 0.9);
1558       legend.SetHeader(("SiStrip Noise " + estimatorType(est) + " by region").c_str(),
1559                        "C");  // option "C" allows to center the header
1560       legend.AddEntry(hfirst.get(), ("IOV: " + std::to_string(std::get<0>(firstiov))).c_str(), "PL");
1561       legend.AddEntry(hlast.get(), ("IOV: " + std::to_string(std::get<0>(lastiov))).c_str(), "PL");
1562       legend.SetTextSize(0.025);
1563       legend.Draw("same");
1564 
1565       std::string fileName(this->m_imageFileName);
1566       canvas.SaveAs(fileName.c_str());
1567 
1568       return true;
1569     }
1570 
1571   private:
1572     TrackerTopology m_trackerTopo;
1573   };
1574 
1575   template <SiStripPI::estimator est>
1576   using SiStripNoiseComparatorByRegionSingleTag = SiStripNoiseComparatorByRegionBase<est, 1, MULTI_IOV>;
1577 
1578   template <SiStripPI::estimator est>
1579   using SiStripNoiseComparatorByRegionTwoTags = SiStripNoiseComparatorByRegionBase<est, 2, SINGLE_IOV>;
1580 
1581   typedef SiStripNoiseComparatorByRegionSingleTag<SiStripPI::mean> SiStripNoiseComparatorMeanByRegionSingleTag;
1582   typedef SiStripNoiseComparatorByRegionSingleTag<SiStripPI::min> SiStripNoiseComparatorMinByRegionSingleTag;
1583   typedef SiStripNoiseComparatorByRegionSingleTag<SiStripPI::max> SiStripNoiseComparatorMaxByRegionSingleTag;
1584   typedef SiStripNoiseComparatorByRegionSingleTag<SiStripPI::rms> SiStripNoiseComparatorRMSByRegionSingleTag;
1585 
1586   typedef SiStripNoiseComparatorByRegionTwoTags<SiStripPI::mean> SiStripNoiseComparatorMeanByRegionTwoTags;
1587   typedef SiStripNoiseComparatorByRegionTwoTags<SiStripPI::min> SiStripNoiseComparatorMinByRegionTwoTags;
1588   typedef SiStripNoiseComparatorByRegionTwoTags<SiStripPI::max> SiStripNoiseComparatorMaxByRegionTwoTags;
1589   typedef SiStripNoiseComparatorByRegionTwoTags<SiStripPI::rms> SiStripNoiseComparatorRMSByRegionTwoTags;
1590 
1591   /************************************************
1592     Noise linearity
1593   *************************************************/
1594   class SiStripNoiseLinearity : public PlotImage<SiStripNoises, SINGLE_IOV> {
1595   public:
1596     SiStripNoiseLinearity()
1597         : PlotImage<SiStripNoises, SINGLE_IOV>("Linearity of Strip Noise as a fuction of strip length") {}
1598 
1599     bool fill() override {
1600       auto tag = PlotBase::getTag<0>();
1601       auto iov = tag.iovs.front();
1602       std::shared_ptr<SiStripNoises> payload = fetchPayload(std::get<1>(iov));
1603 
1604       const auto detInfo =
1605           SiStripDetInfoFileReader::read(edm::FileInPath(SiStripDetInfoFileReader::kDefaultFile).fullPath());
1606 
1607       std::vector<uint32_t> detid;
1608       payload->getDetIds(detid);
1609 
1610       std::map<float, std::tuple<int, float, float>> noisePerStripLength;
1611 
1612       for (const auto& d : detid) {
1613         SiStripNoises::Range range = payload->getRange(d);
1614         for (int it = 0; it < (range.second - range.first) * 8 / 9; ++it) {
1615           auto noise = payload->getNoise(it, range);
1616           //to be used to fill the histogram
1617           float stripL = detInfo.getNumberOfApvsAndStripLength(d).second;
1618           std::get<0>(noisePerStripLength[stripL]) += 1;
1619           std::get<1>(noisePerStripLength[stripL]) += noise;
1620           std::get<2>(noisePerStripLength[stripL]) += (noise * noise);
1621         }  // loop over strips
1622       }  // loop over detIds
1623 
1624       TCanvas canvas("Noise linearity", "noise linearity", 1200, 1000);
1625       canvas.cd();
1626 
1627       std::vector<float> x;
1628       x.reserve(noisePerStripLength.size());
1629       std::vector<float> y;
1630       y.reserve(noisePerStripLength.size());
1631       std::vector<float> ex;
1632       ex.reserve(noisePerStripLength.size());
1633       std::vector<float> ey;
1634       ey.reserve(noisePerStripLength.size());
1635 
1636       for (const auto& element : noisePerStripLength) {
1637         x.push_back(element.first);
1638         ex.push_back(0.);
1639         float sum = std::get<1>(element.second);
1640         float sum2 = std::get<2>(element.second);
1641         float nstrips = std::get<0>(element.second);
1642         float mean = sum / nstrips;
1643         float rms = (sum2 / nstrips - mean * mean) > 0. ? sqrt(sum2 / nstrips - mean * mean) : 0.;
1644         y.push_back(mean);
1645         ey.push_back(rms);
1646         //std::cout<<" strip lenght: " << element.first << " avg noise=" << mean <<" +/-" << rms << std::endl;
1647       }
1648 
1649       auto graph = std::make_unique<TGraphErrors>(noisePerStripLength.size(), &x[0], &y[0], &ex[0], &ey[0]);
1650       graph->SetTitle("SiStrip Noise Linearity");
1651       graph->GetXaxis()->SetTitle("Strip length [cm]");
1652       graph->GetYaxis()->SetTitle("Average Strip Noise [ADC counts]");
1653       graph->SetMarkerColor(kBlue);
1654       graph->SetMarkerStyle(20);
1655       graph->SetMarkerSize(1.5);
1656       canvas.SetBottomMargin(0.13);
1657       canvas.SetLeftMargin(0.17);
1658       canvas.SetTopMargin(0.08);
1659       canvas.SetRightMargin(0.05);
1660       canvas.Modified();
1661       canvas.cd();
1662 
1663       graph->GetXaxis()->CenterTitle(true);
1664       graph->GetYaxis()->CenterTitle(true);
1665       graph->GetXaxis()->SetTitleFont(42);
1666       graph->GetYaxis()->SetTitleFont(42);
1667       graph->GetXaxis()->SetTitleSize(0.05);
1668       graph->GetYaxis()->SetTitleSize(0.05);
1669       graph->GetXaxis()->SetTitleOffset(1.1);
1670       graph->GetYaxis()->SetTitleOffset(1.3);
1671       graph->GetXaxis()->SetLabelFont(42);
1672       graph->GetYaxis()->SetLabelFont(42);
1673       graph->GetYaxis()->SetLabelSize(.05);
1674       graph->GetXaxis()->SetLabelSize(.05);
1675 
1676       graph->Draw("AP");
1677       graph->Fit("pol1");
1678       //Access the fit resuts
1679       TF1* f1 = graph->GetFunction("pol1");
1680       f1->SetLineWidth(2);
1681       f1->SetLineColor(kBlue);
1682       f1->Draw("same");
1683 
1684       auto fits = std::make_unique<TPaveText>(0.2, 0.72, 0.6, 0.9, "NDC");
1685       char buffer[255];
1686       sprintf(buffer, "fit function: p_{0} + p_{1} * l_{strip}");
1687       fits->AddText(buffer);
1688       sprintf(buffer, "p_{0} : %5.2f [ADC counts]", f1->GetParameter(0));
1689       fits->AddText(buffer);
1690       sprintf(buffer, "p_{1} : %5.2f [ADC counts/cm]", f1->GetParameter(1));
1691       fits->AddText(buffer);
1692       sprintf(buffer, "#chi^{2}/ndf = %5.2f / %i ", f1->GetChisquare(), f1->GetNDF());
1693       fits->AddText(buffer);
1694       fits->SetTextFont(42);
1695       fits->SetTextColor(kBlue);
1696       fits->SetFillColor(0);
1697       fits->SetTextSize(0.03);
1698       fits->SetBorderSize(1);
1699       fits->SetLineColor(kBlue);
1700       fits->SetMargin(0.05);
1701       fits->SetTextAlign(12);
1702       fits->Draw();
1703 
1704       std::string fileName(m_imageFileName);
1705       canvas.SaveAs(fileName.c_str());
1706 
1707       delete f1;
1708       return true;
1709     }
1710   };
1711 
1712   /************************************************
1713    template Noise history per subdetector
1714   *************************************************/
1715 
1716   template <StripSubdetector::SubDetector sub>
1717   class NoiseHistory : public HistoryPlot<SiStripNoises, std::pair<double, double>> {
1718   public:
1719     NoiseHistory()
1720         : HistoryPlot<SiStripNoises, std::pair<double, double>>(
1721               "Average " + SiStripPI::getStringFromSubdet(sub) + " noise vs run number",
1722               "average " + SiStripPI::getStringFromSubdet(sub) + " Noise") {}
1723 
1724     std::pair<double, double> getFromPayload(SiStripNoises& payload) override {
1725       std::vector<uint32_t> detid;
1726       payload.getDetIds(detid);
1727 
1728       int nStrips = 0;
1729       float sum = 0., sum2 = 0.;
1730 
1731       for (const auto& d : detid) {
1732         int subid = DetId(d).subdetId();
1733         if (subid != sub)
1734           continue;
1735         SiStripNoises::Range range = payload.getRange(d);
1736         for (int it = 0; it < (range.second - range.first) * 8 / 9; ++it) {
1737           nStrips++;
1738           auto noise = payload.getNoise(it, range);
1739           sum += noise;
1740           sum2 += (noise * noise);
1741         }  // loop on strips
1742       }  // loop on detIds
1743 
1744       float mean = sum / nStrips;
1745       float rms = (sum2 / nStrips - mean * mean) > 0. ? sqrt(sum2 / nStrips - mean * mean) : 0.;
1746 
1747       return std::make_pair(mean, rms);
1748 
1749     }  // close getFromPayload
1750   };
1751 
1752   typedef NoiseHistory<StripSubdetector::TIB> TIBNoiseHistory;
1753   typedef NoiseHistory<StripSubdetector::TOB> TOBNoiseHistory;
1754   typedef NoiseHistory<StripSubdetector::TID> TIDNoiseHistory;
1755   typedef NoiseHistory<StripSubdetector::TEC> TECNoiseHistory;
1756 
1757   /************************************************
1758    template Noise run history  per subdetector
1759   *************************************************/
1760 
1761   template <StripSubdetector::SubDetector sub>
1762   class NoiseRunHistory : public RunHistoryPlot<SiStripNoises, std::pair<double, double>> {
1763   public:
1764     NoiseRunHistory()
1765         : RunHistoryPlot<SiStripNoises, std::pair<double, double>>(
1766               "Average " + SiStripPI::getStringFromSubdet(sub) + " noise vs run number",
1767               "average " + SiStripPI::getStringFromSubdet(sub) + " Noise") {}
1768 
1769     std::pair<double, double> getFromPayload(SiStripNoises& payload) override {
1770       std::vector<uint32_t> detid;
1771       payload.getDetIds(detid);
1772 
1773       int nStrips = 0;
1774       float sum = 0., sum2 = 0.;
1775 
1776       for (const auto& d : detid) {
1777         int subid = DetId(d).subdetId();
1778         if (subid != sub)
1779           continue;
1780         SiStripNoises::Range range = payload.getRange(d);
1781         for (int it = 0; it < (range.second - range.first) * 8 / 9; ++it) {
1782           nStrips++;
1783           auto noise = payload.getNoise(it, range);
1784           sum += noise;
1785           sum2 += (noise * noise);
1786         }  // loop on strips
1787       }  // loop on detIds
1788 
1789       float mean = sum / nStrips;
1790       float rms = (sum2 / nStrips - mean * mean) > 0. ? sqrt(sum2 / nStrips - mean * mean) : 0.;
1791 
1792       return std::make_pair(mean, rms);
1793 
1794     }  // close getFromPayload
1795   };
1796 
1797   typedef NoiseRunHistory<StripSubdetector::TIB> TIBNoiseRunHistory;
1798   typedef NoiseRunHistory<StripSubdetector::TOB> TOBNoiseRunHistory;
1799   typedef NoiseRunHistory<StripSubdetector::TID> TIDNoiseRunHistory;
1800   typedef NoiseRunHistory<StripSubdetector::TEC> TECNoiseRunHistory;
1801 
1802   /************************************************
1803    template Noise Time history per subdetector
1804   *************************************************/
1805 
1806   template <StripSubdetector::SubDetector sub>
1807   class NoiseTimeHistory : public TimeHistoryPlot<SiStripNoises, std::pair<double, double>> {
1808   public:
1809     NoiseTimeHistory()
1810         : TimeHistoryPlot<SiStripNoises, std::pair<double, double>>(
1811               "Average " + SiStripPI::getStringFromSubdet(sub) + " noise vs run number",
1812               "average " + SiStripPI::getStringFromSubdet(sub) + " Noise") {}
1813 
1814     std::pair<double, double> getFromPayload(SiStripNoises& payload) override {
1815       std::vector<uint32_t> detid;
1816       payload.getDetIds(detid);
1817 
1818       int nStrips = 0;
1819       float sum = 0., sum2 = 0.;
1820 
1821       for (const auto& d : detid) {
1822         int subid = DetId(d).subdetId();
1823         if (subid != sub)
1824           continue;
1825         SiStripNoises::Range range = payload.getRange(d);
1826         for (int it = 0; it < (range.second - range.first) * 8 / 9; ++it) {
1827           nStrips++;
1828           auto noise = payload.getNoise(it, range);
1829           sum += noise;
1830           sum2 += (noise * noise);
1831         }  // loop on strips
1832       }  // loop on detIds
1833 
1834       float mean = sum / nStrips;
1835       float rms = (sum2 / nStrips - mean * mean) > 0. ? sqrt(sum2 / nStrips - mean * mean) : 0.;
1836 
1837       return std::make_pair(mean, rms);
1838 
1839     }  // close getFromPayload
1840   };
1841 
1842   typedef NoiseTimeHistory<StripSubdetector::TIB> TIBNoiseTimeHistory;
1843   typedef NoiseTimeHistory<StripSubdetector::TOB> TOBNoiseTimeHistory;
1844   typedef NoiseTimeHistory<StripSubdetector::TID> TIDNoiseTimeHistory;
1845   typedef NoiseTimeHistory<StripSubdetector::TEC> TECNoiseTimeHistory;
1846 
1847   /************************************************
1848    template Noise run history  per layer
1849   *************************************************/
1850   template <StripSubdetector::SubDetector sub>
1851   class NoiseLayerRunHistory : public PlotImage<SiStripNoises, MULTI_IOV> {
1852   public:
1853     NoiseLayerRunHistory()
1854         : PlotImage<SiStripNoises, MULTI_IOV>("SiStrip Noise values comparison"),
1855           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
1856               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath())} {}
1857 
1858     bool fill() override {
1859       auto tag = PlotBase::getTag<0>();
1860       auto sorted_iovs = tag.iovs;
1861 
1862       // make absolute sure the IOVs are sortd by since
1863       std::sort(begin(sorted_iovs), end(sorted_iovs), [](auto const& t1, auto const& t2) {
1864         return std::get<0>(t1) < std::get<0>(t2);
1865       });
1866 
1867       std::unordered_map<int, std::vector<float>> noises_avg;
1868       std::unordered_map<int, std::vector<float>> noises_err;
1869       std::vector<float> runs;
1870       std::vector<float> runs_err;
1871 
1872       for (auto const& iov : sorted_iovs) {
1873         std::unordered_map<int, std::vector<float>> noises;  //map with noises per layer
1874 
1875         std::shared_ptr<SiStripNoises> payload = fetchPayload(std::get<1>(iov));
1876         unsigned int run = std::get<0>(iov);
1877         runs.push_back(run);
1878         runs_err.push_back(0);
1879 
1880         if (payload.get()) {
1881           std::vector<uint32_t> detid;
1882           payload->getDetIds(detid);
1883 
1884           for (const auto& d : detid) {
1885             int subid = DetId(d).subdetId();
1886             int layer = -1;
1887             if (subid != sub)
1888               continue;
1889             if (subid == StripSubdetector::TIB) {
1890               layer = m_trackerTopo.tibLayer(d);
1891             } else if (subid == StripSubdetector::TOB) {
1892               layer = m_trackerTopo.tobLayer(d);
1893             } else if (subid == StripSubdetector::TID) {
1894               layer = m_trackerTopo.tidWheel(d);
1895             } else if (subid == StripSubdetector::TEC) {
1896               layer = m_trackerTopo.tecWheel(d);
1897             }
1898 
1899             SiStripNoises::Range range = payload->getRange(d);
1900             for (int it = 0; it < (range.second - range.first) * 8 / 9; ++it) {
1901               auto noise = payload->getNoise(it, range);
1902               if (noises.find(layer) == noises.end())
1903                 noises.emplace(layer, std::vector<float>{});
1904               noises[layer].push_back(noise);
1905             }  // loop on strips
1906           }  // loop on detIds
1907 
1908           for (auto& entry : noises) {
1909             double sum = std::accumulate(entry.second.begin(), entry.second.end(), 0.0);
1910             double mean = sum / entry.second.size();
1911 
1912             //double sq_sum = std::inner_product(entry.second.begin(), entry.second.end(), entry.second.begin(), 0.0);
1913             //double stdev = std::sqrt(sq_sum / entry.second.size() - mean * mean);
1914 
1915             if (noises_avg.find(entry.first) == noises_avg.end())
1916               noises_avg.emplace(entry.first, std::vector<float>{});
1917             noises_avg[entry.first].push_back(mean);
1918 
1919             if (noises_err.find(entry.first) == noises_err.end())
1920               noises_err.emplace(entry.first, std::vector<float>{});
1921             noises_err[entry.first].push_back(0);
1922           }  //get
1923         }  //run on iov
1924       }
1925       TCanvas canvas("Partition summary", "partition summary", 2000, 1000);
1926       canvas.cd();
1927       canvas.SetBottomMargin(0.11);
1928       canvas.SetLeftMargin(0.13);
1929       canvas.SetRightMargin(0.05);
1930       canvas.Modified();
1931 
1932       TLegend legend = TLegend(0.73, 0.13, 0.89, 0.43);
1933       //legend.SetHeader("Layers","C"); // option "C" allows to center the header
1934       legend.SetTextSize(0.03);
1935 
1936       std::unique_ptr<TGraphErrors> graph[noises_avg.size()];
1937 
1938       int colors[18] = {1, 2, 3, 4, 5, 6, 7, 8, 9, 20, 30, 40, 42, 46, 48, 32, 36, 38};
1939 
1940       int el = 0;
1941 
1942       for (auto& entry : noises_avg) {
1943         graph[el] = std::make_unique<TGraphErrors>(
1944             runs.size(), &runs[0], &(entry.second[0]), &runs_err[0], &(noises_err[entry.first][0]));
1945         char title[100];
1946         char name[100];
1947         snprintf(name, sizeof(name), "gr%d", entry.first);
1948         graph[el]->SetName(name);
1949 
1950         if (sub == StripSubdetector::TIB) {
1951           snprintf(title, sizeof(title), "SiStrip avg noise per layer -- TIB");
1952           graph[el]->GetYaxis()->SetTitle("Average Noise per Layer [ADC counts]");
1953         } else if (sub == StripSubdetector::TOB) {
1954           snprintf(title, sizeof(title), "SiStrip avg noise per layer -- TOB");
1955           graph[el]->GetYaxis()->SetTitle("Average Noise per Layer [ADC counts]");
1956         } else if (sub == StripSubdetector::TID) {
1957           snprintf(title, sizeof(title), "SiStrip avg noise per disk -- TID");
1958           graph[el]->GetYaxis()->SetTitle("Average Noise per Disk [ADC counts]");
1959         } else if (sub == StripSubdetector::TEC) {
1960           snprintf(title, sizeof(title), "SiStrip avg noise per disk -- TEC");
1961           graph[el]->GetYaxis()->SetTitle("Average Noise per Disk [ADC counts]");
1962         }
1963 
1964         graph[el]->SetTitle(title);
1965         graph[el]->GetXaxis()->SetTitle("run");
1966         graph[el]->SetMarkerColor(colors[el]);
1967         graph[el]->SetMarkerStyle(20);
1968         graph[el]->SetMarkerSize(1.5);
1969         graph[el]->GetXaxis()->CenterTitle(true);
1970         graph[el]->GetYaxis()->CenterTitle(true);
1971         graph[el]->GetXaxis()->SetTitleFont(42);
1972         graph[el]->GetYaxis()->SetTitleFont(42);
1973         graph[el]->GetXaxis()->SetTitleSize(0.05);
1974         graph[el]->GetYaxis()->SetTitleSize(0.05);
1975         graph[el]->GetXaxis()->SetTitleOffset(1.1);
1976         graph[el]->GetYaxis()->SetTitleOffset(1.3);
1977         graph[el]->GetXaxis()->SetLabelFont(42);
1978         graph[el]->GetYaxis()->SetLabelFont(42);
1979         graph[el]->GetYaxis()->SetLabelSize(.05);
1980         graph[el]->GetXaxis()->SetLabelSize(.05);
1981         graph[el]->SetMinimum(3);
1982         graph[el]->SetMaximum(7.5);
1983 
1984         if (el == 0)
1985           graph[el]->Draw("AP");
1986         else
1987           graph[el]->Draw("P");
1988 
1989         if (sub == StripSubdetector::TIB) {
1990           legend.AddEntry(name, ("layer " + std::to_string(entry.first)).c_str(), "lep");
1991         } else if (sub == StripSubdetector::TOB) {
1992           legend.AddEntry(name, ("layer " + std::to_string(entry.first)).c_str(), "lep");
1993         } else if (sub == StripSubdetector::TID) {
1994           legend.AddEntry(name, ("disk " + std::to_string(entry.first)).c_str(), "lep");
1995         } else if (sub == StripSubdetector::TEC) {
1996           legend.AddEntry(name, ("disk " + std::to_string(entry.first)).c_str(), "lep");
1997         }
1998 
1999         if (el == 0)
2000           legend.Draw();
2001         else
2002           legend.Draw("same");
2003         el++;
2004       }
2005       //canvas.BuildLegend();
2006       std::string fileName(m_imageFileName);
2007       canvas.SaveAs(fileName.c_str());
2008       return true;
2009     }
2010 
2011   private:
2012     TrackerTopology m_trackerTopo;
2013   };
2014 
2015   typedef NoiseLayerRunHistory<StripSubdetector::TIB> TIBNoiseLayerRunHistory;
2016   typedef NoiseLayerRunHistory<StripSubdetector::TOB> TOBNoiseLayerRunHistory;
2017   typedef NoiseLayerRunHistory<StripSubdetector::TID> TIDNoiseLayerRunHistory;
2018   typedef NoiseLayerRunHistory<StripSubdetector::TEC> TECNoiseLayerRunHistory;
2019 
2020 }  // namespace
2021 
2022 PAYLOAD_INSPECTOR_MODULE(SiStripNoises) {
2023   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseConsistencyCheck);
2024   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseCompareByPartitionSingleTag);
2025   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseDiffByPartitionSingleTag);
2026   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseCompareByPartitionTwoTags);
2027   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseDiffByPartitionTwoTags);
2028   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseCorrelationByPartition);
2029   PAYLOAD_INSPECTOR_CLASS(SiStripNoisesTest);
2030   PAYLOAD_INSPECTOR_CLASS(SiStripNoisePerDetId);
2031   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseValue);
2032   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseValuePerDetId);
2033   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseValuePerStrip);
2034   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseValuePerAPV);
2035   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseValuePerModule);
2036   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseValueComparisonSingleTag);
2037   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseValueComparisonTwoTags);
2038   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseValueComparisonPerStripSingleTag);
2039   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseValueComparisonPerAPVSingleTag);
2040   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseValueComparisonPerModuleSingleTag);
2041   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseValueComparisonPerStripTwoTags);
2042   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseValueComparisonPerAPVTwoTags);
2043   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseValueComparisonPerModuleTwoTags);
2044   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseMin_TrackerMap);
2045   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseMax_TrackerMap);
2046   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseMean_TrackerMap);
2047   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseRMS_TrackerMap);
2048   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseMeanByRegion);
2049   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseMinByRegion);
2050   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseMaxByRegion);
2051   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseRMSByRegion);
2052   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseComparatorMeanByRegionSingleTag);
2053   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseComparatorMinByRegionSingleTag);
2054   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseComparatorMaxByRegionSingleTag);
2055   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseComparatorRMSByRegionSingleTag);
2056   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseComparatorMeanByRegionTwoTags);
2057   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseComparatorMinByRegionTwoTags);
2058   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseComparatorMaxByRegionTwoTags);
2059   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseComparatorRMSByRegionTwoTags);
2060   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseMin_RatioWithPreviousIOVTrackerMapSingleTag);
2061   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseMax_RatioWithPreviousIOVTrackerMapSingleTag);
2062   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseMean_RatioWithPreviousIOVTrackerMapSingleTag);
2063   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseRms_RatioWithPreviousIOVTrackerMapSingleTag);
2064   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseMin_RatioWithPreviousIOVTrackerMapTwoTags);
2065   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseMax_RatioWithPreviousIOVTrackerMapTwoTags);
2066   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseMean_RatioWithPreviousIOVTrackerMapTwoTags);
2067   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseRms_RatioWithPreviousIOVTrackerMapTwoTags);
2068   PAYLOAD_INSPECTOR_CLASS(SiStripNoiseLinearity);
2069   PAYLOAD_INSPECTOR_CLASS(TIBNoiseHistory);
2070   PAYLOAD_INSPECTOR_CLASS(TOBNoiseHistory);
2071   PAYLOAD_INSPECTOR_CLASS(TIDNoiseHistory);
2072   PAYLOAD_INSPECTOR_CLASS(TECNoiseHistory);
2073   PAYLOAD_INSPECTOR_CLASS(TIBNoiseRunHistory);
2074   PAYLOAD_INSPECTOR_CLASS(TOBNoiseRunHistory);
2075   PAYLOAD_INSPECTOR_CLASS(TIDNoiseRunHistory);
2076   PAYLOAD_INSPECTOR_CLASS(TECNoiseRunHistory);
2077   PAYLOAD_INSPECTOR_CLASS(TIBNoiseLayerRunHistory);
2078   PAYLOAD_INSPECTOR_CLASS(TOBNoiseLayerRunHistory);
2079   PAYLOAD_INSPECTOR_CLASS(TIDNoiseLayerRunHistory);
2080   PAYLOAD_INSPECTOR_CLASS(TECNoiseLayerRunHistory);
2081   PAYLOAD_INSPECTOR_CLASS(TIBNoiseTimeHistory);
2082   PAYLOAD_INSPECTOR_CLASS(TOBNoiseTimeHistory);
2083   PAYLOAD_INSPECTOR_CLASS(TIDNoiseTimeHistory);
2084   PAYLOAD_INSPECTOR_CLASS(TECNoiseTimeHistory);
2085 }