Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:07:10

0001 /*
0002  * \file DTScalerInfoTask.cc
0003  *
0004  * \author C. Battilana - CIEMAT
0005  *
0006 */
0007 
0008 #include "DQM/DTMonitorModule/src/DTScalerInfoTask.h"
0009 
0010 // Framework
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 
0013 // DT DQM
0014 #include "DQM/DTMonitorModule/interface/DTTimeEvolutionHisto.h"
0015 
0016 #include <sstream>
0017 #include <iostream>
0018 #include <fstream>
0019 
0020 using namespace edm;
0021 using namespace std;
0022 
0023 DTScalerInfoTask::DTScalerInfoTask(const edm::ParameterSet& ps) : nEvents(0) {
0024   LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask") << "[DTScalerInfoTask]: Constructor" << endl;
0025 
0026   scalerToken_ = consumes<LumiScalersCollection>(ps.getUntrackedParameter<InputTag>("inputTagScaler"));
0027   theParams = ps;
0028 }
0029 
0030 DTScalerInfoTask::~DTScalerInfoTask() {
0031   LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask") << "[DTScalerInfoTask]: analyzed " << nEvents << " events" << endl;
0032 }
0033 
0034 void DTScalerInfoTask::dqmBeginRun(const edm::Run& run, const edm::EventSetup& context) {
0035   LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask") << "[DTScalerInfoTask]: BeginRun" << endl;
0036 }
0037 
0038 void DTScalerInfoTask::beginLuminosityBlock(const LuminosityBlock& lumiSeg, const EventSetup& context) {
0039   nEventsInLS = 0;
0040 
0041   LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask") << "[DTScalerInfoTask]: Begin of LS transition" << endl;
0042 }
0043 
0044 void DTScalerInfoTask::endLuminosityBlock(const LuminosityBlock& lumiSeg, const EventSetup& context) {
0045   LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask") << "[DTScalerInfoTask]: End of LS transition" << endl;
0046 
0047   int block = lumiSeg.luminosityBlock();
0048 
0049   map<string, DTTimeEvolutionHisto*>::const_iterator histoIt = trendHistos.begin();
0050   map<string, DTTimeEvolutionHisto*>::const_iterator histoEnd = trendHistos.end();
0051   for (; histoIt != histoEnd; ++histoIt) {
0052     histoIt->second->updateTimeSlot(block, nEventsInLS);
0053   }
0054 }
0055 
0056 void DTScalerInfoTask::analyze(const edm::Event& e, const edm::EventSetup& c) {
0057   nEvents++;
0058   nEventsInLS++;
0059   nEventMonitor->Fill(nEvents);
0060 
0061   //retrieve the luminosity
0062   edm::Handle<LumiScalersCollection> lumiScalers;
0063   if (e.getByToken(scalerToken_, lumiScalers)) {
0064     if (lumiScalers->begin() != lumiScalers->end()) {
0065       LumiScalersCollection::const_iterator lumiIt = lumiScalers->begin();
0066       trendHistos["AvgLumivsLumiSec"]->accumulateValueTimeSlot(lumiIt->instantLumi());
0067     } else {
0068       LogVerbatim("DTDQM|DTMonitorModule|DTScalerInfoTask")
0069           << "[DTScalerInfoTask]: LumiScalersCollection size == 0" << endl;
0070     }
0071   } else {
0072     LogVerbatim("DTDQM|DTMonitorModule|DTScalerInfoTask")
0073         << "[DTScalerInfoTask]: LumiScalersCollection getByToken call failed" << endl;
0074   }
0075 }
0076 
0077 void DTScalerInfoTask::bookHistograms(DQMStore::IBooker& ibooker,
0078                                       edm::Run const& iRun,
0079                                       edm::EventSetup const& context) {
0080   ibooker.setCurrentFolder("DT/EventInfo/Counters");
0081   nEventMonitor = ibooker.bookFloat("nProcessedEventsScalerInfo");
0082 
0083   ibooker.setCurrentFolder("DT/00-DataIntegrity/ScalerInfo");
0084 
0085   string histoName = "AvgLumivsLumiSec";
0086   string histoTitle = "Average Lumi vs LumiSec";
0087   trendHistos[histoName] = new DTTimeEvolutionHisto(ibooker, histoName, histoTitle, 200, 10, true, 0);
0088 }
0089 
0090 // Local Variables:
0091 // show-trailing-whitespace: t
0092 // truncate-lines: t
0093 // End: