deal.II version 9.7.0
\(\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
polynomials_vector_anisotropic.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2005 - 2025 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
19
20#include <iomanip>
21#include <iostream>
22#include <memory>
23
24
26
27
28namespace
29{
35 std::vector<std::vector<Point<1>>>
36 create_anisotropic_support_points(const unsigned int dim,
37 const unsigned int normal_degree,
38 const unsigned int tangential_degree)
39 {
40 Assert(dim > 0 && dim <= 3, ExcImpossibleInDim(dim));
41
42 std::vector<std::vector<Point<1>>> points_normal_tangential(dim);
43 points_normal_tangential[0] =
44 normal_degree == 0 ? QMidpoint<1>().get_points() :
45 QGaussLobatto<1>(normal_degree + 1).get_points();
46 for (unsigned int d = 1; d < dim; ++d)
47 {
48 points_normal_tangential[d] =
49 tangential_degree == 0 ?
50 QMidpoint<1>().get_points() :
51 QGaussLobatto<1>(tangential_degree + 1).get_points();
52 }
53 return points_normal_tangential;
54 }
55
56
57
58 // Create nodal anisotropic polynomials as the tensor product of Lagrange
59 // polynomials on Gauss-Lobatto points of the given degrees in the normal and
60 // tangential directions, respectively (we could also choose Lagrange
61 // polynomials on Gauss points but those are slightly more expensive to handle
62 // in classes).
63 std::vector<std::vector<Polynomials::Polynomial<double>>>
64 create_aniso_polynomials(const unsigned int dim,
65 const unsigned int normal_degree,
66 const unsigned int tangential_degree)
67 {
68 std::vector<std::vector<Point<1>>> points_aniso;
69 points_aniso =
70 create_anisotropic_support_points(dim, normal_degree, tangential_degree);
71 std::vector<std::vector<Polynomials::Polynomial<double>>> pols(dim);
72 pols[0] = Polynomials::generate_complete_Lagrange_basis(points_aniso[0]);
73 for (unsigned int d = 1; d < dim; ++d)
74 pols[d] = Polynomials::generate_complete_Lagrange_basis(points_aniso[d]);
75 return pols;
76 }
77} // namespace
78
79
80
81template <int dim>
83 const unsigned int normal_degree,
84 const unsigned int tangential_degree,
85 const std::vector<unsigned int> &polynomial_ordering)
91 create_aniso_polynomials(dim, normal_degree, tangential_degree))
92 , lexicographic_to_hierarchic(polynomial_ordering)
94 Utilities::invert_permutation(lexicographic_to_hierarchic))
95{
96 // create renumbering of the unknowns from the lexicographic order to the
97 // actual order required by the finite element class with unknowns on
98 // faces placed first
99 const unsigned int n_pols = polynomial_space.n();
100
101 // since we only store an anisotropic polynomial for the first component,
102 // we set up a second numbering to switch out the actual coordinate
103 // direction
104 renumber_aniso[0].resize(n_pols);
105 for (unsigned int i = 0; i < n_pols; ++i)
106 renumber_aniso[0][i] = i;
107 if (dim == 2)
108 {
109 // switch x and y component (i and j loops)
110 renumber_aniso[1].resize(n_pols);
111 for (unsigned int j = 0; j < normal_degree + 1; ++j)
112 for (unsigned int i = 0; i < tangential_degree + 1; ++i)
113 renumber_aniso[1][j * (tangential_degree + 1) + i] =
114 j + i * (normal_degree + 1);
115 }
116 if (dim == 3)
117 {
118 // switch x, y, and z component (i, j, k) -> (j, k, i)
119 renumber_aniso[1].resize(n_pols);
120 for (unsigned int k = 0; k < tangential_degree + 1; ++k)
121 for (unsigned int j = 0; j < normal_degree + 1; ++j)
122 for (unsigned int i = 0; i < tangential_degree + 1; ++i)
123 renumber_aniso[1][(k * (normal_degree + 1) + j) *
124 (tangential_degree + 1) +
125 i] =
126 j + (normal_degree + 1) * (k + i * (tangential_degree + 1));
127
128 // switch x, y, and z component (i, j, k) -> (k, i, j)
129 renumber_aniso[2].resize(n_pols);
130 for (unsigned int k = 0; k < normal_degree + 1; ++k)
131 for (unsigned int j = 0; j < tangential_degree + 1; ++j)
132 for (unsigned int i = 0; i < tangential_degree + 1; ++i)
133 renumber_aniso[2][(k * (tangential_degree + 1) + j) *
134 (tangential_degree + 1) +
135 i] =
136 k + (normal_degree + 1) * (i + j * (tangential_degree + 1));
137 }
138}
139
140
141
142template <int dim>
143void
145 const Point<dim> &unit_point,
146 std::vector<Tensor<1, dim>> &values,
147 std::vector<Tensor<2, dim>> &grads,
148 std::vector<Tensor<3, dim>> &grad_grads,
149 std::vector<Tensor<4, dim>> &third_derivatives,
150 std::vector<Tensor<5, dim>> &fourth_derivatives) const
151{
152 Assert(values.size() == this->n() || values.empty(),
153 ExcDimensionMismatch(values.size(), this->n()));
154 Assert(grads.size() == this->n() || grads.empty(),
155 ExcDimensionMismatch(grads.size(), this->n()));
156 Assert(grad_grads.size() == this->n() || grad_grads.empty(),
157 ExcDimensionMismatch(grad_grads.size(), this->n()));
158 Assert(third_derivatives.size() == this->n() || third_derivatives.empty(),
159 ExcDimensionMismatch(third_derivatives.size(), this->n()));
160 Assert(fourth_derivatives.size() == this->n() || fourth_derivatives.empty(),
161 ExcDimensionMismatch(fourth_derivatives.size(), this->n()));
162
163 std::vector<double> p_values;
164 std::vector<Tensor<1, dim>> p_grads;
165 std::vector<Tensor<2, dim>> p_grad_grads;
166 std::vector<Tensor<3, dim>> p_third_derivatives;
167 std::vector<Tensor<4, dim>> p_fourth_derivatives;
168
169 const unsigned int n_sub = polynomial_space.n();
170 p_values.resize((values.empty()) ? 0 : n_sub);
171 p_grads.resize((grads.empty()) ? 0 : n_sub);
172 p_grad_grads.resize((grad_grads.empty()) ? 0 : n_sub);
173 p_third_derivatives.resize((third_derivatives.empty()) ? 0 : n_sub);
174 p_fourth_derivatives.resize((fourth_derivatives.empty()) ? 0 : n_sub);
175
176 for (unsigned int d = 0; d < dim; ++d)
177 {
178 // First we copy the point. The polynomial space for component d
179 // consists of polynomials of degree k in x_d and degree k+1 in the
180 // other variables. in order to simplify this, we use the same
181 // AnisotropicPolynomial space and simply rotate the coordinates
182 // through all directions.
183 Point<dim> p;
184 for (unsigned int c = 0; c < dim; ++c)
185 p[c] = unit_point[(c + d) % dim];
186
187 polynomial_space.evaluate(p,
188 p_values,
189 p_grads,
190 p_grad_grads,
191 p_third_derivatives,
192 p_fourth_derivatives);
193
194 for (unsigned int i = 0; i < p_values.size(); ++i)
195 values[lexicographic_to_hierarchic[i + d * n_sub]][d] =
196 p_values[renumber_aniso[d][i]];
197
198 for (unsigned int i = 0; i < p_grads.size(); ++i)
199 for (unsigned int d1 = 0; d1 < dim; ++d1)
200 grads[lexicographic_to_hierarchic[i + d * n_sub]][d][(d1 + d) % dim] =
201 p_grads[renumber_aniso[d][i]][d1];
202
203 for (unsigned int i = 0; i < p_grad_grads.size(); ++i)
204 for (unsigned int d1 = 0; d1 < dim; ++d1)
205 for (unsigned int d2 = 0; d2 < dim; ++d2)
206 grad_grads[lexicographic_to_hierarchic[i + d * n_sub]][d]
207 [(d1 + d) % dim][(d2 + d) % dim] =
208 p_grad_grads[renumber_aniso[d][i]][d1][d2];
209
210 for (unsigned int i = 0; i < p_third_derivatives.size(); ++i)
211 for (unsigned int d1 = 0; d1 < dim; ++d1)
212 for (unsigned int d2 = 0; d2 < dim; ++d2)
213 for (unsigned int d3 = 0; d3 < dim; ++d3)
214 third_derivatives[lexicographic_to_hierarchic[i + d * n_sub]][d]
215 [(d1 + d) % dim][(d2 + d) % dim]
216 [(d3 + d) % dim] =
217 p_third_derivatives[renumber_aniso[d][i]][d1]
218 [d2][d3];
219
220 for (unsigned int i = 0; i < p_fourth_derivatives.size(); ++i)
221 for (unsigned int d1 = 0; d1 < dim; ++d1)
222 for (unsigned int d2 = 0; d2 < dim; ++d2)
223 for (unsigned int d3 = 0; d3 < dim; ++d3)
224 for (unsigned int d4 = 0; d4 < dim; ++d4)
225 fourth_derivatives[lexicographic_to_hierarchic[i + d * n_sub]]
226 [d][(d1 + d) % dim][(d2 + d) % dim]
227 [(d3 + d) % dim][(d4 + d) % dim] =
228 p_fourth_derivatives[renumber_aniso[d][i]]
229 [d1][d2][d3][d4];
230 }
231}
232
233
234
235template <int dim>
236std::string
238{
239 return "VectorAnisotropic";
240}
241
242
243
244template <int dim>
245unsigned int
247 const unsigned int normal_degree,
248 const unsigned int tangential_degree)
249{
250 return dim * (normal_degree + 1) *
251 Utilities::pow(tangential_degree + 1, dim - 1);
252}
253
254
255
256template <int dim>
257unsigned int
262
263
264
265template <int dim>
266unsigned int
271
272
273
274template <int dim>
275std::unique_ptr<TensorPolynomialsBase<dim>>
277{
278 return std::make_unique<PolynomialsVectorAnisotropic<dim>>(*this);
279}
280
281
282
283template <int dim>
284std::vector<Point<dim>>
286{
287 Assert(dim > 0 && dim <= 3, ExcImpossibleInDim(dim));
288 const std::vector<std::vector<Point<1>>> points_aniso =
289 create_anisotropic_support_points(dim, normal_degree, tangential_degree);
290 const unsigned int n_sub = polynomial_space.n();
291 std::vector<Point<dim>> points(dim * n_sub);
293 for (unsigned int d = 0; d < dim; ++d)
294 for (unsigned int i = 0; i < n_sub; ++i)
295 {
296 unsigned int renumbered_index = renumber_aniso[d][i];
297 std::array<unsigned int, dim> indices_points_anisotropic;
298 indices_points_anisotropic[0] = renumbered_index % (normal_degree + 1);
299 if (dim > 1)
300 {
301 renumbered_index /= (normal_degree + 1);
302 indices_points_anisotropic[1] =
303 renumbered_index % (tangential_degree + 1);
304 }
305 if (dim > 2)
306 indices_points_anisotropic[2] =
307 renumbered_index / (tangential_degree + 1);
308 for (unsigned int c = 0; c < dim; ++c)
309 {
310 points[lexicographic_to_hierarchic[i + d * n_sub]][(c + d) % dim] =
311 points_aniso[c][indices_points_anisotropic[c]][0];
312 }
313 }
314 return points;
315}
316
317
318
322
323
Definition point.h:113
std::vector< unsigned int > hierarchic_to_lexicographic
std::array< std::vector< unsigned int >, dim > renumber_aniso
PolynomialsVectorAnisotropic(const unsigned int degree_normal, const unsigned int degree_tangential, const std::vector< unsigned int > &polynomial_ordering)
const AnisotropicPolynomials< dim > polynomial_space
std::vector< Point< dim > > get_polynomial_support_points() const
virtual std::unique_ptr< TensorPolynomialsBase< dim > > clone() const override
void evaluate(const Point< dim > &unit_point, std::vector< Tensor< 1, dim > > &values, std::vector< Tensor< 2, dim > > &grads, std::vector< Tensor< 3, dim > > &grad_grads, std::vector< Tensor< 4, dim > > &third_derivatives, std::vector< Tensor< 5, dim > > &fourth_derivatives) const override
std::vector< unsigned int > lexicographic_to_hierarchic
static unsigned int n_polynomials(const unsigned int normal_degree, const unsigned int tangential_degree)
TensorPolynomialsBase(const unsigned int deg, const unsigned int n_polynomials)
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
std::vector< Polynomial< double > > generate_complete_Lagrange_basis(const std::vector< Point< 1 > > &points)
constexpr T pow(const T base, const int iexp)
Definition utilities.h:967
STL namespace.