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
polynomials_wedge.h
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
15
16#ifndef dealii_base_polynomials_wedge_h
17#define dealii_base_polynomials_wedge_h
18
19#include <deal.II/base/config.h>
20
24
26
27
28namespace internal
29{
35 static const constexpr ::ndarray<unsigned int, 6, 2> wedge_table_1{
36 {{{0, 0}}, {{1, 0}}, {{2, 0}}, {{0, 1}}, {{1, 1}}, {{2, 1}}}};
37
43 static const constexpr ::ndarray<unsigned int, 18, 2> wedge_table_2{
44 {{{0, 0}},
45 {{1, 0}},
46 {{2, 0}},
47 {{0, 1}},
48 {{1, 1}},
49 {{2, 1}},
50 {{3, 0}},
51 {{4, 0}},
52 {{5, 0}},
53 {{3, 1}},
54 {{4, 1}},
55 {{5, 1}},
56 {{0, 2}},
57 {{1, 2}},
58 {{2, 2}},
59 {{3, 2}},
60 {{4, 2}},
61 {{5, 2}}}};
62} // namespace internal
63
64
74template <int dim>
76{
77public:
81 static constexpr unsigned int dimension = dim;
82
83 /*
84 * Constructor taking the polynomial @p degree as input.
85 *
86 * @note Currently, only linear (degree=1) and quadratic polynomials
87 * (degree=2) are implemented.
88 */
89 ScalarLagrangePolynomialWedge(const unsigned int degree);
90
96 void
97 evaluate(const Point<dim> &unit_point,
98 std::vector<double> &values,
99 std::vector<Tensor<1, dim>> &grads,
100 std::vector<Tensor<2, dim>> &grad_grads,
101 std::vector<Tensor<3, dim>> &third_derivatives,
102 std::vector<Tensor<4, dim>> &fourth_derivatives) const override;
103
104 double
105 compute_value(const unsigned int i, const Point<dim> &p) const override;
106
112 template <int order>
114 compute_derivative(const unsigned int i, const Point<dim> &p) const;
115
117 compute_1st_derivative(const unsigned int i,
118 const Point<dim> &p) const override;
119
126 compute_2nd_derivative(const unsigned int i,
127 const Point<dim> &p) const override;
128
135 compute_3rd_derivative(const unsigned int i,
136 const Point<dim> &p) const override;
137
144 compute_4th_derivative(const unsigned int i,
145 const Point<dim> &p) const override;
146
153 compute_grad(const unsigned int i, const Point<dim> &p) const override;
154
161 compute_grad_grad(const unsigned int i, const Point<dim> &p) const override;
162
163 std::string
164 name() const override;
165
166 virtual std::unique_ptr<ScalarPolynomialsBase<dim>>
167 clone() const override;
168
169private:
174
179};
180
181
182
183template <int dim>
184template <int order>
187 const unsigned int i,
188 const Point<dim> &p) const
189{
191
192 AssertDimension(order, 1);
193 const auto grad = compute_grad(i, p);
194
195 for (unsigned int i = 0; i < dim; ++i)
196 der[i] = grad[i];
197
198 return der;
199}
200
202
203#endif
Definition point.h:111
void evaluate(const Point< dim > &unit_point, std::vector< double > &values, std::vector< Tensor< 1, dim > > &grads, std::vector< Tensor< 2, dim > > &grad_grads, std::vector< Tensor< 3, dim > > &third_derivatives, std::vector< Tensor< 4, dim > > &fourth_derivatives) const override
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const override
Tensor< 1, dim > compute_1st_derivative(const unsigned int i, const Point< dim > &p) const override
Tensor< 4, dim > compute_4th_derivative(const unsigned int i, const Point< dim > &p) const override
std::string name() const override
Tensor< 3, dim > compute_3rd_derivative(const unsigned int i, const Point< dim > &p) const override
Tensor< 2, dim > compute_2nd_derivative(const unsigned int i, const Point< dim > &p) const override
Tensor< order, dim > compute_derivative(const unsigned int i, const Point< dim > &p) const
virtual std::unique_ptr< ScalarPolynomialsBase< dim > > clone() const override
static constexpr unsigned int dimension
double compute_value(const unsigned int i, const Point< dim > &p) const override
const BarycentricPolynomials< 1 > poly_line
const BarycentricPolynomials< 2 > poly_tri
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const override
virtual unsigned int degree() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
#define AssertDimension(dim1, dim2)
static const constexpr ::ndarray< unsigned int, 18, 2 > wedge_table_2
static const constexpr ::ndarray< unsigned int, 6, 2 > wedge_table_1