test
CMS 3D CMS Logo

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