CMS 3D CMS Logo

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  nEvents(0) {
25 
26  LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask")
27  << "[DTScalerInfoTask]: Constructor"<<endl;
28 
29  scalerToken_ = consumes<LumiScalersCollection>(
30  ps.getUntrackedParameter<InputTag>("inputTagScaler"));
31  theParams = ps;
32 }
33 
34 
36 
37  LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask")
38  << "[DTScalerInfoTask]: analyzed " << nEvents << " events" << endl;
39 
40 }
41 
42 
44 
45  LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask")
46  << "[DTScalerInfoTask]: BeginRun" << endl;
47 }
48 
49 
51 
52  nEventsInLS=0;
53 
54  LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask")
55  << "[DTScalerInfoTask]: Begin of LS transition" << endl;
56 
57  }
58 
59 void DTScalerInfoTask::endLuminosityBlock(const LuminosityBlock& lumiSeg, const EventSetup& context) {
60 
61  LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask")
62  << "[DTScalerInfoTask]: End of LS transition" << endl;
63 
64 
65  int block = lumiSeg.luminosityBlock();
66 
67  map<string,DTTimeEvolutionHisto* >::const_iterator histoIt = trendHistos.begin();
68  map<string,DTTimeEvolutionHisto* >::const_iterator histoEnd = trendHistos.end();
69  for(;histoIt!=histoEnd;++histoIt) {
70  histoIt->second->updateTimeSlot(block, nEventsInLS);
71  }
72 
73 }
74 
76 
77  nEvents++;
78  nEventsInLS++;
80 
81  //retrieve the luminosity
83  if (e.getByToken(scalerToken_, lumiScalers)) {
84  if (lumiScalers->begin() != lumiScalers->end()) {
85  LumiScalersCollection::const_iterator lumiIt = lumiScalers->begin();
86  trendHistos["AvgLumivsLumiSec"]->accumulateValueTimeSlot(lumiIt->instantLumi());
87  }
88  else {
89  LogVerbatim("DTDQM|DTMonitorModule|DTScalerInfoTask")
90  << "[DTScalerInfoTask]: LumiScalersCollection size == 0" << endl;
91  }
92  }
93  else {
94  LogVerbatim("DTDQM|DTMonitorModule|DTScalerInfoTask")
95  << "[DTScalerInfoTask]: LumiScalersCollection getByToken call failed" << endl;
96  }
97 
98 }
99 
100 void DTScalerInfoTask::bookHistograms(DQMStore::IBooker & ibooker, edm::Run const & iRun, edm::EventSetup const & context) {
101 
102  ibooker.setCurrentFolder("DT/EventInfo/Counters");
103  nEventMonitor = ibooker.bookFloat("nProcessedEventsScalerInfo");
104 
105  ibooker.setCurrentFolder("DT/00-DataIntegrity/ScalerInfo");
106 
107  string histoName = "AvgLumivsLumiSec";
108  string histoTitle = "Average Lumi vs LumiSec";
109  trendHistos[histoName] = new DTTimeEvolutionHisto(ibooker,histoName,histoTitle,200,10,true,0);
110 
111 }
112 
113 // Local Variables:
114 // show-trailing-whitespace: t
115 // truncate-lines: t
116 // End:
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * nEventMonitor
DTScalerInfoTask(const edm::ParameterSet &ps)
Constructor.
edm::EDGetTokenT< LumiScalersCollection > scalerToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
void beginLuminosityBlock(const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &context) override
To reset the MEs.
virtual ~DTScalerInfoTask()
Destructor.
void Fill(long long x)
LuminosityBlockNumber_t luminosityBlock() const
void analyze(const edm::Event &e, const edm::EventSetup &c) override
Analyze.
#define LogTrace(id)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
void endLuminosityBlock(const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &context) override
Perform trend plot operations.
std::map< std::string,DTTimeEvolutionHisto * > trendHistos
HLT enums.
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
MonitorElement * bookFloat(Args &&...args)
Definition: DQMStore.h:109
UInt_t nEvents
Definition: hcalCalib.cc:42
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
Beginrun.
Definition: Run.h:42
edm::ParameterSet theParams