CMS 3D CMS Logo

RPCDqmClient.cc
Go to the documentation of this file.
1 // Package: RPCDqmClient
2 // Original Author: Anna Cimmino
6 //include client headers
12 //Geometry
15 //Framework
18 
19 #include <fmt/format.h>
20 
22  edm::LogVerbatim("rpcdqmclient") << "[RPCDqmClient]: Constructor";
23 
24  offlineDQM_ = pset.getUntrackedParameter<bool>("OfflineDQM", true);
25  useRollInfo_ = pset.getUntrackedParameter<bool>("UseRollInfo", false);
26  //check enabling
27  enableDQMClients_ = pset.getUntrackedParameter<bool>("EnableRPCDqmClient", true);
28  minimumEvents_ = pset.getUntrackedParameter<int>("MinimumRPCEvents", 10000);
29  numberOfDisks_ = pset.getUntrackedParameter<int>("NumberOfEndcapDisks", 4);
30  numberOfRings_ = pset.getUntrackedParameter<int>("NumberOfEndcapRings", 2);
31 
32  std::string subsystemFolder = pset.getUntrackedParameter<std::string>("RPCFolder", "RPC");
33  std::string recHitTypeFolder = pset.getUntrackedParameter<std::string>("RecHitTypeFolder", "AllHits");
34  std::string summaryFolder = pset.getUntrackedParameter<std::string>("SummaryFolder", "SummaryHistograms");
35 
36  prefixDir_ = subsystemFolder + "/" + recHitTypeFolder;
37  globalFolder_ = subsystemFolder + "/" + recHitTypeFolder + "/" + summaryFolder;
38 
39  //get prescale factor
40  prescaleGlobalFactor_ = pset.getUntrackedParameter<int>("DiagnosticGlobalPrescale", 5);
41 
42  //make default client list
43  clientList_ = {{"RPCMultiplicityTest", "RPCDeadChannelTest", "RPCClusterSizeTest"}};
44  clientList_ = pset.getUntrackedParameter<std::vector<std::string> >("RPCDqmClientList", clientList_);
45 
46  //get all the possible RPC DQM clients
47  this->makeClientMap(pset);
48 
49  //clear counters
50  lumiCounter_ = 0;
51 
52  rpcGeomToken_ = esConsumes<edm::Transition::EndLuminosityBlock>();
53 }
54 
56  if (!enableDQMClients_) {
57  return;
58  };
59  edm::LogVerbatim("rpcdqmclient") << "[RPCDqmClient]: Begin Job";
60 
61  //Do whatever the begin jobs of all client modules do
62  for (auto& module : clientModules_) {
63  module->beginJob(globalFolder_);
64  }
65 }
66 
68  DQMStore::IGetter& igetter,
69  edm::LuminosityBlock const& lumiSeg,
70  edm::EventSetup const& c) {
71  if (!enableDQMClients_) {
72  return;
73  }
74  edm::LogVerbatim("rpcdqmclient") << "[RPCDqmClient]: End DQM LB";
75 
76  if (myDetIds_.empty()) {
77  //Get RPCdetId...
78 
79  this->getRPCdetId(c);
80 
81  //...book summary histograms
82  for (auto& module : clientModules_) {
83  module->myBooker(ibooker);
84  }
85  }
86 
87  if (!offlineDQM_) { //Do this only for the online
88 
89  if (lumiCounter_ == 0) { //only for the first lumi section do this...
90  // ...get chamber based histograms and pass them to the client modules
91  this->getMonitorElements(igetter);
92  }
93 
94  //Do not perform client oparations every lumi block
95  ++lumiCounter_;
97  return;
98  }
99 
100  //Check if there's enough statistics
101  float rpcevents = minimumEvents_;
102  if (RPCEvents_) {
103  rpcevents = RPCEvents_->getBinContent(1);
104  }
105  if (rpcevents < minimumEvents_) {
106  return;
107  }
108 
109  edm::LogVerbatim("rpcdqmclient") << "[RPCDqmClient]: Client operations";
110  for (auto& module : clientModules_) {
111  module->clientOperation();
112  }
113  } //end of online operations
114 }
115 
117  if (!enableDQMClients_) {
118  return;
119  }
120 
121  edm::LogVerbatim("rpcdqmclient") << "[RPCDqmClient]: End DQM Job";
122 
123  if (offlineDQM_) { // ...get chamber based histograms and pass them to the client modules
124  this->getMonitorElements(igetter);
125  }
126 
127  float rpcevents = minimumEvents_;
128  if (RPCEvents_) {
129  rpcevents = RPCEvents_->getBinContent(1);
130  }
131  if (rpcevents < minimumEvents_) {
132  return;
133  }
134 
135  edm::LogVerbatim("rpcdqmclient") << "[RPCDqmClient]: Client operations";
136  for (auto& module : clientModules_) {
137  module->clientOperation();
138  }
139 }
140 
142  std::vector<MonitorElement*> myMeVect;
143  std::vector<RPCDetId> myDetIds;
144 
145  //loop on all geometry and get all histos
146  for (auto& detId : myDetIds_) {
147  //Get name
148  const std::string rollName = RPCNameHelper::name(detId, useRollInfo_);
150 
151  //loop on clients
152  for (unsigned int cl = 0, nCL = clientModules_.size(); cl < nCL; ++cl) {
153  if (clientHisto_[cl] == "ClusterSize")
154  continue;
155 
156  MonitorElement* myMe = igetter.get(fmt::format("{}/{}/{}_{}", prefixDir_, folder, clientHisto_[cl], rollName));
157  if (!myMe)
158  continue;
159 
160  myMeVect.push_back(myMe);
161  myDetIds.push_back(detId);
162 
163  } //end loop on clients
164  } //end loop on all geometry and get all histos
165 
166  //Clustersize
167  std::vector<MonitorElement*> myMeVectCl;
168  const std::array<std::string, 4> chNames = {{"CH01-CH09", "CH10-CH18", "CH19-CH27", "CH28-CH36"}};
169 
170  //Retrieve barrel clustersize
171  for (int wheel = -2; wheel <= 2; wheel++) {
172  for (int sector = 1; sector <= 12; sector++) {
173  MonitorElement* myMeCl = igetter.get(fmt::format(
174  "{}/Barrel/Wheel_{}/SummaryBySectors/ClusterSize_Wheel_{}_Sector_{}", prefixDir_, wheel, wheel, sector));
175  myMeVectCl.push_back(myMeCl);
176  }
177  }
178  //Retrieve endcaps clustersize
179  for (int region = -1; region <= 1; region++) {
180  if (region == 0)
181  continue;
182 
183  std::string regionName = "Endcap-";
184  if (region == 1)
185  regionName = "Endcap+";
186 
187  for (int disk = 1; disk <= numberOfDisks_; disk++) {
188  for (int ring = numberOfRings_; ring <= 3; ring++) {
189  for (unsigned int ich = 0; ich < chNames.size(); ich++) {
190  MonitorElement* myMeCl =
191  igetter.get(fmt::format("{}/{}/Disk_{}/SummaryByRings/ClusterSize_Disk_{}_Ring_{}_{}",
192  prefixDir_,
193  regionName,
194  (region * disk),
195  (region * disk),
196  ring,
197  chNames[ich]));
198  myMeVectCl.push_back(myMeCl);
199  }
200  }
201  }
202  }
203 
204  RPCEvents_ = igetter.get(prefixDir_ + "/RPCEvents");
205  for (unsigned int cl = 0; cl < clientModules_.size(); ++cl) {
206  if (clientHisto_[cl] == "ClusterSize")
207  clientModules_[cl]->getMonitorElements(myMeVectCl, myDetIds, clientHisto_[cl]);
208  else
209  clientModules_[cl]->getMonitorElements(myMeVect, myDetIds, clientHisto_[cl]);
210  }
211 }
212 
214  myDetIds_.clear();
215 
216  auto rpcGeo = eventSetup.getHandle(rpcGeomToken_);
217 
218  for (auto& det : rpcGeo->dets()) {
219  const RPCChamber* ch = dynamic_cast<const RPCChamber*>(det);
220  if (!ch)
221  continue;
222 
223  //Loop on rolls in given chamber
224  for (auto& r : ch->rolls()) {
225  myDetIds_.push_back(r->id());
226  }
227  }
228 }
229 
231  for (unsigned int i = 0; i < clientList_.size(); i++) {
232  if (clientList_[i] == "RPCMultiplicityTest") {
233  clientHisto_.push_back("Multiplicity");
234  // clientTag_.push_back(rpcdqm::MULTIPLICITY);
235  clientModules_.emplace_back(new RPCMultiplicityTest(pset));
236  } else if (clientList_[i] == "RPCDeadChannelTest") {
237  clientHisto_.push_back("Occupancy");
238  clientModules_.emplace_back(new RPCDeadChannelTest(pset));
239  // clientTag_.push_back(rpcdqm::OCCUPANCY);
240  } else if (clientList_[i] == "RPCClusterSizeTest") {
241  clientHisto_.push_back("ClusterSize");
242  clientModules_.emplace_back(new RPCClusterSizeTest(pset));
243  // clientTag_.push_back(rpcdqm::CLUSTERSIZE);
244  } else if (clientList_[i] == "RPCOccupancyTest") {
245  clientHisto_.push_back("Occupancy");
246  clientModules_.emplace_back(new RPCOccupancyTest(pset));
247  // clientTag_.push_back(rpcdqm::OCCUPANCY);
248  } else if (clientList_[i] == "RPCNoisyStripTest") {
249  clientHisto_.push_back("Occupancy");
250  clientModules_.emplace_back(new RPCNoisyStripTest(pset));
251  //clientTag_.push_back(rpcdqm::OCCUPANCY);
252  }
253  }
254 
255  return;
256 }
void dqmEndLuminosityBlock(DQMStore::IBooker &, DQMStore::IGetter &, edm::LuminosityBlock const &, edm::EventSetup const &) override
Definition: RPCDqmClient.cc:67
Log< level::Info, true > LogVerbatim
MonitorElement * RPCEvents_
Definition: RPCDqmClient.h:46
void beginJob() override
Definition: RPCDqmClient.cc:55
std::string prefixDir_
Definition: RPCDqmClient.h:42
std::vector< std::unique_ptr< RPCClient > > clientModules_
Definition: RPCDqmClient.h:49
int numberOfRings_
Definition: RPCDqmClient.h:39
void getMonitorElements(DQMStore::IGetter &)
std::vector< std::string > clientList_
Definition: RPCDqmClient.h:44
edm::ESGetToken< RPCGeometry, MuonGeometryRecord > rpcGeomToken_
Definition: RPCDqmClient.h:53
void makeClientMap(const edm::ParameterSet &parameters_)
std::vector< std::string > clientHisto_
Definition: RPCDqmClient.h:48
RPCDeadChannelTest
RPC Client Modules #######################.
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
bool useRollInfo_
Definition: RPCDqmClient.h:41
std::string globalFolder_
Definition: RPCDqmClient.h:43
int numberOfDisks_
Definition: RPCDqmClient.h:39
static std::string name(const RPCDetId &detId, const bool useRoll)
Definition: RPCNameHelper.cc:6
void getRPCdetId(const edm::EventSetup &)
int minimumEvents_
Definition: RPCDqmClient.h:38
static std::string folderStructure(const RPCDetId &detId)
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:673
std::vector< RPCDetId > myDetIds_
Definition: RPCDqmClient.h:47
bool enableDQMClients_
Definition: RPCDqmClient.h:41
RPCDqmClient(const edm::ParameterSet &ps)
Constructor.
Definition: RPCDqmClient.cc:21
const std::vector< const RPCRoll * > & rolls() const
Return the Rolls.
Definition: RPCChamber.cc:40
int prescaleGlobalFactor_
Definition: RPCDqmClient.h:38
bool offlineDQM_
Definition: RPCDqmClient.h:37
virtual double getBinContent(int binx) const
get content of bin (1-D)