CMS 3D CMS Logo

DQMEventInfo.cc
Go to the documentation of this file.
1 /*
2  * \file DQMEventInfo.cc
3  * \author M. Zanetti - CERN PH
4  * Last Update:
5  *
6  */
7 #include "DQMEventInfo.h"
10 #include <TSystem.h>
11 
12 #include <algorithm>
13 #include <cstdio>
14 #include <sstream>
15 #include <cmath>
16 
17 #include <boost/algorithm/string/join.hpp>
18 
19 
20 static inline double stampToReal(edm::Timestamp time)
21 { return (time.value() >> 32) + 1e-6*(time.value() & 0xffffffff); }
22 
23 static inline double stampToReal(const timeval &time)
24 { return time.tv_sec + 1e-6*time.tv_usec; }
25 
26 
28 
29  struct timeval now;
30  gettimeofday(&now, nullptr);
31 
32  parameters_ = ps;
33  pEvent_ = 0;
34  evtRateCount_ = 0;
36 
37  // read config parms
38  std::string folder = parameters_.getUntrackedParameter<std::string>("eventInfoFolder", "EventInfo") ;
39  subsystemname_ = parameters_.getUntrackedParameter<std::string>("subSystemFolder", "YourSubsystem") ;
40 
41  eventInfoFolder_ = subsystemname_ + "/" + folder ;
42  evtRateWindow_ = parameters_.getUntrackedParameter<double>("eventRateWindow", 0.5);
43  if(evtRateWindow_<=0.15) evtRateWindow_=0.15;
44 
45 }
46 
47 DQMEventInfo::~DQMEventInfo() = default;
48 
50  edm::Run const & iRun,
51  edm::EventSetup const & /* iSetup */)
52 {
54 
55  //Event specific contents
56  runId_ = ibooker.bookInt("iRun");
57  runId_->Fill(iRun.id().run());
58  lumisecId_ = ibooker.bookInt("iLumiSection");
59  lumisecId_->Fill(-1);
60  eventId_ = ibooker.bookInt("iEvent");
61  eventId_->Fill(-1);
62  eventTimeStamp_ = ibooker.bookFloat("eventTimeStamp");
63 
65  //Process specific contents
66  processTimeStamp_ = ibooker.bookFloat("processTimeStamp");
68  processLatency_ = ibooker.bookFloat("processLatency");
70  processEvents_ = ibooker.bookInt("processedEvents");
72  processEventRate_ = ibooker.bookFloat("processEventRate");
74  nUpdates_= ibooker.bookInt("processUpdates");
75  nUpdates_->Fill(-1);
76 
77  //Static Contents
78  processId_= ibooker.bookInt("processID");
79  processId_->Fill(getpid());
80  processStartTimeStamp_ = ibooker.bookFloat("processStartTimeStamp");
82  runStartTimeStamp_ = ibooker.bookFloat("runStartTimeStamp");
84  char hostname[65];
85  gethostname(hostname,64);
86  hostname[64] = 0;
87  hostName_= ibooker.bookString("hostName",hostname);
88  processName_= ibooker.bookString("processName",subsystemname_);
89  char* pwd = getcwd(nullptr, 0);
90  workingDir_= ibooker.bookString("workingDir",pwd);
91  free(pwd);
92  cmsswVer_= ibooker.bookString("CMSSW_Version",edm::getReleaseVersion());
93 
94  // Folder to be populated by sub-systems' code
95  std::string subfolder = eventInfoFolder_ + "/reportSummaryContents" ;
96  ibooker.setCurrentFolder(subfolder);
97 
98  //Online static histograms
99  const edm::ParameterSet &sourcePSet =
101  .getParameterSet("@main_input");
102 
103  if (sourcePSet.getParameter<std::string>("@module_type") == "DQMStreamerReader" ){
104  std::string evSelection;
105  std::vector<std::string> evSelectionList;
106  std::string delimiter( ", " );
107  evSelectionList = sourcePSet.getUntrackedParameter<std::vector<std::string> >("SelectEvents");
108  // add single quotes inline in the vector of HLT paths:
109  // we do copy assignment, and getUntrackedParameter returns
110  // a by-value copy of the vector of strings
111  std::for_each( evSelectionList.begin(), evSelectionList.end(),
112  []( std::string & s ){ std::string squote( "'" );
113  s = squote + s + squote;
114  }
115  );
116  evSelection = boost::algorithm::join( evSelectionList, delimiter );
117  // if no HLT paths are specified, no selections are performed:
118  // we mark this with an asterisk.
119  if( evSelection.empty() ) {
120  evSelection = std::string( "'*'" );
121  }
123  ibooker.bookString("eventSelection",evSelection);
124  }
125 
126 
127 }
128 
129 
131 {
133 }
134 
136 
137  eventId_->Fill(e.id().event()); // Handing edm::EventNumber_t to Fill method which will handle further casting
139 
140  pEvent_++;
141  evtRateCount_++;
143 
144  struct timeval now;
145  gettimeofday(&now, nullptr);
147  currentTime_ = stampToReal(now);
148 
151 
152  double delta = currentTime_ - lastAvgTime_;
153  if (delta >= (evtRateWindow_*60.0))
154  {
156  evtRateCount_ = 0;
157  lastAvgTime_ = currentTime_;
158  }
159 
160  return;
161 }
LuminosityBlockID id() const
dbl * delta
Definition: mlp_gen.cc:36
double evtRateWindow_
Definition: DQMEventInfo.h:64
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
double currentTime_
Definition: DQMEventInfo.h:62
MonitorElement * hostName_
of event processed so far
Definition: DQMEventInfo.h:88
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * workingDir_
DQM "name" of the job (eg, Hcal or DT)
Definition: DQMEventInfo.h:90
double lastUpdateTime_
Definition: DQMEventInfo.h:62
RunID const & id() const
Definition: RunBase.h:39
MonitorElement * eventTimeStamp_
Definition: DQMEventInfo.h:76
RunNumber_t run() const
Definition: RunID.h:39
void analyze(const edm::Event &e, const edm::EventSetup &c) override
Analyze.
edm::ParameterSet parameters_
Definition: DQMEventInfo.h:55
MonitorElement * processId_
Number of collector updates (TBD)
Definition: DQMEventInfo.h:82
MonitorElement * processStartTimeStamp_
The PID associated with this job.
Definition: DQMEventInfo.h:83
MonitorElement * bookInt(Args &&...args)
Definition: DQMStore.h:106
MonitorElement * bookString(Args &&...args)
Definition: DQMStore.h:100
int64_t pEvent_
Definition: DQMEventInfo.h:66
ModuleDescription const & moduleDescription() const
ParameterSet const & getProcessParameterSetContainingModule(ModuleDescription const &moduleDescription)
void Fill(long long x)
MonitorElement * cmsswVer_
Current working directory of the job.
Definition: DQMEventInfo.h:91
MonitorElement * eventId_
UTC time of the run start.
Definition: DQMEventInfo.h:74
std::string eventInfoFolder_
Definition: DQMEventInfo.h:56
MonitorElement * lumisecId_
Definition: DQMEventInfo.h:75
MonitorElement * runStartTimeStamp_
Definition: DQMEventInfo.h:73
DQMEventInfo(const edm::ParameterSet &ps)
Constructor.
Definition: DQMEventInfo.cc:27
double lastAvgTime_
Definition: DQMEventInfo.h:62
MonitorElement * processEvents_
Avg # of events in programmable window (default: 5 min)
Definition: DQMEventInfo.h:87
MonitorElement * processLatency_
The UTC time of the last event.
Definition: DQMEventInfo.h:85
MonitorElement * processName_
Hostname of the local machine.
Definition: DQMEventInfo.h:89
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: DQMEventInfo.cc:49
~DQMEventInfo() override
Destructor.
std::string getReleaseVersion()
int64_t evtRateCount_
Definition: DQMEventInfo.h:65
MonitorElement * runId_
Definition: DQMEventInfo.h:72
static std::string join(char **cmd)
Definition: RemoteFile.cc:18
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
ParameterSet const & getParameterSet(std::string const &) const
MonitorElement * nUpdates_
These MEs are either static or updated upon each analyze() call.
Definition: DQMEventInfo.h:81
LuminosityBlockNumber_t luminosityBlock() const
MonitorElement * processEventRate_
Time elapsed since the last event.
Definition: DQMEventInfo.h:86
Timestamp const & beginTime() const
Definition: RunBase.h:41
edm::EventID id() const
Definition: EventBase.h:60
void beginLuminosityBlock(const edm::LuminosityBlock &l, const edm::EventSetup &c) override
MonitorElement * bookFloat(Args &&...args)
Definition: DQMStore.h:112
MonitorElement * processTimeStamp_
The UTC time of the first event processed.
Definition: DQMEventInfo.h:84
static double stampToReal(edm::Timestamp time)
Definition: DQMEventInfo.cc:20
TimeValue_t value() const
Definition: Timestamp.h:56
std::string subsystemname_
Definition: DQMEventInfo.h:57
edm::Timestamp time() const
Definition: EventBase.h:61
Definition: Run.h:44