Flow123d  JS_before_hm-1011-g79d7237
sys_profiler.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file sys_profiler.cc
15  * @ingroup system
16  * @brief Profiler
17  */
18 
19 
20 // Fat header
21 
22 #include <fstream>
23 #include <iomanip>
24 #include <sys/param.h>
25 #include <unordered_map>
26 
27 #ifdef FLOW123D_HAVE_PYTHON
28  #include "Python.h"
29 #endif // FLOW123D_HAVE_PYTHON
30 
31 #include "sys_profiler.hh"
32 #include "system/system.hh"
33 #include "system/python_loader.hh"
34 #include <iostream>
35 #include <boost/format.hpp>
36 
37 #include "system/file_path.hh"
38 #include "mpi.h"
39 #include "time_point.hh"
40 
41 /*
42  * These should be replaced by using boost MPI interface
43  */
44 int MPI_Functions::sum(int* val, MPI_Comm comm) {
45  int total = 0;
46  MPI_Reduce(val, &total, 1, MPI_INT, MPI_SUM, 0, comm);
47  return total;
48  }
49 
50 double MPI_Functions::sum(double* val, MPI_Comm comm) {
51  double total = 0;
52  MPI_Reduce(val, &total, 1, MPI_DOUBLE, MPI_SUM, 0, comm);
53  return total;
54  }
55 
56 long MPI_Functions::sum(long* val, MPI_Comm comm) {
57  long total = 0;
58  MPI_Reduce(val, &total, 1, MPI_LONG, MPI_SUM, 0, comm);
59  return total;
60  }
61 
62 int MPI_Functions::min(int* val, MPI_Comm comm) {
63  int min = 0;
64  MPI_Reduce(val, &min, 1, MPI_INT, MPI_MIN, 0, comm);
65  return min;
66  }
67 
68 double MPI_Functions::min(double* val, MPI_Comm comm) {
69  double min = 0;
70  MPI_Reduce(val, &min, 1, MPI_DOUBLE, MPI_MIN, 0, comm);
71  return min;
72  }
73 
74 long MPI_Functions::min(long* val, MPI_Comm comm) {
75  long min = 0;
76  MPI_Reduce(val, &min, 1, MPI_LONG, MPI_MIN, 0, comm);
77  return min;
78  }
79 
80 int MPI_Functions::max(int* val, MPI_Comm comm) {
81  int max = 0;
82  MPI_Reduce(val, &max, 1, MPI_INT, MPI_MAX, 0, comm);
83  return max;
84  }
85 
86 double MPI_Functions::max(double* val, MPI_Comm comm) {
87  double max = 0;
88  MPI_Reduce(val, &max, 1, MPI_DOUBLE, MPI_MAX, 0, comm);
89  return max;
90  }
91 
92 long MPI_Functions::max(long* val, MPI_Comm comm) {
93  long max = 0;
94  MPI_Reduce(val, &max, 1, MPI_LONG, MPI_MAX, 0, comm);
95  return max;
96  }
97 
98 
99 #ifdef FLOW123D_DEBUG_PROFILER
100 /*********************************************************************************************
101  * Implementation of class Timer
102  */
103 
104 const int timer_no_child=-1;
105 
106 Timer::Timer(const CodePoint &cp, int parent)
107 : start_time(TimePoint()),
108  cumul_time(0.0),
109  call_count(0),
110  start_count(0),
111  code_point_(&cp),
112  full_hash_(cp.hash_),
113  hash_idx_(cp.hash_idx_),
114  parent_timer(parent),
115  total_allocated_(0),
116  total_deallocated_(0),
117  max_allocated_(0),
118  current_allocated_(0),
119  alloc_called(0),
120  dealloc_called(0)
121 #ifdef FLOW123D_HAVE_PETSC
122 , petsc_start_memory(0),
123  petsc_end_memory (0),
124  petsc_memory_difference(0),
125  petsc_peak_memory(0),
126  petsc_local_peak_memory(0)
127 #endif // FLOW123D_HAVE_PETSC
128 {
129  for(unsigned int i=0; i< max_n_childs ;i++) child_timers[i]=timer_no_child;
130 }
131 
132 
133 /**
134  * Debug information of the timer
135  */
136 ostream & operator <<(ostream& os, const Timer& timer) {
137  os << " Timer: " << timer.tag() << endl;
138  os << " malloc: " << timer.total_allocated_ << endl;
139  os << " dalloc: " << timer.total_deallocated_ << endl;
140  #ifdef FLOW123D_HAVE_PETSC
141  os << " start: " << timer.petsc_start_memory << endl;
142  os << " stop : " << timer.petsc_end_memory << endl;
143  os << " diff : " << timer.petsc_memory_difference << " (" << timer.petsc_end_memory - timer.petsc_start_memory << ")" << endl;
144  os << " peak : " << timer.petsc_peak_memory << " (" << timer.petsc_local_peak_memory << ")" << endl;
145  #endif // FLOW123D_HAVE_PETSC
146  os << endl;
147  return os;
148 }
149 
150 
151 double Timer::cumulative_time() const {
152  return cumul_time;
153 }
154 
155 void Profiler::accept_from_child(Timer &parent, Timer &child) {
156  int child_timer = 0;
157  for (unsigned int i = 0; i < Timer::max_n_childs; i++) {
158  child_timer = child.child_timers[i];
159  if (child_timer != timer_no_child) {
160  // propagate metrics from child to parent
161  accept_from_child(child, timers_[child_timer]);
162  }
163  }
164  // compute totals by adding values from child
165  parent.total_allocated_ += child.total_allocated_;
166  parent.total_deallocated_ += child.total_deallocated_;
167  parent.alloc_called += child.alloc_called;
168  parent.dealloc_called += child.dealloc_called;
169 
170 #ifdef FLOW123D_HAVE_PETSC
171  if (petsc_monitor_memory) {
172  // add differences from child
173  parent.petsc_memory_difference += child.petsc_memory_difference;
174  parent.current_allocated_ += child.current_allocated_;
175 
176  // when computing maximum, we take greater value from parent and child
177  // PetscMemoryGetCurrentUsage always return absolute (not relative) value
178  parent.petsc_peak_memory = max(parent.petsc_peak_memory, child.petsc_peak_memory);
179  }
180 #endif // FLOW123D_HAVE_PETSC
181 
182  parent.max_allocated_ = max(parent.max_allocated_, child.max_allocated_);
183 }
184 
185 
186 void Timer::pause() {
187 #ifdef FLOW123D_HAVE_PETSC
188  if (Profiler::get_petsc_memory_monitoring()) {
189  // get the maximum resident set size (memory used) for the program.
190  PetscMemoryGetMaximumUsage(&petsc_local_peak_memory);
191  if (petsc_peak_memory < petsc_local_peak_memory)
192  petsc_peak_memory = petsc_local_peak_memory;
193  }
194 #endif // FLOW123D_HAVE_PETSC
195 }
196 
197 void Timer::resume() {
198 #ifdef FLOW123D_HAVE_PETSC
199  if (Profiler::get_petsc_memory_monitoring()) {
200  // tell PETSc to monitor the maximum memory usage so
201  // that PetscMemoryGetMaximumUsage() will work.
202  PetscMemorySetGetMaximumUsage();
203  }
204 #endif // FLOW123D_HAVE_PETSC
205 }
206 
207 void Timer::start() {
208 #ifdef FLOW123D_HAVE_PETSC
209  if (Profiler::get_petsc_memory_monitoring()) {
210  // Tell PETSc to monitor the maximum memory usage so
211  // that PetscMemoryGetMaximumUsage() will work.
212  PetscMemorySetGetMaximumUsage();
213  PetscMemoryGetCurrentUsage (&petsc_start_memory);
214  }
215 #endif // FLOW123D_HAVE_PETSC
216 
217  if (start_count == 0) {
218  start_time = TimePoint();
219  }
220  call_count++;
221  start_count++;
222 }
223 
224 
225 
226 bool Timer::stop(bool forced) {
227 #ifdef FLOW123D_HAVE_PETSC
228  if (Profiler::get_petsc_memory_monitoring()) {
229  // get current memory usage
230  PetscMemoryGetCurrentUsage (&petsc_end_memory);
231  petsc_memory_difference += petsc_end_memory - petsc_start_memory;
232 
233  // get the maximum resident set size (memory used) for the program.
234  PetscMemoryGetMaximumUsage(&petsc_local_peak_memory);
235  if (petsc_peak_memory < petsc_local_peak_memory)
236  petsc_peak_memory = petsc_local_peak_memory;
237  }
238 #endif // FLOW123D_HAVE_PETSC
239 
240  if (forced) start_count=1;
241 
242  if (start_count == 1) {
243  cumul_time += (TimePoint() - start_time);
244  start_count--;
245  return true;
246  } else {
247  start_count--;
248  }
249  return false;
250 }
251 
252 
253 
254 void Timer::add_child(int child_index, const Timer &child)
255 {
256  unsigned int idx = child.hash_idx_;
257  if (child_timers[idx] != timer_no_child) {
258  // hash collision, find first empty place
259  unsigned int i=idx;
260  do {
261  i=( i < max_n_childs ? i+1 : 0);
262  } while (i!=idx && child_timers[i] != timer_no_child);
263  ASSERT(i!=idx)(tag()).error("Too many children of the timer");
264  idx=i;
265  }
266  //DebugOut().fmt("Adding child {} at index: {}\n", child_index, idx);
267  child_timers[idx] = child_index;
268 }
269 
270 
271 
272 string Timer::code_point_str() const {
273  return fmt::format("{%s}:{%d}, {%s}()",
274  code_point_->file_, code_point_->line_, code_point_->func_ );
275 }
276 
277 
278 /***********************************************************************************************
279  * Implementation of Profiler
280  */
281 
282 
283 Profiler * Profiler::instance(bool clear) {
284  static Profiler * _instance = NULL;
285 
286  if (clear) {
287  if (_instance != NULL) {
288  delete _instance;
289  _instance = NULL;
290  }
291  return _instance;
292  }
293 
294  if (_instance == NULL) {
295  MemoryAlloc::malloc_map().reserve(Profiler::malloc_map_reserve);
296  _instance = new Profiler();
297  }
298 
299  return _instance;
300  }
301 
302 
303 // static CONSTEXPR_ CodePoint main_cp = CODE_POINT("Whole Program");
304 const long Profiler::malloc_map_reserve = 100 * 1000;
305 
307 : actual_node(0),
308  task_size_(1),
309  start_time( time(NULL) ),
310  json_filepath(""),
311  none_timer_(CODE_POINT("NONE TIMER"), 0),
312  calibration_time_(-1)
313 
314 {
315  static CONSTEXPR_ CodePoint main_cp = CODE_POINT("Whole Program");
316  set_memory_monitoring(true, true);
317 #ifdef FLOW123D_DEBUG_PROFILER
318  MemoryAlloc::malloc_map().reserve(Profiler::malloc_map_reserve);
319  timers_.push_back( Timer(main_cp, 0) );
320  timers_[0].start();
321 #endif
322 }
323 
324 
325 void Profiler::calibrate() {
326 
327  uint SIZE = 64 * 1024;
328  uint HALF = SIZE / 2;
329  START_TIMER("UNIT_PAYLOAD");
330  Timer &timer = timers_[actual_node];
331 
332  uint count = 0;
333  while(TimePoint() - timer.start_time < 0.1) {
334  double * block = new double[SIZE];
335  for(uint i=0; i<HALF; i++) {
336  block[HALF+i] = block[i]*block[i] + i;
337  }
338  delete block;
339  count++;
340  }
341 
342 
343 
344  END_TIMER("UNIT_PAYLOAD");
345 
346  const double reference_count = 7730;
347  // reference count is average of 5 measurements on payload free laptop:
348  // DELL Latitude E5570, Intel(R) Core(TM) i7-6820HQ CPU @ 2.70GHz
349  // new/delete makes about 3% of time
350  calibration_time_ = 10 * timer.cumul_time * reference_count / count;
351  LogOut() << "Profiler calibration count: " << count << std::endl;
352  LogOut() << "Profiler calibration time: " << calibration_time_ << std::endl;
353 
354 }
355 
356 
357 
358 void Profiler::propagate_timers() {
359  for (unsigned int i = 0; i < Timer::max_n_childs; i++) {
360  unsigned int child_timer = timers_[0].child_timers[i];
361  if ((signed int)child_timer != timer_no_child) {
362  // propagate metrics from child to Whole-Program time-frame
363  accept_from_child(timers_[0], timers_[child_timer]);
364  }
365  }
366 }
367 
368 
369 
370 void Profiler::set_task_info(string description, int size) {
371  task_description_ = description;
372  task_size_ = size;
373 }
374 
375 
376 
377 void Profiler::set_program_info(string program_name, string program_version, string branch, string revision, string build) {
378  flow_name_ = program_name;
379  flow_version_ = program_version;
380  flow_branch_ = branch;
381  flow_revision_ = revision;
382  flow_build_ = build;
383 }
384 
385 
386 
387 int Profiler::start_timer(const CodePoint &cp) {
388  unsigned int parent_node = actual_node;
389  //DebugOut().fmt("Start timer: {}\n", cp.tag_);
390  int child_idx = find_child(cp);
391  if (child_idx < 0) {
392  //DebugOut().fmt("Adding timer: {}\n", cp.tag_);
393  // tag not present - create new timer
394  child_idx=timers_.size();
395  timers_.push_back( Timer(cp, actual_node) );
396  timers_[actual_node].add_child(child_idx , timers_.back() );
397  }
398  actual_node=child_idx;
399 
400  // pause current timer
401  timers_[parent_node].pause();
402 
403  timers_[actual_node].start();
404 
405  return actual_node;
406 }
407 
408 
409 
410 int Profiler::find_child(const CodePoint &cp) {
411  Timer &timer =timers_[actual_node];
412  unsigned int idx = cp.hash_idx_;
413  unsigned int child_idx;
414  do {
415  if (timer.child_timers[idx] == timer_no_child) break; // tag is not there
416 
417  child_idx=timer.child_timers[idx];
418  ASSERT_LT(child_idx, timers_.size()).error();
419  if (timers_[child_idx].full_hash_ == cp.hash_) return child_idx;
420  idx = ( (unsigned int)(idx)==(Timer::max_n_childs - 1) ? 0 : idx+1 );
421  } while ( (unsigned int)(idx) != cp.hash_idx_ ); // passed through whole array
422  return -1;
423 }
424 
425 
426 
427 void Profiler::stop_timer(const CodePoint &cp) {
428 #ifdef FLOW123D_DEBUG
429  // check that all childrens are closed
430  Timer &timer=timers_[actual_node];
431  for(unsigned int i=0; i < Timer::max_n_childs; i++)
432  if (timer.child_timers[i] != timer_no_child)
433  ASSERT(! timers_[timer.child_timers[i]].running())(timers_[timer.child_timers[i]].tag())(timer.tag())
434  .error("Child timer running while closing timer.");
435 #endif
436  unsigned int child_timer = actual_node;
437  if ( cp.hash_ != timers_[actual_node].full_hash_) {
438  // timer to close is not actual - we search for it above actual
439  for(unsigned int node=actual_node; node != 0; node=timers_[node].parent_timer) {
440  if ( cp.hash_ == timers_[node].full_hash_) {
441  // found above - close all nodes between
442  for(; (unsigned int)(actual_node) != node; actual_node=timers_[actual_node].parent_timer) {
443  WarningOut() << "Timer to close '" << cp.tag_ << "' do not match actual timer '"
444  << timers_[actual_node].tag() << "'. Force closing actual." << std::endl;
445  timers_[actual_node].stop(true);
446  }
447  // close 'node' itself
448  timers_[actual_node].stop(false);
449  actual_node = timers_[actual_node].parent_timer;
450 
451  // actual_node == child_timer indicates this is root
452  if (actual_node == child_timer)
453  return;
454 
455  // resume current timer
456  timers_[actual_node].resume();
457  return;
458  }
459  }
460  // node not found - do nothing
461  return;
462  }
463  // node to close match the actual
464  timers_[actual_node].stop(false);
465  actual_node = timers_[actual_node].parent_timer;
466 
467  // actual_node == child_timer indicates this is root
468  if (actual_node == child_timer)
469  return;
470 
471  // resume current timer
472  timers_[actual_node].resume();
473 }
474 
475 
476 
477 void Profiler::stop_timer(int timer_index) {
478  // stop_timer with CodePoint type
479  // timer which is still running MUST be the same as actual_node index
480  // if timer is not running index will differ
481  if (timers_[timer_index].running()) {
482  ASSERT_EQ(timer_index, (int)actual_node).error();
483  stop_timer(*timers_[timer_index].code_point_);
484  }
485 
486 }
487 
488 
489 
490 void Profiler::add_calls(unsigned int n_calls) {
491  timers_[actual_node].call_count += n_calls-1;
492 }
493 
494 
495 
496 void Profiler::notify_malloc(const size_t size, const long p) {
497  if (!global_monitor_memory)
498  return;
499 
500  MemoryAlloc::malloc_map()[p] = static_cast<int>(size);
501  timers_[actual_node].total_allocated_ += size;
502  timers_[actual_node].current_allocated_ += size;
503  timers_[actual_node].alloc_called++;
504 
505  if (timers_[actual_node].current_allocated_ > timers_[actual_node].max_allocated_)
506  timers_[actual_node].max_allocated_ = timers_[actual_node].current_allocated_;
507 }
508 
509 
510 
511 void Profiler::notify_free(const long p) {
512  if (!global_monitor_memory)
513  return;
514 
515  int size = sizeof(p);
516  if (MemoryAlloc::malloc_map()[(long)p] > 0) {
517  size = MemoryAlloc::malloc_map()[(long)p];
518  MemoryAlloc::malloc_map().erase((long)p);
519  }
520  timers_[actual_node].total_deallocated_ += size;
521  timers_[actual_node].current_allocated_ -= size;
522  timers_[actual_node].dealloc_called++;
523 }
524 
525 
526 double Profiler::get_resolution () {
527  const int measurements = 100;
528  double result = 0;
529 
530  // perform 100 measurements
531  for (unsigned int i = 1; i < measurements; i++) {
532  TimePoint t1 = TimePoint ();
533  TimePoint t2 = TimePoint ();
534 
535  // double comparison should be avoided
536  while ((t2 - t1) == 0) t2 = TimePoint ();
537  // while ((t2.ticks - t1.ticks) == 0) t2 = TimePoint ();
538 
539  result += t2 - t1;
540  }
541 
542  return (result / measurements) * 1000; // ticks to seconds to microseconds conversion
543 }
544 
545 
546 Timer Profiler::find_timer(string tag) {
547  for(auto t : timers_) {
548  if (t.tag() == tag) return t;
549  }
550  return none_timer_;
551 }
552 
553 
554 std::string common_prefix( std::string a, std::string b ) {
555  if( a.size() > b.size() ) std::swap(a,b) ;
556  return std::string( a.begin(), std::mismatch( a.begin(), a.end(), b.begin() ).first ) ;
557 }
558 
559 
560 
561 template<typename ReduceFunctor>
562 void Profiler::add_timer_info(ReduceFunctor reduce, nlohmann::json* holder, int timer_idx, double parent_time) {
563 
564  // get timer and check preconditions
565  Timer &timer = timers_[timer_idx];
566  ASSERT(timer_idx >=0)(timer_idx).error("Wrong timer index.");
567  ASSERT(timer.parent_timer >=0).error("Inconsistent tree.");
568 
569  // fix path
570  string filepath = timer.code_point_->file_;
571 
572  // if constant FLOW123D_SOURCE_DIR is defined, we try to erase it from beginning of each CodePoint's filepath
573  #ifdef FLOW123D_SOURCE_DIR
574  string common_path = common_prefix (string(FLOW123D_SOURCE_DIR), filepath);
575  filepath.erase (0, common_path.size());
576  #endif
577 
578 
579  // generate node representing this timer
580  // add basic information
581  nlohmann::json node;
582  double cumul_time_sum;
583  node["tag"] = timer.tag();
584  node["file-path"] = filepath;
585  node["file-line"] = timer.code_point_->line_;
586  node["function"] = timer.code_point_->func_;
587  cumul_time_sum = reduce(timer, node);
588 
589 
590  // statistical info
591  if (timer_idx == 0) parent_time = cumul_time_sum;
592  double percent = parent_time > 1.0e-10 ? cumul_time_sum / parent_time * 100.0 : 0.0;
593  node["percent"] = percent;
594 
595  // write times children timers using secured child_timers array
596  nlohmann::json children;
597  bool has_children = false;
598  for (unsigned int i = 0; i < Timer::max_n_childs; i++) {
599  if (timer.child_timers[i] != timer_no_child) {
600  add_timer_info (reduce, &children, timer.child_timers[i], cumul_time_sum);
601  /*
602  if (child_timers[i] != timer_no_child) {
603  add_timer_info (reduce, &children, child_timers[i], cumul_time_sum); */
604  has_children = true;
605  }
606  }
607 
608  // add children tag and other info if present
609  if (has_children) {
610  node["children"] = children;
611  }
612 
613  // add to the array
614  holder->push_back(node);
615 }
616 
617 
618 template <class T>
619 void save_nonmpi_metric (nlohmann::json &node, T * ptr, string name) {
620  node[name + "-min"] = *ptr;
621  node[name + "-max"] = *ptr;
622  node[name + "-sum"] = *ptr;
623 }
624 
625 std::shared_ptr<std::ostream> Profiler::get_default_output_stream() {
626  json_filepath = FilePath("profiler_info.log.json", FilePath::output_file);
627 
628  //LogOut() << "output into: " << json_filepath << std::endl;
629  return make_shared<ofstream>(json_filepath.c_str());
630 }
631 
632 
633 #ifdef FLOW123D_HAVE_MPI
634 template <class T>
635 void save_mpi_metric (nlohmann::json &node, MPI_Comm comm, T * ptr, string name) {
636  node[name + "-min"] = MPI_Functions::min(ptr, comm);
637  node[name + "-max"] = MPI_Functions::max(ptr, comm);
638  node[name + "-sum"] = MPI_Functions::sum(ptr, comm);
639 }
640 
641 void Profiler::output(MPI_Comm comm, ostream &os) {
642  int mpi_rank, mpi_size;
643  //wait until profiling on all processors is finished
644  MPI_Barrier(comm);
645  stop_timer(0);
646  propagate_timers();
647 
648  // stop monitoring memory
649  bool temp_memory_monitoring = global_monitor_memory;
650  set_memory_monitoring(false, petsc_monitor_memory);
651 
652  chkerr( MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank) );
653  MPI_Comm_size(comm, &mpi_size);
654 
655  // output header
656  nlohmann::json jsonRoot, jsonChildren;
657 
658  // recursively add all timers info
659  // define lambda function which reduces timer from multiple processors
660  // MPI implementation uses MPI call to reduce values
661  auto reduce = [=] (Timer &timer, nlohmann::json &node) -> double {
662  int call_count = timer.call_count;
663  double cumul_time = timer.cumulative_time ();
664 
665  long memory_allocated = (long)timer.total_allocated_;
666  long memory_deallocated = (long)timer.total_deallocated_;
667  long memory_peak = (long)timer.max_allocated_;
668 
669  int alloc_called = timer.alloc_called;
670  int dealloc_called = timer.dealloc_called;
671 
672 
673  save_mpi_metric<double>(node, comm, &cumul_time, "cumul-time");
674  save_mpi_metric<int>(node, comm, &call_count, "call-count");
675 
676  save_mpi_metric<long>(node, comm, &memory_allocated, "memory-alloc");
677  save_mpi_metric<long>(node, comm, &memory_deallocated, "memory-dealloc");
678  save_mpi_metric<long>(node, comm, &memory_peak, "memory-peak");
679  //
680  save_mpi_metric<int>(node, comm, &alloc_called, "memory-alloc-called");
681  save_mpi_metric<int>(node, comm, &dealloc_called, "memory-dealloc-called");
682 
683 #ifdef FLOW123D_HAVE_PETSC
684  long petsc_memory_difference = (long)timer.petsc_memory_difference;
685  long petsc_peak_memory = (long)timer.petsc_peak_memory;
686  save_mpi_metric<long>(node, comm, &petsc_memory_difference, "memory-petsc-diff");
687  save_mpi_metric<long>(node, comm, &petsc_peak_memory, "memory-petsc-peak");
688 #endif // FLOW123D_HAVE_PETSC
689 
690  return MPI_Functions::sum(&cumul_time, comm);
691  };
692 
693  add_timer_info (reduce, &jsonChildren, 0, 0.0);
694  jsonRoot["children"] = jsonChildren;
695  output_header(jsonRoot, mpi_size);
696 
697 
698  // create profiler output only once (on the first processor)
699  // only active communicator should be the one with mpi_rank 0
700  if (mpi_rank == 0) {
701  try {
702  /**
703  * indent size
704  * results in json human readable format (indents, newlines)
705  */
706  const int FLOW123D_JSON_HUMAN_READABLE = 2;
707  // write result to stream
708  os << jsonRoot.dump(FLOW123D_JSON_HUMAN_READABLE) << endl;
709 
710  } catch (exception & e) {
711  stringstream ss;
712  ss << "nlohmann::json::dump error: " << e.what() << "\n";
713  THROW( ExcMessage() << EI_Message(ss.str()) );
714  }
715  }
716  // restore memory monitoring
717  set_memory_monitoring(temp_memory_monitoring, petsc_monitor_memory);
718 }
719 
720 
721 void Profiler::output(MPI_Comm comm, string profiler_path /* = "" */) {
722  int mpi_rank;
723  chkerr(MPI_Comm_rank(MPI_COMM_WORLD, &mpi_rank));
724 
725  if (mpi_rank == 0) {
726  if (profiler_path == "") {
727  output(comm, *get_default_output_stream());
728  } else {
729  json_filepath = profiler_path;
730  std::shared_ptr<std::ostream> os = make_shared<ofstream>(profiler_path.c_str());
731  output(comm, *os);
732  }
733  } else {
734  ostringstream os;
735  output(comm, os);
736  }
737 }
738 
739 #endif /* FLOW123D_HAVE_MPI */
740 
741 void Profiler::output(ostream &os) {
742  // last update
743  stop_timer(0);
744  propagate_timers();
745 
746  // output header
747  nlohmann::json jsonRoot, jsonChildren;
748  /**
749  * Constant representing number of MPI processes
750  * where there is no MPI to work with (so 1 process)
751  */
752  const int FLOW123D_MPI_SINGLE_PROCESS = 1;
753  output_header(jsonRoot, FLOW123D_MPI_SINGLE_PROCESS);
754 
755 
756  // recursively add all timers info
757  // define lambda function which reduces timer from multiple processors
758  // non-MPI implementation is just dummy repetition of initial value
759  auto reduce = [=] (Timer &timer, nlohmann::json &node) -> double {
760  int call_count = timer.call_count;
761  double cumul_time = timer.cumulative_time ();
762 
763  long memory_allocated = (long)timer.total_allocated_;
764  long memory_deallocated = (long)timer.total_deallocated_;
765  long memory_peak = (long)timer.max_allocated_;
766 
767  int alloc_called = timer.alloc_called;
768  int dealloc_called = timer.dealloc_called;
769 
770  save_nonmpi_metric<double>(node, &cumul_time, "cumul-time");
771  save_nonmpi_metric<int>(node, &call_count, "call-count");
772 
773  save_nonmpi_metric<long>(node, &memory_allocated, "memory-alloc");
774  save_nonmpi_metric<long>(node, &memory_deallocated, "memory-dealloc");
775  save_nonmpi_metric<long>(node, &memory_peak, "memory-peak");
776 
777  save_nonmpi_metric<int>(node, &alloc_called, "memory-alloc-called");
778  save_nonmpi_metric<int>(node, &dealloc_called, "memory-dealloc-called");
779 
780 #ifdef FLOW123D_HAVE_PETSC
781  long petsc_memory_difference = (long)timer.petsc_memory_difference;
782  long petsc_peak_memory = (long)timer.petsc_peak_memory;
783  save_nonmpi_metric<long>(node, &petsc_memory_difference, "memory-petsc-diff");
784  save_nonmpi_metric<long>(node, &petsc_peak_memory, "memory-petsc-peak");
785 #endif // FLOW123D_HAVE_PETSC
786 
787  return cumul_time;
788  };
789 
790  add_timer_info(reduce, &jsonChildren, 0, 0.0);
791  jsonRoot["children"] = jsonChildren;
792 
793  try {
794  /**
795  * indent size
796  * results in json human readable format (indents, newlines)
797  */
798  const int FLOW123D_JSON_HUMAN_READABLE = 2;
799  // write result to stream
800  os << jsonRoot.dump(FLOW123D_JSON_HUMAN_READABLE) << endl;
801 
802  } catch (exception & e) {
803  stringstream ss;
804  ss << "nlohmann::json::dump error: " << e.what() << "\n";
805  THROW( ExcMessage() << EI_Message(ss.str()) );
806  }
807 }
808 
809 
810 void Profiler::output(string profiler_path /* = "" */) {
811  if(profiler_path == "") {
812  output(*get_default_output_stream());
813  } else {
814  json_filepath = profiler_path;
815  std::shared_ptr<std::ostream> os = make_shared<ofstream>(profiler_path.c_str());
816  output(*os);
817  }
818 }
819 
820 void Profiler::output_header (nlohmann::json &root, int mpi_size) {
821  time_t end_time = time(NULL);
822 
823  const char format[] = "%x %X";
824  char start_time_string[BUFSIZ] = {0};
825  strftime(start_time_string, sizeof (start_time_string) - 1, format, localtime(&start_time));
826 
827  char end_time_string[BUFSIZ] = {0};
828  strftime(end_time_string, sizeof (end_time_string) - 1, format, localtime(&end_time));
829 
830  if (timers_[0].cumul_time > 60) {
831  calibrate();
832  }
833  // generate current run details
834  root["program-name"] = flow_name_;
835  root["program-version"] = flow_version_;
836  root["program-branch"] = flow_branch_;
837  root["program-revision"] = flow_revision_;
838  root["program-build"] = flow_build_;
839  root["timer-resolution"] = Profiler::get_resolution();
840  root["timer-calibration"] = calibration_time_;
841 
842  // print some information about the task at the beginning
843  root["task-description"] = task_description_;
844  root["task-size"] = task_size_;
845 
846  //print some information about the task at the beginning
847  root["run-process-count"] = mpi_size;
848  root["run-started-at"] = start_time_string;
849  root["run-finished-at"] = end_time_string;
850 }
851 
852 #ifdef FLOW123D_HAVE_PYTHON
853 void Profiler::transform_profiler_data (const string &output_file_suffix, const string &formatter) {
854 
855  if (json_filepath == "") return;
856 
857  // error under CYGWIN environment : more details in this repo
858  // https://github.com/x3mSpeedy/cygwin-python-issue
859  //
860  // For now we only support profiler report conversion in UNIX environment
861  // Windows users will have to use a python script located in bin folder
862  //
863 
864  #ifndef FLOW123D_HAVE_CYGWIN
865  // grab module and function by importing module profiler_formatter_module.py
866  PyObject * python_module = PythonLoader::load_module_by_name ("profiler.profiler_formatter_module");
867  //
868  // def convert (json_location, output_file, formatter):
869  //
870  PyObject * convert_method = PythonLoader::get_callable (python_module, "convert" );
871 
872  int argument_index = 0;
873  PyObject * arguments = PyTuple_New (3);
874 
875  // set json path location as first argument
876  PyObject * tmp = PyUnicode_FromString (json_filepath.c_str());
877  PyTuple_SetItem (arguments, argument_index++, tmp);
878 
879  // set output path location as second argument
880  tmp = PyUnicode_FromString ((json_filepath + output_file_suffix).c_str());
881  PyTuple_SetItem (arguments, argument_index++, tmp);
882 
883  // set Formatter class as third value
884  tmp = PyUnicode_FromString (formatter.c_str());
885  PyTuple_SetItem (arguments, argument_index++, tmp);
886 
887  // execute method with arguments
888  PyObject_CallObject (convert_method, arguments);
889 
890  PythonLoader::check_error();
891 
892  #else
893 
894  // print information about windows-cygwin issue and offer manual solution
895  MessageOut() << "# Note: converting json profiler reports is not"
896  << " supported under Windows or Cygwin environment for now.\n"
897  << "# You can use python script located in bin/python folder"
898  << " in order to convert json report to txt or csv format.\n"
899  << "python profiler_formatter_script.py --input \"" << json_filepath
900  << "\" --output \"profiler.txt\"" << std::endl;
901  #endif // FLOW123D_HAVE_CYGWIN
902 }
903 #else
904 void Profiler::transform_profiler_data (const string &, const string &) {
905 }
906 
907 #endif // FLOW123D_HAVE_PYTHON
908 
909 
910 void Profiler::uninitialize() {
911  if (Profiler::instance()) {
912  ASSERT(Profiler::instance()->actual_node==0)(Profiler::instance()->timers_[Profiler::instance()->actual_node].tag())
913  .error("Forbidden to uninitialize the Profiler when actual timer is not zero.");
914  Profiler::instance()->stop_timer(0);
915  set_memory_monitoring(false, false);
916  Profiler::instance(true);
917  }
918 }
919 bool Profiler::global_monitor_memory = false;
920 bool Profiler::petsc_monitor_memory = true;
921 void Profiler::set_memory_monitoring(const bool global_monitor, const bool petsc_monitor) {
922  global_monitor_memory = global_monitor;
923  petsc_monitor_memory = petsc_monitor;
924 }
925 
926 bool Profiler::get_global_memory_monitoring() {
927  return global_monitor_memory;
928 }
929 
930 bool Profiler::get_petsc_memory_monitoring() {
931  return petsc_monitor_memory;
932 }
933 
934 unordered_map_with_alloc & MemoryAlloc::malloc_map() {
935  static unordered_map_with_alloc static_malloc_map;
936  return static_malloc_map;
937 }
938 
939 void * Profiler::operator new (size_t size) {
940  return malloc (size);
941 }
942 
943 void Profiler::operator delete (void* p) {
944  free(p);
945 }
946 
947 void *operator new (std::size_t size) OPERATOR_NEW_THROW_EXCEPTION {
948  void * p = malloc(size);
949  Profiler::instance()->notify_malloc(size, (long)p);
950  return p;
951 }
952 
953 void *operator new[] (std::size_t size) OPERATOR_NEW_THROW_EXCEPTION {
954  void * p = malloc(size);
955  Profiler::instance()->notify_malloc(size, (long)p);
956  return p;
957 }
958 
959 void *operator new[] (std::size_t size, const std::nothrow_t&) throw() {
960  void * p = malloc(size);
961  Profiler::instance()->notify_malloc(size, (long)p);
962  return p;
963 }
964 
965 void operator delete( void *p) throw() {
966  Profiler::instance()->notify_free((long)p);
967  free(p);
968 }
969 
970 void operator delete( void *p, std::size_t) throw() {
971  Profiler::instance()->notify_free((long)p);
972  free(p);
973 }
974 
975 void operator delete[]( void *p) throw() {
976  Profiler::instance()->notify_free((long)p);
977  free(p);
978 }
979 
980 void operator delete[]( void *p, std::size_t) throw() {
981  Profiler::instance()->notify_free((long)p);
982  free(p);
983 }
984 
985 #else // def FLOW123D_DEBUG_PROFILER
986 
987 Profiler * Profiler::instance(bool clear) {
988  static Profiler * _instance = new Profiler();
989  return _instance;
990 }
991 
992 
993 // void Profiler::initialize() {
994 // if (_instance == NULL) {
995 // _instance = new Profiler();
996 // set_memory_monitoring(true, true);
997 // }
998 // }
999 
1001 
1002 
1003 #endif // def FLOW123D_DEBUG_PROFILER
#define MPI_LONG
Definition: mpi.h:161
void set_task_info(string, int)
static int min(int *val, MPI_Comm comm)
Definition: sys_profiler.cc:62
#define CONSTEXPR_
Definition: sys_profiler.hh:86
unsigned int uint
int MPI_Comm
Definition: mpi.h:141
#define OPERATOR_NEW_THROW_EXCEPTION
Definition: system.hh:54
static int sum(int *val, MPI_Comm comm)
Definition: sys_profiler.cc:44
void output()
a class to store JSON values
Definition: json.hpp:173
#define MessageOut()
Macro defining &#39;message&#39; record of log.
Definition: logger.hh:267
double get_resolution() const
std::string format(CStringRef format_str, ArgList args)
Definition: format.h:3141
#define MPI_SUM
Definition: mpi.h:196
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
Definition: system.hh:148
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:347
void set_program_info(string, string, string, string, string)
#define LogOut()
Macro defining &#39;log&#39; record of log.
Definition: logger.hh:273
#define MPI_Reduce(sendbuf, recvbuf, count, datatype, op, root, comm)
Definition: mpi.h:608
void calibrate()
static void uninitialize()
void notify_malloc(const size_t)
#define MPI_MIN
Definition: mpi.h:198
#define START_TIMER(tag)
Starts a timer with specified tag.
static int max(int *val, MPI_Comm comm)
Definition: sys_profiler.cc:80
void start()
Definition: memory.cc:39
void swap(nlohmann::json &j1, nlohmann::json &j2) noexcept(is_nothrow_move_constructible< nlohmann::json >::value andis_nothrow_move_assignable< nlohmann::json >::value)
exchanges the values of two JSON objects
Definition: json.hpp:8688
#define MPI_Comm_size
Definition: mpi.h:235
#define MPI_DOUBLE
Definition: mpi.h:156
STREAM & operator<<(STREAM &s, UpdateFlags u)
void push_back(basic_json &&val)
add an object to an array
Definition: json.hpp:4686
#define MPI_Comm_rank
Definition: mpi.h:236
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
void notify_free(const size_t)
#define MPI_INT
Definition: mpi.h:160
#define ASSERT_LT(a, b)
Definition of comparative assert macro (Less Than)
Definition: asserts.hh:296
string_t dump(const int indent=-1) const
serialization
Definition: json.hpp:2079
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:270
#define MPI_COMM_WORLD
Definition: mpi.h:123
#define END_TIMER(tag)
Ends a timer with specified tag.
#define MPI_MAX
Definition: mpi.h:197
static Profiler * instance(bool clear=false)
Definition: memory.cc:33
#define MPI_Barrier(comm)
Definition: mpi.h:531
void stop()
Definition: memory.cc:42
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
Definition: exceptions.hh:53
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual)
Definition: asserts.hh:328
void transform_profiler_data(const string &, const string &)