16 #include <deal.II/base/config.h> 22 # include <winsock2.h> 25 #include <deal.II/base/exceptions.h> 26 #include <deal.II/base/mpi.h> 27 #include <deal.II/base/point.h> 28 #include <deal.II/base/thread_local_storage.h> 29 #include <deal.II/base/utilities.h> 31 #include <boost/lexical_cast.hpp> 32 #include <boost/random.hpp> 48 #if defined(DEAL_II_HAVE_UNISTD_H) && defined(DEAL_II_HAVE_GETHOSTNAME) 57 #ifdef DEAL_II_WITH_TRILINOS 58 # ifdef DEAL_II_WITH_MPI 59 # include <deal.II/lac/trilinos_parallel_block_vector.h> 60 # include <deal.II/lac/trilinos_vector.h> 61 # include <deal.II/lac/vector_memory.h> 63 # include <Epetra_MpiComm.h> 64 # include <Teuchos_DefaultComm.hpp> 66 # include <Epetra_SerialComm.h> 67 # include <Teuchos_RCP.hpp> 70 DEAL_II_NAMESPACE_OPEN
78 <<
"When trying to convert " << arg1 <<
" to a string with " 79 << arg2 <<
" digits");
83 <<
"Can't convert the string " << arg1
84 <<
" to the desired type");
90 return DEAL_II_PACKAGE_NAME
" version " DEAL_II_PACKAGE_VERSION;
101 std::vector<std::array<std::uint64_t, effective_dim>>
102 inverse_Hilbert_space_filling_curve_effective(
105 const std::array<LongDouble, dim> & extents,
106 const std::bitset<dim> & valid_extents,
108 const Integer max_int)
110 std::vector<std::array<Integer, effective_dim>> int_points(points.size());
112 for (
unsigned int i = 0; i < points.size(); ++i)
115 unsigned int eff_d = 0;
116 for (
unsigned int d = 0; d < dim; ++d)
117 if (valid_extents[d])
120 const LongDouble v = (
static_cast<LongDouble
>(points[i][d]) -
121 static_cast<LongDouble>(bl[d])) /
125 int_points[i][eff_d] =
126 static_cast<Integer
>(v *
static_cast<LongDouble
>(max_int));
132 return inverse_Hilbert_space_filling_curve<effective_dim>(int_points,
137 template <
int dim,
typename Number>
138 std::vector<std::array<std::uint64_t, dim>>
141 const int bits_per_dim)
143 using Integer = std::uint64_t;
145 using LongDouble =
long double;
148 if (points.size() == 0)
149 return std::vector<std::array<std::uint64_t, dim>>();
153 for (
const auto &p : points)
154 for (
unsigned int d = 0; d < dim; ++d)
156 const double cid = p[d];
157 bl[d] = std::min(cid, bl[d]);
158 tr[d] = std::max(cid, tr[d]);
161 std::array<LongDouble, dim> extents;
162 std::bitset<dim> valid_extents;
163 for (
unsigned int i = 0; i < dim; ++i)
166 static_cast<LongDouble
>(tr[i]) - static_cast<LongDouble>(bl[i]);
167 valid_extents[i] = (extents[i] > 0.);
173 std::min(bits_per_dim,
174 std::min(std::numeric_limits<Integer>::digits,
175 std::numeric_limits<LongDouble>::digits));
178 const Integer max_int = (min_bits == std::numeric_limits<Integer>::digits ?
179 std::numeric_limits<Integer>::max() :
180 (Integer(1) << min_bits) - 1);
182 const unsigned int effective_dim = valid_extents.count();
183 if (effective_dim == dim)
185 return inverse_Hilbert_space_filling_curve_effective<dim,
190 points, bl, extents, valid_extents, min_bits, max_int);
194 std::array<std::uint64_t, dim> zero_ind;
195 for (
unsigned int d = 0; d < dim; ++d)
198 std::vector<std::array<std::uint64_t, dim>> ind(points.size(), zero_ind);
200 if (dim == 3 && effective_dim == 2)
203 inverse_Hilbert_space_filling_curve_effective<dim,
208 points, bl, extents, valid_extents, min_bits, max_int);
210 for (
unsigned int i = 0; i < ind.size(); ++i)
211 for (
unsigned int d = 0; d < 2; ++d)
212 ind[i][d + 1] = ind2[i][d];
216 else if (effective_dim == 1)
219 inverse_Hilbert_space_filling_curve_effective<dim,
224 points, bl, extents, valid_extents, min_bits, max_int);
226 for (
unsigned int i = 0; i < ind.size(); ++i)
227 ind[i][dim - 1] = ind1[i][0];
238 for (
unsigned int i = 0; i < points.size(); ++i)
247 std::vector<std::array<std::uint64_t, dim>>
249 const std::vector<std::array<std::uint64_t, dim>> &points,
250 const int bits_per_dim)
252 using Integer = std::uint64_t;
254 std::vector<std::array<Integer, dim>> int_points(points);
256 std::vector<std::array<Integer, dim>> res(int_points.size());
279 Assert(bits_per_dim <= std::numeric_limits<Integer>::digits,
280 ExcMessage(
"This integer type can not hold " +
281 std::to_string(bits_per_dim) +
" bits."));
283 const Integer M = Integer(1) << (bits_per_dim - 1);
285 for (
unsigned int index = 0; index < int_points.size(); ++index)
287 auto &X = int_points[index];
288 auto &L = res[index];
291 for (Integer q = M; q > 1; q >>= 1)
293 const Integer p = q - 1;
294 for (
unsigned int i = 0; i < dim; i++)
304 const Integer t = (X[0] ^ X[i]) & p;
312 for (
unsigned int i = 1; i < dim; i++)
316 for (Integer q = M; q > 1; q >>= 1)
319 for (
unsigned int i = 0; i < dim; i++)
334 for (
unsigned int i = 0; i < dim; ++i)
338 for (Integer q = M; q > 0; q >>= 1)
361 const int bits_per_dim)
363 using Integer = std::uint64_t;
368 const Integer mask = (Integer(1) << bits_per_dim) - 1;
371 for (
unsigned int i = 0; i < dim; ++i)
374 const Integer v = (mask & index[dim - 1 - i]) << (bits_per_dim * i);
390 template <
typename number>
392 to_string(
const number value,
const unsigned int digits)
394 std::string lc_string = boost::lexical_cast<std::string>(value);
398 else if (lc_string.size() < digits)
401 const unsigned int padding_position = (lc_string[0] ==
'-') ? 1 : 0;
403 const std::string padding(digits - lc_string.size(),
'0');
404 lc_string.insert(padding_position, padding);
412 const std::string &from,
413 const std::string &to)
418 std::string out = input;
419 std::string::size_type pos = out.find(from);
421 while (pos != std::string::npos)
423 out.replace(pos, from.size(), to);
424 pos = out.find(from, pos + to.size());
432 std::string::size_type left = 0;
433 std::string::size_type right = input.size() > 0 ? input.size() - 1 : 0;
435 for (; left < input.size(); ++left)
437 if (!std::isspace(input[left]))
443 for (; right >= left; --right)
445 if (!std::isspace(input[right]))
451 return std::string(input, left, right - left + 1);
471 if (max_number < 100)
473 if (max_number < 1000)
475 if (max_number < 10000)
477 if (max_number < 100000)
479 if (max_number < 1000000)
492 while ((s.size() > 0) && (s[0] ==
' '))
494 while ((s.size() > 0) && (s[s.size() - 1] ==
' '))
495 s.erase(s.end() - 1);
505 const int i = std::strtol(s.c_str(), &p, 10);
507 ((s.size() > 0) && (*p !=
'\0'))),
508 ExcMessage(
"Can't convert <" + s +
"> to an integer."));
518 std::vector<int> tmp(s.size());
519 for (
unsigned int i = 0; i < s.size(); ++i)
531 while ((s.size() > 0) && (s[0] ==
' '))
533 while ((s.size() > 0) && (s[s.size() - 1] ==
' '))
534 s.erase(s.end() - 1);
544 const double d = std::strtod(s.c_str(), &p);
546 ((s.size() > 0) && (*p !=
'\0'))),
547 ExcMessage(
"Can't convert <" + s +
"> to a double."));
557 std::vector<double> tmp(s.size());
558 for (
unsigned int i = 0; i < s.size(); ++i)
565 std::vector<std::string>
574 while (tmp.length() != 0 && tmp[tmp.length() - 1] ==
' ')
575 tmp.erase(tmp.length() - 1, 1);
583 std::vector<std::string> split_list;
584 while (tmp.length() != 0)
589 if (name.find(delimiter) != std::string::npos)
591 name.erase(name.find(delimiter), std::string::npos);
592 tmp.erase(0, tmp.find(delimiter) + delimiter.size());
598 while ((name.length() != 0) && (name[0] ==
' '))
600 while (name.length() != 0 && name[name.length() - 1] ==
' ')
601 name.erase(name.length() - 1, 1);
603 split_list.push_back(name);
610 std::vector<std::string>
619 std::vector<std::string>
621 const unsigned int width,
622 const char delimiter)
624 std::string text = original_text;
625 std::vector<std::string> lines;
628 while ((text.length() != 0) && (text[text.length() - 1] == delimiter))
629 text.erase(text.length() - 1, 1);
632 while (text.length() != 0)
636 while ((text.length() != 0) && (text[0] == delimiter))
639 std::size_t pos_newline = text.find_first_of(
'\n', 0);
640 if (pos_newline != std::string::npos && pos_newline <= width)
642 std::string line(text, 0, pos_newline);
643 while ((line.length() != 0) &&
644 (line[line.length() - 1] == delimiter))
645 line.erase(line.length() - 1, 1);
646 lines.push_back(line);
647 text.erase(0, pos_newline + 1);
654 if (text.length() < width)
657 while ((text.length() != 0) &&
658 (text[text.length() - 1] == delimiter))
659 text.erase(text.length() - 1, 1);
660 lines.push_back(text);
668 int location = std::min<int>(width, text.length() - 1);
669 for (; location > 0; --location)
670 if (text[location] == delimiter)
676 for (location = std::min<int>(width, text.length() - 1);
677 location < static_cast<int>(text.length());
679 if (text[location] == delimiter)
685 std::string line(text, 0, location);
686 while ((line.length() != 0) &&
687 (line[line.length() - 1] == delimiter))
688 line.erase(line.length() - 1, 1);
689 lines.push_back(line);
690 text.erase(0, location);
702 if (pattern.size() > name.size())
705 for (
unsigned int i = 0; i < pattern.size(); ++i)
706 if (pattern[i] != name[i])
714 std::pair<int, unsigned int>
719 const std::string test_string(name.begin() + position, name.end());
721 std::istringstream str(test_string);
731 return std::make_pair(i, 1U);
733 return std::make_pair(i, 2U);
735 return std::make_pair(i, 3U);
737 return std::make_pair(i, 4U);
739 return std::make_pair(i, 5U);
740 else if (i < 1000000)
741 return std::make_pair(i, 6U);
742 else if (i < 10000000)
743 return std::make_pair(i, 7U);
769 return boost::normal_distribution<>(a,
770 sigma)(random_number_generator.
get());
775 std::vector<unsigned int>
778 const unsigned int n = permutation.size();
780 std::vector<unsigned int> out(n);
781 for (
unsigned int i = 0; i < n; ++i)
782 out[i] = n - 1 - permutation[i];
789 std::vector<unsigned int>
792 const unsigned int n = permutation.size();
796 for (
unsigned int i = 0; i < n; ++i)
799 out[permutation[i]] = i;
804 for (
unsigned int i = 0; i < n; ++i)
806 ExcMessage(
"The given input permutation had duplicate entries!"));
811 std::vector<unsigned long long int>
814 const unsigned long long int n = permutation.size();
816 std::vector<unsigned long long int> out(n);
817 for (
unsigned long long int i = 0; i < n; ++i)
818 out[i] = n - 1 - permutation[i];
825 std::vector<unsigned long long int>
828 const unsigned long long int n = permutation.size();
832 for (
unsigned long long int i = 0; i < n; ++i)
835 out[permutation[i]] = i;
840 for (
unsigned long long int i = 0; i < n; ++i)
842 ExcMessage(
"The given input permutation had duplicate entries!"));
848 template <
typename Integer>
852 const unsigned int n = permutation.size();
854 std::vector<Integer> out(n);
855 for (
unsigned int i = 0; i < n; ++i)
856 out[i] = n - 1 - permutation[i];
863 template <
typename Integer>
867 const unsigned int n = permutation.size();
871 for (
unsigned int i = 0; i < n; ++i)
874 out[permutation[i]] = i;
879 for (
unsigned int i = 0; i < n; ++i)
881 ExcMessage(
"The given input permutation had duplicate entries!"));
890 #if defined(__linux__) 895 std::ifstream cpuinfo;
896 cpuinfo.open(
"/proc/loadavg");
919 switch (DEAL_II_COMPILER_VECTORIZATION_LEVEL)
936 "Invalid DEAL_II_COMPILER_VECTORIZATION_LEVEL."));
950 #if defined(__linux__) 951 std::ifstream file(
"/proc/self/status");
957 if (name ==
"VmPeak:")
959 else if (name ==
"VmSize:")
961 else if (name ==
"VmHWM:")
963 else if (name ==
"VmRSS:")
979 #if defined(DEAL_II_HAVE_UNISTD_H) && defined(DEAL_II_HAVE_GETHOSTNAME) 980 const unsigned int N = 1024;
982 gethostname(&(hostname[0]), N - 1);
984 std::string hostname(
"unknown");
994 std::time_t time1 = std::time(
nullptr);
995 std::tm * time = std::localtime(&time1);
997 std::ostringstream o;
998 o << time->tm_hour <<
":" << (time->tm_min < 10 ?
"0" :
"")
999 << time->tm_min <<
":" << (time->tm_sec < 10 ?
"0" :
"")
1010 std::time_t time1 = std::time(
nullptr);
1011 std::tm * time = std::localtime(&time1);
1013 std::ostringstream o;
1014 o << time->tm_year + 1900 <<
"/" << time->tm_mon + 1 <<
"/" 1025 #ifndef DEAL_II_MSVC 1033 *memptr = malloc(size);
1049 #ifdef DEAL_II_WITH_TRILINOS 1056 # ifdef DEAL_II_WITH_MPI 1057 static Teuchos::RCP<Epetra_MpiComm> communicator =
1058 Teuchos::rcp(
new Epetra_MpiComm(MPI_COMM_WORLD),
true);
1060 static Teuchos::RCP<Epetra_SerialComm> communicator =
1061 Teuchos::rcp(
new Epetra_SerialComm(),
true);
1064 return *communicator;
1069 const Teuchos::RCP<const Teuchos::Comm<int>> &
1072 # ifdef DEAL_II_WITH_MPI 1073 static auto communicator = Teuchos::RCP<const Teuchos::Comm<int>>(
1074 new Teuchos::MpiComm<int>(MPI_COMM_SELF));
1076 static auto communicator =
1077 Teuchos::RCP<const Teuchos::Comm<int>>(
new Teuchos::Comm<int>());
1080 return communicator;
1088 # ifdef DEAL_II_WITH_MPI 1089 static Teuchos::RCP<Epetra_MpiComm> communicator =
1090 Teuchos::rcp(
new Epetra_MpiComm(MPI_COMM_SELF),
true);
1092 static Teuchos::RCP<Epetra_SerialComm> communicator =
1093 Teuchos::rcp(
new Epetra_SerialComm(),
true);
1096 return *communicator;
1104 # ifdef DEAL_II_WITH_MPI 1109 const Epetra_MpiComm *mpi_comm =
1110 dynamic_cast<const Epetra_MpiComm *
>(&communicator);
1111 if (mpi_comm !=
nullptr)
1112 return new Epetra_MpiComm(
1120 Assert(dynamic_cast<const Epetra_SerialComm *>(&communicator) !=
nullptr,
1122 return new Epetra_SerialComm(
1123 dynamic_cast<const Epetra_SerialComm &>(communicator));
1133 # ifdef DEAL_II_WITH_MPI 1134 Epetra_MpiComm *mpi_comm =
dynamic_cast<Epetra_MpiComm *
>(&communicator);
1135 if (mpi_comm !=
nullptr)
1137 MPI_Comm comm = mpi_comm->GetMpiComm();
1138 *mpi_comm = Epetra_MpiComm(MPI_COMM_SELF);
1139 const int ierr = MPI_Comm_free(&comm);
1150 return mpi_communicator.NumProc();
1157 return static_cast<unsigned int>(mpi_communicator.MyPID());
1165 if (map.LinearMap() ==
true)
1172 return Epetra_Map(map.NumGlobalElements(),
1173 map.NumMyElements(),
1181 return Epetra_Map(map.NumGlobalElements(),
1182 map.NumMyElements(),
1183 map.MyGlobalElements(),
1192 template std::string
1193 to_string<int>(int,
unsigned int);
1194 template std::string
1195 to_string<long int>(
long int,
unsigned int);
1196 template std::string
1197 to_string<long long int>(
long long int,
unsigned int);
1198 template std::string
1199 to_string<unsigned int>(
unsigned int,
unsigned int);
1200 template std::string
1201 to_string<unsigned long int>(
unsigned long int,
unsigned int);
1202 template std::string
1203 to_string<unsigned long long int>(
unsigned long long int,
unsigned int);
1204 template std::string
1205 to_string<float>(float,
unsigned int);
1206 template std::string
1207 to_string<double>(double,
unsigned int);
1208 template std::string
1209 to_string<long double>(
long double,
unsigned int);
1211 template std::vector<std::array<std::uint64_t, 1>>
1212 inverse_Hilbert_space_filling_curve<1, double>(
1213 const std::vector<Point<1, double>> &,
1215 template std::vector<std::array<std::uint64_t, 1>>
1216 inverse_Hilbert_space_filling_curve<1>(
1217 const std::vector<std::array<std::uint64_t, 1>> &,
1219 template std::vector<std::array<std::uint64_t, 2>>
1220 inverse_Hilbert_space_filling_curve<2, double>(
1221 const std::vector<Point<2, double>> &,
1223 template std::vector<std::array<std::uint64_t, 2>>
1224 inverse_Hilbert_space_filling_curve<2>(
1225 const std::vector<std::array<std::uint64_t, 2>> &,
1227 template std::vector<std::array<std::uint64_t, 3>>
1228 inverse_Hilbert_space_filling_curve<3, double>(
1229 const std::vector<Point<3, double>> &,
1231 template std::vector<std::array<std::uint64_t, 3>>
1232 inverse_Hilbert_space_filling_curve<3>(
1233 const std::vector<std::array<std::uint64_t, 3>> &,
1236 template std::uint64_t
1237 pack_integers<1>(
const std::array<std::uint64_t, 1> &,
const int);
1238 template std::uint64_t
1239 pack_integers<2>(
const std::array<std::uint64_t, 2> &,
const int);
1240 template std::uint64_t
1241 pack_integers<3>(
const std::array<std::uint64_t, 3> &,
const int);
1244 DEAL_II_NAMESPACE_CLOSE
std::vector< std::string > split_string_list(const std::string &s, const std::string &delimiter=",")
void posix_memalign(void **memptr, std::size_t alignment, std::size_t size)
static const unsigned int invalid_unsigned_int
std::uint64_t pack_integers(const std::array< std::uint64_t, dim > &index, const int bits_per_dim)
#define DeclException2(Exception2, type1, type2, outsequence)
A class that provides a separate storage location on each thread that accesses the object...
static ::ExceptionBase & ExcIO()
std::string dealii_version_string()
unsigned int get_n_mpi_processes(const Epetra_Comm &mpi_communicator)
std::string trim(const std::string &input)
#define AssertIndexRange(index, range)
std::pair< int, unsigned int > get_integer_at_position(const std::string &name, const unsigned int position)
std::vector< std::string > break_text_into_lines(const std::string &original_text, const unsigned int width, const char delimiter=' ')
#define AssertThrow(cond, exc)
static ::ExceptionBase & ExcOutOfMemory()
const Epetra_Comm & comm_self()
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
void get_memory_stats(MemoryStats &stats)
double string_to_double(const std::string &s)
const Teuchos::RCP< const Teuchos::Comm< int > > & tpetra_comm_self()
static ::ExceptionBase & ExcMessage(std::string arg1)
Epetra_Comm * duplicate_communicator(const Epetra_Comm &communicator)
#define DeclException1(Exception1, type1, outsequence)
double generate_normal_random_number(const double a, const double sigma)
#define Assert(cond, exc)
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)
std::vector< unsigned int > invert_permutation(const std::vector< unsigned int > &permutation)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
std::string replace_in_string(const std::string &input, const std::string &from, const std::string &to)
std::string dim_string(const int dim, const int spacedim)
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)
void destroy_communicator(Epetra_Comm &communicator)
std::string get_hostname()
#define AssertThrowMPI(error_code)
MPI_Comm duplicate_communicator(const MPI_Comm &mpi_communicator)
const Epetra_Comm & comm_world()
int string_to_int(const std::string &s)
std::vector< unsigned int > reverse_permutation(const std::vector< unsigned int > &permutation)
static ::ExceptionBase & ExcNotImplemented()
unsigned int get_this_mpi_process(const Epetra_Comm &mpi_communicator)
const std::string get_current_vectorization_level()
Epetra_Map duplicate_map(const Epetra_BlockMap &map, const Epetra_Comm &comm)
unsigned int needed_digits(const unsigned int max_number)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcInvalidNumber2StringConversersion(unsigned int arg1, unsigned int arg2)