Reference documentation for deal.II version GIT relicensing-437-g81ec864850 2024-04-19 07:30:02+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
function_tools.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2021 - 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
16
17#include <cmath>
18
19
21namespace FunctionTools
22{
23 template <int dim>
24 void
26 const Function<dim> &function,
27 const BoundingBox<dim> &box,
28 std::pair<double, double> &value_bounds,
29 std::array<std::pair<double, double>, dim> &gradient_bounds,
30 const unsigned int component)
31 {
32 const Point<dim> center = box.center();
33 const double value = function.value(center, component);
34 const Tensor<1, dim> gradient = function.gradient(center, component);
35 const SymmetricTensor<2, dim> hessian = function.hessian(center, component);
36
37 // Deviation from function value at the center, based on the
38 // Taylor-expansion: |f'| * dx + 1/2 * |f''| * dx^2, (in 1d). dx is half
39 // the side-length of the box.
40 double taylor_bound_f = 0;
41
42 for (unsigned int i = 0; i < dim; ++i)
43 {
44 const double dx_i = .5 * box.side_length(i);
45
46 taylor_bound_f += std::abs(gradient[i]) * dx_i;
47
48 // Deviation from value of df/dx_i at the center,
49 // |f''| * dx, (in 1d).
50 double taylor_bound_dfdxi = 0;
51
52 for (unsigned int j = 0; j < dim; ++j)
53 {
54 const double dx_j = .5 * box.side_length(j);
55
56 taylor_bound_dfdxi += std::abs(hessian[i][j]) * dx_j;
57 taylor_bound_f += .5 * std::abs(hessian[i][j]) * dx_i * dx_j;
58 }
59
60 gradient_bounds[i].first = gradient[i] - taylor_bound_dfdxi;
61 gradient_bounds[i].second = gradient[i] + taylor_bound_dfdxi;
62 }
63
64 value_bounds.first = value - taylor_bound_f;
65 value_bounds.second = value + taylor_bound_f;
66 }
67
68} // namespace FunctionTools
69
70#include "function_tools.inst"
71
Point< spacedim, Number > center() const
Number side_length(const unsigned int direction) const
virtual SymmetricTensor< 2, dim, RangeNumberType > hessian(const Point< dim > &p, const unsigned int component=0) const
virtual Tensor< 1, dim, RangeNumberType > gradient(const Point< dim > &p, const unsigned int component=0) const
virtual RangeNumberType value(const Point< dim > &p, const unsigned int component=0) const
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
Point< 3 > center
void taylor_estimate_function_bounds(const Function< dim > &function, const BoundingBox< dim > &box, std::pair< double, double > &value_bounds, std::array< std::pair< double, double >, dim > &gradient_bounds, const unsigned int component=0)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)