17 #include <deal.II/base/utilities.h> 18 #include <deal.II/base/mpi.h> 19 #include <deal.II/base/exceptions.h> 20 #include <deal.II/base/thread_local_storage.h> 22 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
23 #include <boost/math/special_functions/erf.hpp> 24 #include <boost/lexical_cast.hpp> 25 #include <boost/random.hpp> 26 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
46 # include <winsock2.h> 50 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
51 #ifdef DEAL_II_WITH_TRILINOS 52 # ifdef DEAL_II_WITH_MPI 53 # include <Epetra_MpiComm.h> 54 # include <deal.II/lac/vector_memory.h> 55 # include <deal.II/lac/trilinos_vector.h> 56 # include <deal.II/lac/trilinos_block_vector.h> 58 # include "Teuchos_RCP.hpp" 59 # include "Epetra_SerialComm.h" 61 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
65 DEAL_II_NAMESPACE_OPEN
73 unsigned int,
unsigned int,
74 <<
"When trying to convert " << arg1
75 <<
" to a string with " << arg2 <<
" digits");
78 <<
"Invalid number " << arg1);
81 <<
"Can't convert the string " << arg1
82 <<
" to the desired type");
90 template <
typename number>
92 to_string (
const number value,
const unsigned int digits)
94 std::string lc_string = boost::lexical_cast<std::string>(value);
98 else if (lc_string.size() < digits)
101 const unsigned int padding_position = (lc_string[0] ==
'-')
107 const std::string padding(digits - lc_string.size(),
'0');
108 lc_string.insert(padding_position, padding);
116 const std::string &from,
117 const std::string &to)
122 std::string out = input;
123 std::string::size_type pos = out.find(from);
125 while (pos != std::string::npos)
127 out.replace(pos, from.size(), to);
128 pos = out.find(from, pos + to.size());
136 std::string::size_type left = 0;
137 std::string::size_type right = input.size() > 0 ? input.size() - 1 : 0;
139 for (; left < input.size(); ++left)
141 if (!std::isspace(input[left]))
147 for (; right >= left; --right)
149 if (!std::isspace(input[right]))
155 return std::string(input, left, right - left + 1);
175 if (max_number < 100)
177 if (max_number < 1000)
179 if (max_number < 10000)
181 if (max_number < 100000)
183 if (max_number < 1000000)
196 while ((s.size() > 0) && (s[0] ==
' '))
198 while ((s.size() > 0) && (s[s.size()-1] ==
' '))
209 const int i = std::strtol(s.c_str(), &p, 10);
210 AssertThrow ( !((errno != 0) || (s.size() == 0) || ((s.size()>0) && (*p !=
'\0'))),
211 ExcMessage (
"Can't convert <" + s +
"> to an integer."));
221 std::vector<int> tmp (s.size());
222 for (
unsigned int i=0; i<s.size(); ++i)
234 while ((s.size() > 0) && (s[0] ==
' '))
236 while ((s.size() > 0) && (s[s.size()-1] ==
' '))
247 const double d = std::strtod(s.c_str(), &p);
248 AssertThrow ( !((errno != 0) || (s.size() == 0) || ((s.size()>0) && (*p !=
'\0'))),
249 ExcMessage (
"Can't convert <" + s +
"> to a double."));
259 std::vector<double> tmp (s.size());
260 for (
unsigned int i=0; i<s.size(); ++i)
267 std::vector<std::string>
269 const char delimiter)
277 while (tmp.length() != 0 && tmp[tmp.length()-1] ==
' ')
278 tmp.erase (tmp.length()-1, 1);
286 std::vector<std::string> split_list;
287 split_list.reserve (std::count (tmp.begin(), tmp.end(), delimiter)+1);
288 while (tmp.length() != 0)
293 if (name.find(delimiter) != std::string::npos)
295 name.erase (name.find(delimiter), std::string::npos);
296 tmp.erase (0, tmp.find(delimiter)+1);
302 while ((name.length() != 0) && (name[0] ==
' '))
304 while (name.length() != 0 && name[name.length()-1] ==
' ')
305 name.erase (name.length()-1, 1);
307 split_list.push_back (name);
315 std::vector<std::string>
317 const unsigned int width,
318 const char delimiter)
320 std::string text = original_text;
321 std::vector<std::string> lines;
324 while ((text.length() != 0) && (text[text.length()-1] == delimiter))
325 text.erase(text.length()-1,1);
328 while (text.length() != 0)
332 while ((text.length() != 0) && (text[0] == delimiter))
335 std::size_t pos_newline = text.find_first_of(
"\n", 0);
336 if (pos_newline != std::string::npos && pos_newline <= width)
338 std::string line (text, 0, pos_newline);
339 while ((line.length() != 0) && (line[line.length()-1] == delimiter))
340 line.erase(line.length()-1,1);
341 lines.push_back (line);
342 text.erase (0, pos_newline+1);
349 if (text.length() < width)
352 while ((text.length() != 0) && (text[text.length()-1] == delimiter))
353 text.erase(text.length()-1,1);
354 lines.push_back (text);
362 int location = std::min<int>(width,text.length()-1);
363 for (; location>0; --location)
364 if (text[location] == delimiter)
370 for (location = std::min<int>(width,text.length()-1);
371 location<static_cast<int>(text.length());
373 if (text[location] == delimiter)
379 std::string line (text, 0, location);
380 while ((line.length() != 0) && (line[line.length()-1] == delimiter))
381 line.erase(line.length()-1,1);
382 lines.push_back (line);
383 text.erase (0, location);
394 const std::string &pattern)
396 if (pattern.size() > name.size())
399 for (
unsigned int i=0; i<pattern.size(); ++i)
400 if (pattern[i] != name[i])
408 std::pair<int, unsigned int>
410 const unsigned int position)
414 const std::string test_string (name.begin()+position,
417 std::istringstream str(test_string);
427 return std::make_pair (i, 1U);
429 return std::make_pair (i, 2U);
431 return std::make_pair (i, 3U);
433 return std::make_pair (i, 4U);
435 return std::make_pair (i, 5U);
437 return std::make_pair (i, 6U);
439 return std::make_pair (i, 7U);
466 return boost::normal_distribution<>(a,sigma)(random_number_generator.
get());
471 std::vector<unsigned int>
474 const unsigned int n = permutation.size();
476 std::vector<unsigned int> out (n);
477 for (
unsigned int i=0; i<n; ++i)
478 out[i] = n - 1 - permutation[i];
485 std::vector<unsigned int>
488 const unsigned int n = permutation.size();
492 for (
unsigned int i=0; i<n; ++i)
495 out[permutation[i]] = i;
500 for (
unsigned int i=0; i<n; ++i)
502 ExcMessage (
"The given input permutation had duplicate entries!"));
507 std::vector<unsigned long long int>
510 const unsigned long long int n = permutation.size();
512 std::vector<unsigned long long int> out (n);
513 for (
unsigned long long int i=0; i<n; ++i)
514 out[i] = n - 1 - permutation[i];
521 std::vector<unsigned long long int>
524 const unsigned long long int n = permutation.size();
528 for (
unsigned long long int i=0; i<n; ++i)
531 out[permutation[i]] = i;
536 for (
unsigned long long int i=0; i<n; ++i)
538 ExcMessage (
"The given input permutation had duplicate entries!"));
544 template <
typename Integer>
548 const unsigned int n = permutation.size();
550 std::vector<Integer> out (n);
551 for (
unsigned int i=0; i<n; ++i)
552 out[i] = n - 1 - permutation[i];
559 template <
typename Integer>
563 const unsigned int n = permutation.size();
567 for (
unsigned int i=0; i<n; ++i)
570 out[permutation[i]] = i;
575 for (
unsigned int i=0; i<n; ++i)
577 ExcMessage (
"The given input permutation had duplicate entries!"));
587 #if defined(__linux__) 591 std::ifstream cpuinfo;
592 cpuinfo.open(
"/proc/loadavg");
620 #if defined(__linux__) 621 std::ifstream file(
"/proc/self/status");
627 if (name ==
"VmPeak:")
629 else if (name ==
"VmSize:")
631 else if (name ==
"VmHWM:")
633 else if (name ==
"VmRSS:")
648 #if defined(DEAL_II_HAVE_UNISTD_H) && defined(DEAL_II_HAVE_GETHOSTNAME) 649 const unsigned int N=1024;
651 gethostname (&(hostname[0]), N-1);
653 std::string hostname(
"unknown");
662 std::time_t time1= std::time (0);
663 std::tm *time = std::localtime(&time1);
665 std::ostringstream o;
666 o << time->tm_hour <<
":" 667 << (time->tm_min < 10 ?
"0" :
"") << time->tm_min <<
":" 668 << (time->tm_sec < 10 ?
"0" :
"") << time->tm_sec;
677 std::time_t time1= std::time (0);
678 std::tm *time = std::localtime(&time1);
680 std::ostringstream o;
681 o << time->tm_year + 1900 <<
"/" 682 << time->tm_mon + 1 <<
"/" 700 *memptr = malloc (size);
715 #ifdef DEAL_II_WITH_TRILINOS 722 #ifdef DEAL_II_WITH_MPI 723 static Teuchos::RCP<Epetra_MpiComm>
724 communicator = Teuchos::rcp (
new Epetra_MpiComm (MPI_COMM_WORLD),
true);
726 static Teuchos::RCP<Epetra_SerialComm>
727 communicator = Teuchos::rcp (
new Epetra_SerialComm (),
true);
730 return *communicator;
738 #ifdef DEAL_II_WITH_MPI 739 static Teuchos::RCP<Epetra_MpiComm>
740 communicator = Teuchos::rcp (
new Epetra_MpiComm (MPI_COMM_SELF),
true);
742 static Teuchos::RCP<Epetra_SerialComm>
743 communicator = Teuchos::rcp (
new Epetra_SerialComm (),
true);
746 return *communicator;
754 #ifdef DEAL_II_WITH_MPI 760 *mpi_comm =
dynamic_cast<const Epetra_MpiComm *
>(&communicator);
770 Assert (dynamic_cast<const Epetra_SerialComm *>(&communicator)
773 return new Epetra_SerialComm(dynamic_cast<const Epetra_SerialComm &>(communicator));
782 #ifdef DEAL_II_WITH_MPI 784 *mpi_comm =
dynamic_cast<Epetra_MpiComm *
>(&communicator);
787 MPI_Comm comm = mpi_comm->GetMpiComm();
788 *mpi_comm = Epetra_MpiComm(MPI_COMM_SELF);
789 const int ierr = MPI_Comm_free (&comm);
799 return mpi_communicator.NumProc();
805 return (
unsigned int)mpi_communicator.MyPID();
812 const Epetra_Comm &comm)
814 if (map.LinearMap() ==
true)
821 return Epetra_Map (map.NumGlobalElements(),
830 return Epetra_Map (map.NumGlobalElements(),
832 map.MyGlobalElements (),
841 template std::string to_string<int> (int,
unsigned int);
842 template std::string to_string<long int> (
long int,
unsigned int);
843 template std::string to_string<long long int> (
long long int,
unsigned int);
844 template std::string to_string<unsigned int> (
unsigned int,
unsigned int);
845 template std::string to_string<unsigned long int> (
unsigned long int,
unsigned int);
846 template std::string to_string<unsigned long long int> (
unsigned long long int,
unsigned int);
847 template std::string to_string<float> (float,
unsigned int);
848 template std::string to_string<double> (double,
unsigned int);
849 template std::string to_string<long double> (
long double,
unsigned int);
853 DEAL_II_NAMESPACE_CLOSE
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()
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)
std::vector< std::string > split_string_list(const std::string &s, const char delimiter=',')
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)
bool job_supports_mpi() 1
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)
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)