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