deal.II version GIT relicensing-1989-gd7a2c90e4e 2024-10-14 01:50:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
auto_derivative_function.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2001 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_auto_derivative_function_h
16#define dealii_auto_derivative_function_h
17
18
19#include <deal.II/base/config.h>
20
23
25
80template <int dim>
82{
83public:
117
132 AutoDerivativeFunction(const double h,
133 const unsigned int n_components = 1,
134 const double initial_time = 0.0);
135
139 virtual ~AutoDerivativeFunction() override = default;
140
145 void
147
156 void
157 set_h(const double h);
158
166 virtual Tensor<1, dim>
167 gradient(const Point<dim> &p,
168 const unsigned int component = 0) const override;
169
176 virtual void
178 std::vector<Tensor<1, dim>> &gradients) const override;
179
189 virtual void
190 gradient_list(const std::vector<Point<dim>> &points,
191 std::vector<Tensor<1, dim>> &gradients,
192 const unsigned int component = 0) const override;
193
206 virtual void
208 const std::vector<Point<dim>> &points,
209 std::vector<std::vector<Tensor<1, dim>>> &gradients) const override;
210
214 static DifferenceFormula
215 get_formula_of_order(const unsigned int ord);
216
217
218private:
222 double h;
223
227 std::vector<Tensor<1, dim>> ht;
228
233};
234
235
237
238#endif
virtual ~AutoDerivativeFunction() override=default
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual void vector_gradient(const Point< dim > &p, std::vector< Tensor< 1, dim > > &gradients) const override
virtual void vector_gradient_list(const std::vector< Point< dim > > &points, std::vector< std::vector< Tensor< 1, dim > > > &gradients) const override
std::vector< Tensor< 1, dim > > ht
void set_formula(const DifferenceFormula formula=Euler)
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component=0) const override
static DifferenceFormula get_formula_of_order(const unsigned int ord)
const unsigned int n_components
Definition function.h:163
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499