CMS 3D CMS Logo

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