CMS 3D CMS Logo

HcalFEDIntegrityTask.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HcalTasks
4 // Class: HcalFEDIntegrityTask
5 // Original Author: Long Wang - University of Maryland
6 //
7 
12 
15 
17 
18 #include <vector>
19 #include <string>
20 
21 using namespace std;
22 using namespace edm;
23 
25 public:
27  ~HcalFEDIntegrityTask() override;
28 
29  void dqmBeginRun(const edm::Run &, edm::EventSetup const &) override;
30  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
31  void analyze(const edm::Event &, const edm::EventSetup &) override;
32 
33  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
34 
35 private:
36  void labelBins(MonitorElement *me);
37 
40 
41  int numOfFED_, minFEDNum_, maxFEDNum_;
46 };
47 
49  : tokFEDs_(consumes<FEDRawDataCollection>(
50  ps.getUntrackedParameter<edm::InputTag>("tagFEDs", edm::InputTag("rawDataCollector")))),
51  tokReport_(consumes<HcalUnpackerReport>(
52  ps.getUntrackedParameter<edm::InputTag>("tagReport", edm::InputTag("hcalDigis")))),
53  minFEDNum_(ps.getUntrackedParameter<int>("MinHcalFEDID", FEDNumbering::MINHCALuTCAFEDID)),
54  maxFEDNum_(ps.getUntrackedParameter<int>("MaxHcalFEDID", FEDNumbering::MAXHCALuTCAFEDID)),
55  dirName_(ps.getUntrackedParameter<std::string>("DirName", "Hcal/FEDIntegrity/")) {
56  LogInfo("HcalDQM") << "HcalFEDIntegrityTask::HcalFEDIntegrityTask: Constructor Initialization" << endl;
58 }
59 
61  LogInfo("HcalDQM") << "HcalFEDIntegrityTask::~HcalFEDIntegrityTask: Destructor" << endl;
62 }
63 
64 void HcalFEDIntegrityTask::dqmBeginRun(const edm::Run &iRun, const edm::EventSetup &iSetup) {}
65 
67  edm::Run const &iRun,
68  edm::EventSetup const &iSetup) {
69  iBooker.cd();
70  iBooker.setCurrentFolder(dirName_);
71 
72  meFEDEntries_ = iBooker.book1D("FEDEntries", "FED Entries", numOfFED_, minFEDNum_, maxFEDNum_ + 1);
73  this->labelBins(meFEDEntries_);
74  meFEDFatal_ = iBooker.book1D("FEDFatal", "FED Fatal Errors", numOfFED_, minFEDNum_, maxFEDNum_ + 1);
75  this->labelBins(meFEDFatal_);
76  meFEDNonFatal_ = iBooker.book1D("FEDNonFatal", "FED NON Fatal Errors", numOfFED_, minFEDNum_, maxFEDNum_ + 1);
78 }
79 
83  iEvent.getByToken(tokFEDs_, raw);
84  iEvent.getByToken(tokReport_, report);
85 
86  // FEDs with unpacking errors: https://github.com/cms-sw/cmssw/blob/master/EventFilter/HcalRawToDigi/plugins/HcalRawToDigi.cc#L235-L262
87  const std::vector<int> FedsError = (*report).getFedsError();
88  for (auto &fed : FedsError) {
89  if (fed < 1000)
90  LogWarning("HcalDQM") << "HcalFEDIntegrityTask::analyze: obsoleteVME FEDs from HcalUnpackerReport" << endl;
91  meFEDFatal_->Fill(fed);
92  }
93 
95  // Same checks as the RawTask Summary map
97  for (int fed = FEDNumbering::MINHCALuTCAFEDID; fed <= FEDNumbering::MAXHCALuTCAFEDID; fed++) {
98  const FEDRawData &fedRawData = raw->FEDData(fed);
99  if (fedRawData.size() != 0) {
100  meFEDEntries_->Fill(fed);
101  }
102 
103  hcal::AMC13Header const *amc13 = (hcal::AMC13Header const *)fedRawData.data();
104  if (!amc13) {
105  continue;
106  }
107 
108  uint32_t bcn = amc13->bunchId();
109  uint32_t orn = amc13->orbitNumber() & 0xFFFF; // LS 16bits only
110  uint32_t evn = amc13->l1aNumber();
111  int namc = amc13->NAMC();
112 
113  // looping over AMC in this packet
114  for (int iamc = 0; iamc < namc; iamc++) {
115  if (!amc13->AMCEnabled(iamc) || !amc13->AMCDataPresent(iamc) || !amc13->AMCCRCOk(iamc) ||
116  amc13->AMCSegmented(iamc)) {
117  LogWarning("HcalDQM") << "HcalFEDIntegrityTask::analyze: AMC issue on iamc" << iamc << endl;
118  continue;
119  }
120 
121  HcalUHTRData uhtr(amc13->AMCPayload(iamc), amc13->AMCSize(iamc));
122  uint32_t uhtr_evn = uhtr.l1ANumber();
123  uint32_t uhtr_bcn = uhtr.bunchNumber();
124  uint32_t uhtr_orn = uhtr.orbitNumber();
125 
126  if (uhtr_evn != evn || uhtr_bcn != bcn || uhtr_orn != orn) {
127  if (std::find(FedsError.begin(), FedsError.end(), fed) ==
128  FedsError.end()) // FED not already in the error list from unpacker report
129  meFEDFatal_->Fill(fed);
130  break; // one mismatch is sufficient enough to determine it's a bad data
131  }
132  }
133 
134  } // end of Hcal FED looping
135 }
136 
138  int xbins = me->getNbinsX();
139 
140  if (xbins != numOfFED_)
141  return;
142 
143  for (int i = 0; i < xbins; i++) {
144  const std::string xLabel = fmt::format("{}", minFEDNum_ + i);
145  me->setBinLabel(i + 1, xLabel, 1);
146  }
147 }
148 
151  desc.addUntracked<std::string>("name", "HcalFEDIntegrityTask");
152  desc.addUntracked<int>("debug", 0);
153  desc.addUntracked<edm::InputTag>("tagFEDs", edm::InputTag("rawDataCollector"));
154  desc.addUntracked<edm::InputTag>("tagReport", edm::InputTag("hcalDigis"));
155  desc.addUntracked<int>(
156  "MinHcalFEDID",
157  FEDNumbering::MINHCALuTCAFEDID); // Assuming no more VME FEDs after LS2, according to Hcal Phase1 upgrade.
158  desc.addUntracked<int>("MaxHcalFEDID", FEDNumbering::MAXHCALuTCAFEDID);
159  desc.addUntracked<std::string>("DirName", "Hcal/FEDIntegrity/");
160  descriptions.addWithDefaultLabel(desc);
161 }
162 
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
void labelBins(MonitorElement *me)
void analyze(const edm::Event &, const edm::EventSetup &) override
const double xbins[]
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
edm::EDGetTokenT< FEDRawDataCollection > tokFEDs_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
MonitorElement * meFEDEntries_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
example_stream void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
MonitorElement * meFEDFatal_
void dqmBeginRun(const edm::Run &, edm::EventSetup const &) override
const FEDRawData & FEDData(int fedid) const
retrieve data for fed
Log< level::Info, false > LogInfo
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
uint32_t l1ANumber() const
Get the HTR event number.
Definition: HcalUHTRData.h:87
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
HLT enums.
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Log< level::Warning, false > LogWarning
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
edm::EDGetTokenT< HcalUnpackerReport > tokReport_
MonitorElement * meFEDNonFatal_
Definition: Run.h:45
HcalFEDIntegrityTask(const edm::ParameterSet &ps)