Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-12-01 23:40:01

0001 /*!
0002   \file SiStripApvGains_PayloadInspector
0003   \Payload Inspector Plugin for SiStrip Gain
0004   \author M. Musich
0005   \version $Revision: 1.0 $
0006   \date $Date: 2017/07/02 17: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/SiStripApvGain.h"
0015 #include "DataFormats/DetId/interface/DetId.h"
0016 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
0017 #include "CondFormats/SiStripObjects/interface/SiStripDetSummary.h"
0018 
0019 // needed for the tracker map
0020 #include "CommonTools/TrackerMap/interface/TrackerMap.h"
0021 
0022 // auxilliary functions
0023 #include "CondCore/SiStripPlugins/interface/SiStripPayloadInspectorHelper.h"
0024 #include "CalibTracker/StandaloneTrackerTopology/interface/StandaloneTrackerTopology.h"
0025 #include "SiStripCondObjectRepresent.h"
0026 
0027 #include <memory>
0028 #include <sstream>
0029 #include <iostream>
0030 
0031 // include ROOT
0032 #include "TProfile.h"
0033 #include "TH2F.h"
0034 #include "THStack.h"
0035 #include "TLegend.h"
0036 #include "TCanvas.h"
0037 #include "TLine.h"
0038 #include "TStyle.h"
0039 #include "TLatex.h"
0040 #include "TPave.h"
0041 #include "TPaveStats.h"
0042 
0043 namespace {
0044 
0045   using namespace cond::payloadInspector;
0046 
0047   class SiStripApvGainContainer : public SiStripCondObjectRepresent::SiStripDataContainer<SiStripApvGain, float> {
0048   public:
0049     SiStripApvGainContainer(const std::shared_ptr<SiStripApvGain>& payload,
0050                             const SiStripPI::MetaData& metadata,
0051                             const std::string& tagName)
0052         : SiStripCondObjectRepresent::SiStripDataContainer<SiStripApvGain, float>(payload, metadata, tagName) {
0053       payloadType_ = "SiStripApvGain";
0054       setGranularity(SiStripCondObjectRepresent::PERAPV);
0055     }
0056 
0057     void storeAllValues() override {
0058       std::vector<uint32_t> detid;
0059       payload_->getDetIds(detid);
0060 
0061       for (const auto& d : detid) {
0062         SiStripApvGain::Range range = payload_->getRange(d);
0063         for (int it = 0; it < range.second - range.first; it++) {
0064           // to be used to fill the histogram
0065           SiStripCondData_.fillByPushBack(d, payload_->getApvGain(it, range));
0066         }
0067       }
0068     }
0069   };
0070 
0071   /************************************************
0072     testing the machinery
0073   ************************************************/
0074   class SiStripApvGainTest : public Histogram1D<SiStripApvGain, SINGLE_IOV> {
0075   public:
0076     SiStripApvGainTest()
0077         : Histogram1D<SiStripApvGain, SINGLE_IOV>("SiStrip ApvGain values", "SiStrip ApvGain values", 1, 0.0, 1.) {}
0078 
0079     bool fill() override {
0080       auto tag = PlotBase::getTag<0>();
0081       auto tagname = tag.name;
0082       for (auto const& iov : tag.iovs) {
0083         std::shared_ptr<SiStripApvGain> payload = Base::fetchPayload(std::get<1>(iov));
0084         if (payload.get()) {
0085           SiStripApvGainContainer* objContainer = new SiStripApvGainContainer(payload, iov, tagname);
0086           objContainer->printAll();
0087 
0088         }  // payload
0089       }  // iovs
0090       return true;
0091     }  // fill
0092   };
0093 
0094   class SiStripApvGainByPartition : public PlotImage<SiStripApvGain, SINGLE_IOV> {
0095   public:
0096     SiStripApvGainByPartition() : PlotImage<SiStripApvGain, SINGLE_IOV>("SiStrip ApvGains By Partition") {}
0097 
0098     bool fill() override {
0099       auto tag = PlotBase::getTag<0>();
0100       auto iov = tag.iovs.front();
0101       auto tagname = tag.name;
0102       std::shared_ptr<SiStripApvGain> payload = fetchPayload(std::get<1>(iov));
0103       if (payload.get()) {
0104         SiStripApvGainContainer* objContainer = new SiStripApvGainContainer(payload, iov, tagname);
0105         //objContainer->printAll();
0106 
0107         TCanvas canvas("Partition summary", "partition summary", 1400, 1000);
0108         objContainer->fillByPartition(canvas, 100, 0., 2.);
0109 
0110         std::string fileName(m_imageFileName);
0111         canvas.SaveAs(fileName.c_str());
0112       }  // payload
0113       return true;
0114     }  // fill
0115   };
0116 
0117   template <int ntags, IOVMultiplicity nIOVs>
0118   class SiStripApvGainCompareByPartitionBase : public PlotImage<SiStripApvGain, nIOVs, ntags> {
0119   public:
0120     SiStripApvGainCompareByPartitionBase()
0121         : PlotImage<SiStripApvGain, nIOVs, ntags>("SiStrip Compare ApvGains By Partition") {}
0122 
0123     bool fill() override {
0124       // trick to deal with the multi-ioved tag and two tag case at the same time
0125       auto theIOVs = PlotBase::getTag<0>().iovs;
0126       auto tagname1 = PlotBase::getTag<0>().name;
0127       std::string tagname2{tagname1};
0128       auto firstiov = theIOVs.front();
0129       SiStripPI::MetaData lastiov;
0130 
0131       // we don't support (yet) comparison with more than 2 tags
0132       assert(this->m_plotAnnotations.ntags < 3);
0133 
0134       if (this->m_plotAnnotations.ntags == 2) {
0135         auto tag2iovs = PlotBase::getTag<1>().iovs;
0136         tagname2 = PlotBase::getTag<1>().name;
0137         lastiov = tag2iovs.front();
0138       } else {
0139         lastiov = theIOVs.back();
0140       }
0141 
0142       std::shared_ptr<SiStripApvGain> last_payload = this->fetchPayload(std::get<1>(lastiov));
0143       std::shared_ptr<SiStripApvGain> first_payload = this->fetchPayload(std::get<1>(firstiov));
0144 
0145       SiStripApvGainContainer* l_objContainer = new SiStripApvGainContainer(last_payload, lastiov, tagname1);
0146       SiStripApvGainContainer* f_objContainer = new SiStripApvGainContainer(first_payload, firstiov, tagname2);
0147 
0148       l_objContainer->compare(f_objContainer);
0149 
0150       //l_objContainer->printAll();
0151 
0152       TCanvas canvas("Partition summary", "partition summary", 1400, 1000);
0153       l_objContainer->fillByPartition(canvas, 100, 0.5, 1.5);
0154 
0155       std::string fileName(this->m_imageFileName);
0156       canvas.SaveAs(fileName.c_str());
0157 
0158       return true;
0159     }  // fill
0160   };
0161 
0162   using SiStripApvGainCompareByPartitionSingleTag = SiStripApvGainCompareByPartitionBase<1, MULTI_IOV>;
0163   using SiStripApvGainCompareByPartitionTwoTags = SiStripApvGainCompareByPartitionBase<2, SINGLE_IOV>;
0164 
0165   template <int ntags, IOVMultiplicity nIOVs>
0166   class SiStripApvGainRatioByPartitionBase : public PlotImage<SiStripApvGain, nIOVs, ntags> {
0167   public:
0168     SiStripApvGainRatioByPartitionBase()
0169         : PlotImage<SiStripApvGain, nIOVs, ntags>("SiStrip Ratio ApvGains By Partition") {}
0170 
0171     bool fill() override {
0172       // trick to deal with the multi-ioved tag and two tag case at the same time
0173       auto theIOVs = PlotBase::getTag<0>().iovs;
0174       auto tagname1 = PlotBase::getTag<0>().name;
0175       std::string tagname2{tagname1};
0176       auto firstiov = theIOVs.front();
0177       SiStripPI::MetaData lastiov;
0178 
0179       // we don't support (yet) comparison with more than 2 tags
0180       assert(this->m_plotAnnotations.ntags < 3);
0181 
0182       if (this->m_plotAnnotations.ntags == 2) {
0183         auto tag2iovs = PlotBase::getTag<1>().iovs;
0184         tagname2 = PlotBase::getTag<1>().name;
0185         lastiov = tag2iovs.front();
0186       } else {
0187         lastiov = theIOVs.back();
0188       }
0189 
0190       std::shared_ptr<SiStripApvGain> last_payload = this->fetchPayload(std::get<1>(lastiov));
0191       std::shared_ptr<SiStripApvGain> first_payload = this->fetchPayload(std::get<1>(firstiov));
0192 
0193       SiStripApvGainContainer* l_objContainer = new SiStripApvGainContainer(last_payload, lastiov, tagname1);
0194       SiStripApvGainContainer* f_objContainer = new SiStripApvGainContainer(first_payload, firstiov, tagname2);
0195 
0196       l_objContainer->divide(f_objContainer);
0197 
0198       //l_objContainer->printAll();
0199 
0200       TCanvas canvas("Partition summary", "partition summary", 1400, 1000);
0201       l_objContainer->fillByPartition(canvas, 200, 0.5, 1.5);
0202       //for (int i = 1; i <= 4; i++)
0203       //  canvas.cd(i)->SetLogy();
0204 
0205       std::string fileName(this->m_imageFileName);
0206       canvas.SaveAs(fileName.c_str());
0207 
0208       return true;
0209     }  // fill
0210   };
0211 
0212   using SiStripApvGainRatioByPartitionSingleTag = SiStripApvGainRatioByPartitionBase<1, MULTI_IOV>;
0213   using SiStripApvGainRatioByPartitionTwoTags = SiStripApvGainRatioByPartitionBase<2, SINGLE_IOV>;
0214 
0215   template <int ntags, IOVMultiplicity nIOVs>
0216   class SiStripApvGainDiffByPartitionBase : public PlotImage<SiStripApvGain, nIOVs, ntags> {
0217   public:
0218     SiStripApvGainDiffByPartitionBase()
0219         : PlotImage<SiStripApvGain, nIOVs, ntags>("SiStrip Diff ApvGains By Partition") {}
0220 
0221     bool fill() override {
0222       // trick to deal with the multi-ioved tag and two tag case at the same time
0223       auto theIOVs = PlotBase::getTag<0>().iovs;
0224       auto tagname1 = PlotBase::getTag<0>().name;
0225       std::string tagname2{tagname1};
0226       auto firstiov = theIOVs.front();
0227       SiStripPI::MetaData lastiov;
0228 
0229       // we don't support (yet) comparison with more than 2 tags
0230       assert(this->m_plotAnnotations.ntags < 3);
0231 
0232       if (this->m_plotAnnotations.ntags == 2) {
0233         auto tag2iovs = PlotBase::getTag<1>().iovs;
0234         tagname2 = PlotBase::getTag<1>().name;
0235         lastiov = tag2iovs.front();
0236       } else {
0237         lastiov = theIOVs.back();
0238       }
0239 
0240       std::shared_ptr<SiStripApvGain> last_payload = this->fetchPayload(std::get<1>(lastiov));
0241       std::shared_ptr<SiStripApvGain> first_payload = this->fetchPayload(std::get<1>(firstiov));
0242 
0243       SiStripApvGainContainer* l_objContainer = new SiStripApvGainContainer(last_payload, lastiov, tagname1);
0244       SiStripApvGainContainer* f_objContainer = new SiStripApvGainContainer(first_payload, firstiov, tagname2);
0245 
0246       l_objContainer->subtract(f_objContainer);
0247 
0248       //l_objContainer->printAll();
0249 
0250       TCanvas canvas("Partition summary", "partition summary", 1400, 1000);
0251       l_objContainer->fillByPartition(canvas, 100, -0.1, 0.1);
0252 
0253       std::string fileName(this->m_imageFileName);
0254       canvas.SaveAs(fileName.c_str());
0255 
0256       return true;
0257     }  // fill
0258   };
0259 
0260   using SiStripApvGainDiffByPartitionSingleTag = SiStripApvGainDiffByPartitionBase<1, MULTI_IOV>;
0261   using SiStripApvGainDiffByPartitionTwoTags = SiStripApvGainDiffByPartitionBase<2, SINGLE_IOV>;
0262 
0263   /************************************************
0264     1d histogram of SiStripApvGains of 1 IOV 
0265   *************************************************/
0266 
0267   // inherit from one of the predefined plot class: Histogram1D
0268   class SiStripApvGainsValue : public Histogram1D<SiStripApvGain, SINGLE_IOV> {
0269   public:
0270     SiStripApvGainsValue()
0271         : Histogram1D<SiStripApvGain, SINGLE_IOV>("SiStripApv Gains values", "SiStripApv Gains values", 200, 0.0, 2.0) {
0272     }
0273 
0274     bool fill() override {
0275       auto tag = PlotBase::getTag<0>();
0276       for (auto const& iov : tag.iovs) {
0277         std::shared_ptr<SiStripApvGain> payload = Base::fetchPayload(std::get<1>(iov));
0278         if (payload.get()) {
0279           std::vector<uint32_t> detid;
0280           payload->getDetIds(detid);
0281 
0282           for (const auto& d : detid) {
0283             SiStripApvGain::Range range = payload->getRange(d);
0284             for (int it = 0; it < range.second - range.first; it++) {
0285               // to be used to fill the histogram
0286               fillWithValue(payload->getApvGain(it, range));
0287 
0288             }  // loop over APVs
0289           }  // loop over detIds
0290         }  // payload
0291       }  // iovs
0292       return true;
0293     }  // fill
0294   };
0295 
0296   /************************************************
0297     1d histogram of means of SiStripApvGains
0298     for Tracker Barrel of 1 IOV 
0299   *************************************************/
0300 
0301   // inherit from one of the predefined plot class: Histogram1D
0302   class SiStripApvBarrelGainsByLayer : public Histogram1D<SiStripApvGain, SINGLE_IOV> {
0303   public:
0304     SiStripApvBarrelGainsByLayer()
0305         : Histogram1D<SiStripApvGain, SINGLE_IOV>("SiStripApv Gains averages by Barrel layer",
0306                                                   "Barrel layer (0-3: TIB), (4-9: TOB)",
0307                                                   10,
0308                                                   0,
0309                                                   10,
0310                                                   "average SiStripApv Gain") {}
0311 
0312     bool fill() override {
0313       auto tag = PlotBase::getTag<0>();
0314       for (auto const& iov : tag.iovs) {
0315         std::shared_ptr<SiStripApvGain> payload = Base::fetchPayload(std::get<1>(iov));
0316         if (payload.get()) {
0317           TrackerTopology tTopo = StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0318               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath());
0319 
0320           std::vector<uint32_t> detid;
0321           payload->getDetIds(detid);
0322 
0323           std::map<int, std::pair<float, float>> sumOfGainsByLayer;
0324 
0325           for (const auto& d : detid) {
0326             int subid = DetId(d).subdetId();
0327             int layer(-1);
0328             if (subid != StripSubdetector::TIB && subid != StripSubdetector::TOB)
0329               continue;
0330             if (subid == StripSubdetector::TIB) {
0331               layer = tTopo.tibLayer(d);
0332             } else if (subid == StripSubdetector::TOB) {
0333               // layers of TOB start at 5th bin
0334               layer = tTopo.tobLayer(d);
0335               layer += 4;
0336             }
0337 
0338             SiStripApvGain::Range range = payload->getRange(d);
0339             for (int it = 0; it < range.second - range.first; it++) {
0340               sumOfGainsByLayer[layer].first += payload->getApvGain(it, range);
0341               sumOfGainsByLayer[layer].second += 1.;
0342             }  // loop over APVs
0343           }  // loop over detIds
0344 
0345           // loop on the map to fill the plot
0346           for (auto& data : sumOfGainsByLayer) {
0347             fillWithBinAndValue(data.first - 1, (data.second.first / data.second.second));
0348           }
0349 
0350         }  // payload
0351       }  // iovs
0352       return true;
0353     }  // fill
0354   };
0355 
0356   /************************************************
0357     2d histogram of absolute (i.e. not average)
0358     SiStripApvGains for Tracker Barrel of 1 IOV
0359   *************************************************/
0360 
0361   class SiStripApvAbsoluteBarrelGainsByLayer : public Histogram2D<SiStripApvGain, SINGLE_IOV> {
0362   public:
0363     SiStripApvAbsoluteBarrelGainsByLayer()
0364         : Histogram2D<SiStripApvGain, SINGLE_IOV>("SiStripApv Gains by Barrel layer",
0365                                                   "Barrel layer (0-3: TIB), (4-9: TOB)",
0366                                                   10,
0367                                                   0,
0368                                                   10,
0369                                                   "SiStripApv Gain",
0370                                                   200,
0371                                                   0.0,
0372                                                   2.0) {}
0373     bool fill() override {
0374       auto tag = PlotBase::getTag<0>();
0375       for (auto const& iov : tag.iovs) {
0376         std::shared_ptr<SiStripApvGain> payload = Base::fetchPayload(std::get<1>(iov));
0377         if (payload.get()) {
0378           TrackerTopology tTopo = StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0379               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath());
0380 
0381           std::vector<uint32_t> detid;
0382           payload->getDetIds(detid);
0383           for (const auto& d : detid) {
0384             int subid = DetId(d).subdetId();
0385             if (subid != 3 && subid != 5)
0386               continue;
0387 
0388             SiStripApvGain::Range range = payload->getRange(d);
0389             for (int it = 0; it < range.second - range.first; it++) {
0390               float gain = payload->getApvGain(it, range);
0391               fillWithValue(static_cast<float>((subid == 5) ? tTopo.tobLayer(d) + 3 : tTopo.tibLayer(d) - 1),
0392                             (gain > 2.0) ? 2.0 : gain);
0393             }
0394           }  //loop over detIds
0395         }  // loop over payloads
0396       }  // loop over iovs
0397       return true;
0398     }  // fill
0399   };
0400 
0401   /************************************************
0402     1d histogram of means of SiStripApvGains
0403     for Tracker Endcaps (minus side) of 1 IOV 
0404   *************************************************/
0405 
0406   // inherit from one of the predefined plot class: Histogram1D
0407   class SiStripApvEndcapMinusGainsByDisk : public Histogram1D<SiStripApvGain, SINGLE_IOV> {
0408   public:
0409     SiStripApvEndcapMinusGainsByDisk()
0410         : Histogram1D<SiStripApvGain, SINGLE_IOV>("SiStripApv Gains averages by Endcap (minus) disk",
0411                                                   "Endcap (minus) disk (0-2: TID), (3-11: TEC)",
0412                                                   12,
0413                                                   0,
0414                                                   12,
0415                                                   "average SiStripApv Gain") {}
0416 
0417     bool fill() override {
0418       auto tag = PlotBase::getTag<0>();
0419       for (auto const& iov : tag.iovs) {
0420         std::shared_ptr<SiStripApvGain> payload = Base::fetchPayload(std::get<1>(iov));
0421         if (payload.get()) {
0422           TrackerTopology tTopo = StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0423               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath());
0424 
0425           std::vector<uint32_t> detid;
0426           payload->getDetIds(detid);
0427 
0428           std::map<int, std::pair<float, float>> sumOfGainsByDisk;
0429 
0430           for (const auto& d : detid) {
0431             int disk = -1;
0432             int side = -1;
0433             int subid = DetId(d).subdetId();
0434             if (subid != StripSubdetector::TID && subid != StripSubdetector::TEC)
0435               continue;
0436 
0437             if (subid == StripSubdetector::TID) {
0438               side = tTopo.tidSide(d);
0439               disk = tTopo.tidWheel(d);
0440             } else {
0441               side = tTopo.tecSide(d);
0442               disk = tTopo.tecWheel(d);
0443               // disks of TEC start at 4th bin
0444               disk += 3;
0445             }
0446 
0447             // only negative side
0448             if (side != 1)
0449               continue;
0450 
0451             SiStripApvGain::Range range = payload->getRange(d);
0452             for (int it = 0; it < range.second - range.first; it++) {
0453               sumOfGainsByDisk[disk].first += payload->getApvGain(it, range);
0454               sumOfGainsByDisk[disk].second += 1.;
0455             }  // loop over APVs
0456           }  // loop over detIds
0457 
0458           // loop on the map to fill the plot
0459           for (auto& data : sumOfGainsByDisk) {
0460             fillWithBinAndValue(data.first - 1, (data.second.first / data.second.second));
0461           }
0462 
0463         }  // payload
0464       }  // iovs
0465       return true;
0466     }  // fill
0467   };
0468 
0469   /************************************************
0470     1d histogram of means of SiStripApvGains
0471     for Tracker Endcaps (plus side) of 1 IOV 
0472   *************************************************/
0473 
0474   // inherit from one of the predefined plot class: Histogram1D
0475   class SiStripApvEndcapPlusGainsByDisk : public Histogram1D<SiStripApvGain, SINGLE_IOV> {
0476   public:
0477     SiStripApvEndcapPlusGainsByDisk()
0478         : Histogram1D<SiStripApvGain, SINGLE_IOV>("SiStripApv Gains averages by Endcap (plus) disk",
0479                                                   "Endcap (plus) disk (0-2: TID), (3-11: TEC)",
0480                                                   12,
0481                                                   0,
0482                                                   12,
0483                                                   "average SiStripApv Gain") {}
0484 
0485     bool fill() override {
0486       auto tag = PlotBase::getTag<0>();
0487       for (auto const& iov : tag.iovs) {
0488         std::shared_ptr<SiStripApvGain> payload = Base::fetchPayload(std::get<1>(iov));
0489         if (payload.get()) {
0490           TrackerTopology tTopo = StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0491               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath());
0492 
0493           std::vector<uint32_t> detid;
0494           payload->getDetIds(detid);
0495 
0496           std::map<int, std::pair<float, float>> sumOfGainsByDisk;
0497 
0498           for (const auto& d : detid) {
0499             int disk = -1;
0500             int side = -1;
0501             int subid = DetId(d).subdetId();
0502             if (subid != StripSubdetector::TID && subid != StripSubdetector::TEC)
0503               continue;
0504 
0505             if (subid == StripSubdetector::TID) {
0506               side = tTopo.tidSide(d);
0507               disk = tTopo.tidWheel(d);
0508               ;
0509             } else {
0510               side = tTopo.tecSide(d);
0511               disk = tTopo.tecWheel(d);
0512               // disks of TEC start at 4th bin
0513               disk += 3;
0514             }
0515 
0516             // only positive side
0517             if (side != 2)
0518               continue;
0519 
0520             SiStripApvGain::Range range = payload->getRange(d);
0521             for (int it = 0; it < range.second - range.first; it++) {
0522               sumOfGainsByDisk[disk].first += payload->getApvGain(it, range);
0523               sumOfGainsByDisk[disk].second += 1.;
0524             }  // loop over APVs
0525           }  // loop over detIds
0526 
0527           // loop on the map to fill the plot
0528           for (auto& data : sumOfGainsByDisk) {
0529             fillWithBinAndValue(data.first - 1, (data.second.first / data.second.second));
0530           }
0531 
0532         }  // payload
0533       }  // iovs
0534       return true;
0535     }  // fill
0536   };
0537 
0538   /************************************************
0539     2D histogram of absolute (i.e. not average)
0540     SiStripApv Gains on the Endcap- for 1 IOV
0541    ************************************************/
0542   class SiStripApvAbsoluteEndcapMinusGainsByDisk : public Histogram2D<SiStripApvGain, SINGLE_IOV> {
0543   public:
0544     SiStripApvAbsoluteEndcapMinusGainsByDisk()
0545         : Histogram2D<SiStripApvGain, SINGLE_IOV>("SiStripApv Gains averages by Endcap (minus) disk",
0546                                                   "Endcap (minus) disk (0-2: TID), (3-11: TEC)",
0547                                                   12,
0548                                                   0,
0549                                                   12,
0550                                                   "SiStripApv Gain",
0551                                                   200,
0552                                                   0.0,
0553                                                   2.0) {}
0554 
0555     bool fill() override {
0556       auto tag = PlotBase::getTag<0>();
0557       for (auto const& iov : tag.iovs) {
0558         std::shared_ptr<SiStripApvGain> payload = Base::fetchPayload(std::get<1>(iov));
0559         if (payload.get()) {
0560           TrackerTopology tTopo = StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0561               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath());
0562 
0563           std::vector<uint32_t> detid;
0564           payload->getDetIds(detid);
0565 
0566           for (const auto& d : detid) {
0567             int subid = DetId(d).subdetId(), side = -1, disk = -1;
0568 
0569             switch (subid) {
0570               case 4:
0571                 side = tTopo.tidSide(d);
0572                 disk = tTopo.tidWheel(d);
0573                 break;
0574               case 6:
0575                 side = tTopo.tecSide(d);
0576                 disk = tTopo.tecWheel(d) + 4;
0577                 break;
0578               default:
0579                 continue;
0580             }
0581 
0582             if (side != 1)
0583               continue;
0584             SiStripApvGain::Range range = payload->getRange(d);
0585             for (int it = 0; it < range.second - range.first; it++) {
0586               float gain = payload->getApvGain(it, range);
0587               fillWithValue((float)disk - 1, (gain > 2.0) ? 2.0 : gain);
0588             }  // apvs
0589           }  // detids
0590         }
0591       }  // iovs
0592       return true;
0593     }  // fill
0594   };
0595 
0596   /************************************************
0597     2D histogram of absolute (i.e. not average)
0598     SiStripApv Gains on the Endcap+ for 1 IOV
0599    ************************************************/
0600   class SiStripApvAbsoluteEndcapPlusGainsByDisk : public Histogram2D<SiStripApvGain, SINGLE_IOV> {
0601   public:
0602     SiStripApvAbsoluteEndcapPlusGainsByDisk()
0603         : Histogram2D<SiStripApvGain, SINGLE_IOV>("SiStripApv Gains averages by Endcap (plus) disk",
0604                                                   "Endcap (plus) disk (0-2: TID), (3-11: TEC)",
0605                                                   12,
0606                                                   0,
0607                                                   12,
0608                                                   "SiStripApv Gain",
0609                                                   200,
0610                                                   0.0,
0611                                                   2.0) {}
0612     bool fill() override {
0613       auto tag = PlotBase::getTag<0>();
0614       for (auto const& iov : tag.iovs) {
0615         std::shared_ptr<SiStripApvGain> payload = Base::fetchPayload(std::get<1>(iov));
0616         if (payload.get()) {
0617           TrackerTopology tTopo = StandaloneTrackerTopology::fromTrackerParametersXMLFile(
0618               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath());
0619 
0620           std::vector<uint32_t> detid;
0621           payload->getDetIds(detid);
0622 
0623           for (const auto& d : detid) {
0624             int subid = DetId(d).subdetId(), side = -1, disk = -1;
0625 
0626             switch (subid) {
0627               case 4:
0628                 side = tTopo.tidSide(d);
0629                 disk = tTopo.tidWheel(d);
0630                 break;
0631               case 6:
0632                 side = tTopo.tecSide(d);
0633                 disk = tTopo.tecWheel(d) + 4;
0634                 break;
0635               default:
0636                 continue;
0637             }
0638 
0639             if (side != 2)
0640               continue;
0641             SiStripApvGain::Range range = payload->getRange(d);
0642             for (int it = 0; it < range.second - range.first; it++) {
0643               float gain = payload->getApvGain(it, range);
0644               fillWithValue((float)disk - 1, (gain > 2.0) ? 2.0 : gain);
0645             }  //apvs
0646           }  //detids
0647         }
0648       }  // iovs
0649       return true;
0650     }  // fill
0651   };
0652 
0653   /************************************************
0654     TrackerMap of SiStripApvGains (average gain per detid)
0655   *************************************************/
0656   class SiStripApvGainsAverageTrackerMap : public PlotImage<SiStripApvGain, SINGLE_IOV> {
0657   public:
0658     SiStripApvGainsAverageTrackerMap() : PlotImage<SiStripApvGain, SINGLE_IOV>("Tracker Map of average SiStripGains") {}
0659 
0660     bool fill() override {
0661       auto tag = PlotBase::getTag<0>();
0662       auto iov = tag.iovs.front();
0663       std::shared_ptr<SiStripApvGain> payload = fetchPayload(std::get<1>(iov));
0664 
0665       std::string titleMap = "SiStrip APV Gain average per module (payload : " + std::get<1>(iov) + ")";
0666 
0667       std::unique_ptr<TrackerMap> tmap = std::make_unique<TrackerMap>("SiStripApvGains");
0668       tmap->setTitle(titleMap);
0669       tmap->setPalette(1);
0670 
0671       std::vector<uint32_t> detid;
0672       payload->getDetIds(detid);
0673 
0674       std::map<uint32_t, float> store;
0675 
0676       for (const auto& d : detid) {
0677         SiStripApvGain::Range range = payload->getRange(d);
0678         float sumOfGains = 0;
0679         float nAPVsPerModule = 0.;
0680         for (int it = 0; it < range.second - range.first; it++) {
0681           nAPVsPerModule += 1;
0682           sumOfGains += payload->getApvGain(it, range);
0683         }  // loop over APVs
0684         // fill the tracker map taking the average gain on a single DetId
0685         if (nAPVsPerModule != 0.) {
0686           store[d] = (sumOfGains / nAPVsPerModule);
0687           tmap->fill(d, (sumOfGains / nAPVsPerModule));
0688         }
0689       }  // loop over detIds
0690 
0691       //=========================
0692       // saturate at 2 std deviations
0693       auto range = SiStripPI::getTheRange(store, 2);
0694 
0695       std::string fileName(m_imageFileName);
0696       tmap->save(true, range.first, range.second, fileName);
0697 
0698       return true;
0699     }
0700   };
0701 
0702   /************************************************
0703     TrackerMap of SiStripApvGains (module with default)
0704   *************************************************/
0705   class SiStripApvGainsDefaultTrackerMap : public PlotImage<SiStripApvGain, SINGLE_IOV> {
0706   public:
0707     SiStripApvGainsDefaultTrackerMap()
0708         : PlotImage<SiStripApvGain, SINGLE_IOV>("Tracker Map of SiStripGains to default") {}
0709 
0710     bool fill() override {
0711       auto tag = PlotBase::getTag<0>();
0712       auto iov = tag.iovs.front();
0713       std::shared_ptr<SiStripApvGain> payload = fetchPayload(std::get<1>(iov));
0714 
0715       std::unique_ptr<TrackerMap> tmap = std::make_unique<TrackerMap>("SiStripApvGains");
0716 
0717       tmap->setPalette(1);
0718 
0719       std::vector<uint32_t> detid;
0720       payload->getDetIds(detid);
0721 
0722       /*
0723     the default G1 value comes from the ratio of DefaultTickHeight/GainNormalizationFactor
0724     as defined in the default of the O2O producer: OnlineDB/SiStripESSources/src/SiStripCondObjBuilderFromDb.cc
0725       */
0726 
0727       constexpr float G1default = 690. / 640.;
0728       constexpr float G2default = 1.f;
0729 
0730       int totalG1DefaultAPVs = 0;
0731       int totalG2DefaultAPVs = 0;
0732 
0733       for (const auto& d : detid) {
0734         SiStripApvGain::Range range = payload->getRange(d);
0735         float sumOfGains = 0;
0736         float nAPVsPerModule = 0.;
0737         int countDefaults = 0;
0738         for (int it = 0; it < range.second - range.first; it++) {
0739           nAPVsPerModule += 1;
0740           sumOfGains += payload->getApvGain(it, range);
0741           if ((payload->getApvGain(it, range)) == G1default || (payload->getApvGain(it, range)) == G2default)
0742             countDefaults++;
0743         }  // loop over APVs
0744         // fill the tracker map taking the average gain on a single DetId
0745         if (countDefaults > 0.) {
0746           tmap->fill(d, countDefaults);
0747           // check here if it exists at least a full module with all APVs defaulted
0748           if (countDefaults == nAPVsPerModule) {
0749             if (std::fmod((sumOfGains / countDefaults), G1default) == 0.) {
0750               totalG1DefaultAPVs += countDefaults;
0751             } else if (std::fmod((sumOfGains / countDefaults), G2default) == 0.) {
0752               totalG2DefaultAPVs += countDefaults;
0753             }
0754           }
0755         }
0756       }  // loop over detIds
0757 
0758       //=========================
0759       std::string gainType{""};
0760       if (totalG2DefaultAPVs == 0) {
0761         gainType = "G1 value (=690./640.)";
0762       } else if (totalG1DefaultAPVs == 0) {
0763         gainType = "G2 value (=1)";
0764       }
0765 
0766       std::string titleMap = "# of APVs/module w/ default " + gainType + " (payload : " + std::get<1>(iov) + ")";
0767       tmap->setTitle(titleMap);
0768 
0769       std::string fileName(m_imageFileName);
0770       tmap->save(true, 0, 0, fileName);
0771 
0772       return true;
0773     }
0774   };
0775 
0776   /************************************************
0777     TrackerMap of SiStripApvGains (ratio with previous gain per detid)
0778   *************************************************/
0779 
0780   template <int ntags, IOVMultiplicity nIOVs>
0781   class SiStripApvGainsRatioTrackerMapBase : public PlotImage<SiStripApvGain, nIOVs, ntags> {
0782   public:
0783     SiStripApvGainsRatioTrackerMapBase()
0784         : PlotImage<SiStripApvGain, nIOVs, ntags>("Tracker Map of ratio of SiStripGains with previous IOV") {
0785       PlotBase::addInputParam("nsigma");
0786     }
0787 
0788     bool fill() override {
0789       // determine n. sigmas
0790       unsigned int nsigma(1);
0791 
0792       auto paramValues = PlotBase::inputParamValues();
0793       auto ip = paramValues.find("nsigma");
0794       if (ip != paramValues.end()) {
0795         nsigma = std::stoul(ip->second);
0796       }
0797 
0798       // trick to deal with the multi-ioved tag and two tag case at the same time
0799       auto theIOVs = PlotBase::getTag<0>().iovs;
0800       auto tagname1 = PlotBase::getTag<0>().name;
0801       std::string tagname2 = "";
0802       auto firstiov = theIOVs.front();
0803       SiStripPI::MetaData lastiov;
0804 
0805       // we don't support (yet) comparison with more than 2 tags
0806       assert(this->m_plotAnnotations.ntags < 3);
0807 
0808       if (this->m_plotAnnotations.ntags == 2) {
0809         auto tag2iovs = PlotBase::getTag<1>().iovs;
0810         tagname2 = PlotBase::getTag<1>().name;
0811         lastiov = tag2iovs.front();
0812       } else {
0813         lastiov = theIOVs.back();
0814       }
0815 
0816       std::shared_ptr<SiStripApvGain> last_payload = this->fetchPayload(std::get<1>(lastiov));
0817       std::shared_ptr<SiStripApvGain> first_payload = this->fetchPayload(std::get<1>(firstiov));
0818 
0819       std::string titleMap = "SiStrip APV Gain ratio per module average (IOV: ";
0820 
0821       titleMap += std::to_string(std::get<0>(firstiov));
0822       titleMap += "/ IOV:";
0823       titleMap += std::to_string(std::get<0>(lastiov));
0824       titleMap += ")";
0825 
0826       titleMap += +" " + std::to_string(nsigma) + " std. dev. saturation";
0827 
0828       std::unique_ptr<TrackerMap> tmap = std::make_unique<TrackerMap>("SiStripApvGains");
0829       tmap->setTitle(titleMap);
0830       tmap->setPalette(1);
0831 
0832       std::map<uint32_t, float> lastmap, firstmap;
0833 
0834       std::vector<uint32_t> detid;
0835       last_payload->getDetIds(detid);
0836 
0837       // cache the last IOV
0838       for (const auto& d : detid) {
0839         SiStripApvGain::Range range = last_payload->getRange(d);
0840         float Gain = 0;
0841         float nAPV = 0;
0842         for (int it = 0; it < range.second - range.first; it++) {
0843           nAPV += 1;
0844           Gain += last_payload->getApvGain(it, range);
0845         }  // loop over APVs
0846         lastmap[d] = (Gain / nAPV);
0847       }  // loop over detIds
0848 
0849       detid.clear();
0850 
0851       first_payload->getDetIds(detid);
0852 
0853       // cache the first IOV
0854       for (const auto& d : detid) {
0855         SiStripApvGain::Range range = first_payload->getRange(d);
0856         float Gain = 0;
0857         float nAPV = 0;
0858         for (int it = 0; it < range.second - range.first; it++) {
0859           nAPV += 1;
0860           Gain += first_payload->getApvGain(it, range);
0861         }  // loop over APVs
0862         firstmap[d] = (Gain / nAPV);
0863       }  // loop over detIds
0864 
0865       std::map<uint32_t, float> cachedRatio;
0866       for (const auto& d : detid) {
0867         float ratio = firstmap[d] / lastmap[d];
0868         tmap->fill(d, ratio);
0869         cachedRatio[d] = ratio;
0870       }
0871 
0872       //=========================
0873       auto range = SiStripPI::getTheRange(cachedRatio, nsigma);
0874 
0875       std::string fileName(this->m_imageFileName);
0876       tmap->save(true, range.first, range.second, fileName);
0877 
0878       return true;
0879     }
0880   };
0881 
0882   using SiStripApvGainsAvgDeviationRatioWithPreviousIOVTrackerMap = SiStripApvGainsRatioTrackerMapBase<1, MULTI_IOV>;
0883   using SiStripApvGainsAvgDeviationRatioTrackerMapTwoTags = SiStripApvGainsRatioTrackerMapBase<2, SINGLE_IOV>;
0884 
0885   /************************************************
0886    TrackerMap of SiStripApvGains (ratio for largest deviation with previous gain per detid)
0887   *************************************************/
0888 
0889   template <int ntags, IOVMultiplicity nIOVs>
0890   class SiStripApvGainsRatioMaxDeviationTrackerMapBase : public PlotImage<SiStripApvGain, nIOVs, ntags> {
0891   public:
0892     SiStripApvGainsRatioMaxDeviationTrackerMapBase()
0893         : PlotImage<SiStripApvGain, nIOVs, ntags>(
0894               "Tracker Map of ratio (for largest deviation) of SiStripGains with previous IOV") {
0895       PlotBase::addInputParam("nsigma");
0896     }
0897 
0898     bool fill() override {
0899       unsigned int nsigma(1);
0900       auto paramValues = PlotBase::inputParamValues();
0901       auto ip = paramValues.find("nsigma");
0902       if (ip != paramValues.end()) {
0903         nsigma = std::stoul(ip->second);
0904         edm::LogPrint("SiStripApvGain_PayloadInspector")
0905             << "using custom z-axis saturation: " << nsigma << " sigmas" << std::endl;
0906       } else {
0907         edm::LogPrint("SiStripApvGain_PayloadInspector")
0908             << "using default saturation: " << nsigma << " sigmas" << std::endl;
0909       }
0910 
0911       // trick to deal with the multi-ioved tag and two tag case at the same time
0912       auto theIOVs = PlotBase::getTag<0>().iovs;
0913       auto tagname1 = PlotBase::getTag<0>().name;
0914       std::string tagname2 = "";
0915       auto firstiov = theIOVs.front();
0916       SiStripPI::MetaData lastiov;
0917 
0918       // we don't support (yet) comparison with more than 2 tags
0919       assert(this->m_plotAnnotations.ntags < 3);
0920 
0921       if (this->m_plotAnnotations.ntags == 2) {
0922         auto tag2iovs = PlotBase::getTag<1>().iovs;
0923         tagname2 = PlotBase::getTag<1>().name;
0924         lastiov = tag2iovs.front();
0925       } else {
0926         lastiov = theIOVs.back();
0927       }
0928 
0929       std::shared_ptr<SiStripApvGain> last_payload = this->fetchPayload(std::get<1>(lastiov));
0930       std::shared_ptr<SiStripApvGain> first_payload = this->fetchPayload(std::get<1>(firstiov));
0931 
0932       std::string titleMap = "SiStrip APV Gain ratio for largest deviation per module (IOV: ";
0933 
0934       titleMap += std::to_string(std::get<0>(firstiov));
0935       titleMap += "/ IOV:";
0936       titleMap += std::to_string(std::get<0>(lastiov));
0937       titleMap += ") ";
0938 
0939       titleMap += +" - " + std::to_string(nsigma) + " std. dev. saturation";
0940 
0941       std::unique_ptr<TrackerMap> tmap = std::make_unique<TrackerMap>("SiStripApvGains");
0942       tmap->setTitle(titleMap);
0943       tmap->setPalette(1);
0944 
0945       std::map<std::pair<uint32_t, int>, float> lastmap, firstmap;
0946 
0947       std::vector<uint32_t> detid;
0948       last_payload->getDetIds(detid);
0949 
0950       // cache the last IOV
0951       for (const auto& d : detid) {
0952         SiStripApvGain::Range range = last_payload->getRange(d);
0953         float nAPV = 0;
0954         for (int it = 0; it < range.second - range.first; it++) {
0955           nAPV += 1;
0956           float Gain = last_payload->getApvGain(it, range);
0957           std::pair<uint32_t, int> index = std::make_pair(d, nAPV);
0958           lastmap[index] = Gain;
0959         }  // loop over APVs
0960       }  // loop over detIds
0961 
0962       detid.clear();
0963 
0964       first_payload->getDetIds(detid);
0965 
0966       // cache the first IOV
0967       for (const auto& d : detid) {
0968         SiStripApvGain::Range range = first_payload->getRange(d);
0969         float nAPV = 0;
0970         for (int it = 0; it < range.second - range.first; it++) {
0971           nAPV += 1;
0972           float Gain = first_payload->getApvGain(it, range);
0973           std::pair<uint32_t, int> index = std::make_pair(d, nAPV);
0974           firstmap[index] = Gain;
0975         }  // loop over APVs
0976       }  // loop over detIds
0977 
0978       // find the largest deviation
0979       std::map<uint32_t, float> cachedRatio;
0980 
0981       for (const auto& item : firstmap) {
0982         // packed index (detid,APV)
0983         auto index = item.first;
0984         auto mod = item.first.first;
0985 
0986         float ratio = firstmap[index] / lastmap[index];
0987         // if we have already cached something
0988         if (cachedRatio[mod]) {
0989           // if the discrepancy with 1 of ratio is larger than the cached value
0990           if (std::abs(ratio - 1.) > std::abs(cachedRatio[mod] - 1.)) {
0991             cachedRatio[mod] = ratio;
0992           }
0993         } else {
0994           cachedRatio[mod] = ratio;
0995         }
0996       }
0997 
0998       for (const auto& element : cachedRatio) {
0999         tmap->fill(element.first, element.second);
1000       }
1001 
1002       // get the range of the TrackerMap (saturate at +/-n std deviations)
1003       auto range = SiStripPI::getTheRange(cachedRatio, nsigma);
1004 
1005       //=========================
1006 
1007       std::string fileName(this->m_imageFileName);
1008       tmap->save(true, range.first, range.second, fileName);
1009 
1010       return true;
1011     }
1012   };
1013 
1014   using SiStripApvGainsMaxDeviationRatioWithPreviousIOVTrackerMap =
1015       SiStripApvGainsRatioMaxDeviationTrackerMapBase<1, MULTI_IOV>;
1016 
1017   using SiStripApvGainsMaxDeviationRatioTrackerMapTwoTags =
1018       SiStripApvGainsRatioMaxDeviationTrackerMapBase<2, SINGLE_IOV>;
1019 
1020   /************************************************
1021     TrackerMap of SiStripApvGains (maximum gain per detid)
1022   *************************************************/
1023   class SiStripApvGainsMaximumTrackerMap : public PlotImage<SiStripApvGain, SINGLE_IOV> {
1024   public:
1025     SiStripApvGainsMaximumTrackerMap()
1026         : PlotImage<SiStripApvGain, SINGLE_IOV>("Tracker Map of SiStripAPVGains (maximum per DetId)") {}
1027 
1028     bool fill() override {
1029       auto tag = PlotBase::getTag<0>();
1030       auto iov = tag.iovs.front();
1031       std::shared_ptr<SiStripApvGain> payload = fetchPayload(std::get<1>(iov));
1032 
1033       std::string titleMap = "SiStrip APV Gain maximum per module (payload : " + std::get<1>(iov) + ")";
1034 
1035       std::unique_ptr<TrackerMap> tmap = std::make_unique<TrackerMap>("SiStripApvGains");
1036       tmap->setTitle(titleMap);
1037       tmap->setPalette(1);
1038 
1039       std::vector<uint32_t> detid;
1040       payload->getDetIds(detid);
1041 
1042       for (const auto& d : detid) {
1043         SiStripApvGain::Range range = payload->getRange(d);
1044         float theMaxGain = 0;
1045         for (int it = 0; it < range.second - range.first; it++) {
1046           float currentGain = payload->getApvGain(it, range);
1047           if (currentGain > theMaxGain) {
1048             theMaxGain = currentGain;
1049           }
1050         }  // loop over APVs
1051         // fill the tracker map taking the maximum gain on a single DetId
1052         tmap->fill(d, theMaxGain);
1053       }  // loop over detIds
1054 
1055       //=========================
1056 
1057       std::pair<float, float> extrema = tmap->getAutomaticRange();
1058 
1059       std::string fileName(m_imageFileName);
1060 
1061       // protect against uniform values (gains are defined positive)
1062       if (extrema.first != extrema.second) {
1063         tmap->save(true, 0, 0, fileName);
1064       } else {
1065         tmap->save(true, extrema.first * 0.95, extrema.first * 1.05, fileName);
1066       }
1067 
1068       return true;
1069     }
1070   };
1071 
1072   /************************************************
1073     TrackerMap of SiStripApvGains (minimum gain per detid)
1074   *************************************************/
1075   class SiStripApvGainsMinimumTrackerMap : public PlotImage<SiStripApvGain, SINGLE_IOV> {
1076   public:
1077     SiStripApvGainsMinimumTrackerMap()
1078         : PlotImage<SiStripApvGain, SINGLE_IOV>("Tracker Map of SiStripAPVGains (minimum per DetId)") {}
1079 
1080     bool fill() override {
1081       auto tag = PlotBase::getTag<0>();
1082       auto iov = tag.iovs.front();
1083 
1084       std::shared_ptr<SiStripApvGain> payload = fetchPayload(std::get<1>(iov));
1085 
1086       std::string titleMap = "SiStrip APV Gain minumum per module (payload : " + std::get<1>(iov) + ")";
1087 
1088       std::unique_ptr<TrackerMap> tmap = std::make_unique<TrackerMap>("SiStripApvGains");
1089       tmap->setTitle(titleMap);
1090       tmap->setPalette(1);
1091 
1092       std::vector<uint32_t> detid;
1093       payload->getDetIds(detid);
1094 
1095       for (const auto& d : detid) {
1096         SiStripApvGain::Range range = payload->getRange(d);
1097         float theMinGain = 999.;
1098         for (int it = 0; it < range.second - range.first; it++) {
1099           float currentGain = payload->getApvGain(it, range);
1100           if (currentGain < theMinGain) {
1101             theMinGain = currentGain;
1102           }
1103         }  // loop over APVs
1104         // fill the tracker map taking the minimum gain on a single DetId
1105         tmap->fill(d, theMinGain);
1106       }  // loop over detIds
1107 
1108       //=========================
1109 
1110       std::pair<float, float> extrema = tmap->getAutomaticRange();
1111 
1112       std::string fileName(m_imageFileName);
1113 
1114       // protect against uniform values (gains are defined positive)
1115       if (extrema.first != extrema.second) {
1116         tmap->save(true, 0, 0, fileName);
1117       } else {
1118         tmap->save(true, extrema.first * 0.95, extrema.first * 1.05, fileName);
1119       }
1120 
1121       return true;
1122     }
1123   };
1124 
1125   /************************************************
1126     time history histogram of SiStripApvGains 
1127   *************************************************/
1128 
1129   class SiStripApvGainByRunMeans : public HistoryPlot<SiStripApvGain, float> {
1130   public:
1131     SiStripApvGainByRunMeans()
1132         : HistoryPlot<SiStripApvGain, float>("SiStripApv Gains average", "average Strip APV gain value") {}
1133     ~SiStripApvGainByRunMeans() override = default;
1134 
1135     float getFromPayload(SiStripApvGain& payload) override {
1136       std::vector<uint32_t> detid;
1137       payload.getDetIds(detid);
1138 
1139       float nAPVs = 0;
1140       float sumOfGains = 0;
1141 
1142       for (const auto& d : detid) {
1143         SiStripApvGain::Range range = payload.getRange(d);
1144         for (int it = 0; it < range.second - range.first; it++) {
1145           nAPVs += 1;
1146           sumOfGains += payload.getApvGain(it, range);
1147         }  // loop over APVs
1148       }  // loop over detIds
1149 
1150       return sumOfGains / nAPVs;
1151     }  // payload
1152   };
1153 
1154   /************************************************
1155     time history of SiStripApvGains properties
1156   *************************************************/
1157 
1158   template <SiStripPI::estimator est>
1159   class SiStripApvGainProperties : public HistoryPlot<SiStripApvGain, float> {
1160   public:
1161     SiStripApvGainProperties()
1162         : HistoryPlot<SiStripApvGain, float>("SiStripApv Gains " + estimatorType(est),
1163                                              estimatorType(est) + " Strip APV gain value") {}
1164     ~SiStripApvGainProperties() override = default;
1165 
1166     float getFromPayload(SiStripApvGain& payload) override {
1167       std::vector<uint32_t> detid;
1168       payload.getDetIds(detid);
1169 
1170       float nAPVs = 0;
1171       float sumOfGains = 0;
1172       float meanOfGains = 0;
1173       float rmsOfGains = 0;
1174       float min(0.), max(0.);
1175 
1176       for (const auto& d : detid) {
1177         SiStripApvGain::Range range = payload.getRange(d);
1178         for (int it = 0; it < range.second - range.first; it++) {
1179           nAPVs += 1;
1180           float gain = payload.getApvGain(it, range);
1181           if (gain < min)
1182             min = gain;
1183           if (gain > max)
1184             max = gain;
1185           sumOfGains += gain;
1186           rmsOfGains += (gain * gain);
1187         }  // loop over APVs
1188       }  // loop over detIds
1189 
1190       meanOfGains = sumOfGains / nAPVs;
1191 
1192       switch (est) {
1193         case SiStripPI::min:
1194           return min;
1195           break;
1196         case SiStripPI::max:
1197           return max;
1198           break;
1199         case SiStripPI::mean:
1200           return meanOfGains;
1201           break;
1202         case SiStripPI::rms:
1203           if ((rmsOfGains / nAPVs - meanOfGains * meanOfGains) > 0.) {
1204             return sqrt(rmsOfGains / nAPVs - meanOfGains * meanOfGains);
1205           } else {
1206             return 0.;
1207           }
1208           break;
1209         default:
1210           edm::LogWarning("LogicError") << "Unknown estimator: " << est;
1211           break;
1212       }
1213       return 0.;
1214     }  // payload
1215   };
1216 
1217   typedef SiStripApvGainProperties<SiStripPI::min> SiStripApvGainMin_History;
1218   typedef SiStripApvGainProperties<SiStripPI::max> SiStripApvGainMax_History;
1219   typedef SiStripApvGainProperties<SiStripPI::mean> SiStripApvGainMean_History;
1220   typedef SiStripApvGainProperties<SiStripPI::rms> SiStripApvGainRMS_History;
1221 
1222   /************************************************
1223     time history histogram of TIB SiStripApvGains 
1224   *************************************************/
1225 
1226   class SiStripApvTIBGainByRunMeans : public HistoryPlot<SiStripApvGain, float> {
1227   public:
1228     SiStripApvTIBGainByRunMeans()
1229         : HistoryPlot<SiStripApvGain, float>("SiStripApv Gains average",
1230                                              "average Tracker Inner Barrel APV gain value") {}
1231     ~SiStripApvTIBGainByRunMeans() override = default;
1232 
1233     float getFromPayload(SiStripApvGain& payload) override {
1234       std::vector<uint32_t> detid;
1235       payload.getDetIds(detid);
1236 
1237       float nAPVs = 0;
1238       float sumOfGains = 0;
1239 
1240       for (const auto& d : detid) {
1241         int subid = DetId(d).subdetId();
1242         if (subid != StripSubdetector::TIB)
1243           continue;
1244 
1245         SiStripApvGain::Range range = payload.getRange(d);
1246         for (int it = 0; it < range.second - range.first; it++) {
1247           nAPVs += 1;
1248           sumOfGains += payload.getApvGain(it, range);
1249         }  // loop over APVs
1250       }  // loop over detIds
1251 
1252       return sumOfGains / nAPVs;
1253 
1254     }  // payload
1255   };
1256 
1257   /************************************************
1258     time history histogram of TOB SiStripApvGains 
1259   *************************************************/
1260 
1261   class SiStripApvTOBGainByRunMeans : public HistoryPlot<SiStripApvGain, float> {
1262   public:
1263     SiStripApvTOBGainByRunMeans()
1264         : HistoryPlot<SiStripApvGain, float>("SiStripApv Gains average", "average Tracker Outer Barrel gain value") {}
1265     ~SiStripApvTOBGainByRunMeans() override = default;
1266 
1267     float getFromPayload(SiStripApvGain& payload) override {
1268       std::vector<uint32_t> detid;
1269       payload.getDetIds(detid);
1270 
1271       float nAPVs = 0;
1272       float sumOfGains = 0;
1273 
1274       for (const auto& d : detid) {
1275         int subid = DetId(d).subdetId();
1276         if (subid != StripSubdetector::TOB)
1277           continue;
1278 
1279         SiStripApvGain::Range range = payload.getRange(d);
1280         for (int it = 0; it < range.second - range.first; it++) {
1281           nAPVs += 1;
1282           sumOfGains += payload.getApvGain(it, range);
1283         }  // loop over APVs
1284       }  // loop over detIds
1285 
1286       return sumOfGains / nAPVs;
1287 
1288     }  // payload
1289   };
1290 
1291   /************************************************
1292     time history histogram of TID SiStripApvGains 
1293   *************************************************/
1294 
1295   class SiStripApvTIDGainByRunMeans : public HistoryPlot<SiStripApvGain, float> {
1296   public:
1297     SiStripApvTIDGainByRunMeans()
1298         : HistoryPlot<SiStripApvGain, float>("SiStripApv Gains average", "average Tracker Inner Disks APV gain value") {
1299     }
1300     ~SiStripApvTIDGainByRunMeans() override = default;
1301 
1302     float getFromPayload(SiStripApvGain& payload) override {
1303       std::vector<uint32_t> detid;
1304       payload.getDetIds(detid);
1305 
1306       float nAPVs = 0;
1307       float sumOfGains = 0;
1308       for (const auto& d : detid) {
1309         int subid = DetId(d).subdetId();
1310         if (subid != StripSubdetector::TID)
1311           continue;
1312 
1313         SiStripApvGain::Range range = payload.getRange(d);
1314         for (int it = 0; it < range.second - range.first; it++) {
1315           nAPVs += 1;
1316           sumOfGains += payload.getApvGain(it, range);
1317         }  // loop over APVs
1318       }  // loop over detIds
1319 
1320       return sumOfGains / nAPVs;
1321 
1322     }  // payload
1323   };
1324 
1325   /************************************************
1326     time history histogram of TEC SiStripApvGains 
1327   *************************************************/
1328 
1329   class SiStripApvTECGainByRunMeans : public HistoryPlot<SiStripApvGain, float> {
1330   public:
1331     SiStripApvTECGainByRunMeans()
1332         : HistoryPlot<SiStripApvGain, float>("SiStripApv Gains average in TEC",
1333                                              "average Tracker Endcaps APV gain value") {}
1334     ~SiStripApvTECGainByRunMeans() override = default;
1335 
1336     float getFromPayload(SiStripApvGain& payload) override {
1337       std::vector<uint32_t> detid;
1338       payload.getDetIds(detid);
1339 
1340       float nAPVs = 0;
1341       float sumOfGains = 0;
1342 
1343       for (const auto& d : detid) {
1344         int subid = DetId(d).subdetId();
1345         if (subid != StripSubdetector::TEC)
1346           continue;
1347 
1348         SiStripApvGain::Range range = payload.getRange(d);
1349         for (int it = 0; it < range.second - range.first; it++) {
1350           nAPVs += 1;
1351           sumOfGains += payload.getApvGain(it, range);
1352         }  // loop over APVs
1353       }  // loop over detIds
1354 
1355       return sumOfGains / nAPVs;
1356 
1357     }  // payload
1358   };
1359 
1360   /************************************************
1361     test class
1362   *************************************************/
1363 
1364   class SiStripApvGainsTest : public Histogram1D<SiStripApvGain, SINGLE_IOV> {
1365   public:
1366     SiStripApvGainsTest()
1367         : Histogram1D<SiStripApvGain, SINGLE_IOV>("SiStripApv Gains test", "SiStripApv Gains test", 10, 0.0, 10.0),
1368           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
1369               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath())} {}
1370 
1371     bool fill() override {
1372       auto tag = PlotBase::getTag<0>();
1373       for (auto const& iov : tag.iovs) {
1374         std::shared_ptr<SiStripApvGain> payload = Base::fetchPayload(std::get<1>(iov));
1375         if (payload.get()) {
1376           std::vector<uint32_t> detid;
1377           payload->getDetIds(detid);
1378 
1379           SiStripDetSummary summaryGain{&m_trackerTopo};
1380 
1381           for (const auto& d : detid) {
1382             SiStripApvGain::Range range = payload->getRange(d);
1383             for (int it = 0; it < range.second - range.first; ++it) {
1384               summaryGain.add(d, payload->getApvGain(it, range));
1385               fillWithValue(payload->getApvGain(it, range));
1386             }
1387           }
1388           std::map<unsigned int, SiStripDetSummary::Values> map = summaryGain.getCounts();
1389 
1390           //SiStripPI::printSummary(map);
1391 
1392           std::stringstream ss;
1393           ss << "Summary of gain values:" << std::endl;
1394           summaryGain.print(ss, true);
1395           std::cout << ss.str() << std::endl;
1396 
1397         }  // payload
1398       }  // iovs
1399       return true;
1400     }  // fill
1401   private:
1402     TrackerTopology m_trackerTopo;
1403   };
1404 
1405   /************************************************
1406     Compare Gains from 2 IOVs, 2 pads canvas, firsr for ratio, second for scatter plot
1407   *************************************************/
1408 
1409   template <int ntags, IOVMultiplicity nIOVs>
1410   class SiStripApvGainsComparatorBase : public PlotImage<SiStripApvGain, nIOVs, ntags> {
1411   public:
1412     SiStripApvGainsComparatorBase() : PlotImage<SiStripApvGain, nIOVs, ntags>("SiStripGains Comparison") {}
1413 
1414     bool fill() override {
1415       // trick to deal with the multi-ioved tag and two tag case at the same time
1416       auto theIOVs = PlotBase::getTag<0>().iovs;
1417       auto tagname1 = PlotBase::getTag<0>().name;
1418       std::string tagname2 = "";
1419       auto firstiov = theIOVs.front();
1420       SiStripPI::MetaData lastiov;
1421 
1422       // we don't support (yet) comparison with more than 2 tags
1423       assert(this->m_plotAnnotations.ntags < 3);
1424 
1425       if (this->m_plotAnnotations.ntags == 2) {
1426         auto tag2iovs = PlotBase::getTag<1>().iovs;
1427         tagname2 = PlotBase::getTag<1>().name;
1428         lastiov = tag2iovs.front();
1429       } else {
1430         lastiov = theIOVs.back();
1431       }
1432 
1433       std::shared_ptr<SiStripApvGain> last_payload = this->fetchPayload(std::get<1>(lastiov));
1434       std::shared_ptr<SiStripApvGain> first_payload = this->fetchPayload(std::get<1>(firstiov));
1435 
1436       std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
1437       std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
1438 
1439       std::vector<uint32_t> detid;
1440       last_payload->getDetIds(detid);
1441 
1442       std::map<std::pair<uint32_t, int>, float> lastmap, firstmap;
1443 
1444       // loop on the last payload
1445       for (const auto& d : detid) {
1446         SiStripApvGain::Range range = last_payload->getRange(d);
1447         float Gain = 0;
1448         float nAPV = 0;
1449         for (int it = 0; it < range.second - range.first; ++it) {
1450           nAPV += 1;
1451           Gain = last_payload->getApvGain(it, range);
1452           std::pair<uint32_t, int> index = std::make_pair(d, nAPV);
1453           lastmap[index] = Gain;
1454         }  // end loop on APVs
1455       }  // end loop on detids
1456 
1457       detid.clear();
1458       first_payload->getDetIds(detid);
1459 
1460       // loop on the first payload
1461       for (const auto& d : detid) {
1462         SiStripApvGain::Range range = first_payload->getRange(d);
1463         float Gain = 0;
1464         float nAPV = 0;
1465         for (int it = 0; it < range.second - range.first; ++it) {
1466           nAPV += 1;
1467           Gain = first_payload->getApvGain(it, range);
1468           std::pair<uint32_t, int> index = std::make_pair(d, nAPV);
1469           firstmap[index] = Gain;
1470         }  // end loop on APVs
1471       }  // end loop on detids
1472 
1473       TCanvas canvas("Payload comparison", "payload comparison", 1400, 1000);
1474       canvas.Divide(2, 1);
1475 
1476       std::map<std::string, std::shared_ptr<TH1F>> ratios;
1477       std::map<std::string, std::shared_ptr<TH2F>> scatters;
1478       std::map<std::string, int> colormap;
1479       std::map<std::string, int> markermap;
1480       colormap["TIB"] = kRed;
1481       markermap["TIB"] = kFullCircle;
1482       colormap["TOB"] = kGreen;
1483       markermap["TOB"] = kFullTriangleUp;
1484       colormap["TID"] = kBlack;
1485       markermap["TID"] = kFullSquare;
1486       colormap["TEC"] = kBlue;
1487       markermap["TEC"] = kFullTriangleDown;
1488 
1489       std::vector<std::string> parts = {"TEC", "TOB", "TIB", "TID"};
1490 
1491       for (const auto& part : parts) {
1492         ratios[part] = std::make_shared<TH1F>(
1493             Form("hRatio_%s", part.c_str()),
1494             Form("Gains ratio IOV: %s/ IOV: %s ;Previous Gain (%s) / New Gain (%s);Number of APV",
1495                  firstIOVsince.c_str(),
1496                  lastIOVsince.c_str(),
1497                  firstIOVsince.c_str(),
1498                  lastIOVsince.c_str()),
1499             200,
1500             0.,
1501             2.);
1502         scatters[part] =
1503             std::make_shared<TH2F>(Form("hScatter_%s", part.c_str()),
1504                                    Form("new Gain (%s) vs previous Gain (%s);Previous Gain (%s);New Gain (%s)",
1505                                         lastIOVsince.c_str(),
1506                                         firstIOVsince.c_str(),
1507                                         firstIOVsince.c_str(),
1508                                         lastIOVsince.c_str()),
1509                                    100,
1510                                    0.5,
1511                                    1.8,
1512                                    100,
1513                                    0.5,
1514                                    1.8);
1515       }
1516 
1517       // now loop on the cached maps
1518       for (const auto& item : firstmap) {
1519         // packed index (detid,APV)
1520         auto index = item.first;
1521         auto mod = item.first.first;
1522 
1523         int subid = DetId(mod).subdetId();
1524         float ratio = firstmap[index] / lastmap[index];
1525 
1526         if (subid == StripSubdetector::TIB) {
1527           ratios["TIB"]->Fill(ratio);
1528           scatters["TIB"]->Fill(firstmap[index], lastmap[index]);
1529         }
1530 
1531         if (subid == StripSubdetector::TOB) {
1532           ratios["TOB"]->Fill(ratio);
1533           scatters["TOB"]->Fill(firstmap[index], lastmap[index]);
1534         }
1535 
1536         if (subid == StripSubdetector::TID) {
1537           ratios["TID"]->Fill(ratio);
1538           scatters["TID"]->Fill(firstmap[index], lastmap[index]);
1539         }
1540 
1541         if (subid == StripSubdetector::TEC) {
1542           ratios["TEC"]->Fill(ratio);
1543           scatters["TEC"]->Fill(firstmap[index], lastmap[index]);
1544         }
1545       }
1546 
1547       auto legend = TLegend(0.60, 0.8, 0.92, 0.95);
1548       legend.SetTextSize(0.05);
1549       canvas.cd(1)->SetLogy();
1550       canvas.cd(1)->SetTopMargin(0.05);
1551       canvas.cd(1)->SetLeftMargin(0.13);
1552       canvas.cd(1)->SetRightMargin(0.08);
1553 
1554       for (const auto& part : parts) {
1555         SiStripPI::makeNicePlotStyle(ratios[part].get());
1556         ratios[part]->SetMinimum(1.);
1557         ratios[part]->SetStats(false);
1558         ratios[part]->SetLineWidth(2);
1559         ratios[part]->SetLineColor(colormap[part]);
1560         if (part == "TEC")
1561           ratios[part]->Draw();
1562         else
1563           ratios[part]->Draw("same");
1564         legend.AddEntry(ratios[part].get(), part.c_str(), "L");
1565       }
1566 
1567       legend.Draw("same");
1568       SiStripPI::drawStatBox(ratios, colormap, parts);
1569 
1570       auto legend2 = TLegend(0.60, 0.8, 0.92, 0.95);
1571       legend2.SetTextSize(0.05);
1572       canvas.cd(2);
1573       canvas.cd(2)->SetTopMargin(0.05);
1574       canvas.cd(2)->SetLeftMargin(0.13);
1575       canvas.cd(2)->SetRightMargin(0.08);
1576 
1577       for (const auto& part : parts) {
1578         SiStripPI::makeNicePlotStyle(scatters[part].get());
1579         scatters[part]->SetStats(false);
1580         scatters[part]->SetMarkerColor(colormap[part]);
1581         scatters[part]->SetMarkerStyle(markermap[part]);
1582         scatters[part]->SetMarkerSize(0.5);
1583 
1584         auto temp = (TH2F*)(scatters[part]->Clone());
1585         temp->SetMarkerSize(1.3);
1586 
1587         if (part == "TEC")
1588           scatters[part]->Draw("SCAT");
1589         else
1590           scatters[part]->Draw("SCATsame");
1591 
1592         legend2.AddEntry(temp, part.c_str(), "P");
1593       }
1594 
1595       TLine diagonal(0.5, 0.5, 1.8, 1.8);
1596       diagonal.SetLineWidth(3);
1597       diagonal.SetLineStyle(2);
1598       diagonal.Draw("same");
1599 
1600       legend2.Draw("same");
1601 
1602       std::string fileName(this->m_imageFileName);
1603       canvas.SaveAs(fileName.c_str());
1604 
1605       return true;
1606     }
1607   };
1608 
1609   using SiStripApvGainsComparatorSingleTag = SiStripApvGainsComparatorBase<1, MULTI_IOV>;
1610   using SiStripApvGainsComparatorTwoTags = SiStripApvGainsComparatorBase<2, SINGLE_IOV>;
1611 
1612   /************************************************
1613     Plot stack of gain by region
1614   *************************************************/
1615 
1616   class SiStripApvGainsTHStack : public PlotImage<SiStripApvGain, SINGLE_IOV> {
1617   public:
1618     SiStripApvGainsTHStack()
1619         : PlotImage<SiStripApvGain, SINGLE_IOV>("Stack of SiStrip APV gains values"),
1620           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
1621               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath())} {}
1622 
1623     bool fill() override {
1624       //TH1F::SetDefaultSumw2(true);
1625       auto tag = PlotBase::getTag<0>();
1626       auto iov = tag.iovs.front();
1627 
1628       std::shared_ptr<SiStripApvGain> payload = this->fetchPayload(std::get<1>(iov));
1629       std::string IOVsince = std::to_string(std::get<0>(iov));
1630 
1631       std::map<partition, std::shared_ptr<TH1F>> h_gains;
1632 
1633       //std::vector< SiStripPI::TrackerRegion > regions;
1634       std::vector<partition> regions;
1635 
1636       std::vector<uint32_t> detid;
1637       payload->getDetIds(detid);
1638 
1639       // fill the vector of regions
1640       for (const auto& d : detid) {
1641         //auto region = this->getTheRegion(d);
1642         auto region = this->getThePartition(d);
1643         if (std::find(regions.begin(), regions.end(), region) == regions.end()) {
1644           regions.push_back(region);
1645         }
1646       }
1647 
1648       LogDebug("SiStripApvGainsTHStack") << "regions.size()=" << regions.size() << std::endl;
1649 
1650       for (const auto& r : regions) {
1651         //auto part = std::string(SiStripPI::regionType(r).second);
1652 
1653         auto part = std::string(this->partitionName(r));
1654 
1655         h_gains[r] = std::make_shared<TH1F>(Form("hGains_%s", part.c_str()),
1656                                             Form("Gains values for IOV: %s ;Gain;Number of APV", IOVsince.c_str()),
1657                                             100,
1658                                             0.5,
1659                                             1.5);
1660       }
1661 
1662       // loop on the payload
1663       for (const auto& d : detid) {
1664         //auto region = this->getTheRegion(d);
1665         auto region = this->getThePartition(d);
1666         SiStripApvGain::Range range = payload->getRange(d);
1667         for (int it = 0; it < range.second - range.first; ++it) {
1668           float gain = payload->getApvGain(it, range);
1669           h_gains[region]->Fill(gain);
1670         }  // end loop on APVs
1671       }  // end loop on detids
1672 
1673       TCanvas canvas("Payload breakout", "payload breakout", 1200, 800);
1674       canvas.Divide(2, 1);
1675 
1676       std::array<int, 6> colors = {{kRed, kBlue, kGreen, kCyan, 8, kMagenta}};
1677 
1678       THStack* hs = new THStack("hs", Form("Gains values for IOV: %s;Gain;Number of APV", IOVsince.c_str()));
1679       int colorCounter = 0;
1680       for (const auto& r : regions) {
1681         hs->Add(h_gains[r].get());
1682         SiStripPI::makeNicePlotStyle(h_gains[r].get());
1683         h_gains[r]->SetFillColor(colors[colorCounter]);
1684         //h_gains[r]->SetLineColor(colorCounter);
1685         h_gains[r]->SetLineWidth(2);
1686         colorCounter++;
1687       }
1688 
1689       TLegend legend = TLegend(0.60, 0.65, 0.95, 0.93);
1690       legend.SetTextSize(0.05);
1691       legend.SetHeader("Gain break-out", "C");  // option "C" allows to center the header
1692       for (const auto& r : regions) {
1693         auto part = std::string(this->partitionName(r));
1694         legend.AddEntry(h_gains[r].get(), part.c_str(), "F");
1695       }
1696 
1697       canvas.cd(1)->SetLogy();
1698       canvas.cd(1)->SetTopMargin(0.07);
1699       canvas.cd(1)->SetBottomMargin(0.10);
1700       canvas.cd(1)->SetLeftMargin(0.15);
1701       canvas.cd(1)->SetRightMargin(0.05);
1702       //      hs->Draw("NOSTACKB");
1703 
1704       int count(0);
1705       auto stack = hs->GetHists();
1706       double maximum = hs->GetMaximum("nostack");  //SiStripPI::getMaximum(stack);
1707 
1708       TLegend legend2 = TLegend(0.70, 0.65, 0.95, 0.93);
1709       legend2.SetTextSize(0.05);
1710       legend2.SetHeader("Partition", "C");  // option "C" allows to center the header
1711 
1712       for (const auto&& elem : *stack) {
1713         auto clone = (TH1F*)(elem->Clone(Form("hclone_%s", elem->GetName())));
1714         SiStripPI::makeNicePlotStyle(clone);
1715         clone->SetFillColor(0);
1716         clone->SetMarkerStyle(20);
1717         clone->SetLineColor(colors[count]);
1718         clone->SetMarkerColor(colors[count]);
1719         clone->SetMaximum(maximum * 10);
1720         TString candName = clone->GetName();
1721         legend2.AddEntry(clone, candName.ReplaceAll("hclone_hGains_", ""), "L");
1722         if (count == 0) {
1723           clone->Draw("HIST");
1724         } else {
1725           clone->Draw("HISTsame");
1726         }
1727         count++;
1728       }
1729 
1730       legend2.Draw("same");
1731 
1732       canvas.cd(2);  //->SetLogy();
1733       canvas.cd(2)->SetTopMargin(0.07);
1734       canvas.cd(2)->SetBottomMargin(0.10);
1735       canvas.cd(2)->SetLeftMargin(0.12);
1736       canvas.cd(2)->SetRightMargin(0.05);
1737       hs->Draw();
1738       // all graphics manipulations *after* drawing the stack!
1739       hs->GetYaxis()->SetMaxDigits(2);
1740       SiStripPI::makeNiceStyle<THStack>(hs);
1741       legend.Draw("same");
1742 
1743       std::string fileName(this->m_imageFileName);
1744       canvas.SaveAs(fileName.c_str());
1745 
1746       return true;
1747     }
1748 
1749   private:
1750     TrackerTopology m_trackerTopo;
1751     enum partition { TIB = 30, TIDP = 41, TIDM = 42, TOB = 50, TECP = 61, TECM = 62, END_OF_PARTS };
1752 
1753     const char* partitionName(partition part) {
1754       std::map<partition, const char*> mapping = {{partition::TIB, "TIB"},
1755                                                   {partition::TIDP, "TIPp"},
1756                                                   {partition::TIDM, "TIDm"},
1757                                                   {partition::TOB, "TOB"},
1758                                                   {partition::TECP, "TECp"},
1759                                                   {partition::TECM, "TECm"}};
1760 
1761       if (mapping.find(part) == mapping.end()) {
1762         throw cms::Exception("Invalid Partition passed");
1763       } else {
1764         return mapping[part];
1765       }
1766     }
1767 
1768     partition getThePartition(DetId detid) {
1769       int detNum = 0;
1770       int side = 0;
1771       switch (detid.subdetId()) {
1772         case StripSubdetector::TIB:
1773           detNum = 30;
1774           break;
1775         case StripSubdetector::TOB:
1776           detNum = 50;
1777           break;
1778         case StripSubdetector::TEC:
1779           // is this module in TEC+ or TEC-?
1780           side = m_trackerTopo.tecSide(detid);
1781           detNum = 60;
1782           break;
1783         case StripSubdetector::TID:
1784           // is this module in TID+ or TID-?
1785           side = m_trackerTopo.tidSide(detid);
1786           detNum = 40;
1787           break;
1788       }
1789 
1790       detNum += side;
1791       return static_cast<partition>(detNum);
1792     }
1793 
1794     SiStripPI::TrackerRegion getTheRegion(DetId detid) {
1795       int layer = 0;
1796       int stereo = 0;
1797       int detNum = 0;
1798 
1799       switch (detid.subdetId()) {
1800         case StripSubdetector::TIB:
1801           layer = m_trackerTopo.tibLayer(detid);
1802           stereo = m_trackerTopo.tibStereo(detid);
1803           detNum = 1000;
1804           break;
1805         case StripSubdetector::TOB:
1806           layer = m_trackerTopo.tobLayer(detid);
1807           stereo = m_trackerTopo.tobStereo(detid);
1808           detNum = 2000;
1809           break;
1810         case StripSubdetector::TEC:
1811           // is this module in TEC+ or TEC-?
1812           layer = m_trackerTopo.tecWheel(detid);
1813           stereo = m_trackerTopo.tecStereo(detid);
1814           detNum = 3000;
1815           break;
1816         case StripSubdetector::TID:
1817           // is this module in TID+ or TID-?
1818           layer = m_trackerTopo.tidWheel(detid);
1819           stereo = m_trackerTopo.tidStereo(detid);
1820           detNum = 4000;
1821           break;
1822       }
1823 
1824       detNum += layer * 10 + stereo;
1825       return static_cast<SiStripPI::TrackerRegion>(detNum);
1826     }
1827   };
1828 
1829   //*******************************************//
1830   // Compare Gains from 2 IOVs
1831   //******************************************//
1832 
1833   template <int ntags, IOVMultiplicity nIOVs>
1834   class SiStripApvGainsValuesComparatorBase : public PlotImage<SiStripApvGain, nIOVs, ntags> {
1835   public:
1836     SiStripApvGainsValuesComparatorBase()
1837         : PlotImage<SiStripApvGain, nIOVs, ntags>("Comparison of SiStrip APV gains values") {}
1838 
1839     bool fill() override {
1840       TH1F::SetDefaultSumw2(true);
1841 
1842       // trick to deal with the multi-ioved tag and two tag case at the same time
1843       auto theIOVs = PlotBase::getTag<0>().iovs;
1844       auto tagname1 = PlotBase::getTag<0>().name;
1845       std::string tagname2 = "";
1846       auto firstiov = theIOVs.front();
1847       SiStripPI::MetaData lastiov;
1848 
1849       // we don't support (yet) comparison with more than 2 tags
1850       assert(this->m_plotAnnotations.ntags < 3);
1851 
1852       if (this->m_plotAnnotations.ntags == 2) {
1853         auto tag2iovs = PlotBase::getTag<1>().iovs;
1854         tagname2 = PlotBase::getTag<1>().name;
1855         lastiov = tag2iovs.front();
1856       } else {
1857         lastiov = theIOVs.back();
1858       }
1859 
1860       std::shared_ptr<SiStripApvGain> last_payload = this->fetchPayload(std::get<1>(lastiov));
1861       std::shared_ptr<SiStripApvGain> first_payload = this->fetchPayload(std::get<1>(firstiov));
1862 
1863       std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
1864       std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
1865 
1866       std::vector<uint32_t> detid;
1867       last_payload->getDetIds(detid);
1868 
1869       std::map<std::pair<uint32_t, int>, float> lastmap, firstmap;
1870 
1871       // loop on the last payload
1872       for (const auto& d : detid) {
1873         SiStripApvGain::Range range = last_payload->getRange(d);
1874         float nAPV = 0;
1875         for (int it = 0; it < range.second - range.first; ++it) {
1876           nAPV += 1;
1877           auto index = std::make_pair(d, nAPV);
1878           lastmap[index] = last_payload->getApvGain(it, range);
1879         }  // end loop on APVs
1880       }  // end loop on detids
1881 
1882       detid.clear();
1883       first_payload->getDetIds(detid);
1884 
1885       // loop on the first payload
1886       for (const auto& d : detid) {
1887         SiStripApvGain::Range range = first_payload->getRange(d);
1888         float nAPV = 0;
1889         for (int it = 0; it < range.second - range.first; ++it) {
1890           nAPV += 1;
1891           auto index = std::make_pair(d, nAPV);
1892           firstmap[index] = last_payload->getApvGain(it, range);
1893         }  // end loop on APVs
1894       }  // end loop on detids
1895 
1896       TCanvas canvas("Payload comparison", "payload comparison", 1000, 1000);
1897       canvas.cd();
1898 
1899       TPad pad1("pad1", "pad1", 0, 0.3, 1, 1.0);
1900       pad1.SetBottomMargin(0.02);  // Upper and lower plot are joined
1901       pad1.SetTopMargin(0.1);
1902       pad1.SetRightMargin(0.05);
1903       pad1.SetLeftMargin(0.15);
1904       pad1.Draw();  // Draw the upper pad: pad1
1905       pad1.cd();    // pad1 becomes the current pad
1906 
1907       auto h_firstGains = std::make_shared<TH1F>("hFirstGains", ";SiStrip APV Gains;n. APVs", 200, 0.2, 1.8);
1908       auto h_lastGains = std::make_shared<TH1F>("hLastGains", ";SiStrip APV Gains;n. APVs", 200, 0.2, 1.8);
1909 
1910       for (const auto& item : firstmap) {
1911         h_firstGains->Fill(item.second);
1912       }
1913 
1914       for (const auto& item : lastmap) {
1915         h_lastGains->Fill(item.second);
1916       }
1917 
1918       SiStripPI::makeNicePlotStyle(h_lastGains.get());
1919       SiStripPI::makeNicePlotStyle(h_firstGains.get());
1920 
1921       TH1F* hratio = (TH1F*)h_firstGains->Clone("hratio");
1922 
1923       h_firstGains->SetLineColor(kRed);
1924       h_lastGains->SetLineColor(kBlue);
1925 
1926       h_firstGains->SetMarkerColor(kRed);
1927       h_lastGains->SetMarkerColor(kBlue);
1928 
1929       h_firstGains->SetMarkerSize(1.);
1930       h_lastGains->SetMarkerSize(1.);
1931 
1932       h_firstGains->SetLineWidth(2);
1933       h_lastGains->SetLineWidth(2);
1934 
1935       h_firstGains->SetMarkerStyle(20);
1936       h_lastGains->SetMarkerStyle(21);
1937 
1938       h_firstGains->GetXaxis()->SetLabelOffset(2.);
1939       h_lastGains->GetXaxis()->SetLabelOffset(2.);
1940 
1941       float max = (h_firstGains->GetMaximum() > h_lastGains->GetMaximum()) ? h_firstGains->GetMaximum()
1942                                                                            : h_lastGains->GetMaximum();
1943 
1944       h_firstGains->GetYaxis()->SetRangeUser(0., std::max(0., max * 1.20));
1945 
1946       h_firstGains->Draw("HIST");
1947       h_lastGains->Draw("HISTsame");
1948 
1949       TLegend legend = TLegend(0.30, 0.78, 0.95, 0.9);
1950       legend.SetHeader("#font[22]{SiStrip APV Gains Comparison}", "C");  // option "C" allows to center the header
1951       legend.AddEntry(h_firstGains.get(), ("payload: #color[2]{" + std::get<1>(firstiov) + "}").c_str(), "F");
1952       legend.AddEntry(h_lastGains.get(), ("payload: #color[4]{" + std::get<1>(lastiov) + "}").c_str(), "F");
1953       legend.SetTextSize(0.025);
1954       legend.Draw("same");
1955 
1956       auto ltx = TLatex();
1957       ltx.SetTextFont(62);
1958       ltx.SetTextSize(0.037);
1959       ltx.SetTextAlign(11);
1960       std::string ltxText;
1961       if (this->m_plotAnnotations.ntags == 2) {
1962         ltxText = fmt::sprintf(
1963             "#splitline{#color[2]{%s, %s} vs}{#color[4]{%s, %s}}", tagname1, firstIOVsince, tagname2, lastIOVsince);
1964       } else {
1965         ltxText = fmt::sprintf(
1966             "#splitline{%s}{IOV: #color[2]{%s} vs IOV: #color[4]{%s}}", tagname1, firstIOVsince, lastIOVsince);
1967       }
1968       ltx.DrawLatexNDC(gPad->GetLeftMargin(), 1 - gPad->GetTopMargin() + 0.04, ltxText.c_str());
1969 
1970       // lower plot will be in pad
1971       canvas.cd();  // Go back to the main canvas before defining pad2
1972       TPad pad2("pad2", "pad2", 0, 0.005, 1, 0.3);
1973       pad2.SetTopMargin(0.01);
1974       pad2.SetBottomMargin(0.2);
1975       pad2.SetRightMargin(0.05);
1976       pad2.SetLeftMargin(0.15);
1977       pad2.SetGridy();  // horizontal grid
1978       pad2.Draw();
1979       pad2.cd();  // pad2 becomes the current pad
1980 
1981       // Define the ratio plot
1982       hratio->SetLineColor(kBlack);
1983       hratio->SetMarkerColor(kBlack);
1984       hratio->SetTitle("");
1985       hratio->SetMinimum(0.55);  // Define Y ..
1986       hratio->SetMaximum(1.55);  // .. range
1987       hratio->SetStats(false);   // No statistics on lower plot
1988       hratio->Divide(h_lastGains.get());
1989       hratio->SetMarkerStyle(20);
1990       hratio->Draw("ep");  // Draw the ratio plot
1991 
1992       // Y axis ratio plot settings
1993       hratio->GetYaxis()->SetTitle(
1994           ("ratio " + std::to_string(std::get<0>(firstiov)) + " / " + std::to_string(std::get<0>(lastiov))).c_str());
1995 
1996       hratio->GetYaxis()->SetNdivisions(505);
1997 
1998       SiStripPI::makeNicePlotStyle(hratio);
1999 
2000       hratio->GetYaxis()->SetTitleSize(25);
2001       hratio->GetXaxis()->SetLabelSize(25);
2002 
2003       hratio->GetYaxis()->SetTitleFont(43);
2004       hratio->GetYaxis()->SetTitleOffset(2.5);
2005       hratio->GetYaxis()->SetLabelFont(43);  // Absolute font size in pixel (precision 3)
2006       hratio->GetYaxis()->SetLabelSize(25);
2007 
2008       // X axis ratio plot settings
2009       hratio->GetXaxis()->SetTitleSize(30);
2010       hratio->GetXaxis()->SetTitleFont(43);
2011       hratio->GetXaxis()->SetTitle("SiStrip APV Gains");
2012       hratio->GetXaxis()->SetLabelFont(43);  // Absolute font size in pixel (precision 3)
2013       hratio->GetXaxis()->SetTitleOffset(3.);
2014 
2015       std::string fileName(this->m_imageFileName);
2016       canvas.SaveAs(fileName.c_str());
2017 
2018       return true;
2019     }
2020   };
2021 
2022   using SiStripApvGainsValuesComparatorSingleTag = SiStripApvGainsValuesComparatorBase<1, MULTI_IOV>;
2023   using SiStripApvGainsValuesComparatorTwoTags = SiStripApvGainsValuesComparatorBase<2, SINGLE_IOV>;
2024 
2025   //*******************************************//
2026   // Compare Gains ratio from 2 IOVs, region by region
2027   //******************************************//
2028 
2029   template <int ntags, IOVMultiplicity nIOVs>
2030   class SiStripApvGainsRatioComparatorByRegionBase : public PlotImage<SiStripApvGain, nIOVs, ntags> {
2031   public:
2032     SiStripApvGainsRatioComparatorByRegionBase()
2033         : PlotImage<SiStripApvGain, nIOVs, ntags>("Module by Module Comparison of SiStrip APV gains"),
2034           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
2035               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath())} {}
2036 
2037     bool fill() override {
2038       //gStyle->SetPalette(5);
2039       SiStripPI::setPaletteStyle(SiStripPI::GRAY);
2040 
2041       // trick to deal with the multi-ioved tag and two tag case at the same time
2042       auto theIOVs = PlotBase::getTag<0>().iovs;
2043       auto tagname1 = PlotBase::getTag<0>().name;
2044       std::string tagname2 = "";
2045       auto firstiov = theIOVs.front();
2046       SiStripPI::MetaData lastiov;
2047 
2048       // we don't support (yet) comparison with more than 2 tags
2049       assert(this->m_plotAnnotations.ntags < 3);
2050 
2051       if (this->m_plotAnnotations.ntags == 2) {
2052         auto tag2iovs = PlotBase::getTag<1>().iovs;
2053         tagname2 = PlotBase::getTag<1>().name;
2054         lastiov = tag2iovs.front();
2055       } else {
2056         lastiov = theIOVs.back();
2057       }
2058 
2059       std::shared_ptr<SiStripApvGain> last_payload = this->fetchPayload(std::get<1>(lastiov));
2060       std::shared_ptr<SiStripApvGain> first_payload = this->fetchPayload(std::get<1>(firstiov));
2061 
2062       std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
2063       std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
2064 
2065       std::vector<uint32_t> detid;
2066       last_payload->getDetIds(detid);
2067 
2068       std::map<std::pair<uint32_t, int>, float> lastmap, firstmap;
2069 
2070       // loop on the last payload
2071       for (const auto& d : detid) {
2072         SiStripApvGain::Range range = last_payload->getRange(d);
2073         float Gain = 0;
2074         float nAPV = 0;
2075         for (int it = 0; it < range.second - range.first; ++it) {
2076           nAPV += 1;
2077           Gain = last_payload->getApvGain(it, range);
2078           std::pair<uint32_t, int> index = std::make_pair(d, nAPV);
2079           lastmap[index] = Gain;
2080         }  // end loop on APVs
2081       }  // end loop on detids
2082 
2083       detid.clear();
2084       first_payload->getDetIds(detid);
2085 
2086       // loop on the first payload
2087       for (const auto& d : detid) {
2088         SiStripApvGain::Range range = first_payload->getRange(d);
2089         float Gain = 0;
2090         float nAPV = 0;
2091         for (int it = 0; it < range.second - range.first; ++it) {
2092           nAPV += 1;
2093           Gain = first_payload->getApvGain(it, range);
2094           std::pair<uint32_t, int> index = std::make_pair(d, nAPV);
2095           firstmap[index] = Gain;
2096         }  // end loop on APVs
2097       }  // end loop on detids
2098 
2099       TCanvas canvas("Payload comparison by Tracker Region", "payload comparison by Tracker Region", 1800, 800);
2100       canvas.Divide(2, 1);
2101 
2102       auto h2first = std::make_unique<TH2F>(
2103           "byRegion1", "SiStrip APV Gain values by region;; average SiStrip Gain", 38, 1., 39., 100., 0., 2.);
2104       auto h2last = std::make_unique<TH2F>(
2105           "byRegion2", "SiStrip APV Gain values by region;; average SiStrip Gain", 38, 1., 39., 100., 0., 2.);
2106 
2107       auto h2ratio = std::make_unique<TH2F>("byRegionRatio",
2108                                             Form("SiStrip APV Gains ratio by region;; Gains ratio IOV: %s/ IOV %s",
2109                                                  lastIOVsince.c_str(),
2110                                                  firstIOVsince.c_str()),
2111                                             38,
2112                                             1.,
2113                                             39.,
2114                                             100.,
2115                                             0.85,
2116                                             1.15);
2117 
2118       h2first->SetStats(false);
2119       h2last->SetStats(false);
2120       h2ratio->SetStats(false);
2121 
2122       canvas.cd(1)->SetBottomMargin(0.18);
2123       canvas.cd(1)->SetLeftMargin(0.12);
2124       canvas.cd(1)->SetRightMargin(0.05);
2125       canvas.Modified();
2126 
2127       std::vector<int> boundaries;
2128       std::string detector;
2129       std::string currentDetector;
2130 
2131       for (const auto& element : lastmap) {
2132         auto region = this->getTheRegion(element.first.first);
2133         auto bin = SiStripPI::regionType(region).first;
2134         auto label = SiStripPI::regionType(region).second;
2135 
2136         h2last->Fill(bin, element.second);
2137         h2last->GetXaxis()->SetBinLabel(bin, label);
2138         h2ratio->Fill(bin, element.second / firstmap[element.first]);
2139         h2ratio->GetXaxis()->SetBinLabel(bin, label);
2140       }
2141 
2142       for (const auto& element : firstmap) {
2143         auto region = this->getTheRegion(element.first.first);
2144         auto bin = SiStripPI::regionType(region).first;
2145         auto label = SiStripPI::regionType(region).second;
2146 
2147         h2first->Fill(bin, element.second);
2148         h2first->GetXaxis()->SetBinLabel(bin, label);
2149       }
2150 
2151       h2first->GetXaxis()->LabelsOption("v");
2152       h2last->GetXaxis()->LabelsOption("v");
2153       h2ratio->GetXaxis()->LabelsOption("v");
2154 
2155       h2last->SetLineColor(kBlue);
2156       h2first->SetLineColor(kRed);
2157       h2first->SetFillColor(kRed);
2158 
2159       h2first->SetMarkerStyle(20);
2160       h2last->SetMarkerStyle(21);
2161 
2162       h2first->SetMarkerColor(kRed);
2163       h2last->SetMarkerColor(kBlue);
2164 
2165       canvas.cd(1);
2166       h2first->Draw("BOX");
2167       h2last->Draw("BOXsame");
2168 
2169       TLegend legend = TLegend(0.70, 0.8, 0.95, 0.9);
2170       legend.SetHeader("Gain Comparison", "C");  // option "C" allows to center the header
2171       legend.AddEntry(h2first.get(), ("IOV: " + std::to_string(std::get<0>(firstiov))).c_str(), "F");
2172       legend.AddEntry(h2last.get(), ("IOV: " + std::to_string(std::get<0>(lastiov))).c_str(), "F");
2173       legend.Draw("same");
2174 
2175       canvas.cd(2);
2176       canvas.cd(2)->SetBottomMargin(0.18);
2177       canvas.cd(2)->SetLeftMargin(0.12);
2178       canvas.cd(2)->SetRightMargin(0.12);
2179 
2180       h2ratio->Draw("COLZ");
2181       auto hpfx_tmp = (TProfile*)(h2ratio->ProfileX("_pfx", 1, -1, "o"));
2182       hpfx_tmp->SetStats(kFALSE);
2183       hpfx_tmp->SetMarkerColor(kRed);
2184       hpfx_tmp->SetLineColor(kRed);
2185       hpfx_tmp->SetMarkerSize(1.2);
2186       hpfx_tmp->SetMarkerStyle(20);
2187       hpfx_tmp->Draw("same");
2188 
2189       std::string fileName(this->m_imageFileName);
2190       canvas.SaveAs(fileName.c_str());
2191 
2192       delete hpfx_tmp;
2193       return true;
2194     }
2195 
2196   private:
2197     TrackerTopology m_trackerTopo;
2198 
2199     SiStripPI::TrackerRegion getTheRegion(DetId detid) {
2200       int layer = 0;
2201       int stereo = 0;
2202       int detNum = 0;
2203 
2204       switch (detid.subdetId()) {
2205         case StripSubdetector::TIB:
2206           layer = m_trackerTopo.tibLayer(detid);
2207           stereo = m_trackerTopo.tibStereo(detid);
2208           detNum = 1000;
2209           break;
2210         case StripSubdetector::TOB:
2211           layer = m_trackerTopo.tobLayer(detid);
2212           stereo = m_trackerTopo.tobStereo(detid);
2213           detNum = 2000;
2214           break;
2215         case StripSubdetector::TEC:
2216           // is this module in TEC+ or TEC-?
2217           layer = m_trackerTopo.tecWheel(detid);
2218           stereo = m_trackerTopo.tecStereo(detid);
2219           detNum = 3000;
2220           break;
2221         case StripSubdetector::TID:
2222           // is this module in TID+ or TID-?
2223           layer = m_trackerTopo.tidWheel(detid);
2224           stereo = m_trackerTopo.tidStereo(detid);
2225           detNum = 4000;
2226           break;
2227       }
2228 
2229       detNum += layer * 10 + stereo;
2230       return static_cast<SiStripPI::TrackerRegion>(detNum);
2231     }
2232   };
2233 
2234   using SiStripApvGainsRatioComparatorByRegionSingleTag = SiStripApvGainsRatioComparatorByRegionBase<1, MULTI_IOV>;
2235   using SiStripApvGainsRatioComparatorByRegionTwoTags = SiStripApvGainsRatioComparatorByRegionBase<2, SINGLE_IOV>;
2236 
2237   /************************************************
2238     Compare Gains for each tracker region
2239   *************************************************/
2240 
2241   template <int ntags, IOVMultiplicity nIOVs>
2242   class SiStripApvGainsComparatorByRegionBase : public PlotImage<SiStripApvGain, nIOVs, ntags> {
2243   public:
2244     SiStripApvGainsComparatorByRegionBase()
2245         : PlotImage<SiStripApvGain, nIOVs, ntags>("SiStripGains Comparison By Region"),
2246           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
2247               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath())} {}
2248 
2249     bool fill() override {
2250       // trick to deal with the multi-ioved tag and two tag case at the same time
2251       auto theIOVs = PlotBase::getTag<0>().iovs;
2252       auto tagname1 = PlotBase::getTag<0>().name;
2253       std::string tagname2 = "";
2254       auto firstiov = theIOVs.front();
2255       SiStripPI::MetaData lastiov;
2256 
2257       // we don't support (yet) comparison with more than 2 tags
2258       assert(this->m_plotAnnotations.ntags < 3);
2259 
2260       if (this->m_plotAnnotations.ntags == 2) {
2261         auto tag2iovs = PlotBase::getTag<1>().iovs;
2262         tagname2 = PlotBase::getTag<1>().name;
2263         lastiov = tag2iovs.front();
2264       } else {
2265         lastiov = theIOVs.back();
2266       }
2267 
2268       std::shared_ptr<SiStripApvGain> last_payload = this->fetchPayload(std::get<1>(lastiov));
2269       std::shared_ptr<SiStripApvGain> first_payload = this->fetchPayload(std::get<1>(firstiov));
2270 
2271       std::vector<uint32_t> detid;
2272       last_payload->getDetIds(detid);
2273 
2274       SiStripDetSummary summaryLastGain{&m_trackerTopo};
2275 
2276       for (const auto& d : detid) {
2277         SiStripApvGain::Range range = last_payload->getRange(d);
2278         for (int it = 0; it < range.second - range.first; ++it) {
2279           summaryLastGain.add(d, last_payload->getApvGain(it, range));
2280         }
2281       }
2282 
2283       SiStripDetSummary summaryFirstGain{&m_trackerTopo};
2284 
2285       for (const auto& d : detid) {
2286         SiStripApvGain::Range range = first_payload->getRange(d);
2287         for (int it = 0; it < range.second - range.first; ++it) {
2288           summaryFirstGain.add(d, first_payload->getApvGain(it, range));
2289         }
2290       }
2291 
2292       std::map<unsigned int, SiStripDetSummary::Values> firstmap = summaryFirstGain.getCounts();
2293       std::map<unsigned int, SiStripDetSummary::Values> lastmap = summaryLastGain.getCounts();
2294       //=========================
2295 
2296       TCanvas canvas("Region summary", "region summary", 1200, 1000);
2297       canvas.cd();
2298 
2299       auto hfirst = std::make_unique<TH1F>("byRegion1",
2300                                            "SiStrip APV Gain average by region;; average SiStrip Gain",
2301                                            firstmap.size(),
2302                                            0.,
2303                                            firstmap.size());
2304       auto hlast = std::make_unique<TH1F>(
2305           "byRegion2", "SiStrip APV Gain average by region;; average SiStrip Gain", lastmap.size(), 0., lastmap.size());
2306 
2307       hfirst->SetStats(false);
2308       hlast->SetStats(false);
2309 
2310       canvas.SetBottomMargin(0.18);
2311       canvas.SetLeftMargin(0.12);
2312       canvas.SetRightMargin(0.05);
2313       canvas.Modified();
2314 
2315       std::vector<int> boundaries;
2316       unsigned int iBin = 0;
2317 
2318       std::string detector;
2319       std::string currentDetector;
2320 
2321       for (const auto& element : lastmap) {
2322         iBin++;
2323         int count = element.second.count;
2324         double mean = (element.second.mean) / count;
2325 
2326         if (currentDetector.empty())
2327           currentDetector = "TIB";
2328 
2329         switch ((element.first) / 1000) {
2330           case 1:
2331             detector = "TIB";
2332             break;
2333           case 2:
2334             detector = "TOB";
2335             break;
2336           case 3:
2337             detector = "TEC";
2338             break;
2339           case 4:
2340             detector = "TID";
2341             break;
2342         }
2343 
2344         hlast->SetBinContent(iBin, mean);
2345         hlast->SetBinError(iBin, mean / 10000.);
2346         hlast->GetXaxis()->SetBinLabel(iBin, SiStripPI::regionType(element.first).second);
2347         hlast->GetXaxis()->LabelsOption("v");
2348 
2349         if (detector != currentDetector) {
2350           boundaries.push_back(iBin);
2351           currentDetector = detector;
2352         }
2353       }
2354 
2355       // reset the count
2356       iBin = 0;
2357 
2358       for (const auto& element : firstmap) {
2359         iBin++;
2360         int count = element.second.count;
2361         double mean = (element.second.mean) / count;
2362 
2363         hfirst->SetBinContent(iBin, mean);
2364         hfirst->SetBinError(iBin, mean / 10000.);
2365         hfirst->GetXaxis()->SetBinLabel(iBin, SiStripPI::regionType(element.first).second);
2366         hfirst->GetXaxis()->LabelsOption("v");
2367       }
2368 
2369       auto extrema = SiStripPI::getExtrema(hfirst.get(), hlast.get());
2370       hlast->GetYaxis()->SetRangeUser(extrema.first, extrema.second);
2371 
2372       hlast->SetMarkerStyle(20);
2373       hlast->SetMarkerSize(1);
2374       hlast->Draw("E1");
2375       hlast->Draw("Psame");
2376 
2377       hfirst->SetMarkerStyle(18);
2378       hfirst->SetMarkerSize(1);
2379       hfirst->SetLineColor(kBlue);
2380       hfirst->SetMarkerColor(kBlue);
2381       hfirst->Draw("E1same");
2382       hfirst->Draw("Psame");
2383 
2384       canvas.Update();
2385       canvas.cd();
2386 
2387       TLine l[boundaries.size()];
2388       unsigned int i = 0;
2389       for (const auto& line : boundaries) {
2390         l[i] = TLine(
2391             hfirst->GetBinLowEdge(line), canvas.cd()->GetUymin(), hfirst->GetBinLowEdge(line), canvas.cd()->GetUymax());
2392         l[i].SetLineWidth(1);
2393         l[i].SetLineStyle(9);
2394         l[i].SetLineColor(2);
2395         l[i].Draw("same");
2396         i++;
2397       }
2398 
2399       TLegend legend = TLegend(0.70, 0.8, 0.95, 0.9);
2400       legend.SetHeader("Gain Comparison", "C");  // option "C" allows to center the header
2401       legend.AddEntry(hfirst.get(), ("IOV: " + std::to_string(std::get<0>(firstiov))).c_str(), "PL");
2402       legend.AddEntry(hlast.get(), ("IOV: " + std::to_string(std::get<0>(lastiov))).c_str(), "PL");
2403       legend.Draw("same");
2404 
2405       std::string fileName(this->m_imageFileName);
2406       canvas.SaveAs(fileName.c_str());
2407 
2408       return true;
2409     }
2410 
2411   private:
2412     TrackerTopology m_trackerTopo;
2413   };
2414 
2415   using SiStripApvGainsComparatorByRegionSingleTag = SiStripApvGainsComparatorByRegionBase<1, MULTI_IOV>;
2416   using SiStripApvGainsComparatorByRegionTwoTags = SiStripApvGainsComparatorByRegionBase<2, SINGLE_IOV>;
2417 
2418   /************************************************
2419     Plot gain averages by region 
2420   *************************************************/
2421 
2422   class SiStripApvGainsByRegion : public PlotImage<SiStripApvGain, SINGLE_IOV> {
2423   public:
2424     SiStripApvGainsByRegion()
2425         : PlotImage<SiStripApvGain, SINGLE_IOV>("SiStripGains By Region"),
2426           m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(
2427               edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath())} {}
2428 
2429     bool fill() override {
2430       auto tag = PlotBase::getTag<0>();
2431       auto iov = tag.iovs.front();
2432       std::shared_ptr<SiStripApvGain> payload = fetchPayload(std::get<1>(iov));
2433 
2434       std::vector<uint32_t> detid;
2435       payload->getDetIds(detid);
2436 
2437       SiStripDetSummary summaryGain{&m_trackerTopo};
2438 
2439       for (const auto& d : detid) {
2440         SiStripApvGain::Range range = payload->getRange(d);
2441         for (int it = 0; it < range.second - range.first; ++it) {
2442           summaryGain.add(d, payload->getApvGain(it, range));
2443         }
2444       }
2445 
2446       std::map<unsigned int, SiStripDetSummary::Values> map = summaryGain.getCounts();
2447       //=========================
2448 
2449       TCanvas canvas("Region summary", "region summary", 1200, 1000);
2450       canvas.cd();
2451       auto h1 = std::make_unique<TH1F>(
2452           "byRegion", "SiStrip Gain average by region;; average SiStrip Gain", map.size(), 0., map.size());
2453       h1->SetStats(false);
2454       canvas.SetBottomMargin(0.18);
2455       canvas.SetLeftMargin(0.12);
2456       canvas.SetRightMargin(0.05);
2457       canvas.Modified();
2458 
2459       std::vector<int> boundaries;
2460       unsigned int iBin = 0;
2461 
2462       std::string detector;
2463       std::string currentDetector;
2464 
2465       for (const auto& element : map) {
2466         iBin++;
2467         int count = element.second.count;
2468         double mean = (element.second.mean) / count;
2469 
2470         if (currentDetector.empty())
2471           currentDetector = "TIB";
2472 
2473         switch ((element.first) / 1000) {
2474           case 1:
2475             detector = "TIB";
2476             break;
2477           case 2:
2478             detector = "TOB";
2479             break;
2480           case 3:
2481             detector = "TEC";
2482             break;
2483           case 4:
2484             detector = "TID";
2485             break;
2486         }
2487 
2488         h1->SetBinContent(iBin, mean);
2489         h1->GetXaxis()->SetBinLabel(iBin, SiStripPI::regionType(element.first).second);
2490         h1->GetXaxis()->LabelsOption("v");
2491 
2492         if (detector != currentDetector) {
2493           boundaries.push_back(iBin);
2494           currentDetector = detector;
2495         }
2496       }
2497 
2498       h1->SetMarkerStyle(20);
2499       h1->SetMarkerSize(1);
2500       h1->Draw("HIST");
2501       h1->Draw("Psame");
2502 
2503       canvas.Update();
2504 
2505       TLine l[boundaries.size()];
2506       unsigned int i = 0;
2507       for (const auto& line : boundaries) {
2508         l[i] = TLine(h1->GetBinLowEdge(line), canvas.GetUymin(), h1->GetBinLowEdge(line), canvas.GetUymax());
2509         l[i].SetLineWidth(1);
2510         l[i].SetLineStyle(9);
2511         l[i].SetLineColor(2);
2512         l[i].Draw("same");
2513         i++;
2514       }
2515 
2516       TLegend legend = TLegend(0.52, 0.82, 0.95, 0.9);
2517       legend.SetHeader((std::get<1>(iov)).c_str(), "C");  // option "C" allows to center the header
2518       legend.AddEntry(h1.get(), ("IOV: " + std::to_string(std::get<0>(iov))).c_str(), "PL");
2519       legend.SetTextSize(0.025);
2520       legend.Draw("same");
2521 
2522       std::string fileName(m_imageFileName);
2523       canvas.SaveAs(fileName.c_str());
2524 
2525       return true;
2526     }
2527 
2528   private:
2529     TrackerTopology m_trackerTopo;
2530   };
2531 
2532 }  // namespace
2533 
2534 // Register the classes as boost python plugin
2535 PAYLOAD_INSPECTOR_MODULE(SiStripApvGain) {
2536   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsValue);
2537   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainTest);
2538   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainByPartition);
2539   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainCompareByPartitionSingleTag);
2540   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainRatioByPartitionSingleTag);
2541   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainDiffByPartitionSingleTag);
2542   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainCompareByPartitionTwoTags);
2543   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainRatioByPartitionTwoTags);
2544   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainDiffByPartitionTwoTags);
2545   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsTest);
2546   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsByRegion);
2547   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsComparatorSingleTag);
2548   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsComparatorTwoTags);
2549   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsValuesComparatorSingleTag);
2550   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsValuesComparatorTwoTags);
2551   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsComparatorByRegionSingleTag);
2552   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsComparatorByRegionTwoTags);
2553   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsRatioComparatorByRegionSingleTag);
2554   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsRatioComparatorByRegionTwoTags);
2555   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsTHStack);
2556   PAYLOAD_INSPECTOR_CLASS(SiStripApvBarrelGainsByLayer);
2557   PAYLOAD_INSPECTOR_CLASS(SiStripApvAbsoluteBarrelGainsByLayer);
2558   PAYLOAD_INSPECTOR_CLASS(SiStripApvEndcapMinusGainsByDisk);
2559   PAYLOAD_INSPECTOR_CLASS(SiStripApvEndcapPlusGainsByDisk);
2560   PAYLOAD_INSPECTOR_CLASS(SiStripApvAbsoluteEndcapMinusGainsByDisk);
2561   PAYLOAD_INSPECTOR_CLASS(SiStripApvAbsoluteEndcapPlusGainsByDisk);
2562   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsAverageTrackerMap);
2563   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsDefaultTrackerMap);
2564   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsMaximumTrackerMap);
2565   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsMinimumTrackerMap);
2566   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsAvgDeviationRatioWithPreviousIOVTrackerMap);
2567   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsAvgDeviationRatioTrackerMapTwoTags);
2568   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsMaxDeviationRatioWithPreviousIOVTrackerMap);
2569   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainsMaxDeviationRatioTrackerMapTwoTags);
2570   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainByRunMeans);
2571   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainMin_History);
2572   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainMax_History);
2573   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainMean_History);
2574   PAYLOAD_INSPECTOR_CLASS(SiStripApvGainRMS_History);
2575   PAYLOAD_INSPECTOR_CLASS(SiStripApvTIBGainByRunMeans);
2576   PAYLOAD_INSPECTOR_CLASS(SiStripApvTIDGainByRunMeans);
2577   PAYLOAD_INSPECTOR_CLASS(SiStripApvTOBGainByRunMeans);
2578   PAYLOAD_INSPECTOR_CLASS(SiStripApvTECGainByRunMeans);
2579 }