Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // Package:    RPCDqmClient
0002 // Original Author:  Anna Cimmino
0003 #include "DQM/RPCMonitorClient/interface/RPCDqmClient.h"
0004 #include "DQM/RPCMonitorClient/interface/RPCNameHelper.h"
0005 #include "DQM/RPCMonitorClient/interface/RPCBookFolderStructure.h"
0006 //include client headers
0007 #include "DQM/RPCMonitorClient/interface/RPCDeadChannelTest.h"
0008 #include "DQM/RPCMonitorClient/interface/RPCMultiplicityTest.h"
0009 #include "DQM/RPCMonitorClient/interface/RPCClusterSizeTest.h"
0010 #include "DQM/RPCMonitorClient/interface/RPCOccupancyTest.h"
0011 #include "DQM/RPCMonitorClient/interface/RPCNoisyStripTest.h"
0012 //Geometry
0013 #include "Geometry/RPCGeometry/interface/RPCGeomServ.h"
0014 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0015 //Framework
0016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0017 #include "FWCore/Framework/interface/ESHandle.h"
0018 
0019 #include <fmt/format.h>
0020 
0021 RPCDqmClient::RPCDqmClient(const edm::ParameterSet& pset) {
0022   edm::LogVerbatim("rpcdqmclient") << "[RPCDqmClient]: Constructor";
0023 
0024   offlineDQM_ = pset.getUntrackedParameter<bool>("OfflineDQM", true);
0025   useRollInfo_ = pset.getUntrackedParameter<bool>("UseRollInfo", false);
0026   //check enabling
0027   enableDQMClients_ = pset.getUntrackedParameter<bool>("EnableRPCDqmClient", true);
0028   minimumEvents_ = pset.getUntrackedParameter<int>("MinimumRPCEvents", 10000);
0029   numberOfDisks_ = pset.getUntrackedParameter<int>("NumberOfEndcapDisks", 4);
0030   numberOfRings_ = pset.getUntrackedParameter<int>("NumberOfEndcapRings", 2);
0031 
0032   std::string subsystemFolder = pset.getUntrackedParameter<std::string>("RPCFolder", "RPC");
0033   std::string recHitTypeFolder = pset.getUntrackedParameter<std::string>("RecHitTypeFolder", "AllHits");
0034   std::string summaryFolder = pset.getUntrackedParameter<std::string>("SummaryFolder", "SummaryHistograms");
0035 
0036   prefixDir_ = subsystemFolder + "/" + recHitTypeFolder;
0037   globalFolder_ = subsystemFolder + "/" + recHitTypeFolder + "/" + summaryFolder;
0038 
0039   //get prescale factor
0040   prescaleGlobalFactor_ = pset.getUntrackedParameter<int>("DiagnosticGlobalPrescale", 5);
0041 
0042   //make default client list
0043   clientList_ = {{"RPCMultiplicityTest", "RPCDeadChannelTest", "RPCClusterSizeTest"}};
0044   clientList_ = pset.getUntrackedParameter<std::vector<std::string> >("RPCDqmClientList", clientList_);
0045 
0046   //get all the possible RPC DQM clients
0047   this->makeClientMap(pset);
0048 
0049   //clear counters
0050   lumiCounter_ = 0;
0051   RPCEvents_ = nullptr;
0052 
0053   rpcGeomToken_ = esConsumes<edm::Transition::EndLuminosityBlock>();
0054 }
0055 
0056 void RPCDqmClient::beginJob() {
0057   if (!enableDQMClients_) {
0058     return;
0059   };
0060   edm::LogVerbatim("rpcdqmclient") << "[RPCDqmClient]: Begin Job";
0061 
0062   //Do whatever the begin jobs of all client modules do
0063   for (auto& module : clientModules_) {
0064     module->beginJob(globalFolder_);
0065   }
0066 }
0067 
0068 void RPCDqmClient::dqmEndLuminosityBlock(DQMStore::IBooker& ibooker,
0069                                          DQMStore::IGetter& igetter,
0070                                          edm::LuminosityBlock const& lumiSeg,
0071                                          edm::EventSetup const& c) {
0072   if (!enableDQMClients_) {
0073     return;
0074   }
0075   edm::LogVerbatim("rpcdqmclient") << "[RPCDqmClient]: End DQM LB";
0076 
0077   if (myDetIds_.empty()) {
0078     //Get RPCdetId...
0079 
0080     this->getRPCdetId(c);
0081 
0082     //...book summary histograms
0083     for (auto& module : clientModules_) {
0084       module->myBooker(ibooker);
0085     }
0086   }
0087 
0088   if (!offlineDQM_) {  //Do this only for the online
0089 
0090     if (lumiCounter_ == 0) {  //only for the first lumi section do this...
0091       // ...get chamber based histograms and pass them to the client modules
0092       this->getMonitorElements(igetter);
0093     }
0094 
0095     //Do not perform client oparations every lumi block
0096     ++lumiCounter_;
0097     if (lumiCounter_ % prescaleGlobalFactor_ != 0) {
0098       return;
0099     }
0100 
0101     //Check if there's enough statistics
0102     float rpcevents = minimumEvents_;
0103     if (RPCEvents_) {
0104       rpcevents = RPCEvents_->getBinContent(1);
0105     }
0106     if (rpcevents < minimumEvents_) {
0107       return;
0108     }
0109 
0110     edm::LogVerbatim("rpcdqmclient") << "[RPCDqmClient]: Client operations";
0111     for (auto& module : clientModules_) {
0112       module->clientOperation();
0113     }
0114   }  //end of online operations
0115 }
0116 
0117 void RPCDqmClient::dqmEndJob(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
0118   if (!enableDQMClients_) {
0119     return;
0120   }
0121 
0122   edm::LogVerbatim("rpcdqmclient") << "[RPCDqmClient]: End DQM Job";
0123 
0124   if (offlineDQM_) {  // ...get chamber based histograms and pass them to the client modules
0125     this->getMonitorElements(igetter);
0126   }
0127 
0128   float rpcevents = minimumEvents_;
0129   if (RPCEvents_) {
0130     rpcevents = RPCEvents_->getBinContent(1);
0131   }
0132   if (rpcevents < minimumEvents_) {
0133     return;
0134   }
0135 
0136   edm::LogVerbatim("rpcdqmclient") << "[RPCDqmClient]: Client operations";
0137   for (auto& module : clientModules_) {
0138     module->clientOperation();
0139   }
0140 }
0141 
0142 void RPCDqmClient::getMonitorElements(DQMStore::IGetter& igetter) {
0143   std::vector<MonitorElement*> myMeVect;
0144   std::vector<RPCDetId> myDetIds;
0145 
0146   //loop on all geometry and get all histos
0147   for (auto& detId : myDetIds_) {
0148     //Get name
0149     const std::string rollName = RPCNameHelper::name(detId, useRollInfo_);
0150     const std::string folder = RPCBookFolderStructure::folderStructure(detId);
0151 
0152     //loop on clients
0153     for (unsigned int cl = 0, nCL = clientModules_.size(); cl < nCL; ++cl) {
0154       if (clientHisto_[cl] == "ClusterSize")
0155         continue;
0156 
0157       MonitorElement* myMe = igetter.get(fmt::format("{}/{}/{}_{}", prefixDir_, folder, clientHisto_[cl], rollName));
0158       if (!myMe)
0159         continue;
0160 
0161       myMeVect.push_back(myMe);
0162       myDetIds.push_back(detId);
0163 
0164     }  //end loop on clients
0165   }    //end loop on all geometry and get all histos
0166 
0167   //Clustersize
0168   std::vector<MonitorElement*> myMeVectCl;
0169   const std::array<std::string, 4> chNames = {{"CH01-CH09", "CH10-CH18", "CH19-CH27", "CH28-CH36"}};
0170 
0171   //Retrieve barrel clustersize
0172   for (int wheel = -2; wheel <= 2; wheel++) {
0173     for (int sector = 1; sector <= 12; sector++) {
0174       MonitorElement* myMeCl = igetter.get(fmt::format(
0175           "{}/Barrel/Wheel_{}/SummaryBySectors/ClusterSize_Wheel_{}_Sector_{}", prefixDir_, wheel, wheel, sector));
0176       myMeVectCl.push_back(myMeCl);
0177     }
0178   }
0179   //Retrieve endcaps clustersize
0180   for (int region = -1; region <= 1; region++) {
0181     if (region == 0)
0182       continue;
0183 
0184     std::string regionName = "Endcap-";
0185     if (region == 1)
0186       regionName = "Endcap+";
0187 
0188     for (int disk = 1; disk <= numberOfDisks_; disk++) {
0189       for (int ring = numberOfRings_; ring <= 3; ring++) {
0190         for (unsigned int ich = 0; ich < chNames.size(); ich++) {
0191           MonitorElement* myMeCl =
0192               igetter.get(fmt::format("{}/{}/Disk_{}/SummaryByRings/ClusterSize_Disk_{}_Ring_{}_{}",
0193                                       prefixDir_,
0194                                       regionName,
0195                                       (region * disk),
0196                                       (region * disk),
0197                                       ring,
0198                                       chNames[ich]));
0199           myMeVectCl.push_back(myMeCl);
0200         }
0201       }
0202     }
0203   }
0204 
0205   RPCEvents_ = igetter.get(prefixDir_ + "/RPCEvents");
0206   for (unsigned int cl = 0; cl < clientModules_.size(); ++cl) {
0207     if (clientHisto_[cl] == "ClusterSize")
0208       clientModules_[cl]->getMonitorElements(myMeVectCl, myDetIds, clientHisto_[cl]);
0209     else
0210       clientModules_[cl]->getMonitorElements(myMeVect, myDetIds, clientHisto_[cl]);
0211   }
0212 }
0213 
0214 void RPCDqmClient::getRPCdetId(const edm::EventSetup& eventSetup) {
0215   myDetIds_.clear();
0216 
0217   auto rpcGeo = eventSetup.getHandle(rpcGeomToken_);
0218 
0219   for (auto& det : rpcGeo->dets()) {
0220     const RPCChamber* ch = dynamic_cast<const RPCChamber*>(det);
0221     if (!ch)
0222       continue;
0223 
0224     //Loop on rolls in given chamber
0225     for (auto& r : ch->rolls()) {
0226       myDetIds_.push_back(r->id());
0227     }
0228   }
0229 }
0230 
0231 void RPCDqmClient::makeClientMap(const edm::ParameterSet& pset) {
0232   for (unsigned int i = 0; i < clientList_.size(); i++) {
0233     if (clientList_[i] == "RPCMultiplicityTest") {
0234       clientHisto_.push_back("Multiplicity");
0235       // clientTag_.push_back(rpcdqm::MULTIPLICITY);
0236       clientModules_.emplace_back(new RPCMultiplicityTest(pset));
0237     } else if (clientList_[i] == "RPCDeadChannelTest") {
0238       clientHisto_.push_back("Occupancy");
0239       clientModules_.emplace_back(new RPCDeadChannelTest(pset));
0240       // clientTag_.push_back(rpcdqm::OCCUPANCY);
0241     } else if (clientList_[i] == "RPCClusterSizeTest") {
0242       clientHisto_.push_back("ClusterSize");
0243       clientModules_.emplace_back(new RPCClusterSizeTest(pset));
0244       // clientTag_.push_back(rpcdqm::CLUSTERSIZE);
0245     } else if (clientList_[i] == "RPCOccupancyTest") {
0246       clientHisto_.push_back("Occupancy");
0247       clientModules_.emplace_back(new RPCOccupancyTest(pset));
0248       // clientTag_.push_back(rpcdqm::OCCUPANCY);
0249     } else if (clientList_[i] == "RPCNoisyStripTest") {
0250       clientHisto_.push_back("Occupancy");
0251       clientModules_.emplace_back(new RPCNoisyStripTest(pset));
0252       //clientTag_.push_back(rpcdqm::OCCUPANCY);
0253     }
0254   }
0255 
0256   return;
0257 }