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
15 {
18  gainToME_(),
19  pnGainToME_(),
20  minChannelEntries_(0),
21  expectedMean_(0.),
22  toleranceMean_(0.),
23  toleranceRMSEB_(0),
24  toleranceRMSEE_(0),
25  expectedPNMean_(0.),
26  tolerancePNMean_(0.),
27  tolerancePNRMS_(0)
28  {
29  }
30 
31  void
33  {
34  minChannelEntries_ = _params.getUntrackedParameter<int>("minChannelEntries");
35  expectedMean_ = _params.getUntrackedParameter<double>("expectedMean");
36  toleranceMean_ = _params.getUntrackedParameter<double>("toleranceMean");
37  expectedPNMean_ = _params.getUntrackedParameter<double>("expectedPNMean");
38  tolerancePNMean_ = _params.getUntrackedParameter<double>("tolerancePNMean");
39 
40  std::vector<int> MGPAGains(_params.getUntrackedParameter<std::vector<int> >("MGPAGains"));
41  std::vector<int> MGPAGainsPN(_params.getUntrackedParameter<std::vector<int> >("MGPAGainsPN"));
42 
44 
45  MESetMulti const& pedestal(static_cast<MESetMulti const&>(sources_.at("Pedestal")));
46  unsigned nG(MGPAGains.size());
47  for(unsigned iG(0); iG != nG; ++iG){
48  int gain(MGPAGains[iG]);
49  if(gain != 1 && gain != 6 && gain != 12) throw cms::Exception("InvalidConfiguration") << "MGPA gain";
50  repl["gain"] = std::to_string(gain);
51  gainToME_[gain] = pedestal.getIndex(repl);
52  }
53 
54  repl.clear();
55 
56  MESetMulti const& pnPedestal(static_cast<MESetMulti const&>(sources_.at("PNPedestal")));
57  unsigned nGPN(MGPAGainsPN.size());
58  for(unsigned iG(0); iG != nGPN; ++iG){
59  int gain(MGPAGainsPN[iG]);
60  if(gain != 1 && gain != 16) throw cms::Exception("InvalidConfiguration") << "PN MGPA gain";
61  repl["pngain"] = std::to_string(gain);
62  pnGainToME_[gain] = pnPedestal.getIndex(repl);
63  }
64 
65  toleranceRMSEB_.resize(nG);
66  toleranceRMSEE_.resize(nG);
67 
68  std::vector<double> inToleranceRMSEB(_params.getUntrackedParameter<std::vector<double> >("toleranceRMSEB"));
69  std::vector<double> inToleranceRMSEE(_params.getUntrackedParameter<std::vector<double> >("toleranceRMSEE"));
70 
71  for(std::map<int, unsigned>::iterator gainItr(gainToME_.begin()); gainItr != gainToME_.end(); ++gainItr){
72  unsigned iME(gainItr->second);
73  unsigned iGain(0);
74  switch(gainItr->first){
75  case 1:
76  iGain = 0; break;
77  case 6:
78  iGain = 1; break;
79  case 12:
80  iGain = 2; 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; break;
97  case 16:
98  iGain = 1; break;
99  }
100 
101  tolerancePNRMS_[iME] = inTolerancePNRMS[iGain];
102  }
103 
104  qualitySummaries_.insert("Quality");
105  qualitySummaries_.insert("QualitySummary");
106  qualitySummaries_.insert("PNQualitySummary");
107  }
108 
109  void
111  {
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());
151  MESet::const_iterator pItr(sPedestal);
152  for(MESet::iterator qItr(meQuality.beginChannel()); qItr != qEnd; qItr.toNextChannel()){
153 
154  DetId id(qItr->getId());
155 
156  bool doMask(meQuality.maskMatches(id, mask, statusManager_));
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(id, mean);
171  meRMS.fill(id, rms);
172 
173  float toleranceRMS_ = ( id.subdetId() == EcalBarrel ) ? toleranceRMSEB_[gainItr->second] : toleranceRMSEE_[gainItr->second];
174 
175  if(abs(mean - expectedMean_) > toleranceMean_ || rms > toleranceRMS_)
176  qItr->setBinContent(doMask ? kMBad : kBad);
177  else
178  qItr->setBinContent(doMask ? kMGood : kGood);
179  }
180 
181  towerAverage_(meQualitySummary, meQuality, 0.2);
182  }
183 
184  for(map<int, unsigned>::iterator gainItr(pnGainToME_.begin()); gainItr != pnGainToME_.end(); ++gainItr){
185  mePNQualitySummary.use(gainItr->second);
186  mePNRMS.use(gainItr->second);
187 
188  sPNPedestal.use(gainItr->second);
189 
190  uint32_t mask(0);
191  switch(gainItr->first){
192  case 1:
195  break;
196  case 16:
199  break;
200  default:
201  break;
202  }
203 
204  for(unsigned iDCC(0); iDCC < nDCC; ++iDCC){
205 
206  if(memDCCIndex(iDCC + 1) == unsigned(-1)) continue;
207 
208  for(unsigned iPN(0); iPN < 10; ++iPN){
209  int subdet(0);
210  if(iDCC >= kEBmLow && iDCC <= kEBpHigh) subdet = EcalBarrel;
211  else subdet = EcalEndcap;
212 
213  EcalPnDiodeDetId id(subdet, iDCC + 1, iPN + 1);
214 
215  bool doMask(mePNQualitySummary.maskMatches(id, mask, statusManager_));
216 
217  float entries(sPNPedestal.getBinEntries(id));
218 
219  if(entries < minChannelEntries_){
220  mePNQualitySummary.setBinContent(id, doMask ? kMUnknown : kUnknown);
221  continue;
222  }
223 
224  float mean(sPNPedestal.getBinContent(id));
225  float rms(sPNPedestal.getBinError(id) * sqrt(entries));
226 
227  mePNRMS.fill(id, rms);
228 
229  if(abs(mean - expectedPNMean_) > tolerancePNMean_ || rms > tolerancePNRMS_[gainItr->second])
230  mePNQualitySummary.setBinContent(id, doMask ? kMBad : kBad);
231  else
232  mePNQualitySummary.setBinContent(id, doMask ? kMGood : kGood);
233  }
234  }
235  }
236  }
237 
239 }
unsigned memDCCIndex(unsigned)
T getUntrackedParameter(std::string const &, T const &) const
static const int PEDESTAL_MIDDLE_GAIN_RMS_ERROR
void setParams(edm::ParameterSet const &) override
void towerAverage_(MESet &, MESet const &, float)
#define DEFINE_ECALDQM_WORKER(TYPE)
Definition: DQWorker.h:109
static const int PEDESTAL_HIGH_GAIN_MEAN_ERROR
double getBinContent(DetId const &_id, int _bin=0) const override
Definition: MESetMulti.h:81
void setBinContent(DetId const &_id, double _content) override
Definition: MESetMulti.h:40
double getBinError(DetId const &_id, int _bin=0) const override
Definition: MESetMulti.h:89
static const int PEDESTAL_LOW_GAIN_RMS_ERROR
const_iterator & toNextChannel()
Definition: MESet.h:290
void fill(DetId const &_id, double _xyw=1., double _yw=1., double _w=1.) override
Definition: MESetMulti.h:29
static const int PEDESTAL_HIGH_GAIN_RMS_ERROR
const_iterator beginChannel() const override
Definition: MESetMulti.h:123
void use(unsigned) const
Definition: MESetMulti.cc:130
std::set< std::string > qualitySummaries_
T sqrt(T t)
Definition: SSEVec.h:18
StatusManager const * statusManager_
std::map< int, unsigned > pnGainToME_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool maskMatches(DetId const &_id, uint32_t _mask, StatusManager const *_statusManager) const override
Definition: MESetMulti.h:106
MESetCollection sources_
const_iterator end() const override
Definition: MESetMulti.h:122
static const int PEDESTAL_LOW_GAIN_MEAN_ERROR
Definition: DetId.h:18
std::vector< float > toleranceRMSEB_
std::vector< float > tolerancePNRMS_
MESetCollection MEs_
Definition: DQWorker.h:75
void producePlots(ProcessType) override
double getBinEntries(DetId const &_id, int _bin=0) const override
Definition: MESetMulti.h:95
std::map< int, unsigned > gainToME_
static const int PEDESTAL_MIDDLE_GAIN_MEAN_ERROR
std::map< std::string, std::string > PathReplacements
Definition: MESet.h:30
unsigned getIndex(PathReplacements const &) const
Definition: MESetMulti.cc:137
std::vector< float > toleranceRMSEE_