16 #ifndef dealii_utilities_h
17 #define dealii_utilities_h
26 #include <type_traits>
33 #ifdef DEAL_II_WITH_TRILINOS
34 # include <Epetra_Comm.h>
35 # include <Epetra_Map.h>
36 # include <Teuchos_Comm.hpp>
37 # include <Teuchos_RCP.hpp>
38 # ifdef DEAL_II_WITH_MPI
39 # include <Epetra_MpiComm.h>
41 # include <Epetra_SerialComm.h>
45 #include <boost/archive/binary_iarchive.hpp>
46 #include <boost/archive/binary_oarchive.hpp>
47 #include <boost/core/demangle.hpp>
48 #include <boost/iostreams/device/array.hpp>
49 #include <boost/iostreams/device/back_inserter.hpp>
50 #include <boost/iostreams/filtering_streambuf.hpp>
51 #include <boost/serialization/array.hpp>
52 #include <boost/serialization/complex.hpp>
53 #include <boost/serialization/vector.hpp>
55 #ifdef DEAL_II_WITH_ZLIB
56 # include <boost/iostreams/filter/gzip.hpp>
65 template <
int dim,
typename Number>
103 template <
int dim,
typename Number>
104 std::vector<std::array<std::uint64_t, dim>>
107 const int bits_per_dim = 64);
113 std::vector<std::array<std::uint64_t, dim>>
115 const std::vector<std::array<std::uint64_t, dim>> &points,
116 const int bits_per_dim = 64);
136 const int bits_per_dim);
167 decompress(
const std::string &compressed_input);
183 encode_base64(
const std::vector<unsigned char> &binary_input);
193 std::vector<unsigned char>
232 template <
typename number>
252 template <
typename Number>
275 dim_string(
const int dim,
const int spacedim);
342 std::vector<std::string>
350 std::vector<std::string>
363 std::vector<std::string>
365 const unsigned int width,
366 const char delimiter =
' ');
383 std::pair<int, unsigned int>
392 const std::string &from,
393 const std::string &to);
401 trim(
const std::string &input);
450 template <
int N,
typename T>
459 template <
typename T>
461 pow(
const T base,
const int iexp)
463 #if defined(DEBUG) && !defined(DEAL_II_CXX14_CONSTEXPR_BUG)
469 abort_or_throw_on_exception,
474 "ExcMessage(\"The exponent must not be negative!\")",
475 ExcMessage(
"The exponent must not be negative!"));
498 static_assert(std::is_integral<T>::value,
"Only integral types supported");
503 (((iexp % 2 == 1) ? base : 1) *
528 template <
typename Iterator,
typename T>
538 template <
typename Iterator,
typename T,
typename Comp>
547 template <
typename Integer>
556 template <
typename Integer>
575 template <
typename T>
577 pack(
const T &
object,
578 std::vector<char> &dest_buffer,
579 const bool allow_compression =
true);
589 template <
typename T>
591 pack(
const T &
object,
const bool allow_compression =
true);
623 template <
typename T>
625 unpack(
const std::vector<char> &buffer,
const bool allow_compression =
true);
635 template <
typename T>
637 unpack(
const std::vector<char>::const_iterator &cbegin,
638 const std::vector<char>::const_iterator &cend,
639 const bool allow_compression =
true);
673 template <
typename T,
int N>
675 unpack(
const std::vector<char> &buffer,
676 T (&unpacked_object)[
N],
677 const bool allow_compression =
true);
687 template <
typename T,
int N>
689 unpack(
const std::vector<char>::const_iterator &cbegin,
690 const std::vector<char>::const_iterator &cend,
691 T (&unpacked_object)[
N],
692 const bool allow_compression =
true);
698 get_bit(
const unsigned char number,
const unsigned int n);
705 set_bit(
unsigned char &number,
const unsigned int n,
const bool x);
765 template <
typename To,
typename From>
772 template <
typename T>
779 template <
typename T>
786 template <
typename T>
793 template <
typename T>
800 template <
typename T>
929 posix_memalign(
void **memptr, std::size_t alignment, std::size_t size);
933 #ifdef DEAL_II_WITH_TRILINOS
973 const Teuchos::RCP<const Teuchos::Comm<int>> &
1080 template <
int N,
typename T>
1085 !std::is_integral<T>::value || (
N >= 0),
1087 "The non-type template parameter N must be a non-negative integer for integral type T"));
1096 fixed_power<N / 2>(x * x));
1105 return boost::core::demangle(
typeid(t).name());
1110 template <
typename Iterator,
typename T>
1119 template <
typename Iterator,
typename T,
typename Comp>
1128 "The given iterators do not satisfy the proper ordering."));
1130 unsigned int len =
static_cast<unsigned int>(last -
first);
1147 if (!comp(*
first, val))
1152 if (!comp(*
first, val))
1157 if (!comp(*
first, val))
1162 if (!comp(*
first, val))
1167 if (!comp(*
first, val))
1172 if (!comp(*
first, val))
1177 if (!comp(*
first, val))
1196 const unsigned int half = len >> 1;
1197 const Iterator middle =
first + half;
1204 if (comp(*middle, val))
1217 template <
typename T>
1220 std::vector<char> &dest_buffer,
1221 const bool allow_compression)
1223 std::size_t size = 0;
1228 #ifdef DEAL_II_HAVE_CXX17
1229 if constexpr (std::is_trivially_copyable<T>() &&
sizeof(
T) < 256)
1231 if (std::is_trivially_copyable<T>() &&
sizeof(
T) < 256)
1234 (void)allow_compression;
1235 const std::size_t previous_size = dest_buffer.size();
1236 dest_buffer.resize(previous_size +
sizeof(
T));
1238 std::memcpy(dest_buffer.data() + previous_size, &
object,
sizeof(
T));
1246 const std::size_t previous_size = dest_buffer.size();
1248 boost::iostreams::filtering_ostreambuf fosb;
1249 #ifdef DEAL_II_WITH_ZLIB
1250 if (allow_compression)
1251 fosb.push(boost::iostreams::gzip_compressor());
1253 (void)allow_compression;
1255 fosb.push(boost::iostreams::back_inserter(dest_buffer));
1257 boost::archive::binary_oarchive boa(fosb);
1262 size = dest_buffer.size() - previous_size;
1269 template <
typename T>
1271 pack(
const T &
object,
const bool allow_compression)
1273 std::vector<char> buffer;
1274 pack<T>(
object, buffer, allow_compression);
1279 template <
typename T>
1281 unpack(
const std::vector<char>::const_iterator &cbegin,
1282 const std::vector<char>::const_iterator &cend,
1283 const bool allow_compression)
1290 #ifdef DEAL_II_HAVE_CXX17
1291 if constexpr (std::is_trivially_copyable<T>() &&
sizeof(
T) < 256)
1293 if (std::is_trivially_copyable<T>() &&
sizeof(
T) < 256)
1296 (void)allow_compression;
1298 std::memcpy(&
object, &*cbegin,
sizeof(
T));
1303 boost::iostreams::filtering_istreambuf fisb;
1304 #ifdef DEAL_II_WITH_ZLIB
1305 if (allow_compression)
1306 fisb.push(boost::iostreams::gzip_decompressor());
1308 (void)allow_compression;
1310 fisb.push(boost::iostreams::array_source(&*cbegin, cend - cbegin));
1312 boost::archive::binary_iarchive bia(fisb);
1320 template <
typename T>
1322 unpack(
const std::vector<char> &buffer,
const bool allow_compression)
1324 return unpack<T>(buffer.cbegin(), buffer.cend(), allow_compression);
1328 template <
typename T,
int N>
1330 unpack(
const std::vector<char>::const_iterator &cbegin,
1331 const std::vector<char>::const_iterator &cend,
1332 T (&unpacked_object)[
N],
1333 const bool allow_compression)
1338 #ifdef DEAL_II_HAVE_CXX17
1339 if constexpr (std::is_trivially_copyable<T>() &&
sizeof(
T) *
N < 256)
1341 if (std::is_trivially_copyable<T>() &&
sizeof(
T) *
N < 256)
1344 Assert(std::distance(cbegin, cend) ==
sizeof(
T) *
N,
1346 std::memcpy(unpacked_object, &*cbegin,
sizeof(
T) *
N);
1351 boost::iostreams::filtering_istreambuf fisb;
1352 #ifdef DEAL_II_WITH_ZLIB
1353 if (allow_compression)
1354 fisb.push(boost::iostreams::gzip_decompressor());
1356 (void)allow_compression;
1358 fisb.push(boost::iostreams::array_source(&*cbegin, cend - cbegin));
1360 boost::archive::binary_iarchive bia(fisb);
1361 bia >> unpacked_object;
1366 template <
typename T,
int N>
1369 T (&unpacked_object)[
N],
1370 const bool allow_compression)
1372 unpack<T, N>(buffer.cbegin(),
1381 get_bit(
const unsigned char number,
const unsigned int n)
1388 return ((number >> n) & 1U) != 0u;
1394 set_bit(
unsigned char &number,
const unsigned int n,
const bool x)
1401 number ^= (-
static_cast<unsigned char>(x) ^ number) & (1U << n);
1406 template <
typename To,
typename From>
1407 inline std::unique_ptr<To>
1413 if (To *cast =
dynamic_cast<To *
>(p.get()))
1415 std::unique_ptr<To> result(cast);
1420 throw std::bad_cast();
1425 template <
typename T>
1434 template <
typename T>
1443 template <
typename T>
1452 template <
typename T>
1461 template <
typename T>
1470 template <
typename Integer>
1471 std::vector<Integer>
1474 const std::size_t n = permutation.size();
1476 std::vector<Integer> out(n);
1477 for (std::size_t i = 0; i < n; ++i)
1478 out[i] = n - 1 - permutation[i];
1485 template <
typename Integer>
1486 std::vector<Integer>
1489 const std::size_t n = permutation.size();
1493 for (std::size_t i = 0; i < n; ++i)
1496 out[permutation[i]] = i;
1501 for (std::size_t i = 0; i < n; ++i)
1503 ExcMessage(
"The given input permutation had duplicate entries!"));
1515 namespace serialization
1522 template <
class Archive,
typename... Args>
1524 serialize(Archive &ar, std::tuple<Args...> &t,
const unsigned int version)
1526 ar &std::get<
N - 1>(t);
1527 Serialize<N - 1>::serialize(ar, t, version);
1534 template <
class Archive,
typename... Args>
1536 serialize(Archive &ar, std::tuple<Args...> &t,
const unsigned int version)
1544 template <
class Archive,
typename... Args>
1546 serialize(Archive &ar, std::tuple<Args...> &t,
const unsigned int version)
1548 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
static ::ExceptionBase & ExcInternalError()
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
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()
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::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::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::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