CMS 3D CMS Logo

ECALpedestalPCLHarvester.cc
Go to the documentation of this file.
1 
2 // system include files
3 #include <memory>
4 
5 // user include files
11 #include <string>
12 
14  : minEntries_(ps.getParameter<int>("MinEntries")),
15  checkAnomalies_(ps.getParameter<bool>("checkAnomalies")),
16  nSigma_(ps.getParameter<double>("nSigma")),
17  thresholdAnomalies_(ps.getParameter<double>("thresholdAnomalies")),
18  dqmDir_(ps.getParameter<std::string>("dqmDir")),
19  labelG6G1_(ps.getParameter<std::string>("labelG6G1")),
20  threshDiffEB_(ps.getParameter<double>("threshDiffEB")),
21  threshDiffEE_(ps.getParameter<double>("threshDiffEE")),
22  threshChannelsAnalyzed_(ps.getParameter<double>("threshChannelsAnalyzed")),
23  channelsStatusToken_(esConsumes<edm::Transition::EndRun>()),
24  pedestalsToken_(esConsumes<edm::Transition::EndRun>()),
25  g6g1PedestalsToken_(esConsumes<edm::Transition::EndRun>(edm::ESInputTag("", labelG6G1_))),
26  currentPedestals_(nullptr),
27  g6g1Pedestals_(nullptr),
28  channelStatus_(nullptr) {
29  chStatusToExclude_ = StringToEnumValue<EcalChannelStatusCode::Code>(
30  ps.getParameter<std::vector<std::string> >("ChannelStatusToExclude"));
31 }
32 
34  // calculate pedestals and fill db record
35  EcalPedestals pedestals;
36  float kBarrelSize = 61200;
37  float kEndcapSize = 2 * 7324;
38  float skipped_channels_EB = 0;
39  float skipped_channels_EE = 0;
40 
41  for (uint16_t i = 0; i < EBDetId::kSizeForDenseIndexing; ++i) {
42  std::string hname = dqmDir_ + "/EB/" + std::to_string(int(i / 100)) + "/eb_" + std::to_string(i);
43  MonitorElement* ch = igetter_.get(hname);
44  if (ch == nullptr) {
45  edm::LogWarning("MissingMonitorElement") << "failed to find MonitorElement " << hname;
46  entriesEB_[i] = 0;
47  continue;
48  }
49  double mean = ch->getMean();
50  double rms = ch->getRMS();
51  entriesEB_[i] = ch->getEntries();
52 
54  EcalPedestal ped;
55  EcalPedestal oldped = *currentPedestals_->find(id.rawId());
56  EcalPedestal g6g1ped = *g6g1Pedestals_->find(id.rawId());
57 
58  ped.mean_x12 = mean;
59  ped.rms_x12 = rms;
60 
61  float diff = std::abs(mean - oldped.mean_x12);
62 
63  // if bad channel or low stat skip or the difference is too large wrt to previous record
64  if (ch->getEntries() < minEntries_ || !checkStatusCode(id) || diff > threshDiffEB_) {
65  ped.mean_x12 = oldped.mean_x12;
66  ped.rms_x12 = oldped.rms_x12;
67 
68  skipped_channels_EB++;
69  }
70 
71  // copy g6 and g1 from the corressponding record
72  ped.mean_x6 = g6g1ped.mean_x6;
73  ped.rms_x6 = g6g1ped.rms_x6;
74  ped.mean_x1 = g6g1ped.mean_x1;
75  ped.rms_x1 = g6g1ped.rms_x1;
76 
77  pedestals.setValue(id.rawId(), ped);
78  }
79 
80  for (uint16_t i = 0; i < EEDetId::kSizeForDenseIndexing; ++i) {
81  std::string hname = dqmDir_ + "/EE/" + std::to_string(int(i / 100)) + "/ee_" + std::to_string(i);
82 
83  MonitorElement* ch = igetter_.get(hname);
84  if (ch == nullptr) {
85  edm::LogWarning("MissingMonitorElement") << "failed to find MonitorElement " << hname;
86  entriesEE_[i] = 0;
87  continue;
88  }
89  double mean = ch->getMean();
90  double rms = ch->getRMS();
91  entriesEE_[i] = ch->getEntries();
92 
94  EcalPedestal ped;
95  EcalPedestal oldped = *currentPedestals_->find(id.rawId());
96  EcalPedestal g6g1ped = *g6g1Pedestals_->find(id.rawId());
97 
98  ped.mean_x12 = mean;
99  ped.rms_x12 = rms;
100 
101  float diff = std::abs(mean - oldped.mean_x12);
102 
103  // if bad channel or low stat skip or the difference is too large wrt to previous record
104  if (ch->getEntries() < minEntries_ || !checkStatusCode(id) || diff > threshDiffEE_) {
105  ped.mean_x12 = oldped.mean_x12;
106  ped.rms_x12 = oldped.rms_x12;
107 
108  skipped_channels_EE++;
109  }
110 
111  // copy g6 and g1 pedestals from corresponding record
112  ped.mean_x6 = g6g1ped.mean_x6;
113  ped.rms_x6 = g6g1ped.rms_x6;
114  ped.mean_x1 = g6g1ped.mean_x1;
115  ped.rms_x1 = g6g1ped.rms_x1;
116 
117  pedestals.setValue(id.rawId(), ped);
118  }
119 
120  bool enough_stat = false;
121  if ((kBarrelSize - skipped_channels_EB) / kBarrelSize > threshChannelsAnalyzed_ &&
122  (kEndcapSize - skipped_channels_EE) / kEndcapSize > threshChannelsAnalyzed_) {
123  enough_stat = true;
124  }
125 
126  dqmPlots(pedestals, ibooker_);
127 
128  // check if there are large variations wrt exisiting pedstals
129 
130  if (checkAnomalies_) {
131  if (checkVariation(*currentPedestals_, pedestals)) {
132  edm::LogError("Large Variations found wrt to old pedestals, no file created");
133  return;
134  }
135  }
136 
137  if (!enough_stat)
138  return;
139 
140  // write out pedestal record
142 
143  if (!poolDbService.isAvailable()) {
144  throw std::runtime_error("PoolDBService required.");
145  }
146 
147  poolDbService->writeOneIOV(pedestals, poolDbService->currentTime(), "EcalPedestalsRcd");
148 }
149 
150 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
153  desc.setUnknown();
154  descriptions.addDefault(desc);
155 }
156 
161 }
162 
165  dbstatusPtr = channelStatus_->getMap().find(id.rawId());
166  EcalChannelStatusCode::Code dbstatus = dbstatusPtr->getStatusCode();
167 
168  std::vector<int>::const_iterator res = std::find(chStatusToExclude_.begin(), chStatusToExclude_.end(), dbstatus);
169  if (res != chStatusToExclude_.end())
170  return false;
171 
172  return true;
173 }
174 
177  dbstatusPtr = channelStatus_->getMap().find(id.rawId());
178  if (dbstatusPtr == channelStatus_->getMap().end())
179  edm::LogError("Invalid DetId supplied");
180  EcalChannelStatusCode::Code dbstatus = dbstatusPtr->getStatusCode();
181  if (dbstatus == 0)
182  return true;
183  return false;
184 }
185 
187  const EcalPedestalsMap& newPedestals) {
188  uint32_t nAnomaliesEB = 0;
189  uint32_t nAnomaliesEE = 0;
190 
191  for (uint16_t i = 0; i < EBDetId::kSizeForDenseIndexing; ++i) {
193  const EcalPedestal& newped = *newPedestals.find(id.rawId());
194  const EcalPedestal& oldped = *oldPedestals.find(id.rawId());
195 
196  if (std::abs(newped.mean_x12 - oldped.mean_x12) > nSigma_ * oldped.rms_x12)
197  nAnomaliesEB++;
198  }
199 
200  for (uint16_t i = 0; i < EEDetId::kSizeForDenseIndexing; ++i) {
202  const EcalPedestal& newped = *newPedestals.find(id.rawId());
203  const EcalPedestal& oldped = *oldPedestals.find(id.rawId());
204 
205  if (std::abs(newped.mean_x12 - oldped.mean_x12) > nSigma_ * oldped.rms_x12)
206  nAnomaliesEE++;
207  }
208 
211  return true;
212 
213  return false;
214 }
215 
217  ibooker.cd();
218  ibooker.setCurrentFolder(dqmDir_ + "/Summary");
219 
220  MonitorElement* pmeb = ibooker.book2D("meaneb", "Pedestal Means EB", 360, 1., 361., 171, -85., 86.);
221  MonitorElement* preb = ibooker.book2D("rmseb", "Pedestal RMS EB ", 360, 1., 361., 171, -85., 86.);
222 
223  MonitorElement* pmeep = ibooker.book2D("meaneep", "Pedestal Means EEP", 100, 1, 101, 100, 1, 101);
224  MonitorElement* preep = ibooker.book2D("rmseep", "Pedestal RMS EEP", 100, 1, 101, 100, 1, 101);
225 
226  MonitorElement* pmeem = ibooker.book2D("meaneem", "Pedestal Means EEM", 100, 1, 101, 100, 1, 101);
227  MonitorElement* preem = ibooker.book2D("rmseem", "Pedestal RMS EEM", 100, 1, 101, 100, 1, 101);
228 
229  MonitorElement* pmebd = ibooker.book2D("meanebdiff", "Abs Rel Pedestal Means Diff EB", 360, 1., 361., 171, -85., 86.);
230  MonitorElement* prebd = ibooker.book2D("rmsebdiff", "Abs Rel Pedestal RMS Diff E ", 360, 1., 361., 171, -85., 86.);
231 
232  MonitorElement* pmeepd = ibooker.book2D("meaneepdiff", "Abs Rel Pedestal Means Diff EEP", 100, 1, 101, 100, 1, 101);
233  MonitorElement* preepd = ibooker.book2D("rmseepdiff", "Abs Rel Pedestal RMS Diff EEP", 100, 1, 101, 100, 1, 101);
234 
235  MonitorElement* pmeemd = ibooker.book2D("meaneemdiff", "Abs Rel Pedestal Means Diff EEM", 100, 1, 101, 100, 1, 101);
236  MonitorElement* preemd = ibooker.book2D("rmseemdiff", "Abs RelPedestal RMS Diff EEM", 100, 1, 101, 100, 1, 101);
237 
238  MonitorElement* poeb = ibooker.book2D("occeb", "Occupancy EB", 360, 1., 361., 171, -85., 86.);
239  MonitorElement* poeep = ibooker.book2D("occeep", "Occupancy EEP", 100, 1, 101, 100, 1, 101);
240  MonitorElement* poeem = ibooker.book2D("occeem", "Occupancy EEM", 100, 1, 101, 100, 1, 101);
241 
242  MonitorElement* hdiffeb = ibooker.book1D("diffeb", "Pedestal Differences EB", 100, -2.5, 2.5);
243  MonitorElement* hdiffee = ibooker.book1D("diffee", "Pedestal Differences EE", 100, -2.5, 2.5);
244 
245  for (int hash = 0; hash < EBDetId::kSizeForDenseIndexing; ++hash) {
247  float mean = newpeds[di].mean_x12;
248  float rms = newpeds[di].rms_x12;
249 
250  float cmean = (*currentPedestals_)[di].mean_x12;
251  float crms = (*currentPedestals_)[di].rms_x12;
252 
253  if (!isGood(di))
254  continue; // only good channels are plotted
255 
256  pmeb->Fill(di.iphi(), di.ieta(), mean);
257  preb->Fill(di.iphi(), di.ieta(), rms);
258  if (cmean)
259  pmebd->Fill(di.iphi(), di.ieta(), std::abs(mean - cmean) / cmean);
260  if (crms)
261  prebd->Fill(di.iphi(), di.ieta(), std::abs(rms - crms) / crms);
262  poeb->Fill(di.iphi(), di.ieta(), entriesEB_[hash]);
263  hdiffeb->Fill(mean - cmean);
264  }
265 
266  for (int hash = 0; hash < EEDetId::kSizeForDenseIndexing; ++hash) {
268  float mean = newpeds[di].mean_x12;
269  float rms = newpeds[di].rms_x12;
270  float cmean = (*currentPedestals_)[di].mean_x12;
271  float crms = (*currentPedestals_)[di].rms_x12;
272 
273  if (!isGood(di))
274  continue; // only good channels are plotted
275 
276  if (di.zside() > 0) {
277  pmeep->Fill(di.ix(), di.iy(), mean);
278  preep->Fill(di.ix(), di.iy(), rms);
279  poeep->Fill(di.ix(), di.iy(), entriesEE_[hash]);
280  if (cmean)
281  pmeepd->Fill(di.ix(), di.iy(), std::abs(mean - cmean) / cmean);
282  if (crms)
283  preepd->Fill(di.ix(), di.iy(), std::abs(rms - crms) / crms);
284  } else {
285  pmeem->Fill(di.ix(), di.iy(), mean);
286  preem->Fill(di.ix(), di.iy(), rms);
287  if (cmean)
288  pmeemd->Fill(di.ix(), di.iy(), std::abs(mean - cmean) / cmean);
289  poeem->Fill(di.ix(), di.iy(), entriesEE_[hash]);
290  if (crms)
291  preemd->Fill(di.ix(), di.iy(), std::abs(rms - crms) / crms);
292  }
293  hdiffee->Fill(mean - cmean);
294  }
295 }
static EEDetId detIdFromDenseIndex(uint32_t din)
Definition: EEDetId.h:220
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
const edm::ESGetToken< EcalPedestals, EcalPedestalsRcd > g6g1PedestalsToken_
int entriesEE_[EEDetId::kSizeForDenseIndexing]
int iphi() const
get the crystal iphi
Definition: EBDetId.h:51
bool checkVariation(const EcalPedestalsMap &oldPedestals, const EcalPedestalsMap &newPedestals)
int ix() const
Definition: EEDetId.h:77
std::string to_string(const V &value)
Definition: OMSAccess.h:71
Log< level::Error, false > LogError
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
void dqmEndJob(DQMStore::IBooker &ibooker_, DQMStore::IGetter &igetter_) override
Definition: Electron.h:6
int ieta() const
get the crystal ieta
Definition: EBDetId.h:49
static EBDetId detIdFromDenseIndex(uint32_t di)
Definition: EBDetId.h:107
void Fill(long long x)
void setValue(const uint32_t id, const Item &item)
ECALpedestalPCLHarvester(const edm::ParameterSet &ps)
virtual double getRMS(int axis=1) const
get RMS of histogram along x, y or z axis (axis=1, 2, 3 respectively)
void addDefault(ParameterSetDescription const &psetDescription)
const EcalPedestals * g6g1Pedestals_
std::vector< int > chStatusToExclude_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Hash writeOneIOV(const T &payload, Time_t time, const std::string &recordName)
Transition
Definition: Transition.h:12
virtual double getEntries() const
get # of entries
bool getData(T &iHolder) const
Definition: EventSetup.h:122
void dqmPlots(const EcalPedestals &newpeds, DQMStore::IBooker &ibooker)
const_iterator find(uint32_t rawId) const
Definition: DetId.h:17
int zside() const
Definition: EEDetId.h:71
const edm::ESGetToken< EcalPedestals, EcalPedestalsRcd > pedestalsToken_
const EcalChannelStatus * channelStatus_
std::vector< Item >::const_iterator const_iterator
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:212
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:690
virtual double getMean(int axis=1) const
get mean value of histogram along x, y or z axis (axis=1, 2, 3 respectively)
bool checkStatusCode(const DetId &id)
const EcalPedestals * currentPedestals_
HLT enums.
const_iterator end() const
bool isAvailable() const
Definition: Service.h:40
void endRun(edm::Run const &run, edm::EventSetup const &isetup) override
Log< level::Warning, false > LogWarning
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
int entriesEB_[EBDetId::kSizeForDenseIndexing]
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const edm::ESGetToken< EcalChannelStatus, EcalChannelStatusRcd > channelsStatusToken_
Definition: Run.h:45
int iy() const
Definition: EEDetId.h:83