Reference documentation for deal.II version 9.3.3
\(\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\}}\)
polynomials_rannacher_turek.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2015 - 2020 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16
19
20#include <memory>
21
23
24
25template <int dim>
27 : ScalarPolynomialsBase<dim>(2, ::GeometryInfo<dim>::faces_per_cell)
28{
29 Assert(dim == 2, ExcNotImplemented());
30}
31
32
33
34template <int dim>
35double
37 const Point<dim> & p) const
38{
39 Assert(dim == 2, ExcNotImplemented());
40 if (i == 0)
41 {
42 return (0.75 - 2.5 * p(0) + 1.5 * p(1) +
43 1.5 * (p(0) * p(0) - p(1) * p(1)));
44 }
45 else if (i == 1)
46 {
47 return (-0.25 - 0.5 * p(0) + 1.5 * p(1) +
48 1.5 * (p(0) * p(0) - p(1) * p(1)));
49 }
50 else if (i == 2)
51 {
52 return (0.75 + 1.5 * p(0) - 2.5 * p(1) -
53 1.5 * (p(0) * p(0) - p(1) * p(1)));
54 }
55 else if (i == 3)
56 {
57 return (-0.25 + 1.5 * p(0) - 0.5 * p(1) -
58 1.5 * (p(0) * p(0) - p(1) * p(1)));
59 }
60
61 Assert(false, ExcNotImplemented());
62 return 0;
63}
64
65
66
67template <int dim>
70 const Point<dim> & p) const
71{
72 Assert(dim == 2, ExcNotImplemented());
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 {
96 Assert(false, ExcNotImplemented());
97 }
98
99 return grad;
100}
101
102
103
104template <int dim>
107 const unsigned int i,
108 const Point<dim> & /*p*/) const
109{
110 Assert(dim == 2, ExcNotImplemented());
111 Tensor<2, dim> grad_grad;
112 if (i == 0)
113 {
114 grad_grad[0][0] = 3;
115 grad_grad[0][1] = 0;
116 grad_grad[1][0] = 0;
117 grad_grad[1][1] = -3;
118 }
119 else if (i == 1)
120 {
121 grad_grad[0][0] = 3;
122 grad_grad[0][1] = 0;
123 grad_grad[1][0] = 0;
124 grad_grad[1][1] = -3;
125 }
126 else if (i == 2)
127 {
128 grad_grad[0][0] = -3;
129 grad_grad[0][1] = 0;
130 grad_grad[1][0] = 0;
131 grad_grad[1][1] = 3;
132 }
133 else if (i == 3)
134 {
135 grad_grad[0][0] = -3;
136 grad_grad[0][1] = 0;
137 grad_grad[1][0] = 0;
138 grad_grad[1][1] = 3;
139 }
140 return grad_grad;
141}
142
143
144
145template <int dim>
146void
148 const Point<dim> & unit_point,
149 std::vector<double> & values,
150 std::vector<Tensor<1, dim>> &grads,
151 std::vector<Tensor<2, dim>> &grad_grads,
152 std::vector<Tensor<3, dim>> &third_derivatives,
153 std::vector<Tensor<4, dim>> &fourth_derivatives) const
154{
155 const unsigned int n_pols = this->n();
156 Assert(values.size() == n_pols || values.size() == 0,
157 ExcDimensionMismatch(values.size(), n_pols));
158 Assert(grads.size() == n_pols || grads.size() == 0,
159 ExcDimensionMismatch(grads.size(), n_pols));
160 Assert(grad_grads.size() == n_pols || grad_grads.size() == 0,
161 ExcDimensionMismatch(grad_grads.size(), n_pols));
162 Assert(third_derivatives.size() == n_pols || third_derivatives.size() == 0,
163 ExcDimensionMismatch(third_derivatives.size(), n_pols));
164 Assert(fourth_derivatives.size() == n_pols || fourth_derivatives.size() == 0,
165 ExcDimensionMismatch(fourth_derivatives.size(), n_pols));
166
167 for (unsigned int i = 0; i < n_pols; ++i)
168 {
169 if (values.size() != 0)
170 {
171 values[i] = compute_value(i, unit_point);
172 }
173 if (grads.size() != 0)
174 {
175 grads[i] = compute_grad(i, unit_point);
176 }
177 if (grad_grads.size() != 0)
178 {
179 grad_grads[i] = compute_grad_grad(i, unit_point);
180 }
181 if (third_derivatives.size() != 0)
182 {
183 third_derivatives[i] = compute_derivative<3>(i, unit_point);
184 }
185 if (fourth_derivatives.size() != 0)
186 {
187 fourth_derivatives[i] = compute_derivative<4>(i, unit_point);
188 }
189 }
190}
191
192
193
194template <int dim>
195std::unique_ptr<ScalarPolynomialsBase<dim>>
197{
198 return std::make_unique<PolynomialsRannacherTurek<dim>>(*this);
199}
200
201
202// explicit instantiations
203#include "polynomials_rannacher_turek.inst"
204
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:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)