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::endJob(void) { LogInfo("ESTrendTask") << "analyzed " << ievt_ << " events"; }
92 
93 void ESTrendTask::analyze(const Event& e, const EventSetup& c) {
94  ievt_++;
95 
96  // Collect time information
97  updateTime(e);
98 
99  long int diff_current_start = current_time_ - start_time_;
100  long int diff_last_start = last_time_ - start_time_;
101  //LogInfo("ESTrendTask") << "time difference is negative in " << ievt_ << " events\n"
102  //<< "\tcurrent - start time = " << diff_current_start
103  //<< ", \tlast - start time = " << diff_last_start << endl;
104 
105  // std::cout << "current_time : " << current_time_ << ", diff : " << diff_current_start << std::endl;
106 
107  // Calculate time interval and bin width
108  // int minuteBinWidth = int(nBasicClusterMinutely_->getTProfile()->GetXaxis()->GetBinWidth(1));
109  int minuteBinWidth = 5;
110  long int minuteBinDiff = diff_current_start / 60 / minuteBinWidth - diff_last_start / 60 / minuteBinWidth;
111  long int minuteDiff = (current_time_ - last_time_) / 60;
112 
113  // int hourBinWidth = int(nBasicClusterHourly_->getTProfile()->GetXaxis()->GetBinWidth(1));
114  int hourBinWidth = 1;
115  long int hourBinDiff = diff_current_start / 3600 / hourBinWidth - diff_last_start / 3600 / hourBinWidth;
116  long int hourDiff = (current_time_ - last_time_) / 3600;
117 
118  if (minuteDiff >= minuteBinWidth) {
119  while (minuteDiff >= minuteBinWidth)
120  minuteDiff -= minuteBinWidth;
121  }
122  if (hourDiff >= hourBinWidth) {
123  while (hourDiff >= hourBinWidth)
124  hourDiff -= hourBinWidth;
125  }
126 
127  // ES DCC
128  Int_t slinkCRCErr = 0;
129  Int_t fiberErr = 0;
130  vector<int> fiberStatus;
132  if (e.getByToken(dccCollections_, dccs)) {
133  for (ESRawDataCollection::const_iterator dccItr = dccs->begin(); dccItr != dccs->end(); ++dccItr) {
134  ESDCCHeaderBlock dcc = (*dccItr);
135 
136  if (dcc.getDCCErrors() == 101)
137  slinkCRCErr++;
138 
139  fiberStatus = dcc.getFEChannelStatus();
140  for (unsigned int i = 0; i < fiberStatus.size(); ++i) {
141  if (fiberStatus[i] == 4 || fiberStatus[i] == 8 || fiberStatus[i] == 10 || fiberStatus[i] == 11 ||
142  fiberStatus[i] == 12)
143  fiberErr++;
144  }
145  }
146  }
147 
148  shift2Right(hESSLinkErrTrend_->getTProfile(), minuteBinDiff);
149  hESSLinkErrTrend_->Fill(minuteDiff, slinkCRCErr);
150 
151  shift2Right(hESFiberErrTrend_->getTProfile(), minuteBinDiff);
152  hESFiberErrTrend_->Fill(minuteDiff, fiberErr);
153 
154  shift2Right(hESSLinkErrTrendHr_->getTProfile(), hourBinDiff);
155  hESSLinkErrTrendHr_->Fill(hourDiff, slinkCRCErr);
156 
157  shift2Right(hESFiberErrTrendHr_->getTProfile(), hourBinDiff);
158  hESFiberErrTrendHr_->Fill(hourDiff, fiberErr);
159 
160  // ES RecHits
161  int zside, plane;
162  int nrh[2][2];
163  for (int i = 0; i < 2; i++)
164  for (int j = 0; j < 2; j++) {
165  nrh[i][j] = 0;
166  }
167 
169  if (e.getByToken(rechittoken_, ESRecHit)) {
170  for (ESRecHitCollection::const_iterator hitItr = ESRecHit->begin(); hitItr != ESRecHit->end(); ++hitItr) {
171  ESDetId id = ESDetId(hitItr->id());
172 
173  zside = id.zside();
174  plane = id.plane();
175 
176  int i = (zside == 1) ? 0 : 1;
177  int j = plane - 1;
178 
179  nrh[i][j]++;
180  }
181  } else {
182  LogWarning("ESTrendTask") << "RecHitCollection not available";
183  }
184 
185  for (int i = 0; i < 2; ++i)
186  for (int j = 0; j < 2; ++j) {
187  shift2Right(hESRecHitTrend_[i][j]->getTProfile(), minuteBinDiff);
188  hESRecHitTrend_[i][j]->Fill(minuteDiff, nrh[i][j] / (1072 * 32.));
189 
190  shift2Right(hESRecHitTrendHr_[i][j]->getTProfile(), hourBinDiff);
191  hESRecHitTrendHr_[i][j]->Fill(hourDiff, nrh[i][j] / (1072 * 32.));
192  }
193 }
194 
196  last_time_ = current_time_;
197  current_time_ = e.time().unixTime();
198 }
199 
200 void ESTrendTask::shift2Right(TProfile* p, int bins) {
201  if (bins <= 0)
202  return;
203 
204  if (!p->GetSumw2())
205  p->Sumw2();
206  int nBins = p->GetXaxis()->GetNbins();
207 
208  // by shifting n bin to the right, the number of entries are
209  // reduced by the number in n bins including the overflow bin.
210  double nentries = p->GetEntries();
211  for (int i = 0; i < bins; i++)
212  nentries -= p->GetBinEntries(nBins + 1 - bins);
213  p->SetEntries(nentries);
214 
215  // the last bin goes to overflow
216  // each bin moves to the right
217 
218  TArrayD* sumw2 = p->GetSumw2();
219 
220  for (int i = nBins + 1; i > bins; i--) {
221  // GetBinContent return binContent/binEntries
222  p->SetBinContent(i, p->GetBinContent(i - bins) * p->GetBinEntries(i - bins));
223  p->SetBinEntries(i, p->GetBinEntries(i - bins));
224  sumw2->SetAt(sumw2->GetAt(i - bins), i);
225  }
226 }
227 
228 void ESTrendTask::shift2Left(TProfile* p, int bins) {
229  if (bins <= 0)
230  return;
231 
232  if (!p->GetSumw2())
233  p->Sumw2();
234  int nBins = p->GetXaxis()->GetNbins();
235 
236  // by shifting n bin to the left, the number of entries are
237  // reduced by the number in n bins including the underflow bin.
238  double nentries = p->GetEntries();
239  for (int i = 0; i < bins; i++)
240  nentries -= p->GetBinEntries(i);
241  p->SetEntries(nentries);
242 
243  // the first bin goes to underflow
244  // each bin moves to the right
245 
246  TArrayD* sumw2 = p->GetSumw2();
247 
248  for (int i = 0; i <= nBins + 1 - bins; i++) {
249  // GetBinContent return binContent/binEntries
250  p->SetBinContent(i, p->GetBinContent(i + bins) * p->GetBinEntries(i + bins));
251  p->SetBinEntries(i, p->GetBinEntries(i + bins));
252  sumw2->SetAt(sumw2->GetAt(i + bins), i);
253  }
254 }
255 
void dqmBeginRun(const edm::Run &r, const edm::EventSetup &c) override
Definition: ESTrendTask.cc:45
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:418
std::vector< T >::const_iterator const_iterator
void shift2Right(TProfile *p, int bins=1)
Definition: ESTrendTask.cc:200
int zside(DetId const &)
const std::vector< int > & getFEChannelStatus() const
void shift2Left(TProfile *p, int bins=1)
Definition: ESTrendTask.cc:228
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: ESTrendTask.cc:51
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
MonitorElement * bookProfile(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, char const *option="s")
Definition: DQMStore.cc:333
unsigned int unixTime() const
Time in seconds since January 1, 1970.
Definition: Timestamp.h:40
int zside() const
Definition: ESDetId.h:39
Namespace of DDCMS conversion namespace.
const_iterator end() const
ESTrendTask(const edm::ParameterSet &ps)
Definition: ESTrendTask.cc:22
Timestamp const & beginTime() const
Definition: RunBase.h:41
HLT enums.
int getDCCErrors() const
void analyze(const edm::Event &e, const edm::EventSetup &c) override
Definition: ESTrendTask.cc:93
void updateTime(const edm::Event &)
Definition: ESTrendTask.cc:195
edm::Timestamp time() const
Definition: EventBase.h:60
const_iterator begin() const
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)
void endJob(void) override
Definition: ESTrendTask.cc:91