22 # include <winsock2.h>
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 ?
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 " +
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(
402 bio::gzip_params(boost::iostreams::gzip::default_compression)));
406 return compressed.str();
417 #ifdef DEAL_II_WITH_ZLIB
418 namespace bio = boost::iostreams;
420 std::stringstream compressed(compressed_input);
421 std::stringstream decompressed;
423 bio::filtering_streambuf<bio::input> out;
424 out.push(bio::gzip_decompressor());
425 out.push(compressed);
428 return decompressed.str();
430 return compressed_input;
439 using namespace boost::archive::iterators;
440 using It = base64_from_binary<
441 transform_width<std::vector<unsigned char>::const_iterator, 6, 8>>;
442 auto base64 = std::string(It(binary_input.begin()), It(binary_input.end()));
444 return base64.append((3 - binary_input.size() % 3) % 3,
'=');
449 std::vector<unsigned char>
452 using namespace boost::archive::iterators;
454 transform_width<binary_from_base64<std::string::const_iterator>, 8, 6>;
455 auto binary = std::vector<unsigned char>(It(base64_input.begin()),
456 It(base64_input.end()));
458 auto length = base64_input.size();
459 if (binary.size() > 2 && base64_input[length - 1] ==
'=' &&
460 base64_input[length - 2] ==
'=')
462 binary.erase(binary.end() - 2, binary.end());
464 else if (binary.size() > 1 && base64_input[length - 1] ==
'=')
466 binary.erase(binary.end() - 1, binary.end());
481 template <
typename number>
495 boost::lexical_cast<std::string>(
value));
498 (lc_string.size() < digits))
501 const unsigned int padding_position = (lc_string[0] ==
'-') ? 1 : 0;
503 const std::string padding(digits - lc_string.size(),
'0');
504 lc_string.insert(padding_position, padding);
514 const std::string &from,
515 const std::string &to)
520 std::string out = input;
523 while (pos != std::string::npos)
525 out.replace(pos, from.size(), to);
526 pos = out.find(from, pos + to.size());
537 for (; left < input.size(); ++left)
539 if (!std::isspace(input[left]))
545 for (; right >= left; --right)
547 if (!std::isspace(input[right]))
553 return std::string(input, left, right - left + 1);
572 return static_cast<int>(
580 template <
typename Number>
592 const int shift = -order +
static_cast<int>(n_digits) - 1;
597 "Overflow. Use a smaller value for n_digits and/or make sure "
598 "that the absolute value of 'number' does not become too small."));
600 const Number factor = std::pow(10.0,
static_cast<Number
>(
shift));
602 const Number number_cutoff = std::trunc(number * factor) / factor;
604 return number_cutoff;
613 while ((s.size() > 0) && (s[0] ==
' '))
615 while ((s.size() > 0) && (s[s.size() - 1] ==
' '))
616 s.erase(s.end() - 1);
626 const int i = std::strtol(s.c_str(), &p, 10);
637 ((s.size() > 0) && (*p !=
'\0'))),
638 ExcMessage(
"Can't convert <" + s +
"> to an integer."));
648 std::vector<int> tmp(s.size());
649 for (
unsigned int i = 0; i < s.size(); ++i)
661 while ((s.size() > 0) && (s[0] ==
' '))
663 while ((s.size() > 0) && (s[s.size() - 1] ==
' '))
664 s.erase(s.end() - 1);
674 const double d = std::strtod(s.c_str(), &p);
685 ((s.size() > 0) && (*p !=
'\0'))),
686 ExcMessage(
"Can't convert <" + s +
"> to a double."));
696 std::vector<double> tmp(s.size());
697 for (
unsigned int i = 0; i < s.size(); ++i)
704 std::vector<std::string>
713 while (tmp.length() != 0 && tmp[tmp.length() - 1] ==
' ')
714 tmp.erase(tmp.length() - 1, 1);
722 std::vector<std::string> split_list;
723 while (tmp.length() != 0)
728 if (name.find(delimiter) != std::string::npos)
730 name.erase(name.find(delimiter), std::string::npos);
731 tmp.erase(0, tmp.find(delimiter) + delimiter.size());
737 while ((name.length() != 0) && (name[0] ==
' '))
739 while (name.length() != 0 && name[name.length() - 1] ==
' ')
740 name.erase(name.length() - 1, 1);
742 split_list.push_back(name);
749 std::vector<std::string>
758 std::vector<std::string>
760 const unsigned int width,
761 const char delimiter)
763 std::string text = original_text;
764 std::vector<std::string> lines;
767 while ((text.length() != 0) && (text[text.length() - 1] == delimiter))
768 text.erase(text.length() - 1, 1);
771 while (text.length() != 0)
775 while ((text.length() != 0) && (text[0] == delimiter))
778 std::size_t pos_newline = text.find_first_of(
'\n', 0);
779 if (pos_newline != std::string::npos && pos_newline <= width)
781 std::string line(text, 0, pos_newline);
782 while ((line.length() != 0) &&
783 (line[line.length() - 1] == delimiter))
784 line.erase(line.length() - 1, 1);
785 lines.push_back(line);
786 text.erase(0, pos_newline + 1);
793 if (text.length() < width)
796 while ((text.length() != 0) &&
797 (text[text.length() - 1] == delimiter))
798 text.erase(text.length() - 1, 1);
799 lines.push_back(text);
807 int location = std::min<int>(width, text.length() - 1);
808 for (; location > 0; --location)
809 if (text[location] == delimiter)
815 for (location = std::min<int>(width, text.length() - 1);
816 location <
static_cast<int>(text.length());
818 if (text[location] == delimiter)
824 std::string line(text, 0, location);
825 while ((line.length() != 0) &&
826 (line[line.length() - 1] == delimiter))
827 line.erase(line.length() - 1, 1);
828 lines.push_back(line);
829 text.erase(0, location);
841 if (pattern.size() > name.size())
844 for (
unsigned int i = 0; i < pattern.size(); ++i)
845 if (pattern[i] != name[i])
853 std::pair<int, unsigned int>
858 const std::string test_string(name.begin() + position, name.end());
860 std::istringstream str(test_string);
870 return std::make_pair(i, 1
U);
872 return std::make_pair(i, 2
U);
874 return std::make_pair(i, 3
U);
876 return std::make_pair(i, 4
U);
878 return std::make_pair(i, 5
U);
879 else if (i < 1000000)
880 return std::make_pair(i, 6
U);
881 else if (i < 10000000)
882 return std::make_pair(i, 7
U);
910 return boost::normal_distribution<>(a,
911 sigma)(random_number_generator.
get());
918 #if defined(__linux__)
923 std::ifstream cpuinfo;
924 cpuinfo.open(
"/proc/loadavg");
964 "Invalid DEAL_II_VECTORIZATION_WIDTH_IN_BITS."));
978 #if defined(__linux__)
979 std::ifstream file(
"/proc/self/status");
985 if (name ==
"VmPeak:")
987 else if (name ==
"VmSize:")
989 else if (name ==
"VmHWM:")
991 else if (name ==
"VmRSS:")
1007 #if defined(DEAL_II_HAVE_UNISTD_H) && defined(DEAL_II_HAVE_GETHOSTNAME)
1008 const unsigned int N = 1024;
1010 gethostname(&(hostname[0]),
N - 1);
1012 std::string hostname(
"unknown");
1022 std::time_t time1 = std::time(
nullptr);
1023 std::tm * time = std::localtime(&time1);
1025 std::ostringstream o;
1026 o << time->tm_hour <<
":" << (time->tm_min < 10 ?
"0" :
"")
1027 << time->tm_min <<
":" << (time->tm_sec < 10 ?
"0" :
"")
1038 std::time_t time1 = std::time(
nullptr);
1039 std::tm * time = std::localtime(&time1);
1041 std::ostringstream o;
1042 o << time->tm_year + 1900 <<
"/" << time->tm_mon + 1 <<
"/"
1053 #ifndef DEAL_II_MSVC
1077 #ifdef DEAL_II_WITH_TRILINOS
1084 # ifdef DEAL_II_WITH_MPI
1085 static Teuchos::RCP<Epetra_MpiComm> communicator =
1086 Teuchos::rcp(
new Epetra_MpiComm(MPI_COMM_WORLD),
true);
1088 static Teuchos::RCP<Epetra_SerialComm> communicator =
1089 Teuchos::rcp(
new Epetra_SerialComm(),
true);
1092 return *communicator;
1097 const Teuchos::RCP<const Teuchos::Comm<int>> &
1100 # ifdef DEAL_II_WITH_MPI
1101 static auto communicator = Teuchos::RCP<const Teuchos::Comm<int>>(
1102 new Teuchos::MpiComm<int>(MPI_COMM_SELF));
1104 static auto communicator =
1105 Teuchos::RCP<const Teuchos::Comm<int>>(
new Teuchos::Comm<int>());
1108 return communicator;
1116 # ifdef DEAL_II_WITH_MPI
1117 static Teuchos::RCP<Epetra_MpiComm> communicator =
1118 Teuchos::rcp(
new Epetra_MpiComm(MPI_COMM_SELF),
true);
1120 static Teuchos::RCP<Epetra_SerialComm> communicator =
1121 Teuchos::rcp(
new Epetra_SerialComm(),
true);
1124 return *communicator;
1132 # ifdef DEAL_II_WITH_MPI
1137 const Epetra_MpiComm *mpi_comm =
1138 dynamic_cast<const Epetra_MpiComm *
>(&communicator);
1139 if (mpi_comm !=
nullptr)
1140 return new Epetra_MpiComm(
1148 Assert(
dynamic_cast<const Epetra_SerialComm *
>(&communicator) !=
nullptr,
1150 return new Epetra_SerialComm(
1151 dynamic_cast<const Epetra_SerialComm &
>(communicator));
1161 # ifdef DEAL_II_WITH_MPI
1162 Epetra_MpiComm *mpi_comm =
dynamic_cast<Epetra_MpiComm *
>(&communicator);
1163 if (mpi_comm !=
nullptr)
1165 MPI_Comm comm = mpi_comm->GetMpiComm();
1166 *mpi_comm = Epetra_MpiComm(MPI_COMM_SELF);
1167 const int ierr = MPI_Comm_free(&comm);
1178 return mpi_communicator.NumProc();
1185 return static_cast<unsigned int>(mpi_communicator.MyPID());
1193 if (map.LinearMap() ==
true)
1200 return Epetra_Map(map.NumGlobalElements(),
1201 map.NumMyElements(),
1209 return Epetra_Map(map.NumGlobalElements(),
1210 map.NumMyElements(),
1211 map.MyGlobalElements(),
1220 template std::string
1222 template std::string
1224 template std::string
1226 template std::string
1228 template std::string
1230 template std::string
1232 template std::string
1234 template std::string
1236 template std::string
1244 template std::vector<std::array<std::uint64_t, 1>>
1248 template std::vector<std::array<std::uint64_t, 1>>
1250 const std::vector<std::array<std::uint64_t, 1>> &,
1252 template std::vector<std::array<std::uint64_t, 2>>
1256 template std::vector<std::array<std::uint64_t, 2>>
1258 const std::vector<std::array<std::uint64_t, 2>> &,
1260 template std::vector<std::array<std::uint64_t, 3>>
1264 template std::vector<std::array<std::uint64_t, 3>>
1266 const std::vector<std::array<std::uint64_t, 3>> &,
1269 template std::uint64_t
1271 template std::uint64_t
1273 template std::uint64_t