Reference documentation for deal.II version 8.5.1
utilities.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2005 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #include <deal.II/base/utilities.h>
18 #include <deal.II/base/mpi.h>
19 #include <deal.II/base/exceptions.h>
20 #include <deal.II/base/thread_local_storage.h>
21 
22 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
23 #include <boost/math/special_functions/erf.hpp>
24 #include <boost/lexical_cast.hpp>
25 #include <boost/random.hpp>
26 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
27 
28 #include <algorithm>
29 #include <cctype>
30 #include <cerrno>
31 #include <cmath>
32 #include <cstddef>
33 #include <cstdio>
34 #include <ctime>
35 #include <fstream>
36 #include <iomanip>
37 #include <iostream>
38 #include <limits>
39 #include <sstream>
40 
41 #ifndef DEAL_II_MSVC
42 # include <stdlib.h>
43 #endif
44 
45 #ifdef DEAL_II_MSVC
46 # include <winsock2.h>
47 #endif
48 
49 
50 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
51 #ifdef DEAL_II_WITH_TRILINOS
52 # ifdef DEAL_II_WITH_MPI
53 # include <Epetra_MpiComm.h>
54 # include <deal.II/lac/vector_memory.h>
55 # include <deal.II/lac/trilinos_vector.h>
56 # include <deal.II/lac/trilinos_block_vector.h>
57 # endif
58 # include "Teuchos_RCP.hpp"
59 # include "Epetra_SerialComm.h"
60 #endif
61 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
62 
63 
64 
65 DEAL_II_NAMESPACE_OPEN
66 
67 
68 namespace Utilities
69 {
70 
71 
73  unsigned int, unsigned int,
74  << "When trying to convert " << arg1
75  << " to a string with " << arg2 << " digits");
77  unsigned int,
78  << "Invalid number " << arg1);
80  std::string,
81  << "Can't convert the string " << arg1
82  << " to the desired type");
83 
84  std::string
85  int_to_string (const unsigned int value, const unsigned int digits)
86  {
87  return to_string(value,digits);
88  }
89 
90  template <typename number>
91  std::string
92  to_string (const number value, const unsigned int digits)
93  {
94  std::string lc_string = boost::lexical_cast<std::string>(value);
95 
96  if (digits == numbers::invalid_unsigned_int)
97  return lc_string;
98  else if (lc_string.size() < digits)
99  {
100  // We have to add the padding zeroes in front of the number
101  const unsigned int padding_position = (lc_string[0] == '-')
102  ?
103  1
104  :
105  0;
106 
107  const std::string padding(digits - lc_string.size(), '0');
108  lc_string.insert(padding_position, padding);
109  }
110  return lc_string;
111  }
112 
113 
114  std::string
115  replace_in_string(const std::string &input,
116  const std::string &from,
117  const std::string &to)
118  {
119  if (from.empty())
120  return input;
121 
122  std::string out = input;
123  std::string::size_type pos = out.find(from);
124 
125  while (pos != std::string::npos)
126  {
127  out.replace(pos, from.size(), to);
128  pos = out.find(from, pos + to.size());
129  }
130  return out;
131  }
132 
133  std::string
134  trim(const std::string &input)
135  {
136  std::string::size_type left = 0;
137  std::string::size_type right = input.size() > 0 ? input.size() - 1 : 0;
138 
139  for (; left < input.size(); ++left)
140  {
141  if (!std::isspace(input[left]))
142  {
143  break;
144  }
145  }
146 
147  for (; right >= left; --right)
148  {
149  if (!std::isspace(input[right]))
150  {
151  break;
152  }
153  }
154 
155  return std::string(input, left, right - left + 1);
156  }
157 
158 
159 
160  std::string
161  dim_string(const int dim, const int spacedim)
162  {
163  if (dim == spacedim)
164  return int_to_string(dim);
165  else
166  return int_to_string(dim)+","+int_to_string(spacedim);
167  }
168 
169 
170  unsigned int
171  needed_digits (const unsigned int max_number)
172  {
173  if (max_number < 10)
174  return 1;
175  if (max_number < 100)
176  return 2;
177  if (max_number < 1000)
178  return 3;
179  if (max_number < 10000)
180  return 4;
181  if (max_number < 100000)
182  return 5;
183  if (max_number < 1000000)
184  return 6;
185  AssertThrow (false, ExcInvalidNumber(max_number));
186  return 0;
187  }
188 
189 
190 
191  int
192  string_to_int (const std::string &s_)
193  {
194  // trim whitespace on either side of the text if necessary
195  std::string s = s_;
196  while ((s.size() > 0) && (s[0] == ' '))
197  s.erase (s.begin());
198  while ((s.size() > 0) && (s[s.size()-1] == ' '))
199  s.erase (s.end()-1);
200 
201  // now convert and see whether we succeed. note that strtol only
202  // touches errno if an error occurred, so if we want to check
203  // whether an error happened, we need to make sure that errno==0
204  // before calling strtol since otherwise it may be that the
205  // conversion succeeds and that errno remains at the value it
206  // was before, whatever that was
207  char *p;
208  errno = 0;
209  const int i = std::strtol(s.c_str(), &p, 10);
210  AssertThrow ( !((errno != 0) || (s.size() == 0) || ((s.size()>0) && (*p != '\0'))),
211  ExcMessage ("Can't convert <" + s + "> to an integer."));
212 
213  return i;
214  }
215 
216 
217 
218  std::vector<int>
219  string_to_int (const std::vector<std::string> &s)
220  {
221  std::vector<int> tmp (s.size());
222  for (unsigned int i=0; i<s.size(); ++i)
223  tmp[i] = string_to_int (s[i]);
224  return tmp;
225  }
226 
227 
228 
229  double
230  string_to_double (const std::string &s_)
231  {
232  // trim whitespace on either side of the text if necessary
233  std::string s = s_;
234  while ((s.size() > 0) && (s[0] == ' '))
235  s.erase (s.begin());
236  while ((s.size() > 0) && (s[s.size()-1] == ' '))
237  s.erase (s.end()-1);
238 
239  // now convert and see whether we succeed. note that strtol only
240  // touches errno if an error occurred, so if we want to check
241  // whether an error happened, we need to make sure that errno==0
242  // before calling strtol since otherwise it may be that the
243  // conversion succeeds and that errno remains at the value it
244  // was before, whatever that was
245  char *p;
246  errno = 0;
247  const double d = std::strtod(s.c_str(), &p);
248  AssertThrow ( !((errno != 0) || (s.size() == 0) || ((s.size()>0) && (*p != '\0'))),
249  ExcMessage ("Can't convert <" + s + "> to a double."));
250 
251  return d;
252  }
253 
254 
255 
256  std::vector<double>
257  string_to_double (const std::vector<std::string> &s)
258  {
259  std::vector<double> tmp (s.size());
260  for (unsigned int i=0; i<s.size(); ++i)
261  tmp[i] = string_to_double (s[i]);
262  return tmp;
263  }
264 
265 
266 
267  std::vector<std::string>
268  split_string_list (const std::string &s,
269  const char delimiter)
270  {
271  // keep the currently remaining part of the input string in 'tmp' and
272  // keep chopping elements of the list off the front
273  std::string tmp = s;
274 
275  // as discussed in the documentation, eat whitespace from the end
276  // of the string
277  while (tmp.length() != 0 && tmp[tmp.length()-1] == ' ')
278  tmp.erase (tmp.length()-1, 1);
279 
280  // split the input list until it is empty. since in every iteration
281  // 'tmp' is what's left of the string after the next delimiter,
282  // and since we've stripped trailing space already, 'tmp' will
283  // be empty at one point if 's' ended in a delimiter, even if
284  // there was space after the last delimiter. this matches what's
285  // discussed in the documentation
286  std::vector<std::string> split_list;
287  split_list.reserve (std::count (tmp.begin(), tmp.end(), delimiter)+1);
288  while (tmp.length() != 0)
289  {
290  std::string name;
291  name = tmp;
292 
293  if (name.find(delimiter) != std::string::npos)
294  {
295  name.erase (name.find(delimiter), std::string::npos);
296  tmp.erase (0, tmp.find(delimiter)+1);
297  }
298  else
299  tmp = "";
300 
301  // strip spaces from this element's front and end
302  while ((name.length() != 0) && (name[0] == ' '))
303  name.erase (0,1);
304  while (name.length() != 0 && name[name.length()-1] == ' ')
305  name.erase (name.length()-1, 1);
306 
307  split_list.push_back (name);
308  }
309 
310  return split_list;
311  }
312 
313 
314 
315  std::vector<std::string>
316  break_text_into_lines (const std::string &original_text,
317  const unsigned int width,
318  const char delimiter)
319  {
320  std::string text = original_text;
321  std::vector<std::string> lines;
322 
323  // remove trailing spaces
324  while ((text.length() != 0) && (text[text.length()-1] == delimiter))
325  text.erase(text.length()-1,1);
326 
327  // then split the text into lines
328  while (text.length() != 0)
329  {
330  // in each iteration, first remove
331  // leading spaces
332  while ((text.length() != 0) && (text[0] == delimiter))
333  text.erase(0, 1);
334 
335  std::size_t pos_newline = text.find_first_of("\n", 0);
336  if (pos_newline != std::string::npos && pos_newline <= width)
337  {
338  std::string line (text, 0, pos_newline);
339  while ((line.length() != 0) && (line[line.length()-1] == delimiter))
340  line.erase(line.length()-1,1);
341  lines.push_back (line);
342  text.erase (0, pos_newline+1);
343  continue;
344  }
345 
346  // if we can fit everything into one
347  // line, then do so. otherwise, we have
348  // to keep breaking
349  if (text.length() < width)
350  {
351  // remove trailing spaces
352  while ((text.length() != 0) && (text[text.length()-1] == delimiter))
353  text.erase(text.length()-1,1);
354  lines.push_back (text);
355  text = "";
356  }
357  else
358  {
359  // starting at position width, find the
360  // location of the previous space, so
361  // that we can break around there
362  int location = std::min<int>(width,text.length()-1);
363  for (; location>0; --location)
364  if (text[location] == delimiter)
365  break;
366 
367  // if there are no spaces, then try if
368  // there are spaces coming up
369  if (location == 0)
370  for (location = std::min<int>(width,text.length()-1);
371  location<static_cast<int>(text.length());
372  ++location)
373  if (text[location] == delimiter)
374  break;
375 
376  // now take the text up to the found
377  // location and put it into a single
378  // line, and remove it from 'text'
379  std::string line (text, 0, location);
380  while ((line.length() != 0) && (line[line.length()-1] == delimiter))
381  line.erase(line.length()-1,1);
382  lines.push_back (line);
383  text.erase (0, location);
384  }
385  }
386 
387  return lines;
388  }
389 
390 
391 
392  bool
393  match_at_string_start (const std::string &name,
394  const std::string &pattern)
395  {
396  if (pattern.size() > name.size())
397  return false;
398 
399  for (unsigned int i=0; i<pattern.size(); ++i)
400  if (pattern[i] != name[i])
401  return false;
402 
403  return true;
404  }
405 
406 
407 
408  std::pair<int, unsigned int>
409  get_integer_at_position (const std::string &name,
410  const unsigned int position)
411  {
412  Assert (position < name.size(), ExcInternalError());
413 
414  const std::string test_string (name.begin()+position,
415  name.end());
416 
417  std::istringstream str(test_string);
418 
419  int i;
420  if (str >> i)
421  {
422  // compute the number of
423  // digits of i. assuming it
424  // is less than 8 is likely
425  // ok
426  if (i<10)
427  return std::make_pair (i, 1U);
428  else if (i<100)
429  return std::make_pair (i, 2U);
430  else if (i<1000)
431  return std::make_pair (i, 3U);
432  else if (i<10000)
433  return std::make_pair (i, 4U);
434  else if (i<100000)
435  return std::make_pair (i, 5U);
436  else if (i<1000000)
437  return std::make_pair (i, 6U);
438  else if (i<10000000)
439  return std::make_pair (i, 7U);
440  else
441  {
442  Assert (false, ExcNotImplemented());
443  return std::make_pair (-1, numbers::invalid_unsigned_int);
444  }
445  }
446  else
447  return std::make_pair (-1, numbers::invalid_unsigned_int);
448  }
449 
450 
451 
452  double
454  const double sigma)
455  {
456  // if no noise: return now
457  if (sigma == 0)
458  return a;
459 
460  // we would want to use rand(), but that function is not reentrant
461  // in a thread context. one could use rand_r, but this does not
462  // produce reproducible results between threads either (though at
463  // least it is reentrant). these two approaches being
464  // non-workable, use a thread-local random number generator here
465  static Threads::ThreadLocalStorage<boost::mt19937> random_number_generator;
466  return boost::normal_distribution<>(a,sigma)(random_number_generator.get());
467  }
468 
469 
470 
471  std::vector<unsigned int>
472  reverse_permutation (const std::vector<unsigned int> &permutation)
473  {
474  const unsigned int n = permutation.size();
475 
476  std::vector<unsigned int> out (n);
477  for (unsigned int i=0; i<n; ++i)
478  out[i] = n - 1 - permutation[i];
479 
480  return out;
481  }
482 
483 
484 
485  std::vector<unsigned int>
486  invert_permutation (const std::vector<unsigned int> &permutation)
487  {
488  const unsigned int n = permutation.size();
489 
490  std::vector<unsigned int> out (n, numbers::invalid_unsigned_int);
491 
492  for (unsigned int i=0; i<n; ++i)
493  {
494  Assert (permutation[i] < n, ExcIndexRange (permutation[i], 0, n));
495  out[permutation[i]] = i;
496  }
497 
498  // check that we have actually reached
499  // all indices
500  for (unsigned int i=0; i<n; ++i)
502  ExcMessage ("The given input permutation had duplicate entries!"));
503 
504  return out;
505  }
506 
507  std::vector<unsigned long long int>
508  reverse_permutation (const std::vector<unsigned long long int> &permutation)
509  {
510  const unsigned long long int n = permutation.size();
511 
512  std::vector<unsigned long long int> out (n);
513  for (unsigned long long int i=0; i<n; ++i)
514  out[i] = n - 1 - permutation[i];
515 
516  return out;
517  }
518 
519 
520 
521  std::vector<unsigned long long int>
522  invert_permutation (const std::vector<unsigned long long int> &permutation)
523  {
524  const unsigned long long int n = permutation.size();
525 
526  std::vector<unsigned long long int> out (n, numbers::invalid_unsigned_int);
527 
528  for (unsigned long long int i=0; i<n; ++i)
529  {
530  Assert (permutation[i] < n, ExcIndexRange (permutation[i], 0, n));
531  out[permutation[i]] = i;
532  }
533 
534  // check that we have actually reached
535  // all indices
536  for (unsigned long long int i=0; i<n; ++i)
538  ExcMessage ("The given input permutation had duplicate entries!"));
539 
540  return out;
541  }
542 
543 
544  template <typename Integer>
545  std::vector<Integer>
546  reverse_permutation (const std::vector<Integer> &permutation)
547  {
548  const unsigned int n = permutation.size();
549 
550  std::vector<Integer> out (n);
551  for (unsigned int i=0; i<n; ++i)
552  out[i] = n - 1 - permutation[i];
553 
554  return out;
555  }
556 
557 
558 
559  template <typename Integer>
560  std::vector<Integer>
561  invert_permutation (const std::vector<Integer> &permutation)
562  {
563  const unsigned int n = permutation.size();
564 
565  std::vector<Integer> out (n, numbers::invalid_unsigned_int);
566 
567  for (unsigned int i=0; i<n; ++i)
568  {
569  Assert (permutation[i] < n, ExcIndexRange (permutation[i], 0, n));
570  out[permutation[i]] = i;
571  }
572 
573  // check that we have actually reached
574  // all indices
575  for (unsigned int i=0; i<n; ++i)
577  ExcMessage ("The given input permutation had duplicate entries!"));
578 
579  return out;
580  }
581 
582 
583 
584 
585  namespace System
586  {
587 #if defined(__linux__)
588 
589  double get_cpu_load ()
590  {
591  std::ifstream cpuinfo;
592  cpuinfo.open("/proc/loadavg");
593 
594  AssertThrow(cpuinfo, ExcIO());
595 
596  double load;
597  cpuinfo >> load;
598 
599  return load;
600  }
601 
602 #else
603 
604  double get_cpu_load ()
605  {
606  return 0.;
607  }
608 
609 #endif
610 
611 
612 
614  {
615  stats.VmPeak = stats.VmSize = stats.VmHWM = stats.VmRSS = 0;
616 
617  // parsing /proc/self/stat would be a
618  // lot easier, but it does not contain
619  // VmHWM, so we use /status instead.
620 #if defined(__linux__)
621  std::ifstream file("/proc/self/status");
622  std::string line;
623  std::string name;
624  while (!file.eof())
625  {
626  file >> name;
627  if (name == "VmPeak:")
628  file >> stats.VmPeak;
629  else if (name == "VmSize:")
630  file >> stats.VmSize;
631  else if (name == "VmHWM:")
632  file >> stats.VmHWM;
633  else if (name == "VmRSS:")
634  {
635  file >> stats.VmRSS;
636  break; //this is always the last entry
637  }
638 
639  getline(file, line);
640  }
641 #endif
642  }
643 
644 
645 
646  std::string get_hostname ()
647  {
648 #if defined(DEAL_II_HAVE_UNISTD_H) && defined(DEAL_II_HAVE_GETHOSTNAME)
649  const unsigned int N=1024;
650  char hostname[N];
651  gethostname (&(hostname[0]), N-1);
652 #else
653  std::string hostname("unknown");
654 #endif
655  return hostname;
656  }
657 
658 
659 
660  std::string get_time ()
661  {
662  std::time_t time1= std::time (0);
663  std::tm *time = std::localtime(&time1);
664 
665  std::ostringstream o;
666  o << time->tm_hour << ":"
667  << (time->tm_min < 10 ? "0" : "") << time->tm_min << ":"
668  << (time->tm_sec < 10 ? "0" : "") << time->tm_sec;
669 
670  return o.str();
671  }
672 
673 
674 
675  std::string get_date ()
676  {
677  std::time_t time1= std::time (0);
678  std::tm *time = std::localtime(&time1);
679 
680  std::ostringstream o;
681  o << time->tm_year + 1900 << "/"
682  << time->tm_mon + 1 << "/"
683  << time->tm_mday;
684 
685  return o.str();
686  }
687 
688 
689 
690  void posix_memalign (void **memptr, size_t alignment, size_t size)
691  {
692 #ifndef DEAL_II_MSVC
693  const int ierr = ::posix_memalign (memptr, alignment, size);
694 
695  AssertThrow (ierr == 0, ExcOutOfMemory());
696  AssertThrow (*memptr != 0, ExcOutOfMemory());
697 #else
698  // Windows does not appear to have posix_memalign. just use the
699  // regular malloc in that case
700  *memptr = malloc (size);
701  (void)alignment;
702  AssertThrow (*memptr != 0, ExcOutOfMemory());
703 #endif
704  }
705 
706 
707 
709  {
711  }
712  }
713 
714 
715 #ifdef DEAL_II_WITH_TRILINOS
716 
717  namespace Trilinos
718  {
719  const Epetra_Comm &
721  {
722 #ifdef DEAL_II_WITH_MPI
723  static Teuchos::RCP<Epetra_MpiComm>
724  communicator = Teuchos::rcp (new Epetra_MpiComm (MPI_COMM_WORLD), true);
725 #else
726  static Teuchos::RCP<Epetra_SerialComm>
727  communicator = Teuchos::rcp (new Epetra_SerialComm (), true);
728 #endif
729 
730  return *communicator;
731  }
732 
733 
734 
735  const Epetra_Comm &
737  {
738 #ifdef DEAL_II_WITH_MPI
739  static Teuchos::RCP<Epetra_MpiComm>
740  communicator = Teuchos::rcp (new Epetra_MpiComm (MPI_COMM_SELF), true);
741 #else
742  static Teuchos::RCP<Epetra_SerialComm>
743  communicator = Teuchos::rcp (new Epetra_SerialComm (), true);
744 #endif
745 
746  return *communicator;
747  }
748 
749 
750 
751  Epetra_Comm *
752  duplicate_communicator (const Epetra_Comm &communicator)
753  {
754 #ifdef DEAL_II_WITH_MPI
755 
756  // see if the communicator is in fact a
757  // parallel MPI communicator; if so,
758  // return a duplicate of it
759  const Epetra_MpiComm
760  *mpi_comm = dynamic_cast<const Epetra_MpiComm *>(&communicator);
761  if (mpi_comm != 0)
762  return new Epetra_MpiComm(Utilities::MPI::
763  duplicate_communicator(mpi_comm->GetMpiComm()));
764 #endif
765 
766  // if we don't support MPI, or if the
767  // communicator in question was in fact
768  // not an MPI communicator, return a
769  // copy of the same object again
770  Assert (dynamic_cast<const Epetra_SerialComm *>(&communicator)
771  != 0,
772  ExcInternalError());
773  return new Epetra_SerialComm(dynamic_cast<const Epetra_SerialComm &>(communicator));
774  }
775 
776 
777 
778  void destroy_communicator (Epetra_Comm &communicator)
779  {
780  // save the communicator, reset the map, and delete the communicator if
781  // this whole thing was created as an MPI communicator
782 #ifdef DEAL_II_WITH_MPI
783  Epetra_MpiComm
784  *mpi_comm = dynamic_cast<Epetra_MpiComm *>(&communicator);
785  if (mpi_comm != 0)
786  {
787  MPI_Comm comm = mpi_comm->GetMpiComm();
788  *mpi_comm = Epetra_MpiComm(MPI_COMM_SELF);
789  const int ierr = MPI_Comm_free (&comm);
790  AssertThrowMPI(ierr);
791  }
792 #endif
793  }
794 
795 
796 
797  unsigned int get_n_mpi_processes (const Epetra_Comm &mpi_communicator)
798  {
799  return mpi_communicator.NumProc();
800  }
801 
802 
803  unsigned int get_this_mpi_process (const Epetra_Comm &mpi_communicator)
804  {
805  return (unsigned int)mpi_communicator.MyPID();
806  }
807 
808 
809 
810  Epetra_Map
811  duplicate_map (const Epetra_BlockMap &map,
812  const Epetra_Comm &comm)
813  {
814  if (map.LinearMap() == true)
815  {
816  // each processor stores a
817  // contiguous range of
818  // elements in the
819  // following constructor
820  // call
821  return Epetra_Map (map.NumGlobalElements(),
822  map.NumMyElements(),
823  map.IndexBase(),
824  comm);
825  }
826  else
827  {
828  // the range is not
829  // contiguous
830  return Epetra_Map (map.NumGlobalElements(),
831  map.NumMyElements(),
832  map.MyGlobalElements (),
833  0,
834  comm);
835  }
836  }
837  }
838 
839 #endif
840 
841  template std::string to_string<int> (int, unsigned int);
842  template std::string to_string<long int> (long int, unsigned int);
843  template std::string to_string<long long int> (long long int, unsigned int);
844  template std::string to_string<unsigned int> (unsigned int, unsigned int);
845  template std::string to_string<unsigned long int> (unsigned long int, unsigned int);
846  template std::string to_string<unsigned long long int> (unsigned long long int, unsigned int);
847  template std::string to_string<float> (float, unsigned int);
848  template std::string to_string<double> (double, unsigned int);
849  template std::string to_string<long double> (long double, unsigned int);
850 
851 }
852 
853 DEAL_II_NAMESPACE_CLOSE
static const unsigned int invalid_unsigned_int
Definition: types.h:170
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:576
A class that provides a separate storage location on each thread that accesses the object...
static ::ExceptionBase & ExcIO()
unsigned int get_n_mpi_processes(const Epetra_Comm &mpi_communicator)
Definition: utilities.cc:797
std::string trim(const std::string &input)
Definition: utilities.cc:134
std::pair< int, unsigned int > get_integer_at_position(const std::string &name, const unsigned int position)
Definition: utilities.cc:409
std::vector< std::string > break_text_into_lines(const std::string &original_text, const unsigned int width, const char delimiter=' ')
Definition: utilities.cc:316
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
static ::ExceptionBase & ExcOutOfMemory()
const Epetra_Comm & comm_self()
Definition: utilities.cc:736
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
std::string get_date()
Definition: utilities.cc:675
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:92
void get_memory_stats(MemoryStats &stats)
Definition: utilities.cc:613
double string_to_double(const std::string &s)
Definition: utilities.cc:230
double get_cpu_load()
Definition: utilities.cc:604
static ::ExceptionBase & ExcMessage(std::string arg1)
Epetra_Comm * duplicate_communicator(const Epetra_Comm &communicator)
Definition: utilities.cc:752
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:564
double generate_normal_random_number(const double a, const double sigma)
Definition: utilities.cc:453
#define Assert(cond, exc)
Definition: exceptions.h:313
std::vector< std::string > split_string_list(const std::string &s, const char delimiter=',')
Definition: utilities.cc:268
unsigned long int VmSize
Definition: utilities.h:400
static ::ExceptionBase & ExcCantConvertString(std::string arg1)
static ::ExceptionBase & ExcInvalidNumber(unsigned int arg1)
bool match_at_string_start(const std::string &name, const std::string &pattern)
Definition: utilities.cc:393
void posix_memalign(void **memptr, size_t alignment, size_t size)
Definition: utilities.cc:690
std::vector< unsigned int > invert_permutation(const std::vector< unsigned int > &permutation)
Definition: utilities.cc:486
bool job_supports_mpi() 1
Definition: utilities.cc:708
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:85
std::string replace_in_string(const std::string &input, const std::string &from, const std::string &to)
Definition: utilities.cc:115
std::string dim_string(const int dim, const int spacedim)
Definition: utilities.cc:161
void destroy_communicator(Epetra_Comm &communicator)
Definition: utilities.cc:778
Definition: mpi.h:48
std::string get_hostname()
Definition: utilities.cc:646
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1211
unsigned long int VmHWM
Definition: utilities.h:405
const Epetra_Comm & comm_world()
Definition: utilities.cc:720
int string_to_int(const std::string &s)
Definition: utilities.cc:192
std::vector< unsigned int > reverse_permutation(const std::vector< unsigned int > &permutation)
Definition: utilities.cc:472
static ::ExceptionBase & ExcNotImplemented()
unsigned long int VmRSS
Definition: utilities.h:410
unsigned int get_this_mpi_process(const Epetra_Comm &mpi_communicator)
Definition: utilities.cc:803
Epetra_Map duplicate_map(const Epetra_BlockMap &map, const Epetra_Comm &comm)
Definition: utilities.cc:811
std::string get_time()
Definition: utilities.cc:660
bool job_supports_mpi()
Definition: mpi.cc:519
unsigned long int VmPeak
Definition: utilities.h:395
unsigned int needed_digits(const unsigned int max_number)
Definition: utilities.cc:171
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcInvalidNumber2StringConversersion(unsigned int arg1, unsigned int arg2)