CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
DTScalerInfoTask.cc
Go to the documentation of this file.
1 /*
2  * \file DTScalerInfoTask.cc
3  *
4  * \author C. Battilana - CIEMAT
5  *
6 */
7 
9 
10 // Framework
12 
13 // DT DQM
15 
16 #include <sstream>
17 #include <iostream>
18 #include <fstream>
19 
20 using namespace edm;
21 using namespace std;
22 
24  LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask") << "[DTScalerInfoTask]: Constructor" << endl;
25 
26  scalerToken_ = consumes<LumiScalersCollection>(ps.getUntrackedParameter<InputTag>("inputTagScaler"));
27  theParams = ps;
28 }
29 
31  LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask") << "[DTScalerInfoTask]: analyzed " << nEvents << " events" << endl;
32 }
33 
35  LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask") << "[DTScalerInfoTask]: BeginRun" << endl;
36 }
37 
39  nEventsInLS = 0;
40 
41  LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask") << "[DTScalerInfoTask]: Begin of LS transition" << endl;
42 }
43 
45  LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask") << "[DTScalerInfoTask]: End of LS transition" << endl;
46 
47  int block = lumiSeg.luminosityBlock();
48 
49  map<string, DTTimeEvolutionHisto*>::const_iterator histoIt = trendHistos.begin();
50  map<string, DTTimeEvolutionHisto*>::const_iterator histoEnd = trendHistos.end();
51  for (; histoIt != histoEnd; ++histoIt) {
52  histoIt->second->updateTimeSlot(block, nEventsInLS);
53  }
54 }
55 
57  nEvents++;
58  nEventsInLS++;
60 
61  //retrieve the luminosity
63  if (e.getByToken(scalerToken_, lumiScalers)) {
64  if (lumiScalers->begin() != lumiScalers->end()) {
65  LumiScalersCollection::const_iterator lumiIt = lumiScalers->begin();
66  trendHistos["AvgLumivsLumiSec"]->accumulateValueTimeSlot(lumiIt->instantLumi());
67  } else {
68  LogVerbatim("DTDQM|DTMonitorModule|DTScalerInfoTask")
69  << "[DTScalerInfoTask]: LumiScalersCollection size == 0" << endl;
70  }
71  } else {
72  LogVerbatim("DTDQM|DTMonitorModule|DTScalerInfoTask")
73  << "[DTScalerInfoTask]: LumiScalersCollection getByToken call failed" << endl;
74  }
75 }
76 
78  edm::Run const& iRun,
79  edm::EventSetup const& context) {
80  ibooker.setCurrentFolder("DT/EventInfo/Counters");
81  nEventMonitor = ibooker.bookFloat("nProcessedEventsScalerInfo");
82 
83  ibooker.setCurrentFolder("DT/00-DataIntegrity/ScalerInfo");
84 
85  string histoName = "AvgLumivsLumiSec";
86  string histoTitle = "Average Lumi vs LumiSec";
87  trendHistos[histoName] = new DTTimeEvolutionHisto(ibooker, histoName, histoTitle, 200, 10, true, 0);
88 }
89 
90 // Local Variables:
91 // show-trailing-whitespace: t
92 // truncate-lines: t
93 // End:
Log< level::Info, true > LogVerbatim
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * bookFloat(TString const &name, FUNC onbooking=NOOP())
Definition: DQMStore.h:80
const edm::EventSetup & c
MonitorElement * nEventMonitor
DTScalerInfoTask(const edm::ParameterSet &ps)
Constructor.
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
edm::EDGetTokenT< LumiScalersCollection > scalerToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
void beginLuminosityBlock(const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &context) override
To reset the MEs.
#define LogTrace(id)
LuminosityBlockNumber_t luminosityBlock() const
void Fill(long long x)
~DTScalerInfoTask() override
Destructor.
void analyze(const edm::Event &e, const edm::EventSetup &c) override
Analyze.
std::map< std::string, DTTimeEvolutionHisto * > trendHistos
void endLuminosityBlock(const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &context) override
Perform trend plot operations.
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
Beginrun.
Definition: Run.h:45
edm::ParameterSet theParams