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