Reference documentation for deal.II version 8.5.1
polynomials_rannacher_turek.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2015 - 2016 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 #ifndef dealii__polynomials_rannacher_turek_h
18 #define dealii__polynomials_rannacher_turek_h
19 
20 #include <deal.II/base/point.h>
21 #include <deal.II/base/tensor.h>
22 #include <vector>
23 
24 DEAL_II_NAMESPACE_OPEN
25 
26 
39 template <int dim>
41 {
42 public:
46  static const unsigned int dimension = dim;
47 
52 
56  double compute_value(const unsigned int i,
57  const Point<dim> &p) const;
58 
64  template <int order>
65  Tensor<order,dim> compute_derivative (const unsigned int i,
66  const Point<dim> &p) const;
67 
71  Tensor<1, dim> compute_grad(const unsigned int i,
72  const Point<dim> &p) const;
73 
77  Tensor<2, dim> compute_grad_grad(const unsigned int i,
78  const Point<dim> &p) const;
79 
86  void compute(const Point<dim> &unit_point,
87  std::vector<double> &values,
88  std::vector<Tensor<1, dim> > &grads,
89  std::vector<Tensor<2,dim> > &grad_grads,
90  std::vector<Tensor<3,dim> > &third_derivatives,
91  std::vector<Tensor<4,dim> > &fourth_derivatives) const;
92 };
93 
94 
95 namespace internal
96 {
98  {
99  template <int order, int dim>
100  inline
102  compute_derivative (const unsigned int,
103  const Point<dim> &)
104  {
105  Assert (dim == 2, ExcNotImplemented());
106  return Tensor<order,dim>();
107  }
108 
109 
110  template <int order>
111  inline
113  compute_derivative (const unsigned int i,
114  const Point<2> &p)
115  {
116  const unsigned int dim = 2;
117 
118  Tensor<order,dim> derivative;
119  switch (order)
120  {
121  case 1:
122  {
123  Tensor<1,dim> &grad = *reinterpret_cast<Tensor<1,dim>*>(&derivative);
124  if (i == 0)
125  {
126  grad[0] = -2.5 + 3*p(0);
127  grad[1] = 1.5 - 3*p(1);
128  }
129  else if (i == 1)
130  {
131  grad[0] = -0.5 + 3.0*p(0);
132  grad[1] = 1.5 - 3.0*p(1);
133  }
134  else if (i == 2)
135  {
136  grad[0] = 1.5 - 3.0*p(0);
137  grad[1] = -2.5 + 3.0*p(1);
138  }
139  else if (i == 3)
140  {
141  grad[0] = 1.5 - 3.0*p(0);
142  grad[1] = -0.5 + 3.0*p(1);
143  }
144  else
145  {
146  Assert(false, ExcNotImplemented());
147  }
148  return derivative;
149  }
150  case 2:
151  {
152  Tensor<2,dim> &grad_grad = *reinterpret_cast<Tensor<2,dim>*>(&derivative);
153  if (i == 0)
154  {
155  grad_grad[0][0] = 3;
156  grad_grad[0][1] = 0;
157  grad_grad[1][0] = 0;
158  grad_grad[1][1] = -3;
159  }
160  else if (i == 1)
161  {
162  grad_grad[0][0] = 3;
163  grad_grad[0][1] = 0;
164  grad_grad[1][0] = 0;
165  grad_grad[1][1] = -3;
166  }
167  else if (i == 2)
168  {
169  grad_grad[0][0] = -3;
170  grad_grad[0][1] = 0;
171  grad_grad[1][0] = 0;
172  grad_grad[1][1] = 3;
173  }
174  else if (i == 3)
175  {
176  grad_grad[0][0] = -3;
177  grad_grad[0][1] = 0;
178  grad_grad[1][0] = 0;
179  grad_grad[1][1] = 3;
180  }
181  return derivative;
182  }
183  default:
184  {
185  // higher derivatives are all zero
186  return Tensor<order,dim>();
187  }
188  }
189  }
190  }
191 }
192 
193 
194 
195 // template functions
196 template <int dim>
197 template <int order>
200  const Point<dim> &p) const
201 {
202  return internal::PolynomialsRannacherTurek::compute_derivative<order> (i, p);
203 }
204 
205 
206 DEAL_II_NAMESPACE_CLOSE
207 
208 #endif
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const
static const unsigned int dimension
double compute_value(const unsigned int i, const Point< dim > &p) const
#define Assert(cond, exc)
Definition: exceptions.h:313
Tensor< order, dim > compute_derivative(const unsigned int i, const Point< dim > &p) const
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const
Definition: mpi.h:41
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()