CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
ecaldqm::LaserClient Class Reference

#include <LaserClient.h>

Inheritance diagram for ecaldqm::LaserClient:
ecaldqm::DQWorkerClient ecaldqm::DQWorker

Public Member Functions

 LaserClient ()
 
void producePlots (ProcessType) override
 
 ~LaserClient () override
 
- Public Member Functions inherited from ecaldqm::DQWorkerClient
void bookMEs (DQMStore::IBooker &) override
 
 DQWorkerClient ()
 
void endLuminosityBlock (edm::LuminosityBlock const &, edm::EventSetup const &) override
 
void releaseMEs () override
 
void releaseSource ()
 
virtual void resetMEs ()
 
void resetPerLumi ()
 
bool retrieveSource (DQMStore::IGetter &, ProcessType)
 
bool runsOn (ProcessType _type) const
 
void setStatusManager (StatusManager const &_manager)
 
 ~DQWorkerClient () override
 
- Public Member Functions inherited from ecaldqm::DQWorker
virtual void beginLuminosityBlock (edm::LuminosityBlock const &, edm::EventSetup const &)
 
virtual void beginRun (edm::Run const &, edm::EventSetup const &)
 
 DQWorker ()
 
virtual void endRun (edm::Run const &, edm::EventSetup const &)
 
const EcalDQMSetupObjects getEcalDQMSetupObjects ()
 
const EcalElectronicsMappingGetElectronicsMap ()
 
const CaloGeometryGetGeometry ()
 
const std::string & getName () const
 
const CaloTopologyGetTopology ()
 
const EcalTrigTowerConstituentsMapGetTrigTowerMap ()
 
bool onlineMode () const
 
void setEventNumber (edm::EventNumber_t _e)
 
void setLumiNumber (edm::LuminosityBlockNumber_t _l)
 
void setRunNumber (edm::RunNumber_t _r)
 
void setSetupObjects (edm::EventSetup const &)
 
void setTime (time_t _t)
 
virtual ~DQWorker () noexcept(false)
 

Private Member Functions

void setParams (edm::ParameterSet const &) override
 

Private Attributes

std::vector< float > expectedAmplitude_
 
std::vector< float > expectedPNAmplitude_
 
std::vector< float > expectedTiming_
 
float forwardFactor_
 
int minChannelEntries_
 
float toleranceAmplitudeHi_
 
float toleranceAmplitudeLo_
 
float toleranceAmpRMSRatio_
 
float tolerancePNAmp_
 
float tolerancePNRMSRatio_
 
float toleranceTiming_
 
float toleranceTimRMS_
 
std::map< int, unsigned > wlToME_
 

Additional Inherited Members

- Public Types inherited from ecaldqm::DQWorkerClient
enum  ProcessType { kLumi, kJob, nProcessType }
 
enum  Quality {
  kBad = 0, kGood = 1, kUnknown = 2, kMBad = 3,
  kMGood = 4, kMUnknown = 5
}
 
- Static Public Member Functions inherited from ecaldqm::DQWorkerClient
static void fillDescriptions (edm::ParameterSetDescription &)
 
- Static Public Member Functions inherited from ecaldqm::DQWorker
static void fillDescriptions (edm::ParameterSetDescription &_desc)
 
- Protected Types inherited from ecaldqm::DQWorker
typedef dqm::legacy::DQMStore DQMStore
 
typedef dqm::legacy::MonitorElement MonitorElement
 
- Protected Member Functions inherited from ecaldqm::DQWorkerClient
void setME (edm::ParameterSet const &_ps) final
 
void setSource (edm::ParameterSet const &) override
 
void towerAverage_ (MESet &, MESet const &, float)
 
bool using_ (std::string const &_name, ProcessType _type=kJob) const
 
- Protected Member Functions inherited from ecaldqm::DQWorker
void initialize (std::string const &_name, edm::ParameterSet const &)
 
void print_ (std::string const &, int=0) const
 
void setVerbosity (int _verbosity)
 
- Protected Attributes inherited from ecaldqm::DQWorkerClient
bool hasLumiPlots_
 
std::set< std::string > qualitySummaries_
 
MESetCollection sources_
 
StatusManager const * statusManager_
 
- Protected Attributes inherited from ecaldqm::DQWorker
bool booked_
 
MESetCollection MEs_
 
std::string name_
 
bool onlineMode_
 
Timestamp timestamp_
 
int verbosity_
 
bool willConvertToEDM_
 

Detailed Description

Definition at line 7 of file LaserClient.h.

Constructor & Destructor Documentation

◆ LaserClient()

ecaldqm::LaserClient::LaserClient ( )

Definition at line 15 of file LaserClient.cc.

◆ ~LaserClient()

ecaldqm::LaserClient::~LaserClient ( )
inlineoverride

Definition at line 10 of file LaserClient.h.

10 {}

Member Function Documentation

◆ producePlots()

void ecaldqm::LaserClient::producePlots ( ProcessType  )
overridevirtual

Implements ecaldqm::DQWorkerClient.

Definition at line 82 of file LaserClient.cc.

82  {
85 
86  MESetMulti& meQuality(static_cast<MESetMulti&>(MEs_.at("Quality")));
87  MESetMulti& meQualitySummary(static_cast<MESetMulti&>(MEs_.at("QualitySummary")));
88  MESetMulti& meAmplitudeMean(static_cast<MESetMulti&>(MEs_.at("AmplitudeMean")));
89  MESetMulti& meAmplitudeRMS(static_cast<MESetMulti&>(MEs_.at("AmplitudeRMS")));
90  MESetMulti& meTimingMean(static_cast<MESetMulti&>(MEs_.at("TimingMean")));
91  MESetMulti& meTimingRMSMap(static_cast<MESetMulti&>(MEs_.at("TimingRMSMap")));
92  MESetMulti& meTimingRMS(static_cast<MESetMulti&>(MEs_.at("TimingRMS")));
93  MESetMulti& mePNQualitySummary(static_cast<MESetMulti&>(MEs_.at("PNQualitySummary")));
94 
95  MESetMulti const& sAmplitude(static_cast<MESetMulti const&>(sources_.at("Amplitude")));
96  MESetMulti const& sTiming(static_cast<MESetMulti const&>(sources_.at("Timing")));
97  MESetMulti const& sPNAmplitude(static_cast<MESetMulti const&>(sources_.at("PNAmplitude")));
98  MESet const& sCalibStatus(static_cast<MESet const&>(sources_.at("CalibStatus")));
99 
100  for (std::map<int, unsigned>::iterator wlItr(wlToME_.begin()); wlItr != wlToME_.end(); ++wlItr) {
101  meQuality.use(wlItr->second);
102  meQualitySummary.use(wlItr->second);
103  meAmplitudeMean.use(wlItr->second);
104  meAmplitudeRMS.use(wlItr->second);
105  meTimingMean.use(wlItr->second);
106  meTimingRMSMap.use(wlItr->second);
107  meTimingRMS.use(wlItr->second);
108  mePNQualitySummary.use(wlItr->second);
109 
110  sAmplitude.use(wlItr->second);
111  sTiming.use(wlItr->second);
112  sPNAmplitude.use(wlItr->second);
113 
114  MESet::iterator qEnd(meQuality.end(GetElectronicsMap()));
115 
116  MESet::const_iterator tItr(GetElectronicsMap(), sTiming);
117  MESet::const_iterator aItr(GetElectronicsMap(), sAmplitude);
118 
119  int wl(wlItr->first - 1);
120  bool enabled(wl < 0 ? false : sCalibStatus.getBinContent(getEcalDQMSetupObjects(), wl) > 0 ? true : false);
121  for (MESet::iterator qItr(meQuality.beginChannel(GetElectronicsMap())); qItr != qEnd;
122  qItr.toNextChannel(GetElectronicsMap())) {
123  DetId id(qItr->getId());
124 
125  bool doMask(meQuality.maskMatches(id, mask, statusManager_, GetTrigTowerMap()));
126 
127  aItr = qItr;
128 
129  float aEntries(aItr->getBinEntries());
130 
131  if (aEntries < minChannelEntries_) {
132  qItr->setBinContent(enabled ? (doMask ? kMUnknown : kUnknown) : kMUnknown);
133  continue;
134  }
135 
136  float aMean(aItr->getBinContent());
137  float aRms(aItr->getBinError() * sqrt(aEntries));
138 
139  meAmplitudeMean.fill(getEcalDQMSetupObjects(), id, aMean);
140  meAmplitudeRMS.setBinContent(getEcalDQMSetupObjects(), id, aRms);
141 
142  tItr = qItr;
143 
144  float tEntries(tItr->getBinEntries());
145 
146  if (tEntries < minChannelEntries_)
147  continue;
148 
149  float tMean(tItr->getBinContent());
150  float tRms(tItr->getBinError() * sqrt(tEntries));
151 
152  meTimingMean.fill(getEcalDQMSetupObjects(), id, tMean);
153  meTimingRMS.fill(getEcalDQMSetupObjects(), id, tRms);
154  meTimingRMSMap.setBinContent(getEcalDQMSetupObjects(), id, tRms);
155 
156  float intensity(aMean / expectedAmplitude_[wlItr->second]);
157  if (isForward(id))
158  intensity /= forwardFactor_;
159 
160  if (intensity < toleranceAmplitudeLo_ || intensity > toleranceAmplitudeHi_ ||
161  aRms > aMean * toleranceAmpRMSRatio_ ||
162  std::abs(tMean - expectedTiming_[wlItr->second]) > toleranceTiming_ /*|| tRms > toleranceTimRMS_*/)
163  qItr->setBinContent(doMask ? kMBad : kBad);
164  else
165  qItr->setBinContent(doMask ? kMGood : kGood);
166  }
167 
168  towerAverage_(meQualitySummary, meQuality, 0.2);
169 
170  for (unsigned iDCC(0); iDCC < nDCC; ++iDCC) {
171  if (memDCCIndex(iDCC + 1) == unsigned(-1))
172  continue;
173  int subdet(0);
174  if (iDCC >= kEBmLow && iDCC <= kEBpHigh)
175  subdet = EcalBarrel;
176  else
177  subdet = EcalEndcap;
178 
179  for (unsigned iPN(0); iPN < 10; ++iPN) {
180  EcalPnDiodeDetId id(subdet, iDCC + 1, iPN + 1);
181 
182  bool doMask(mePNQualitySummary.maskMatches(id, mask, statusManager_, GetTrigTowerMap()));
183 
184  float pEntries(sPNAmplitude.getBinEntries(getEcalDQMSetupObjects(), id));
185 
186  if (pEntries < minChannelEntries_) {
187  mePNQualitySummary.setBinContent(getEcalDQMSetupObjects(), id, doMask ? kMUnknown : kUnknown);
188  continue;
189  }
190 
191  float pMean(sPNAmplitude.getBinContent(getEcalDQMSetupObjects(), id));
192  float pRms(sPNAmplitude.getBinError(getEcalDQMSetupObjects(), id) * sqrt(pEntries));
193  float intensity(pMean / expectedPNAmplitude_[wlItr->second]);
194 
195  if (intensity < tolerancePNAmp_ || pRms > pMean * tolerancePNRMSRatio_)
196  mePNQualitySummary.setBinContent(getEcalDQMSetupObjects(), id, doMask ? kMBad : kBad);
197  else
198  mePNQualitySummary.setBinContent(getEcalDQMSetupObjects(), id, doMask ? kMGood : kGood);
199  }
200  }
201  }
202  }

References funct::abs(), ecaldqm::MESetCollection::at(), ecaldqm::MESetMulti::beginChannel(), EcalBarrel, EcalEndcap, pixel_dqm_sourceclient-live_cfg::enabled, ecaldqm::MESetMulti::end(), expectedAmplitude_, expectedPNAmplitude_, expectedTiming_, ecaldqm::MESetMulti::fill(), forwardFactor_, ecaldqm::MESet::getBinContent(), ecaldqm::MESetMulti::getBinContent(), ecaldqm::MESetMulti::getBinEntries(), ecaldqm::MESetMulti::getBinError(), ecaldqm::DQWorker::getEcalDQMSetupObjects(), ecaldqm::DQWorker::GetElectronicsMap(), ecaldqm::DQWorker::GetTrigTowerMap(), triggerObjects_cff::id, ecaldqm::isForward(), ecaldqm::DQWorkerClient::kBad, ecaldqm::kEBmLow, ecaldqm::kEBpHigh, ecaldqm::DQWorkerClient::kGood, ecaldqm::DQWorkerClient::kMBad, ecaldqm::DQWorkerClient::kMGood, ecaldqm::DQWorkerClient::kMUnknown, ecaldqm::DQWorkerClient::kUnknown, EcalDQMStatusHelper::LASER_MEAN_ERROR, EcalDQMStatusHelper::LASER_RMS_ERROR, EcalDQMStatusHelper::LASER_TIMING_MEAN_ERROR, EcalDQMStatusHelper::LASER_TIMING_RMS_ERROR, ecaldqm::MESetMulti::maskMatches(), ecaldqm::memDCCIndex(), ecaldqm::DQWorker::MEs_, minChannelEntries_, ecaldqm::nDCC, ecaldqm::MESetMulti::setBinContent(), ecaldqm::DQWorkerClient::sources_, mathSSE::sqrt(), ecaldqm::DQWorkerClient::statusManager_, toleranceAmplitudeHi_, toleranceAmpRMSRatio_, tolerancePNRMSRatio_, toleranceTiming_, ecaldqm::MESet::iterator::toNextChannel(), ecaldqm::DQWorkerClient::towerAverage_(), funct::true, ecaldqm::MESetMulti::use(), LaserClient_cfi::wl, and wlToME_.

◆ setParams()

void ecaldqm::LaserClient::setParams ( edm::ParameterSet const &  _params)
overrideprivatevirtual

Reimplemented from ecaldqm::DQWorker.

Definition at line 31 of file LaserClient.cc.

31  {
32  minChannelEntries_ = _params.getUntrackedParameter<int>("minChannelEntries");
33  toleranceAmplitudeLo_ = _params.getUntrackedParameter<double>("toleranceAmplitudeLo");
34  toleranceAmplitudeHi_ = _params.getUntrackedParameter<double>("toleranceAmplitudeHi");
35  toleranceAmpRMSRatio_ = _params.getUntrackedParameter<double>("toleranceAmpRMSRatio");
36  toleranceTiming_ = _params.getUntrackedParameter<double>("toleranceTiming");
37  toleranceTimRMS_ = _params.getUntrackedParameter<double>("toleranceTimRMS");
38  tolerancePNAmp_ = _params.getUntrackedParameter<double>("tolerancePNAmp");
39  tolerancePNRMSRatio_ = _params.getUntrackedParameter<double>("tolerancePNRMSRatio");
40  forwardFactor_ = _params.getUntrackedParameter<double>("forwardFactor");
41 
42  std::vector<int> laserWavelengths(_params.getUntrackedParameter<std::vector<int> >("laserWavelengths"));
43 
44  // wavelengths are not necessarily ordered
45  // create a map wl -> MESet index
46  // using Amplitude here but any multi-wavelength plot is fine
47 
49 
50  MESetMulti const& amplitude(static_cast<MESetMulti const&>(sources_.at("Amplitude")));
51  unsigned nWL(laserWavelengths.size());
52  for (unsigned iWL(0); iWL != nWL; ++iWL) {
53  int wl(laserWavelengths[iWL]);
54  if (wl <= 0 || wl >= 5)
55  throw cms::Exception("InvalidConfiguration") << "Laser Wavelength";
56  repl["wl"] = std::to_string(wl);
57  wlToME_[wl] = amplitude.getIndex(repl);
58  }
59 
60  expectedAmplitude_.resize(nWL);
61  expectedTiming_.resize(nWL);
62  expectedPNAmplitude_.resize(nWL);
63 
64  std::vector<double> inExpectedAmplitude(_params.getUntrackedParameter<std::vector<double> >("expectedAmplitude"));
65  std::vector<double> inExpectedTiming(_params.getUntrackedParameter<std::vector<double> >("expectedTiming"));
66  std::vector<double> inExpectedPNAmplitude(
67  _params.getUntrackedParameter<std::vector<double> >("expectedPNAmplitude"));
68 
69  for (std::map<int, unsigned>::iterator wlItr(wlToME_.begin()); wlItr != wlToME_.end(); ++wlItr) {
70  unsigned iME(wlItr->second);
71  int iWL(wlItr->first - 1);
72  expectedAmplitude_[iME] = inExpectedAmplitude[iWL];
73  expectedTiming_[iME] = inExpectedTiming[iWL];
74  expectedPNAmplitude_[iME] = inExpectedPNAmplitude[iWL];
75  }
76 
77  qualitySummaries_.insert("Quality");
78  qualitySummaries_.insert("QualitySummary");
79  qualitySummaries_.insert("PNQualitySummary");
80  }

References l1extraParticles_cfi::_params, CustomPhysics_cfi::amplitude, ecaldqm::MESetCollection::at(), Exception, expectedAmplitude_, expectedPNAmplitude_, expectedTiming_, forwardFactor_, CalibrationSummaryClient_cfi::laserWavelengths, minChannelEntries_, ecaldqm::DQWorkerClient::qualitySummaries_, ecaldqm::DQWorkerClient::sources_, toleranceAmplitudeHi_, toleranceAmplitudeLo_, toleranceAmpRMSRatio_, tolerancePNAmp_, tolerancePNRMSRatio_, toleranceTiming_, toleranceTimRMS_, LaserClient_cfi::wl, and wlToME_.

Member Data Documentation

◆ expectedAmplitude_

std::vector<float> ecaldqm::LaserClient::expectedAmplitude_
private

Definition at line 20 of file LaserClient.h.

Referenced by producePlots(), and setParams().

◆ expectedPNAmplitude_

std::vector<float> ecaldqm::LaserClient::expectedPNAmplitude_
private

Definition at line 27 of file LaserClient.h.

Referenced by producePlots(), and setParams().

◆ expectedTiming_

std::vector<float> ecaldqm::LaserClient::expectedTiming_
private

Definition at line 24 of file LaserClient.h.

Referenced by producePlots(), and setParams().

◆ forwardFactor_

float ecaldqm::LaserClient::forwardFactor_
private

Definition at line 30 of file LaserClient.h.

Referenced by producePlots(), and setParams().

◆ minChannelEntries_

int ecaldqm::LaserClient::minChannelEntries_
private

Definition at line 19 of file LaserClient.h.

Referenced by producePlots(), and setParams().

◆ toleranceAmplitudeHi_

float ecaldqm::LaserClient::toleranceAmplitudeHi_
private

Definition at line 22 of file LaserClient.h.

Referenced by producePlots(), and setParams().

◆ toleranceAmplitudeLo_

float ecaldqm::LaserClient::toleranceAmplitudeLo_
private

Definition at line 21 of file LaserClient.h.

Referenced by setParams().

◆ toleranceAmpRMSRatio_

float ecaldqm::LaserClient::toleranceAmpRMSRatio_
private

Definition at line 23 of file LaserClient.h.

Referenced by producePlots(), and setParams().

◆ tolerancePNAmp_

float ecaldqm::LaserClient::tolerancePNAmp_
private

Definition at line 28 of file LaserClient.h.

Referenced by setParams().

◆ tolerancePNRMSRatio_

float ecaldqm::LaserClient::tolerancePNRMSRatio_
private

Definition at line 29 of file LaserClient.h.

Referenced by producePlots(), and setParams().

◆ toleranceTiming_

float ecaldqm::LaserClient::toleranceTiming_
private

Definition at line 25 of file LaserClient.h.

Referenced by producePlots(), and setParams().

◆ toleranceTimRMS_

float ecaldqm::LaserClient::toleranceTimRMS_
private

Definition at line 26 of file LaserClient.h.

Referenced by setParams().

◆ wlToME_

std::map<int, unsigned> ecaldqm::LaserClient::wlToME_
private

Definition at line 17 of file LaserClient.h.

Referenced by producePlots(), and setParams().

CustomPhysics_cfi.amplitude
amplitude
Definition: CustomPhysics_cfi.py:12
ecaldqm::LaserClient::wlToME_
std::map< int, unsigned > wlToME_
Definition: LaserClient.h:17
ecaldqm::LaserClient::forwardFactor_
float forwardFactor_
Definition: LaserClient.h:30
EcalDQMStatusHelper::LASER_MEAN_ERROR
static const int LASER_MEAN_ERROR
Definition: EcalDQMStatusHelper.h:34
ecaldqm::LaserClient::toleranceAmplitudeHi_
float toleranceAmplitudeHi_
Definition: LaserClient.h:22
LaserClient_cfi.wl
wl
Definition: LaserClient_cfi.py:46
ecaldqm::LaserClient::toleranceTiming_
float toleranceTiming_
Definition: LaserClient.h:25
ecaldqm::MESet::PathReplacements
std::map< std::string, std::string > PathReplacements
Definition: MESet.h:46
ecaldqm::LaserClient::expectedTiming_
std::vector< float > expectedTiming_
Definition: LaserClient.h:24
l1extraParticles_cfi._params
_params
Definition: l1extraParticles_cfi.py:29
ecaldqm::memDCCIndex
unsigned memDCCIndex(unsigned)
Definition: EcalDQMCommonUtils.cc:46
ecaldqm::DQWorkerClient::statusManager_
StatusManager const * statusManager_
Definition: DQWorkerClient.h:60
ecaldqm::LaserClient::toleranceAmpRMSRatio_
float toleranceAmpRMSRatio_
Definition: LaserClient.h:23
EcalPnDiodeDetId
Definition: EcalPnDiodeDetId.h:22
ecaldqm::LaserClient::toleranceAmplitudeLo_
float toleranceAmplitudeLo_
Definition: LaserClient.h:21
EcalDQMStatusHelper::LASER_TIMING_MEAN_ERROR
static const int LASER_TIMING_MEAN_ERROR
Definition: EcalDQMStatusHelper.h:36
ecaldqm::isForward
bool isForward(DetId const &)
Definition: EcalDQMCommonUtils.cc:243
ecaldqm::DQWorker::GetElectronicsMap
const EcalElectronicsMapping * GetElectronicsMap()
Definition: DQWorker.cc:104
EcalBarrel
Definition: EcalSubdetector.h:10
DetId
Definition: DetId.h:17
ecaldqm::DQWorkerClient::qualitySummaries_
std::set< std::string > qualitySummaries_
Definition: DQWorkerClient.h:56
EcalDQMStatusHelper::LASER_TIMING_RMS_ERROR
static const int LASER_TIMING_RMS_ERROR
Definition: EcalDQMStatusHelper.h:37
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
ecaldqm::kEBpHigh
Definition: EcalDQMCommonUtils.h:87
ecaldqm::LaserClient::tolerancePNRMSRatio_
float tolerancePNRMSRatio_
Definition: LaserClient.h:29
ecaldqm::DQWorker::GetTrigTowerMap
const EcalTrigTowerConstituentsMap * GetTrigTowerMap()
Definition: DQWorker.cc:110
ecaldqm::DQWorker::MEs_
MESetCollection MEs_
Definition: DQWorker.h:104
ecaldqm::DQWorkerClient::kMBad
Definition: DQWorkerClient.h:37
EcalEndcap
Definition: EcalSubdetector.h:10
ecaldqm::DQWorkerClient::kUnknown
Definition: DQWorkerClient.h:37
ecaldqm::DQWorkerClient::kBad
Definition: DQWorkerClient.h:37
CalibrationSummaryClient_cfi.laserWavelengths
laserWavelengths
Definition: CalibrationSummaryClient_cfi.py:15
funct::true
true
Definition: Factorize.h:173
ecaldqm::DQWorkerClient::sources_
MESetCollection sources_
Definition: DQWorkerClient.h:55
ecaldqm::DQWorkerClient::kMGood
Definition: DQWorkerClient.h:37
ecaldqm::LaserClient::minChannelEntries_
int minChannelEntries_
Definition: LaserClient.h:19
ecaldqm::MESetCollection::at
MESet & at(const std::string &key)
Definition: MESet.h:399
ecaldqm::DQWorkerClient::kMUnknown
Definition: DQWorkerClient.h:37
ecaldqm::LaserClient::toleranceTimRMS_
float toleranceTimRMS_
Definition: LaserClient.h:26
ecaldqm::kEBmLow
Definition: EcalDQMCommonUtils.h:84
ecaldqm::DQWorkerClient::DQWorkerClient
DQWorkerClient()
Definition: DQWorkerClient.cc:17
ecaldqm::DQWorker::getEcalDQMSetupObjects
const EcalDQMSetupObjects getEcalDQMSetupObjects()
Definition: DQWorker.cc:128
triggerObjects_cff.id
id
Definition: triggerObjects_cff.py:29
ecaldqm::DQWorkerClient::towerAverage_
void towerAverage_(MESet &, MESet const &, float)
Definition: DQWorkerClient.cc:159
ecaldqm::DQWorkerClient::kGood
Definition: DQWorkerClient.h:37
Exception
Definition: hltDiff.cc:245
ecaldqm::nDCC
Definition: EcalDQMCommonUtils.h:91
ecaldqm::LaserClient::expectedPNAmplitude_
std::vector< float > expectedPNAmplitude_
Definition: LaserClient.h:27
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ecaldqm::LaserClient::expectedAmplitude_
std::vector< float > expectedAmplitude_
Definition: LaserClient.h:20
EcalDQMStatusHelper::LASER_RMS_ERROR
static const int LASER_RMS_ERROR
Definition: EcalDQMStatusHelper.h:35
ecaldqm::LaserClient::tolerancePNAmp_
float tolerancePNAmp_
Definition: LaserClient.h:28
pixel_dqm_sourceclient-live_cfg.enabled
enabled
Definition: pixel_dqm_sourceclient-live_cfg.py:136