Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-11 04:32:27

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               std::map<std::string, MonitorElement*>* innerME = &(secME[sector_id]);
0204 
0205               if (whME[wh].find(fullName("CorrectBXPhiIn")) == whME[wh].end()) {
0206                 bookWheelHistos(ibooker, wh, "ResidualBXPhiIn");
0207                 bookWheelHistos(ibooker, wh, "CorrectBXPhiIn");
0208                 bookWheelHistos(ibooker, wh, "CorrFractionPhiIn");
0209                 bookWheelHistos(ibooker, wh, "2ndFractionPhiIn");
0210                 bookWheelHistos(ibooker, wh, "TriggerInclusivePhiIn");
0211               }
0212 
0213               innerME = &(whME[wh]);
0214               innerME->find(fullName("CorrectBXPhiIn"))->second->setBinContent(sect, stat, BX_OK + 0.00001);
0215               innerME->find(fullName("ResidualBXPhiIn"))
0216                   ->second->setBinContent(sect, stat, round(25. * (BXMean - BX_OK)) + 0.00001);
0217               innerME->find(fullName("CorrFractionPhiIn"))->second->setBinContent(sect, stat, corrFrac);
0218               innerME->find(fullName("TriggerInclusivePhiIn"))->second->setBinContent(sect, stat, besttrigs);
0219               innerME->find(fullName("2ndFractionPhiIn"))->second->setBinContent(sect, stat, secondFrac);
0220 
0221               whME[wh].find(fullName("CorrFractionSummaryIn"))->second->setBinContent(sect, stat, corrSummary);
0222               whME[wh].find(fullName("2ndFractionSummaryIn"))->second->setBinContent(sect, stat, secondSummary);
0223 
0224             }  // closes BXvsQual && Flag1stvsQual && BestQual
0225 
0226             if (hwSource == "TM") {
0227               //Out part
0228 
0229               TH2F* BXvsQual = getHisto<TH2F>(igetter.get(getMEName("BXvsQual_Out", "LocalTriggerPhiOut", chId)));
0230               TH1F* BestQual = getHisto<TH1F>(igetter.get(getMEName("BestQual_Out", "LocalTriggerPhiOut", chId)));
0231               TH2F* Flag1stvsQual =
0232                   getHisto<TH2F>(igetter.get(getMEName("Flag1stvsQual_Out", "LocalTriggerPhiOut", chId)));
0233               if (BXvsQual && Flag1stvsQual && BestQual) {
0234                 int corrSummary = 1;
0235                 int secondSummary = 1;
0236                 //default values for histograms
0237                 double BX_OK = 51.;
0238                 double BXMean = 51.;
0239                 double corrFrac = 0.;
0240                 double secondFrac = 0.;
0241                 double besttrigs = 0.;
0242 
0243                 if (BestQual->GetEntries() > 1) {
0244                   TH1D* BXHH = BXvsQual->ProjectionY("", 6, 7, "");
0245                   TH1D* Flag1st = Flag1stvsQual->ProjectionY();
0246                   int BXOK_bin = BXHH->GetEntries() >= 1 ? BXHH->GetMaximumBin() : 51;
0247                   BXMean = BXHH->GetEntries() >= 1 ? BXHH->GetMean() : 51;
0248                   BX_OK = BXvsQual->GetYaxis()->GetBinCenter(BXOK_bin);
0249                   double trigsFlag2nd = Flag1st->GetBinContent(2);
0250                   double trigs = Flag1st->GetEntries();
0251                   besttrigs = BestQual->GetEntries();
0252                   double besttrigsCorr = BestQual->Integral(5, 7, "");
0253                   delete BXHH;
0254                   delete Flag1st;
0255 
0256                   corrFrac = besttrigsCorr / besttrigs;
0257                   secondFrac = trigsFlag2nd / trigs;
0258                   if (corrFrac < parameters.getUntrackedParameter<double>("corrFracError", .5)) {
0259                     corrSummary = 2;
0260                   } else if (corrFrac < parameters.getUntrackedParameter<double>("corrFracWarning", .6)) {
0261                     corrSummary = 3;
0262                   } else {
0263                     corrSummary = 0;
0264                   }
0265                   if (secondFrac > parameters.getUntrackedParameter<double>("secondFracError", .2)) {
0266                     secondSummary = 2;
0267                   } else if (secondFrac > parameters.getUntrackedParameter<double>("secondFracWarning", .1)) {
0268                     secondSummary = 3;
0269                   } else {
0270                     secondSummary = 0;
0271                   }
0272 
0273                   if (secME[sector_id].find(fullName("BXDistribPhiOut")) == secME[sector_id].end()) {
0274                     bookSectorHistos(ibooker, wh, sect, "QualDistribPhiOut");
0275                     bookSectorHistos(ibooker, wh, sect, "BXDistribPhiOut");
0276                   }
0277 
0278                   TH1D* BXDistr = BXvsQual->ProjectionY();
0279                   TH1D* QualDistr = BXvsQual->ProjectionX();
0280                   std::map<std::string, MonitorElement*>* innerME = &(secME[sector_id]);
0281 
0282                   int nbinsBX = BXDistr->GetNbinsX();
0283                   int firstBinCenter = static_cast<int>(BXDistr->GetBinCenter(1));
0284                   int lastBinCenter = static_cast<int>(BXDistr->GetBinCenter(nbinsBX));
0285                   int iMin = firstBinCenter > -4 ? firstBinCenter : -4;
0286                   int iMax = lastBinCenter < 20 ? lastBinCenter : 20;
0287                   for (int ibin = iMin + 5; ibin <= iMax + 5; ++ibin) {
0288                     innerME->find(fullName("BXDistribPhiOut"))
0289                         ->second->setBinContent(ibin, stat, BXDistr->GetBinContent(ibin - 5 - firstBinCenter + 1));
0290                   }
0291                   for (int ibin = 1; ibin <= 7; ++ibin) {
0292                     innerME->find(fullName("QualDistribPhiOut"))
0293                         ->second->setBinContent(ibin, stat, QualDistr->GetBinContent(ibin));
0294                   }
0295 
0296                   delete BXDistr;
0297                   delete QualDistr;
0298                 }
0299 
0300                 std::map<std::string, MonitorElement*>* innerME = &(secME[sector_id]);
0301 
0302                 if (whME[wh].find(fullName("CorrectBXPhiOut")) == whME[wh].end()) {
0303                   bookWheelHistos(ibooker, wh, "ResidualBXPhiOut");
0304                   bookWheelHistos(ibooker, wh, "CorrectBXPhiOut");
0305                   bookWheelHistos(ibooker, wh, "CorrFractionPhiOut");
0306                   bookWheelHistos(ibooker, wh, "2ndFractionPhiOut");
0307                   bookWheelHistos(ibooker, wh, "TriggerInclusivePhiOut");
0308                 }
0309 
0310                 innerME = &(whME[wh]);
0311                 innerME->find(fullName("CorrectBXPhiOut"))->second->setBinContent(sect, stat, BX_OK + 0.00001);
0312                 innerME->find(fullName("ResidualBXPhiOut"))
0313                     ->second->setBinContent(sect, stat, round(25. * (BXMean - BX_OK)) + 0.00001);
0314                 innerME->find(fullName("CorrFractionPhiOut"))->second->setBinContent(sect, stat, corrFrac);
0315                 innerME->find(fullName("TriggerInclusivePhiOut"))->second->setBinContent(sect, stat, besttrigs);
0316                 innerME->find(fullName("2ndFractionPhiOut"))->second->setBinContent(sect, stat, secondFrac);
0317 
0318                 whME[wh].find(fullName("CorrFractionSummaryOut"))->second->setBinContent(sect, stat, corrSummary);
0319                 whME[wh].find(fullName("2ndFractionSummaryOut"))->second->setBinContent(sect, stat, secondSummary);
0320 
0321               }  // closes BXvsQual && Flag1stvsQual && BestQual
0322 
0323             }  // Check on TM source
0324             //Theta part
0325             if (hwSource == "TM") {
0326               // Perform TM plot analysis (Theta ones)
0327               TH2F* ThetaPosvsBX = getHisto<TH2F>(igetter.get(getMEName("PositionvsBX", "LocalTriggerTheta", chId)));
0328               double BX_OK = 48;
0329               // no theta triggers in stat 4!
0330               if (ThetaPosvsBX && stat < 4 && ThetaPosvsBX->GetEntries() > 1) {
0331                 TH1D* BX = ThetaPosvsBX->ProjectionX();
0332                 int BXOK_bin = BX->GetEffectiveEntries() >= 1 ? BX->GetMaximumBin() : 10;
0333                 BX_OK = ThetaPosvsBX->GetXaxis()->GetBinCenter(BXOK_bin);
0334                 delete BX;
0335 
0336                 if (whME[wh].find(fullName("CorrectBXTheta")) == whME[wh].end()) {
0337                   bookWheelHistos(ibooker, wh, "CorrectBXTheta");
0338                 }
0339                 std::map<std::string, MonitorElement*>* innerME = &(whME.find(wh)->second);
0340                 innerME->find(fullName("CorrectBXTheta"))->second->setBinContent(sect, stat, BX_OK + 0.00001);
0341               }
0342               // Adding trigger info to compute H fraction (11/10/2016) M.C.Fouz
0343               TH2F* ThetaBXvsQual = getHisto<TH2F>(igetter.get(getMEName("ThetaBXvsQual", "LocalTriggerTheta", chId)));
0344               TH1F* ThetaBestQual = getHisto<TH1F>(igetter.get(getMEName("ThetaBestQual", "LocalTriggerTheta", chId)));
0345               if (ThetaBXvsQual && ThetaBestQual && stat < 4 && ThetaBestQual->GetEntries() > 1) {
0346                 double trigs = ThetaBestQual->GetEntries();
0347                 double trigsH = ThetaBestQual->GetBinContent(
0348                     2);  // Note that for the new plots H is at bin=2 and not 4 as in DDU!!!!
0349                 if (whME[wh].find(fullName("HFractionTheta")) == whME[wh].end()) {
0350                   bookWheelHistos(ibooker, wh, "HFractionTheta");
0351                 }
0352                 std::map<std::string, MonitorElement*>* innerME = &(whME.find(wh)->second);
0353                 innerME->find(fullName("HFractionTheta"))->second->setBinContent(sect, stat, trigsH / trigs);
0354               }
0355               // END ADDING H Fraction info
0356             }
0357           }
0358         }
0359       }
0360     }
0361   }
0362 
0363   for (vector<string>::const_iterator iTr = trigSources.begin(); iTr != trigSources.end(); ++iTr) {
0364     trigSource = (*iTr);
0365     for (vector<string>::const_iterator iHw = hwSources.begin(); iHw != hwSources.end(); ++iHw) {
0366       hwSource = (*iHw);
0367       for (int wh = -2; wh <= 2; ++wh) {
0368         std::map<std::string, MonitorElement*>* innerME = &(whME[wh]);
0369         // In part
0370         TH2F* corrWhSummaryIn = getHisto<TH2F>(innerME->find(fullName("CorrFractionSummaryIn"))->second);
0371         TH2F* secondWhSummaryIn = getHisto<TH2F>(innerME->find(fullName("2ndFractionSummaryIn"))->second);
0372         for (int sect = 1; sect <= 12; ++sect) {
0373           int corrErr = 0;
0374           int secondErr = 0;
0375           int corrNoData = 0;
0376           int secondNoData = 0;
0377           for (int stat = 1; stat <= 4; ++stat) {
0378             switch (static_cast<int>(corrWhSummaryIn->GetBinContent(sect, stat))) {
0379               case 1:
0380                 corrNoData++;
0381                 [[fallthrough]];
0382               case 2:
0383                 corrErr++;
0384             }
0385             switch (static_cast<int>(secondWhSummaryIn->GetBinContent(sect, stat))) {
0386               case 1:
0387                 secondNoData++;
0388                 [[fallthrough]];
0389               case 2:
0390                 secondErr++;
0391             }
0392           }
0393           if (corrNoData == 4)
0394             corrErr = 5;
0395           if (secondNoData == 4)
0396             secondErr = 5;
0397           cmsME.find(fullName("CorrFractionSummaryIn"))->second->setBinContent(sect, wh + wheelArrayShift, corrErr);
0398           cmsME.find(fullName("2ndFractionSummaryIn"))->second->setBinContent(sect, wh + wheelArrayShift, secondErr);
0399         }
0400         // Out part
0401         TH2F* corrWhSummaryOut = getHisto<TH2F>(innerME->find(fullName("CorrFractionSummaryOut"))->second);
0402         TH2F* secondWhSummaryOut = getHisto<TH2F>(innerME->find(fullName("2ndFractionSummaryOut"))->second);
0403         for (int sect = 1; sect <= 12; ++sect) {
0404           int corrErr = 0;
0405           int secondErr = 0;
0406           int corrNoData = 0;
0407           int secondNoData = 0;
0408           for (int stat = 1; stat <= 4; ++stat) {
0409             switch (static_cast<int>(corrWhSummaryOut->GetBinContent(sect, stat))) {
0410               case 1:
0411                 corrNoData++;
0412                 [[fallthrough]];
0413               case 2:
0414                 corrErr++;
0415             }
0416             switch (static_cast<int>(secondWhSummaryOut->GetBinContent(sect, stat))) {
0417               case 1:
0418                 secondNoData++;
0419                 [[fallthrough]];
0420               case 2:
0421                 secondErr++;
0422             }
0423           }
0424           if (corrNoData == 4)
0425             corrErr = 5;
0426           if (secondNoData == 4)
0427             secondErr = 5;
0428           cmsME.find(fullName("CorrFractionSummaryOut"))->second->setBinContent(sect, wh + wheelArrayShift, corrErr);
0429           cmsME.find(fullName("2ndFractionSummaryOut"))->second->setBinContent(sect, wh + wheelArrayShift, secondErr);
0430         }
0431       }
0432     }
0433   }
0434   fillGlobalSummary(igetter);
0435 }
0436 
0437 void DTLocalTriggerTest::fillGlobalSummary(DQMStore::IGetter& igetter) {
0438   float glbPerc[5] = {1., 0.9, 0.6, 0.3, 0.01};
0439   trigSource = "";
0440   hwSource = "TM";
0441 
0442   int nSecReadout = 0;
0443 
0444   for (int wh = -2; wh <= 2; ++wh) {
0445     for (int sect = 1; sect <= 12; ++sect) {
0446       float maxErr = 8.;
0447       int corr = cmsME.find(fullName("CorrFractionSummaryIn"))->second->getBinContent(sect, wh + wheelArrayShift);
0448       int second = cmsME.find(fullName("2ndFractionSummaryIn"))->second->getBinContent(sect, wh + wheelArrayShift);
0449       int lut = 0;
0450       MonitorElement* lutsME = igetter.get(topFolder() + "Summaries/TrigLutSummary");
0451       if (lutsME) {
0452         lut = lutsME->getBinContent(sect, wh + wheelArrayShift);
0453         maxErr += 4;
0454       } else {
0455         LogTrace(category()) << "[" << testName << "Test]: TM Lut test Summary histo not found." << endl;
0456       }
0457       (corr < 5 || second < 5) && nSecReadout++;
0458       int errcode = ((corr < 5 ? corr : 4) + (second < 5 ? second : 4) + (lut < 5 ? lut : 4));
0459       errcode = min(int((errcode / maxErr + 0.01) * 4), 4);
0460       cmsME.find("TrigGlbSummary")->second->setBinContent(sect, wh + wheelArrayShift, glbPerc[errcode]);
0461     }
0462   }
0463 
0464   if (!nSecReadout)
0465     cmsME.find("TrigGlbSummary")->second->Reset();  // white histo id TM is not RO
0466 
0467   string nEvtsName = "DT/EventInfo/Counters/nProcessedEventsTrigger";
0468   MonitorElement* meProcEvts = igetter.get(nEvtsName);
0469 
0470   if (meProcEvts) {
0471     int nProcEvts = meProcEvts->getFloatValue();
0472     cmsME.find("TrigGlbSummary")->second->setEntries(nProcEvts < nMinEvts ? 10. : nProcEvts);
0473   } else {
0474     cmsME.find("TrigGlbSummary")->second->setEntries(nMinEvts + 1);
0475     LogVerbatim(category()) << "[" << testName << "Test]: ME: " << nEvtsName << " not found!" << endl;
0476   }
0477 }