Main Page
Namespaces
Classes
Package Documentation
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Groups
Pages
DQM
RPCMonitorDigi
src
RPCRawDataCountsHistoMaker.cc
Go to the documentation of this file.
1
#include "
DQM/RPCMonitorDigi/interface/RPCRawDataCountsHistoMaker.h
"
2
#include "
EventFilter/RPCRawToDigi/interface/RPCRawDataCounts.h
"
3
#include "
EventFilter/RPCRawToDigi/interface/DataRecord.h
"
4
#include "
EventFilter/RPCRawToDigi/interface/ReadoutError.h
"
5
6
#include <vector>
7
#include <fmt/format.h>
8
9
TH1F*
RPCRawDataCountsHistoMaker::emptyReadoutErrorHisto
(
int
fedId
) {
10
const
std::string
str
=
fmt::format
(
"readoutErrors_{}"
, fedId);
11
TH1F*
result
=
new
TH1F(str.c_str(), str.c_str(), 9, 0.5, 9.5);
12
for
(
unsigned
int
i
= 1;
i
<= 9; ++
i
) {
13
rpcrawtodigi::ReadoutError::ReadoutErrorType
code =
static_cast<
rpcrawtodigi::ReadoutError::ReadoutErrorType
>
(
i
);
14
result->GetXaxis()->SetBinLabel(
i
,
rpcrawtodigi::ReadoutError::name
(code).c_str());
15
}
16
return
result
;
17
}
18
19
TH1F*
RPCRawDataCountsHistoMaker::emptyRecordTypeHisto
(
int
fedId
) {
20
const
std::string
str
=
fmt::format
(
"recordType_{}"
, fedId);
21
TH1F*
result
=
new
TH1F(str.c_str(), str.c_str(), 9, 0.5, 9.5);
22
result->SetTitleOffset(1.4,
"x"
);
23
for
(
unsigned
int
i
= 1;
i
<= 9; ++
i
) {
24
rpcrawtodigi::DataRecord::DataRecordType
code =
static_cast<
rpcrawtodigi::DataRecord::DataRecordType
>
(
i
);
25
result->GetXaxis()->SetBinLabel(
i
,
rpcrawtodigi::DataRecord::name
(code).c_str());
26
}
27
return
result
;
28
}
29
30
TH2F*
RPCRawDataCountsHistoMaker::emptyReadoutErrorMapHisto
(
int
fedId
,
int
type
) {
31
rpcrawtodigi::ReadoutError::ReadoutErrorType
code =
static_cast<
rpcrawtodigi::ReadoutError::ReadoutErrorType
>
(
type
);
32
const
std::string
str
=
fmt::format
(
"errors_{}_{}"
,
rpcrawtodigi::ReadoutError::name
(code), fedId);
33
TH2F*
result
=
new
TH2F(str.c_str(), str.c_str(), 36, -0.5, 35.5, 18, -0.5, 17.5);
34
result->GetXaxis()->SetNdivisions(512);
35
result->GetYaxis()->SetNdivisions(505);
36
result->SetXTitle(
"rmb"
);
37
result->SetYTitle(
"link"
);
38
result->SetStats(
false
);
39
return
result
;
40
}
RPCRawDataCountsHistoMaker.h
diffTreeTool.format
def format
Definition:
diffTreeTool.py:28
RPCRawDataCounts.h
mps_fire.i
i
Definition:
mps_fire.py:428
rpcrawtodigi::ReadoutError::name
std::string name() const
Definition:
ReadoutError.h:32
type
type
Definition:
SiPixelVCal_PayloadInspector.cc:39
RPCRawDataCountsHistoMaker::emptyReadoutErrorHisto
static TH1F * emptyReadoutErrorHisto(int fedId)
Definition:
RPCRawDataCountsHistoMaker.cc:9
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:256
mps_fire.result
tuple result
Definition:
mps_fire.py:311
rpcrawtodigi::DataRecord::DataRecordType
DataRecordType
Definition:
DataRecord.h:13
gainCalibHelper::gainCalibPI::type
type
Definition:
SiPixelGainCalibHelper.h:40
rpcrawtodigi::DataRecord::name
static std::string name(const DataRecordType &code)
Definition:
DataRecord.cc:58
RPCRawDataCountsHistoMaker::emptyReadoutErrorMapHisto
static TH2F * emptyReadoutErrorMapHisto(int fedId, int type)
Definition:
RPCRawDataCountsHistoMaker.cc:30
DataRecord.h
l1tstage2_dqm_sourceclient-live_cfg.fedId
tuple fedId
Definition:
l1tstage2_dqm_sourceclient-live_cfg.py:89
RPCRawDataCountsHistoMaker::emptyRecordTypeHisto
static TH1F * emptyRecordTypeHisto(int fedId)
Definition:
RPCRawDataCountsHistoMaker.cc:19
rpcrawtodigi::ReadoutError::ReadoutErrorType
ReadoutErrorType
Definition:
ReadoutError.h:10
str
#define str(s)
Definition:
TestProcessor.cc:53
ReadoutError.h
Generated for CMSSW Reference Manual by
1.8.5