CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes
ECALpedestalPCLworker Class Reference

#include <ECALpedestalPCLworker.h>

Inheritance diagram for ECALpedestalPCLworker:
one::DQMEDAnalyzer< T > one::dqmimplementation::DQMBaseClass< T... >

Public Member Functions

 ECALpedestalPCLworker (const edm::ParameterSet &)
 
- Public Member Functions inherited from one::DQMEDAnalyzer< T >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void beginJob () override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
void endJob () override
 

Static Private Member Functions

static bool adc_compare (uint16_t a, uint16_t b)
 

Private Attributes

bool checkSignal_
 
edm::EDGetTokenT< EBDigiCollectiondigiTokenEB_
 
edm::EDGetTokenT< EEDigiCollectiondigiTokenEE_
 
std::string dqmDir_
 
bool dynamicBooking_
 
int fixedBookingCenterBin_
 
std::vector< MonitorElement * > meEB_
 
std::vector< MonitorElement * > meEE_
 
int nBins_
 
uint32_t pedestalSamples_
 
bool requireStableBeam_
 
uint32_t sThresholdEB_
 
uint32_t sThresholdEE_
 
edm::EDGetTokenT< TCDSRecordtcdsToken_
 

Detailed Description

Description: Fill DQM histograms with pedestals. Intended to be used on laser data from the TestEnablesEcalHcal dataset

Definition at line 35 of file ECALpedestalPCLworker.h.

Constructor & Destructor Documentation

ECALpedestalPCLworker::ECALpedestalPCLworker ( const edm::ParameterSet iConfig)
explicit

Definition at line 14 of file ECALpedestalPCLworker.cc.

References checkSignal_, digiTokenEB_, digiTokenEE_, dqmDir_, dynamicBooking_, fixedBookingCenterBin_, edm::ParameterSet::getParameter(), nBins_, pedestalSamples_, requireStableBeam_, sThresholdEB_, sThresholdEE_, AlCaHLTBitMon_QueryRunRegistry::string, and tcdsToken_.

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 }
T getParameter(std::string const &) const
edm::EDGetTokenT< EEDigiCollection > digiTokenEE_
edm::EDGetTokenT< EBDigiCollection > digiTokenEB_
edm::EDGetTokenT< TCDSRecord > tcdsToken_

Member Function Documentation

static bool ECALpedestalPCLworker::adc_compare ( uint16_t  a,
uint16_t  b 
)
inlinestaticprivate

Definition at line 72 of file ECALpedestalPCLworker.h.

Referenced by analyze().

72 { return ( a&0xFFF ) < (b&0xFFF);}
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
void ECALpedestalPCLworker::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 42 of file ECALpedestalPCLworker.cc.

References adc_compare(), edm::DataFrame::begin(), edm::DataFrameContainer::begin(), checkSignal_, digiTokenEB_, digiTokenEE_, edm::DataFrame::end(), edm::DataFrameContainer::end(), HcalObjRepresent::Fill(), EcalDataFrame::frame(), BSTRecord::getBeamMode(), TCDSRecord::getBST(), edm::Event::getByToken(), EBDetId::hashedIndex(), EEDetId::hashedIndex(), meEB_, meEE_, pedestalSamples_, requireStableBeam_, sThresholdEB_, sThresholdEE_, and tcdsToken_.

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 }
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
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
std::vector< MonitorElement * > meEE_
static bool adc_compare(uint16_t a, uint16_t b)
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
int hashedIndex() const
Definition: EEDetId.h:182
edm::EDGetTokenT< EBDigiCollection > digiTokenEB_
const_iterator end() const
std::vector< MonitorElement * > meEB_
HLT enums.
data_type * iterator
Definition: DataFrame.h:21
edm::EDGetTokenT< TCDSRecord > tcdsToken_
void ECALpedestalPCLworker::beginJob ( void  )
overrideprivate

Definition at line 113 of file ECALpedestalPCLworker.cc.

114 {
115 }
void ECALpedestalPCLworker::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  run,
edm::EventSetup const &  es 
)
overrideprivate

Definition at line 133 of file ECALpedestalPCLworker.cc.

References DQMStore::IBooker::book1D(), DQMStore::IBooker::cd(), EBDetId::detIdFromDenseIndex(), EEDetId::detIdFromDenseIndex(), dqmDir_, dynamicBooking_, EcalCondObjectContainer< T >::find(), fixedBookingCenterBin_, edm::EventSetup::get(), mps_fire::i, createfilelist::int, EBDetId::kSizeForDenseIndexing, EEDetId::kSizeForDenseIndexing, SiStripPI::max, meEB_, meEE_, min(), nBins_, DQMStore::IBooker::setCurrentFolder(), and AlCaHLTBitMon_QueryRunRegistry::string.

133  {
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
static EBDetId detIdFromDenseIndex(uint32_t di)
Definition: EBDetId.h:111
std::vector< MonitorElement * > meEE_
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
T min(T a, T b)
Definition: MathUtil.h:58
Definition: DetId.h:18
std::vector< MonitorElement * > meEB_
const_iterator find(uint32_t rawId) const
void ECALpedestalPCLworker::endJob ( void  )
overrideprivate

Definition at line 119 of file ECALpedestalPCLworker.cc.

Referenced by o2olib.O2ORunMgr::executeJob().

120 {
121 }
void ECALpedestalPCLworker::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Member Data Documentation

bool ECALpedestalPCLworker::checkSignal_
private

Definition at line 58 of file ECALpedestalPCLworker.h.

Referenced by analyze(), and ECALpedestalPCLworker().

edm::EDGetTokenT<EBDigiCollection> ECALpedestalPCLworker::digiTokenEB_
private

Definition at line 50 of file ECALpedestalPCLworker.h.

Referenced by analyze(), and ECALpedestalPCLworker().

edm::EDGetTokenT<EEDigiCollection> ECALpedestalPCLworker::digiTokenEE_
private

Definition at line 51 of file ECALpedestalPCLworker.h.

Referenced by analyze(), and ECALpedestalPCLworker().

std::string ECALpedestalPCLworker::dqmDir_
private

Definition at line 67 of file ECALpedestalPCLworker.h.

Referenced by bookHistograms(), and ECALpedestalPCLworker().

bool ECALpedestalPCLworker::dynamicBooking_
private

Definition at line 64 of file ECALpedestalPCLworker.h.

Referenced by bookHistograms(), and ECALpedestalPCLworker().

int ECALpedestalPCLworker::fixedBookingCenterBin_
private

Definition at line 65 of file ECALpedestalPCLworker.h.

Referenced by bookHistograms(), and ECALpedestalPCLworker().

std::vector<MonitorElement *> ECALpedestalPCLworker::meEB_
private

Definition at line 54 of file ECALpedestalPCLworker.h.

Referenced by analyze(), and bookHistograms().

std::vector<MonitorElement *> ECALpedestalPCLworker::meEE_
private

Definition at line 55 of file ECALpedestalPCLworker.h.

Referenced by analyze(), and bookHistograms().

int ECALpedestalPCLworker::nBins_
private

Definition at line 66 of file ECALpedestalPCLworker.h.

Referenced by bookHistograms(), and ECALpedestalPCLworker().

uint32_t ECALpedestalPCLworker::pedestalSamples_
private

Definition at line 57 of file ECALpedestalPCLworker.h.

Referenced by analyze(), and ECALpedestalPCLworker().

bool ECALpedestalPCLworker::requireStableBeam_
private

Definition at line 69 of file ECALpedestalPCLworker.h.

Referenced by analyze(), and ECALpedestalPCLworker().

uint32_t ECALpedestalPCLworker::sThresholdEB_
private

Definition at line 59 of file ECALpedestalPCLworker.h.

Referenced by analyze(), and ECALpedestalPCLworker().

uint32_t ECALpedestalPCLworker::sThresholdEE_
private

Definition at line 62 of file ECALpedestalPCLworker.h.

Referenced by analyze(), and ECALpedestalPCLworker().

edm::EDGetTokenT<TCDSRecord> ECALpedestalPCLworker::tcdsToken_
private

Definition at line 52 of file ECALpedestalPCLworker.h.

Referenced by analyze(), and ECALpedestalPCLworker().