CMS 3D CMS Logo

ECALpedestalPCLworker.cc
Go to the documentation of this file.
1 
6 #include <sstream>
7 
9  : pedestalToken_(esConsumes<edm::Transition::BeginRun>()) {
10  edm::InputTag digiTagEB = iConfig.getParameter<edm::InputTag>("BarrelDigis");
11  edm::InputTag digiTagEE = iConfig.getParameter<edm::InputTag>("EndcapDigis");
12 
13  digiTokenEB_ = consumes<EBDigiCollection>(digiTagEB);
14  digiTokenEE_ = consumes<EEDigiCollection>(digiTagEE);
15 
16  pedestalSamples_ = iConfig.getParameter<uint32_t>("pedestalSamples");
17  checkSignal_ = iConfig.getParameter<bool>("checkSignal");
18  sThresholdEB_ = iConfig.getParameter<uint32_t>("sThresholdEB");
19  sThresholdEE_ = iConfig.getParameter<uint32_t>("sThresholdEE");
20 
21  dynamicBooking_ = iConfig.getParameter<bool>("dynamicBooking");
22  fixedBookingCenterBin_ = iConfig.getParameter<int>("fixedBookingCenterBin");
23  nBins_ = iConfig.getParameter<int>("nBins");
24  dqmDir_ = iConfig.getParameter<std::string>("dqmDir");
25 
26  edm::InputTag tcdsRecord = iConfig.getParameter<edm::InputTag>("tcdsRecord");
27  tcdsToken_ = consumes<TCDSRecord>(tcdsRecord);
28  requireStableBeam_ = iConfig.getParameter<bool>("requireStableBeam");
29 }
30 
31 // ------------ method called for each event ------------
33  using namespace edm;
34 
36  iEvent.getByToken(digiTokenEB_, pDigiEB);
37 
39  iEvent.getByToken(digiTokenEE_, pDigiEE);
40 
41  // Only Events with stable beam
42 
43  if (requireStableBeam_) {
44  edm::Handle<TCDSRecord> tcdsData;
45  iEvent.getByToken(tcdsToken_, tcdsData);
46  int beamMode = tcdsData->getBST().getBeamMode();
47  if (beamMode != BSTRecord::BeamMode::STABLE)
48  return;
49  }
50 
51  for (EBDigiCollection::const_iterator pDigi = pDigiEB->begin(); pDigi != pDigiEB->end(); ++pDigi) {
52  EBDetId id = pDigi->id();
53  uint32_t hashedId = id.hashedIndex();
54 
55  EBDataFrame digi(*pDigi);
56 
57  if (checkSignal_) {
58  uint16_t maxdiff = *std::max_element(digi.frame().begin(), digi.frame().end(), adc_compare) -
59  *std::min_element(digi.frame().begin(), digi.frame().end(), adc_compare);
60  if (maxdiff > sThresholdEB_)
61  continue; // assume there is signal in this frame
62  }
63 
64  //for (auto& mgpasample : digi.frame()) meEB_[hashedId]->Fill(mgpasample&0xFFF);
65  for (edm::DataFrame::iterator mgpasample = digi.frame().begin();
66  mgpasample != digi.frame().begin() + pedestalSamples_;
67  ++mgpasample)
68  meEB_[hashedId]->Fill(*mgpasample & 0xFFF);
69 
70  } // eb digis
71 
72  for (EEDigiCollection::const_iterator pDigi = pDigiEE->begin(); pDigi != pDigiEE->end(); ++pDigi) {
73  EEDetId id = pDigi->id();
74  uint32_t hashedId = id.hashedIndex();
75 
76  EEDataFrame digi(*pDigi);
77 
78  if (checkSignal_) {
79  uint16_t maxdiff = *std::max_element(digi.frame().begin(), digi.frame().end(), adc_compare) -
80  *std::min_element(digi.frame().begin(), digi.frame().end(), adc_compare);
81  if (maxdiff > sThresholdEE_)
82  continue; // assume there is signal in this frame
83  }
84 
85  //for (auto& mgpasample : digi.frame()) meEE_[hashedId]->Fill(mgpasample&0xFFF);
86  for (edm::DataFrame::iterator mgpasample = digi.frame().begin();
87  mgpasample != digi.frame().begin() + pedestalSamples_;
88  ++mgpasample)
89  meEE_[hashedId]->Fill(*mgpasample & 0xFFF);
90 
91  } // ee digis
92 }
93 
96  desc.setUnknown();
97  descriptions.addDefault(desc);
98 }
99 
101  ibooker.cd();
102  ibooker.setCurrentFolder(dqmDir_);
103 
104  const auto& peds = es.getData(pedestalToken_);
105 
106  for (uint32_t i = 0; i < EBDetId::kSizeForDenseIndexing; ++i) {
107  ibooker.setCurrentFolder(dqmDir_ + "/EB/" + std::to_string(int(i / 100)));
108 
109  std::string hname = "eb_" + std::to_string(i);
111  int centralBin = fixedBookingCenterBin_;
112 
113  if (dynamicBooking_) {
114  centralBin = int((peds.find(id))->mean_x12);
115  }
116 
117  int min = centralBin - nBins_ / 2;
118  int max = centralBin + nBins_ / 2;
119 
120  meEB_.push_back(ibooker.book1D(hname, hname, nBins_, min, max));
121  }
122 
123  for (uint32_t i = 0; i < EEDetId::kSizeForDenseIndexing; ++i) {
124  ibooker.setCurrentFolder(dqmDir_ + "/EE/" + std::to_string(int(i / 100)));
125 
126  std::string hname = "ee_" + std::to_string(i);
127 
129  int centralBin = fixedBookingCenterBin_;
130 
131  if (dynamicBooking_) {
132  centralBin = int((peds.find(id))->mean_x12);
133  }
134 
135  int min = centralBin - nBins_ / 2;
136  int max = centralBin + nBins_ / 2;
137 
138  meEE_.push_back(ibooker.book1D(hname, hname, nBins_, min, max));
139  }
140 }
static EEDetId detIdFromDenseIndex(uint32_t din)
Definition: EEDetId.h:220
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
constexpr iterator end()
Definition: DataFrame.h:35
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
edm::EDGetTokenT< EEDigiCollection > digiTokenEE_
const edm::ESGetToken< EcalPedestals, EcalPedestalsRcd > pedestalToken_
static std::string to_string(const XMLCh *ch)
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)
Transition
Definition: Transition.h:12
const BSTRecord & getBST() const
Definition: TCDSRecord.h:100
const_iterator end() const
uint16_t const getBeamMode() const
Definition: BSTRecord.h:70
Definition: DetId.h:17
ECALpedestalPCLworker(const edm::ParameterSet &)
const_iterator begin() const
The iterator returned can not safely be used across threads.
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_
std::vector< MonitorElement * > meEB_
HLT enums.
data_type * iterator
Definition: DataFrame.h:20
edm::EDGetTokenT< TCDSRecord > tcdsToken_
int hashedIndex() const
get a compact index for arrays
Definition: EBDetId.h:82
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
int hashedIndex() const
Definition: EEDetId.h:183
Definition: Run.h:45