Reference documentation for deal.II version 9.2.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\}}\)
polynomials_rannacher_turek.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2015 - 2019 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 
20 
22 
23 
24 template <int dim>
26  : ScalarPolynomialsBase<dim>(2, ::GeometryInfo<dim>::faces_per_cell)
27 {
28  Assert(dim == 2, ExcNotImplemented());
29 }
30 
31 
32 
33 template <int dim>
34 double
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 
60  Assert(false, ExcNotImplemented());
61  return 0;
62 }
63 
64 
65 
66 template <int dim>
69  const Point<dim> & p) const
70 {
71  Assert(dim == 2, ExcNotImplemented());
72  Tensor<1, dim> grad;
73  if (i == 0)
74  {
75  grad[0] = -2.5 + 3 * p(0);
76  grad[1] = 1.5 - 3 * p(1);
77  }
78  else if (i == 1)
79  {
80  grad[0] = -0.5 + 3.0 * p(0);
81  grad[1] = 1.5 - 3.0 * p(1);
82  }
83  else if (i == 2)
84  {
85  grad[0] = 1.5 - 3.0 * p(0);
86  grad[1] = -2.5 + 3.0 * p(1);
87  }
88  else if (i == 3)
89  {
90  grad[0] = 1.5 - 3.0 * p(0);
91  grad[1] = -0.5 + 3.0 * p(1);
92  }
93  else
94  {
95  Assert(false, ExcNotImplemented());
96  }
97 
98  return grad;
99 }
100 
101 
102 
103 template <int dim>
106  const unsigned int i,
107  const Point<dim> & /*p*/) const
108 {
109  Assert(dim == 2, ExcNotImplemented());
110  Tensor<2, dim> grad_grad;
111  if (i == 0)
112  {
113  grad_grad[0][0] = 3;
114  grad_grad[0][1] = 0;
115  grad_grad[1][0] = 0;
116  grad_grad[1][1] = -3;
117  }
118  else if (i == 1)
119  {
120  grad_grad[0][0] = 3;
121  grad_grad[0][1] = 0;
122  grad_grad[1][0] = 0;
123  grad_grad[1][1] = -3;
124  }
125  else if (i == 2)
126  {
127  grad_grad[0][0] = -3;
128  grad_grad[0][1] = 0;
129  grad_grad[1][0] = 0;
130  grad_grad[1][1] = 3;
131  }
132  else if (i == 3)
133  {
134  grad_grad[0][0] = -3;
135  grad_grad[0][1] = 0;
136  grad_grad[1][0] = 0;
137  grad_grad[1][1] = 3;
138  }
139  return grad_grad;
140 }
141 
142 
143 
144 template <int dim>
145 void
147  const Point<dim> & unit_point,
148  std::vector<double> & values,
149  std::vector<Tensor<1, dim>> &grads,
150  std::vector<Tensor<2, dim>> &grad_grads,
151  std::vector<Tensor<3, dim>> &third_derivatives,
152  std::vector<Tensor<4, dim>> &fourth_derivatives) const
153 {
154  const unsigned int n_pols = this->n();
155  Assert(values.size() == n_pols || values.size() == 0,
156  ExcDimensionMismatch(values.size(), n_pols));
157  Assert(grads.size() == n_pols || grads.size() == 0,
158  ExcDimensionMismatch(grads.size(), n_pols));
159  Assert(grad_grads.size() == n_pols || grad_grads.size() == 0,
160  ExcDimensionMismatch(grad_grads.size(), n_pols));
161  Assert(third_derivatives.size() == n_pols || third_derivatives.size() == 0,
162  ExcDimensionMismatch(third_derivatives.size(), n_pols));
163  Assert(fourth_derivatives.size() == n_pols || fourth_derivatives.size() == 0,
164  ExcDimensionMismatch(fourth_derivatives.size(), n_pols));
165 
166  for (unsigned int i = 0; i < n_pols; ++i)
167  {
168  if (values.size() != 0)
169  {
170  values[i] = compute_value(i, unit_point);
171  }
172  if (grads.size() != 0)
173  {
174  grads[i] = compute_grad(i, unit_point);
175  }
176  if (grad_grads.size() != 0)
177  {
178  grad_grads[i] = compute_grad_grad(i, unit_point);
179  }
180  if (third_derivatives.size() != 0)
181  {
182  third_derivatives[i] = compute_derivative<3>(i, unit_point);
183  }
184  if (fourth_derivatives.size() != 0)
185  {
186  fourth_derivatives[i] = compute_derivative<4>(i, unit_point);
187  }
188  }
189 }
190 
191 
192 
193 template <int dim>
194 std::unique_ptr<ScalarPolynomialsBase<dim>>
196 {
197  return std_cxx14::make_unique<PolynomialsRannacherTurek<dim>>(*this);
198 }
199 
200 
201 // explicit instantiations
202 #include "polynomials_rannacher_turek.inst"
203 
polynomials_rannacher_turek.h
PolynomialsRannacherTurek::compute_value
double compute_value(const unsigned int i, const Point< dim > &p) const
Definition: polynomials_rannacher_turek.cc:35
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
GeometryInfo
Definition: geometry_info.h:1224
PolynomialsRannacherTurek::compute_grad
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const
Definition: polynomials_rannacher_turek.cc:68
PolynomialsRannacherTurek::PolynomialsRannacherTurek
PolynomialsRannacherTurek()
Definition: polynomials_rannacher_turek.cc:25
ScalarPolynomialsBase::n
unsigned int n() const
Definition: scalar_polynomials_base.h:164
PolynomialsRannacherTurek::compute_grad_grad
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const
Definition: polynomials_rannacher_turek.cc:105
geometry_info.h
Tensor< 1, dim >
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
PolynomialsRannacherTurek::clone
virtual std::unique_ptr< ScalarPolynomialsBase< dim > > clone() const override
Definition: polynomials_rannacher_turek.cc:195
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
ScalarPolynomialsBase::n_pols
const unsigned int n_pols
Definition: scalar_polynomials_base.h:157
StandardExceptions::ExcDimensionMismatch
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Point< dim >
memory.h
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
PolynomialsRannacherTurek::evaluate
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
Definition: polynomials_rannacher_turek.cc:146
ScalarPolynomialsBase
Definition: scalar_polynomials_base.h:63