CMS 3D CMS Logo

TestPulseClient.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  amplitudeThreshold_(0),
21  toleranceRMS_(0),
22  PNAmplitudeThreshold_(0),
23  tolerancePNRMS_(0) {}
24 
26  minChannelEntries_ = _params.getUntrackedParameter<int>("minChannelEntries");
27 
28  std::vector<int> MGPAGains(_params.getUntrackedParameter<std::vector<int> >("MGPAGains"));
29  std::vector<int> MGPAGainsPN(_params.getUntrackedParameter<std::vector<int> >("MGPAGainsPN"));
30 
32 
33  MESetMulti const& amplitude(static_cast<MESetMulti const&>(sources_.at("Amplitude")));
34  unsigned nG(MGPAGains.size());
35  for (unsigned iG(0); iG != nG; ++iG) {
36  int gain(MGPAGains[iG]);
37  if (gain != 1 && gain != 6 && gain != 12)
38  throw cms::Exception("InvalidConfiguration") << "MGPA gain";
39  repl["gain"] = std::to_string(gain);
40  gainToME_[gain] = amplitude.getIndex(repl);
41  }
42 
43  repl.clear();
44 
45  MESetMulti const& pnAmplitude(static_cast<MESetMulti const&>(sources_.at("PNAmplitude")));
46  unsigned nGPN(MGPAGainsPN.size());
47  for (unsigned iG(0); iG != nGPN; ++iG) {
48  int gain(MGPAGainsPN[iG]);
49  if (gain != 1 && gain != 16)
50  throw cms::Exception("InvalidConfiguration") << "PN MGPA gain";
51  repl["pngain"] = std::to_string(gain);
52  pnGainToME_[gain] = pnAmplitude.getIndex(repl);
53  }
54 
55  amplitudeThreshold_.resize(nG);
56  toleranceRMS_.resize(nG);
57 
58  std::vector<double> inAmplitudeThreshold(_params.getUntrackedParameter<std::vector<double> >("amplitudeThreshold"));
59  std::vector<double> inToleranceRMS(_params.getUntrackedParameter<std::vector<double> >("toleranceRMS"));
60 
61  for (std::map<int, unsigned>::iterator gainItr(gainToME_.begin()); gainItr != gainToME_.end(); ++gainItr) {
62  unsigned iME(gainItr->second);
63  unsigned iGain(0);
64  switch (gainItr->first) {
65  case 1:
66  iGain = 0;
67  break;
68  case 6:
69  iGain = 1;
70  break;
71  case 12:
72  iGain = 2;
73  break;
74  }
75 
76  amplitudeThreshold_[iME] = inAmplitudeThreshold[iGain];
77  toleranceRMS_[iME] = inToleranceRMS[iGain];
78  }
79 
80  PNAmplitudeThreshold_.resize(nGPN);
81  tolerancePNRMS_.resize(nGPN);
82 
83  std::vector<double> inPNAmplitudeThreshold(
84  _params.getUntrackedParameter<std::vector<double> >("PNAmplitudeThreshold"));
85  std::vector<double> inTolerancePNRMS(_params.getUntrackedParameter<std::vector<double> >("tolerancePNRMS"));
86 
87  for (std::map<int, unsigned>::iterator gainItr(pnGainToME_.begin()); gainItr != pnGainToME_.end(); ++gainItr) {
88  unsigned iME(gainItr->second);
89  unsigned iGain(0);
90  switch (gainItr->first) {
91  case 1:
92  iGain = 0;
93  break;
94  case 16:
95  iGain = 1;
96  break;
97  }
98 
99  PNAmplitudeThreshold_[iME] = inPNAmplitudeThreshold[iGain];
100  tolerancePNRMS_[iME] = inTolerancePNRMS[iGain];
101  }
102 
103  qualitySummaries_.insert("Quality");
104  qualitySummaries_.insert("QualitySummary");
105  qualitySummaries_.insert("PNQualitySummary");
106  }
107 
109  using namespace std;
110 
111  MESetMulti& meQuality(static_cast<MESetMulti&>(MEs_.at("Quality")));
112  MESetMulti& meAmplitudeRMS(static_cast<MESetMulti&>(MEs_.at("AmplitudeRMS")));
113  MESetMulti& meQualitySummary(static_cast<MESetMulti&>(MEs_.at("QualitySummary")));
114  MESetMulti& mePNQualitySummary(static_cast<MESetMulti&>(MEs_.at("PNQualitySummary")));
115 
116  MESetMulti const& sAmplitude(static_cast<MESetMulti const&>(sources_.at("Amplitude")));
117  MESetMulti const& sPNAmplitude(static_cast<MESetMulti const&>(sources_.at("PNAmplitude")));
118 
119  for (map<int, unsigned>::iterator gainItr(gainToME_.begin()); gainItr != gainToME_.end(); ++gainItr) {
120  meQuality.use(gainItr->second);
121  meQualitySummary.use(gainItr->second);
122  meAmplitudeRMS.use(gainItr->second);
123 
124  sAmplitude.use(gainItr->second);
125 
126  uint32_t mask(0);
127  switch (gainItr->first) {
128  case 1:
131  break;
132  case 6:
135  break;
136  case 12:
139  break;
140  default:
141  break;
142  }
143 
144  MESet::iterator qEnd(meQuality.end(GetElectronicsMap()));
145  MESet::iterator rItr(GetElectronicsMap(), meAmplitudeRMS);
146  MESet::const_iterator aItr(GetElectronicsMap(), sAmplitude);
147  for (MESet::iterator qItr(meQuality.beginChannel(GetElectronicsMap())); qItr != qEnd;
149  DetId id(qItr->getId());
150 
151  bool doMask(meQuality.maskMatches(id, mask, statusManager_, GetTrigTowerMap()));
152 
153  aItr = qItr;
154  rItr = qItr;
155 
156  float entries(aItr->getBinEntries());
157 
158  if (entries < minChannelEntries_) {
159  qItr->setBinContent(doMask ? kMUnknown : kUnknown);
160  continue;
161  }
162 
163  float amp(aItr->getBinContent());
164  float rms(aItr->getBinError() * sqrt(entries));
165 
166  rItr->setBinContent(rms);
167 
168  if (amp < amplitudeThreshold_[gainItr->second] || rms > toleranceRMS_[gainItr->second])
169  qItr->setBinContent(doMask ? kMBad : kBad);
170  else
171  qItr->setBinContent(doMask ? kMGood : kGood);
172  }
173 
174  towerAverage_(meQualitySummary, meQuality, 0.2);
175  }
176 
177  for (map<int, unsigned>::iterator gainItr(pnGainToME_.begin()); gainItr != pnGainToME_.end(); ++gainItr) {
178  mePNQualitySummary.use(gainItr->second);
179 
180  sPNAmplitude.use(gainItr->second);
181 
182  uint32_t mask(0);
183  switch (gainItr->first) {
184  case 1:
187  break;
188  case 16:
191  break;
192  default:
193  break;
194  }
195 
196  for (unsigned iDCC(0); iDCC < nDCC; ++iDCC) {
197  if (memDCCIndex(iDCC + 1) == unsigned(-1))
198  continue;
199 
200  for (unsigned iPN(0); iPN < 10; ++iPN) {
201  int subdet(0);
202  if (iDCC >= kEBmLow && iDCC <= kEBpHigh)
203  subdet = EcalBarrel;
204  else
205  subdet = EcalEndcap;
206 
207  EcalPnDiodeDetId id(subdet, iDCC + 1, iPN + 1);
208 
209  bool doMask(mePNQualitySummary.maskMatches(id, mask, statusManager_, GetTrigTowerMap()));
210 
211  float amp(sPNAmplitude.getBinContent(getEcalDQMSetupObjects(), id));
212  float entries(sPNAmplitude.getBinEntries(getEcalDQMSetupObjects(), id));
213  float rms(sPNAmplitude.getBinError(getEcalDQMSetupObjects(), id) * sqrt(entries));
214 
215  if (entries < minChannelEntries_) {
216  mePNQualitySummary.setBinContent(getEcalDQMSetupObjects(), id, doMask ? kMUnknown : kUnknown);
217  continue;
218  }
219 
220  if (amp < PNAmplitudeThreshold_[gainItr->second] || rms > tolerancePNRMS_[gainItr->second])
221  mePNQualitySummary.setBinContent(getEcalDQMSetupObjects(), id, doMask ? kMBad : kBad);
222  else
223  mePNQualitySummary.setBinContent(getEcalDQMSetupObjects(), id, doMask ? kMGood : kGood);
224  }
225  }
226  }
227  }
228 
230 } // namespace ecaldqm
unsigned memDCCIndex(unsigned)
void use(unsigned) const
Definition: MESetMulti.cc:133
void towerAverage_(MESet &, MESet const &, float)
static const int TESTPULSE_LOW_GAIN_RMS_ERROR
#define DEFINE_ECALDQM_WORKER(TYPE)
Definition: DQWorker.h:168
MESet & at(const std::string &key)
Definition: MESet.h:401
static constexpr int kGood
const_iterator & toNextChannel(EcalElectronicsMapping const *electronicsMap)
Definition: MESet.h:322
double getBinContent(EcalDQMSetupObjects const edso, DetId const &_id, int _bin=0) const override
Definition: MESetMulti.h:111
std::vector< float > PNAmplitudeThreshold_
static constexpr int kMUnknown
const_iterator beginChannel(EcalElectronicsMapping const *electronicsMap) const override
Definition: MESetMulti.h:168
static const int TESTPULSE_MIDDLE_GAIN_MEAN_ERROR
std::vector< float > tolerancePNRMS_
static constexpr int kUnknown
static std::string to_string(const XMLCh *ch)
static constexpr int kMBad
std::map< int, unsigned > gainToME_
static const int TESTPULSE_MIDDLE_GAIN_RMS_ERROR
std::set< std::string > qualitySummaries_
T sqrt(T t)
Definition: SSEVec.h:23
StatusManager const * statusManager_
void setBinContent(EcalDQMSetupObjects const edso, DetId const &_id, double _content) override
Definition: MESetMulti.h:48
static const int TESTPULSE_HIGH_GAIN_RMS_ERROR
static constexpr int kBad
MESetCollection sources_
double getBinError(EcalDQMSetupObjects const edso, DetId const &_id, int _bin=0) const override
Definition: MESetMulti.h:121
void producePlots(ProcessType) override
std::map< int, unsigned > pnGainToME_
EcalElectronicsMapping const * GetElectronicsMap()
Definition: DQWorker.cc:150
EcalDQMSetupObjects const getEcalDQMSetupObjects()
Definition: DQWorker.cc:170
Definition: DetId.h:17
unsigned getIndex(PathReplacements const &) const
Definition: MESetMulti.cc:140
static const int TESTPULSE_LOW_GAIN_MEAN_ERROR
static constexpr int nDCC
double getBinEntries(EcalDQMSetupObjects const edso, DetId const &_id, int _bin=0) const override
Definition: MESetMulti.h:131
MESetCollection MEs_
Definition: DQWorker.h:131
std::vector< float > toleranceRMS_
std::vector< float > amplitudeThreshold_
static constexpr int kMGood
EcalTrigTowerConstituentsMap const * GetTrigTowerMap()
Definition: DQWorker.cc:155
static const int TESTPULSE_HIGH_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
void setParams(edm::ParameterSet const &) override