16#ifndef dealii_utilities_h
17#define dealii_utilities_h
34#ifdef DEAL_II_WITH_TRILINOS
35# include <Epetra_Comm.h>
36# include <Epetra_Map.h>
37# include <Teuchos_Comm.hpp>
38# include <Teuchos_RCP.hpp>
39# ifdef DEAL_II_WITH_MPI
40# include <Epetra_MpiComm.h>
42# include <Epetra_SerialComm.h>
46#include <boost/archive/binary_iarchive.hpp>
47#include <boost/archive/binary_oarchive.hpp>
48#include <boost/core/demangle.hpp>
49#include <boost/iostreams/device/array.hpp>
50#include <boost/iostreams/device/back_inserter.hpp>
51#include <boost/iostreams/filtering_streambuf.hpp>
52#include <boost/serialization/array.hpp>
53#include <boost/serialization/complex.hpp>
54#include <boost/serialization/vector.hpp>
56#ifdef DEAL_II_WITH_ZLIB
57# include <boost/iostreams/filter/gzip.hpp>
66template <
int dim,
typename Number>
104 template <
int dim,
typename Number>
105 std::vector<std::array<std::uint64_t, dim>>
108 const int bits_per_dim = 64);
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);
137 const int bits_per_dim);
168 decompress(
const std::string &compressed_input);
184 encode_base64(
const std::vector<unsigned char> &binary_input);
194 std::vector<unsigned char>
233 template <
typename number>
253 template <
typename Number>
276 dim_string(
const int dim,
const int spacedim);
343 std::vector<std::string>
351 std::vector<std::string>
364 std::vector<std::string>
366 const unsigned int width,
367 const char delimiter =
' ');
384 std::pair<int, unsigned int>
393 const std::string &from,
394 const std::string &to);
402 trim(
const std::string &input);
451 template <
int N,
typename T>
460 template <
typename T>
462 pow(
const T base,
const int iexp)
464#if defined(DEBUG) && !defined(DEAL_II_CXX14_CONSTEXPR_BUG)
470 abort_or_throw_on_exception,
475 "ExcMessage(\"The exponent must not be negative!\")",
476 ExcMessage(
"The exponent must not be negative!"));
499 static_assert(std::is_integral<T>::value,
"Only integral types supported");
504 (((iexp % 2 == 1) ? base : 1) *
529 template <
typename Iterator,
typename T>
539 template <
typename Iterator,
typename T,
typename Comp>
548 template <
typename Integer>
557 template <
typename Integer>
615 template <
typename T>
617 pack(
const T &
object,
618 std::vector<char> &dest_buffer,
619 const bool allow_compression =
true);
629 template <
typename T>
631 pack(
const T &
object,
const bool allow_compression =
true);
664 template <
typename T>
666 unpack(
const std::vector<char> &buffer,
const bool allow_compression =
true);
676 template <
typename T>
678 unpack(
const std::vector<char>::const_iterator &cbegin,
679 const std::vector<char>::const_iterator &cend,
680 const bool allow_compression =
true);
714 template <
typename T,
int N>
716 unpack(
const std::vector<char> &buffer,
717 T (&unpacked_object)[N],
718 const bool allow_compression =
true);
728 template <
typename T,
int N>
730 unpack(
const std::vector<char>::const_iterator &cbegin,
731 const std::vector<char>::const_iterator &cend,
732 T (&unpacked_object)[N],
733 const bool allow_compression =
true);
739 get_bit(
const unsigned char number,
const unsigned int n);
746 set_bit(
unsigned char &number,
const unsigned int n,
const bool x);
806 template <
typename To,
typename From>
813 template <
typename T>
820 template <
typename T>
827 template <
typename T>
834 template <
typename T>
841 template <
typename T>
970 posix_memalign(
void **memptr, std::size_t alignment, std::size_t size);
974#ifdef DEAL_II_WITH_TRILINOS
1014 const Teuchos::RCP<const Teuchos::Comm<int>> &
1121 template <
int N,
typename T>
1126 !std::is_integral<T>::value || (N >= 0),
1128 "The non-type template parameter N must be a non-negative integer for integral type T"));
1136 return ((N % 2 == 1) ? x *
fixed_power<N / 2>(x * x) :
1137 fixed_power<N / 2>(x * x));
1146 return boost::core::demangle(
typeid(t).name());
1151 template <
typename Iterator,
typename T>
1160 template <
typename Iterator,
typename T,
typename Comp>
1169 "The given iterators do not satisfy the proper ordering."));
1171 unsigned int len =
static_cast<unsigned int>(last -
first);
1188 if (!comp(*
first, val))
1193 if (!comp(*
first, val))
1198 if (!comp(*
first, val))
1203 if (!comp(*
first, val))
1208 if (!comp(*
first, val))
1213 if (!comp(*
first, val))
1218 if (!comp(*
first, val))
1237 const unsigned int half = len >> 1;
1238 const Iterator middle =
first + half;
1245 if (comp(*middle, val))
1265 template <
typename T>
1273 template <
typename T>
1277 std::is_trivially_copyable<T>::value && !std::is_same<T, bool>::value;
1282 template <
typename T>
1286 std::is_trivially_copyable<T>::value && !std::is_same<T, bool>::value;
1299 template <
typename T>
1302 std::vector<char> &)
1310 template <
typename T,
1311 typename = std::enable_if_t<!std::is_same<T, bool>::value &&
1312 std::is_trivially_copyable<T>::value>>
1315 const std::vector<T> &
object,
1316 std::vector<char> & dest_buffer)
1318 const typename std::vector<T>::size_type vector_size =
object.size();
1322 dest_buffer.reserve(dest_buffer.size() +
sizeof(vector_size) +
1323 vector_size *
sizeof(T));
1326 dest_buffer.insert(dest_buffer.end(),
1327 reinterpret_cast<const char *
>(&vector_size),
1328 reinterpret_cast<const char *
>(&vector_size + 1));
1331 if (vector_size > 0)
1332 dest_buffer.insert(dest_buffer.end(),
1333 reinterpret_cast<const char *
>(
object.data()),
1334 reinterpret_cast<const char *
>(
object.data() +
1340 template <
typename T,
1341 typename = std::enable_if_t<!std::is_same<T, bool>::value &&
1342 std::is_trivially_copyable<T>::value>>
1345 const std::vector<std::vector<T>> &
object,
1346 std::vector<char> & dest_buffer)
1348 using size_type =
typename std::vector<T>::size_type;
1349 const size_type vector_size =
object.size();
1351 typename std::vector<T>::size_type aggregated_size = 0;
1352 std::vector<size_type> sizes;
1353 sizes.reserve(vector_size);
1354 for (
const auto &a :
object)
1356 aggregated_size += a.size();
1357 sizes.push_back(a.size());
1362 dest_buffer.reserve(dest_buffer.size() +
1363 sizeof(vector_size) * (1 + vector_size) +
1364 aggregated_size *
sizeof(T));
1367 dest_buffer.insert(dest_buffer.end(),
1368 reinterpret_cast<const char *
>(&vector_size),
1369 reinterpret_cast<const char *
>(&vector_size + 1));
1372 if (vector_size > 0)
1373 dest_buffer.insert(dest_buffer.end(),
1374 reinterpret_cast<const char *
>(sizes.data()),
1375 reinterpret_cast<const char *
>(sizes.data() +
1379 for (
const auto &a :
object)
1380 dest_buffer.insert(dest_buffer.end(),
1381 reinterpret_cast<const char *
>(a.data()),
1382 reinterpret_cast<const char *
>(a.data() + a.size()));
1387 template <
typename T>
1390 const std::vector<char>::const_iterator &,
1391 const std::vector<char>::const_iterator &,
1400 template <
typename T,
1401 typename = std::enable_if_t<!std::is_same<T, bool>::value &&
1402 std::is_trivially_copyable<T>::value>>
1405 const std::vector<char>::const_iterator &cbegin,
1406 const std::vector<char>::const_iterator &cend,
1407 std::vector<T> &
object)
1410 typename std::vector<T>::size_type vector_size;
1411 memcpy(&vector_size, &*cbegin,
sizeof(vector_size));
1413 Assert(
static_cast<std::ptrdiff_t
>(cend - cbegin) ==
1414 static_cast<std::ptrdiff_t
>(
sizeof(vector_size) +
1415 vector_size *
sizeof(T)),
1416 ExcMessage(
"The given buffer has the wrong size."));
1421 if (vector_size > 0)
1422 object.insert(
object.end(),
1423 reinterpret_cast<const T *
>(&*cbegin +
1424 sizeof(vector_size)),
1425 reinterpret_cast<const T *
>(&*cend));
1430 template <
typename T,
1431 typename = std::enable_if_t<!std::is_same<T, bool>::value &&
1432 std::is_trivially_copyable<T>::value>>
1435 const std::vector<char>::const_iterator &cbegin,
1436 const std::vector<char>::const_iterator &cend,
1437 std::vector<std::vector<T>> &
object)
1440 using size_type =
typename std::vector<T>::size_type;
1441 std::vector<char>::const_iterator iterator = cbegin;
1442 size_type vector_size;
1443 memcpy(&vector_size, &*iterator,
sizeof(vector_size));
1445 object.resize(vector_size);
1446 std::vector<size_type> sizes(vector_size);
1447 if (vector_size > 0)
1448 memcpy(sizes.data(),
1449 &*iterator +
sizeof(vector_size),
1450 vector_size *
sizeof(size_type));
1452 iterator +=
sizeof(vector_size) * (1 + vector_size);
1453 size_type aggregated_size = 0;
1454 for (
const auto a : sizes)
1455 aggregated_size += a;
1457 Assert(
static_cast<std::ptrdiff_t
>(cend - iterator) ==
1458 static_cast<std::ptrdiff_t
>(aggregated_size *
sizeof(T)),
1459 ExcMessage(
"The given buffer has the wrong size."));
1463 for (
unsigned int i = 0; i < vector_size; ++i)
1466 object[i].insert(
object[i].end(),
1467 reinterpret_cast<const T *
>(&*iterator),
1468 reinterpret_cast<const T *
>(&*iterator +
1469 sizeof(T) * sizes[i]));
1470 iterator +=
sizeof(T) * sizes[i];
1474 ExcMessage(
"The given buffer has the wrong size."));
1481 template <
typename T>
1484 std::vector<char> &dest_buffer,
1485 const bool allow_compression)
1487 std::size_t size = 0;
1493#ifdef DEAL_II_HAVE_CXX17
1494 if constexpr (std::is_trivially_copyable<T>() &&
sizeof(T) < 256)
1496 if (std::is_trivially_copyable<T>() &&
sizeof(T) < 256)
1505 size = (std::is_same<T, std::tuple<>>::value ? 0 :
sizeof(T));
1507 (void)allow_compression;
1508 const std::size_t previous_size = dest_buffer.size();
1509 dest_buffer.resize(previous_size + size);
1512 std::memcpy(dest_buffer.data() + previous_size, &
object, size);
1520 (allow_compression ==
false))
1522 const std::size_t previous_size = dest_buffer.size();
1531 size = dest_buffer.size() - previous_size;
1537 const std::size_t previous_size = dest_buffer.size();
1539 boost::iostreams::filtering_ostreambuf fosb;
1540#ifdef DEAL_II_WITH_ZLIB
1541 if (allow_compression)
1542 fosb.push(boost::iostreams::gzip_compressor());
1544 (void)allow_compression;
1546 fosb.push(boost::iostreams::back_inserter(dest_buffer));
1548 boost::archive::binary_oarchive boa(fosb);
1553 size = dest_buffer.size() - previous_size;
1560 template <
typename T>
1562 pack(
const T &
object,
const bool allow_compression)
1564 std::vector<char> buffer;
1565 pack<T>(
object, buffer, allow_compression);
1571 template <
typename T>
1573 unpack(
const std::vector<char>::const_iterator &cbegin,
1574 const std::vector<char>::const_iterator &cend,
1575 const bool allow_compression)
1580#ifdef DEAL_II_HAVE_CXX17
1581 if constexpr (std::is_trivially_copyable<T>() &&
sizeof(T) < 256)
1583 if (std::is_trivially_copyable<T>() &&
sizeof(T) < 256)
1592 const std::size_t size =
1593 (std::is_same<T, std::tuple<>>::value ? 0 :
sizeof(T));
1597 (void)allow_compression;
1601 std::memcpy(&
object, &*cbegin, size);
1611 (allow_compression ==
false))
1626 boost::iostreams::filtering_istreambuf fisb;
1627#ifdef DEAL_II_WITH_ZLIB
1628 if (allow_compression)
1629 fisb.push(boost::iostreams::gzip_decompressor());
1631 (void)allow_compression;
1633 fisb.push(boost::iostreams::array_source(&*cbegin, cend - cbegin));
1635 boost::archive::binary_iarchive bia(fisb);
1646 template <
typename T>
1648 unpack(
const std::vector<char> &buffer,
const bool allow_compression)
1650 return unpack<T>(buffer.cbegin(), buffer.cend(), allow_compression);
1654 template <
typename T,
int N>
1656 unpack(
const std::vector<char>::const_iterator &cbegin,
1657 const std::vector<char>::const_iterator &cend,
1658 T (&unpacked_object)[N],
1659 const bool allow_compression)
1664#ifdef DEAL_II_HAVE_CXX17
1665 if constexpr (std::is_trivially_copyable<T>() &&
sizeof(T) * N < 256)
1667 if (std::is_trivially_copyable<T>() &&
sizeof(T) * N < 256)
1670 Assert(std::distance(cbegin, cend) ==
sizeof(T) * N,
1672 std::memcpy(unpacked_object, &*cbegin,
sizeof(T) * N);
1677 boost::iostreams::filtering_istreambuf fisb;
1678#ifdef DEAL_II_WITH_ZLIB
1679 if (allow_compression)
1680 fisb.push(boost::iostreams::gzip_decompressor());
1682 (void)allow_compression;
1684 fisb.push(boost::iostreams::array_source(&*cbegin, cend - cbegin));
1686 boost::archive::binary_iarchive bia(fisb);
1687 bia >> unpacked_object;
1692 template <
typename T,
int N>
1695 T (&unpacked_object)[N],
1696 const bool allow_compression)
1698 unpack<T, N>(buffer.cbegin(),
1707 get_bit(
const unsigned char number,
const unsigned int n)
1714 return ((number >> n) & 1U) != 0u;
1720 set_bit(
unsigned char &number,
const unsigned int n,
const bool x)
1727 number ^= (-
static_cast<unsigned char>(x) ^ number) & (1U << n);
1732 template <
typename To,
typename From>
1733 inline std::unique_ptr<To>
1739 if (To *cast =
dynamic_cast<To *
>(p.get()))
1741 std::unique_ptr<To> result(cast);
1746 throw std::bad_cast();
1751 template <
typename T>
1760 template <
typename T>
1769 template <
typename T>
1778 template <
typename T>
1787 template <
typename T>
1796 template <
typename Integer>
1797 std::vector<Integer>
1800 const std::size_t n = permutation.size();
1802 std::vector<Integer> out(n);
1803 for (std::size_t i = 0; i < n; ++i)
1804 out[i] = n - 1 - permutation[i];
1811 template <
typename Integer>
1812 std::vector<Integer>
1815 const std::size_t n = permutation.size();
1819 for (std::size_t i = 0; i < n; ++i)
1822 out[permutation[i]] = i;
1827 for (std::size_t i = 0; i < n; ++i)
1829 ExcMessage(
"The given input permutation had duplicate entries!"));
1841 namespace serialization
1848 template <
class Archive,
typename... Args>
1850 serialize(Archive &ar, std::tuple<Args...> &t,
const unsigned int version)
1852 ar &std::get<N - 1>(t);
1853 Serialize<N - 1>::serialize(ar, t, version);
1860 template <
class Archive,
typename... Args>
1862 serialize(Archive &ar, std::tuple<Args...> &t,
const unsigned int version)
1870 template <
class Archive,
typename... Args>
1872 serialize(Archive &ar, std::tuple<Args...> &t,
const unsigned int version)
1874 Serialize<
sizeof...(Args)>::serialize(ar, t, version);
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_FALLTHROUGH
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
void get_memory_stats(MemoryStats &stats)
std::string get_hostname()
void posix_memalign(void **memptr, std::size_t alignment, std::size_t size)
const std::string get_current_vectorization_level()
void destroy_communicator(Epetra_Comm &communicator)
Epetra_Map duplicate_map(const Epetra_BlockMap &map, const Epetra_Comm &comm)
unsigned int get_this_mpi_process(const Epetra_Comm &mpi_communicator)
unsigned int get_n_mpi_processes(const Epetra_Comm &mpi_communicator)
const Epetra_Comm & comm_self()
const Teuchos::RCP< const Teuchos::Comm< int > > & tpetra_comm_self()
Epetra_Comm * duplicate_communicator(const Epetra_Comm &communicator)
const Epetra_Comm & comm_world()
void create_vector_of_trivially_copyable_from_buffer(const std::vector< char >::const_iterator &, const std::vector< char >::const_iterator &, T &)
void append_vector_of_trivially_copyable_to_buffer(const T &, std::vector< char > &)
constexpr T pow(const T base, const int iexp)
std::vector< std::string > split_string_list(const std::string &s, const std::string &delimiter=",")
T & get_underlying_value(T &p)
Number truncate_to_n_digits(const Number number, const unsigned int n_digits)
std::string type_to_string(const T &t)
size_t pack(const T &object, std::vector< char > &dest_buffer, const bool allow_compression=true)
std::string dim_string(const int dim, const int spacedim)
std::uint64_t pack_integers(const std::array< std::uint64_t, dim > &index, const int bits_per_dim)
bool get_bit(const unsigned char number, const unsigned int n)
std::pair< int, unsigned int > get_integer_at_position(const std::string &name, const unsigned int position)
std::string encode_base64(const std::vector< unsigned char > &binary_input)
std::unique_ptr< To > dynamic_unique_cast(std::unique_ptr< From > &&p)
std::string replace_in_string(const std::string &input, const std::string &from, const std::string &to)
std::vector< unsigned char > decode_base64(const std::string &base64_input)
std::vector< std::string > break_text_into_lines(const std::string &original_text, const unsigned int width, const char delimiter=' ')
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
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)
std::string compress(const std::string &input)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
bool match_at_string_start(const std::string &name, const std::string &pattern)
T unpack(const std::vector< char > &buffer, const bool allow_compression=true)
std::string decompress(const std::string &compressed_input)
unsigned int needed_digits(const unsigned int max_number)
Iterator lower_bound(Iterator first, Iterator last, const T &val)
double string_to_double(const std::string &s)
std::string dealii_version_string()
void set_bit(unsigned char &number, const unsigned int n, const bool x)
std::string trim(const std::string &input)
double generate_normal_random_number(const double a, const double sigma)
std::vector< Integer > reverse_permutation(const std::vector< Integer > &permutation)
std::vector< Integer > invert_permutation(const std::vector< Integer > &permutation)
int string_to_int(const std::string &s)
void issue_error_noreturn(ExceptionHandling handling, const char *file, int line, const char *function, const char *cond, const char *exc_name, ExceptionType e)
static const unsigned int invalid_unsigned_int
static constexpr bool value