CMS 3D CMS Logo

RPCDaqInfo.cc
Go to the documentation of this file.
3 
8 
9 #include <fmt/format.h>
10 
12  FEDRange_.first = ps.getUntrackedParameter<unsigned int>("MinimumRPCFEDId", 790);
13  FEDRange_.second = ps.getUntrackedParameter<unsigned int>("MaximumRPCFEDId", 792);
14 
15  NumberOfFeds_ = FEDRange_.second - FEDRange_.first + 1;
16 
17  numberOfDisks_ = ps.getUntrackedParameter<int>("NumberOfEndcapDisks", 4);
18 
19  runInfoToken_ = esConsumes<edm::Transition::EndLuminosityBlock>();
20 
21  init_ = false;
22 }
23 
26  DQMStore::IGetter& igetter,
27  edm::LuminosityBlock const& LB,
28  edm::EventSetup const& iSetup) {
29  if (!init_) {
30  this->myBooker(ibooker);
31  }
32 
33  if (auto runInfoRec = iSetup.tryToGet<RunInfoRcd>()) {
34  //get fed summary information
35  auto sumFED = runInfoRec->get(runInfoToken_);
36  const std::vector<int> FedsInIds = sumFED.m_fed_in;
37 
38  int FedCount = 0;
39 
40  //loop on all active feds
41  for (const int fedID : FedsInIds) {
42  //make sure fed id is in allowed range
43  if (fedID >= FEDRange_.first && fedID <= FEDRange_.second)
44  ++FedCount;
45  }
46 
47  //Fill active fed fraction ME
48  if (NumberOfFeds_ > 0)
49  DaqFraction_->Fill(FedCount / NumberOfFeds_);
50  else
51  DaqFraction_->Fill(-1);
52 
53  } else {
54  DaqFraction_->Fill(-1);
55  return;
56  }
57 }
58 
60 
62  //fraction of alive FEDs
63  ibooker.setCurrentFolder("RPC/EventInfo/DAQContents");
64 
65  const int limit = std::max(2, numberOfDisks_);
66 
67  for (int i = -limit; i <= limit; ++i) { //loop on wheels and disks
68  if (i > -3 && i < nWheels_ - 2) { //wheels
69  const std::string meName = fmt::format("RPC_Wheel{}", i);
70  daqWheelFractions[i + 2] = ibooker.bookFloat(meName);
71  daqWheelFractions[i + 2]->Fill(-1);
72  }
73 
74  if (i == 0 || i > numberOfDisks_ || i < -numberOfDisks_)
75  continue;
76 
77  if (i > -3 && i < nDisks_ - 2) {
78  const std::string meName = fmt::format("RPC_Disk{}", i);
79  daqDiskFractions[i + 2] = ibooker.bookFloat(meName);
80  daqDiskFractions[i + 2]->Fill(-1);
81  }
82  }
83 
84  //daq summary for RPCs
85  ibooker.setCurrentFolder("RPC/EventInfo");
86 
87  DaqFraction_ = ibooker.bookFloat("DAQSummary");
88  DaqMap_ = RPCSummaryMapHisto::book(ibooker, "DAQSummaryMap", "RPC DAQ Summary Map");
91 
92  init_ = true;
93 }
MonitorElement * bookFloat(TString const &name, FUNC onbooking=NOOP())
Definition: DQMStore.h:80
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
Definition: RPCDaqInfo.cc:59
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
bool init_
Definition: RPCDaqInfo.h:29
MonitorElement * daqWheelFractions[nWheels_]
Definition: RPCDaqInfo.h:34
std::optional< T > tryToGet() const
Definition: EventSetup.h:100
static void setBinsEndcap(MonitorElement *me, const double value)
MonitorElement * DaqFraction_
Definition: RPCDaqInfo.h:31
edm::ESGetToken< RunInfo, RunInfoRcd > runInfoToken_
Definition: RPCDaqInfo.h:27
static constexpr int nDisks_
Definition: RPCDaqInfo.h:35
void beginJob() override
Definition: RPCDaqInfo.cc:24
MonitorElement * daqDiskFractions[nDisks_]
Definition: RPCDaqInfo.h:36
T getUntrackedParameter(std::string const &, T const &) const
void Fill(long long x)
RPCDaqInfo(const edm::ParameterSet &)
Definition: RPCDaqInfo.cc:11
static MonitorElement * book(IBooker &booker, const std::string &name, const std::string &title)
std::pair< int, int > FEDRange_
Definition: RPCDaqInfo.h:38
void myBooker(DQMStore::IBooker &)
Definition: RPCDaqInfo.cc:61
void dqmEndLuminosityBlock(DQMStore::IBooker &, DQMStore::IGetter &, edm::LuminosityBlock const &, edm::EventSetup const &) override
Definition: RPCDaqInfo.cc:25
static constexpr int nWheels_
Definition: RPCDaqInfo.h:33
MonitorElement * DaqMap_
Definition: RPCDaqInfo.h:32
static void setBinsBarrel(MonitorElement *me, const double value)
int numberOfDisks_
Definition: RPCDaqInfo.h:40
int NumberOfFeds_
Definition: RPCDaqInfo.h:40