]> https://gitweb.dealii.org/ - dealii.git/commitdiff
ScalarFunctionFromFunctionObject for time-dependent functions 17731/head
authorPeter Munch <peterrmuench@gmail.com>
Thu, 26 Sep 2024 16:09:29 +0000 (18:09 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Thu, 26 Sep 2024 17:47:32 +0000 (19:47 +0200)
doc/news/changes/minor/20240926Munch [new file with mode: 0644]
include/deal.II/base/function.h
include/deal.II/base/function.templates.h
tests/base/functions_15.cc [new file with mode: 0644]
tests/base/functions_15.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20240926Munch b/doc/news/changes/minor/20240926Munch
new file mode 100644 (file)
index 0000000..5f243ad
--- /dev/null
@@ -0,0 +1,4 @@
+Improved: The class ScalarFunctionFromFunctionObject can now
+also handle time-dependent functions.
+<br>
+(Peter Munch, 2024/09/26)
index 36cd270d6ed43e6fc20a2abd7b05dc015481176b..cb5b913238ab4972e505db3b02440e01407c685b 100644 (file)
@@ -810,6 +810,15 @@ public:
   explicit ScalarFunctionFromFunctionObject(
     const std::function<RangeNumberType(const Point<dim> &)> &function_object);
 
+  /**
+   * Given a function object that takes  a time and a Point and returns a
+   * RangeNumberType value, convert this into an object that matches the
+   * Function<dim, RangeNumberType> interface.
+   */
+  explicit ScalarFunctionFromFunctionObject(
+    const std::function<RangeNumberType(const double, const Point<dim> &)>
+      &function_object_t);
+
   /**
    * Return the value of the function at the given point. Returns the value
    * the function given to the constructor produces for this point.
@@ -822,7 +831,8 @@ private:
    * The function object which we call when this class's value() or
    * value_list() functions are called.
    */
-  const std::function<RangeNumberType(const Point<dim> &)> function_object;
+  const std::function<RangeNumberType(const double, const Point<dim> &)>
+    function_object;
 };
 
 
index 7405c291ad4078807322831febb1ce65e5a75a40..1103012fca5768934d7c8d2904b3e72336fcb934 100644 (file)
@@ -691,7 +691,21 @@ ComponentSelectFunction<dim, RangeNumberType>::memory_consumption() const
 template <int dim, typename RangeNumberType>
 ScalarFunctionFromFunctionObject<dim, RangeNumberType>::
   ScalarFunctionFromFunctionObject(
-    const std::function<RangeNumberType(const Point<dim> &)> &function_object)
+    const std::function<RangeNumberType(const Point<dim> &)> &fu)
+  : ScalarFunctionFromFunctionObject<dim, RangeNumberType>(
+      [fu](const double t, const Point<dim> &x) {
+        (void)t; // we got a function object that only takes 'x', so ignore 't'
+        return fu(x);
+      })
+{}
+
+
+
+template <int dim, typename RangeNumberType>
+ScalarFunctionFromFunctionObject<dim, RangeNumberType>::
+  ScalarFunctionFromFunctionObject(
+    const std::function<RangeNumberType(const double, const Point<dim> &)>
+      &function_object)
   : Function<dim, RangeNumberType>(1)
   , function_object(function_object)
 {}
@@ -707,7 +721,8 @@ ScalarFunctionFromFunctionObject<dim, RangeNumberType>::value(
   (void)component;
   Assert(component == 0,
          ExcMessage("This object represents only scalar functions"));
-  return function_object(p);
+
+  return function_object(this->get_time(), p);
 }
 
 
diff --git a/tests/base/functions_15.cc b/tests/base/functions_15.cc
new file mode 100644 (file)
index 0000000..b7bf939
--- /dev/null
@@ -0,0 +1,79 @@
+// ------------------------------------------------------------------------
+//
+// SPDX-License-Identifier: LGPL-2.1-or-later
+// Copyright (C) 2024 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// Part of the source code is dual licensed under Apache-2.0 WITH
+// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
+// governing the source code and code contributions can be found in
+// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
+//
+// ------------------------------------------------------------------------
+
+
+// Test ScalarFunctionFromFunctionObject with functions with and without
+// time dependency.
+
+#include <deal.II/base/function.h>
+
+#include "../tests.h"
+
+
+template <int dim>
+void
+check1()
+{
+  ScalarFunctionFromFunctionObject<dim> object(
+    [](const auto &p) { return p.norm(); });
+
+  for (unsigned int i = 0; i < 10; ++i)
+    {
+      Point<dim> p;
+      for (unsigned int d = 0; d < dim; ++d)
+        p[d] = i + d;
+
+      AssertThrow(object.value(p) == p.norm(), ExcInternalError());
+    }
+
+  deallog << "OK" << std::endl;
+}
+
+
+template <int dim>
+void
+check2()
+{
+  ScalarFunctionFromFunctionObject<dim> object(
+    [](const auto t, const auto &p) { return p.norm() * t; });
+
+  for (unsigned int i = 0; i < 10; ++i)
+    {
+      object.set_time(i);
+
+      Point<dim> p;
+      for (unsigned int d = 0; d < dim; ++d)
+        p[d] = i + d;
+
+      AssertThrow(object.value(p) == (p.norm() * static_cast<double>(i)),
+                  ExcInternalError());
+    }
+
+  deallog << "OK" << std::endl;
+}
+
+
+int
+main()
+{
+  initlog();
+
+  check1<1>();
+  check1<2>();
+  check1<3>();
+
+  check2<1>();
+  check2<2>();
+  check2<3>();
+}
diff --git a/tests/base/functions_15.output b/tests/base/functions_15.output
new file mode 100644 (file)
index 0000000..0aa61ff
--- /dev/null
@@ -0,0 +1,7 @@
+
+DEAL::OK
+DEAL::OK
+DEAL::OK
+DEAL::OK
+DEAL::OK
+DEAL::OK

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.