CMS 3D CMS Logo

PedestalClient.cc
Go to the documentation of this file.
2 
5 
7 
9 
11 
12 #include <iomanip>
13 
14 namespace ecaldqm {
16  : DQWorkerClient(),
17  gainToME_(),
18  pnGainToME_(),
19  minChannelEntries_(0),
20  expectedMean_(0.),
21  toleranceMean_(0.),
22  toleranceRMSEB_(0),
23  toleranceRMSEE_(0),
24  expectedPNMean_(0.),
25  tolerancePNMean_(0.),
26  tolerancePNRMS_(0) {}
27 
29  minChannelEntries_ = _params.getUntrackedParameter<int>("minChannelEntries");
30  expectedMean_ = _params.getUntrackedParameter<double>("expectedMean");
31  toleranceMean_ = _params.getUntrackedParameter<double>("toleranceMean");
32  expectedPNMean_ = _params.getUntrackedParameter<double>("expectedPNMean");
33  tolerancePNMean_ = _params.getUntrackedParameter<double>("tolerancePNMean");
34 
35  std::vector<int> MGPAGains(_params.getUntrackedParameter<std::vector<int> >("MGPAGains"));
36  std::vector<int> MGPAGainsPN(_params.getUntrackedParameter<std::vector<int> >("MGPAGainsPN"));
37 
39 
40  MESetMulti const& pedestal(static_cast<MESetMulti const&>(sources_.at("Pedestal")));
41  unsigned nG(MGPAGains.size());
42  for (unsigned iG(0); iG != nG; ++iG) {
43  int gain(MGPAGains[iG]);
44  if (gain != 1 && gain != 6 && gain != 12)
45  throw cms::Exception("InvalidConfiguration") << "MGPA gain";
46  repl["gain"] = std::to_string(gain);
47  gainToME_[gain] = pedestal.getIndex(repl);
48  }
49 
50  repl.clear();
51 
52  MESetMulti const& pnPedestal(static_cast<MESetMulti const&>(sources_.at("PNPedestal")));
53  unsigned nGPN(MGPAGainsPN.size());
54  for (unsigned iG(0); iG != nGPN; ++iG) {
55  int gain(MGPAGainsPN[iG]);
56  if (gain != 1 && gain != 16)
57  throw cms::Exception("InvalidConfiguration") << "PN MGPA gain";
58  repl["pngain"] = std::to_string(gain);
59  pnGainToME_[gain] = pnPedestal.getIndex(repl);
60  }
61 
62  toleranceRMSEB_.resize(nG);
63  toleranceRMSEE_.resize(nG);
64 
65  std::vector<double> inToleranceRMSEB(_params.getUntrackedParameter<std::vector<double> >("toleranceRMSEB"));
66  std::vector<double> inToleranceRMSEE(_params.getUntrackedParameter<std::vector<double> >("toleranceRMSEE"));
67 
68  for (std::map<int, unsigned>::iterator gainItr(gainToME_.begin()); gainItr != gainToME_.end(); ++gainItr) {
69  unsigned iME(gainItr->second);
70  unsigned iGain(0);
71  switch (gainItr->first) {
72  case 1:
73  iGain = 0;
74  break;
75  case 6:
76  iGain = 1;
77  break;
78  case 12:
79  iGain = 2;
80  break;
81  }
82 
83  toleranceRMSEB_[iME] = inToleranceRMSEB[iGain];
84  toleranceRMSEE_[iME] = inToleranceRMSEE[iGain];
85  }
86 
87  tolerancePNRMS_.resize(nGPN);
88 
89  std::vector<double> inTolerancePNRMS(_params.getUntrackedParameter<std::vector<double> >("tolerancePNRMS"));
90 
91  for (std::map<int, unsigned>::iterator gainItr(pnGainToME_.begin()); gainItr != pnGainToME_.end(); ++gainItr) {
92  unsigned iME(gainItr->second);
93  unsigned iGain(0);
94  switch (gainItr->first) {
95  case 1:
96  iGain = 0;
97  break;
98  case 16:
99  iGain = 1;
100  break;
101  }
102 
103  tolerancePNRMS_[iME] = inTolerancePNRMS[iGain];
104  }
105 
106  qualitySummaries_.insert("Quality");
107  qualitySummaries_.insert("QualitySummary");
108  qualitySummaries_.insert("PNQualitySummary");
109  }
110 
112  using namespace std;
113 
114  MESetMulti& meQuality(static_cast<MESetMulti&>(MEs_.at("Quality")));
115  MESetMulti& meQualitySummary(static_cast<MESetMulti&>(MEs_.at("QualitySummary")));
116  MESetMulti& meMean(static_cast<MESetMulti&>(MEs_.at("Mean")));
117  MESetMulti& meRMS(static_cast<MESetMulti&>(MEs_.at("RMS")));
118  MESetMulti& mePNQualitySummary(static_cast<MESetMulti&>(MEs_.at("PNQualitySummary")));
119  MESetMulti& mePNRMS(static_cast<MESetMulti&>(MEs_.at("PNRMS")));
120 
121  MESetMulti const& sPedestal(static_cast<MESetMulti const&>(sources_.at("Pedestal")));
122  MESetMulti const& sPNPedestal(static_cast<MESetMulti const&>(sources_.at("PNPedestal")));
123 
124  for (map<int, unsigned>::iterator gainItr(gainToME_.begin()); gainItr != gainToME_.end(); ++gainItr) {
125  meQuality.use(gainItr->second);
126  meQualitySummary.use(gainItr->second);
127  meMean.use(gainItr->second);
128  meRMS.use(gainItr->second);
129 
130  sPedestal.use(gainItr->second);
131 
132  uint32_t mask(0);
133  switch (gainItr->first) {
134  case 1:
137  break;
138  case 6:
141  break;
142  case 12:
145  break;
146  default:
147  break;
148  }
149 
150  MESet::iterator qEnd(meQuality.end(GetElectronicsMap()));
151  MESet::const_iterator pItr(GetElectronicsMap(), sPedestal);
152  for (MESet::iterator qItr(meQuality.beginChannel(GetElectronicsMap())); qItr != qEnd;
154  DetId id(qItr->getId());
155 
156  bool doMask(meQuality.maskMatches(id, mask, statusManager_, GetTrigTowerMap()));
157 
158  pItr = qItr;
159 
160  float entries(pItr->getBinEntries());
161 
162  if (entries < minChannelEntries_) {
163  qItr->setBinContent(doMask ? kMUnknown : kUnknown);
164  continue;
165  }
166 
167  float mean(pItr->getBinContent());
168  float rms(pItr->getBinError() * sqrt(entries));
169 
170  meMean.fill(getEcalDQMSetupObjects(), id, mean);
171  meRMS.fill(getEcalDQMSetupObjects(), id, rms);
172 
173  float toleranceRMS_ =
174  (id.subdetId() == EcalBarrel) ? toleranceRMSEB_[gainItr->second] : toleranceRMSEE_[gainItr->second];
175 
176  if (abs(mean - expectedMean_) > toleranceMean_ || rms > toleranceRMS_)
177  qItr->setBinContent(doMask ? kMBad : kBad);
178  else
179  qItr->setBinContent(doMask ? kMGood : kGood);
180  }
181 
182  towerAverage_(meQualitySummary, meQuality, 0.2);
183  }
184 
185  for (map<int, unsigned>::iterator gainItr(pnGainToME_.begin()); gainItr != pnGainToME_.end(); ++gainItr) {
186  mePNQualitySummary.use(gainItr->second);
187  mePNRMS.use(gainItr->second);
188 
189  sPNPedestal.use(gainItr->second);
190 
191  uint32_t mask(0);
192  switch (gainItr->first) {
193  case 1:
196  break;
197  case 16:
200  break;
201  default:
202  break;
203  }
204 
205  for (unsigned iDCC(0); iDCC < nDCC; ++iDCC) {
206  if (memDCCIndex(iDCC + 1) == unsigned(-1))
207  continue;
208 
209  for (unsigned iPN(0); iPN < 10; ++iPN) {
210  int subdet(0);
211  if (iDCC >= kEBmLow && iDCC <= kEBpHigh)
212  subdet = EcalBarrel;
213  else
214  subdet = EcalEndcap;
215 
216  EcalPnDiodeDetId id(subdet, iDCC + 1, iPN + 1);
217 
218  bool doMask(mePNQualitySummary.maskMatches(id, mask, statusManager_, GetTrigTowerMap()));
219 
220  float entries(sPNPedestal.getBinEntries(getEcalDQMSetupObjects(), id));
221 
222  if (entries < minChannelEntries_) {
223  mePNQualitySummary.setBinContent(getEcalDQMSetupObjects(), id, doMask ? kMUnknown : kUnknown);
224  continue;
225  }
226 
227  float mean(sPNPedestal.getBinContent(getEcalDQMSetupObjects(), id));
228  float rms(sPNPedestal.getBinError(getEcalDQMSetupObjects(), id) * sqrt(entries));
229 
230  mePNRMS.fill(getEcalDQMSetupObjects(), id, rms);
231 
232  if (abs(mean - expectedPNMean_) > tolerancePNMean_ || rms > tolerancePNRMS_[gainItr->second])
233  mePNQualitySummary.setBinContent(getEcalDQMSetupObjects(), id, doMask ? kMBad : kBad);
234  else
235  mePNQualitySummary.setBinContent(getEcalDQMSetupObjects(), id, doMask ? kMGood : kGood);
236  }
237  }
238  }
239  }
240 
242 } // namespace ecaldqm
unsigned memDCCIndex(unsigned)
static const int PEDESTAL_MIDDLE_GAIN_RMS_ERROR
void setParams(edm::ParameterSet const &) override
void use(unsigned) const
Definition: MESetMulti.cc:135
void towerAverage_(MESet &, MESet const &, float)
#define DEFINE_ECALDQM_WORKER(TYPE)
Definition: DQWorker.h:168
MESet & at(const std::string &key)
Definition: MESet.h:399
const_iterator & toNextChannel(EcalElectronicsMapping const *electronicsMap)
Definition: MESet.h:320
double getBinContent(EcalDQMSetupObjects const edso, DetId const &_id, int _bin=0) const override
Definition: MESetMulti.h:111
static const int PEDESTAL_HIGH_GAIN_MEAN_ERROR
std::string to_string(const V &value)
Definition: OMSAccess.h:71
const_iterator beginChannel(EcalElectronicsMapping const *electronicsMap) const override
Definition: MESetMulti.h:168
static const int PEDESTAL_LOW_GAIN_RMS_ERROR
static const int PEDESTAL_HIGH_GAIN_RMS_ERROR
constexpr uint32_t mask
Definition: gpuClustering.h:26
void fill(EcalDQMSetupObjects const edso, DetId const &_id, double _xyw=1., double _yw=1., double _w=1.) override
Definition: MESetMulti.h:29
std::set< std::string > qualitySummaries_
T sqrt(T t)
Definition: SSEVec.h:19
StatusManager const * statusManager_
void setBinContent(EcalDQMSetupObjects const edso, DetId const &_id, double _content) override
Definition: MESetMulti.h:48
std::map< int, unsigned > pnGainToME_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
MESetCollection sources_
double getBinError(EcalDQMSetupObjects const edso, DetId const &_id, int _bin=0) const override
Definition: MESetMulti.h:121
EcalElectronicsMapping const * GetElectronicsMap()
Definition: DQWorker.cc:150
static const int PEDESTAL_LOW_GAIN_MEAN_ERROR
EcalDQMSetupObjects const getEcalDQMSetupObjects()
Definition: DQWorker.cc:170
Definition: DetId.h:17
unsigned getIndex(PathReplacements const &) const
Definition: MESetMulti.cc:142
std::vector< float > toleranceRMSEB_
std::vector< float > tolerancePNRMS_
double getBinEntries(EcalDQMSetupObjects const edso, DetId const &_id, int _bin=0) const override
Definition: MESetMulti.h:131
MESetCollection MEs_
Definition: DQWorker.h:131
void producePlots(ProcessType) override
std::map< int, unsigned > gainToME_
EcalTrigTowerConstituentsMap const * GetTrigTowerMap()
Definition: DQWorker.cc:155
static const int PEDESTAL_MIDDLE_GAIN_MEAN_ERROR
bool maskMatches(DetId const &_id, uint32_t _mask, StatusManager const *_statusManager, EcalTrigTowerConstituentsMap const *trigTowerMap) const override
Definition: MESetMulti.h:144
std::map< std::string, std::string > PathReplacements
Definition: MESet.h:46
const_iterator end(EcalElectronicsMapping const *electronicsMap) const override
Definition: MESetMulti.h:165
std::vector< float > toleranceRMSEE_