CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PresampleTask.cc
Go to the documentation of this file.
3 
5 
7 
9 
10 namespace ecaldqm {
12  : DQWorkerTask(), doPulseMaxCheck_(true), pulseMaxPosition_(0), nSamples_(0), mePedestalByLS(nullptr) {}
13 
15  doPulseMaxCheck_ = _params.getUntrackedParameter<bool>("doPulseMaxCheck");
16  pulseMaxPosition_ = _params.getUntrackedParameter<int>("pulseMaxPosition");
17  nSamples_ = _params.getUntrackedParameter<int>("nSamples");
18  }
20 
21  bool PresampleTask::filterRunType(short const* _runType) {
22  for (int iFED(0); iFED < nDCC; iFED++) {
23  if (_runType[iFED] == EcalDCCHeaderBlock::COSMIC || _runType[iFED] == EcalDCCHeaderBlock::MTCC ||
24  _runType[iFED] == EcalDCCHeaderBlock::COSMICS_GLOBAL ||
25  _runType[iFED] == EcalDCCHeaderBlock::PHYSICS_GLOBAL || _runType[iFED] == EcalDCCHeaderBlock::COSMICS_LOCAL ||
26  _runType[iFED] == EcalDCCHeaderBlock::PHYSICS_LOCAL)
27  return true;
28  }
29 
30  return false;
31  }
32 
33  void PresampleTask::beginRun(edm::Run const&, edm::EventSetup const& _es) { FillPedestal = true; }
34 
36  edm::EventSetup const& _es,
37  bool const& ByLumiResetSwitch,
38  bool&) {
39  if (ByLumiResetSwitch) {
40  // Fill separate MEs with only 10 LSs worth of stats
41  // Used to correctly fill Presample Trend plots:
42  // 1 pt:10 LS in Trend plots
43  mePedestalByLS = &MEs_.at("PedestalByLS");
44  if (timestamp_.iLumi % 10 == 0)
46  }
47 
48  MESet& mePedestalProjEtaG1(MEs_.at("PedestalProjEtaG1"));
49  MESet& mePedestalProjEtaG6(MEs_.at("PedestalProjEtaG6"));
50  MESet& mePedestalProjEtaG12(MEs_.at("PedestalProjEtaG12"));
51 
52  if (FillPedestal) {
53  const EcalPedestals* myped = &_es.getData(Pedtoken_);
54 
55  for (int i = 0; i < EBDetId::kSizeForDenseIndexing; i++) {
57  continue;
59  EcalPedestals::const_iterator it = myped->find(ebid.rawId());
60  if (it != myped->end()) {
61  mePedestalProjEtaG1.fill(getEcalDQMSetupObjects(), ebid, (*it).rms_x1);
62  mePedestalProjEtaG6.fill(getEcalDQMSetupObjects(), ebid, (*it).rms_x6);
63  mePedestalProjEtaG12.fill(getEcalDQMSetupObjects(), ebid, (*it).rms_x12);
64  }
65  }
66  for (int i = 0; i < EEDetId::kSizeForDenseIndexing; i++) {
68  continue;
70  EcalPedestals::const_iterator it = myped->find(eeid.rawId());
71  if (it != myped->end()) {
72  mePedestalProjEtaG1.fill(getEcalDQMSetupObjects(), eeid, (*it).rms_x1);
73  mePedestalProjEtaG6.fill(getEcalDQMSetupObjects(), eeid, (*it).rms_x6);
74  mePedestalProjEtaG12.fill(getEcalDQMSetupObjects(), eeid, (*it).rms_x12);
75  }
76  }
77 
78  FillPedestal = false;
79  }
80  }
81 
82  template <typename DigiCollection>
83  void PresampleTask::runOnDigis(DigiCollection const& _digis) {
84  MESet& mePedestal(MEs_.at("Pedestal")); // contains cumulative run stats => not suitable for Trend plots
85 
86  for (typename DigiCollection::const_iterator digiItr(_digis.begin()); digiItr != _digis.end(); ++digiItr) {
87  DetId id(digiItr->id());
88 
89  // EcalDataFrame is not a derived class of edm::DataFrame, but can take edm::DataFrame in the constructor
90  EcalDataFrame dataFrame(*digiItr);
91 
92  // Check that the digi pulse maximum occurs on the 6th sample
93  // For cosmics: disable this check to preserve statistics
94  if (doPulseMaxCheck_) {
95  bool gainSwitch(false);
96  int iMax(-1);
97  int maxADC(0);
98  for (int iSample(0); iSample < EcalDataFrame::MAXSAMPLES; ++iSample) {
99  int adc(dataFrame.sample(iSample).adc());
100  if (adc > maxADC) {
101  iMax = iSample;
102  maxADC = adc;
103  }
104  if (iSample < nSamples_ && dataFrame.sample(iSample).gainId() != 1) {
105  gainSwitch = true;
106  break;
107  }
108  } // iSample
109  if (iMax != pulseMaxPosition_ || gainSwitch)
110  continue;
111  } // PulseMaxCheck
112 
113  for (int iSample(0); iSample < nSamples_; ++iSample) {
114  mePedestal.fill(getEcalDQMSetupObjects(), id, double(dataFrame.sample(iSample).adc()));
115  mePedestalByLS->fill(getEcalDQMSetupObjects(), id, double(dataFrame.sample(iSample).adc()));
116  }
117 
118  } // _digis loop
119  } // runOnDigis
120 
122 } // namespace ecaldqm
T getUntrackedParameter(std::string const &, T const &) const
#define DEFINE_ECALDQM_WORKER(TYPE)
Definition: DQWorker.h:162
uint16_t *__restrict__ id
edm::LuminosityBlockNumber_t iLumi
Definition: DQWorker.h:48
MESet & at(const std::string &key)
Definition: MESet.h:399
void beginRun(edm::Run const &, edm::EventSetup const &) override
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
void beginEvent(edm::Event const &, edm::EventSetup const &, bool const &, bool &) override
void runOnDigis(DigiCollection const &)
static EEDetId unhashIndex(int hi)
Definition: EEDetId.cc:65
void setParams(edm::ParameterSet const &) override
bool getData(T &iHolder) const
Definition: EventSetup.h:128
void setTokens(edm::ConsumesCollector &) override
virtual void fill(EcalDQMSetupObjects const, DetId const &, double=1., double=1., double=1.)
Definition: MESet.h:74
EcalDQMSetupObjects const getEcalDQMSetupObjects()
Definition: DQWorker.cc:142
Definition: DetId.h:17
Timestamp timestamp_
Definition: DQWorker.h:128
std::vector< Item >::const_iterator const_iterator
MESetCollection MEs_
Definition: DQWorker.h:125
static EBDetId unhashIndex(int hi)
get a DetId from a compact index for arrays
Definition: EBDetId.h:110
static bool validDenseIndex(uint32_t din)
Definition: EEDetId.h:213
const_iterator find(uint32_t rawId) const
edm::ESGetToken< EcalPedestals, EcalPedestalsRcd > Pedtoken_
Definition: PresampleTask.h:31
EcalElectronicsMapping const * GetElectronicsMap()
Definition: DQWorker.cc:118
const_iterator end() const
static constexpr int MAXSAMPLES
Definition: EcalDataFrame.h:48
virtual void reset(EcalElectronicsMapping const *, double=0., double=0., double=0.)
Definition: MESet.cc:98
bool filterRunType(short const *) override
static bool validDenseIndex(uint32_t din)
Definition: EBDetId.h:105
Definition: Run.h:45
uint16_t *__restrict__ uint16_t const *__restrict__ adc