CMS 3D CMS Logo

List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
RPCFEDIntegrity Class Reference

#include <RPCFEDIntegrity.h>

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

Public Member Functions

void analyze (const edm::Event &iEvent, const edm::EventSetup &c) override
 Analyze. More...
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 Begin Lumi block. More...
 
 RPCFEDIntegrity (const edm::ParameterSet &ps)
 Constructor. More...
 
 ~RPCFEDIntegrity () override
 Destructor. More...
 
- Public Member Functions inherited from one::DQMEDAnalyzer< T >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Private Types

enum  fedHisto { Entries, Fatal, NonFatal }
 

Private Member Functions

void bookFEDMe (DQMStore::IBooker &)
 
void labelBins (MonitorElement *myMe)
 

Private Attributes

int FATAL_LIMIT
 
MonitorElementfedMe_ [3]
 
std::vector< std::string > histoName_
 
bool init_
 
int maxFEDNum_
 
bool merge_
 
int minFEDNum_
 
int numOfFED_
 
std::string prefixDir_
 
edm::EDGetTokenT< RPCRawDataCountsrawCountsLabel_
 

Detailed Description

Definition at line 21 of file RPCFEDIntegrity.h.

Member Enumeration Documentation

Enumerator
Entries 
Fatal 
NonFatal 

Definition at line 52 of file RPCFEDIntegrity.h.

Constructor & Destructor Documentation

RPCFEDIntegrity::RPCFEDIntegrity ( const edm::ParameterSet ps)

Constructor.

Definition at line 20 of file RPCFEDIntegrity.cc.

References FATAL_LIMIT, edm::ParameterSet::getUntrackedParameter(), init_, maxFEDNum_, merge_, minFEDNum_, numOfFED_, prefixDir_, rawCountsLabel_, and AlCaHLTBitMon_QueryRunRegistry::string.

20  {
21  edm::LogVerbatim ("rpcfedintegrity") << "[RPCFEDIntegrity]: Constructor";
22 
23  rawCountsLabel_ = consumes<RPCRawDataCounts>(ps.getUntrackedParameter<edm::InputTag>("RPCRawCountsInputTag"));
24  prefixDir_ = ps.getUntrackedParameter<std::string>("RPCPrefixDir", "RPC/FEDIntegrity");
25  merge_ = ps.getUntrackedParameter<bool>("MergeRuns", false);
26  minFEDNum_ = ps.getUntrackedParameter<int>("MinimumFEDID", 790);
27  maxFEDNum_ = ps.getUntrackedParameter<int>("MaximumFEDID", 792);
28 
29  init_ = false;
31  FATAL_LIMIT = 5;
32 }
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< RPCRawDataCounts > rawCountsLabel_
std::string prefixDir_
RPCFEDIntegrity::~RPCFEDIntegrity ( )
override

Destructor.

Definition at line 34 of file RPCFEDIntegrity.cc.

34  {
35  edm::LogVerbatim ("rpcfedintegrity") << "[RPCFEDIntegrity]: Destructor ";
36  // dbe_=0;
37 }

Member Function Documentation

void RPCFEDIntegrity::analyze ( const edm::Event iEvent,
const edm::EventSetup c 
)
override

Analyze.

Definition at line 51 of file RPCFEDIntegrity.cc.

References Entries, Fatal, RPCRawDataCounts::fedBxRecords(), RPCRawDataCounts::fedErrorRecords(), RPCRawDataCounts::fedFormatErrors(), fedMe_, HcalObjRepresent::Fill(), edm::Event::getByToken(), edm::HandleBase::isValid(), maxFEDNum_, minFEDNum_, NonFatal, edm::Handle< T >::product(), and rawCountsLabel_.

51  {
52 
53  //get hold of raw data counts
55  iEvent.getByToken (rawCountsLabel_, rawCounts);
56  if(!rawCounts.isValid()) return;
57 
58 
59  const RPCRawDataCounts & counts = *rawCounts.product();
60 
61  for (int fed=minFEDNum_; fed <=maxFEDNum_; ++fed) {
62  if (counts.fedBxRecords(fed) ) fedMe_[Entries]->Fill(fed);
63  if (counts.fedFormatErrors(fed)) fedMe_[Fatal]->Fill(fed);
64  if (counts.fedErrorRecords(fed)) fedMe_[NonFatal]->Fill(fed);
65  }
66 }
int fedErrorRecords(int fedId) const
int fedBxRecords(int fedId) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
int fedFormatErrors(int fedId) const
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
bool isValid() const
Definition: HandleBase.h:74
MonitorElement * fedMe_[3]
T const * product() const
Definition: Handle.h:81
edm::EDGetTokenT< RPCRawDataCounts > rawCountsLabel_
void RPCFEDIntegrity::bookFEDMe ( DQMStore::IBooker ibooker)
private

Definition at line 70 of file RPCFEDIntegrity.cc.

References DQMStore::IBooker::book1D(), DQMStore::IBooker::cd(), Entries, Fatal, fedMe_, init_, labelBins(), maxFEDNum_, minFEDNum_, NonFatal, numOfFED_, prefixDir_, and DQMStore::IBooker::setCurrentFolder().

Referenced by bookHistograms().

70  {
71 
72 
73  ibooker.cd();
75 
76  fedMe_[Entries] = ibooker.book1D("FEDEntries","FED Entries",numOfFED_, minFEDNum_, maxFEDNum_ +1);
77  this->labelBins(fedMe_[Entries]);
78  fedMe_[Fatal] = ibooker.book1D("FEDFatal","FED Fatal Errors",numOfFED_, minFEDNum_, maxFEDNum_ +1);
79  this->labelBins(fedMe_[Fatal]);
80  fedMe_[NonFatal] = ibooker.book1D("FEDNonFatal","FED NON Fatal Errors",numOfFED_, minFEDNum_, maxFEDNum_ +1);
81  this->labelBins(fedMe_[NonFatal]);
82 
83  init_ = true;
84 }
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
MonitorElement * fedMe_[3]
void labelBins(MonitorElement *myMe)
std::string prefixDir_
void RPCFEDIntegrity::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  iRun,
edm::EventSetup const &  iSetup 
)
override

Begin Lumi block.

Definition at line 40 of file RPCFEDIntegrity.cc.

References bookFEDMe().

42  {
43 
44 
45  edm::LogVerbatim ("rpcfedintegrity") << "[RPCFEDIntegrity]: Begin booking histograms ";
46 
47  this->bookFEDMe(ibooker);
48 }
void bookFEDMe(DQMStore::IBooker &)
void RPCFEDIntegrity::labelBins ( MonitorElement myMe)
private

Definition at line 86 of file RPCFEDIntegrity.cc.

References MonitorElement::getNbinsX(), mps_fire::i, minFEDNum_, numOfFED_, MonitorElement::setBinLabel(), and fw3dlego::xbins.

Referenced by bookFEDMe().

86  {
87 
88  int xbins = myMe->getNbinsX();
89 
90  if (xbins!= numOfFED_ ) return;
91  std::stringstream xLabel;
92 
93  for (int i = 0; i<xbins; i++){
94  xLabel.str("");
95  int fedNum = minFEDNum_ + i;
96  xLabel<<fedNum;
97  myMe->setBinLabel(i+1, xLabel.str(),1);
98  }
99 }
const double xbins[]
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
int getNbinsX() const
get # of bins in X-axis

Member Data Documentation

int RPCFEDIntegrity::FATAL_LIMIT
private

Definition at line 50 of file RPCFEDIntegrity.h.

Referenced by RPCFEDIntegrity().

MonitorElement* RPCFEDIntegrity::fedMe_[3]
private

Definition at line 54 of file RPCFEDIntegrity.h.

Referenced by analyze(), and bookFEDMe().

std::vector<std::string> RPCFEDIntegrity::histoName_
private

Definition at line 57 of file RPCFEDIntegrity.h.

bool RPCFEDIntegrity::init_
private

Definition at line 48 of file RPCFEDIntegrity.h.

Referenced by bookFEDMe(), and RPCFEDIntegrity().

int RPCFEDIntegrity::maxFEDNum_
private

Definition at line 56 of file RPCFEDIntegrity.h.

Referenced by analyze(), bookFEDMe(), and RPCFEDIntegrity().

bool RPCFEDIntegrity::merge_
private

Definition at line 48 of file RPCFEDIntegrity.h.

Referenced by RPCFEDIntegrity().

int RPCFEDIntegrity::minFEDNum_
private

Definition at line 56 of file RPCFEDIntegrity.h.

Referenced by analyze(), bookFEDMe(), labelBins(), and RPCFEDIntegrity().

int RPCFEDIntegrity::numOfFED_
private

Definition at line 56 of file RPCFEDIntegrity.h.

Referenced by bookFEDMe(), labelBins(), and RPCFEDIntegrity().

std::string RPCFEDIntegrity::prefixDir_
private

Definition at line 46 of file RPCFEDIntegrity.h.

Referenced by bookFEDMe(), and RPCFEDIntegrity().

edm::EDGetTokenT<RPCRawDataCounts> RPCFEDIntegrity::rawCountsLabel_
private

Definition at line 43 of file RPCFEDIntegrity.h.

Referenced by analyze(), and RPCFEDIntegrity().