Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:08:12

0001 #include "DQM/RPCMonitorClient/interface/RPCClusterSizeTest.h"
0002 #include <DQM/RPCMonitorClient/interface/RPCRollMapHisto.h>
0003 #include "DQM/RPCMonitorClient/interface/utils.h"
0004 
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 #include "Geometry/RPCGeometry/interface/RPCGeomServ.h"
0007 
0008 RPCClusterSizeTest::RPCClusterSizeTest(const edm::ParameterSet& ps) {
0009   edm::LogVerbatim("rpceventsummary") << "[RPCClusterSizeTest]: Constructor";
0010 
0011   prescaleFactor_ = ps.getUntrackedParameter<int>("DiagnosticPrescale", 1);
0012 
0013   numberOfDisks_ = ps.getUntrackedParameter<int>("NumberOfEndcapDisks", 4);
0014   numberOfRings_ = ps.getUntrackedParameter<int>("NumberOfEndcapRings", 2);
0015   testMode_ = ps.getUntrackedParameter<bool>("testMode", false);
0016   useRollInfo_ = ps.getUntrackedParameter<bool>("useRollInfo", true);
0017 
0018   resetMEArrays();
0019 }
0020 
0021 void RPCClusterSizeTest::beginJob(std::string& workingFolder) {
0022   edm::LogVerbatim("rpceventsummary") << "[RPCClusterSizeTest]: Begin job ";
0023 
0024   globalFolder_ = workingFolder;
0025 }
0026 
0027 void RPCClusterSizeTest::getMonitorElements(std::vector<MonitorElement*>& meVector,
0028                                             std::vector<RPCDetId>& detIdVector,
0029                                             std::string& clientHistoName) {
0030   //Get  ME for each roll
0031   for (unsigned int i = 0; i < meVector.size(); i++)
0032     myClusterMe_.push_back(meVector[i]);
0033 }
0034 
0035 void RPCClusterSizeTest::clientOperation() {
0036   edm::LogVerbatim("rpceventsummary") << "[RPCClusterSizeTest]:Client Operation";
0037 
0038   //check some statements and prescale Factor
0039   if (myClusterMe_.empty())
0040     return;
0041 
0042   MonitorElement* MEAN = nullptr;   // Mean ClusterSize, Roll vs Sector
0043   MonitorElement* MEAND = nullptr;  // Mean ClusterSize, Distribution
0044 
0045   //Loop on summary histograms
0046   for (unsigned int i = 0; i < myClusterMe_.size(); ++i) {
0047     MonitorElement* myMe = myClusterMe_[i];
0048     if (!myMe || myMe->getEntries() == 0)
0049       continue;
0050     std::string meName = myMe->getName();
0051 
0052     int xBin = 0, yBin = 0;
0053     float meanCLS = 0.;
0054 
0055     //Barrel
0056     for (int wheel = -2; wheel <= 2; wheel++) {
0057       for (int sector = 1; sector <= 12; sector++) {
0058         std::string tmpName = fmt::format("ClusterSize_Wheel_{}_Sector_{}", wheel, sector);
0059 
0060         if (tmpName == meName) {
0061           MEAN = MEANWheel[wheel + 2];
0062           if (testMode_) {
0063             MEAND = MEANDWheel[wheel + 2];
0064           }
0065 
0066           xBin = sector;
0067           for (int yBin_ = 1; yBin_ <= myMe->getNbinsY(); yBin_++) {
0068             int nbinsX = myMe->getNbinsX() - 1;  //Exclude overflow
0069             TH1F* h = new TH1F("h", "h", nbinsX, 0.5, 0.5 + nbinsX);
0070             for (int xBin_ = 1; xBin_ <= nbinsX; xBin_++) {
0071               int cl = myMe->getBinContent(xBin_, yBin_);
0072               h->SetBinContent(xBin_, cl);
0073             }
0074             yBin = yBin_;
0075 
0076             meanCLS = h->GetMean();
0077 
0078             if (MEAN)
0079               MEAN->setBinContent(xBin, yBin, meanCLS);
0080 
0081             if (testMode_) {
0082               if (MEAND)
0083                 MEAND->Fill(meanCLS);
0084             }
0085           }
0086         }
0087       }
0088     }
0089 
0090     //Endcap
0091     const std::array<std::string, 4> chNames = {{"CH01-CH09", "CH10-CH18", "CH19-CH27", "CH28-CH36"}};
0092 
0093     for (int region = -1; region <= 1; region++) {
0094       if (region == 0)
0095         continue;
0096 
0097       for (int disk = 1; disk <= numberOfDisks_; disk++) {
0098         for (int ring = 2; ring < numberOfRings_ + 2; ring++) {
0099           for (unsigned int ich = 0; ich < chNames.size(); ich++) {
0100             std::string tmpName = fmt::format("ClusterSize_Disk_{}_Ring_{}_{}", (region * disk), ring, chNames[ich]);
0101 
0102             if (tmpName == meName) {
0103               xBin = (ich * 9);
0104 
0105               if (((disk * region) + numberOfDisks_) >= 0) {
0106                 if (region < 0) {
0107                   MEAN = MEANDisk[(disk * region) + numberOfDisks_];
0108                   if (testMode_) {
0109                     MEAND = MEANDDisk[(disk * region) + numberOfDisks_];
0110                   }
0111                 } else {
0112                   MEAN = MEANDisk[(disk * region) + numberOfDisks_ - 1];
0113                   if (testMode_) {
0114                     MEAND = MEANDDisk[(disk * region) + numberOfDisks_ - 1];
0115                   }
0116                 }
0117               }
0118 
0119               for (int yBin_ = 1; yBin_ <= myMe->getNbinsY(); yBin_++) {
0120                 int nbinsX = myMe->getNbinsX() - 1;  //Exclude overflow
0121                 TH1F* h = new TH1F("h", "h", nbinsX, 0.5, 0.5 + nbinsX);
0122                 for (int xBin_ = 1; xBin_ <= nbinsX; xBin_++) {
0123                   int cl = myMe->getBinContent(xBin_, yBin_);
0124                   h->SetBinContent(xBin_, cl);
0125                 }
0126 
0127                 yBin = yBin_ - (3 * ((yBin_ - 1) / 3)) + (3 * (ring - 2));
0128                 if (yBin % 3 == 1)
0129                   xBin++;
0130 
0131                 meanCLS = h->GetMean();
0132 
0133                 if (MEAN)
0134                   MEAN->setBinContent(xBin, yBin, meanCLS);
0135 
0136                 if (testMode_) {
0137                   if (MEAND)
0138                     MEAND->Fill(meanCLS);
0139                 }
0140               }
0141             }  //name check
0142           }
0143         }
0144       }
0145     }
0146   }  //End loop on chambers
0147 }
0148 
0149 void RPCClusterSizeTest::resetMEArrays(void) {
0150   memset((void*)MEANWheel, 0, sizeof(MonitorElement*) * kWheels);
0151   memset((void*)MEANDWheel, 0, sizeof(MonitorElement*) * kWheels);
0152 
0153   memset((void*)MEANDisk, 0, sizeof(MonitorElement*) * kDisks);
0154   memset((void*)MEANDDisk, 0, sizeof(MonitorElement*) * kDisks);
0155 }
0156 
0157 void RPCClusterSizeTest::myBooker(DQMStore::IBooker& ibooker) {
0158   resetMEArrays();
0159 
0160   ibooker.setCurrentFolder(globalFolder_);
0161 
0162   std::stringstream histoName;
0163 
0164   rpcdqm::utils rpcUtils;
0165 
0166   // Loop over wheels
0167   for (int w = -2; w <= 2; w++) {
0168     histoName.str("");
0169     histoName << "ClusterSizeMean_Roll_vs_Sector_Wheel" << w;  // Avarage ClusterSize (2D Roll vs Sector)
0170     auto me = RPCRollMapHisto::bookBarrel(ibooker, w, histoName.str(), histoName.str(), useRollInfo_);
0171     MEANWheel[w + 2] = dynamic_cast<MonitorElement*>(me);
0172 
0173     if (testMode_) {
0174       histoName.str("");
0175       histoName << "ClusterSizeMean_Distribution_Wheel" << w;  //  Avarage ClusterSize Distribution
0176       MEANDWheel[w + 2] = ibooker.book1D(histoName.str().c_str(), histoName.str().c_str(), 100, 1, 11);
0177     }
0178   }  //end loop on wheels
0179 
0180   for (int d = -numberOfDisks_; d <= numberOfDisks_; d++) {
0181     if (d == 0)
0182       continue;
0183     //Endcap
0184     int offset = numberOfDisks_;
0185     if (d > 0)
0186       offset--;
0187 
0188     histoName.str("");
0189     histoName << "ClusterSizeMean_Ring_vs_Segment_Disk" << d;  // Avarage ClusterSize (2D Roll vs Sector)
0190     auto me = RPCRollMapHisto::bookEndcap(ibooker, d, histoName.str(), histoName.str(), useRollInfo_);
0191     MEANDisk[d + offset] = dynamic_cast<MonitorElement*>(me);
0192 
0193     if (testMode_) {
0194       histoName.str("");
0195       histoName << "ClusterSizeMean_Distribution_Disk" << d;  //  Avarage ClusterSize Distribution
0196       MEANDDisk[d + offset] = ibooker.book1D(histoName.str().c_str(), histoName.str().c_str(), 100, 1, 11);
0197     }
0198   }
0199 }