CMS 3D CMS Logo

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