CMS 3D CMS Logo

PresampleTask.cc
Go to the documentation of this file.
2 
4 
6 
8 
9 namespace ecaldqm {
11  : DQWorkerTask(), doPulseMaxCheck_(true), pulseMaxPosition_(0), nSamples_(0), mePedestalByLS(nullptr) {}
12 
14  doPulseMaxCheck_ = _params.getUntrackedParameter<bool>("doPulseMaxCheck");
15  pulseMaxPosition_ = _params.getUntrackedParameter<int>("pulseMaxPosition");
16  nSamples_ = _params.getUntrackedParameter<int>("nSamples");
17  }
18 
19  bool PresampleTask::filterRunType(short const* _runType) {
20  for (int iFED(0); iFED < nDCC; iFED++) {
21  if (_runType[iFED] == EcalDCCHeaderBlock::COSMIC || _runType[iFED] == EcalDCCHeaderBlock::MTCC ||
22  _runType[iFED] == EcalDCCHeaderBlock::COSMICS_GLOBAL ||
23  _runType[iFED] == EcalDCCHeaderBlock::PHYSICS_GLOBAL || _runType[iFED] == EcalDCCHeaderBlock::COSMICS_LOCAL ||
24  _runType[iFED] == EcalDCCHeaderBlock::PHYSICS_LOCAL)
25  return true;
26  }
27 
28  return false;
29  }
30 
32  // Fill separate MEs with only 10 LSs worth of stats
33  // Used to correctly fill Presample Trend plots:
34  // 1 pt:10 LS in Trend plots
35  mePedestalByLS = &MEs_.at("PedestalByLS");
36  if (timestamp_.iLumi % 10 == 0)
38  }
39 
40  template <typename DigiCollection>
42  MESet& mePedestal(MEs_.at("Pedestal")); // contains cumulative run stats => not suitable for Trend plots
43 
44  for (typename DigiCollection::const_iterator digiItr(_digis.begin()); digiItr != _digis.end(); ++digiItr) {
45  DetId id(digiItr->id());
46 
47  // EcalDataFrame is not a derived class of edm::DataFrame, but can take edm::DataFrame in the constructor
48  EcalDataFrame dataFrame(*digiItr);
49 
50  // Check that the digi pulse maximum occurs on the 6th sample
51  // For cosmics: disable this check to preserve statistics
52  if (doPulseMaxCheck_) {
53  bool gainSwitch(false);
54  int iMax(-1);
55  int maxADC(0);
56  for (int iSample(0); iSample < EcalDataFrame::MAXSAMPLES; ++iSample) {
57  int adc(dataFrame.sample(iSample).adc());
58  if (adc > maxADC) {
59  iMax = iSample;
60  maxADC = adc;
61  }
62  if (iSample < nSamples_ && dataFrame.sample(iSample).gainId() != 1) {
63  gainSwitch = true;
64  break;
65  }
66  } // iSample
67  if (iMax != pulseMaxPosition_ || gainSwitch)
68  continue;
69  } // PulseMaxCheck
70 
71  for (int iSample(0); iSample < nSamples_; ++iSample) {
72  mePedestal.fill(id, double(dataFrame.sample(iSample).adc()));
73  mePedestalByLS->fill(id, double(dataFrame.sample(iSample).adc()));
74  }
75 
76  } // _digis loop
77  } // runOnDigis
78 
80 } // namespace ecaldqm
T getUntrackedParameter(std::string const &, T const &) const
#define DEFINE_ECALDQM_WORKER(TYPE)
Definition: DQWorker.h:112
edm::LuminosityBlockNumber_t iLumi
Definition: DQWorker.h:35
void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) override
#define nullptr
void runOnDigis(DigiCollection const &)
void setParams(edm::ParameterSet const &) override
virtual void reset(double=0., double=0., double=0.)
Definition: MESet.cc:98
constexpr int adc(sample_type sample)
get the ADC sample (12 bits)
Definition: DetId.h:17
Timestamp timestamp_
Definition: DQWorker.h:81
virtual void fill(DetId const &, double=1., double=1., double=1.)
Definition: MESet.h:46
MESetCollection MEs_
Definition: DQWorker.h:78
static constexpr int MAXSAMPLES
Definition: EcalDataFrame.h:48
bool filterRunType(short const *) override