32#define BOOST_BIND_GLOBAL_PLACEHOLDERS
33#include <boost/archive/iterators/base64_from_binary.hpp>
34#include <boost/archive/iterators/binary_from_base64.hpp>
35#include <boost/archive/iterators/transform_width.hpp>
36#include <boost/iostreams/copy.hpp>
37#include <boost/lexical_cast.hpp>
38#include <boost/random.hpp>
39#undef BOOST_BIND_GLOBAL_PLACEHOLDERS
57#if defined(DEAL_II_HAVE_UNISTD_H) && defined(DEAL_II_HAVE_GETHOSTNAME)
66#ifdef DEAL_II_WITH_TRILINOS
67# ifdef DEAL_II_WITH_MPI
72# include <Epetra_MpiComm.h>
73# include <Teuchos_DefaultComm.hpp>
75# include <Epetra_SerialComm.h>
76# include <Teuchos_RCP.hpp>
87 <<
"When trying to convert " << arg1 <<
" to a string with "
88 << arg2 <<
" digits");
92 <<
"Can't convert the string " << arg1
93 <<
" to the desired type");
110 std::vector<std::array<std::uint64_t, effective_dim>>
111 inverse_Hilbert_space_filling_curve_effective(
114 const std::array<LongDouble, dim> & extents,
115 const std::bitset<dim> & valid_extents,
117 const Integer max_int)
119 std::vector<std::array<Integer, effective_dim>> int_points(points.size());
121 for (
unsigned int i = 0; i < points.size(); ++i)
124 unsigned int eff_d = 0;
125 for (
unsigned int d = 0; d < dim; ++d)
126 if (valid_extents[d])
129 const LongDouble v = (
static_cast<LongDouble
>(points[i][d]) -
130 static_cast<LongDouble
>(bl[d])) /
134 int_points[i][eff_d] =
135 static_cast<Integer
>(v *
static_cast<LongDouble
>(max_int));
141 return inverse_Hilbert_space_filling_curve<effective_dim>(int_points,
146 template <
int dim,
typename Number>
147 std::vector<std::array<std::uint64_t, dim>>
150 const int bits_per_dim)
152 using Integer = std::uint64_t;
154 using LongDouble =
long double;
157 if (points.size() == 0)
158 return std::vector<std::array<std::uint64_t, dim>>();
162 for (
const auto &p : points)
163 for (
unsigned int d = 0; d < dim; ++d)
165 const double cid = p[d];
170 std::array<LongDouble, dim> extents;
171 std::bitset<dim> valid_extents;
172 for (
unsigned int i = 0; i < dim; ++i)
175 static_cast<LongDouble
>(tr[i]) -
static_cast<LongDouble
>(bl[i]);
176 valid_extents[i] = (extents[i] > 0.);
183 std::min(std::numeric_limits<Integer>::digits,
184 std::numeric_limits<LongDouble>::digits));
187 const Integer max_int = (min_bits == std::numeric_limits<Integer>::digits ?
188 std::numeric_limits<Integer>::max() :
189 (Integer(1) << min_bits) - 1);
191 const unsigned int effective_dim = valid_extents.count();
192 if (effective_dim == dim)
194 return inverse_Hilbert_space_filling_curve_effective<dim,
199 points, bl, extents, valid_extents, min_bits, max_int);
203 std::array<std::uint64_t, dim> zero_ind;
204 for (
unsigned int d = 0; d < dim; ++d)
207 std::vector<std::array<std::uint64_t, dim>> ind(points.size(), zero_ind);
209 if (dim == 3 && effective_dim == 2)
212 inverse_Hilbert_space_filling_curve_effective<dim,
217 points, bl, extents, valid_extents, min_bits, max_int);
219 for (
unsigned int i = 0; i < ind.size(); ++i)
220 for (
unsigned int d = 0; d < 2; ++d)
221 ind[i][d + 1] = ind2[i][d];
225 else if (effective_dim == 1)
228 inverse_Hilbert_space_filling_curve_effective<dim,
233 points, bl, extents, valid_extents, min_bits, max_int);
235 for (
unsigned int i = 0; i < ind.size(); ++i)
236 ind[i][dim - 1] = ind1[i][0];
247 for (
unsigned int i = 0; i < points.size(); ++i)
256 std::vector<std::array<std::uint64_t, dim>>
258 const std::vector<std::array<std::uint64_t, dim>> &points,
259 const int bits_per_dim)
261 using Integer = std::uint64_t;
263 std::vector<std::array<Integer, dim>> int_points(points);
265 std::vector<std::array<Integer, dim>> res(int_points.size());
288 Assert(bits_per_dim <= std::numeric_limits<Integer>::digits,
289 ExcMessage(
"This integer type can not hold " +
290 std::to_string(bits_per_dim) +
" bits."));
292 const Integer M = Integer(1) << (bits_per_dim - 1);
294 for (
unsigned int index = 0; index < int_points.size(); ++index)
296 auto &X = int_points[index];
297 auto &L = res[index];
300 for (Integer q = M; q > 1; q >>= 1)
302 const Integer p = q - 1;
303 for (
unsigned int i = 0; i < dim; ++i)
313 const Integer t = (X[0] ^ X[i]) & p;
321 for (
unsigned int i = 1; i < dim; ++i)
325 for (Integer q = M; q > 1; q >>= 1)
328 for (
unsigned int i = 0; i < dim; ++i)
343 for (
unsigned int i = 0; i < dim; ++i)
347 for (Integer q = M; q > 0; q >>= 1)
370 const int bits_per_dim)
372 using Integer = std::uint64_t;
377 const Integer mask = (Integer(1) << bits_per_dim) - 1;
380 for (
unsigned int i = 0; i < dim; ++i)
383 const Integer v = (mask & index[dim - 1 - i]) << (bits_per_dim * i);
394#ifdef DEAL_II_WITH_ZLIB
395 namespace bio = boost::iostreams;
397 std::stringstream compressed;
398 std::stringstream origin(input);
400 bio::filtering_streambuf<bio::input> out;
401 out.push(bio::gzip_compressor());
403 bio::copy(out, compressed);
405 return compressed.str();
416#ifdef DEAL_II_WITH_ZLIB
417 namespace bio = boost::iostreams;
419 std::stringstream compressed(compressed_input);
420 std::stringstream decompressed;
422 bio::filtering_streambuf<bio::input> out;
423 out.push(bio::gzip_decompressor());
424 out.push(compressed);
425 bio::copy(out, decompressed);
427 return decompressed.str();
429 return compressed_input;
438 using namespace boost::archive::iterators;
439 using It = base64_from_binary<
440 transform_width<std::vector<unsigned char>::const_iterator, 6, 8>>;
441 auto base64 = std::string(It(binary_input.begin()), It(binary_input.end()));
443 return base64.append((3 - binary_input.size() % 3) % 3,
'=');
448 std::vector<unsigned char>
451 using namespace boost::archive::iterators;
453 transform_width<binary_from_base64<std::string::const_iterator>, 8, 6>;
454 auto binary = std::vector<unsigned char>(It(base64_input.begin()),
455 It(base64_input.end()));
457 auto length = base64_input.size();
458 if (binary.size() > 2 && base64_input[length - 1] ==
'=' &&
459 base64_input[length - 2] ==
'=')
461 binary.erase(binary.end() - 2, binary.end());
463 else if (binary.size() > 1 && base64_input[length - 1] ==
'=')
465 binary.erase(binary.end() - 1, binary.end());
480 template <
typename number>
482 to_string(
const number value,
const unsigned int digits)
492 std::string lc_string = (std::is_integral<number>::value ?
493 std::to_string(value) :
494 boost::lexical_cast<std::string>(value));
497 (lc_string.size() < digits))
500 const unsigned int padding_position = (lc_string[0] ==
'-') ? 1 : 0;
502 const std::string padding(digits - lc_string.size(),
'0');
503 lc_string.insert(padding_position, padding);
513 const std::string &from,
514 const std::string &to)
519 std::string out = input;
520 std::string::size_type pos = out.find(from);
522 while (pos != std::string::npos)
524 out.replace(pos, from.size(), to);
525 pos = out.find(from, pos + to.size());
533 std::string::size_type left = 0;
534 std::string::size_type right = input.size() > 0 ? input.size() - 1 : 0;
536 for (; left < input.size(); ++left)
538 if (std::isspace(input[left]) == 0)
544 for (; right >= left; --right)
546 if (std::isspace(input[right]) == 0)
552 return std::string(input, left, right - left + 1);
571 return static_cast<int>(
572 std::ceil(std::log10(std::fabs(max_number + 0.1))));
579 template <
typename Number>
585 if (!(std::fabs(number) > std::numeric_limits<Number>::min()))
589 static_cast<int>(std::floor(std::log10(std::fabs(number))));
591 const int shift = -order +
static_cast<int>(n_digits) - 1;
593 Assert(shift <=
static_cast<int>(std::floor(
594 std::log10(std::numeric_limits<Number>::max()))),
596 "Overflow. Use a smaller value for n_digits and/or make sure "
597 "that the absolute value of 'number' does not become too small."));
599 const Number factor =
std::pow(10.0,
static_cast<Number
>(shift));
601 const Number number_cutoff = std::trunc(number * factor) / factor;
603 return number_cutoff;
612 while ((s.size() > 0) && (s[0] ==
' '))
614 while ((s.size() > 0) && (s.back() ==
' '))
615 s.erase(s.end() - 1);
625 const int i = std::strtol(s.c_str(), &p, 10);
636 ((s.size() > 0) && (*p !=
'\0'))),
637 ExcMessage(
"Can't convert <" + s +
"> to an integer."));
647 std::vector<int> tmp(s.size());
648 for (
unsigned int i = 0; i < s.size(); ++i)
660 while ((s.size() > 0) && (s[0] ==
' '))
662 while ((s.size() > 0) && (s.back() ==
' '))
663 s.erase(s.end() - 1);
673 const double d = std::strtod(s.c_str(), &p);
684 ((s.size() > 0) && (*p !=
'\0'))),
685 ExcMessage(
"Can't convert <" + s +
"> to a double."));
695 std::vector<double> tmp(s.size());
696 for (
unsigned int i = 0; i < s.size(); ++i)
703 std::vector<std::string>
712 while (tmp.size() != 0 && tmp.back() ==
' ')
713 tmp.erase(tmp.size() - 1, 1);
721 std::vector<std::string> split_list;
722 while (tmp.size() != 0)
727 if (name.find(delimiter) != std::string::npos)
729 name.erase(name.find(delimiter), std::string::npos);
730 tmp.erase(0, tmp.find(delimiter) + delimiter.size());
736 while ((name.size() != 0) && (name[0] ==
' '))
738 while (name.size() != 0 && name.back() ==
' ')
739 name.erase(name.size() - 1, 1);
741 split_list.push_back(name);
748 std::vector<std::string>
757 std::vector<std::string>
759 const unsigned int width,
760 const char delimiter)
762 std::string text = original_text;
763 std::vector<std::string> lines;
766 while ((text.size() != 0) && (text.back() == delimiter))
767 text.erase(text.size() - 1, 1);
770 while (text.size() != 0)
774 while ((text.size() != 0) && (text[0] == delimiter))
777 std::size_t pos_newline = text.find_first_of(
'\n', 0);
778 if (pos_newline != std::string::npos && pos_newline <= width)
780 std::string line(text, 0, pos_newline);
781 while ((line.size() != 0) && (line.back() == delimiter))
782 line.erase(line.size() - 1, 1);
783 lines.push_back(line);
784 text.erase(0, pos_newline + 1);
791 if (text.size() < width)
794 while ((text.size() != 0) && (text.back() == delimiter))
795 text.erase(text.size() - 1, 1);
796 lines.push_back(text);
804 int location = std::min<int>(width, text.size() - 1);
805 for (; location > 0; --location)
806 if (text[location] == delimiter)
812 for (location = std::min<int>(width, text.size() - 1);
813 location <
static_cast<int>(text.size());
815 if (text[location] == delimiter)
821 std::string line(text, 0, location);
822 while ((line.size() != 0) && (line.back() == delimiter))
823 line.erase(line.size() - 1, 1);
824 lines.push_back(line);
825 text.erase(0, location);
837 if (pattern.size() > name.size())
840 for (
unsigned int i = 0; i < pattern.size(); ++i)
841 if (pattern[i] != name[i])
849 std::pair<int, unsigned int>
854 const std::string test_string(name.begin() + position, name.end());
856 std::istringstream str(test_string);
866 return std::make_pair(i, 1U);
868 return std::make_pair(i, 2U);
870 return std::make_pair(i, 3U);
872 return std::make_pair(i, 4U);
874 return std::make_pair(i, 5U);
875 else if (i < 1000000)
876 return std::make_pair(i, 6U);
877 else if (i < 10000000)
878 return std::make_pair(i, 7U);
906 return boost::normal_distribution<>(a,
907 sigma)(random_number_generator.
get());
919 std::ifstream cpuinfo;
920 cpuinfo.open(
"/proc/loadavg");
960 "Invalid DEAL_II_VECTORIZATION_WIDTH_IN_BITS."));
975 std::ifstream file(
"/proc/self/status");
981 if (name ==
"VmPeak:")
983 else if (name ==
"VmSize:")
985 else if (name ==
"VmHWM:")
987 else if (name ==
"VmRSS:")
1003#if defined(DEAL_II_HAVE_UNISTD_H) && defined(DEAL_II_HAVE_GETHOSTNAME)
1004 const unsigned int N = 1024;
1006 gethostname(&(hostname[0]), N - 1);
1008 std::string hostname(
"unknown");
1018 std::time_t time1 = std::time(
nullptr);
1019 std::tm * time = std::localtime(&time1);
1021 std::ostringstream o;
1022 o << time->tm_hour <<
":" << (time->tm_min < 10 ?
"0" :
"")
1023 << time->tm_min <<
":" << (time->tm_sec < 10 ?
"0" :
"")
1034 std::time_t time1 = std::time(
nullptr);
1035 std::tm * time = std::localtime(&time1);
1037 std::ostringstream o;
1038 o << time->tm_year + 1900 <<
"/" << time->tm_mon + 1 <<
"/"
1057 *memptr = malloc(size);
1073#ifdef DEAL_II_WITH_TRILINOS
1080# ifdef DEAL_II_WITH_MPI
1081 static Teuchos::RCP<Epetra_MpiComm> communicator =
1082 Teuchos::rcp(
new Epetra_MpiComm(MPI_COMM_WORLD),
true);
1084 static Teuchos::RCP<Epetra_SerialComm> communicator =
1085 Teuchos::rcp(
new Epetra_SerialComm(),
true);
1088 return *communicator;
1093 const Teuchos::RCP<const Teuchos::Comm<int>> &
1096# ifdef DEAL_II_WITH_MPI
1097 static auto communicator = Teuchos::RCP<const Teuchos::Comm<int>>(
1098 new Teuchos::MpiComm<int>(MPI_COMM_SELF));
1100 static auto communicator =
1101 Teuchos::RCP<const Teuchos::Comm<int>>(
new Teuchos::Comm<int>());
1104 return communicator;
1112# ifdef DEAL_II_WITH_MPI
1113 static Teuchos::RCP<Epetra_MpiComm> communicator =
1114 Teuchos::rcp(
new Epetra_MpiComm(MPI_COMM_SELF),
true);
1116 static Teuchos::RCP<Epetra_SerialComm> communicator =
1117 Teuchos::rcp(
new Epetra_SerialComm(),
true);
1120 return *communicator;
1128# ifdef DEAL_II_WITH_MPI
1133 const Epetra_MpiComm *mpi_comm =
1134 dynamic_cast<const Epetra_MpiComm *
>(&communicator);
1135 if (mpi_comm !=
nullptr)
1136 return new Epetra_MpiComm(
1144 Assert(
dynamic_cast<const Epetra_SerialComm *
>(&communicator) !=
nullptr,
1146 return new Epetra_SerialComm(
1147 dynamic_cast<const Epetra_SerialComm &
>(communicator));
1157# ifdef DEAL_II_WITH_MPI
1158 Epetra_MpiComm *mpi_comm =
dynamic_cast<Epetra_MpiComm *
>(&communicator);
1159 if (mpi_comm !=
nullptr)
1162 *mpi_comm = Epetra_MpiComm(MPI_COMM_SELF);
1174 return mpi_communicator.NumProc();
1181 return static_cast<unsigned int>(mpi_communicator.MyPID());
1189 if (map.LinearMap() ==
true)
1196 return Epetra_Map(map.NumGlobalElements(),
1197 map.NumMyElements(),
1205 return Epetra_Map(map.NumGlobalElements(),
1206 map.NumMyElements(),
1207 map.MyGlobalElements(),
1217 template std::string
1218 to_string<int>(
int,
unsigned int);
1219 template std::string
1220 to_string<long int>(
long int,
unsigned int);
1221 template std::string
1222 to_string<long long int>(
long long int,
unsigned int);
1223 template std::string
1224 to_string<unsigned int>(
unsigned int,
unsigned int);
1225 template std::string
1226 to_string<unsigned long int>(
unsigned long int,
unsigned int);
1227 template std::string
1228 to_string<unsigned long long int>(
unsigned long long int,
unsigned int);
1229 template std::string
1230 to_string<float>(
float,
unsigned int);
1231 template std::string
1232 to_string<double>(
double,
unsigned int);
1233 template std::string
1234 to_string<long double>(
long double,
unsigned int);
1241 template std::vector<std::array<std::uint64_t, 1>>
1242 inverse_Hilbert_space_filling_curve<1, double>(
1245 template std::vector<std::array<std::uint64_t, 1>>
1246 inverse_Hilbert_space_filling_curve<1>(
1247 const std::vector<std::array<std::uint64_t, 1>> &,
1249 template std::vector<std::array<std::uint64_t, 2>>
1250 inverse_Hilbert_space_filling_curve<2, double>(
1253 template std::vector<std::array<std::uint64_t, 2>>
1254 inverse_Hilbert_space_filling_curve<2>(
1255 const std::vector<std::array<std::uint64_t, 2>> &,
1257 template std::vector<std::array<std::uint64_t, 3>>
1258 inverse_Hilbert_space_filling_curve<3, double>(
1261 template std::vector<std::array<std::uint64_t, 3>>
1262 inverse_Hilbert_space_filling_curve<3>(
1263 const std::vector<std::array<std::uint64_t, 3>> &,
1266 template std::uint64_t
1267 pack_integers<1>(
const std::array<std::uint64_t, 1> &,
const int);
1268 template std::uint64_t
1269 pack_integers<2>(
const std::array<std::uint64_t, 2> &,
const int);
1270 template std::uint64_t
1271 pack_integers<3>(
const std::array<std::uint64_t, 3> &,
const int);
A class that provides a separate storage location on each thread that accesses the object.
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_PACKAGE_VERSION
#define DEAL_II_VECTORIZATION_WIDTH_IN_BITS
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_PACKAGE_NAME
static ::ExceptionBase & ExcInvalidNumber2StringConversersion(unsigned int arg1, unsigned int arg2)
static ::ExceptionBase & ExcOutOfMemory(std::size_t arg1)
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcInvalidNumber(unsigned int arg1)
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define DeclException2(Exception2, type1, type2, outsequence)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcCantConvertString(std::string arg1)
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
void free_communicator(MPI_Comm &mpi_communicator)
MPI_Comm duplicate_communicator(const MPI_Comm &mpi_communicator)
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()
std::vector< std::string > split_string_list(const std::string &s, const std::string &delimiter=",")
Number truncate_to_n_digits(const Number number, const unsigned int n_digits)
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)
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::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)
std::string decompress(const std::string &compressed_input)
unsigned int needed_digits(const unsigned int max_number)
double string_to_double(const std::string &s)
std::string dealii_version_string()
std::string trim(const std::string &input)
double generate_normal_random_number(const double a, const double sigma)
int string_to_int(const std::string &s)
static const unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)