Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-03-13 02:31:34

0001 /*
0002  *  See header file for a description of this class.
0003  *
0004  *  \author C. Battilana S. Marcellini - INFN Bologna
0005  *
0006  *  threadsafe version (//-) oct/nov 2014 - WATWanAbdullah -ncpp-um-my
0007  *
0008  */
0009 
0010 // This class header
0011 #include "DQM/DTMonitorClient/src/DTLocalTriggerTest.h"
0012 
0013 // Framework headers
0014 #include "FWCore/Framework/interface/EventSetup.h"
0015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0016 #include "DQMServices/Core/interface/DQMStore.h"
0017 
0018 // Geometry
0019 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0020 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0021 
0022 // Root
0023 #include "TF1.h"
0024 #include "TProfile.h"
0025 
0026 //C++ headers
0027 #include <iostream>
0028 #include <sstream>
0029 
0030 using namespace edm;
0031 using namespace std;
0032 
0033 DTLocalTriggerTest::DTLocalTriggerTest(const edm::ParameterSet& ps) {
0034   setConfig(ps, "DTLocalTrigger");
0035   baseFolderTM = "DT/03-LocalTrigger-TM/";
0036   nMinEvts = ps.getUntrackedParameter<int>("nEventsCert", 5000);
0037 
0038   bookingdone = false;
0039 }
0040 
0041 DTLocalTriggerTest::~DTLocalTriggerTest() {}
0042 
0043 void DTLocalTriggerTest::Bookings(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
0044   vector<string>::const_iterator iTr = trigSources.begin();
0045   vector<string>::const_iterator trEnd = trigSources.end();
0046   vector<string>::const_iterator iHw = hwSources.begin();
0047   vector<string>::const_iterator hwEnd = hwSources.end();
0048 
0049   //Booking
0050   if (parameters.getUntrackedParameter<bool>("staticBooking", true)) {
0051     for (; iTr != trEnd; ++iTr) {
0052       trigSource = (*iTr);
0053       for (; iHw != hwEnd; ++iHw) {
0054         hwSource = (*iHw);
0055         // Loop over the TriggerUnits
0056         for (int wh = -2; wh <= 2; ++wh) {
0057           for (int sect = 1; sect <= 12; ++sect) {
0058             bookSectorHistos(ibooker, wh, sect, "BXDistribPhiIn");
0059             bookSectorHistos(ibooker, wh, sect, "QualDistribPhiIn");
0060             bookSectorHistos(ibooker, wh, sect, "BXDistribPhiOut");
0061             bookSectorHistos(ibooker, wh, sect, "QualDistribPhiOut");
0062           }
0063 
0064           bookWheelHistos(ibooker, wh, "CorrectBXPhiIn");
0065           bookWheelHistos(ibooker, wh, "ResidualBXPhiIn");
0066           bookWheelHistos(ibooker, wh, "CorrFractionPhiIn");
0067           bookWheelHistos(ibooker, wh, "2ndFractionPhiIn");
0068           bookWheelHistos(ibooker, wh, "TriggerInclusivePhiIn");
0069 
0070           bookWheelHistos(ibooker, wh, "CorrectBXPhiOut");
0071           bookWheelHistos(ibooker, wh, "ResidualBXPhiOut");
0072           bookWheelHistos(ibooker, wh, "CorrFractionPhiOut");
0073           bookWheelHistos(ibooker, wh, "2ndFractionPhiOut");
0074           bookWheelHistos(ibooker, wh, "TriggerInclusivePhiOut");
0075         }
0076       }
0077     }
0078   }
0079   // Summary test histo booking (only static)
0080   for (iTr = trigSources.begin(); iTr != trEnd; ++iTr) {
0081     trigSource = (*iTr);
0082     for (iHw = hwSources.begin(); iHw != hwSources.end(); ++iHw) {
0083       hwSource = (*iHw);
0084       // Loop over the TriggerUnits
0085       for (int wh = -2; wh <= 2; ++wh) {
0086         bookWheelHistos(ibooker, wh, "CorrFractionSummaryIn", "Summaries");
0087         bookWheelHistos(ibooker, wh, "2ndFractionSummaryIn", "Summaries");
0088         bookWheelHistos(ibooker, wh, "CorrFractionSummaryOut", "Summaries");
0089         bookWheelHistos(ibooker, wh, "2ndFractionSummaryOut", "Summaries");
0090       }
0091       bookCmsHistos(ibooker, "CorrFractionSummaryIn");
0092       bookCmsHistos(ibooker, "2ndFractionSummaryIn");
0093       bookCmsHistos(ibooker, "CorrFractionSummaryOut");
0094       bookCmsHistos(ibooker, "2ndFractionSummaryOut");
0095 
0096       if (hwSource == "TM") {
0097         bookCmsHistos(ibooker, "TrigGlbSummary", "", true);
0098       }
0099     }
0100   }
0101 
0102   bookingdone = true;
0103 }
0104 
0105 void DTLocalTriggerTest::beginRun(const edm::Run& r, const edm::EventSetup& c) {
0106   DTLocalTriggerBaseTest::beginRun(r, c);
0107 }
0108 
0109 void DTLocalTriggerTest::runClientDiagnostic(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
0110   if (!bookingdone)
0111     Bookings(ibooker, igetter);
0112 
0113   // Loop over Trig & Hw sources
0114   for (vector<string>::const_iterator iTr = trigSources.begin(); iTr != trigSources.end(); ++iTr) {
0115     trigSource = (*iTr);
0116 
0117     for (vector<string>::const_iterator iHw = hwSources.begin(); iHw != hwSources.end(); ++iHw) {
0118       hwSource = (*iHw);
0119       // Loop over the TriggerUnits
0120       for (int stat = 1; stat <= 4; ++stat) {
0121         for (int wh = -2; wh <= 2; ++wh) {
0122           for (int sect = 1; sect <= 12; ++sect) {
0123             DTChamberId chId(wh, stat, sect);
0124             int sector_id = (wh + wheelArrayShift) + (sect - 1) * 5;
0125 
0126             // IN part
0127             TH2F* BXvsQual = getHisto<TH2F>(igetter.get(getMEName("BXvsQual_In", "LocalTriggerPhiIn", chId)));
0128             TH1F* BestQual = getHisto<TH1F>(igetter.get(getMEName("BestQual_In", "LocalTriggerPhiIn", chId)));
0129             TH2F* Flag1stvsQual = getHisto<TH2F>(igetter.get(getMEName("Flag1stvsQual_In", "LocalTriggerPhiIn", chId)));
0130             if (BXvsQual && Flag1stvsQual && BestQual) {
0131               int corrSummary = 1;
0132               int secondSummary = 1;
0133               //default values for histograms
0134               double BX_OK = 51.;
0135               double BXMean = 51.;
0136               double corrFrac = 0.;
0137               double secondFrac = 0.;
0138               double besttrigs = 0.;
0139               if (BestQual->GetEntries() > 1) {
0140                 TH1D* BXHH = BXvsQual->ProjectionY("", 6, 7, "");
0141                 TH1D* Flag1st = Flag1stvsQual->ProjectionY();
0142                 int BXOK_bin = BXHH->GetEntries() >= 1 ? BXHH->GetMaximumBin() : 51;
0143                 BXMean = BXHH->GetEntries() >= 1 ? BXHH->GetMean() : 51;
0144                 BX_OK = BXvsQual->GetYaxis()->GetBinCenter(BXOK_bin);
0145                 double trigsFlag2nd = Flag1st->GetBinContent(2);
0146                 double trigs = Flag1st->GetEntries();
0147                 besttrigs = BestQual->GetEntries();
0148                 double besttrigsCorr = BestQual->Integral(5, 7, "");
0149                 delete BXHH;
0150                 delete Flag1st;
0151 
0152                 if (besttrigs != 0)
0153                   corrFrac = besttrigsCorr / besttrigs;
0154                 else
0155                   corrFrac = 1;
0156                 if (trigs != 0)
0157                   secondFrac = trigsFlag2nd / trigs;
0158                 else
0159                   secondFrac = 0;
0160 
0161                 if (corrFrac < parameters.getUntrackedParameter<double>("corrFracError", .5)) {
0162                   corrSummary = 2;
0163                 } else if (corrFrac < parameters.getUntrackedParameter<double>("corrFracWarning", .6)) {
0164                   corrSummary = 3;
0165                 } else {
0166                   corrSummary = 0;
0167                 }
0168                 if (secondFrac > parameters.getUntrackedParameter<double>("secondFracError", .2)) {
0169                   secondSummary = 2;
0170                 } else if (secondFrac > parameters.getUntrackedParameter<double>("secondFracWarning", .1)) {
0171                   secondSummary = 3;
0172                 } else {
0173                   secondSummary = 0;
0174                 }
0175 
0176                 if (secME[sector_id].find(fullName("BXDistribPhiIn")) == secME[sector_id].end()) {
0177                   bookSectorHistos(ibooker, wh, sect, "QualDistribPhiIn");
0178                   bookSectorHistos(ibooker, wh, sect, "BXDistribPhiIn");
0179                 }
0180 
0181                 TH1D* BXDistr = BXvsQual->ProjectionY();
0182                 TH1D* QualDistr = BXvsQual->ProjectionX();
0183                 std::map<std::string, MonitorElement*>* innerME = &(secME[sector_id]);
0184 
0185                 int nbinsBX = BXDistr->GetNbinsX();
0186                 int firstBinCenter = static_cast<int>(BXDistr->GetBinCenter(1));
0187                 int lastBinCenter = static_cast<int>(BXDistr->GetBinCenter(nbinsBX));
0188                 int iMin = firstBinCenter > -4 ? firstBinCenter : -4;
0189                 int iMax = lastBinCenter < 20 ? lastBinCenter : 20;
0190                 for (int ibin = iMin + 5; ibin <= iMax + 5; ++ibin) {
0191                   innerME->find(fullName("BXDistribPhiIn"))
0192                       ->second->setBinContent(ibin, stat, BXDistr->GetBinContent(ibin - 5 - firstBinCenter + 1));
0193                 }
0194                 for (int ibin = 1; ibin <= 7; ++ibin) {
0195                   innerME->find(fullName("QualDistribPhiIn"))
0196                       ->second->setBinContent(ibin, stat, QualDistr->GetBinContent(ibin));
0197                 }
0198 
0199                 delete BXDistr;
0200                 delete QualDistr;
0201               }
0202 
0203               if (whME[wh].find(fullName("CorrectBXPhiIn")) == whME[wh].end()) {
0204                 bookWheelHistos(ibooker, wh, "ResidualBXPhiIn");
0205                 bookWheelHistos(ibooker, wh, "CorrectBXPhiIn");
0206                 bookWheelHistos(ibooker, wh, "CorrFractionPhiIn");
0207                 bookWheelHistos(ibooker, wh, "2ndFractionPhiIn");
0208                 bookWheelHistos(ibooker, wh, "TriggerInclusivePhiIn");
0209               }
0210 
0211               std::map<std::string, MonitorElement*>* innerME = &(whME[wh]);
0212               innerME->find(fullName("CorrectBXPhiIn"))->second->setBinContent(sect, stat, BX_OK + 0.00001);
0213               innerME->find(fullName("ResidualBXPhiIn"))
0214                   ->second->setBinContent(sect, stat, round(25. * (BXMean - BX_OK)) + 0.00001);
0215               innerME->find(fullName("CorrFractionPhiIn"))->second->setBinContent(sect, stat, corrFrac);
0216               innerME->find(fullName("TriggerInclusivePhiIn"))->second->setBinContent(sect, stat, besttrigs);
0217               innerME->find(fullName("2ndFractionPhiIn"))->second->setBinContent(sect, stat, secondFrac);
0218 
0219               whME[wh].find(fullName("CorrFractionSummaryIn"))->second->setBinContent(sect, stat, corrSummary);
0220               whME[wh].find(fullName("2ndFractionSummaryIn"))->second->setBinContent(sect, stat, secondSummary);
0221 
0222             }  // closes BXvsQual && Flag1stvsQual && BestQual
0223 
0224             if (hwSource == "TM") {
0225               //Out part
0226 
0227               TH2F* BXvsQual = getHisto<TH2F>(igetter.get(getMEName("BXvsQual_Out", "LocalTriggerPhiOut", chId)));
0228               TH1F* BestQual = getHisto<TH1F>(igetter.get(getMEName("BestQual_Out", "LocalTriggerPhiOut", chId)));
0229               TH2F* Flag1stvsQual =
0230                   getHisto<TH2F>(igetter.get(getMEName("Flag1stvsQual_Out", "LocalTriggerPhiOut", chId)));
0231               if (BXvsQual && Flag1stvsQual && BestQual) {
0232                 int corrSummary = 1;
0233                 int secondSummary = 1;
0234                 //default values for histograms
0235                 double BX_OK = 51.;
0236                 double BXMean = 51.;
0237                 double corrFrac = 0.;
0238                 double secondFrac = 0.;
0239                 double besttrigs = 0.;
0240 
0241                 if (BestQual->GetEntries() > 1) {
0242                   TH1D* BXHH = BXvsQual->ProjectionY("", 6, 7, "");
0243                   TH1D* Flag1st = Flag1stvsQual->ProjectionY();
0244                   int BXOK_bin = BXHH->GetEntries() >= 1 ? BXHH->GetMaximumBin() : 51;
0245                   BXMean = BXHH->GetEntries() >= 1 ? BXHH->GetMean() : 51;
0246                   BX_OK = BXvsQual->GetYaxis()->GetBinCenter(BXOK_bin);
0247                   double trigsFlag2nd = Flag1st->GetBinContent(2);
0248                   double trigs = Flag1st->GetEntries();
0249                   besttrigs = BestQual->GetEntries();
0250                   double besttrigsCorr = BestQual->Integral(5, 7, "");
0251                   delete BXHH;
0252                   delete Flag1st;
0253 
0254                   corrFrac = besttrigsCorr / besttrigs;
0255                   secondFrac = trigsFlag2nd / trigs;
0256                   if (corrFrac < parameters.getUntrackedParameter<double>("corrFracError", .5)) {
0257                     corrSummary = 2;
0258                   } else if (corrFrac < parameters.getUntrackedParameter<double>("corrFracWarning", .6)) {
0259                     corrSummary = 3;
0260                   } else {
0261                     corrSummary = 0;
0262                   }
0263                   if (secondFrac > parameters.getUntrackedParameter<double>("secondFracError", .2)) {
0264                     secondSummary = 2;
0265                   } else if (secondFrac > parameters.getUntrackedParameter<double>("secondFracWarning", .1)) {
0266                     secondSummary = 3;
0267                   } else {
0268                     secondSummary = 0;
0269                   }
0270 
0271                   if (secME[sector_id].find(fullName("BXDistribPhiOut")) == secME[sector_id].end()) {
0272                     bookSectorHistos(ibooker, wh, sect, "QualDistribPhiOut");
0273                     bookSectorHistos(ibooker, wh, sect, "BXDistribPhiOut");
0274                   }
0275 
0276                   TH1D* BXDistr = BXvsQual->ProjectionY();
0277                   TH1D* QualDistr = BXvsQual->ProjectionX();
0278                   std::map<std::string, MonitorElement*>* innerME = &(secME[sector_id]);
0279 
0280                   int nbinsBX = BXDistr->GetNbinsX();
0281                   int firstBinCenter = static_cast<int>(BXDistr->GetBinCenter(1));
0282                   int lastBinCenter = static_cast<int>(BXDistr->GetBinCenter(nbinsBX));
0283                   int iMin = firstBinCenter > -4 ? firstBinCenter : -4;
0284                   int iMax = lastBinCenter < 20 ? lastBinCenter : 20;
0285                   for (int ibin = iMin + 5; ibin <= iMax + 5; ++ibin) {
0286                     innerME->find(fullName("BXDistribPhiOut"))
0287                         ->second->setBinContent(ibin, stat, BXDistr->GetBinContent(ibin - 5 - firstBinCenter + 1));
0288                   }
0289                   for (int ibin = 1; ibin <= 7; ++ibin) {
0290                     innerME->find(fullName("QualDistribPhiOut"))
0291                         ->second->setBinContent(ibin, stat, QualDistr->GetBinContent(ibin));
0292                   }
0293 
0294                   delete BXDistr;
0295                   delete QualDistr;
0296                 }
0297 
0298                 if (whME[wh].find(fullName("CorrectBXPhiOut")) == whME[wh].end()) {
0299                   bookWheelHistos(ibooker, wh, "ResidualBXPhiOut");
0300                   bookWheelHistos(ibooker, wh, "CorrectBXPhiOut");
0301                   bookWheelHistos(ibooker, wh, "CorrFractionPhiOut");
0302                   bookWheelHistos(ibooker, wh, "2ndFractionPhiOut");
0303                   bookWheelHistos(ibooker, wh, "TriggerInclusivePhiOut");
0304                 }
0305 
0306                 std::map<std::string, MonitorElement*>* innerME = &(whME[wh]);
0307                 innerME->find(fullName("CorrectBXPhiOut"))->second->setBinContent(sect, stat, BX_OK + 0.00001);
0308                 innerME->find(fullName("ResidualBXPhiOut"))
0309                     ->second->setBinContent(sect, stat, round(25. * (BXMean - BX_OK)) + 0.00001);
0310                 innerME->find(fullName("CorrFractionPhiOut"))->second->setBinContent(sect, stat, corrFrac);
0311                 innerME->find(fullName("TriggerInclusivePhiOut"))->second->setBinContent(sect, stat, besttrigs);
0312                 innerME->find(fullName("2ndFractionPhiOut"))->second->setBinContent(sect, stat, secondFrac);
0313 
0314                 whME[wh].find(fullName("CorrFractionSummaryOut"))->second->setBinContent(sect, stat, corrSummary);
0315                 whME[wh].find(fullName("2ndFractionSummaryOut"))->second->setBinContent(sect, stat, secondSummary);
0316 
0317               }  // closes BXvsQual && Flag1stvsQual && BestQual
0318 
0319             }  // Check on TM source
0320             //Theta part
0321             if (hwSource == "TM") {
0322               // Perform TM plot analysis (Theta ones)
0323               TH2F* ThetaPosvsBX = getHisto<TH2F>(igetter.get(getMEName("PositionvsBX", "LocalTriggerTheta", chId)));
0324               double BX_OK = 48;
0325               // no theta triggers in stat 4!
0326               if (ThetaPosvsBX && stat < 4 && ThetaPosvsBX->GetEntries() > 1) {
0327                 TH1D* BX = ThetaPosvsBX->ProjectionX();
0328                 int BXOK_bin = BX->GetEffectiveEntries() >= 1 ? BX->GetMaximumBin() : 10;
0329                 BX_OK = ThetaPosvsBX->GetXaxis()->GetBinCenter(BXOK_bin);
0330                 delete BX;
0331 
0332                 if (whME[wh].find(fullName("CorrectBXTheta")) == whME[wh].end()) {
0333                   bookWheelHistos(ibooker, wh, "CorrectBXTheta");
0334                 }
0335                 std::map<std::string, MonitorElement*>* innerME = &(whME.find(wh)->second);
0336                 innerME->find(fullName("CorrectBXTheta"))->second->setBinContent(sect, stat, BX_OK + 0.00001);
0337               }
0338               // Adding trigger info to compute H fraction (11/10/2016) M.C.Fouz
0339               TH2F* ThetaBXvsQual = getHisto<TH2F>(igetter.get(getMEName("ThetaBXvsQual", "LocalTriggerTheta", chId)));
0340               TH1F* ThetaBestQual = getHisto<TH1F>(igetter.get(getMEName("ThetaBestQual", "LocalTriggerTheta", chId)));
0341               if (ThetaBXvsQual && ThetaBestQual && stat < 4 && ThetaBestQual->GetEntries() > 1) {
0342                 double trigs = ThetaBestQual->GetEntries();
0343                 double trigsH = ThetaBestQual->GetBinContent(
0344                     2);  // Note that for the new plots H is at bin=2 and not 4 as in DDU!!!!
0345                 if (whME[wh].find(fullName("HFractionTheta")) == whME[wh].end()) {
0346                   bookWheelHistos(ibooker, wh, "HFractionTheta");
0347                 }
0348                 std::map<std::string, MonitorElement*>* innerME = &(whME.find(wh)->second);
0349                 innerME->find(fullName("HFractionTheta"))->second->setBinContent(sect, stat, trigsH / trigs);
0350               }
0351               // END ADDING H Fraction info
0352             }
0353           }
0354         }
0355       }
0356     }
0357   }
0358 
0359   for (vector<string>::const_iterator iTr = trigSources.begin(); iTr != trigSources.end(); ++iTr) {
0360     trigSource = (*iTr);
0361     for (vector<string>::const_iterator iHw = hwSources.begin(); iHw != hwSources.end(); ++iHw) {
0362       hwSource = (*iHw);
0363       for (int wh = -2; wh <= 2; ++wh) {
0364         std::map<std::string, MonitorElement*>* innerME = &(whME[wh]);
0365         // In part
0366         TH2F* corrWhSummaryIn = getHisto<TH2F>(innerME->find(fullName("CorrFractionSummaryIn"))->second);
0367         TH2F* secondWhSummaryIn = getHisto<TH2F>(innerME->find(fullName("2ndFractionSummaryIn"))->second);
0368         for (int sect = 1; sect <= 12; ++sect) {
0369           int corrErr = 0;
0370           int secondErr = 0;
0371           int corrNoData = 0;
0372           int secondNoData = 0;
0373           for (int stat = 1; stat <= 4; ++stat) {
0374             switch (static_cast<int>(corrWhSummaryIn->GetBinContent(sect, stat))) {
0375               case 1:
0376                 corrNoData++;
0377                 [[fallthrough]];
0378               case 2:
0379                 corrErr++;
0380             }
0381             switch (static_cast<int>(secondWhSummaryIn->GetBinContent(sect, stat))) {
0382               case 1:
0383                 secondNoData++;
0384                 [[fallthrough]];
0385               case 2:
0386                 secondErr++;
0387             }
0388           }
0389           if (corrNoData == 4)
0390             corrErr = 5;
0391           if (secondNoData == 4)
0392             secondErr = 5;
0393           cmsME.find(fullName("CorrFractionSummaryIn"))->second->setBinContent(sect, wh + wheelArrayShift, corrErr);
0394           cmsME.find(fullName("2ndFractionSummaryIn"))->second->setBinContent(sect, wh + wheelArrayShift, secondErr);
0395         }
0396         // Out part
0397         TH2F* corrWhSummaryOut = getHisto<TH2F>(innerME->find(fullName("CorrFractionSummaryOut"))->second);
0398         TH2F* secondWhSummaryOut = getHisto<TH2F>(innerME->find(fullName("2ndFractionSummaryOut"))->second);
0399         for (int sect = 1; sect <= 12; ++sect) {
0400           int corrErr = 0;
0401           int secondErr = 0;
0402           int corrNoData = 0;
0403           int secondNoData = 0;
0404           for (int stat = 1; stat <= 4; ++stat) {
0405             switch (static_cast<int>(corrWhSummaryOut->GetBinContent(sect, stat))) {
0406               case 1:
0407                 corrNoData++;
0408                 [[fallthrough]];
0409               case 2:
0410                 corrErr++;
0411             }
0412             switch (static_cast<int>(secondWhSummaryOut->GetBinContent(sect, stat))) {
0413               case 1:
0414                 secondNoData++;
0415                 [[fallthrough]];
0416               case 2:
0417                 secondErr++;
0418             }
0419           }
0420           if (corrNoData == 4)
0421             corrErr = 5;
0422           if (secondNoData == 4)
0423             secondErr = 5;
0424           cmsME.find(fullName("CorrFractionSummaryOut"))->second->setBinContent(sect, wh + wheelArrayShift, corrErr);
0425           cmsME.find(fullName("2ndFractionSummaryOut"))->second->setBinContent(sect, wh + wheelArrayShift, secondErr);
0426         }
0427       }
0428     }
0429   }
0430   fillGlobalSummary(igetter);
0431 }
0432 
0433 void DTLocalTriggerTest::fillGlobalSummary(DQMStore::IGetter& igetter) {
0434   float glbPerc[5] = {1., 0.9, 0.6, 0.3, 0.01};
0435   trigSource = "";
0436   hwSource = "TM";
0437 
0438   int nSecReadout = 0;
0439 
0440   for (int wh = -2; wh <= 2; ++wh) {
0441     for (int sect = 1; sect <= 12; ++sect) {
0442       float maxErr = 8.;
0443       int corr = cmsME.find(fullName("CorrFractionSummaryIn"))->second->getBinContent(sect, wh + wheelArrayShift);
0444       int second = cmsME.find(fullName("2ndFractionSummaryIn"))->second->getBinContent(sect, wh + wheelArrayShift);
0445       int lut = 0;
0446       MonitorElement* lutsME = igetter.get(topFolder() + "Summaries/TrigLutSummary");
0447       if (lutsME) {
0448         lut = lutsME->getBinContent(sect, wh + wheelArrayShift);
0449         maxErr += 4;
0450       } else {
0451         LogTrace(category()) << "[" << testName << "Test]: TM Lut test Summary histo not found." << endl;
0452       }
0453       (corr < 5 || second < 5) && nSecReadout++;
0454       int errcode = ((corr < 5 ? corr : 4) + (second < 5 ? second : 4) + (lut < 5 ? lut : 4));
0455       errcode = min(int((errcode / maxErr + 0.01) * 4), 4);
0456       cmsME.find("TrigGlbSummary")->second->setBinContent(sect, wh + wheelArrayShift, glbPerc[errcode]);
0457     }
0458   }
0459 
0460   if (!nSecReadout)
0461     cmsME.find("TrigGlbSummary")->second->Reset();  // white histo id TM is not RO
0462 
0463   string nEvtsName = "DT/EventInfo/Counters/nProcessedEventsTrigger";
0464   MonitorElement* meProcEvts = igetter.get(nEvtsName);
0465 
0466   if (meProcEvts) {
0467     int nProcEvts = meProcEvts->getFloatValue();
0468     cmsME.find("TrigGlbSummary")->second->setEntries(nProcEvts < nMinEvts ? 10. : nProcEvts);
0469   } else {
0470     cmsME.find("TrigGlbSummary")->second->setEntries(nMinEvts + 1);
0471     LogVerbatim(category()) << "[" << testName << "Test]: ME: " << nEvtsName << " not found!" << endl;
0472   }
0473 }