CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Private Attributes | Friends
DTScalerInfoTask Class Reference

#include <DTScalerInfoTask.h>

Inheritance diagram for DTScalerInfoTask:
one::DQMEDAnalyzer< edm::one::WatchLuminosityBlocks > one::dqmimplementation::DQMBaseClass< T... >

Public Member Functions

 DTScalerInfoTask (const edm::ParameterSet &ps)
 Constructor. More...
 
 ~DTScalerInfoTask () override
 Destructor. More...
 
- Public Member Functions inherited from one::DQMEDAnalyzer< edm::one::WatchLuminosityBlocks >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Protected Member Functions

void analyze (const edm::Event &e, const edm::EventSetup &c) override
 Analyze. More...
 
void beginLuminosityBlock (const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &context) override
 To reset the MEs. More...
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
void dqmBeginRun (const edm::Run &, const edm::EventSetup &) override
 Beginrun. More...
 
void endLuminosityBlock (const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &context) override
 Perform trend plot operations. More...
 

Private Attributes

MonitorElementnEventMonitor
 
int nEvents
 
int nEventsInLS
 
edm::EDGetTokenT< LumiScalersCollectionscalerToken_
 
edm::ParameterSet theParams
 
std::map< std::string,DTTimeEvolutionHisto * > trendHistos
 

Friends

class DTMonitorModule
 

Detailed Description

Definition at line 38 of file DTScalerInfoTask.h.

Constructor & Destructor Documentation

DTScalerInfoTask::DTScalerInfoTask ( const edm::ParameterSet ps)

Constructor.

Definition at line 23 of file DTScalerInfoTask.cc.

References edm::ParameterSet::getUntrackedParameter(), LogTrace, scalerToken_, and theParams.

23  :
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 }
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< LumiScalersCollection > scalerToken_
#define LogTrace(id)
edm::ParameterSet theParams
DTScalerInfoTask::~DTScalerInfoTask ( )
override

Destructor.

Definition at line 35 of file DTScalerInfoTask.cc.

References LogTrace, and nEvents.

35  {
36 
37  LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask")
38  << "[DTScalerInfoTask]: analyzed " << nEvents << " events" << endl;
39 
40 }
#define LogTrace(id)

Member Function Documentation

void DTScalerInfoTask::analyze ( const edm::Event e,
const edm::EventSetup c 
)
overrideprotected

Analyze.

Definition at line 75 of file DTScalerInfoTask.cc.

References MonitorElement::Fill(), edm::Event::getByToken(), nEventMonitor, nEvents, nEventsInLS, scalerToken_, and trendHistos.

75  {
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 }
MonitorElement * nEventMonitor
edm::EDGetTokenT< LumiScalersCollection > scalerToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
void Fill(long long x)
std::map< std::string,DTTimeEvolutionHisto * > trendHistos
void DTScalerInfoTask::beginLuminosityBlock ( const edm::LuminosityBlock lumiSeg,
const edm::EventSetup context 
)
overrideprotected

To reset the MEs.

Definition at line 50 of file DTScalerInfoTask.cc.

References LogTrace, and nEventsInLS.

50  {
51 
52  nEventsInLS=0;
53 
54  LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask")
55  << "[DTScalerInfoTask]: Begin of LS transition" << endl;
56 
57  }
#define LogTrace(id)
void DTScalerInfoTask::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  iRun,
edm::EventSetup const &  context 
)
overrideprotected

Definition at line 100 of file DTScalerInfoTask.cc.

References DQMStore::IBooker::bookFloat(), nEventMonitor, DQMStore::IBooker::setCurrentFolder(), and trendHistos.

100  {
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 }
MonitorElement * nEventMonitor
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
std::map< std::string,DTTimeEvolutionHisto * > trendHistos
MonitorElement * bookFloat(Args &&...args)
Definition: DQMStore.h:105
void DTScalerInfoTask::dqmBeginRun ( const edm::Run run,
const edm::EventSetup context 
)
overrideprotected

Beginrun.

Definition at line 43 of file DTScalerInfoTask.cc.

References LogTrace.

43  {
44 
45  LogTrace("DTDQM|DTMonitorModule|DTScalerInfoTask")
46  << "[DTScalerInfoTask]: BeginRun" << endl;
47 }
#define LogTrace(id)
void DTScalerInfoTask::endLuminosityBlock ( const edm::LuminosityBlock lumiSeg,
const edm::EventSetup context 
)
overrideprotected

Perform trend plot operations.

Definition at line 59 of file DTScalerInfoTask.cc.

References groupFilesInBlocks::block, LogTrace, edm::LuminosityBlockBase::luminosityBlock(), nEventsInLS, and trendHistos.

59  {
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 }
LuminosityBlockNumber_t luminosityBlock() const
#define LogTrace(id)
std::map< std::string,DTTimeEvolutionHisto * > trendHistos

Friends And Related Function Documentation

friend class DTMonitorModule
friend

Definition at line 40 of file DTScalerInfoTask.h.

Member Data Documentation

MonitorElement* DTScalerInfoTask::nEventMonitor
private

Definition at line 77 of file DTScalerInfoTask.h.

Referenced by analyze(), and bookHistograms().

int DTScalerInfoTask::nEvents
private

Definition at line 69 of file DTScalerInfoTask.h.

Referenced by analyze(), looper.Looper::loop(), and ~DTScalerInfoTask().

int DTScalerInfoTask::nEventsInLS
private

Definition at line 70 of file DTScalerInfoTask.h.

Referenced by analyze(), beginLuminosityBlock(), and endLuminosityBlock().

edm::EDGetTokenT<LumiScalersCollection> DTScalerInfoTask::scalerToken_
private

Definition at line 74 of file DTScalerInfoTask.h.

Referenced by analyze(), and DTScalerInfoTask().

edm::ParameterSet DTScalerInfoTask::theParams
private

Definition at line 72 of file DTScalerInfoTask.h.

Referenced by DTScalerInfoTask().

std::map<std::string ,DTTimeEvolutionHisto* > DTScalerInfoTask::trendHistos
private

Definition at line 76 of file DTScalerInfoTask.h.

Referenced by analyze(), bookHistograms(), and endLuminosityBlock().