Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
utilities.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2005 - 2019 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/config.h>
17 
18 // It's necessary to include winsock2.h before thread_local_storage.h,
19 // because Intel implementation of TBB includes winsock.h,
20 // and we'll get a conflict between winsock.h and winsock2.h otherwise.
21 #ifdef DEAL_II_MSVC
22 # include <winsock2.h>
23 #endif
24 
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>
30 
31 #include <boost/lexical_cast.hpp>
32 #include <boost/random.hpp>
33 
34 #include <algorithm>
35 #include <bitset>
36 #include <cctype>
37 #include <cerrno>
38 #include <cmath>
39 #include <cstddef>
40 #include <cstdio>
41 #include <ctime>
42 #include <fstream>
43 #include <iomanip>
44 #include <iostream>
45 #include <limits>
46 #include <sstream>
47 
48 #if defined(DEAL_II_HAVE_UNISTD_H) && defined(DEAL_II_HAVE_GETHOSTNAME)
49 # include <unistd.h>
50 #endif
51 
52 #ifndef DEAL_II_MSVC
53 # include <cstdlib>
54 #endif
55 
56 
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>
62 
63 # include <Epetra_MpiComm.h>
64 # include <Teuchos_DefaultComm.hpp>
65 # endif
66 # include <Epetra_SerialComm.h>
67 # include <Teuchos_RCP.hpp>
68 #endif
69 
70 DEAL_II_NAMESPACE_OPEN
71 
72 
73 namespace Utilities
74 {
76  unsigned int,
77  unsigned int,
78  << "When trying to convert " << arg1 << " to a string with "
79  << arg2 << " digits");
80  DeclException1(ExcInvalidNumber, unsigned int, << "Invalid number " << arg1);
82  std::string,
83  << "Can't convert the string " << arg1
84  << " to the desired type");
85 
86 
87  std::string
89  {
90  return DEAL_II_PACKAGE_NAME " version " DEAL_II_PACKAGE_VERSION;
91  }
92 
93 
94  namespace
95  {
96  template <int dim,
97  typename Number,
98  int effective_dim,
99  typename LongDouble,
100  typename Integer>
101  std::vector<std::array<std::uint64_t, effective_dim>>
102  inverse_Hilbert_space_filling_curve_effective(
103  const std::vector<Point<dim, Number>> &points,
104  const Point<dim, Number> & bl,
105  const std::array<LongDouble, dim> & extents,
106  const std::bitset<dim> & valid_extents,
107  const int min_bits,
108  const Integer max_int)
109  {
110  std::vector<std::array<Integer, effective_dim>> int_points(points.size());
111 
112  for (unsigned int i = 0; i < points.size(); ++i)
113  {
114  // convert into integers:
115  unsigned int eff_d = 0;
116  for (unsigned int d = 0; d < dim; ++d)
117  if (valid_extents[d])
118  {
119  Assert(extents[d] > 0, ExcInternalError());
120  const LongDouble v = (static_cast<LongDouble>(points[i][d]) -
121  static_cast<LongDouble>(bl[d])) /
122  extents[d];
123  Assert(v >= 0. && v <= 1., ExcInternalError());
124  AssertIndexRange(eff_d, effective_dim);
125  int_points[i][eff_d] =
126  static_cast<Integer>(v * static_cast<LongDouble>(max_int));
127  ++eff_d;
128  }
129  }
130 
131  // note that we call this with "min_bits"
132  return inverse_Hilbert_space_filling_curve<effective_dim>(int_points,
133  min_bits);
134  }
135  } // namespace
136 
137  template <int dim, typename Number>
138  std::vector<std::array<std::uint64_t, dim>>
140  const std::vector<Point<dim, Number>> &points,
141  const int bits_per_dim)
142  {
143  using Integer = std::uint64_t;
144  // take floating point number hopefully with mantissa >= 64bit
145  using LongDouble = long double;
146 
147  // return if there is nothing to do
148  if (points.size() == 0)
149  return std::vector<std::array<std::uint64_t, dim>>();
150 
151  // get bounding box:
152  Point<dim, Number> bl = points[0], tr = points[0];
153  for (const auto &p : points)
154  for (unsigned int d = 0; d < dim; ++d)
155  {
156  const double cid = p[d];
157  bl[d] = std::min(cid, bl[d]);
158  tr[d] = std::max(cid, tr[d]);
159  }
160 
161  std::array<LongDouble, dim> extents;
162  std::bitset<dim> valid_extents;
163  for (unsigned int i = 0; i < dim; ++i)
164  {
165  extents[i] =
166  static_cast<LongDouble>(tr[i]) - static_cast<LongDouble>(bl[i]);
167  valid_extents[i] = (extents[i] > 0.);
168  }
169 
170  // make sure our conversion from fractional coordinates to
171  // Integers work as expected, namely our cast (LongDouble)max_int
172  const int min_bits =
173  std::min(bits_per_dim,
174  std::min(std::numeric_limits<Integer>::digits,
175  std::numeric_limits<LongDouble>::digits));
176 
177  // based on that get the maximum integer:
178  const Integer max_int = (min_bits == std::numeric_limits<Integer>::digits ?
179  std::numeric_limits<Integer>::max() :
180  (Integer(1) << min_bits) - 1);
181 
182  const unsigned int effective_dim = valid_extents.count();
183  if (effective_dim == dim)
184  {
185  return inverse_Hilbert_space_filling_curve_effective<dim,
186  Number,
187  dim,
188  LongDouble,
189  Integer>(
190  points, bl, extents, valid_extents, min_bits, max_int);
191  }
192 
193  // various degenerate cases
194  std::array<std::uint64_t, dim> zero_ind;
195  for (unsigned int d = 0; d < dim; ++d)
196  zero_ind[d] = 0;
197 
198  std::vector<std::array<std::uint64_t, dim>> ind(points.size(), zero_ind);
199  // manually check effective_dim == 1 and effective_dim == 2
200  if (dim == 3 && effective_dim == 2)
201  {
202  const auto ind2 =
203  inverse_Hilbert_space_filling_curve_effective<dim,
204  Number,
205  2,
206  LongDouble,
207  Integer>(
208  points, bl, extents, valid_extents, min_bits, max_int);
209 
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];
213 
214  return ind;
215  }
216  else if (effective_dim == 1)
217  {
218  const auto ind1 =
219  inverse_Hilbert_space_filling_curve_effective<dim,
220  Number,
221  1,
222  LongDouble,
223  Integer>(
224  points, bl, extents, valid_extents, min_bits, max_int);
225 
226  for (unsigned int i = 0; i < ind.size(); ++i)
227  ind[i][dim - 1] = ind1[i][0];
228 
229  return ind;
230  }
231 
232  // we should get here only if effective_dim == 0
233  Assert(effective_dim == 0, ExcInternalError());
234 
235  // if the bounding box is degenerate in all dimensions,
236  // can't do much but exit gracefully by setting index according
237  // to the index of each point so that there is no re-ordering
238  for (unsigned int i = 0; i < points.size(); ++i)
239  ind[i][dim - 1] = i;
240 
241  return ind;
242  }
243 
244 
245 
246  template <int dim>
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)
251  {
252  using Integer = std::uint64_t;
253 
254  std::vector<std::array<Integer, dim>> int_points(points);
255 
256  std::vector<std::array<Integer, dim>> res(int_points.size());
257 
258  // follow
259  // J. Skilling, Programming the Hilbert curve, AIP Conf. Proc. 707, 381
260  // (2004); http://dx.doi.org/10.1063/1.1751381 also see
261  // https://stackoverflow.com/questions/499166/mapping-n-dimensional-value-to-a-point-on-hilbert-curve
262  // https://gitlab.com/octopus-code/octopus/blob/develop/src/grid/hilbert.c
263  // https://github.com/trilinos/Trilinos/blob/master/packages/zoltan/src/hsfc/hsfc_hilbert.c
264  // (Zoltan_HSFC_InvHilbertXd)
265  // https://github.com/aditi137/Hilbert/blob/master/Hilbert/hilbert.cpp
266 
267  // now we can map to 1D coordinate stored in Transpose format
268  // adopt AxestoTranspose function from the paper, that
269  // transforms in-place between geometrical axes and Hilbert transpose.
270  // Example: b=5 bits for each of n=3 coordinates.
271  // 15-bit Hilbert integer = A B C D E F G H I J K L M N O is
272  // stored as its Transpose
273  // X[0] = A D G J M X[2]|
274  // X[1] = B E H K N <-------> | /X[1]
275  // X[2] = C F I L O axes |/
276  // high low 0------ X[0]
277 
278  // Depth of the Hilbert curve
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."));
282 
283  const Integer M = Integer(1) << (bits_per_dim - 1); // largest bit
284 
285  for (unsigned int index = 0; index < int_points.size(); ++index)
286  {
287  auto &X = int_points[index];
288  auto &L = res[index];
289 
290  // Inverse undo
291  for (Integer q = M; q > 1; q >>= 1)
292  {
293  const Integer p = q - 1;
294  for (unsigned int i = 0; i < dim; i++)
295  {
296  // invert
297  if (X[i] & q)
298  {
299  X[0] ^= p;
300  }
301  // exchange
302  else
303  {
304  const Integer t = (X[0] ^ X[i]) & p;
305  X[0] ^= t;
306  X[i] ^= t;
307  }
308  }
309  }
310 
311  // Gray encode (inverse of decode)
312  for (unsigned int i = 1; i < dim; i++)
313  X[i] ^= X[i - 1];
314 
315  Integer t = 0;
316  for (Integer q = M; q > 1; q >>= 1)
317  if (X[dim - 1] & q)
318  t ^= q - 1;
319  for (unsigned int i = 0; i < dim; i++)
320  X[i] ^= t;
321 
322  // now we need to go from index stored in transpose format to
323  // consecutive format, which is better suited for comparators.
324  // we could interleave into some big unsigned int...
325  // https://www.forceflow.be/2013/10/07/morton-encodingdecoding-through-bit-interleaving-implementations/
326  // https://stackoverflow.com/questions/4431522/given-2-16-bit-ints-can-i-interleave-those-bits-to-form-a-single-32-bit-int
327  // ...but we would loose spatial resolution!
328 
329  // interleave using brute force, follow TransposetoLine from
330  // https://github.com/aditi137/Hilbert/blob/master/Hilbert/hilbert.cpp
331  {
332  Integer p = M;
333  unsigned int j = 0;
334  for (unsigned int i = 0; i < dim; ++i)
335  {
336  L[i] = 0;
337  // go through bits using a mask q
338  for (Integer q = M; q > 0; q >>= 1)
339  {
340  if (X[j] & p)
341  L[i] |= q;
342  if (++j == dim)
343  {
344  j = 0;
345  p >>= 1;
346  }
347  }
348  }
349  }
350 
351  } // end of the loop over points
352 
353  return res;
354  }
355 
356 
357 
358  template <int dim>
359  std::uint64_t
360  pack_integers(const std::array<std::uint64_t, dim> &index,
361  const int bits_per_dim)
362  {
363  using Integer = std::uint64_t;
364 
365  AssertIndexRange(bits_per_dim * dim, 65);
366  Assert(bits_per_dim > 0, ExcMessage("bits_per_dim should be positive"));
367 
368  const Integer mask = (Integer(1) << bits_per_dim) - 1;
369 
370  Integer res = 0;
371  for (unsigned int i = 0; i < dim; ++i)
372  {
373  // take bits_per_dim from each integer and shift them
374  const Integer v = (mask & index[dim - 1 - i]) << (bits_per_dim * i);
375  res |= v;
376  }
377  return res;
378  }
379 
380 
381 
382  std::string
383  int_to_string(const unsigned int value, const unsigned int digits)
384  {
385  return to_string(value, digits);
386  }
387 
388 
389 
390  template <typename number>
391  std::string
392  to_string(const number value, const unsigned int digits)
393  {
394  std::string lc_string = boost::lexical_cast<std::string>(value);
395 
396  if (digits == numbers::invalid_unsigned_int)
397  return lc_string;
398  else if (lc_string.size() < digits)
399  {
400  // We have to add the padding zeroes in front of the number
401  const unsigned int padding_position = (lc_string[0] == '-') ? 1 : 0;
402 
403  const std::string padding(digits - lc_string.size(), '0');
404  lc_string.insert(padding_position, padding);
405  }
406  return lc_string;
407  }
408 
409 
410  std::string
411  replace_in_string(const std::string &input,
412  const std::string &from,
413  const std::string &to)
414  {
415  if (from.empty())
416  return input;
417 
418  std::string out = input;
419  std::string::size_type pos = out.find(from);
420 
421  while (pos != std::string::npos)
422  {
423  out.replace(pos, from.size(), to);
424  pos = out.find(from, pos + to.size());
425  }
426  return out;
427  }
428 
429  std::string
430  trim(const std::string &input)
431  {
432  std::string::size_type left = 0;
433  std::string::size_type right = input.size() > 0 ? input.size() - 1 : 0;
434 
435  for (; left < input.size(); ++left)
436  {
437  if (!std::isspace(input[left]))
438  {
439  break;
440  }
441  }
442 
443  for (; right >= left; --right)
444  {
445  if (!std::isspace(input[right]))
446  {
447  break;
448  }
449  }
450 
451  return std::string(input, left, right - left + 1);
452  }
453 
454 
455 
456  std::string
457  dim_string(const int dim, const int spacedim)
458  {
459  if (dim == spacedim)
460  return int_to_string(dim);
461  else
462  return int_to_string(dim) + "," + int_to_string(spacedim);
463  }
464 
465 
466  unsigned int
467  needed_digits(const unsigned int max_number)
468  {
469  if (max_number < 10)
470  return 1;
471  if (max_number < 100)
472  return 2;
473  if (max_number < 1000)
474  return 3;
475  if (max_number < 10000)
476  return 4;
477  if (max_number < 100000)
478  return 5;
479  if (max_number < 1000000)
480  return 6;
481  AssertThrow(false, ExcInvalidNumber(max_number));
482  return 0;
483  }
484 
485 
486 
487  int
488  string_to_int(const std::string &s_)
489  {
490  // trim whitespace on either side of the text if necessary
491  std::string s = s_;
492  while ((s.size() > 0) && (s[0] == ' '))
493  s.erase(s.begin());
494  while ((s.size() > 0) && (s[s.size() - 1] == ' '))
495  s.erase(s.end() - 1);
496 
497  // now convert and see whether we succeed. note that strtol only
498  // touches errno if an error occurred, so if we want to check
499  // whether an error happened, we need to make sure that errno==0
500  // before calling strtol since otherwise it may be that the
501  // conversion succeeds and that errno remains at the value it
502  // was before, whatever that was
503  char *p;
504  errno = 0;
505  const int i = std::strtol(s.c_str(), &p, 10);
506  AssertThrow(!((errno != 0) || (s.size() == 0) ||
507  ((s.size() > 0) && (*p != '\0'))),
508  ExcMessage("Can't convert <" + s + "> to an integer."));
509 
510  return i;
511  }
512 
513 
514 
515  std::vector<int>
516  string_to_int(const std::vector<std::string> &s)
517  {
518  std::vector<int> tmp(s.size());
519  for (unsigned int i = 0; i < s.size(); ++i)
520  tmp[i] = string_to_int(s[i]);
521  return tmp;
522  }
523 
524 
525 
526  double
527  string_to_double(const std::string &s_)
528  {
529  // trim whitespace on either side of the text if necessary
530  std::string s = s_;
531  while ((s.size() > 0) && (s[0] == ' '))
532  s.erase(s.begin());
533  while ((s.size() > 0) && (s[s.size() - 1] == ' '))
534  s.erase(s.end() - 1);
535 
536  // now convert and see whether we succeed. note that strtol only
537  // touches errno if an error occurred, so if we want to check
538  // whether an error happened, we need to make sure that errno==0
539  // before calling strtol since otherwise it may be that the
540  // conversion succeeds and that errno remains at the value it
541  // was before, whatever that was
542  char *p;
543  errno = 0;
544  const double d = std::strtod(s.c_str(), &p);
545  AssertThrow(!((errno != 0) || (s.size() == 0) ||
546  ((s.size() > 0) && (*p != '\0'))),
547  ExcMessage("Can't convert <" + s + "> to a double."));
548 
549  return d;
550  }
551 
552 
553 
554  std::vector<double>
555  string_to_double(const std::vector<std::string> &s)
556  {
557  std::vector<double> tmp(s.size());
558  for (unsigned int i = 0; i < s.size(); ++i)
559  tmp[i] = string_to_double(s[i]);
560  return tmp;
561  }
562 
563 
564 
565  std::vector<std::string>
566  split_string_list(const std::string &s, const std::string &delimiter)
567  {
568  // keep the currently remaining part of the input string in 'tmp' and
569  // keep chopping elements of the list off the front
570  std::string tmp = s;
571 
572  // as discussed in the documentation, eat whitespace from the end
573  // of the string
574  while (tmp.length() != 0 && tmp[tmp.length() - 1] == ' ')
575  tmp.erase(tmp.length() - 1, 1);
576 
577  // split the input list until it is empty. since in every iteration
578  // 'tmp' is what's left of the string after the next delimiter,
579  // and since we've stripped trailing space already, 'tmp' will
580  // be empty at one point if 's' ended in a delimiter, even if
581  // there was space after the last delimiter. this matches what's
582  // discussed in the documentation
583  std::vector<std::string> split_list;
584  while (tmp.length() != 0)
585  {
586  std::string name;
587  name = tmp;
588 
589  if (name.find(delimiter) != std::string::npos)
590  {
591  name.erase(name.find(delimiter), std::string::npos);
592  tmp.erase(0, tmp.find(delimiter) + delimiter.size());
593  }
594  else
595  tmp = "";
596 
597  // strip spaces from this element's front and end
598  while ((name.length() != 0) && (name[0] == ' '))
599  name.erase(0, 1);
600  while (name.length() != 0 && name[name.length() - 1] == ' ')
601  name.erase(name.length() - 1, 1);
602 
603  split_list.push_back(name);
604  }
605 
606  return split_list;
607  }
608 
609 
610  std::vector<std::string>
611  split_string_list(const std::string &s, const char delimiter)
612  {
613  std::string d = ",";
614  d[0] = delimiter;
615  return split_string_list(s, d);
616  }
617 
618 
619  std::vector<std::string>
620  break_text_into_lines(const std::string &original_text,
621  const unsigned int width,
622  const char delimiter)
623  {
624  std::string text = original_text;
625  std::vector<std::string> lines;
626 
627  // remove trailing spaces
628  while ((text.length() != 0) && (text[text.length() - 1] == delimiter))
629  text.erase(text.length() - 1, 1);
630 
631  // then split the text into lines
632  while (text.length() != 0)
633  {
634  // in each iteration, first remove
635  // leading spaces
636  while ((text.length() != 0) && (text[0] == delimiter))
637  text.erase(0, 1);
638 
639  std::size_t pos_newline = text.find_first_of('\n', 0);
640  if (pos_newline != std::string::npos && pos_newline <= width)
641  {
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);
648  continue;
649  }
650 
651  // if we can fit everything into one
652  // line, then do so. otherwise, we have
653  // to keep breaking
654  if (text.length() < width)
655  {
656  // remove trailing spaces
657  while ((text.length() != 0) &&
658  (text[text.length() - 1] == delimiter))
659  text.erase(text.length() - 1, 1);
660  lines.push_back(text);
661  text = "";
662  }
663  else
664  {
665  // starting at position width, find the
666  // location of the previous space, so
667  // that we can break around there
668  int location = std::min<int>(width, text.length() - 1);
669  for (; location > 0; --location)
670  if (text[location] == delimiter)
671  break;
672 
673  // if there are no spaces, then try if
674  // there are spaces coming up
675  if (location == 0)
676  for (location = std::min<int>(width, text.length() - 1);
677  location < static_cast<int>(text.length());
678  ++location)
679  if (text[location] == delimiter)
680  break;
681 
682  // now take the text up to the found
683  // location and put it into a single
684  // line, and remove it from 'text'
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);
691  }
692  }
693 
694  return lines;
695  }
696 
697 
698 
699  bool
700  match_at_string_start(const std::string &name, const std::string &pattern)
701  {
702  if (pattern.size() > name.size())
703  return false;
704 
705  for (unsigned int i = 0; i < pattern.size(); ++i)
706  if (pattern[i] != name[i])
707  return false;
708 
709  return true;
710  }
711 
712 
713 
714  std::pair<int, unsigned int>
715  get_integer_at_position(const std::string &name, const unsigned int position)
716  {
717  Assert(position < name.size(), ExcInternalError());
718 
719  const std::string test_string(name.begin() + position, name.end());
720 
721  std::istringstream str(test_string);
722 
723  int i;
724  if (str >> i)
725  {
726  // compute the number of
727  // digits of i. assuming it
728  // is less than 8 is likely
729  // ok
730  if (i < 10)
731  return std::make_pair(i, 1U);
732  else if (i < 100)
733  return std::make_pair(i, 2U);
734  else if (i < 1000)
735  return std::make_pair(i, 3U);
736  else if (i < 10000)
737  return std::make_pair(i, 4U);
738  else if (i < 100000)
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);
744  else
745  {
746  Assert(false, ExcNotImplemented());
747  return std::make_pair(-1, numbers::invalid_unsigned_int);
748  }
749  }
750  else
751  return std::make_pair(-1, numbers::invalid_unsigned_int);
752  }
753 
754 
755 
756  double
757  generate_normal_random_number(const double a, const double sigma)
758  {
759  // if no noise: return now
760  if (sigma == 0)
761  return a;
762 
763  // we would want to use rand(), but that function is not reentrant
764  // in a thread context. one could use rand_r, but this does not
765  // produce reproducible results between threads either (though at
766  // least it is reentrant). these two approaches being
767  // non-workable, use a thread-local random number generator here
768  static Threads::ThreadLocalStorage<boost::mt19937> random_number_generator;
769  return boost::normal_distribution<>(a,
770  sigma)(random_number_generator.get());
771  }
772 
773 
774 
775  std::vector<unsigned int>
776  reverse_permutation(const std::vector<unsigned int> &permutation)
777  {
778  const unsigned int n = permutation.size();
779 
780  std::vector<unsigned int> out(n);
781  for (unsigned int i = 0; i < n; ++i)
782  out[i] = n - 1 - permutation[i];
783 
784  return out;
785  }
786 
787 
788 
789  std::vector<unsigned int>
790  invert_permutation(const std::vector<unsigned int> &permutation)
791  {
792  const unsigned int n = permutation.size();
793 
794  std::vector<unsigned int> out(n, numbers::invalid_unsigned_int);
795 
796  for (unsigned int i = 0; i < n; ++i)
797  {
798  Assert(permutation[i] < n, ExcIndexRange(permutation[i], 0, n));
799  out[permutation[i]] = i;
800  }
801 
802  // check that we have actually reached
803  // all indices
804  for (unsigned int i = 0; i < n; ++i)
806  ExcMessage("The given input permutation had duplicate entries!"));
807 
808  return out;
809  }
810 
811  std::vector<unsigned long long int>
812  reverse_permutation(const std::vector<unsigned long long int> &permutation)
813  {
814  const unsigned long long int n = permutation.size();
815 
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];
819 
820  return out;
821  }
822 
823 
824 
825  std::vector<unsigned long long int>
826  invert_permutation(const std::vector<unsigned long long int> &permutation)
827  {
828  const unsigned long long int n = permutation.size();
829 
830  std::vector<unsigned long long int> out(n, numbers::invalid_unsigned_int);
831 
832  for (unsigned long long int i = 0; i < n; ++i)
833  {
834  Assert(permutation[i] < n, ExcIndexRange(permutation[i], 0, n));
835  out[permutation[i]] = i;
836  }
837 
838  // check that we have actually reached
839  // all indices
840  for (unsigned long long int i = 0; i < n; ++i)
842  ExcMessage("The given input permutation had duplicate entries!"));
843 
844  return out;
845  }
846 
847 
848  template <typename Integer>
849  std::vector<Integer>
850  reverse_permutation(const std::vector<Integer> &permutation)
851  {
852  const unsigned int n = permutation.size();
853 
854  std::vector<Integer> out(n);
855  for (unsigned int i = 0; i < n; ++i)
856  out[i] = n - 1 - permutation[i];
857 
858  return out;
859  }
860 
861 
862 
863  template <typename Integer>
864  std::vector<Integer>
865  invert_permutation(const std::vector<Integer> &permutation)
866  {
867  const unsigned int n = permutation.size();
868 
869  std::vector<Integer> out(n, numbers::invalid_unsigned_int);
870 
871  for (unsigned int i = 0; i < n; ++i)
872  {
873  Assert(permutation[i] < n, ExcIndexRange(permutation[i], 0, n));
874  out[permutation[i]] = i;
875  }
876 
877  // check that we have actually reached
878  // all indices
879  for (unsigned int i = 0; i < n; ++i)
881  ExcMessage("The given input permutation had duplicate entries!"));
882 
883  return out;
884  }
885 
886 
887 
888  namespace System
889  {
890 #if defined(__linux__)
891 
892  double
893  get_cpu_load()
894  {
895  std::ifstream cpuinfo;
896  cpuinfo.open("/proc/loadavg");
897 
898  AssertThrow(cpuinfo, ExcIO());
899 
900  double load;
901  cpuinfo >> load;
902 
903  return load;
904  }
905 
906 #else
907 
908  double
910  {
911  return 0.;
912  }
913 
914 #endif
915 
916  const std::string
918  {
919  switch (DEAL_II_COMPILER_VECTORIZATION_LEVEL)
920  {
921  case 0:
922  return "disabled";
923  case 1:
924 #ifdef __ALTIVEC__
925  return "AltiVec";
926 #else
927  return "SSE2";
928 #endif
929  case 2:
930  return "AVX";
931  case 3:
932  return "AVX512";
933  default:
934  AssertThrow(false,
936  "Invalid DEAL_II_COMPILER_VECTORIZATION_LEVEL."));
937  return "ERROR";
938  }
939  }
940 
941 
942  void
944  {
945  stats.VmPeak = stats.VmSize = stats.VmHWM = stats.VmRSS = 0;
946 
947  // parsing /proc/self/stat would be a
948  // lot easier, but it does not contain
949  // VmHWM, so we use /status instead.
950 #if defined(__linux__)
951  std::ifstream file("/proc/self/status");
952  std::string line;
953  std::string name;
954  while (!file.eof())
955  {
956  file >> name;
957  if (name == "VmPeak:")
958  file >> stats.VmPeak;
959  else if (name == "VmSize:")
960  file >> stats.VmSize;
961  else if (name == "VmHWM:")
962  file >> stats.VmHWM;
963  else if (name == "VmRSS:")
964  {
965  file >> stats.VmRSS;
966  break; // this is always the last entry
967  }
968 
969  getline(file, line);
970  }
971 #endif
972  }
973 
974 
975 
976  std::string
978  {
979 #if defined(DEAL_II_HAVE_UNISTD_H) && defined(DEAL_II_HAVE_GETHOSTNAME)
980  const unsigned int N = 1024;
981  char hostname[N];
982  gethostname(&(hostname[0]), N - 1);
983 #else
984  std::string hostname("unknown");
985 #endif
986  return hostname;
987  }
988 
989 
990 
991  std::string
993  {
994  std::time_t time1 = std::time(nullptr);
995  std::tm * time = std::localtime(&time1);
996 
997  std::ostringstream o;
998  o << time->tm_hour << ":" << (time->tm_min < 10 ? "0" : "")
999  << time->tm_min << ":" << (time->tm_sec < 10 ? "0" : "")
1000  << time->tm_sec;
1001 
1002  return o.str();
1003  }
1004 
1005 
1006 
1007  std::string
1009  {
1010  std::time_t time1 = std::time(nullptr);
1011  std::tm * time = std::localtime(&time1);
1012 
1013  std::ostringstream o;
1014  o << time->tm_year + 1900 << "/" << time->tm_mon + 1 << "/"
1015  << time->tm_mday;
1016 
1017  return o.str();
1018  }
1019 
1020 
1021 
1022  void
1023  posix_memalign(void **memptr, std::size_t alignment, std::size_t size)
1024  {
1025 #ifndef DEAL_II_MSVC
1026  const int ierr = ::posix_memalign(memptr, alignment, size);
1027 
1028  AssertThrow(ierr == 0, ExcOutOfMemory());
1029  AssertThrow(*memptr != nullptr, ExcOutOfMemory());
1030 #else
1031  // Windows does not appear to have posix_memalign. just use the
1032  // regular malloc in that case
1033  *memptr = malloc(size);
1034  (void)alignment;
1035  AssertThrow(*memptr != 0, ExcOutOfMemory());
1036 #endif
1037  }
1038 
1039 
1040 
1041  bool
1042  job_supports_mpi()
1043  {
1045  }
1046  } // namespace System
1047 
1048 
1049 #ifdef DEAL_II_WITH_TRILINOS
1050 
1051  namespace Trilinos
1052  {
1053  const Epetra_Comm &
1055  {
1056 # ifdef DEAL_II_WITH_MPI
1057  static Teuchos::RCP<Epetra_MpiComm> communicator =
1058  Teuchos::rcp(new Epetra_MpiComm(MPI_COMM_WORLD), true);
1059 # else
1060  static Teuchos::RCP<Epetra_SerialComm> communicator =
1061  Teuchos::rcp(new Epetra_SerialComm(), true);
1062 # endif
1063 
1064  return *communicator;
1065  }
1066 
1067 
1068 
1069  const Teuchos::RCP<const Teuchos::Comm<int>> &
1071  {
1072 # ifdef DEAL_II_WITH_MPI
1073  static auto communicator = Teuchos::RCP<const Teuchos::Comm<int>>(
1074  new Teuchos::MpiComm<int>(MPI_COMM_SELF));
1075 # else
1076  static auto communicator =
1077  Teuchos::RCP<const Teuchos::Comm<int>>(new Teuchos::Comm<int>());
1078 # endif
1079 
1080  return communicator;
1081  }
1082 
1083 
1084 
1085  const Epetra_Comm &
1087  {
1088 # ifdef DEAL_II_WITH_MPI
1089  static Teuchos::RCP<Epetra_MpiComm> communicator =
1090  Teuchos::rcp(new Epetra_MpiComm(MPI_COMM_SELF), true);
1091 # else
1092  static Teuchos::RCP<Epetra_SerialComm> communicator =
1093  Teuchos::rcp(new Epetra_SerialComm(), true);
1094 # endif
1095 
1096  return *communicator;
1097  }
1098 
1099 
1100 
1101  Epetra_Comm *
1102  duplicate_communicator(const Epetra_Comm &communicator)
1103  {
1104 # ifdef DEAL_II_WITH_MPI
1105 
1106  // see if the communicator is in fact a
1107  // parallel MPI communicator; if so,
1108  // return a duplicate of it
1109  const Epetra_MpiComm *mpi_comm =
1110  dynamic_cast<const Epetra_MpiComm *>(&communicator);
1111  if (mpi_comm != nullptr)
1112  return new Epetra_MpiComm(
1113  Utilities::MPI::duplicate_communicator(mpi_comm->GetMpiComm()));
1114 # endif
1115 
1116  // if we don't support MPI, or if the
1117  // communicator in question was in fact
1118  // not an MPI communicator, return a
1119  // copy of the same object again
1120  Assert(dynamic_cast<const Epetra_SerialComm *>(&communicator) != nullptr,
1121  ExcInternalError());
1122  return new Epetra_SerialComm(
1123  dynamic_cast<const Epetra_SerialComm &>(communicator));
1124  }
1125 
1126 
1127 
1128  void
1129  destroy_communicator(Epetra_Comm &communicator)
1130  {
1131  // save the communicator, reset the map, and delete the communicator if
1132  // this whole thing was created as an MPI communicator
1133 # ifdef DEAL_II_WITH_MPI
1134  Epetra_MpiComm *mpi_comm = dynamic_cast<Epetra_MpiComm *>(&communicator);
1135  if (mpi_comm != nullptr)
1136  {
1137  MPI_Comm comm = mpi_comm->GetMpiComm();
1138  *mpi_comm = Epetra_MpiComm(MPI_COMM_SELF);
1139  const int ierr = MPI_Comm_free(&comm);
1140  AssertThrowMPI(ierr);
1141  }
1142 # endif
1143  }
1144 
1145 
1146 
1147  unsigned int
1148  get_n_mpi_processes(const Epetra_Comm &mpi_communicator)
1149  {
1150  return mpi_communicator.NumProc();
1151  }
1152 
1153 
1154  unsigned int
1155  get_this_mpi_process(const Epetra_Comm &mpi_communicator)
1156  {
1157  return static_cast<unsigned int>(mpi_communicator.MyPID());
1158  }
1159 
1160 
1161 
1162  Epetra_Map
1163  duplicate_map(const Epetra_BlockMap &map, const Epetra_Comm &comm)
1164  {
1165  if (map.LinearMap() == true)
1166  {
1167  // each processor stores a
1168  // contiguous range of
1169  // elements in the
1170  // following constructor
1171  // call
1172  return Epetra_Map(map.NumGlobalElements(),
1173  map.NumMyElements(),
1174  map.IndexBase(),
1175  comm);
1176  }
1177  else
1178  {
1179  // the range is not
1180  // contiguous
1181  return Epetra_Map(map.NumGlobalElements(),
1182  map.NumMyElements(),
1183  map.MyGlobalElements(),
1184  0,
1185  comm);
1186  }
1187  }
1188  } // namespace Trilinos
1189 
1190 #endif
1191 
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);
1210 
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>> &,
1214  const int);
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>> &,
1218  const int);
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>> &,
1222  const int);
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>> &,
1226  const int);
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>> &,
1230  const int);
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>> &,
1234  const int);
1235 
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);
1242 } // namespace Utilities
1243 
1244 DEAL_II_NAMESPACE_CLOSE
std::vector< std::string > split_string_list(const std::string &s, const std::string &delimiter=",")
Definition: utilities.cc:566
void posix_memalign(void **memptr, std::size_t alignment, std::size_t size)
Definition: utilities.cc:1023
static const unsigned int invalid_unsigned_int
Definition: types.h:173
std::uint64_t pack_integers(const std::array< std::uint64_t, dim > &index, const int bits_per_dim)
Definition: utilities.cc:360
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:541
A class that provides a separate storage location on each thread that accesses the object...
static ::ExceptionBase & ExcIO()
std::string dealii_version_string()
Definition: utilities.cc:88
unsigned int get_n_mpi_processes(const Epetra_Comm &mpi_communicator)
Definition: utilities.cc:1148
std::string trim(const std::string &input)
Definition: utilities.cc:430
#define AssertIndexRange(index, range)
Definition: exceptions.h:1637
std::pair< int, unsigned int > get_integer_at_position(const std::string &name, const unsigned int position)
Definition: utilities.cc:715
std::vector< std::string > break_text_into_lines(const std::string &original_text, const unsigned int width, const char delimiter=' ')
Definition: utilities.cc:620
#define AssertThrow(cond, exc)
Definition: exceptions.h:1519
static ::ExceptionBase & ExcOutOfMemory()
const Epetra_Comm & comm_self()
Definition: utilities.cc:1086
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
Definition: point.h:110
std::string get_date()
Definition: utilities.cc:1008
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:392
void get_memory_stats(MemoryStats &stats)
Definition: utilities.cc:943
double string_to_double(const std::string &s)
Definition: utilities.cc:527
double get_cpu_load()
Definition: utilities.cc:909
const Teuchos::RCP< const Teuchos::Comm< int > > & tpetra_comm_self()
Definition: utilities.cc:1070
static ::ExceptionBase & ExcMessage(std::string arg1)
Epetra_Comm * duplicate_communicator(const Epetra_Comm &communicator)
Definition: utilities.cc:1102
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:518
double generate_normal_random_number(const double a, const double sigma)
Definition: utilities.cc:757
#define Assert(cond, exc)
Definition: exceptions.h:1407
unsigned long int VmSize
Definition: utilities.h:820
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)
Definition: utilities.cc:700
std::vector< unsigned int > invert_permutation(const std::vector< unsigned int > &permutation)
Definition: utilities.cc:790
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:383
std::string replace_in_string(const std::string &input, const std::string &from, const std::string &to)
Definition: utilities.cc:411
std::string dim_string(const int dim, const int spacedim)
Definition: utilities.cc:457
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)
Definition: utilities.cc:139
void destroy_communicator(Epetra_Comm &communicator)
Definition: utilities.cc:1129
Definition: cuda.h:31
std::string get_hostname()
Definition: utilities.cc:977
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1695
MPI_Comm duplicate_communicator(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:93
unsigned long int VmHWM
Definition: utilities.h:825
const Epetra_Comm & comm_world()
Definition: utilities.cc:1054
int string_to_int(const std::string &s)
Definition: utilities.cc:488
std::vector< unsigned int > reverse_permutation(const std::vector< unsigned int > &permutation)
Definition: utilities.cc:776
static ::ExceptionBase & ExcNotImplemented()
unsigned long int VmRSS
Definition: utilities.h:831
unsigned int get_this_mpi_process(const Epetra_Comm &mpi_communicator)
Definition: utilities.cc:1155
const std::string get_current_vectorization_level()
Definition: utilities.cc:917
Epetra_Map duplicate_map(const Epetra_BlockMap &map, const Epetra_Comm &comm)
Definition: utilities.cc:1163
std::string get_time()
Definition: utilities.cc:992
bool job_supports_mpi()
Definition: mpi.cc:802
unsigned long int VmPeak
Definition: utilities.h:815
unsigned int needed_digits(const unsigned int max_number)
Definition: utilities.cc:467
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcInvalidNumber2StringConversersion(unsigned int arg1, unsigned int arg2)