Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
utilities.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2005 - 2020 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_utilities_h
17 #define dealii_utilities_h
18 
19 #include <deal.II/base/config.h>
20 
22 
23 #include <functional>
24 #include <string>
25 #include <tuple>
26 #include <type_traits>
27 #include <typeinfo>
28 #include <utility>
29 #include <vector>
30 
31 #ifdef DEAL_II_WITH_TRILINOS
32 # include <Epetra_Comm.h>
33 # include <Epetra_Map.h>
34 # include <Teuchos_Comm.hpp>
35 # include <Teuchos_RCP.hpp>
36 # ifdef DEAL_II_WITH_MPI
37 # include <Epetra_MpiComm.h>
38 # else
39 # include <Epetra_SerialComm.h>
40 # endif
41 #endif
42 
44 
45 #include <boost/archive/binary_iarchive.hpp>
46 #include <boost/archive/binary_oarchive.hpp>
47 #include <boost/core/demangle.hpp>
48 #include <boost/serialization/array.hpp>
49 #include <boost/serialization/complex.hpp>
50 #include <boost/serialization/vector.hpp>
51 
52 #ifdef DEAL_II_WITH_ZLIB
53 # include <boost/iostreams/device/back_inserter.hpp>
54 # include <boost/iostreams/filter/gzip.hpp>
55 # include <boost/iostreams/filtering_stream.hpp>
56 # include <boost/iostreams/stream.hpp>
57 #endif
58 
60 
62 
63 // forward declare Point
64 #ifndef DOXYGEN
65 template <int dim, typename Number>
66 class Point;
67 #endif
68 
77 namespace Utilities
78 {
85  std::string
87 
104  template <int dim, typename Number>
105  std::vector<std::array<std::uint64_t, dim>>
107  const std::vector<Point<dim, Number>> &points,
108  const int bits_per_dim = 64);
109 
113  template <int dim>
114  std::vector<std::array<std::uint64_t, dim>>
116  const std::vector<std::array<std::uint64_t, dim>> &points,
117  const int bits_per_dim = 64);
118 
134  template <int dim>
135  std::uint64_t
136  pack_integers(const std::array<std::uint64_t, dim> &index,
137  const int bits_per_dim);
138 
153  std::string
154  compress(const std::string &input);
155 
171  std::string
172  decompress(const std::string &compressed_input);
173 
189  std::string
190  encode_base64(const std::vector<unsigned char> &binary_input);
191 
202  std::vector<unsigned char>
203  decode_base64(const std::string &base64_input);
204 
226  std::string
227  int_to_string(const unsigned int value,
228  const unsigned int digits = numbers::invalid_unsigned_int);
229 
241  template <typename number>
242  std::string
243  to_string(const number value,
244  const unsigned int digits = numbers::invalid_unsigned_int);
245 
252  unsigned int
253  needed_digits(const unsigned int max_number);
254 
265  template <typename Number>
266  Number
267  truncate_to_n_digits(const Number number, const unsigned int n_digits);
268 
273  int
274  string_to_int(const std::string &s);
275 
287  std::string
288  dim_string(const int dim, const int spacedim);
289 
294  std::vector<int>
295  string_to_int(const std::vector<std::string> &s);
296 
301  double
302  string_to_double(const std::string &s);
303 
304 
309  std::vector<double>
310  string_to_double(const std::vector<std::string> &s);
311 
312 
355  std::vector<std::string>
356  split_string_list(const std::string &s, const std::string &delimiter = ",");
357 
358 
363  std::vector<std::string>
364  split_string_list(const std::string &s, const char delimiter);
365 
366 
376  std::vector<std::string>
377  break_text_into_lines(const std::string &original_text,
378  const unsigned int width,
379  const char delimiter = ' ');
380 
385  bool
386  match_at_string_start(const std::string &name, const std::string &pattern);
387 
396  std::pair<int, unsigned int>
397  get_integer_at_position(const std::string &name, const unsigned int position);
398 
403  std::string
404  replace_in_string(const std::string &input,
405  const std::string &from,
406  const std::string &to);
407 
413  std::string
414  trim(const std::string &input);
415 
441  double
442  generate_normal_random_number(const double a, const double sigma);
443 
453  template <class T>
454  std::string
455  type_to_string(const T &t);
456 
465  template <int N, typename T>
466  T
467  fixed_power(const T t);
468 
474  template <typename T>
475  constexpr T
476  pow(const T base, const int iexp)
477  {
478 #if defined(DEBUG) && defined(DEAL_II_HAVE_CXX14_CONSTEXPR)
479  // Up to __builtin_expect this is the same code as in the 'Assert' macro.
480  // The call to __builtin_expect turns out to be problematic.
481  if (!(iexp >= 0))
484  __FILE__,
485  __LINE__,
486  __PRETTY_FUNCTION__,
487  "iexp>=0",
488  "ExcMessage(\"The exponent must not be negative!\")",
489  ExcMessage("The exponent must not be negative!"));
490 #endif
491  // The "exponentiation by squaring" algorithm used below has to be
492  // compressed to one statement due to C++11's restrictions on constexpr
493  // functions. A more descriptive version would be:
494  //
495  // <code>
496  // if (iexp <= 0)
497  // return 1;
498  //
499  // // avoid overflow of one additional recursion with pow(base * base, 0)
500  // if (iexp == 1)
501  // return base;
502  //
503  // // if the current exponent is not divisible by two,
504  // // we need to account for that.
505  // const unsigned int prefactor = (iexp % 2 == 1) ? base : 1;
506  //
507  // // a^b = (a*a)^(b/2) for b even
508  // // a^b = a*(a*a)^((b-1)/2 for b odd
509  // return prefactor * ::Utilities::pow(base*base, iexp/2);
510  // </code>
511 
512  static_assert(std::is_integral<T>::value, "Only integral types supported");
513 
514  return iexp <= 0 ?
515  1 :
516  (iexp == 1 ? base :
517  (((iexp % 2 == 1) ? base : 1) *
518  ::Utilities::pow(base * base, iexp / 2)));
519  }
520 
542  template <typename Iterator, typename T>
543  Iterator
544  lower_bound(Iterator first, Iterator last, const T &val);
545 
546 
552  template <typename Iterator, typename T, typename Comp>
553  Iterator
554  lower_bound(Iterator first, Iterator last, const T &val, const Comp comp);
555 
561  template <typename Integer>
562  std::vector<Integer>
563  reverse_permutation(const std::vector<Integer> &permutation);
564 
570  template <typename Integer>
571  std::vector<Integer>
572  invert_permutation(const std::vector<Integer> &permutation);
573 
591  template <typename T>
592  size_t
593  pack(const T & object,
594  std::vector<char> &dest_buffer,
595  const bool allow_compression = true);
596 
607  template <typename T>
608  std::vector<char>
609  pack(const T &object, const bool allow_compression = true);
610 
643  template <typename T>
644  T
645  unpack(const std::vector<char> &buffer, const bool allow_compression = true);
646 
657  template <typename T>
658  T
659  unpack(const std::vector<char>::const_iterator &cbegin,
660  const std::vector<char>::const_iterator &cend,
661  const bool allow_compression = true);
662 
697  template <typename T, int N>
698  void
699  unpack(const std::vector<char> &buffer,
700  T (&unpacked_object)[N],
701  const bool allow_compression = true);
702 
713  template <typename T, int N>
714  void
715  unpack(const std::vector<char>::const_iterator &cbegin,
716  const std::vector<char>::const_iterator &cend,
717  T (&unpacked_object)[N],
718  const bool allow_compression = true);
719 
777  template <typename To, typename From>
778  std::unique_ptr<To>
779  dynamic_unique_cast(std::unique_ptr<From> &&p)
780  {
781  // Let's see if we can cast from 'From' to 'To'. If so, do the cast,
782  // and then release the pointer from the old
783  // owner
784  if (To *cast = dynamic_cast<To *>(p.get()))
785  {
786  std::unique_ptr<To> result(cast);
787  p.release();
788  return result;
789  }
790  else
791  throw std::bad_cast();
792  }
793 
799  namespace System
800  {
808  double
809  get_cpu_load();
810 
844  const std::string
846 
851  struct MemoryStats
852  {
856  unsigned long int VmPeak;
857 
861  unsigned long int VmSize;
862 
866  unsigned long int VmHWM;
867 
872  unsigned long int VmRSS;
873  };
874 
875 
880  void
882 
883 
887  std::string
888  get_hostname();
889 
890 
894  std::string
895  get_time();
896 
901  std::string
902  get_date();
903 
918  void
919  posix_memalign(void **memptr, std::size_t alignment, std::size_t size);
920  } // namespace System
921 
922 
923 #ifdef DEAL_II_WITH_TRILINOS
924 
929  namespace Trilinos
930  {
940  const Epetra_Comm &
941  comm_world();
942 
952  const Epetra_Comm &
953  comm_self();
954 
963  const Teuchos::RCP<const Teuchos::Comm<int>> &
965 
998  Epetra_Comm *
999  duplicate_communicator(const Epetra_Comm &communicator);
1000 
1023  void
1024  destroy_communicator(Epetra_Comm &communicator);
1025 
1034  unsigned int
1035  get_n_mpi_processes(const Epetra_Comm &mpi_communicator);
1036 
1043  unsigned int
1044  get_this_mpi_process(const Epetra_Comm &mpi_communicator);
1045 
1056  Epetra_Map
1057  duplicate_map(const Epetra_BlockMap &map, const Epetra_Comm &comm);
1058  } // namespace Trilinos
1059 
1060 #endif
1061 
1062 
1063 } // namespace Utilities
1064 
1065 
1066 // --------------------- inline functions
1067 
1068 namespace Utilities
1069 {
1070  template <int N, typename T>
1071  inline T
1072  fixed_power(const T x)
1073  {
1074  Assert(
1075  !std::is_integral<T>::value || (N >= 0),
1076  ExcMessage(
1077  "The non-type template parameter N must be a non-negative integer for integral type T"));
1078 
1079  if (N == 0)
1080  return T(1.);
1081  else if (N < 0)
1082  return T(1.) / fixed_power<-N>(x);
1083  else
1084  // Use exponentiation by squaring:
1085  return ((N % 2 == 1) ? x * fixed_power<N / 2>(x * x) :
1086  fixed_power<N / 2>(x * x));
1087  }
1088 
1089 
1090 
1091  template <class T>
1092  inline std::string
1093  type_to_string(const T &t)
1094  {
1095  return boost::core::demangle(typeid(t).name());
1096  }
1097 
1098 
1099 
1100  template <typename Iterator, typename T>
1101  inline Iterator
1102  lower_bound(Iterator first, Iterator last, const T &val)
1103  {
1104  return Utilities::lower_bound(first, last, val, std::less<T>());
1105  }
1106 
1107 
1108 
1109  template <typename Iterator, typename T, typename Comp>
1110  inline Iterator
1111  lower_bound(Iterator first, Iterator last, const T &val, const Comp comp)
1112  {
1113  // verify that the two iterators are properly ordered. since
1114  // we need operator- for the iterator type anyway, do the
1115  // test as follows, rather than via 'last >= first'
1116  Assert(last - first >= 0,
1117  ExcMessage(
1118  "The given iterators do not satisfy the proper ordering."));
1119 
1120  unsigned int len = static_cast<unsigned int>(last - first);
1121 
1122  if (len == 0)
1123  return first;
1124 
1125  while (true)
1126  {
1127  // if length equals 8 or less,
1128  // then do a rolled out
1129  // search. use a switch without
1130  // breaks for that and roll-out
1131  // the loop somehow
1132  if (len < 8)
1133  {
1134  switch (len)
1135  {
1136  case 7:
1137  if (!comp(*first, val))
1138  return first;
1139  ++first;
1141  case 6:
1142  if (!comp(*first, val))
1143  return first;
1144  ++first;
1146  case 5:
1147  if (!comp(*first, val))
1148  return first;
1149  ++first;
1151  case 4:
1152  if (!comp(*first, val))
1153  return first;
1154  ++first;
1156  case 3:
1157  if (!comp(*first, val))
1158  return first;
1159  ++first;
1161  case 2:
1162  if (!comp(*first, val))
1163  return first;
1164  ++first;
1166  case 1:
1167  if (!comp(*first, val))
1168  return first;
1169  return first + 1;
1170  default:
1171  // indices seem
1172  // to not be
1173  // sorted
1174  // correctly!? or
1175  // did len
1176  // become==0
1177  // somehow? that
1178  // shouldn't have
1179  // happened
1180  Assert(false, ExcInternalError());
1181  }
1182  }
1183 
1184 
1185 
1186  const unsigned int half = len >> 1;
1187  const Iterator middle = first + half;
1188 
1189  // if the value is larger than
1190  // that pointed to by the
1191  // middle pointer, then the
1192  // insertion point must be
1193  // right of it
1194  if (comp(*middle, val))
1195  {
1196  first = middle + 1;
1197  len -= half + 1;
1198  }
1199  else
1200  len = half;
1201  }
1202  }
1203 
1204 
1205  // --------------------- non-inline functions
1206 
1207  template <typename T>
1208  size_t
1209  pack(const T & object,
1210  std::vector<char> &dest_buffer,
1211  const bool allow_compression)
1212  {
1213  // the data is never compressed when we can't use zlib.
1214  (void)allow_compression;
1215 
1216  std::size_t size = 0;
1217 
1218  // see if the object is small and copyable via memcpy. if so, use
1219  // this fast path. otherwise, we have to go through the BOOST
1220  // serialization machinery
1221  //
1222  // we have to work around the fact that GCC 4.8.x claims to be C++
1223  // conforming, but is not actually as it does not implement
1224  // std::is_trivially_copyable.
1225 #if __GNUG__ && __GNUC__ < 5
1226  if (__has_trivial_copy(T) && sizeof(T) < 256)
1227 #else
1228 # ifdef DEAL_II_WITH_CXX17
1229  if constexpr (std::is_trivially_copyable<T>() && sizeof(T) < 256)
1230 # else
1231  if (std::is_trivially_copyable<T>() && sizeof(T) < 256)
1232 # endif
1233 #endif
1234  {
1235  const std::size_t previous_size = dest_buffer.size();
1236  dest_buffer.resize(previous_size + sizeof(T));
1237 
1238  std::memcpy(dest_buffer.data() + previous_size, &object, sizeof(T));
1239 
1240  size = sizeof(T);
1241  }
1242  else
1243  {
1244  // use buffer as the target of a compressing
1245  // stream into which we serialize the current object
1246  const std::size_t previous_size = dest_buffer.size();
1247 #ifdef DEAL_II_WITH_ZLIB
1248  if (allow_compression)
1249  {
1250  boost::iostreams::filtering_ostream out;
1251  out.push(
1252  boost::iostreams::gzip_compressor(boost::iostreams::gzip_params(
1253  boost::iostreams::gzip::default_compression)));
1254  out.push(boost::iostreams::back_inserter(dest_buffer));
1255 
1256  boost::archive::binary_oarchive archive(out);
1257  archive << object;
1258  out.flush();
1259  }
1260  else
1261 #endif
1262  {
1263  std::ostringstream out;
1264  boost::archive::binary_oarchive archive(out);
1265  archive << object;
1266 
1267  const std::string s = out.str();
1268  dest_buffer.reserve(dest_buffer.size() + s.size());
1269  std::move(s.begin(), s.end(), std::back_inserter(dest_buffer));
1270  }
1271 
1272  size = dest_buffer.size() - previous_size;
1273  }
1274 
1275  return size;
1276  }
1277 
1278 
1279  template <typename T>
1280  std::vector<char>
1281  pack(const T &object, const bool allow_compression)
1282  {
1283  std::vector<char> buffer;
1284  pack<T>(object, buffer, allow_compression);
1285  return buffer;
1286  }
1287 
1288 
1289  template <typename T>
1290  T
1291  unpack(const std::vector<char>::const_iterator &cbegin,
1292  const std::vector<char>::const_iterator &cend,
1293  const bool allow_compression)
1294  {
1295  T object;
1296 
1297  // the data is never compressed when we can't use zlib.
1298  (void)allow_compression;
1299 
1300  // see if the object is small and copyable via memcpy. if so, use
1301  // this fast path. otherwise, we have to go through the BOOST
1302  // serialization machinery
1303  //
1304  // we have to work around the fact that GCC 4.8.x claims to be C++
1305  // conforming, but is not actually as it does not implement
1306  // std::is_trivially_copyable.
1307 #if __GNUG__ && __GNUC__ < 5
1308  if (__has_trivial_copy(T) && sizeof(T) < 256)
1309 #else
1310 # ifdef DEAL_II_WITH_CXX17
1311  if constexpr (std::is_trivially_copyable<T>() && sizeof(T) < 256)
1312 # else
1313  if (std::is_trivially_copyable<T>() && sizeof(T) < 256)
1314 # endif
1315 #endif
1316  {
1317  Assert(std::distance(cbegin, cend) == sizeof(T), ExcInternalError());
1318  std::memcpy(&object, &*cbegin, sizeof(T));
1319  }
1320  else
1321  {
1322  std::string decompressed_buffer;
1323 
1324  // first decompress the buffer
1325 #ifdef DEAL_II_WITH_ZLIB
1326  if (allow_compression)
1327  {
1328  boost::iostreams::filtering_ostream decompressing_stream;
1329  decompressing_stream.push(boost::iostreams::gzip_decompressor());
1330  decompressing_stream.push(
1331  boost::iostreams::back_inserter(decompressed_buffer));
1332  decompressing_stream.write(&*cbegin, std::distance(cbegin, cend));
1333  }
1334  else
1335 #endif
1336  {
1337  decompressed_buffer.assign(cbegin, cend);
1338  }
1339 
1340  // then restore the object from the buffer
1341  std::istringstream in(decompressed_buffer);
1342  boost::archive::binary_iarchive archive(in);
1343 
1344  archive >> object;
1345  }
1346 
1347  return object;
1348  }
1349 
1350 
1351  template <typename T>
1352  T
1353  unpack(const std::vector<char> &buffer, const bool allow_compression)
1354  {
1355  return unpack<T>(buffer.cbegin(), buffer.cend(), allow_compression);
1356  }
1357 
1358 
1359  template <typename T, int N>
1360  void
1361  unpack(const std::vector<char>::const_iterator &cbegin,
1362  const std::vector<char>::const_iterator &cend,
1363  T (&unpacked_object)[N],
1364  const bool allow_compression)
1365  {
1366  // see if the object is small and copyable via memcpy. if so, use
1367  // this fast path. otherwise, we have to go through the BOOST
1368  // serialization machinery
1369  //
1370  // we have to work around the fact that GCC 4.8.x claims to be C++
1371  // conforming, but is not actually as it does not implement
1372  // std::is_trivially_copyable.
1373  if (
1374 #if __GNUG__ && __GNUC__ < 5
1375  __has_trivial_copy(T)
1376 #else
1377  std::is_trivially_copyable<T>()
1378 #endif
1379  && sizeof(T) * N < 256)
1380  {
1381  Assert(std::distance(cbegin, cend) == sizeof(T) * N,
1382  ExcInternalError());
1383  std::memcpy(unpacked_object, &*cbegin, sizeof(T) * N);
1384  }
1385  else
1386  {
1387  std::string decompressed_buffer;
1388 
1389  // first decompress the buffer
1390  (void)allow_compression;
1391 #ifdef DEAL_II_WITH_ZLIB
1392  if (allow_compression)
1393  {
1394  boost::iostreams::filtering_ostream decompressing_stream;
1395  decompressing_stream.push(boost::iostreams::gzip_decompressor());
1396  decompressing_stream.push(
1397  boost::iostreams::back_inserter(decompressed_buffer));
1398  decompressing_stream.write(&*cbegin, std::distance(cbegin, cend));
1399  }
1400  else
1401 #endif
1402  {
1403  decompressed_buffer.assign(cbegin, cend);
1404  }
1405 
1406  // then restore the object from the buffer
1407  std::istringstream in(decompressed_buffer);
1408  boost::archive::binary_iarchive archive(in);
1409 
1410  archive >> unpacked_object;
1411  }
1412  }
1413 
1414 
1415  template <typename T, int N>
1416  void
1417  unpack(const std::vector<char> &buffer,
1418  T (&unpacked_object)[N],
1419  const bool allow_compression)
1420  {
1421  unpack<T, N>(buffer.cbegin(),
1422  buffer.cend(),
1423  unpacked_object,
1424  allow_compression);
1425  }
1426 
1427 
1428 
1429  template <typename Integer>
1430  std::vector<Integer>
1431  reverse_permutation(const std::vector<Integer> &permutation)
1432  {
1433  const std::size_t n = permutation.size();
1434 
1435  std::vector<Integer> out(n);
1436  for (std::size_t i = 0; i < n; ++i)
1437  out[i] = n - 1 - permutation[i];
1438 
1439  return out;
1440  }
1441 
1442 
1443 
1444  template <typename Integer>
1445  std::vector<Integer>
1446  invert_permutation(const std::vector<Integer> &permutation)
1447  {
1448  const std::size_t n = permutation.size();
1449 
1450  std::vector<Integer> out(n, numbers::invalid_unsigned_int);
1451 
1452  for (std::size_t i = 0; i < n; ++i)
1453  {
1454  AssertIndexRange(permutation[i], n);
1455  out[permutation[i]] = i;
1456  }
1457 
1458  // check that we have actually reached
1459  // all indices
1460  for (std::size_t i = 0; i < n; ++i)
1462  ExcMessage("The given input permutation had duplicate entries!"));
1463 
1464  return out;
1465  }
1466 } // namespace Utilities
1467 
1468 
1470 
1471 #ifndef DOXYGEN
1472 namespace boost
1473 {
1474  namespace serialization
1475  {
1476  // Provides boost and c++11 with a way to serialize tuples and pairs
1477  // automatically
1478  template <int N>
1479  struct Serialize
1480  {
1481  template <class Archive, typename... Args>
1482  static void
1483  serialize(Archive &ar, std::tuple<Args...> &t, const unsigned int version)
1484  {
1485  ar &std::get<N - 1>(t);
1486  Serialize<N - 1>::serialize(ar, t, version);
1487  }
1488  };
1489 
1490  template <>
1491  struct Serialize<0>
1492  {
1493  template <class Archive, typename... Args>
1494  static void
1495  serialize(Archive &ar, std::tuple<Args...> &t, const unsigned int version)
1496  {
1497  (void)ar;
1498  (void)t;
1499  (void)version;
1500  }
1501  };
1502 
1503  template <class Archive, typename... Args>
1504  void
1505  serialize(Archive &ar, std::tuple<Args...> &t, const unsigned int version)
1506  {
1507  Serialize<sizeof...(Args)>::serialize(ar, t, version);
1508  }
1509  } // namespace serialization
1510 } // namespace boost
1511 #endif
1512 
1513 #endif
Utilities::to_string
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:483
Utilities::Trilinos::duplicate_map
Epetra_Map duplicate_map(const Epetra_BlockMap &map, const Epetra_Comm &comm)
Definition: utilities.cc:1191
Utilities::Trilinos::tpetra_comm_self
const Teuchos::RCP< const Teuchos::Comm< int > > & tpetra_comm_self()
Definition: utilities.cc:1098
Utilities::dealii_version_string
std::string dealii_version_string()
Definition: utilities.cc:97
Utilities::pow
constexpr T pow(const T base, const int iexp)
Definition: utilities.h:476
deal_II_exceptions::internals::issue_error_noreturn
void issue_error_noreturn(ExceptionHandling handling, const char *file, int line, const char *function, const char *cond, const char *exc_name, ExceptionType e)
Definition: exceptions.h:1305
Utilities::System::get_time
std::string get_time()
Definition: utilities.cc:1020
Utilities::int_to_string
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:474
Utilities::replace_in_string
std::string replace_in_string(const std::string &input, const std::string &from, const std::string &to)
Definition: utilities.cc:513
Utilities::System::MemoryStats::VmPeak
unsigned long int VmPeak
Definition: utilities.h:856
DEAL_II_FALLTHROUGH
#define DEAL_II_FALLTHROUGH
Definition: config.h:153
Utilities::System::get_memory_stats
void get_memory_stats(MemoryStats &stats)
Definition: utilities.cc:971
Utilities::decode_base64
std::vector< unsigned char > decode_base64(const std::string &base64_input)
Definition: utilities.cc:450
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
Utilities::System::get_date
std::string get_date()
Definition: utilities.cc:1036
boost
Definition: bounding_box.h:28
Utilities::pack
size_t pack(const T &object, std::vector< char > &dest_buffer, const bool allow_compression=true)
Definition: utilities.h:1209
Utilities::System::get_cpu_load
double get_cpu_load()
Definition: utilities.cc:937
deal_II_exceptions::internals::abort_or_throw_on_exception
@ abort_or_throw_on_exception
Definition: exceptions.h:1279
Utilities::get_integer_at_position
std::pair< int, unsigned int > get_integer_at_position(const std::string &name, const unsigned int position)
Definition: utilities.cc:854
Utilities::reverse_permutation
std::vector< Integer > reverse_permutation(const std::vector< Integer > &permutation)
Definition: utilities.h:1431
Utilities::inverse_Hilbert_space_filling_curve
std::vector< std::array< std::uint64_t, dim > > inverse_Hilbert_space_filling_curve(const std::vector< Point< dim, Number >> &points, const int bits_per_dim=64)
Definition: utilities.cc:148
Utilities::Trilinos::comm_world
const Epetra_Comm & comm_world()
Definition: utilities.cc:1082
Utilities::pack_integers
std::uint64_t pack_integers(const std::array< std::uint64_t, dim > &index, const int bits_per_dim)
Definition: utilities.cc:369
DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition: config.h:409
LAPACKSupport::T
static const char T
Definition: lapack_support.h:163
TransposeTableIterators::Iterator
MatrixTableIterators::Iterator< TransposeTable< T >, Constness, MatrixTableIterators::Storage::column_major > Iterator
Definition: table.h:1913
Utilities::Trilinos::destroy_communicator
void destroy_communicator(Epetra_Comm &communicator)
Definition: utilities.cc:1157
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
Utilities::dim_string
std::string dim_string(const int dim, const int spacedim)
Definition: utilities.cc:559
Utilities::trim
std::string trim(const std::string &input)
Definition: utilities.cc:532
Utilities::fixed_power
T fixed_power(const T t)
Definition: utilities.h:1072
Utilities::lower_bound
Iterator lower_bound(Iterator first, Iterator last, const T &val)
Definition: utilities.h:1102
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
Utilities::unpack
T unpack(const std::vector< char > &buffer, const bool allow_compression=true)
Definition: utilities.h:1353
Utilities::string_to_int
int string_to_int(const std::string &s)
Definition: utilities.cc:609
Utilities::encode_base64
std::string encode_base64(const std::vector< unsigned char > &binary_input)
Definition: utilities.cc:437
Utilities::Trilinos::duplicate_communicator
Epetra_Comm * duplicate_communicator(const Epetra_Comm &communicator)
Definition: utilities.cc:1130
Utilities::System::get_hostname
std::string get_hostname()
Definition: utilities.cc:1005
Utilities::System::MemoryStats::VmSize
unsigned long int VmSize
Definition: utilities.h:861
Utilities::Trilinos::comm_self
const Epetra_Comm & comm_self()
Definition: utilities.cc:1114
exceptions.h
Utilities::needed_digits
unsigned int needed_digits(const unsigned int max_number)
Definition: utilities.cc:569
value
static const bool value
Definition: dof_tools_constraints.cc:433
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
Utilities::Trilinos::get_this_mpi_process
unsigned int get_this_mpi_process(const Epetra_Comm &mpi_communicator)
Definition: utilities.cc:1183
Utilities::split_string_list
std::vector< std::string > split_string_list(const std::string &s, const std::string &delimiter=",")
Definition: utilities.cc:705
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
Utilities::truncate_to_n_digits
Number truncate_to_n_digits(const Number number, const unsigned int n_digits)
Definition: utilities.cc:582
Utilities::type_to_string
std::string type_to_string(const T &t)
Definition: utilities.h:1093
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
Utilities::generate_normal_random_number
double generate_normal_random_number(const double a, const double sigma)
Definition: utilities.cc:896
Utilities::System::MemoryStats::VmHWM
unsigned long int VmHWM
Definition: utilities.h:866
Point
Definition: point.h:111
Utilities::System::get_current_vectorization_level
const std::string get_current_vectorization_level()
Definition: utilities.cc:945
config.h
Utilities::dynamic_unique_cast
std::unique_ptr< To > dynamic_unique_cast(std::unique_ptr< From > &&p)
Definition: utilities.h:779
Utilities::Trilinos::get_n_mpi_processes
unsigned int get_n_mpi_processes(const Epetra_Comm &mpi_communicator)
Definition: utilities.cc:1176
Utilities::string_to_double
double string_to_double(const std::string &s)
Definition: utilities.cc:657
LAPACKSupport::N
static const char N
Definition: lapack_support.h:159
Utilities::match_at_string_start
bool match_at_string_start(const std::string &name, const std::string &pattern)
Definition: utilities.cc:839
Utilities::invert_permutation
std::vector< Integer > invert_permutation(const std::vector< Integer > &permutation)
Definition: utilities.h:1446
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
Utilities::System::posix_memalign
void posix_memalign(void **memptr, std::size_t alignment, std::size_t size)
Definition: utilities.cc:1051
Utilities::break_text_into_lines
std::vector< std::string > break_text_into_lines(const std::string &original_text, const unsigned int width, const char delimiter=' ')
Definition: utilities.cc:759
first
Point< 2 > first
Definition: grid_out.cc:4352
Utilities::compress
std::string compress(const std::string &input)
Definition: utilities.cc:392
Utilities
Definition: cuda.h:31
Utilities::decompress
std::string decompress(const std::string &compressed_input)
Definition: utilities.cc:415
Utilities::System::MemoryStats::VmRSS
unsigned long int VmRSS
Definition: utilities.h:872
DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition: config.h:372
Utilities::System::MemoryStats
Definition: utilities.h:851