CMS 3D CMS Logo

ECALpedestalPCLworker.cc
Go to the documentation of this file.
1 
9 #include <iostream>
10 #include <sstream>
11 
13 
14 {
15  edm::InputTag digiTagEB = iConfig.getParameter<edm::InputTag>("BarrelDigis");
16  edm::InputTag digiTagEE = iConfig.getParameter<edm::InputTag>("EndcapDigis");
17 
18  digiTokenEB_ = consumes<EBDigiCollection>(digiTagEB);
19  digiTokenEE_ = consumes<EEDigiCollection>(digiTagEE);
20 
21  pedestalSamples_ = iConfig.getParameter<uint32_t>("pedestalSamples");
22  checkSignal_ = iConfig.getParameter<bool>("checkSignal");
23  sThresholdEB_ = iConfig.getParameter<uint32_t>("sThresholdEB");
24  sThresholdEE_ = iConfig.getParameter<uint32_t>("sThresholdEE");
25 
26  dynamicBooking_ = iConfig.getParameter<bool>("dynamicBooking");
27  fixedBookingCenterBin_ = iConfig.getParameter<int>("fixedBookingCenterBin");
28  nBins_ = iConfig.getParameter<int>("nBins");
29  dqmDir_ = iConfig.getParameter<std::string>("dqmDir");
30 
31  edm::InputTag tcdsRecord = iConfig.getParameter<edm::InputTag>("tcdsRecord");
32  tcdsToken_ = consumes<TCDSRecord>(tcdsRecord);
33  requireStableBeam_ = iConfig.getParameter<bool>("requireStableBeam");
34 }
35 
36 // ------------ method called for each event ------------
38  using namespace edm;
39 
41  iEvent.getByToken(digiTokenEB_, pDigiEB);
42 
44  iEvent.getByToken(digiTokenEE_, pDigiEE);
45 
46  // Only Events with stable beam
47 
48  if (requireStableBeam_) {
49  edm::Handle<TCDSRecord> tcdsData;
50  iEvent.getByToken(tcdsToken_, tcdsData);
51  int beamMode = tcdsData->getBST().getBeamMode();
52  if (beamMode != BSTRecord::BeamMode::STABLE)
53  return;
54  }
55 
56  for (EBDigiCollection::const_iterator pDigi = pDigiEB->begin(); pDigi != pDigiEB->end(); ++pDigi) {
57  EBDetId id = pDigi->id();
58  uint32_t hashedId = id.hashedIndex();
59 
60  EBDataFrame digi(*pDigi);
61 
62  if (checkSignal_) {
63  uint16_t maxdiff = *std::max_element(digi.frame().begin(), digi.frame().end(), adc_compare) -
64  *std::min_element(digi.frame().begin(), digi.frame().end(), adc_compare);
65  if (maxdiff > sThresholdEB_)
66  continue; // assume there is signal in this frame
67  }
68 
69  //for (auto& mgpasample : digi.frame()) meEB_[hashedId]->Fill(mgpasample&0xFFF);
70  for (edm::DataFrame::iterator mgpasample = digi.frame().begin();
71  mgpasample != digi.frame().begin() + pedestalSamples_;
72  ++mgpasample)
73  meEB_[hashedId]->Fill(*mgpasample & 0xFFF);
74 
75  } // eb digis
76 
77  for (EEDigiCollection::const_iterator pDigi = pDigiEE->begin(); pDigi != pDigiEE->end(); ++pDigi) {
78  EEDetId id = pDigi->id();
79  uint32_t hashedId = id.hashedIndex();
80 
81  EEDataFrame digi(*pDigi);
82 
83  if (checkSignal_) {
84  uint16_t maxdiff = *std::max_element(digi.frame().begin(), digi.frame().end(), adc_compare) -
85  *std::min_element(digi.frame().begin(), digi.frame().end(), adc_compare);
86  if (maxdiff > sThresholdEE_)
87  continue; // assume there is signal in this frame
88  }
89 
90  //for (auto& mgpasample : digi.frame()) meEE_[hashedId]->Fill(mgpasample&0xFFF);
91  for (edm::DataFrame::iterator mgpasample = digi.frame().begin();
92  mgpasample != digi.frame().begin() + pedestalSamples_;
93  ++mgpasample)
94  meEE_[hashedId]->Fill(*mgpasample & 0xFFF);
95 
96  } // ee digis
97 }
98 
100 
102 
105  desc.setUnknown();
106  descriptions.addDefault(desc);
107 }
108 
110  ibooker.cd();
111  ibooker.setCurrentFolder(dqmDir_);
112 
114  es.get<EcalPedestalsRcd>().get(peds);
115 
116  for (uint32_t i = 0; i < EBDetId::kSizeForDenseIndexing; ++i) {
117  ibooker.setCurrentFolder(dqmDir_ + "/EB/" + std::to_string(int(i / 100)));
118 
119  std::string hname = "eb_" + std::to_string(i);
121  int centralBin = fixedBookingCenterBin_;
122 
123  if (dynamicBooking_) {
124  centralBin = int((peds->find(id))->mean_x12);
125  }
126 
127  int min = centralBin - nBins_ / 2;
128  int max = centralBin + nBins_ / 2;
129 
130  meEB_.push_back(ibooker.book1D(hname, hname, nBins_, min, max));
131  }
132 
133  for (uint32_t i = 0; i < EEDetId::kSizeForDenseIndexing; ++i) {
134  ibooker.setCurrentFolder(dqmDir_ + "/EE/" + std::to_string(int(i / 100)));
135 
136  std::string hname = "ee_" + std::to_string(i);
137 
139  int centralBin = fixedBookingCenterBin_;
140 
141  if (dynamicBooking_) {
142  centralBin = int((peds->find(id))->mean_x12);
143  }
144 
145  int min = centralBin - nBins_ / 2;
146  int max = centralBin + nBins_ / 2;
147 
148  meEE_.push_back(ibooker.book1D(hname, hname, nBins_, min, max));
149  }
150 }
static EEDetId detIdFromDenseIndex(uint32_t din)
Definition: EEDetId.h:220
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX)
Definition: DQMStore.cc:239
T getParameter(std::string const &) const
int hashedIndex() const
get a compact index for arrays
Definition: EBDetId.h:82
const BSTRecord & getBST() const
Definition: TCDSRecord.h:100
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
constexpr iterator end()
Definition: DataFrame.h:35
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
edm::EDGetTokenT< EEDigiCollection > digiTokenEE_
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:418
const_iterator begin() const
uint16_t const getBeamMode() const
Definition: BSTRecord.h:70
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
static EBDetId detIdFromDenseIndex(uint32_t di)
Definition: EBDetId.h:107
std::vector< MonitorElement * > meEE_
void analyze(const edm::Event &, const edm::EventSetup &) override
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
static bool adc_compare(uint16_t a, uint16_t b)
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
T min(T a, T b)
Definition: MathUtil.h:58
Definition: DetId.h:17
ECALpedestalPCLworker(const edm::ParameterSet &)
int hashedIndex() const
Definition: EEDetId.h:183
constexpr iterator begin()
Definition: DataFrame.h:33
edm::DataFrame const & frame() const
Definition: EcalDataFrame.h:50
boost::transform_iterator< IterHelp, boost::counting_iterator< int > > const_iterator
edm::EDGetTokenT< EBDigiCollection > digiTokenEB_
const_iterator end() const
std::vector< MonitorElement * > meEB_
HLT enums.
data_type * iterator
Definition: DataFrame.h:20
T get() const
Definition: EventSetup.h:73
const_iterator find(uint32_t rawId) const
edm::EDGetTokenT< TCDSRecord > tcdsToken_
Definition: Run.h:45