1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203
|
// ------------------------------------------------------------------------
//
// SPDX-License-Identifier: LGPL-2.1-or-later
// Copyright (C) 2021 - 2025 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.
//
// ------------------------------------------------------------------------
#include <deal.II/base/bounding_box.h>
#include <deal.II/base/function_restriction.h>
DEAL_II_NAMESPACE_OPEN
namespace internal
{
template <int dim>
Point<dim + 1>
create_higher_dim_point(const Point<dim> &point,
const unsigned int component_in_dim_plus_1,
const double coordinate_value)
{
AssertIndexRange(component_in_dim_plus_1, dim + 1);
Point<dim + 1> output;
output[component_in_dim_plus_1] = coordinate_value;
for (int d = 0; d < dim; ++d)
{
const unsigned int component_to_write_to =
dealii::internal::coordinate_to_one_dim_higher<dim>(
component_in_dim_plus_1, d);
output[component_to_write_to] = point[d];
}
return output;
}
} // namespace internal
namespace Functions
{
template <int dim>
CoordinateRestriction<dim>::CoordinateRestriction(
const Function<dim + 1> &function,
const unsigned int direction,
const double coordinate_value)
: function(&function)
, restricted_direction(direction)
, coordinate_value(coordinate_value)
{
AssertIndexRange(restricted_direction, dim + 1);
}
template <int dim>
double
CoordinateRestriction<dim>::value(const Point<dim> &point,
const unsigned int component) const
{
const Point<dim + 1> full_point =
internal::create_higher_dim_point(point,
restricted_direction,
coordinate_value);
return function->value(full_point, component);
}
template <int dim>
Tensor<1, dim>
CoordinateRestriction<dim>::gradient(const Point<dim> &point,
const unsigned int component) const
{
const Point<dim + 1> full_point =
internal::create_higher_dim_point(point,
restricted_direction,
coordinate_value);
const Tensor<1, dim + 1> full_gradient =
function->gradient(full_point, component);
// CoordinateRestriction is constant in restricted direction. Through away
// the derivatives with respect to this direction and copy the other
// values.
Tensor<1, dim> grad;
for (unsigned int d = 0; d < dim; ++d)
{
const unsigned int index_to_write_from =
internal::coordinate_to_one_dim_higher<dim>(restricted_direction, d);
grad[d] = full_gradient[index_to_write_from];
}
return grad;
}
template <int dim>
SymmetricTensor<2, dim>
CoordinateRestriction<dim>::hessian(const Point<dim> &point,
const unsigned int component) const
{
const Point<dim + 1> full_point =
internal::create_higher_dim_point(point,
restricted_direction,
coordinate_value);
const Tensor<2, dim + 1> full_hessian =
function->hessian(full_point, component);
// CoordinateRestriction is constant in restricted direction. Through away
// the derivatives with respect to this direction and copy the other
// values.
SymmetricTensor<2, dim> hess;
for (unsigned int i = 0; i < dim; ++i)
{
const unsigned int i_to_write_from =
internal::coordinate_to_one_dim_higher<dim>(restricted_direction, i);
for (unsigned int j = 0; j < dim; ++j)
{
const unsigned int j_to_write_from =
internal::coordinate_to_one_dim_higher<dim>(restricted_direction,
j);
hess[i][j] = full_hessian[i_to_write_from][j_to_write_from];
}
}
return hess;
}
template <int dim>
PointRestriction<dim>::PointRestriction(const Function<dim + 1> &function,
const unsigned int open_direction,
const Point<dim> &point)
: function(&function)
, open_direction(open_direction)
, point(point)
{
AssertIndexRange(open_direction, dim + 1);
}
template <int dim>
double
PointRestriction<dim>::value(const Point<1> &point_1D,
const unsigned int component) const
{
const Point<dim + 1> full_point =
internal::create_higher_dim_point(point, open_direction, point_1D[0]);
return function->value(full_point, component);
}
template <int dim>
Tensor<1, 1>
PointRestriction<dim>::gradient(const Point<1> &point_1D,
const unsigned int component) const
{
const Point<dim + 1> full_point =
internal::create_higher_dim_point(point, open_direction, point_1D[0]);
const Tensor<1, dim + 1> full_gradient =
function->gradient(full_point, component);
// The PointRestrictions is constant in all but the open direction. Throw
// away the derivatives in all but this direction.
Tensor<1, 1> grad;
grad[0] = full_gradient[open_direction];
return grad;
}
template <int dim>
SymmetricTensor<2, 1>
PointRestriction<dim>::hessian(const Point<1> &point_1D,
const unsigned int component) const
{
const Point<dim + 1> full_point =
internal::create_higher_dim_point(point, open_direction, point_1D[0]);
const Tensor<2, dim + 1> full_hessian =
function->hessian(full_point, component);
// The PointRestrictions is constant in all but the open direction. Throw
// away the derivatives in all but this direction.
SymmetricTensor<2, 1> hess;
hess[0][0] = full_hessian[open_direction][open_direction];
return hess;
}
} // namespace Functions
#include "base/function_restriction.inst"
DEAL_II_NAMESPACE_CLOSE
|