--- /dev/null
+New: Add Functions::IncrementalFunction which returns the difference in values of
+a given function at two time steps.
+<br>
+(Denis Davydov, Jean-Paul Pelteret, 2018/10/05)
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2018 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.
+//
+// ---------------------------------------------------------------------
+
+#ifndef dealii_incremental_function_h
+#define dealii_incremental_function_h
+
+
+#include <deal.II/base/function.h>
+#include <deal.II/base/thread_management.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+template <typename number>
+class Vector;
+
+namespace Functions
+{
+ /**
+ * This class represents an incremental function. That is, given arbitrary
+ * function <code>func</code>, this class will return
+ * <code>f(t)-f(t-delta)</code>. The decrement is set by the method
+ * set_decrement(). The main application is to transform a given Dirichlet
+ * boundary condition function into the incremental form.
+ *
+ * @ingroup functions
+ * @author Denis Davydov, Jean-Paul Pelteret, 2018
+ */
+ template <int dim, typename RangeNumberType = double>
+ class IncrementalFunction : public Function<dim, RangeNumberType>
+ {
+ public:
+ /**
+ * Export the value of the template parameter as a static member constant.
+ * Sometimes useful for some expression template programming.
+ */
+ static const unsigned int dimension = dim;
+
+ /**
+ * The scalar-valued real type used for representing time.
+ */
+ using time_type = typename Function<dim, RangeNumberType>::time_type;
+
+ /**
+ * Constructor which wraps a given function @p base.
+ *
+ * @note this class stores a non-constant reference to @p base
+ * and will call <code>base.set_time()</code> during evaluation.
+ */
+ IncrementalFunction(Function<dim, RangeNumberType> &base);
+
+ /**
+ * Virtual destructor
+ */
+ virtual ~IncrementalFunction() = default;
+
+ /**
+ * Return the value of the function at the given point. Unless there is only
+ * one component (i.e. the function is scalar), you should state the
+ * component you want to have evaluated; it defaults to zero, i.e. the first
+ * component.
+ */
+ virtual RangeNumberType
+ value(const Point<dim> &p, const unsigned int component = 0) const override;
+
+ /**
+ * Return all components of a vector-valued function at a given point.
+ *
+ * <tt>values</tt> shall have the right size beforehand.
+ */
+ virtual void
+ vector_value(const Point<dim> & p,
+ Vector<RangeNumberType> &values) const override;
+
+ /**
+ * Set (positive) time decrement.
+ */
+ void
+ set_decrement(const time_type delta_t);
+
+ private:
+ /**
+ * Reference to the function being wrapped.
+ */
+ Function<dim, RangeNumberType> &base;
+
+ /**
+ * Time decrement.
+ */
+ time_type delta_t;
+
+ /**
+ * Auxiliary vector to store values
+ */
+ mutable Vector<RangeNumberType> values_old;
+
+ /**
+ * Thread mutex
+ */
+ mutable Threads::Mutex mutex;
+ };
+
+} // namespace Functions
+
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
geometry_info.cc
geometric_utilities.cc
graph_coloring.cc
+ incremental_function.cc
index_set.cc
job_identifier.cc
logstream.cc
data_out_base.inst.in
function.inst.in
function_time.inst.in
+ incremental_function.inst.in
geometric_utilities.inst.in
mpi.inst.in
partitioner.inst.in
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2007 - 2018 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.
+//
+// ---------------------------------------------------------------------
+
+#include <deal.II/base/incremental_function.h>
+
+#include <deal.II/lac/vector.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace Functions
+{
+ template <int dim, typename RangeNumberType>
+ IncrementalFunction<dim, RangeNumberType>::IncrementalFunction(
+ Function<dim, RangeNumberType> &base)
+ : base(base)
+ , delta_t(time_type())
+ , values_old(base.n_components)
+ {}
+
+
+
+ template <int dim, typename RangeNumberType>
+ void
+ IncrementalFunction<dim, RangeNumberType>::set_decrement(
+ const time_type delta)
+ {
+ delta_t = delta;
+ }
+
+
+
+ template <int dim, typename RangeNumberType>
+ void
+ IncrementalFunction<dim, RangeNumberType>::vector_value(
+ const Point<dim> & p,
+ Vector<RangeNumberType> &values) const
+ {
+ // since we modify a mutable member variable, lock the
+ // the data via a mutex
+ std::lock_guard<std::mutex> lock(mutex);
+
+ base.set_time(this->get_time());
+ base.vector_value(p, values);
+
+ base.set_time(this->get_time() - delta_t);
+ base.vector_value(p, values_old);
+
+ values -= values_old;
+ }
+
+
+
+ template <int dim, typename RangeNumberType>
+ RangeNumberType
+ IncrementalFunction<dim, RangeNumberType>::value(const Point<dim> &p,
+ unsigned int comp) const
+ {
+ base.set_time(this->get_time());
+ const RangeNumberType current = base.value(p, comp);
+
+ base.set_time(this->get_time() - delta_t);
+ const RangeNumberType old = base.value(p, comp);
+
+ return current - old;
+ }
+
+
+// Explicit instantiations
+#include "incremental_function.inst"
+} // namespace Functions
+DEAL_II_NAMESPACE_CLOSE
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2015 - 2017 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.
+//
+// ---------------------------------------------------------------------
+
+
+for (S : REAL_SCALARS; dim : SPACE_DIMENSIONS)
+ {
+ template class Functions::IncrementalFunction<dim, S>;
+ }
+
+for (S : COMPLEX_SCALARS; dim : SPACE_DIMENSIONS)
+ {
+ template class Functions::IncrementalFunction<dim, S>;
+ }
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 - 2017 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 IncrementalFunction
+
+#include <deal.II/base/incremental_function.h>
+
+#include "../tests.h"
+
+// homogeneous linear in time function
+template <int dim>
+class MyFunc : public Function<dim>
+{
+public:
+ MyFunc(const double a, const double b)
+ : Function<dim>(dim)
+ , a(a)
+ , b(b)
+ {}
+
+ virtual double
+ value(const Point<dim> /*p*/, const unsigned int component = 0) const
+ {
+ (void)component;
+ return a * this->get_time() + b;
+ }
+
+ virtual void
+ vector_value(const Point<dim> & /*p*/, Vector<double> &values) const override
+ {
+ for (unsigned int i = 0; i < dim; ++i)
+ values[i] = a * this->get_time() + b;
+ }
+
+private:
+ const double a;
+ const double b;
+};
+
+template <int dim>
+void
+check()
+{
+ Vector<double> v1(dim);
+ Vector<double> v2(dim);
+ Vector<double> v3(dim);
+
+ const Point<dim> point;
+
+ MyFunc<dim> func(0.5, 10);
+ Functions::IncrementalFunction<dim> inc(func);
+
+ // default delta is zero:
+ inc.set_time(1.0);
+
+ std::vector<std::pair<double, double>> time_delta = {{1.0, 0.}, {1.0, 0.2}};
+
+ for (auto &el : time_delta)
+ {
+ deallog << "time = " << el.first << " delta = " << el.second << std::endl;
+ func.set_time(el.first);
+ func.vector_value(point, v1);
+ v1.print(deallog);
+ deallog << "-" << std::endl;
+
+ func.set_time(el.first - el.second);
+ func.vector_value(point, v2);
+ v2.print(deallog);
+ deallog << "=" << std::endl;
+
+ inc.set_time(el.first);
+ inc.set_decrement(el.second);
+ inc.vector_value(point, v3);
+ v3.print(deallog);
+
+ v1.add(-1, v2);
+ v1.add(-1, v3);
+ AssertThrow(v1.l2_norm() == 0., ExcInternalError());
+ }
+
+ deallog << "Ok" << std::endl;
+}
+
+int
+main()
+{
+ initlog();
+
+ check<2>();
+}
--- /dev/null
+
+DEAL::time = 1.00000 delta = 0.00000
+DEAL::10.5000 10.5000
+DEAL::-
+DEAL::10.5000 10.5000
+DEAL::=
+DEAL::0.00000 0.00000
+DEAL::time = 1.00000 delta = 0.200000
+DEAL::10.5000 10.5000
+DEAL::-
+DEAL::10.4000 10.4000
+DEAL::=
+DEAL::0.100000 0.100000
+DEAL::Ok