CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
CastorMonitorModule.cc
Go to the documentation of this file.
7 //#include "CondFormats/CastorObjects/interface/CastorChannelQuality.h"
8 #include <string>
9 
10 //**************************************************************//
11 //***************** CastorMonitorModule ******************//
12 //***************** Author: Dmytro Volyanskyy ******************//
13 //***************** Date : 22.11.2008 (first version) *********//
19 //**************************************************************//
20 
21 //---- critical revision 26.06.2014 (Vladimir Popov)
22 
23 //**************************************************************//
24 
25 using namespace std;
26 using namespace edm;
27 
29  : castorDbServiceToken_{esConsumes<CastorDbService, CastorDbRecord>()} {
30  fVerbosity = ps.getUntrackedParameter<int>("debug", 0);
31  subsystemname_ = ps.getUntrackedParameter<std::string>("subSystemFolder", "Castor");
32  inputTokenRaw_ = consumes<FEDRawDataCollection>(ps.getParameter<edm::InputTag>("rawLabel"));
33  inputTokenReport_ = consumes<HcalUnpackerReport>(ps.getParameter<edm::InputTag>("unpackerReportLabel"));
34  inputTokenDigi_ = consumes<CastorDigiCollection>(ps.getParameter<edm::InputTag>("digiLabel"));
35  inputTokenRecHitCASTOR_ = consumes<CastorRecHitCollection>(ps.getParameter<edm::InputTag>("CastorRecHitLabel"));
36  inputTokenCastorTowers_ = consumes<CastorTowerCollection>(ps.getParameter<edm::InputTag>("CastorTowerLabel"));
37  JetAlgorithm = consumes<BasicJetCollection>(ps.getParameter<edm::InputTag>("CastorBasicJetsLabel"));
38  tokenTriggerResults = consumes<edm::TriggerResults>(ps.getParameter<edm::InputTag>("tagTriggerResults"));
39 
40  showTiming_ = ps.getUntrackedParameter<bool>("showTiming", false);
41 
42  if (ps.getUntrackedParameter<bool>("DigiMonitor", false))
43  DigiMon_ = std::make_unique<CastorDigiMonitor>(ps, consumesCollector());
44 
45  if (ps.getUntrackedParameter<bool>("RecHitMonitor", false))
46  RecHitMon_ = std::make_unique<CastorRecHitMonitor>(ps);
47 
48  if (ps.getUntrackedParameter<bool>("LEDMonitor", false))
49  LedMon_ = std::make_unique<CastorLEDMonitor>(ps);
50 
51  ievt_ = 0;
52 }
53 
55 
57  if (fVerbosity > 0)
58  LogPrint("CastorMonitorModule") << "dqmBeginRun(start)";
59 }
60 
62  const edm::Run &iRun,
63  const edm::EventSetup &iSetup) {
64  if (DigiMon_) {
65  // Run histos only since there is endRun processing.
66  auto scope = DQMStore::IBooker::UseRunScope(ibooker);
67  DigiMon_->bookHistograms(ibooker, iRun, iSetup);
68  }
69  if (RecHitMon_) {
70  RecHitMon_->bookHistograms(ibooker, iRun);
71  }
72  if (LedMon_) {
73  LedMon_->bookHistograms(ibooker, iRun);
74  }
75 
77  char s[60];
78  sprintf(s, "CastorEventProducts");
79  CastorEventProduct = ibooker.book1DD(s, s, 6, -0.5, 5.5);
80  CastorEventProduct->setAxisTitle("Events", /* axis */ 2);
81  CastorEventProduct->setBinLabel(1, "FEDs/3");
82  CastorEventProduct->setBinLabel(2, "RawData");
83  CastorEventProduct->setBinLabel(3, "Digi");
84  CastorEventProduct->setBinLabel(4, "RecHits");
85  CastorEventProduct->setBinLabel(5, "Towers");
86  CastorEventProduct->setBinLabel(6, "Jets");
87 
88  sprintf(s, "CASTORUnpackReport");
89  hunpkrep = ibooker.bookProfile(s, s, 6, -0.5, 5.5, 100, 0, 1.e10, "");
90  hunpkrep->setBinLabel(1, "N_FEDs");
91  hunpkrep->setBinLabel(2, "SPIGOT_Err");
92  hunpkrep->setBinLabel(3, "empty");
93  hunpkrep->setBinLabel(4, "busy");
94  hunpkrep->setBinLabel(5, "OvF");
95  hunpkrep->setBinLabel(6, "BadDigis");
96  return;
97 }
98 
100  if (DigiMon_) {
101  DigiMon_->endRun();
102  }
103 }
104 
106  if (fVerbosity > 1)
107  LogPrint("CastorMonitorModule") << "analyze (start)";
108 
109  ievt_++;
110 
111  bool rawOK_ = true;
112  bool digiOK_ = true;
113  bool rechitOK_ = true, towerOK_ = true, jetsOK_ = true;
114  int nDigi = 0, nrecHits = 0, nTowers = 0, nJets = 0;
115 
117  iEvent.getByToken(tokenTriggerResults, TrigResults);
118 
120  iEvent.getByToken(inputTokenRaw_, RawData);
121  if (!RawData.isValid())
122  rawOK_ = false;
123 
124  float fedsUnpacked = 0.;
126  iEvent.getByToken(inputTokenReport_, report);
127  if (!report.isValid())
128  rawOK_ = false;
129  else {
130  const std::vector<int> feds = (*report).getFedsUnpacked();
131  fedsUnpacked = float(feds.size());
132  hunpkrep->Fill(0, fedsUnpacked);
133  hunpkrep->Fill(1, report->spigotFormatErrors());
134  hunpkrep->Fill(2, report->emptyEventSpigots());
135  hunpkrep->Fill(3, report->busySpigots());
136  hunpkrep->Fill(4, report->OFWSpigots());
137  hunpkrep->Fill(5, report->badQualityDigis());
138  }
139 
141  iEvent.getByToken(inputTokenDigi_, CastorDigi);
142  if (CastorDigi.isValid())
143  nDigi = CastorDigi->size();
144  else
145  digiOK_ = false;
146 
148  iEvent.getByToken(inputTokenRecHitCASTOR_, CastorHits);
149  if (CastorHits.isValid())
150  nrecHits = CastorHits->size();
151  else
152  rechitOK_ = false;
153 
155  iEvent.getByToken(inputTokenCastorTowers_, castorTowers);
156  if (castorTowers.isValid())
157  nTowers = castorTowers->size();
158  else
159  towerOK_ = false;
160 
162  iEvent.getByToken(JetAlgorithm, jets);
163  if (jets.isValid())
164  nJets = jets->size();
165  else
166  jetsOK_ = false;
167 
168  if (fVerbosity > 0)
169  LogPrint("CastorMonitorModule") << "CastorProductValid(size): RawDataValid=" << RawData.isValid()
170  << " Digi=" << digiOK_ << "(" << nDigi << ") Hits=" << rechitOK_ << "(" << nrecHits
171  << ")"
172  << " Towers=" << towerOK_ << "(" << nTowers << ")"
173  << " Jets=" << jetsOK_ << "(" << nJets << ")";
174 
175  CastorEventProduct->Fill(0, fedsUnpacked / 3.);
176  CastorEventProduct->Fill(1, rawOK_);
177  CastorEventProduct->Fill(2, digiOK_);
178  CastorEventProduct->Fill(3, rechitOK_);
179  CastorEventProduct->Fill(4, towerOK_);
180  CastorEventProduct->Fill(5, jetsOK_);
181 
182  if (digiOK_) {
184  DigiMon_->processEvent(iEvent, *CastorDigi, *TrigResults, conditions);
185  }
186  if (showTiming_) {
187  cpu_timer.stop();
188  if (DigiMon_ != nullptr)
189  std::cout << "TIMER:: DIGI MONITOR ->" << cpu_timer.cpuTime() << std::endl;
190  cpu_timer.reset();
191  cpu_timer.start();
192  }
193 
194  if (rechitOK_)
195  RecHitMon_->processEvent(*CastorHits);
196  if (showTiming_) {
197  cpu_timer.stop();
198  if (RecHitMon_ != nullptr)
199  std::cout << "TIMER:: RECHIT MONITOR->" << cpu_timer.cpuTime() << std::endl;
200  cpu_timer.reset();
201  cpu_timer.start();
202  }
203 
204  if (towerOK_)
205  RecHitMon_->processEventTowers(*castorTowers);
206  if (jetsOK_)
207  RecHitMon_->processEventJets(*jets);
208 
209  if (fVerbosity > 0 && ievt_ % 100 == 0)
210  LogPrint("CastorMonitorModule") << "processed " << ievt_ << " events";
211  return;
212 }
213 
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< edm::TriggerResults > tokenTriggerResults
std::unique_ptr< CastorLEDMonitor > LedMon_
void start()
Definition: CPUTimer.cc:68
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
edm::EDGetTokenT< FEDRawDataCollection > inputTokenRaw_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
MonitorElement * hunpkrep
void reset()
Definition: CPUTimer.cc:99
edm::ESGetToken< CastorDbService, CastorDbRecord > castorDbServiceToken_
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
edm::EDGetTokenT< CastorDigiCollection > inputTokenDigi_
void Fill(long long x)
edm::EDGetTokenT< CastorTowerCollection > inputTokenCastorTowers_
bool getData(T &iHolder) const
Definition: EventSetup.h:122
tuple report
Definition: zeeHLT_cff.py:9
MonitorElement * book1DD(TString const &name, TString const &title, int nchX, double lowX, double highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:155
int iEvent
Definition: GenABIO.cc:224
MonitorElement * bookProfile(TString const &name, TString const &title, int nchX, double lowX, double highX, int, double lowY, double highY, char const *option="s", FUNC onbooking=NOOP())
Definition: DQMStore.h:399
CastorMonitorModule(const edm::ParameterSet &ps)
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
std::unique_ptr< CastorDigiMonitor > DigiMon_
vector< PseudoJet > jets
Times stop()
Definition: CPUTimer.cc:87
bool isValid() const
Definition: HandleBase.h:70
virtual 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)
Log< level::Warning, true > LogPrint
double cpuTime() const
Definition: CPUTimer.cc:146
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void analyze(const edm::Event &iEvent, const edm::EventSetup &) override
edm::EDGetTokenT< BasicJetCollection > JetAlgorithm
MonitorElement * CastorEventProduct
void bookHistograms(DQMStore::IBooker &, edm::Run const &, const edm::EventSetup &) override
void dqmEndRun(const edm::Run &run, const edm::EventSetup &) override
edm::EDGetTokenT< HcalUnpackerReport > inputTokenReport_
std::unique_ptr< CastorRecHitMonitor > RecHitMon_
tuple cout
Definition: gather_cfg.py:144
UseScope< MonitorElementData::Scope::RUN > UseRunScope
Definition: DQMStore.h:541
edm::EDGetTokenT< CastorRecHitCollection > inputTokenRecHitCASTOR_
Definition: Run.h:45
virtual void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)