CMS 3D CMS Logo

FastTimerServiceClient.cc
Go to the documentation of this file.
1 // C++ headers
2 #include <string>
3 #include <cstring>
4 
5 // boost headers
6 #include <boost/regex.hpp>
7 
8 // Root headers
9 #include <TH1F.h>
10 
11 // CMSSW headers
25 
26 struct MEPSet {
29  int nbins;
30  double xmin;
31  double xmax;
32 };
33 
35 public:
36  explicit FastTimerServiceClient(edm::ParameterSet const &);
37  ~FastTimerServiceClient() override;
38 
39  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
40  static void fillLumiMePSetDescription(edm::ParameterSetDescription & pset);
41  static void fillPUMePSetDescription(edm::ParameterSetDescription & pset);
42 
43 private:
45 
46  void dqmEndLuminosityBlock(DQMStore::IBooker & booker, DQMStore::IGetter & getter, edm::LuminosityBlock const &, edm::EventSetup const&) override;
47  void dqmEndJob(DQMStore::IBooker & booker, DQMStore::IGetter & getter) override;
48 
49  void fillSummaryPlots( DQMStore::IBooker & booker, DQMStore::IGetter & getter);
50  void fillProcessSummaryPlots( DQMStore::IBooker & booker, DQMStore::IGetter & getter, std::string const & path);
51  void fillPathSummaryPlots( DQMStore::IBooker & booker, DQMStore::IGetter & getter, double events, std::string const & path);
52  void fillPlotsVsLumi(DQMStore::IBooker & booker, DQMStore::IGetter & getter, std::string const & current_path, std::string const & suffix, MEPSet pset);
53 
54  static MEPSet getHistoPSet (edm::ParameterSet pset);
55 
59 
63 
64 };
65 
66 
68  m_dqm_path( config.getUntrackedParameter<std::string>( "dqmPath" ) )
69  , doPlotsVsScalLumi_ ( config.getParameter<bool>( "doPlotsVsScalLumi" ) )
70  , doPlotsVsPixelLumi_( config.getParameter<bool>( "doPlotsVsPixelLumi" ) )
71  , doPlotsVsPU_ ( config.getParameter<bool>( "doPlotsVsPU" ) )
72  , scalLumiMEPSet_ ( doPlotsVsScalLumi_ ? getHistoPSet(config.getParameter<edm::ParameterSet>("scalLumiME") ) : MEPSet{} )
73  , pixelLumiMEPSet_( doPlotsVsPixelLumi_ ? getHistoPSet(config.getParameter<edm::ParameterSet>("pixelLumiME")) : MEPSet{} )
74  , puMEPSet_ ( doPlotsVsPU_ ? getHistoPSet(config.getParameter<edm::ParameterSet>("puME") ) : MEPSet{} )
75 {
76 }
77 
79 {
80 }
81 
82 void
84 {
85  fillSummaryPlots(booker, getter);
86 }
87 
88 void
90 {
91  fillSummaryPlots(booker, getter);
92 }
93 
94 void
96 {
97  if (getter.get(m_dqm_path + "/event time_real")) {
98  // the plots are directly in the configured folder
99  fillProcessSummaryPlots(booker, getter, m_dqm_path);
100  } else {
101  static const boost::regex running_n_processes(".*/Running [0-9]+ processes");
102 
104  std::vector<std::string> subdirs = getter.getSubdirs();
105  for (auto const & subdir: subdirs) {
106 
107  // the plots are in a per-number-of-processes folder
108  if (boost::regex_match(subdir, running_n_processes)) {
109 
110  booker.setCurrentFolder(subdir);
111  if ( getter.get(subdir + "/event time_real") )
112  fillProcessSummaryPlots(booker, getter, subdir);
113 
114  std::vector<std::string> subsubdirs = getter.getSubdirs();
115  for (auto const & subsubdir: subsubdirs) {
116  if ( getter.get(subsubdir + "/event time_real") )
117  fillProcessSummaryPlots(booker, getter, subsubdir);
118  }
119 
120  }
121  } // loop on subdirs
122  }
123 }
124 
125 
126 
127 void
129 
130  MonitorElement * me = getter.get(current_path + "/event time_real");
131  if (me == nullptr)
132  // no FastTimerService DQM information
133  return;
134 
135  if ( doPlotsVsScalLumi_ )
136  fillPlotsVsLumi( booker,getter, current_path, "VsScalLumi", scalLumiMEPSet_ );
137  if ( doPlotsVsPixelLumi_ )
138  fillPlotsVsLumi( booker,getter, current_path, "VsPixelLumi", pixelLumiMEPSet_ );
139  if ( doPlotsVsPU_ )
140  fillPlotsVsLumi( booker,getter, current_path, "VsPU", puMEPSet_ );
141 
142  // getter.setCurrentFolder(current_path);
143 
144  double events = me->getTH1F()->GetEntries();
145 
146  // look for per-process directories
147  static const boost::regex process_name(".*/process .*");
148 
149  booker.setCurrentFolder(current_path); // ?!?!?
150  std::vector<std::string> subdirs = getter.getSubdirs();
151  for (auto const & subdir: subdirs) {
152 
153  if (boost::regex_match(subdir, process_name)) {
154 
155  getter.setCurrentFolder(subdir);
156  // look for per-path plots inside each per-process directory
157  std::vector<std::string> subsubdirs = getter.getSubdirs();
158  for (auto const & subsubdir: subsubdirs) {
159  if ( getter.get(subsubdir + "/path time_real") ) {
160  fillPathSummaryPlots(booker, getter, events, subdir);
161  break;
162  }
163  }
164 
165  }
166  } // loop on subdir
167 
168 }
169 
170 void
172  // note: the following checks need to be kept separate, as any of these histograms might be missing
173 
174  booker.setCurrentFolder(current_path);
175  std::vector<std::string> subsubdirs = getter.getSubdirs();
176  size_t npaths = subsubdirs.size();
177 
178  MonitorElement* paths_time = booker.book1D("paths_time_real", "Total (real) time spent in each path", npaths, -0.5, double(npaths)-0.5);
179  MonitorElement* paths_thread = booker.book1D("paths_time_thread","Total (thread) time spent in each path", npaths, -0.5, double(npaths)-0.5);
180  MonitorElement* paths_allocated = booker.book1D("paths_allocated", "Total allocated memory in each path", npaths, -0.5, double(npaths)-0.5);
181  MonitorElement* paths_deallocated = booker.book1D("paths_deallocated","Total deallocated in each path", npaths, -0.5, double(npaths)-0.5);
182 
183  MonitorElement * me;
184  double mean = -1.;
185 
186  // extract the list of Paths and EndPaths from the summary plots
187  int ibin = 1;
188  for (auto const & subsubdir: subsubdirs) {
189 
190  std::string test = "/path ";
191  if ( subsubdir.find(test)==std::string::npos ) continue;
192 
193  static const boost::regex prefix(current_path + "/path ");
194  std::string path = boost::regex_replace(subsubdir,prefix,"");
195 
196  paths_time ->getTH1F()->GetXaxis()->SetBinLabel(ibin,path.c_str());
197  paths_thread ->getTH1F()->GetXaxis()->SetBinLabel(ibin,path.c_str());
198  paths_allocated ->getTH1F()->GetXaxis()->SetBinLabel(ibin,path.c_str());
199  paths_deallocated->getTH1F()->GetXaxis()->SetBinLabel(ibin,path.c_str());
200 
201 
202  if ( (me = getter.get(subsubdir + "/path time_real")) ) {
203  mean = me->getMean();
204  paths_time->getTH1F()->SetBinContent(ibin,mean);
205  }
206  if ( (me = getter.get(subsubdir + "/path time_thread")) ) {
207  mean = me->getMean();
208  paths_thread->getTH1F()->SetBinContent(ibin,mean);
209  }
210  if ( (me = getter.get(subsubdir + "/path allocated")) ) {
211  mean = me->getMean();
212  paths_allocated->getTH1F()->SetBinContent(ibin,mean);
213  }
214 
215  if ( (me = getter.get(subsubdir + "/path deallocated")) ) {
216  mean = me->getMean();
217  paths_deallocated->getTH1F()->SetBinContent(ibin,mean);
218  }
219 
220  ibin++;
221  }
222 
223 
224  for (auto const & subsubdir: subsubdirs) {
225  // for each path, fill histograms with
226  // - the average time spent in each module (total time spent in that module, averaged over all events)
227  // - the running time spent in each module (total time spent in that module, averaged over the events where that module actually ran)
228  // - the "efficiency" of each module (number of time a module succeded divided by the number of times the has run)
229 
230  getter.setCurrentFolder(subsubdir);
231  std::vector<std::string> allmenames = getter.getMEs();
232  if ( allmenames.empty() ) continue;
233 
234  MonitorElement * me_counter = getter.get( subsubdir + "/module_counter" );
235  MonitorElement * me_real_total = getter.get( subsubdir + "/module_time_real_total" );
236  MonitorElement * me_thread_total = getter.get( subsubdir + "/module_time_thread_total" );
237 
238  if (me_counter == nullptr or me_real_total == nullptr)
239  continue;
240 
241  TH1D * counter = me_counter ->getTH1D();
242  TH1D * real_total = me_real_total->getTH1D();
243  TH1D * thread_total = me_thread_total->getTH1D();
244  uint32_t bins = counter->GetXaxis()->GetNbins() - 1;
245  double min = counter->GetXaxis()->GetXmin();
246  double max = counter->GetXaxis()->GetXmax() - 1;
247 
248  TH1F * real_average;
249  TH1F * real_running;
250  TH1F * thread_average;
251  TH1F * thread_running;
252  TH1F * efficiency;
253  MonitorElement * me;
254 
255  booker.setCurrentFolder(subsubdir);
256  me = getter.get( "module_time_real_average" );
257  if (me) {
258  real_average = me->getTH1F();
259  assert( me->getTH1F()->GetXaxis()->GetXmin() == min );
260  assert( me->getTH1F()->GetXaxis()->GetXmax() == max );
261  real_average->Reset();
262  } else {
263  real_average = booker.book1D("module_time_real_average", "module real average timing", bins, min, max)->getTH1F();
264  real_average->SetYTitle("average processing (real) time [ms]");
265  for (uint32_t i = 1; i <= bins; ++i) {
266  const char * module = counter->GetXaxis()->GetBinLabel(i);
267  real_average->GetXaxis()->SetBinLabel(i, module);
268  }
269  }
270 
271  me = getter.get( "module_time_thread_average" );
272  if (me) {
273  thread_average = me->getTH1F();
274  assert( me->getTH1F()->GetXaxis()->GetXmin() == min );
275  assert( me->getTH1F()->GetXaxis()->GetXmax() == max );
276  thread_average->Reset();
277  } else {
278  thread_average = booker.book1D("module_time_thread_average", "module thread average timing", bins, min, max)->getTH1F();
279  thread_average->SetYTitle("average processing (thread) time [ms]");
280  for (uint32_t i = 1; i <= bins; ++i) {
281  const char * module = counter->GetXaxis()->GetBinLabel(i);
282  thread_average->GetXaxis()->SetBinLabel(i, module);
283  }
284  }
285 
286  me = getter.get( "module_time_real_running" );
287  if (me) {
288  real_running = me->getTH1F();
289  assert( me->getTH1F()->GetXaxis()->GetXmin() == min );
290  assert( me->getTH1F()->GetXaxis()->GetXmax() == max );
291  real_running->Reset();
292  } else {
293  real_running = booker.book1D("module_time_real_running", "module real running timing", bins, min, max)->getTH1F();
294  real_running->SetYTitle("running processing (real) time [ms]");
295  for (uint32_t i = 1; i <= bins; ++i) {
296  const char * module = counter->GetXaxis()->GetBinLabel(i);
297  real_running->GetXaxis()->SetBinLabel(i, module);
298  }
299  }
300 
301  me = getter.get( "module_time_thread_running" );
302  if (me) {
303  thread_running = me->getTH1F();
304  assert( me->getTH1F()->GetXaxis()->GetXmin() == min );
305  assert( me->getTH1F()->GetXaxis()->GetXmax() == max );
306  thread_running->Reset();
307  } else {
308  thread_running = booker.book1D("module_time_thread_running", "module thread running timing", bins, min, max)->getTH1F();
309  thread_running->SetYTitle("running processing (thread) time [ms]");
310  for (uint32_t i = 1; i <= bins; ++i) {
311  const char * module = counter->GetXaxis()->GetBinLabel(i);
312  thread_running->GetXaxis()->SetBinLabel(i, module);
313  }
314  }
315 
316  me = getter.get( "module_efficiency" );
317  if (me) {
318  efficiency = me->getTH1F();
319  assert( me->getTH1F()->GetXaxis()->GetXmin() == min );
320  assert( me->getTH1F()->GetXaxis()->GetXmax() == max );
321  efficiency->Reset();
322  } else {
323  efficiency = booker.book1D("module_efficiency", "module efficiency", bins, min, max)->getTH1F();
324  efficiency->SetYTitle("filter efficiency");
325  efficiency->SetMaximum(1.05);
326  for (uint32_t i = 1; i <= bins; ++i) {
327  const char * module = counter->GetXaxis()->GetBinLabel(i);
328  efficiency->GetXaxis()->SetBinLabel(i, module);
329  }
330  }
331 
332 
333  for (uint32_t i = 1; i <= bins; ++i) {
334  double n = counter ->GetBinContent(i);
335  double p = counter ->GetBinContent(i+1);
336  if (n)
337  efficiency->SetBinContent(i, p / n);
338 
339  // real timing
340  double t = real_total->GetBinContent(i);
341  real_average->SetBinContent(i, t / events);
342  if (n)
343  real_running->SetBinContent(i, t / n);
344 
345  // thread timing
346  t = thread_total->GetBinContent(i);
347  thread_average->SetBinContent(i, t / events);
348  if (n)
349  thread_running->SetBinContent(i, t / n);
350  }
351 
352  // vs lumi
353  if ( doPlotsVsScalLumi_ )
354  fillPlotsVsLumi( booker,getter, subsubdir, "VsScalLumi", scalLumiMEPSet_ );
355  if ( doPlotsVsPixelLumi_ )
356  fillPlotsVsLumi( booker,getter, subsubdir, "VsPixelLumi", pixelLumiMEPSet_ );
357  if ( doPlotsVsPU_ )
358  fillPlotsVsLumi( booker,getter, subsubdir, "VsPU", puMEPSet_ );
359 
360  }
361 
362 }
363 
365 void
367 
368  std::vector<std::string> menames;
369 
370  static const boost::regex byls(".*byls");
371  static const boost::regex test(suffix);
372  // get all MEs in the current_path
373  getter.setCurrentFolder(current_path);
374  std::vector<std::string> allmenames = getter.getMEs();
375  for ( const auto & m : allmenames ) {
376  // get only MEs vs LS
377  if (boost::regex_match(m, byls))
378  menames.push_back(m);
379  }
380  // if no MEs available, return
381  if ( menames.empty() )
382  return;
383 
384  // get info for getting the lumi VS LS histogram
385  std::string folder = pset.folder;
386  std::string name = pset.name;
387  int nbins = pset.nbins;
388  double xmin = pset.xmin;
389  double xmax = pset.xmax;
390 
391  // get lumi/PU VS LS ME
392  getter.setCurrentFolder(folder);
393  MonitorElement* lumiVsLS = getter.get(folder+"/"+name);
394  // if no ME available, return
395  if ( !lumiVsLS ) {
396  edm::LogWarning("FastTimerServiceClient") << "no " << name << " ME is available in " << folder << std::endl;
397  return;
398  }
399 
400  // get range and binning for new MEs x-axis
401  size_t size = lumiVsLS->getTProfile()->GetXaxis()->GetNbins();
402  std::string xtitle = lumiVsLS->getTProfile()->GetYaxis()->GetTitle();
403 
404  std::vector<double> lumi;
405  std::vector<int> LS;
406  for ( size_t ibin=1; ibin <= size; ++ibin ) {
407  // avoid to store points w/ no info
408  if ( lumiVsLS->getTProfile()->GetBinContent(ibin) == 0. ) continue;
409 
410  lumi.push_back( lumiVsLS->getTProfile()->GetBinContent(ibin) );
411  LS.push_back ( lumiVsLS->getTProfile()->GetXaxis()->GetBinCenter(ibin) );
412  }
413 
414  booker.setCurrentFolder(current_path);
415  getter.setCurrentFolder(current_path);
416  for ( auto m : menames ) {
417  std::string label = m;
418  label.erase(label.find("_byls"));
419 
420  MonitorElement* me = getter.get(current_path + "/" + m);
421  float ymin = 0.;
423  std::string ytitle = me->getTProfile()->GetYaxis()->GetTitle();
424 
425  MonitorElement* meVsLumi = getter.get( current_path + "/" + label + "_" + suffix );
426  if (meVsLumi) {
427  assert( meVsLumi->getTProfile()->GetXaxis()->GetXmin() == xmin );
428  assert( meVsLumi->getTProfile()->GetXaxis()->GetXmax() == xmax );
429  meVsLumi->Reset(); // do I have to do it ?!?!?
430  } else {
431  meVsLumi = booker.bookProfile(label + "_" + suffix, label + "_" + suffix, nbins, xmin, xmax, ymin, ymax);
432  // TProfile* meVsLumi_p = meVsLumi->getTProfile();
433  meVsLumi->getTProfile()->GetXaxis()->SetTitle(xtitle.c_str());
434  meVsLumi->getTProfile()->GetYaxis()->SetTitle(ytitle.c_str());
435  }
436  for ( size_t ils=0; ils < LS.size(); ++ils ) {
437  int ibin = me->getTProfile()->GetXaxis()->FindBin(LS[ils]);
438  double y = me->getTProfile()->GetBinContent(ibin);
439 
440  meVsLumi->Fill(lumi[ils],y);
441  }
442  }
443 
444 
445 }
446 
447 void
449  pset.add<std::string>("folder", "HLT/LumiMonitoring");
450  pset.add<std::string>("name" , "lumiVsLS");
451  pset.add<int> ("nbins", 440 );
452  pset.add<double>("xmin", 0.);
453  pset.add<double>("xmax", 22000.);
454 }
455 
456 
457 void
459  pset.add<std::string>("folder", "HLT/LumiMonitoring");
460  pset.add<std::string>("name" , "puVsLS");
461  pset.add<int> ("nbins", 260 );
462  pset.add<double>("xmin", 0.);
463  pset.add<double>("xmax", 130.);
464 }
465 
466 
468 {
469  return MEPSet{
470  pset.getParameter<std::string>("folder"),
471  pset.getParameter<std::string>("name"),
472  pset.getParameter<int>("nbins"),
473  pset.getParameter<double>("xmin"),
474  pset.getParameter<double>("xmax"),
475  };
476 }
477 
478 void
480  // The following says we do not know what parameters are allowed so do no validation
481  // Please change this to state exactly what you do use, even if it is no parameters
483  desc.addUntracked<std::string>( "dqmPath", "HLT/TimerService");
484  desc.add<bool>( "doPlotsVsScalLumi", true );
485  desc.add<bool>( "doPlotsVsPixelLumi", false );
486  desc.add<bool>( "doPlotsVsPU", true );
487 
488  edm::ParameterSetDescription scalLumiMEPSet;
489  fillLumiMePSetDescription(scalLumiMEPSet);
490  desc.add<edm::ParameterSetDescription>("scalLumiME", scalLumiMEPSet);
491 
492  edm::ParameterSetDescription pixelLumiMEPSet;
493  fillLumiMePSetDescription(pixelLumiMEPSet);
494  desc.add<edm::ParameterSetDescription>("pixelLumiME", pixelLumiMEPSet);
495 
497  fillPUMePSetDescription(puMEPSet);
498  desc.add<edm::ParameterSetDescription>("puME", puMEPSet);
499 
500  descriptions.add("fastTimerServiceClient", desc);
501 }
502 
503 //define this as a plug-in
size
Write out results.
T getParameter(std::string const &) const
static void fillPUMePSetDescription(edm::ParameterSetDescription &pset)
static MEPSet getHistoPSet(edm::ParameterSet pset)
void dqmEndLuminosityBlock(DQMStore::IBooker &booker, DQMStore::IGetter &getter, edm::LuminosityBlock const &, edm::EventSetup const &) override
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:157
void fillPlotsVsLumi(DQMStore::IBooker &booker, DQMStore::IGetter &getter, std::string const &current_path, std::string const &suffix, MEPSet pset)
MonitorElement * get(const std::string &path)
Definition: DQMStore.cc:305
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::string folder
void fillProcessSummaryPlots(DQMStore::IBooker &booker, DQMStore::IGetter &getter, std::string const &path)
TH1D * getTH1D(void) const
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
Definition: config.py:1
double getMean(int axis=1) const
get mean value of histogram along x, y or z axis (axis=1, 2, 3 respectively)
void Fill(long long x)
void dqmEndJob(DQMStore::IBooker &booker, DQMStore::IGetter &getter) override
void fillPathSummaryPlots(DQMStore::IBooker &booker, DQMStore::IGetter &getter, double events, std::string const &path)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
static void fillLumiMePSetDescription(edm::ParameterSetDescription &pset)
T min(T a, T b)
Definition: MathUtil.h:58
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::vector< std::string > getMEs(void)
Definition: DQMStore.cc:327
FastTimerServiceClient(edm::ParameterSet const &)
std::string name
void fillSummaryPlots(DQMStore::IBooker &booker, DQMStore::IGetter &getter)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
TH1F * getTH1F(void) const
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< std::string > getSubdirs(void)
Definition: DQMStore.cc:323
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:347
HLT enums.
TProfile * getTProfile(void) const
void Reset(void)
reset ME (ie. contents, errors, etc)
Definition: vlib.h:208