CMS 3D CMS Logo

ESTrendTask.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <fstream>
3 #include <vector>
4 #include <cmath>
5 
17 
18 using namespace cms;
19 using namespace edm;
20 using namespace std;
21 
23  prefixME_ = ps.getUntrackedParameter<string>("prefixME", "EcalPreshower");
24  rechittoken_ = consumes<ESRecHitCollection>(ps.getParameter<InputTag>("RecHitLabel"));
25  dccCollections_ = consumes<ESRawDataCollection>(ps.getParameter<InputTag>("ESDCCCollections"));
26 
27  for (int i = 0; i < 2; ++i)
28  for (int j = 0; j < 2; ++j) {
29  hESRecHitTrend_[i][j] = nullptr;
30  hESRecHitTrendHr_[i][j] = nullptr;
31  }
32 
33  hESSLinkErrTrend_ = nullptr;
34  hESFiberErrTrend_ = nullptr;
35  hESSLinkErrTrendHr_ = nullptr;
36  hESFiberErrTrendHr_ = nullptr;
37 
38  start_time_ = 0;
39  current_time_ = 0;
40  last_time_ = 0;
41 
42  ievt_ = 0;
43 }
44 
45 void ESTrendTask::dqmBeginRun(const Run& r, const EventSetup& c) {
46  start_time_ = (r.beginTime()).unixTime();
47 
48  //std::cout << "start time : " << start_time_ << std::endl;
49 }
50 
52  char histo[200];
53 
54  iBooker.setCurrentFolder(prefixME_ + "/ESTrendTask");
55 
56  for (int i = 0; i < 2; ++i)
57  for (int j = 0; j < 2; ++j) {
58  int iz = (i == 0) ? 1 : -1;
59  sprintf(histo, "ES Trending RH Occ per 5 mins Z %d P %d", iz, j + 1);
60  hESRecHitTrend_[i][j] = iBooker.bookProfile(histo, histo, 36, 0.0, 180.0, 100, 0.0, 1.0e6, "s");
61  hESRecHitTrend_[i][j]->setAxisTitle("Elapse time (Minutes)", 1);
62  hESRecHitTrend_[i][j]->setAxisTitle("ES RecHit Occupancy / 5 minutes", 2);
63 
64  sprintf(histo, "ES Trending RH Occ per hour Z %d P %d", iz, j + 1);
65  hESRecHitTrendHr_[i][j] = iBooker.bookProfile(histo, histo, 24, 0.0, 24.0, 100, 0.0, 1.0e6, "s");
66  hESRecHitTrendHr_[i][j]->setAxisTitle("Elapse time (Hours)", 1);
67  hESRecHitTrendHr_[i][j]->setAxisTitle("ES RecHit Occupancy / hour", 2);
68  }
69 
70  sprintf(histo, "ES Trending SLink CRC Error per 5 mins");
71  hESSLinkErrTrend_ = iBooker.bookProfile(histo, histo, 36, 0.0, 180.0, 100, 0.0, 1.0e6, "s");
72  hESSLinkErrTrend_->setAxisTitle("Elapse time (Minutes)", 1);
73  hESSLinkErrTrend_->setAxisTitle("ES SLink CRC Err / 5 minutes", 2);
74 
75  sprintf(histo, "ES Trending Fiber Error per 5 mins");
76  hESFiberErrTrend_ = iBooker.bookProfile(histo, histo, 36, 0.0, 180.0, 100, 0.0, 1.0e6, "s");
77  hESFiberErrTrend_->setAxisTitle("Elapse time (Minutes)", 1);
78  hESFiberErrTrend_->setAxisTitle("ES Fiber Err / 5 minutes", 2);
79 
80  sprintf(histo, "ES Trending SLink CRC Error per hour");
81  hESSLinkErrTrendHr_ = iBooker.bookProfile(histo, histo, 24, 0.0, 24.0, 100, 0.0, 1.0e6, "s");
82  hESSLinkErrTrendHr_->setAxisTitle("Elapse time (Hours)", 1);
83  hESSLinkErrTrendHr_->setAxisTitle("ES SLink CRC Err / hour", 2);
84 
85  sprintf(histo, "ES Trending Fiber Error per hour");
86  hESFiberErrTrendHr_ = iBooker.bookProfile(histo, histo, 24, 0.0, 24.0, 100, 0.0, 1.0e6, "s");
87  hESFiberErrTrendHr_->setAxisTitle("Elapse time (Hours)", 1);
88  hESFiberErrTrendHr_->setAxisTitle("ES Fiber Err / hour", 2);
89 }
90 
91 void ESTrendTask::analyze(const Event& e, const EventSetup& c) {
92  ievt_++;
93 
94  // Collect time information
95  updateTime(e);
96 
97  long int diff_current_start = current_time_ - start_time_;
98  long int diff_last_start = last_time_ - start_time_;
99  //LogInfo("ESTrendTask") << "time difference is negative in " << ievt_ << " events\n"
100  //<< "\tcurrent - start time = " << diff_current_start
101  //<< ", \tlast - start time = " << diff_last_start << endl;
102 
103  // std::cout << "current_time : " << current_time_ << ", diff : " << diff_current_start << std::endl;
104 
105  // Calculate time interval and bin width
106  // int minuteBinWidth = int(nBasicClusterMinutely_->getTProfile()->GetXaxis()->GetBinWidth(1));
107  int minuteBinWidth = 5;
108  long int minuteBinDiff = diff_current_start / 60 / minuteBinWidth - diff_last_start / 60 / minuteBinWidth;
109  long int minuteDiff = (current_time_ - last_time_) / 60;
110 
111  // int hourBinWidth = int(nBasicClusterHourly_->getTProfile()->GetXaxis()->GetBinWidth(1));
112  int hourBinWidth = 1;
113  long int hourBinDiff = diff_current_start / 3600 / hourBinWidth - diff_last_start / 3600 / hourBinWidth;
114  long int hourDiff = (current_time_ - last_time_) / 3600;
115 
116  if (minuteDiff >= minuteBinWidth) {
117  while (minuteDiff >= minuteBinWidth)
118  minuteDiff -= minuteBinWidth;
119  }
120  if (hourDiff >= hourBinWidth) {
121  while (hourDiff >= hourBinWidth)
122  hourDiff -= hourBinWidth;
123  }
124 
125  // ES DCC
126  Int_t slinkCRCErr = 0;
127  Int_t fiberErr = 0;
128  vector<int> fiberStatus;
130  if (e.getByToken(dccCollections_, dccs)) {
131  for (ESRawDataCollection::const_iterator dccItr = dccs->begin(); dccItr != dccs->end(); ++dccItr) {
132  const ESDCCHeaderBlock& dcc = (*dccItr);
133 
134  if (dcc.getDCCErrors() == 101)
135  slinkCRCErr++;
136 
137  fiberStatus = dcc.getFEChannelStatus();
138  for (unsigned int i = 0; i < fiberStatus.size(); ++i) {
139  if (fiberStatus[i] == 4 || fiberStatus[i] == 8 || fiberStatus[i] == 10 || fiberStatus[i] == 11 ||
140  fiberStatus[i] == 12)
141  fiberErr++;
142  }
143  }
144  }
145 
146  shift2Right(hESSLinkErrTrend_->getTProfile(), minuteBinDiff);
147  hESSLinkErrTrend_->Fill(minuteDiff, slinkCRCErr);
148 
149  shift2Right(hESFiberErrTrend_->getTProfile(), minuteBinDiff);
150  hESFiberErrTrend_->Fill(minuteDiff, fiberErr);
151 
152  shift2Right(hESSLinkErrTrendHr_->getTProfile(), hourBinDiff);
153  hESSLinkErrTrendHr_->Fill(hourDiff, slinkCRCErr);
154 
155  shift2Right(hESFiberErrTrendHr_->getTProfile(), hourBinDiff);
156  hESFiberErrTrendHr_->Fill(hourDiff, fiberErr);
157 
158  // ES RecHits
159  int zside, plane;
160  int nrh[2][2];
161  for (int i = 0; i < 2; i++)
162  for (int j = 0; j < 2; j++) {
163  nrh[i][j] = 0;
164  }
165 
167  if (e.getByToken(rechittoken_, ESRecHit)) {
168  for (ESRecHitCollection::const_iterator hitItr = ESRecHit->begin(); hitItr != ESRecHit->end(); ++hitItr) {
169  ESDetId id = ESDetId(hitItr->id());
170 
171  zside = id.zside();
172  plane = id.plane();
173 
174  int i = (zside == 1) ? 0 : 1;
175  int j = plane - 1;
176 
177  nrh[i][j]++;
178  }
179  } else {
180  LogWarning("ESTrendTask") << "RecHitCollection not available";
181  }
182 
183  for (int i = 0; i < 2; ++i)
184  for (int j = 0; j < 2; ++j) {
185  shift2Right(hESRecHitTrend_[i][j]->getTProfile(), minuteBinDiff);
186  hESRecHitTrend_[i][j]->Fill(minuteDiff, nrh[i][j] / (1072 * 32.));
187 
188  shift2Right(hESRecHitTrendHr_[i][j]->getTProfile(), hourBinDiff);
189  hESRecHitTrendHr_[i][j]->Fill(hourDiff, nrh[i][j] / (1072 * 32.));
190  }
191 }
192 
194  last_time_ = current_time_;
195  current_time_ = e.time().unixTime();
196 }
197 
198 void ESTrendTask::shift2Right(TProfile* p, int bins) {
199  if (bins <= 0)
200  return;
201 
202  if (!p->GetSumw2())
203  p->Sumw2();
204  int nBins = p->GetXaxis()->GetNbins();
205 
206  // by shifting n bin to the right, the number of entries are
207  // reduced by the number in n bins including the overflow bin.
208  double nentries = p->GetEntries();
209  for (int i = 0; i < bins; i++)
210  nentries -= p->GetBinEntries(nBins + 1 - bins);
211  p->SetEntries(nentries);
212 
213  // the last bin goes to overflow
214  // each bin moves to the right
215 
216  TArrayD* sumw2 = p->GetSumw2();
217 
218  for (int i = nBins + 1; i > bins; i--) {
219  // GetBinContent return binContent/binEntries
220  p->SetBinContent(i, p->GetBinContent(i - bins) * p->GetBinEntries(i - bins));
221  p->SetBinEntries(i, p->GetBinEntries(i - bins));
222  sumw2->SetAt(sumw2->GetAt(i - bins), i);
223  }
224 }
225 
226 void ESTrendTask::shift2Left(TProfile* p, int bins) {
227  if (bins <= 0)
228  return;
229 
230  if (!p->GetSumw2())
231  p->Sumw2();
232  int nBins = p->GetXaxis()->GetNbins();
233 
234  // by shifting n bin to the left, the number of entries are
235  // reduced by the number in n bins including the underflow bin.
236  double nentries = p->GetEntries();
237  for (int i = 0; i < bins; i++)
238  nentries -= p->GetBinEntries(i);
239  p->SetEntries(nentries);
240 
241  // the first bin goes to underflow
242  // each bin moves to the right
243 
244  TArrayD* sumw2 = p->GetSumw2();
245 
246  for (int i = 0; i <= nBins + 1 - bins; i++) {
247  // GetBinContent return binContent/binEntries
248  p->SetBinContent(i, p->GetBinContent(i + bins) * p->GetBinEntries(i + bins));
249  p->SetBinEntries(i, p->GetBinEntries(i + bins));
250  sumw2->SetAt(sumw2->GetAt(i + bins), i);
251  }
252 }
253 
void dqmBeginRun(const edm::Run &r, const edm::EventSetup &c) override
Definition: ESTrendTask.cc:45
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
std::vector< T >::const_iterator const_iterator
void shift2Right(TProfile *p, int bins=1)
Definition: ESTrendTask.cc:198
int zside(DetId const &)
T getUntrackedParameter(std::string const &, T const &) const
void shift2Left(TProfile *p, int bins=1)
Definition: ESTrendTask.cc:226
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: ESTrendTask.cc:51
MonitorElement * bookProfile(TString const &name, TString const &title, int nchX, double lowX, double highX, int, double lowY, double highY, char const *option="s", FUNC onbooking=NOOP())
Definition: DQMStore.h:399
int getDCCErrors() const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const_iterator begin() const
Namespace of DDCMS conversion namespace.
const_iterator end() const
ESTrendTask(const edm::ParameterSet &ps)
Definition: ESTrendTask.cc:22
HLT enums.
const std::vector< int > & getFEChannelStatus() const
Log< level::Warning, false > LogWarning
void analyze(const edm::Event &e, const edm::EventSetup &c) override
Definition: ESTrendTask.cc:91
void updateTime(const edm::Event &)
Definition: ESTrendTask.cc:193
Definition: Run.h:45
virtual void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)