16 #include <deal.II/base/config.h> 22 # include <winsock2.h> 25 #include <deal.II/base/utilities.h> 26 #include <deal.II/base/mpi.h> 27 #include <deal.II/base/exceptions.h> 28 #include <deal.II/base/thread_local_storage.h> 30 #include <boost/math/special_functions/erf.hpp> 31 #include <boost/lexical_cast.hpp> 32 #include <boost/random.hpp> 47 #if defined(DEAL_II_HAVE_UNISTD_H) && defined(DEAL_II_HAVE_GETHOSTNAME) 56 #ifdef DEAL_II_WITH_TRILINOS 57 # ifdef DEAL_II_WITH_MPI 58 # include <Epetra_MpiComm.h> 59 # include <deal.II/lac/vector_memory.h> 60 # include <deal.II/lac/trilinos_vector.h> 61 # include <deal.II/lac/trilinos_parallel_block_vector.h> 63 # include <Teuchos_RCP.hpp> 64 # include <Epetra_SerialComm.h> 69 DEAL_II_NAMESPACE_OPEN
77 unsigned int,
unsigned int,
78 <<
"When trying to convert " << arg1
79 <<
" to a string with " << arg2 <<
" digits");
82 <<
"Invalid number " << arg1);
85 <<
"Can't convert the string " << arg1
86 <<
" to the desired type");
92 return DEAL_II_PACKAGE_NAME
" version " DEAL_II_PACKAGE_VERSION;
105 template <
typename number>
107 to_string (
const number value,
const unsigned int digits)
109 std::string lc_string = boost::lexical_cast<std::string>(value);
113 else if (lc_string.size() < digits)
116 const unsigned int padding_position = (lc_string[0] ==
'-')
122 const std::string padding(digits - lc_string.size(),
'0');
123 lc_string.insert(padding_position, padding);
131 const std::string &from,
132 const std::string &to)
137 std::string out = input;
138 std::string::size_type pos = out.find(from);
140 while (pos != std::string::npos)
142 out.replace(pos, from.size(), to);
143 pos = out.find(from, pos + to.size());
151 std::string::size_type left = 0;
152 std::string::size_type right = input.size() > 0 ? input.size() - 1 : 0;
154 for (; left < input.size(); ++left)
156 if (!std::isspace(input[left]))
162 for (; right >= left; --right)
164 if (!std::isspace(input[right]))
170 return std::string(input, left, right - left + 1);
190 if (max_number < 100)
192 if (max_number < 1000)
194 if (max_number < 10000)
196 if (max_number < 100000)
198 if (max_number < 1000000)
211 while ((s.size() > 0) && (s[0] ==
' '))
213 while ((s.size() > 0) && (s[s.size()-1] ==
' '))
224 const int i = std::strtol(s.c_str(), &p, 10);
225 AssertThrow ( !((errno != 0) || (s.size() == 0) || ((s.size()>0) && (*p !=
'\0'))),
226 ExcMessage (
"Can't convert <" + s +
"> to an integer."));
236 std::vector<int> tmp (s.size());
237 for (
unsigned int i=0; i<s.size(); ++i)
249 while ((s.size() > 0) && (s[0] ==
' '))
251 while ((s.size() > 0) && (s[s.size()-1] ==
' '))
262 const double d = std::strtod(s.c_str(), &p);
263 AssertThrow ( !((errno != 0) || (s.size() == 0) || ((s.size()>0) && (*p !=
'\0'))),
264 ExcMessage (
"Can't convert <" + s +
"> to a double."));
274 std::vector<double> tmp (s.size());
275 for (
unsigned int i=0; i<s.size(); ++i)
282 std::vector<std::string>
284 const std::string &delimiter)
292 while (tmp.length() != 0 && tmp[tmp.length()-1] ==
' ')
293 tmp.erase (tmp.length()-1, 1);
301 std::vector<std::string> split_list;
302 while (tmp.length() != 0)
307 if (name.find(delimiter) != std::string::npos)
309 name.erase (name.find(delimiter), std::string::npos);
310 tmp.erase (0, tmp.find(delimiter)+delimiter.size());
316 while ((name.length() != 0) && (name[0] ==
' '))
318 while (name.length() != 0 && name[name.length()-1] ==
' ')
319 name.erase (name.length()-1, 1);
321 split_list.push_back (name);
328 std::vector<std::string>
330 const char delimiter)
338 std::vector<std::string>
340 const unsigned int width,
341 const char delimiter)
343 std::string text = original_text;
344 std::vector<std::string> lines;
347 while ((text.length() != 0) && (text[text.length()-1] == delimiter))
348 text.erase(text.length()-1,1);
351 while (text.length() != 0)
355 while ((text.length() != 0) && (text[0] == delimiter))
358 std::size_t pos_newline = text.find_first_of(
'\n', 0);
359 if (pos_newline != std::string::npos && pos_newline <= width)
361 std::string line (text, 0, pos_newline);
362 while ((line.length() != 0) && (line[line.length()-1] == delimiter))
363 line.erase(line.length()-1,1);
364 lines.push_back (line);
365 text.erase (0, pos_newline+1);
372 if (text.length() < width)
375 while ((text.length() != 0) && (text[text.length()-1] == delimiter))
376 text.erase(text.length()-1,1);
377 lines.push_back (text);
385 int location = std::min<int>(width,text.length()-1);
386 for (; location>0; --location)
387 if (text[location] == delimiter)
393 for (location = std::min<int>(width,text.length()-1);
394 location<static_cast<int>(text.length());
396 if (text[location] == delimiter)
402 std::string line (text, 0, location);
403 while ((line.length() != 0) && (line[line.length()-1] == delimiter))
404 line.erase(line.length()-1,1);
405 lines.push_back (line);
406 text.erase (0, location);
417 const std::string &pattern)
419 if (pattern.size() > name.size())
422 for (
unsigned int i=0; i<pattern.size(); ++i)
423 if (pattern[i] != name[i])
431 std::pair<int, unsigned int>
433 const unsigned int position)
437 const std::string test_string (name.begin()+position,
440 std::istringstream str(test_string);
450 return std::make_pair (i, 1U);
452 return std::make_pair (i, 2U);
454 return std::make_pair (i, 3U);
456 return std::make_pair (i, 4U);
458 return std::make_pair (i, 5U);
460 return std::make_pair (i, 6U);
462 return std::make_pair (i, 7U);
489 return boost::normal_distribution<>(a,sigma)(random_number_generator.
get());
494 std::vector<unsigned int>
497 const unsigned int n = permutation.size();
499 std::vector<unsigned int> out (n);
500 for (
unsigned int i=0; i<n; ++i)
501 out[i] = n - 1 - permutation[i];
508 std::vector<unsigned int>
511 const unsigned int n = permutation.size();
515 for (
unsigned int i=0; i<n; ++i)
518 out[permutation[i]] = i;
523 for (
unsigned int i=0; i<n; ++i)
525 ExcMessage (
"The given input permutation had duplicate entries!"));
530 std::vector<unsigned long long int>
533 const unsigned long long int n = permutation.size();
535 std::vector<unsigned long long int> out (n);
536 for (
unsigned long long int i=0; i<n; ++i)
537 out[i] = n - 1 - permutation[i];
544 std::vector<unsigned long long int>
547 const unsigned long long int n = permutation.size();
551 for (
unsigned long long int i=0; i<n; ++i)
554 out[permutation[i]] = i;
559 for (
unsigned long long int i=0; i<n; ++i)
561 ExcMessage (
"The given input permutation had duplicate entries!"));
567 template <
typename Integer>
571 const unsigned int n = permutation.size();
573 std::vector<Integer> out (n);
574 for (
unsigned int i=0; i<n; ++i)
575 out[i] = n - 1 - permutation[i];
582 template <
typename Integer>
586 const unsigned int n = permutation.size();
590 for (
unsigned int i=0; i<n; ++i)
593 out[permutation[i]] = i;
598 for (
unsigned int i=0; i<n; ++i)
600 ExcMessage (
"The given input permutation had duplicate entries!"));
610 #if defined(__linux__) 614 std::ifstream cpuinfo;
615 cpuinfo.open(
"/proc/loadavg");
637 switch (DEAL_II_COMPILER_VECTORIZATION_LEVEL)
661 #if defined(__linux__) 662 std::ifstream file(
"/proc/self/status");
668 if (name ==
"VmPeak:")
670 else if (name ==
"VmSize:")
672 else if (name ==
"VmHWM:")
674 else if (name ==
"VmRSS:")
689 #if defined(DEAL_II_HAVE_UNISTD_H) && defined(DEAL_II_HAVE_GETHOSTNAME) 690 const unsigned int N=1024;
692 gethostname (&(hostname[0]), N-1);
694 std::string hostname(
"unknown");
703 std::time_t time1= std::time (
nullptr);
704 std::tm *time = std::localtime(&time1);
706 std::ostringstream o;
707 o << time->tm_hour <<
":" 708 << (time->tm_min < 10 ?
"0" :
"") << time->tm_min <<
":" 709 << (time->tm_sec < 10 ?
"0" :
"") << time->tm_sec;
718 std::time_t time1= std::time (
nullptr);
719 std::tm *time = std::localtime(&time1);
721 std::ostringstream o;
722 o << time->tm_year + 1900 <<
"/" 723 << time->tm_mon + 1 <<
"/" 741 *memptr = malloc (size);
749 bool job_supports_mpi ()
756 #ifdef DEAL_II_WITH_TRILINOS 763 #ifdef DEAL_II_WITH_MPI 764 static Teuchos::RCP<Epetra_MpiComm>
765 communicator = Teuchos::rcp (
new Epetra_MpiComm (MPI_COMM_WORLD),
true);
767 static Teuchos::RCP<Epetra_SerialComm>
768 communicator = Teuchos::rcp (
new Epetra_SerialComm (),
true);
771 return *communicator;
779 #ifdef DEAL_II_WITH_MPI 780 static Teuchos::RCP<Epetra_MpiComm>
781 communicator = Teuchos::rcp (
new Epetra_MpiComm (MPI_COMM_SELF),
true);
783 static Teuchos::RCP<Epetra_SerialComm>
784 communicator = Teuchos::rcp (
new Epetra_SerialComm (),
true);
787 return *communicator;
795 #ifdef DEAL_II_WITH_MPI 801 *mpi_comm =
dynamic_cast<const Epetra_MpiComm *
>(&communicator);
802 if (mpi_comm !=
nullptr)
811 Assert (dynamic_cast<const Epetra_SerialComm *>(&communicator)
814 return new Epetra_SerialComm(dynamic_cast<const Epetra_SerialComm &>(communicator));
823 #ifdef DEAL_II_WITH_MPI 825 *mpi_comm =
dynamic_cast<Epetra_MpiComm *
>(&communicator);
826 if (mpi_comm !=
nullptr)
828 MPI_Comm comm = mpi_comm->GetMpiComm();
829 *mpi_comm = Epetra_MpiComm(MPI_COMM_SELF);
830 const int ierr = MPI_Comm_free (&comm);
840 return mpi_communicator.NumProc();
846 return (
unsigned int)mpi_communicator.MyPID();
853 const Epetra_Comm &comm)
855 if (map.LinearMap() ==
true)
862 return Epetra_Map (map.NumGlobalElements(),
871 return Epetra_Map (map.NumGlobalElements(),
873 map.MyGlobalElements (),
882 template std::string to_string<int> (int,
unsigned int);
883 template std::string to_string<long int> (
long int,
unsigned int);
884 template std::string to_string<long long int> (
long long int,
unsigned int);
885 template std::string to_string<unsigned int> (
unsigned int,
unsigned int);
886 template std::string to_string<unsigned long int> (
unsigned long int,
unsigned int);
887 template std::string to_string<unsigned long long int> (
unsigned long long int,
unsigned int);
888 template std::string to_string<float> (float,
unsigned int);
889 template std::string to_string<double> (double,
unsigned int);
890 template std::string to_string<long double> (
long double,
unsigned int);
894 DEAL_II_NAMESPACE_CLOSE
std::vector< std::string > split_string_list(const std::string &s, const std::string &delimiter=",")
static const unsigned int invalid_unsigned_int
#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)
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)
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)
void posix_memalign(void **memptr, size_t alignment, size_t size)
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)
void destroy_communicator(Epetra_Comm &communicator)
std::string get_hostname()
#define AssertThrowMPI(error_code)
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)