]> https://gitweb.dealii.org/ - dealii.git/commitdiff
new function Utilities::cutoff_digits() 8927/head
authorNiklas Fehn <fehn@lnm.mw.tum.de>
Fri, 18 Oct 2019 15:54:48 +0000 (17:54 +0200)
committerNiklas Fehn <fehn@lnm.mw.tum.de>
Thu, 24 Oct 2019 17:56:29 +0000 (19:56 +0200)
doc/news/changes/minor/20191018NiklasFehn-2 [new file with mode: 0644]
include/deal.II/base/utilities.h
source/base/utilities.cc
tests/base/utilities_20.cc [new file with mode: 0644]
tests/base/utilities_20.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20191018NiklasFehn-2 b/doc/news/changes/minor/20191018NiklasFehn-2
new file mode 100644 (file)
index 0000000..cbb30ee
--- /dev/null
@@ -0,0 +1,5 @@
+New: new function Utilities::truncate_to_n_digits() that can be used to
+cutoff floating point numbers after a specified number of digits of
+accuracy (decimal places) in scientific floating point notation.
+<br>
+(Niklas Fehn, 2019/10/18)
index 56b93386150d86a9af5562c43f1dcbf1b1045da9..bc5ae83a5f82a16bf420aa72c8eca649b26a014e 100644 (file)
@@ -174,10 +174,26 @@ namespace Utilities
   /**
    * Determine how many digits are needed to represent numbers at most as
    * large as the given number.
+   *
+   * @author Niklas Fehn, 2019
    */
   unsigned int
   needed_digits(const unsigned int max_number);
 
+  /**
+   * This function allows to cut off a floating point number @p number
+   * after @p n_digits of accuracy, i.e., after @p n_digits decimal places
+   * in scientific floating point notation. When interpreted as rounding
+   * operation, this function reduces the absolute value of a floating point
+   * number and always rounds towards zero, since decimal places are simply
+   * cut off.
+   *
+   * @author Niklas Fehn, 2019
+   */
+  template <typename Number>
+  Number
+  truncate_to_n_digits(const Number number, const unsigned int n_digits);
+
   /**
    * Given a string, convert it to an integer. Throw an assertion if that is
    * not possible.
index 1af3fe351fa86f2d1aa8b802af913b962120fcdc..0319a90167b3c4b45c96b6cd1df796d960d99296 100644 (file)
@@ -484,6 +484,34 @@ namespace Utilities
 
 
 
+  template <typename Number>
+  Number
+  truncate_to_n_digits(const Number number, const unsigned int n_digits)
+  {
+    AssertThrow(n_digits >= 1, ExcMessage("invalid parameter."));
+
+    if (!(std::fabs(number) > std::numeric_limits<Number>::min()))
+      return number;
+
+    const int order =
+      static_cast<int>(std::floor(std::log10(std::fabs(number))));
+
+    const int shift = -order + static_cast<int>(n_digits) - 1;
+
+    Assert(shift <= static_cast<int>(std::floor(
+                      std::log10(std::numeric_limits<Number>::max()))),
+           ExcMessage(
+             "Overflow. Use a smaller value for n_digits and/or make sure "
+             "that the absolute value of 'number' does not become too small."));
+
+    const Number factor = std::pow(10.0, static_cast<Number>(shift));
+
+    const Number number_cutoff = std::trunc(number * factor) / factor;
+
+    return number_cutoff;
+  }
+
+
   int
   string_to_int(const std::string &s_)
   {
@@ -1228,6 +1256,11 @@ namespace Utilities
   template std::string
   to_string<long double>(long double, unsigned int);
 
+  template double
+  truncate_to_n_digits(const double, const unsigned int);
+  template float
+  truncate_to_n_digits(const float, const unsigned int);
+
   template std::vector<std::array<std::uint64_t, 1>>
   inverse_Hilbert_space_filling_curve<1, double>(
     const std::vector<Point<1, double>> &,
diff --git a/tests/base/utilities_20.cc b/tests/base/utilities_20.cc
new file mode 100644 (file)
index 0000000..74deb90
--- /dev/null
@@ -0,0 +1,92 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2019 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+// test functions in namespace Utilities
+
+#include <deal.II/base/utilities.h>
+
+#include "../tests.h"
+
+
+void
+test()
+{
+  double number        = 1.23456789e5;
+  double number_1digit = Utilities::truncate_to_n_digits(number, 1);
+  double number_2digit = Utilities::truncate_to_n_digits(number, 2);
+  double number_3digit = Utilities::truncate_to_n_digits(number, 3);
+
+  deallog << std::endl;
+  deallog << "number = " << std::scientific << number << std::endl;
+  deallog << "number_1digit = " << std::scientific << number_1digit
+          << std::endl;
+  deallog << "number_2digit = " << std::scientific << number_2digit
+          << std::endl;
+  deallog << "number_3digit = " << std::scientific << number_3digit
+          << std::endl;
+
+  number        = 0.999999999;
+  number_1digit = Utilities::truncate_to_n_digits(number, 1);
+  number_2digit = Utilities::truncate_to_n_digits(number, 2);
+  number_3digit = Utilities::truncate_to_n_digits(number, 3);
+
+  deallog << std::endl;
+  deallog << "number = " << std::scientific << number << std::endl;
+  deallog << "number_1digit = " << std::scientific << number_1digit
+          << std::endl;
+  deallog << "number_2digit = " << std::scientific << number_2digit
+          << std::endl;
+  deallog << "number_3digit = " << std::scientific << number_3digit
+          << std::endl;
+
+  number        = 0.0;
+  number_1digit = Utilities::truncate_to_n_digits(number, 1);
+  number_2digit = Utilities::truncate_to_n_digits(number, 2);
+  number_3digit = Utilities::truncate_to_n_digits(number, 3);
+
+  deallog << std::endl;
+  deallog << "number = " << std::scientific << number << std::endl;
+  deallog << "number_1digit = " << std::scientific << number_1digit
+          << std::endl;
+  deallog << "number_2digit = " << std::scientific << number_2digit
+          << std::endl;
+  deallog << "number_3digit = " << std::scientific << number_3digit
+          << std::endl;
+
+  float number_float        = -9.87654321e-6;
+  float number_1digit_float = Utilities::truncate_to_n_digits(number_float, 1);
+  float number_2digit_float = Utilities::truncate_to_n_digits(number_float, 2);
+  float number_3digit_float = Utilities::truncate_to_n_digits(number_float, 3);
+
+  deallog << std::endl;
+  deallog << "number = " << std::scientific << number_float << std::endl;
+  deallog << "number_1digit = " << std::scientific << number_1digit_float
+          << std::endl;
+  deallog << "number_2digit = " << std::scientific << number_2digit_float
+          << std::endl;
+  deallog << "number_3digit = " << std::scientific << number_3digit_float
+          << std::endl;
+}
+
+
+
+int
+main()
+{
+  initlog();
+
+  test();
+}
diff --git a/tests/base/utilities_20.output b/tests/base/utilities_20.output
new file mode 100644 (file)
index 0000000..afc9b26
--- /dev/null
@@ -0,0 +1,21 @@
+
+DEAL::
+DEAL::number = 1.234568e+05
+DEAL::number_1digit = 1.000000e+05
+DEAL::number_2digit = 1.200000e+05
+DEAL::number_3digit = 1.230000e+05
+DEAL::
+DEAL::number = 1.000000e+00
+DEAL::number_1digit = 9.000000e-01
+DEAL::number_2digit = 9.900000e-01
+DEAL::number_3digit = 9.990000e-01
+DEAL::
+DEAL::number = 0.000000e+00
+DEAL::number_1digit = 0.000000e+00
+DEAL::number_2digit = 0.000000e+00
+DEAL::number_3digit = 0.000000e+00
+DEAL::
+DEAL::number = -9.876543e-06
+DEAL::number_1digit = -9.000000e-06
+DEAL::number_2digit = -9.800000e-06
+DEAL::number_3digit = -9.870000e-06

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.