CMS 3D CMS Logo

DQMDaqInfo.cc
Go to the documentation of this file.
3 
5 {
6 }
7 
8 DQMDaqInfo::~DQMDaqInfo() = default;
9 
11 
13 
14  if( nullptr != iSetup.find( recordKey ) ) {
15 
17  iSetup.get<RunInfoRcd>().get(sumFED);
18 
19  //const RunInfo* summaryFED=sumFED.product();
20 
21  std::vector<int> FedsInIds= sumFED->m_fed_in;
22 
23  float FedCount[9]={0., 0., 0., 0., 0., 0., 0., 0., 0.};
24 
25  for(int fedID : FedsInIds) {
26  if(fedID>=PixelRange.first && fedID<=PixelRange.second) ++FedCount[Pixel];
27  if(fedID>=TrackerRange.first && fedID<=TrackerRange.second) ++FedCount[SiStrip];
28  if(fedID>=CSCRange.first && fedID<=CSCRange.second) ++FedCount[CSC];
29  if(fedID>=RPCRange.first && fedID<=RPCRange.second) ++FedCount[RPC];
30  if(fedID>=DTRange.first && fedID<=DTRange.second) ++FedCount[DT];
31  if(fedID>=HcalRange.first && fedID<=HcalRange.second) ++FedCount[Hcal];
32  if(fedID>=ECALBarrRange.first && fedID<=ECALBarrRange.second) ++FedCount[EcalBarrel];
33  if((fedID>=ECALEndcapRangeLow.first && fedID<=ECALEndcapRangeLow.second) ||
34  (fedID>=ECALEndcapRangeHigh.first && fedID<=ECALEndcapRangeHigh.second)) ++FedCount[EcalEndcap];
35  if(fedID>=L1TRange.first && fedID<=L1TRange.second) ++FedCount[L1T];
36 
37  }
38 
39  for(int detIndex=0; detIndex<9; ++detIndex) {
40  DaqFraction[detIndex]->Fill( FedCount[detIndex]/NumberOfFeds[detIndex]);
41  }
42 
43  }else{
44 
45  for(auto & detIndex : DaqFraction) detIndex->Fill(-1);
46  return;
47  }
48 
49 
50 }
51 
52 
54 }
55 
56 
57 void
59 {
60  dbe_ = nullptr;
62 
63  std::string commonFolder = "/EventInfo/DAQContents";
64  std::string subsystFolder;
65  std::string curentFolder;
66 
67  subsystFolder="Pixel";
68  curentFolder= subsystFolder+commonFolder;
69  dbe_->setCurrentFolder(curentFolder);
70  DaqFraction[Pixel] = dbe_->bookFloat("PixelDaqFraction");
71 
72 
73  subsystFolder="SiStrip";
74  curentFolder=subsystFolder+commonFolder;
75  dbe_->setCurrentFolder(curentFolder);
76  DaqFraction[SiStrip] = dbe_->bookFloat("SiStripDaqFraction");
77 
78  subsystFolder="RPC";
79  curentFolder=subsystFolder+commonFolder;
80  dbe_->setCurrentFolder(curentFolder);
81  DaqFraction[RPC] = dbe_->bookFloat("RPCDaqFraction");
82 
83  subsystFolder="CSC";
84  curentFolder=subsystFolder+commonFolder;
85  dbe_->setCurrentFolder(curentFolder);
86  DaqFraction[CSC] = dbe_->bookFloat("CSCDaqFraction");
87 
88  subsystFolder="DT";
89  curentFolder=subsystFolder+commonFolder;
90  dbe_->setCurrentFolder(curentFolder);
91  DaqFraction[DT] = dbe_->bookFloat("DTDaqFraction");
92 
93  subsystFolder="Hcal";
94  curentFolder=subsystFolder+commonFolder;
95  dbe_->setCurrentFolder(curentFolder);
96  DaqFraction[Hcal] = dbe_->bookFloat("HcalDaqFraction");
97 
98  subsystFolder="EcalBarrel";
99  curentFolder=subsystFolder+commonFolder;
100  dbe_->setCurrentFolder(curentFolder);
101  DaqFraction[EcalBarrel] = dbe_->bookFloat("EcalBarrDaqFraction");
102 
103  subsystFolder="EcalEndcap";
104  curentFolder=subsystFolder+commonFolder;
105  dbe_->setCurrentFolder(curentFolder);
106  DaqFraction[EcalEndcap] = dbe_->bookFloat("EcalEndDaqFraction");
107 
108  subsystFolder="L1T";
109  curentFolder=subsystFolder+commonFolder;
110  dbe_->setCurrentFolder(curentFolder);
111  DaqFraction[L1T] = dbe_->bookFloat("L1TDaqFraction");
112 
113 
120  RPCRange.first = 790;
121  RPCRange.second = 792;
122  DTRange.first = 770;
123  DTRange.second = 774;
128  ECALBarrRange.first = 610;
129  ECALBarrRange.second = 645;
130  ECALEndcapRangeLow.first = 601;
131  ECALEndcapRangeLow.second = 609;
132  ECALEndcapRangeHigh.first = 646;
133  ECALEndcapRangeHigh.second = 654;
134 
135  NumberOfFeds[Pixel] = PixelRange.second-PixelRange.first +1;
136  NumberOfFeds[SiStrip] = TrackerRange.second-TrackerRange.first +1;
137  NumberOfFeds[CSC] = CSCRange.second-CSCRange.first +1;
138  NumberOfFeds[RPC] = RPCRange.second-RPCRange.first +1;
139  NumberOfFeds[DT] = DTRange.second-DTRange.first +1;
140  NumberOfFeds[Hcal] = HcalRange.second-HcalRange.first +1;
143  NumberOfFeds[L1T] = L1TRange.second-L1TRange.first +1;
144 
145 }
146 
147 
148 void
150 }
151 
152 
153 
154 void
156 {
157 
158 
159 }
160 
std::pair< int, int > ECALBarrRange
Definition: DQMDaqInfo.h:73
std::pair< int, int > TrackerRange
Definition: DQMDaqInfo.h:68
void beginJob() override
Definition: DQMDaqInfo.cc:58
static HCTypeTag findType(char const *iTypeName)
find a type based on the types name, if not found will return default HCTypeTag
float NumberOfFeds[9]
Definition: DQMDaqInfo.h:78
MonitorElement * DaqFraction[9]
Definition: DQMDaqInfo.h:65
std::pair< int, int > PixelRange
Definition: DQMDaqInfo.h:67
void endJob() override
Definition: DQMDaqInfo.cc:149
MonitorElement * bookFloat(const char *name)
Book float.
Definition: DQMStore.cc:977
std::pair< int, int > DTRange
Definition: DQMDaqInfo.h:71
const eventsetup::EventSetupRecord * find(const eventsetup::EventSetupRecordKey &) const
Definition: EventSetup.cc:85
void Fill(long long x)
std::pair< int, int > RPCRange
Definition: DQMDaqInfo.h:70
DQMDaqInfo(const edm::ParameterSet &)
Definition: DQMDaqInfo.cc:4
DQMStore * dbe_
Definition: DQMDaqInfo.h:61
int iEvent
Definition: GenABIO.cc:230
void endLuminosityBlock(const edm::LuminosityBlock &, const edm::EventSetup &) override
Definition: DQMDaqInfo.cc:53
std::pair< int, int > CSCRange
Definition: DQMDaqInfo.h:69
std::pair< int, int > ECALEndcapRangeHigh
Definition: DQMDaqInfo.h:75
std::vector< int > m_fed_in
Definition: RunInfo.h:26
void beginLuminosityBlock(const edm::LuminosityBlock &, const edm::EventSetup &) override
Definition: DQMDaqInfo.cc:10
std::pair< int, int > L1TRange
Definition: DQMDaqInfo.h:76
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: DQMDaqInfo.cc:155
const T & get() const
Definition: EventSetup.h:58
~DQMDaqInfo() override
std::pair< int, int > HcalRange
Definition: DQMDaqInfo.h:72
std::pair< int, int > ECALEndcapRangeLow
Definition: DQMDaqInfo.h:74
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:746