Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef CONDCORE_SISTRIPPLUGINS_SISTRIPPAYLOADINSPECTORHELPER_H
0002 #define CONDCORE_SISTRIPPLUGINS_SISTRIPPAYLOADINSPECTORHELPER_H
0003 
0004 // system includes
0005 #include <numeric>
0006 #include <string>
0007 #include <vector>
0008 
0009 // user includes
0010 #include "CalibFormats/SiStripObjects/interface/SiStripQuality.h"
0011 #include "CalibTracker/SiStripCommon/interface/SiStripDetInfoFileReader.h"
0012 #include "CalibTracker/StandaloneTrackerTopology/interface/StandaloneTrackerTopology.h"
0013 #include "CondCore/Utilities/interface/PayloadInspector.h"
0014 #include "CondFormats/SiStripObjects/interface/SiStripDetSummary.h"
0015 #include "CondFormats/SiStripObjects/interface/SiStripNoises.h"
0016 #include "CondFormats/SiStripObjects/interface/SiStripSummary.h"
0017 #include "DataFormats/SiStripCommon/interface/ConstantsForHardwareSystems.h" /* for STRIPS_PER_APV*/
0018 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
0019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0020 #include "FWCore/ParameterSet/interface/FileInPath.h"
0021 
0022 // ROOT includes
0023 #include "TH1.h"
0024 #include "TH2.h"
0025 #include "TObjArray.h"
0026 #include "TPaveText.h"
0027 #include "TStyle.h"
0028 
0029 namespace SiStripPI {
0030 
0031   //##### for metadata
0032   using MetaData = std::tuple<cond::Time_t, cond::Hash>;
0033 
0034   //##### for plotting
0035 
0036   enum OpMode { STRIP_BASED, APV_BASED, MODULE_BASED };
0037 
0038   class Entry {
0039   public:
0040     Entry() : entries(0), sum(0), sq_sum(0) {}
0041 
0042     double mean() { return sum / entries; }
0043     double std_dev() {
0044       double tmean = mean();
0045       return sqrt((sq_sum - entries * tmean * tmean) / (entries - 1));
0046     }
0047     double mean_rms() { return std_dev() / sqrt(entries); }
0048 
0049     void add(double val) {
0050       entries++;
0051       sum += val;
0052       sq_sum += val * val;
0053     }
0054 
0055     void reset() {
0056       entries = 0;
0057       sum = 0;
0058       sq_sum = 0;
0059     }
0060 
0061   private:
0062     long int entries;
0063     double sum, sq_sum;
0064   };
0065 
0066   // class Monitor1D
0067 
0068   class Monitor1D {
0069   public:
0070     Monitor1D(OpMode mode, const char* name, const char* title, int nbinsx, double xmin, double xmax)
0071         : entry_(), mode_(mode), obj_(name, title, nbinsx, xmin, xmax) {}
0072 
0073     Monitor1D() : entry_(), mode_(OpMode::STRIP_BASED), obj_() {}
0074 
0075     ~Monitor1D() {}
0076 
0077     void Fill(int apv, int det, double vx) {
0078       switch (mode_) {
0079         case (OpMode::APV_BASED):
0080           if (!((apv == prev_apv_ && det == prev_det_) || prev_apv_ == 0)) {
0081             flush();
0082           }
0083           prev_apv_ = apv;
0084           prev_det_ = det;
0085           break;
0086         case (OpMode::MODULE_BASED):
0087           if (!(det == prev_det_ || prev_det_ == 0)) {
0088             flush();
0089           }
0090           prev_det_ = det;
0091           break;
0092         case (OpMode::STRIP_BASED):
0093           flush();
0094           break;
0095       }
0096       entry_.add(vx);
0097     }
0098 
0099     void flush() {
0100       obj_.Fill(entry_.mean());
0101       entry_.reset();
0102     }
0103 
0104     TH1F& hist() {
0105       flush();
0106       return obj_;
0107     }
0108 
0109     TH1F& getHist() { return obj_; }
0110 
0111   private:
0112     int prev_apv_ = 0, prev_det_ = 0;
0113     Entry entry_;
0114     OpMode mode_;
0115     TH1F obj_;
0116   };
0117 
0118   // class monitor 2D
0119 
0120   class Monitor2D {
0121   public:
0122     Monitor2D(OpMode mode,
0123               const char* name,
0124               const char* title,
0125               int nbinsx,
0126               double xmin,
0127               double xmax,
0128               int nbinsy,
0129               double ymin,
0130               double ymax)
0131         : entryx_(), entryy_(), mode_(mode), obj_(name, title, nbinsx, xmin, xmax, nbinsy, ymin, ymax) {}
0132 
0133     Monitor2D() : entryx_(), entryy_(), mode_(OpMode::STRIP_BASED), obj_() {}
0134 
0135     ~Monitor2D() {}
0136 
0137     void Fill(int apv, int det, double vx, double vy) {
0138       switch (mode_) {
0139         case (OpMode::APV_BASED):
0140           if (!((apv == prev_apv_ && det == prev_det_) || prev_apv_ == 0)) {
0141             flush();
0142           }
0143           prev_apv_ = apv;
0144           prev_det_ = det;
0145           break;
0146         case (OpMode::MODULE_BASED):
0147           if (!(det == prev_det_ || prev_det_ == 0)) {
0148             flush();
0149           }
0150           prev_det_ = det;
0151           break;
0152         case (OpMode::STRIP_BASED):
0153           flush();
0154           break;
0155       }
0156       entryx_.add(vx);
0157       entryy_.add(vy);
0158     }
0159 
0160     void flush() {
0161       obj_.Fill(entryx_.mean(), entryy_.mean());
0162       entryx_.reset();
0163       entryy_.reset();
0164     }
0165 
0166     TH2F& hist() {
0167       flush();
0168       return obj_;
0169     }
0170 
0171   private:
0172     int prev_apv_ = 0, prev_det_ = 0;
0173     Entry entryx_, entryy_;
0174     OpMode mode_;
0175     TH2F obj_;
0176   };
0177 
0178   enum estimator { min, max, mean, rms };
0179 
0180   /*--------------------------------------------------------------------*/
0181   inline std::string estimatorType(SiStripPI::estimator e)
0182   /*--------------------------------------------------------------------*/
0183   {
0184     switch (e) {
0185       case SiStripPI::min:
0186         return "minimum";
0187       case SiStripPI::max:
0188         return "maximum";
0189       case SiStripPI::mean:
0190         return "mean";
0191       case SiStripPI::rms:
0192         return "RMS";
0193       default:
0194         return "should never be here";
0195     }
0196   }
0197 
0198   /*--------------------------------------------------------------------*/
0199   inline std::string getStringFromSubdet(StripSubdetector::SubDetector sub)
0200   /*-------------------------------------------------------------------*/
0201   {
0202     switch (sub) {
0203       case StripSubdetector::TIB:
0204         return "TIB";
0205       case StripSubdetector::TOB:
0206         return "TOB";
0207       case StripSubdetector::TID:
0208         return "TID";
0209       case StripSubdetector::TEC:
0210         return "TEC";
0211       default:
0212         return "should never be here";
0213     }
0214   }
0215 
0216   enum TrackerRegion {
0217     TIB1r = 1010,
0218     TIB1s = 1011,
0219     TIB2r = 1020,
0220     TIB2s = 1021,
0221     TIB3r = 1030,
0222     TIB4r = 1040,
0223     TOB1r = 2010,
0224     TOB1s = 2011,
0225     TOB2r = 2020,
0226     TOB2s = 2021,
0227     TOB3r = 2030,
0228     TOB4r = 2040,
0229     TOB5r = 2050,
0230     TOB6r = 2060,
0231     TEC1r = 3010,
0232     TEC1s = 3011,
0233     TEC2r = 3020,
0234     TEC2s = 3021,
0235     TEC3r = 3030,
0236     TEC3s = 3031,
0237     TEC4r = 3040,
0238     TEC4s = 3041,
0239     TEC5r = 3050,
0240     TEC5s = 3051,
0241     TEC6r = 3060,
0242     TEC6s = 3061,
0243     TEC7r = 3070,
0244     TEC7s = 3071,
0245     TEC8r = 3080,
0246     TEC8s = 3081,
0247     TEC9r = 3090,
0248     TEC9s = 3091,
0249     TID1r = 4010,
0250     TID1s = 4011,
0251     TID2r = 4020,
0252     TID2s = 4021,
0253     TID3r = 4030,
0254     TID3s = 4031,
0255     END_OF_REGIONS
0256   };
0257 
0258   /*--------------------------------------------------------------------*/
0259   inline std::pair<int, const char*> regionType(int index)
0260   /*--------------------------------------------------------------------*/
0261   {
0262     auto region = static_cast<std::underlying_type_t<SiStripPI::TrackerRegion>>(index);
0263 
0264     switch (region) {
0265       case SiStripPI::TIB1r:
0266         return std::make_pair(1, "TIB L1 r-#varphi");
0267       case SiStripPI::TIB1s:
0268         return std::make_pair(2, "TIB L1 stereo");
0269       case SiStripPI::TIB2r:
0270         return std::make_pair(3, "TIB L2 r-#varphi");
0271       case SiStripPI::TIB2s:
0272         return std::make_pair(4, "TIB L2 stereo");
0273       case SiStripPI::TIB3r:
0274         return std::make_pair(5, "TIB L3");
0275       case SiStripPI::TIB4r:
0276         return std::make_pair(6, "TIB L4");
0277       case SiStripPI::TOB1r:
0278         return std::make_pair(7, "TOB L1 r-#varphi");
0279       case SiStripPI::TOB1s:
0280         return std::make_pair(8, "TOB L1 stereo");
0281       case SiStripPI::TOB2r:
0282         return std::make_pair(9, "TOB L2 r-#varphi");
0283       case SiStripPI::TOB2s:
0284         return std::make_pair(10, "TOB L2 stereo");
0285       case SiStripPI::TOB3r:
0286         return std::make_pair(11, "TOB L3 r-#varphi");
0287       case SiStripPI::TOB4r:
0288         return std::make_pair(12, "TOB L4");
0289       case SiStripPI::TOB5r:
0290         return std::make_pair(13, "TOB L5");
0291       case SiStripPI::TOB6r:
0292         return std::make_pair(14, "TOB L6");
0293       case SiStripPI::TEC1r:
0294         return std::make_pair(15, "TEC D1 r-#varphi");
0295       case SiStripPI::TEC1s:
0296         return std::make_pair(16, "TEC D1 stereo");
0297       case SiStripPI::TEC2r:
0298         return std::make_pair(17, "TEC D2 r-#varphi");
0299       case SiStripPI::TEC2s:
0300         return std::make_pair(18, "TEC D2 stereo");
0301       case SiStripPI::TEC3r:
0302         return std::make_pair(19, "TEC D3 r-#varphi");
0303       case SiStripPI::TEC3s:
0304         return std::make_pair(20, "TEC D3 stereo");
0305       case SiStripPI::TEC4r:
0306         return std::make_pair(21, "TEC D4 r-#varphi");
0307       case SiStripPI::TEC4s:
0308         return std::make_pair(22, "TEC D4 stereo");
0309       case SiStripPI::TEC5r:
0310         return std::make_pair(23, "TEC D5 r-#varphi");
0311       case SiStripPI::TEC5s:
0312         return std::make_pair(24, "TEC D5 stereo");
0313       case SiStripPI::TEC6r:
0314         return std::make_pair(25, "TEC D6 r-#varphi");
0315       case SiStripPI::TEC6s:
0316         return std::make_pair(26, "TEC D6 stereo");
0317       case SiStripPI::TEC7r:
0318         return std::make_pair(27, "TEC D7 r-#varphi");
0319       case SiStripPI::TEC7s:
0320         return std::make_pair(28, "TEC D7 stereo");
0321       case SiStripPI::TEC8r:
0322         return std::make_pair(29, "TEC D8 r-#varphi");
0323       case SiStripPI::TEC8s:
0324         return std::make_pair(30, "TEC D8 stereo");
0325       case SiStripPI::TEC9r:
0326         return std::make_pair(31, "TEC D9 r-#varphi");
0327       case SiStripPI::TEC9s:
0328         return std::make_pair(32, "TEC D9 stereo");
0329       case SiStripPI::TID1r:
0330         return std::make_pair(33, "TID D1 r-#varphi");
0331       case SiStripPI::TID1s:
0332         return std::make_pair(34, "TID D1 stereo");
0333       case SiStripPI::TID2r:
0334         return std::make_pair(35, "TID D2 r-#varphi");
0335       case SiStripPI::TID2s:
0336         return std::make_pair(36, "TID D2 stereo");
0337       case SiStripPI::TID3r:
0338         return std::make_pair(37, "TID D3 r-#varphi");
0339       case SiStripPI::TID3s:
0340         return std::make_pair(38, "TID D3 stereo");
0341       case SiStripPI::END_OF_REGIONS:
0342         return std::make_pair(-1, "undefined");
0343       default:
0344         return std::make_pair(999, "should never be here");
0345     }
0346   }
0347 
0348   /*--------------------------------------------------------------------*/
0349   inline std::pair<float, float> getTheRange(std::map<uint32_t, float> values, const float nsigma)
0350   /*--------------------------------------------------------------------*/
0351   {
0352     float sum = std::accumulate(
0353         std::begin(values), std::end(values), 0.0, [](float value, const std::map<uint32_t, float>::value_type& p) {
0354           return value + p.second;
0355         });
0356 
0357     float m = sum / values.size();
0358 
0359     float accum = 0.0;
0360     std::for_each(std::begin(values), std::end(values), [&](const std::map<uint32_t, float>::value_type& p) {
0361       accum += (p.second - m) * (p.second - m);
0362     });
0363 
0364     float stdev = sqrt(accum / (values.size() - 1));
0365 
0366     if (stdev != 0.) {
0367       return std::make_pair(m - nsigma * stdev, m + nsigma * stdev);
0368     } else {
0369       return std::make_pair(m > 0. ? 0.95 * m : 1.05 * m, m > 0 ? 1.05 * m : 0.95 * m);
0370     }
0371   }
0372 
0373   /*--------------------------------------------------------------------*/
0374   inline void drawStatBox(std::map<std::string, std::shared_ptr<TH1F>> histos,
0375                           std::map<std::string, int> colormap,
0376                           std::vector<std::string> legend,
0377                           double X = 0.15,
0378                           double Y = 0.93,
0379                           double W = 0.15,
0380                           double H = 0.10)
0381   /*--------------------------------------------------------------------*/
0382   {
0383     char buffer[255];
0384 
0385     int i = 0;
0386     for (const auto& element : legend) {
0387       TPaveText* stat = new TPaveText(X, Y - (i * H), X + W, Y - (i + 1) * H, "NDC");
0388       i++;
0389       auto Histo = histos[element];
0390       sprintf(buffer, "Entries : %i\n", (int)Histo->GetEntries());
0391       stat->AddText(buffer);
0392 
0393       sprintf(buffer, "Mean    : %6.2f\n", Histo->GetMean());
0394       stat->AddText(buffer);
0395 
0396       sprintf(buffer, "RMS     : %6.2f\n", Histo->GetRMS());
0397       stat->AddText(buffer);
0398 
0399       stat->SetFillColor(0);
0400       stat->SetLineColor(colormap[element]);
0401       stat->SetTextColor(colormap[element]);
0402       stat->SetTextSize(0.03);
0403       stat->SetBorderSize(0);
0404       stat->SetMargin(0.05);
0405       stat->SetTextAlign(12);
0406       stat->Draw();
0407     }
0408   }
0409 
0410   /*--------------------------------------------------------------------*/
0411   inline std::pair<float, float> getExtrema(TH1* h1, TH1* h2)
0412   /*--------------------------------------------------------------------*/
0413   {
0414     float theMax(-9999.);
0415     float theMin(9999.);
0416     theMax = h1->GetMaximum() > h2->GetMaximum() ? h1->GetMaximum() : h2->GetMaximum();
0417     theMin = h1->GetMinimum() < h2->GetMaximum() ? h1->GetMinimum() : h2->GetMinimum();
0418 
0419     float add_min = theMin > 0. ? -0.05 : 0.05;
0420     float add_max = theMax > 0. ? 0.05 : -0.05;
0421 
0422     auto result = std::make_pair(theMin * (1 + add_min), theMax * (1 + add_max));
0423     return result;
0424   }
0425 
0426   /*--------------------------------------------------------------------*/
0427   inline double getMaximum(TObjArray* array)
0428   /*--------------------------------------------------------------------*/
0429   {
0430     double theMaximum = (static_cast<TH1*>(array->At(0)))->GetMaximum();
0431     for (int i = 0; i < array->GetSize(); i++) {
0432       if ((static_cast<TH1*>(array->At(i)))->GetMaximum() > theMaximum) {
0433         theMaximum = (static_cast<TH1*>(array->At(i)))->GetMaximum();
0434       }
0435     }
0436     return theMaximum;
0437   }
0438 
0439   /*--------------------------------------------------------------------*/
0440   inline void makeNicePlotStyle(TH1* hist)
0441   /*--------------------------------------------------------------------*/
0442   {
0443     hist->SetStats(kFALSE);
0444     hist->SetLineWidth(2);
0445     hist->GetXaxis()->CenterTitle(true);
0446     hist->GetYaxis()->CenterTitle(true);
0447     hist->GetXaxis()->SetTitleFont(42);
0448     hist->GetYaxis()->SetTitleFont(42);
0449     hist->GetXaxis()->SetTitleSize(0.05);
0450     hist->GetYaxis()->SetTitleSize(0.05);
0451     hist->GetXaxis()->SetTitleOffset(0.9);
0452     hist->GetYaxis()->SetTitleOffset(1.3);
0453     hist->GetXaxis()->SetLabelFont(42);
0454     hist->GetYaxis()->SetLabelFont(42);
0455     hist->GetYaxis()->SetLabelSize(.05);
0456     hist->GetXaxis()->SetLabelSize(.05);
0457   }
0458 
0459   /*--------------------------------------------------------------------*/
0460   template <class T>
0461   inline void makeNiceStyle(T* hist)
0462   /*--------------------------------------------------------------------*/
0463   {
0464     // only for TH1s and inherited classes
0465     if constexpr (std::is_base_of<T, TH1>::value) {
0466       hist->SetStats(kFALSE);
0467       hist->SetLineWidth(2);
0468     }
0469     hist->GetXaxis()->CenterTitle(true);
0470     hist->GetYaxis()->CenterTitle(true);
0471     hist->GetXaxis()->SetTitleFont(42);
0472     hist->GetYaxis()->SetTitleFont(42);
0473     hist->GetXaxis()->SetTitleSize(0.05);
0474     hist->GetYaxis()->SetTitleSize(0.05);
0475     hist->GetXaxis()->SetTitleOffset(0.9);
0476     hist->GetYaxis()->SetTitleOffset(1.3);
0477     hist->GetXaxis()->SetLabelFont(42);
0478     hist->GetYaxis()->SetLabelFont(42);
0479     hist->GetYaxis()->SetLabelSize(.05);
0480     hist->GetXaxis()->SetLabelSize(.05);
0481   }
0482 
0483   /*--------------------------------------------------------------------*/
0484   inline void printSummary(const std::map<unsigned int, SiStripDetSummary::Values>& map)
0485   /*--------------------------------------------------------------------*/
0486   {
0487     for (const auto& element : map) {
0488       int count = element.second.count;
0489       double mean = count > 0 ? (element.second.mean) / count : 0.;
0490       double rms = count > 0 ? (element.second.rms) / count - mean * mean : 0.;
0491       if (rms <= 0)
0492         rms = 0;
0493       else
0494         rms = sqrt(rms);
0495 
0496       std::string detector;
0497 
0498       switch ((element.first) / 1000) {
0499         case 1:
0500           detector = "TIB ";
0501           break;
0502         case 2:
0503           detector = "TOB ";
0504           break;
0505         case 3:
0506           detector = "TEC ";
0507           break;
0508         case 4:
0509           detector = "TID ";
0510           break;
0511       }
0512 
0513       int layer = (element.first) / 10 - (element.first) / 1000 * 100;
0514       int stereo = (element.first) - (layer * 10) - (element.first) / 1000 * 1000;
0515 
0516       std::cout << "key of the map:" << element.first << " ( region: " << regionType(element.first).second << " ) "
0517                 << detector << " layer: " << layer << " stereo:" << stereo << "| count:" << count << " mean: " << mean
0518                 << " rms: " << rms << std::endl;
0519     }
0520   }
0521 
0522   // code is mutuated from CalibTracker/SiStripQuality/plugins/SiStripQualityStatistics
0523 
0524   /*--------------------------------------------------------------------*/
0525   inline void setBadComponents(int i, int component, const SiStripQuality::BadComponent& BC, int NBadComponent[4][19][4])
0526   /*--------------------------------------------------------------------*/
0527   {
0528     if (BC.BadApvs) {
0529       NBadComponent[i][0][2] += std::bitset<16>(BC.BadApvs & 0x3f).count();
0530       NBadComponent[i][component][2] += std::bitset<16>(BC.BadApvs & 0x3f).count();
0531     }
0532 
0533     if (BC.BadFibers) {
0534       NBadComponent[i][0][1] += std::bitset<4>(BC.BadFibers & 0x7).count();
0535       NBadComponent[i][component][1] += std::bitset<4>(BC.BadFibers & 0x7).count();
0536     }
0537 
0538     if (BC.BadModule) {
0539       NBadComponent[i][0][0]++;
0540       NBadComponent[i][component][0]++;
0541     }
0542   }
0543 
0544   // generic code to fill a SiStripDetSummary with Noise payload info
0545   /*--------------------------------------------------------------------*/
0546   inline void fillNoiseDetSummary(SiStripDetSummary& summaryNoise,
0547                                   std::shared_ptr<SiStripNoises> payload,
0548                                   SiStripPI::estimator est)
0549   /*--------------------------------------------------------------------*/
0550   {
0551     SiStripNoises::RegistryIterator rit = payload->getRegistryVectorBegin(), erit = payload->getRegistryVectorEnd();
0552     uint16_t Nstrips;
0553     std::vector<float> vstripnoise;
0554     double mean, rms, min, max;
0555     for (; rit != erit; ++rit) {
0556       Nstrips = (rit->iend - rit->ibegin) * 8 / 9;  //number of strips = number of chars * char size / strip noise size
0557       vstripnoise.resize(Nstrips);
0558       payload->allNoises(
0559           vstripnoise,
0560           make_pair(payload->getDataVectorBegin() + rit->ibegin, payload->getDataVectorBegin() + rit->iend));
0561 
0562       mean = 0;
0563       rms = 0;
0564       min = 10000;
0565       max = 0;
0566 
0567       DetId detId(rit->detid);
0568 
0569       for (size_t i = 0; i < Nstrips; ++i) {
0570         mean += vstripnoise[i];
0571         rms += vstripnoise[i] * vstripnoise[i];
0572         if (vstripnoise[i] < min)
0573           min = vstripnoise[i];
0574         if (vstripnoise[i] > max)
0575           max = vstripnoise[i];
0576       }
0577 
0578       mean /= Nstrips;
0579       if ((rms / Nstrips - mean * mean) > 0.) {
0580         rms = sqrt(rms / Nstrips - mean * mean);
0581       } else {
0582         rms = 0.;
0583       }
0584 
0585       switch (est) {
0586         case SiStripPI::min:
0587           summaryNoise.add(detId, min);
0588           break;
0589         case SiStripPI::max:
0590           summaryNoise.add(detId, max);
0591           break;
0592         case SiStripPI::mean:
0593           summaryNoise.add(detId, mean);
0594           break;
0595         case SiStripPI::rms:
0596           summaryNoise.add(detId, rms);
0597           break;
0598         default:
0599           edm::LogWarning("LogicError") << "Unknown estimator: " << est;
0600           break;
0601       }
0602     }
0603   }
0604 
0605   /*--------------------------------------------------------------------*/
0606   inline void fillTotalComponents(int NTkComponents[4], int NComponents[4][19][4], const TrackerTopology m_trackerTopo)
0607   /*--------------------------------------------------------------------*/
0608   {
0609     const auto detInfo =
0610         SiStripDetInfoFileReader::read(edm::FileInPath(SiStripDetInfoFileReader::kDefaultFile).fullPath());
0611     for (const auto& det : detInfo.getAllData()) {
0612       int nAPVs = detInfo.getNumberOfApvsAndStripLength(det.first).first;
0613       // one fiber connects to 2 APVs
0614       int nFibers = nAPVs / 2;
0615       int nStrips = (sistrip::STRIPS_PER_APV * detInfo.getNumberOfApvsAndStripLength(det.first).first);
0616       NTkComponents[0]++;
0617       NTkComponents[1] += nFibers;
0618       NTkComponents[2] += nAPVs;
0619       NTkComponents[3] += nStrips;
0620 
0621       DetId detectorId = DetId(det.first);
0622       int subDet = detectorId.subdetId();
0623 
0624       int subDetIndex = -1;
0625       int component = -1;
0626       if (subDet == StripSubdetector::TIB) {
0627         subDetIndex = 0;
0628         component = m_trackerTopo.tibLayer(det.first);
0629       } else if (subDet == StripSubdetector::TID) {
0630         subDetIndex = 1;
0631         component = m_trackerTopo.tidSide(det.first) == 2 ? m_trackerTopo.tidWheel(det.first)
0632                                                           : m_trackerTopo.tidWheel(det.first) + 3;
0633       } else if (subDet == StripSubdetector::TOB) {
0634         subDetIndex = 2;
0635         component = m_trackerTopo.tobLayer(det.first);
0636       } else if (subDet == StripSubdetector::TEC) {
0637         subDetIndex = 3;
0638         component = m_trackerTopo.tecSide(det.first) == 2 ? m_trackerTopo.tecWheel(det.first)
0639                                                           : m_trackerTopo.tecWheel(det.first) + 9;
0640       }
0641 
0642       NComponents[subDetIndex][0][0]++;
0643       NComponents[subDetIndex][0][1] += nFibers;
0644       NComponents[subDetIndex][0][2] += nAPVs;
0645       NComponents[subDetIndex][0][3] += nStrips;
0646 
0647       NComponents[subDetIndex][component][0]++;
0648       NComponents[subDetIndex][component][1] += nFibers;
0649       NComponents[subDetIndex][component][2] += nAPVs;
0650       NComponents[subDetIndex][component][3] += nStrips;
0651     }
0652   }
0653 
0654   // generic code to fill the vectors of bad components
0655   /*--------------------------------------------------------------------*/
0656   inline void fillBCArrays(const SiStripQuality* siStripQuality_,
0657                            int NTkBadComponent[4],
0658                            int NBadComponent[4][19][4],
0659                            const TrackerTopology m_trackerTopo)
0660   /*--------------------------------------------------------------------*/
0661   {
0662     std::vector<SiStripQuality::BadComponent> BC = siStripQuality_->getBadComponentList();
0663 
0664     for (size_t i = 0; i < BC.size(); ++i) {
0665       //&&&&&&&&&&&&&
0666       //Full Tk
0667       //&&&&&&&&&&&&&
0668 
0669       if (BC.at(i).BadModule)
0670         NTkBadComponent[0]++;
0671       if (BC.at(i).BadFibers)
0672         NTkBadComponent[1] +=
0673             ((BC.at(i).BadFibers >> 2) & 0x1) + ((BC.at(i).BadFibers >> 1) & 0x1) + ((BC.at(i).BadFibers) & 0x1);
0674       if (BC.at(i).BadApvs)
0675         NTkBadComponent[2] += ((BC.at(i).BadApvs >> 5) & 0x1) + ((BC.at(i).BadApvs >> 4) & 0x1) +
0676                               ((BC.at(i).BadApvs >> 3) & 0x1) + ((BC.at(i).BadApvs >> 2) & 0x1) +
0677                               ((BC.at(i).BadApvs >> 1) & 0x1) + ((BC.at(i).BadApvs) & 0x1);
0678 
0679       //&&&&&&&&&&&&&&&&&
0680       //Single SubSyste
0681       //&&&&&&&&&&&&&&&&&
0682       int component;
0683       DetId detectorId = DetId(BC.at(i).detid);
0684       int subDet = detectorId.subdetId();
0685       if (subDet == StripSubdetector::TIB) {
0686         //&&&&&&&&&&&&&&&&&
0687         //TIB
0688         //&&&&&&&&&&&&&&&&&
0689 
0690         component = m_trackerTopo.tibLayer(BC.at(i).detid);
0691         SiStripPI::setBadComponents(0, component, BC.at(i), NBadComponent);
0692 
0693       } else if (subDet == StripSubdetector::TID) {
0694         //&&&&&&&&&&&&&&&&&
0695         //TID
0696         //&&&&&&&&&&&&&&&&&
0697 
0698         component = m_trackerTopo.tidSide(BC.at(i).detid) == 2 ? m_trackerTopo.tidWheel(BC.at(i).detid)
0699                                                                : m_trackerTopo.tidWheel(BC.at(i).detid) + 3;
0700         SiStripPI::setBadComponents(1, component, BC.at(i), NBadComponent);
0701 
0702       } else if (subDet == StripSubdetector::TOB) {
0703         //&&&&&&&&&&&&&&&&&
0704         //TOB
0705         //&&&&&&&&&&&&&&&&&
0706 
0707         component = m_trackerTopo.tobLayer(BC.at(i).detid);
0708         SiStripPI::setBadComponents(2, component, BC.at(i), NBadComponent);
0709 
0710       } else if (subDet == StripSubdetector::TEC) {
0711         //&&&&&&&&&&&&&&&&&
0712         //TEC
0713         //&&&&&&&&&&&&&&&&&
0714 
0715         component = m_trackerTopo.tecSide(BC.at(i).detid) == 2 ? m_trackerTopo.tecWheel(BC.at(i).detid)
0716                                                                : m_trackerTopo.tecWheel(BC.at(i).detid) + 9;
0717         SiStripPI::setBadComponents(3, component, BC.at(i), NBadComponent);
0718       }
0719     }
0720 
0721     //&&&&&&&&&&&&&&&&&&
0722     // Single Strip Info
0723     //&&&&&&&&&&&&&&&&&&
0724 
0725     const auto detInfo =
0726         SiStripDetInfoFileReader::read(edm::FileInPath(SiStripDetInfoFileReader::kDefaultFile).fullPath());
0727 
0728     float percentage = 0;
0729 
0730     SiStripQuality::RegistryIterator rbegin = siStripQuality_->getRegistryVectorBegin();
0731     SiStripQuality::RegistryIterator rend = siStripQuality_->getRegistryVectorEnd();
0732 
0733     for (SiStripBadStrip::RegistryIterator rp = rbegin; rp != rend; ++rp) {
0734       uint32_t detid = rp->detid;
0735 
0736       int subdet = -999;
0737       int component = -999;
0738       DetId detectorId = DetId(detid);
0739       int subDet = detectorId.subdetId();
0740       if (subDet == StripSubdetector::TIB) {
0741         subdet = 0;
0742         component = m_trackerTopo.tibLayer(detid);
0743       } else if (subDet == StripSubdetector::TID) {
0744         subdet = 1;
0745         component =
0746             m_trackerTopo.tidSide(detid) == 2 ? m_trackerTopo.tidWheel(detid) : m_trackerTopo.tidWheel(detid) + 3;
0747       } else if (subDet == StripSubdetector::TOB) {
0748         subdet = 2;
0749         component = m_trackerTopo.tobLayer(detid);
0750       } else if (subDet == StripSubdetector::TEC) {
0751         subdet = 3;
0752         component =
0753             m_trackerTopo.tecSide(detid) == 2 ? m_trackerTopo.tecWheel(detid) : m_trackerTopo.tecWheel(detid) + 9;
0754       }
0755 
0756       SiStripQuality::Range sqrange = SiStripQuality::Range(siStripQuality_->getDataVectorBegin() + rp->ibegin,
0757                                                             siStripQuality_->getDataVectorBegin() + rp->iend);
0758 
0759       percentage = 0;
0760       for (int it = 0; it < sqrange.second - sqrange.first; it++) {
0761         unsigned int range = siStripQuality_->decode(*(sqrange.first + it)).range;
0762         NTkBadComponent[3] += range;
0763         NBadComponent[subdet][0][3] += range;
0764         NBadComponent[subdet][component][3] += range;
0765         percentage += range;
0766       }
0767       if (percentage != 0)
0768         percentage /= sistrip::STRIPS_PER_APV * detInfo.getNumberOfApvsAndStripLength(detid).first;
0769       if (percentage > 1)
0770         edm::LogError("SiStripBadStrip_PayloadInspector")
0771             << "PROBLEM detid " << detid << " value " << percentage << std::endl;
0772     }
0773   }
0774 
0775   /*--------------------------------------------------------------------*/
0776   inline void printBCDebug(int NTkBadComponent[4], int NBadComponent[4][19][4])
0777   /*--------------------------------------------------------------------*/
0778   {
0779     //&&&&&&&&&&&&&&&&&&
0780     // printout
0781     //&&&&&&&&&&&&&&&&&&
0782 
0783     std::stringstream ss;
0784     ss.str("");
0785     ss << "\n-----------------\nGlobal Info\n-----------------";
0786     ss << "\nBadComponent \t   Modules \tFibers "
0787           "\tApvs\tStrips\n----------------------------------------------------------------";
0788     ss << "\nTracker:\t\t" << NTkBadComponent[0] << "\t" << NTkBadComponent[1] << "\t" << NTkBadComponent[2] << "\t"
0789        << NTkBadComponent[3];
0790     ss << "\n";
0791     ss << "\nTIB:\t\t\t" << NBadComponent[0][0][0] << "\t" << NBadComponent[0][0][1] << "\t" << NBadComponent[0][0][2]
0792        << "\t" << NBadComponent[0][0][3];
0793     ss << "\nTID:\t\t\t" << NBadComponent[1][0][0] << "\t" << NBadComponent[1][0][1] << "\t" << NBadComponent[1][0][2]
0794        << "\t" << NBadComponent[1][0][3];
0795     ss << "\nTOB:\t\t\t" << NBadComponent[2][0][0] << "\t" << NBadComponent[2][0][1] << "\t" << NBadComponent[2][0][2]
0796        << "\t" << NBadComponent[2][0][3];
0797     ss << "\nTEC:\t\t\t" << NBadComponent[3][0][0] << "\t" << NBadComponent[3][0][1] << "\t" << NBadComponent[3][0][2]
0798        << "\t" << NBadComponent[3][0][3];
0799     ss << "\n";
0800 
0801     for (int i = 1; i < 5; ++i)
0802       ss << "\nTIB Layer " << i << " :\t\t" << NBadComponent[0][i][0] << "\t" << NBadComponent[0][i][1] << "\t"
0803          << NBadComponent[0][i][2] << "\t" << NBadComponent[0][i][3];
0804     ss << "\n";
0805     for (int i = 1; i < 4; ++i)
0806       ss << "\nTID+ Disk " << i << " :\t\t" << NBadComponent[1][i][0] << "\t" << NBadComponent[1][i][1] << "\t"
0807          << NBadComponent[1][i][2] << "\t" << NBadComponent[1][i][3];
0808     for (int i = 4; i < 7; ++i)
0809       ss << "\nTID- Disk " << i - 3 << " :\t\t" << NBadComponent[1][i][0] << "\t" << NBadComponent[1][i][1] << "\t"
0810          << NBadComponent[1][i][2] << "\t" << NBadComponent[1][i][3];
0811     ss << "\n";
0812     for (int i = 1; i < 7; ++i)
0813       ss << "\nTOB Layer " << i << " :\t\t" << NBadComponent[2][i][0] << "\t" << NBadComponent[2][i][1] << "\t"
0814          << NBadComponent[2][i][2] << "\t" << NBadComponent[2][i][3];
0815     ss << "\n";
0816     for (int i = 1; i < 10; ++i)
0817       ss << "\nTEC+ Disk " << i << " :\t\t" << NBadComponent[3][i][0] << "\t" << NBadComponent[3][i][1] << "\t"
0818          << NBadComponent[3][i][2] << "\t" << NBadComponent[3][i][3];
0819     for (int i = 10; i < 19; ++i)
0820       ss << "\nTEC- Disk " << i - 9 << " :\t\t" << NBadComponent[3][i][0] << "\t" << NBadComponent[3][i][1] << "\t"
0821          << NBadComponent[3][i][2] << "\t" << NBadComponent[3][i][3];
0822     ss << "\n";
0823 
0824     //edm::LogInfo("SiStripBadStrip_PayloadInspector") << ss.str() << std::endl;
0825     std::cout << ss.str() << std::endl;
0826   }
0827 
0828   enum palette { HALFGRAY, GRAY, BLUES, REDS, ANTIGRAY, FIRE, ANTIFIRE, LOGREDBLUE, BLUERED, LOGBLUERED, DEFAULT };
0829 
0830   /*--------------------------------------------------------------------*/
0831   inline void setPaletteStyle(SiStripPI::palette palette)
0832   /*--------------------------------------------------------------------*/
0833   {
0834     TStyle* palettestyle = new TStyle("palettestyle", "Style for P-TDR");
0835 
0836     const int NRGBs = 5;
0837     const int NCont = 255;
0838 
0839     switch (palette) {
0840       case HALFGRAY: {
0841         double stops[NRGBs] = {0.00, 0.34, 0.61, 0.84, 1.00};
0842         double red[NRGBs] = {1.00, 0.91, 0.80, 0.67, 1.00};
0843         double green[NRGBs] = {1.00, 0.91, 0.80, 0.67, 1.00};
0844         double blue[NRGBs] = {1.00, 0.91, 0.80, 0.67, 1.00};
0845         TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
0846       } break;
0847 
0848       case GRAY: {
0849         double stops[NRGBs] = {0.00, 0.01, 0.05, 0.09, 0.1};
0850         double red[NRGBs] = {1.00, 0.84, 0.61, 0.34, 0.00};
0851         double green[NRGBs] = {1.00, 0.84, 0.61, 0.34, 0.00};
0852         double blue[NRGBs] = {1.00, 0.84, 0.61, 0.34, 0.00};
0853         TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
0854       } break;
0855 
0856       case BLUES: {
0857         double stops[NRGBs] = {0.00, 0.34, 0.61, 0.84, 1.00};
0858         double red[NRGBs] = {1.00, 0.84, 0.61, 0.34, 0.00};
0859         double green[NRGBs] = {1.00, 0.84, 0.61, 0.34, 0.00};
0860         double blue[NRGBs] = {1.00, 1.00, 1.00, 1.00, 1.00};
0861         TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
0862 
0863       } break;
0864 
0865       case REDS: {
0866         double stops[NRGBs] = {0.00, 0.34, 0.61, 0.84, 1.00};
0867         double red[NRGBs] = {1.00, 1.00, 1.00, 1.00, 1.00};
0868         double green[NRGBs] = {1.00, 0.84, 0.61, 0.34, 0.00};
0869         double blue[NRGBs] = {1.00, 0.84, 0.61, 0.34, 0.00};
0870         TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
0871       } break;
0872 
0873       case ANTIGRAY: {
0874         double stops[NRGBs] = {0.00, 0.34, 0.61, 0.84, 1.00};
0875         double red[NRGBs] = {0.00, 0.34, 0.61, 0.84, 1.00};
0876         double green[NRGBs] = {0.00, 0.34, 0.61, 0.84, 1.00};
0877         double blue[NRGBs] = {0.00, 0.34, 0.61, 0.84, 1.00};
0878         TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
0879       } break;
0880 
0881       case FIRE: {
0882         const int NCOLs = 4;
0883         double stops[NCOLs] = {0.00, 0.20, 0.80, 1.00};
0884         double red[NCOLs] = {1.00, 1.00, 1.00, 0.50};
0885         double green[NCOLs] = {1.00, 1.00, 0.00, 0.00};
0886         double blue[NCOLs] = {0.20, 0.00, 0.00, 0.00};
0887         TColor::CreateGradientColorTable(NCOLs, stops, red, green, blue, NCont);
0888       } break;
0889 
0890       case ANTIFIRE: {
0891         const int NCOLs = 4;
0892         double stops[NCOLs] = {0.00, 0.20, 0.80, 1.00};
0893         double red[NCOLs] = {0.50, 1.00, 1.00, 1.00};
0894         double green[NCOLs] = {0.00, 0.00, 1.00, 1.00};
0895         double blue[NCOLs] = {0.00, 0.00, 0.00, 0.20};
0896         TColor::CreateGradientColorTable(NCOLs, stops, red, green, blue, NCont);
0897       } break;
0898 
0899       case LOGREDBLUE: {
0900         double stops[NRGBs] = {0.0001, 0.0010, 0.0100, 0.1000, 1.0000};
0901         double red[NRGBs] = {1.00, 0.75, 0.50, 0.25, 0.00};
0902         double green[NRGBs] = {0.00, 0.00, 0.00, 0.00, 0.00};
0903         double blue[NRGBs] = {0.00, 0.25, 0.50, 0.75, 1.00};
0904         TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
0905       } break;
0906 
0907       case LOGBLUERED: {
0908         double stops[NRGBs] = {0.0001, 0.0010, 0.0100, 0.1000, 1.0000};
0909         double red[NRGBs] = {0.00, 0.25, 0.50, 0.75, 1.00};
0910         double green[NRGBs] = {0.00, 0.00, 0.00, 0.00, 0.00};
0911         double blue[NRGBs] = {1.00, 0.75, 0.50, 0.25, 0.00};
0912         TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
0913       } break;
0914 
0915       case BLUERED: {
0916         double stops[NRGBs] = {0.00, 0.34, 0.61, 0.84, 1.00};
0917         double red[NRGBs] = {0.00, 0.25, 0.50, 0.75, 1.00};
0918         double green[NRGBs] = {0.00, 0.00, 0.00, 0.00, 0.00};
0919         double blue[NRGBs] = {1.00, 0.75, 0.50, 0.25, 0.00};
0920         TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
0921       } break;
0922 
0923       case DEFAULT: {
0924         double stops[NRGBs] = {0.00, 0.34, 0.61, 0.84, 1.00};
0925         double red[NRGBs] = {0.00, 0.00, 0.87, 1.00, 0.51};
0926         double green[NRGBs] = {0.00, 0.81, 1.00, 0.20, 0.00};
0927         double blue[NRGBs] = {0.51, 1.00, 0.12, 0.00, 0.00};
0928         TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
0929       } break;
0930       default:
0931         std::cout << "should nevere be here" << std::endl;
0932         break;
0933     }
0934 
0935     palettestyle->SetNumberContours(NCont);
0936   }
0937 
0938 };  // namespace SiStripPI
0939 #endif