deal.II version GIT relicensing-2167-g9622207b8f 2024-11-21 12:40: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
polynomials_rannacher_turek.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2015 - 2024 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
18
19#include <memory>
20
22
23
24template <int dim>
30
31
32
33template <int dim>
34double
36 const Point<dim> &p) const
37{
38 Assert(dim == 2, ExcNotImplemented());
39 if (i == 0)
40 {
41 return (0.75 - 2.5 * p[0] + 1.5 * p[1] +
42 1.5 * (p[0] * p[0] - p[1] * p[1]));
43 }
44 else if (i == 1)
45 {
46 return (-0.25 - 0.5 * p[0] + 1.5 * p[1] +
47 1.5 * (p[0] * p[0] - p[1] * p[1]));
48 }
49 else if (i == 2)
50 {
51 return (0.75 + 1.5 * p[0] - 2.5 * p[1] -
52 1.5 * (p[0] * p[0] - p[1] * p[1]));
53 }
54 else if (i == 3)
55 {
56 return (-0.25 + 1.5 * p[0] - 0.5 * p[1] -
57 1.5 * (p[0] * p[0] - p[1] * p[1]));
58 }
59
61 return 0;
62}
63
64
65
66template <int dim>
69 const Point<dim> &p) const
70{
71 if constexpr (dim == 2)
72 {
73 Tensor<1, dim> grad;
74 if (i == 0)
75 {
76 grad[0] = -2.5 + 3 * p[0];
77 grad[1] = 1.5 - 3 * p[1];
78 }
79 else if (i == 1)
80 {
81 grad[0] = -0.5 + 3.0 * p[0];
82 grad[1] = 1.5 - 3.0 * p[1];
83 }
84 else if (i == 2)
85 {
86 grad[0] = 1.5 - 3.0 * p[0];
87 grad[1] = -2.5 + 3.0 * p[1];
88 }
89 else if (i == 3)
90 {
91 grad[0] = 1.5 - 3.0 * p[0];
92 grad[1] = -0.5 + 3.0 * p[1];
93 }
94 else
95 {
97 }
98
99 return grad;
100 }
101
102 else
103 {
105 return {};
106 }
107}
108
109
110
111template <int dim>
114 const unsigned int i,
115 const Point<dim> & /*p*/) const
116{
117 if constexpr (dim == 2)
118 {
119 Tensor<2, dim> grad_grad;
120 if (i == 0)
121 {
122 grad_grad[0][0] = 3;
123 grad_grad[0][1] = 0;
124 grad_grad[1][0] = 0;
125 grad_grad[1][1] = -3;
126 }
127 else if (i == 1)
128 {
129 grad_grad[0][0] = 3;
130 grad_grad[0][1] = 0;
131 grad_grad[1][0] = 0;
132 grad_grad[1][1] = -3;
133 }
134 else if (i == 2)
135 {
136 grad_grad[0][0] = -3;
137 grad_grad[0][1] = 0;
138 grad_grad[1][0] = 0;
139 grad_grad[1][1] = 3;
140 }
141 else if (i == 3)
142 {
143 grad_grad[0][0] = -3;
144 grad_grad[0][1] = 0;
145 grad_grad[1][0] = 0;
146 grad_grad[1][1] = 3;
147 }
148 return grad_grad;
149 }
150
151 else
152 {
154 return {};
155 }
156}
157
158
159
160template <int dim>
161void
163 const Point<dim> &unit_point,
164 std::vector<double> &values,
165 std::vector<Tensor<1, dim>> &grads,
166 std::vector<Tensor<2, dim>> &grad_grads,
167 std::vector<Tensor<3, dim>> &third_derivatives,
168 std::vector<Tensor<4, dim>> &fourth_derivatives) const
169{
170 const unsigned int n_pols = this->n();
171 Assert(values.size() == n_pols || values.empty(),
172 ExcDimensionMismatch(values.size(), n_pols));
173 Assert(grads.size() == n_pols || grads.empty(),
174 ExcDimensionMismatch(grads.size(), n_pols));
175 Assert(grad_grads.size() == n_pols || grad_grads.empty(),
176 ExcDimensionMismatch(grad_grads.size(), n_pols));
177 Assert(third_derivatives.size() == n_pols || third_derivatives.empty(),
178 ExcDimensionMismatch(third_derivatives.size(), n_pols));
179 Assert(fourth_derivatives.size() == n_pols || fourth_derivatives.empty(),
180 ExcDimensionMismatch(fourth_derivatives.size(), n_pols));
181
182 for (unsigned int i = 0; i < n_pols; ++i)
183 {
184 if (values.size() != 0)
185 {
186 values[i] = compute_value(i, unit_point);
187 }
188 if (grads.size() != 0)
189 {
190 grads[i] = compute_grad(i, unit_point);
191 }
192 if (grad_grads.size() != 0)
193 {
194 grad_grads[i] = compute_grad_grad(i, unit_point);
195 }
196 if (third_derivatives.size() != 0)
197 {
198 third_derivatives[i] = compute_derivative<3>(i, unit_point);
199 }
200 if (fourth_derivatives.size() != 0)
201 {
202 fourth_derivatives[i] = compute_derivative<4>(i, unit_point);
203 }
204 }
205}
206
207
208
209template <int dim>
210std::unique_ptr<ScalarPolynomialsBase<dim>>
212{
213 return std::make_unique<PolynomialsRannacherTurek<dim>>(*this);
214}
215
216
217// explicit instantiations
218#include "polynomials_rannacher_turek.inst"
219
Definition point.h:111
virtual std::unique_ptr< ScalarPolynomialsBase< dim > > clone() const override
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< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const override
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const override
double compute_value(const unsigned int i, const Point< dim > &p) const override
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define DEAL_II_NOT_IMPLEMENTED()