--- /dev/null
+Improved: The class ScalarFunctionFromFunctionObject can now
+also handle time-dependent functions.
+<br>
+(Peter Munch, 2024/09/26)
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.
* 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;
};
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)
{}
(void)component;
Assert(component == 0,
ExcMessage("This object represents only scalar functions"));
- return function_object(p);
+
+ return function_object(this->get_time(), p);
}
--- /dev/null
+// ------------------------------------------------------------------------
+//
+// 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>();
+}
--- /dev/null
+
+DEAL::OK
+DEAL::OK
+DEAL::OK
+DEAL::OK
+DEAL::OK
+DEAL::OK