CMS 3D CMS Logo

RPCDaqInfo.cc
Go to the documentation of this file.
8 
10  FEDRange_.first = ps.getUntrackedParameter<unsigned int>("MinimumRPCFEDId", 790);
11  FEDRange_.second = ps.getUntrackedParameter<unsigned int>("MaximumRPCFEDId", 792);
12 
13  NumberOfFeds_ = FEDRange_.second - FEDRange_.first + 1;
14 
15  numberOfDisks_ = ps.getUntrackedParameter<int>("NumberOfEndcapDisks", 4);
16 
17  init_ = false;
18 }
19 
23  DQMStore::IGetter& igetter,
24  edm::LuminosityBlock const& LB,
25  edm::EventSetup const& iSetup) {
26  if (!init_) {
27  this->myBooker(ibooker);
28  }
29 
30  if (auto runInfoRec = iSetup.tryToGet<RunInfoRcd>()) {
31  //get fed summary information
33  runInfoRec->get(sumFED);
34  std::vector<int> FedsInIds = sumFED->m_fed_in;
35 
36  int FedCount = 0;
37 
38  //loop on all active feds
39  for (unsigned int fedItr = 0; fedItr < FedsInIds.size(); ++fedItr) {
40  int fedID = FedsInIds[fedItr];
41  //make sure fed id is in allowed range
42 
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  int limit = numberOfDisks_;
66  if (numberOfDisks_ < 2)
67  limit = 2;
68 
69  for (int i = -1 * limit; i <= limit; i++) { //loop on wheels and disks
70  if (i > -3 && i < 3) { //wheels
71  std::stringstream streams;
72  streams << "RPC_Wheel" << i;
73  daqWheelFractions[i + 2] = ibooker.bookFloat(streams.str());
74  daqWheelFractions[i + 2]->Fill(-1);
75  }
76 
77  if (i == 0 || i > numberOfDisks_ || i < (-1 * numberOfDisks_))
78  continue;
79 
80  int offset = numberOfDisks_;
81  if (i > 0)
82  offset--; //used to skip case equale to zero
83 
84  std::stringstream streams;
85  streams << "RPC_Disk" << i;
86  daqDiskFractions[i + 2] = ibooker.bookFloat(streams.str());
87  daqDiskFractions[i + 2]->Fill(-1);
88  }
89 
90  //daq summary for RPCs
91  ibooker.setCurrentFolder("RPC/EventInfo");
92 
93  DaqFraction_ = ibooker.bookFloat("DAQSummary");
94 
95  DaqMap_ = ibooker.book2D("DAQSummaryMap", "RPC DAQ Summary Map", 15, -7.5, 7.5, 12, 0.5, 12.5);
96 
97  //customize the 2d histo
98  std::stringstream BinLabel;
99  for (int i = 1; i <= 15; i++) {
100  BinLabel.str("");
101  if (i < 13) {
102  BinLabel << "Sec" << i;
103  DaqMap_->setBinLabel(i, BinLabel.str(), 2);
104  }
105 
106  BinLabel.str("");
107  if (i < 5)
108  BinLabel << "Disk" << i - 5;
109  else if (i > 11)
110  BinLabel << "Disk" << i - 11;
111  else if (i == 11 || i == 5)
112  BinLabel.str("");
113  else
114  BinLabel << "Wheel" << i - 8;
115 
116  DaqMap_->setBinLabel(i, BinLabel.str(), 1);
117  }
118 
119  init_ = true;
120 }
FEDNumbering.h
mps_fire.i
i
Definition: mps_fire.py:428
RPCDaqInfo::FEDRange_
std::pair< int, int > FEDRange_
Definition: RPCDaqInfo.h:35
RPCDaqInfo::myBooker
void myBooker(DQMStore::IBooker &)
Definition: RPCDaqInfo.cc:61
RunSummaryRcd.h
dqm::implementation::IBooker::bookFloat
MonitorElement * bookFloat(TString const &name, FUNC onbooking=NOOP())
Definition: DQMStore.h:80
ESHandle.h
RPCDaqInfo::init_
bool init_
Definition: RPCDaqInfo.h:28
edm::LuminosityBlock
Definition: LuminosityBlock.h:50
RPCDaqInfo::daqWheelFractions
MonitorElement * daqWheelFractions[5]
Definition: RPCDaqInfo.h:32
RPCDaqInfo::DaqMap_
MonitorElement * DaqMap_
Definition: RPCDaqInfo.h:31
dqm::implementation::NavigatorBase::setCurrentFolder
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
RPCDaqInfo.h
RunInfo::m_fed_in
std::vector< int > m_fed_in
Definition: RunInfo.h:25
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
edm::EventSetup::tryToGet
std::optional< T > tryToGet() const
Definition: EventSetup.h:101
RPCDaqInfo::RPCDaqInfo
RPCDaqInfo(const edm::ParameterSet &)
Definition: RPCDaqInfo.cc:9
dqm::impl::MonitorElement::Fill
void Fill(long long x)
Definition: MonitorElement.h:290
edm::ESHandle
Definition: DTSurvey.h:22
CentralityFilter_cfi.BinLabel
BinLabel
Definition: CentralityFilter_cfi.py:5
HLT_Fake1_cff.streams
streams
Definition: HLT_Fake1_cff.py:13
RunInfoRcd
Definition: RunSummaryRcd.h:26
edm::ParameterSet
Definition: ParameterSet.h:47
RPCDaqInfo::NumberOfFeds_
int NumberOfFeds_
Definition: RPCDaqInfo.h:37
dqm::impl::MonitorElement::setBinLabel
virtual void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
Definition: MonitorElement.cc:771
edm::EventSetup
Definition: EventSetup.h:57
RPCDaqInfo::beginJob
void beginJob() override
Definition: RPCDaqInfo.cc:21
RPCDaqInfo::~RPCDaqInfo
~RPCDaqInfo() override
Definition: RPCDaqInfo.cc:20
RunSummary.h
RunInfo.h
remoteMonitoring_LED_IterMethod_cfg.limit
limit
Definition: remoteMonitoring_LED_IterMethod_cfg.py:427
RPCDaqInfo::daqDiskFractions
MonitorElement * daqDiskFractions[10]
Definition: RPCDaqInfo.h:33
RPCDaqInfo::dqmEndLuminosityBlock
void dqmEndLuminosityBlock(DQMStore::IBooker &, DQMStore::IGetter &, edm::LuminosityBlock const &, edm::EventSetup const &) override
Definition: RPCDaqInfo.cc:22
dqm::implementation::IGetter
Definition: DQMStore.h:484
dqm::implementation::IBooker::book2D
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:177
EventSetup.h
dqm::implementation::IBooker
Definition: DQMStore.h:43
RPCDaqInfo::dqmEndJob
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
Definition: RPCDaqInfo.cc:59
RPCDaqInfo::numberOfDisks_
int numberOfDisks_
Definition: RPCDaqInfo.h:37
hltrates_dqm_sourceclient-live_cfg.offset
offset
Definition: hltrates_dqm_sourceclient-live_cfg.py:82
RPCDaqInfo::DaqFraction_
MonitorElement * DaqFraction_
Definition: RPCDaqInfo.h:30