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 
53 
54 
55 void
57 {
58  dbe_ = nullptr;
60 
61  std::string commonFolder = "/EventInfo/DAQContents";
62  std::string subsystFolder;
63  std::string curentFolder;
64 
65  subsystFolder="Pixel";
66  curentFolder= subsystFolder+commonFolder;
67  dbe_->setCurrentFolder(curentFolder);
68  DaqFraction[Pixel] = dbe_->bookFloat("PixelDaqFraction");
69 
70 
71  subsystFolder="SiStrip";
72  curentFolder=subsystFolder+commonFolder;
73  dbe_->setCurrentFolder(curentFolder);
74  DaqFraction[SiStrip] = dbe_->bookFloat("SiStripDaqFraction");
75 
76  subsystFolder="RPC";
77  curentFolder=subsystFolder+commonFolder;
78  dbe_->setCurrentFolder(curentFolder);
79  DaqFraction[RPC] = dbe_->bookFloat("RPCDaqFraction");
80 
81  subsystFolder="CSC";
82  curentFolder=subsystFolder+commonFolder;
83  dbe_->setCurrentFolder(curentFolder);
84  DaqFraction[CSC] = dbe_->bookFloat("CSCDaqFraction");
85 
86  subsystFolder="DT";
87  curentFolder=subsystFolder+commonFolder;
88  dbe_->setCurrentFolder(curentFolder);
89  DaqFraction[DT] = dbe_->bookFloat("DTDaqFraction");
90 
91  subsystFolder="Hcal";
92  curentFolder=subsystFolder+commonFolder;
93  dbe_->setCurrentFolder(curentFolder);
94  DaqFraction[Hcal] = dbe_->bookFloat("HcalDaqFraction");
95 
96  subsystFolder="EcalBarrel";
97  curentFolder=subsystFolder+commonFolder;
98  dbe_->setCurrentFolder(curentFolder);
99  DaqFraction[EcalBarrel] = dbe_->bookFloat("EcalBarrDaqFraction");
100 
101  subsystFolder="EcalEndcap";
102  curentFolder=subsystFolder+commonFolder;
103  dbe_->setCurrentFolder(curentFolder);
104  DaqFraction[EcalEndcap] = dbe_->bookFloat("EcalEndDaqFraction");
105 
106  subsystFolder="L1T";
107  curentFolder=subsystFolder+commonFolder;
108  dbe_->setCurrentFolder(curentFolder);
109  DaqFraction[L1T] = dbe_->bookFloat("L1TDaqFraction");
110 
111 
118  RPCRange.first = 790;
119  RPCRange.second = 792;
120  DTRange.first = 770;
121  DTRange.second = 774;
126  ECALBarrRange.first = 610;
127  ECALBarrRange.second = 645;
128  ECALEndcapRangeLow.first = 601;
129  ECALEndcapRangeLow.second = 609;
130  ECALEndcapRangeHigh.first = 646;
131  ECALEndcapRangeHigh.second = 654;
132 
133  NumberOfFeds[Pixel] = PixelRange.second-PixelRange.first +1;
134  NumberOfFeds[SiStrip] = TrackerRange.second-TrackerRange.first +1;
135  NumberOfFeds[CSC] = CSCRange.second-CSCRange.first +1;
136  NumberOfFeds[RPC] = RPCRange.second-RPCRange.first +1;
137  NumberOfFeds[DT] = DTRange.second-DTRange.first +1;
138  NumberOfFeds[Hcal] = HcalRange.second-HcalRange.first +1;
141  NumberOfFeds[L1T] = L1TRange.second-L1TRange.first +1;
142 
143 }
144 
145 
146 void
148 }
149 
150 
151 
152 void
154 {
155 
156 
157 }
158 
std::pair< int, int > ECALBarrRange
Definition: DQMDaqInfo.h:72
std::pair< int, int > TrackerRange
Definition: DQMDaqInfo.h:67
void beginJob() override
Definition: DQMDaqInfo.cc:56
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:77
MonitorElement * DaqFraction[9]
Definition: DQMDaqInfo.h:64
std::pair< int, int > PixelRange
Definition: DQMDaqInfo.h:66
void endJob() override
Definition: DQMDaqInfo.cc:147
MonitorElement * bookFloat(const char *name)
Book float.
Definition: DQMStore.cc:774
std::pair< int, int > DTRange
Definition: DQMDaqInfo.h:70
void Fill(long long x)
std::pair< int, int > RPCRange
Definition: DQMDaqInfo.h:69
DQMDaqInfo(const edm::ParameterSet &)
Definition: DQMDaqInfo.cc:4
DQMStore * dbe_
Definition: DQMDaqInfo.h:60
int iEvent
Definition: GenABIO.cc:230
std::pair< int, int > CSCRange
Definition: DQMDaqInfo.h:68
std::pair< int, int > ECALEndcapRangeHigh
Definition: DQMDaqInfo.h:74
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:75
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: DQMDaqInfo.cc:153
T get() const
Definition: EventSetup.h:63
boost::optional< eventsetup::EventSetupRecordGeneric > find(const eventsetup::EventSetupRecordKey &) const
Definition: EventSetup.cc:88
~DQMDaqInfo() override
std::pair< int, int > HcalRange
Definition: DQMDaqInfo.h:71
std::pair< int, int > ECALEndcapRangeLow
Definition: DQMDaqInfo.h:73
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:545