Reference documentation for deal.II version 9.0.0
polynomials_rannacher_turek.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2015 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #include <deal.II/base/polynomials_rannacher_turek.h>
18 #include <deal.II/base/geometry_info.h>
19 
20 DEAL_II_NAMESPACE_OPEN
21 
22 
23 template <int dim>
25 {
26  Assert(dim == 2, ExcNotImplemented());
27 }
28 
29 
30 
31 template <int dim>
33  const Point<dim> &p) const
34 {
35  Assert(dim == 2, ExcNotImplemented());
36  if (i == 0)
37  {
38  return (0.75 - 2.5*p(0) + 1.5*p(1) + 1.5*(p(0)*p(0) - p(1)*p(1)));
39  }
40  else if (i == 1)
41  {
42  return (-0.25 - 0.5*p(0) + 1.5*p(1) + 1.5*(p(0)*p(0) - p(1)*p(1)));
43  }
44  else if (i == 2)
45  {
46  return (0.75 + 1.5*p(0) - 2.5*p(1) - 1.5*(p(0)*p(0) - p(1)*p(1)));
47  }
48  else if (i == 3)
49  {
50  return (-0.25 + 1.5*p(0) - 0.5*p(1) - 1.5*(p(0)*p(0) - p(1)*p(1)));
51  }
52 
53  Assert(false, ExcNotImplemented());
54  return 0;
55 }
56 
57 
58 
59 template <int dim>
61  const unsigned int i,
62  const Point<dim> &p) const
63 {
64  Assert(dim == 2, ExcNotImplemented());
65  Tensor<1, dim> grad;
66  if (i == 0)
67  {
68  grad[0] = -2.5 + 3*p(0);
69  grad[1] = 1.5 - 3*p(1);
70  }
71  else if (i == 1)
72  {
73  grad[0] = -0.5 + 3.0*p(0);
74  grad[1] = 1.5 - 3.0*p(1);
75  }
76  else if (i == 2)
77  {
78  grad[0] = 1.5 - 3.0*p(0);
79  grad[1] = -2.5 + 3.0*p(1);
80  }
81  else if (i == 3)
82  {
83  grad[0] = 1.5 - 3.0*p(0);
84  grad[1] = -0.5 + 3.0*p(1);
85  }
86  else
87  {
88  Assert(false, ExcNotImplemented());
89  }
90 
91  return grad;
92 }
93 
94 
95 
96 template <int dim> Tensor<2, dim>
98  const Point<dim> & /*p*/)
99 const
100 {
101  Assert(dim == 2, ExcNotImplemented());
102  Tensor<2, dim> grad_grad;
103  if (i == 0)
104  {
105  grad_grad[0][0] = 3;
106  grad_grad[0][1] = 0;
107  grad_grad[1][0] = 0;
108  grad_grad[1][1] = -3;
109  }
110  else if (i == 1)
111  {
112  grad_grad[0][0] = 3;
113  grad_grad[0][1] = 0;
114  grad_grad[1][0] = 0;
115  grad_grad[1][1] = -3;
116  }
117  else if (i == 2)
118  {
119  grad_grad[0][0] = -3;
120  grad_grad[0][1] = 0;
121  grad_grad[1][0] = 0;
122  grad_grad[1][1] = 3;
123  }
124  else if (i == 3)
125  {
126  grad_grad[0][0] = -3;
127  grad_grad[0][1] = 0;
128  grad_grad[1][0] = 0;
129  grad_grad[1][1] = 3;
130  }
131  return grad_grad;
132 }
133 
134 
135 
136 template <int dim>
138  const Point<dim> &unit_point,
139  std::vector<double> &values,
140  std::vector<Tensor<1, dim> > &grads,
141  std::vector<Tensor<2, dim> > &grad_grads,
142  std::vector<Tensor<3, dim> > &third_derivatives,
143  std::vector<Tensor<4, dim> > &fourth_derivatives) const
144 {
145  const unsigned int n_pols = ::GeometryInfo<dim>::faces_per_cell;
146  Assert(values.size() == n_pols || values.size() == 0,
147  ExcDimensionMismatch(values.size(), n_pols));
148  Assert(grads.size() == n_pols || grads.size() == 0,
149  ExcDimensionMismatch(grads.size(), n_pols));
150  Assert(grad_grads.size() == n_pols || grad_grads.size() == 0,
151  ExcDimensionMismatch(grad_grads.size(), n_pols));
152  Assert(third_derivatives.size() == n_pols || third_derivatives.size() == 0,
153  ExcDimensionMismatch(third_derivatives.size(), n_pols));
154  Assert(fourth_derivatives.size() == n_pols || fourth_derivatives.size() == 0,
155  ExcDimensionMismatch(fourth_derivatives.size(), n_pols));
156 
157  for (unsigned int i = 0; i < n_pols; ++i)
158  {
159  if (values.size() != 0)
160  {
161  values[i] = compute_value(i, unit_point);
162  }
163  if (grads.size() != 0)
164  {
165  grads[i] = compute_grad(i, unit_point);
166  }
167  if (grad_grads.size() != 0)
168  {
169  grad_grads[i] = compute_grad_grad(i, unit_point);
170  }
171  if (third_derivatives.size() != 0)
172  {
173  third_derivatives[i] = compute_derivative<3>(i, unit_point);
174  }
175  if (fourth_derivatives.size() != 0)
176  {
177  fourth_derivatives[i] = compute_derivative<4>(i, unit_point);
178  }
179  }
180 }
181 
182 
183 // explicit instantiations
184 #include "polynomials_rannacher_turek.inst"
185 
186 DEAL_II_NAMESPACE_CLOSE
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const
double compute_value(const unsigned int i, const Point< dim > &p) const
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const
void compute(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
static ::ExceptionBase & ExcNotImplemented()