Reference documentation for deal.II version 8.5.1
tensor_product_polynomials_bubbles.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2012 - 2017 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 #ifndef dealii__tensor_product_polynomials_bubbles_h
17 #define dealii__tensor_product_polynomials_bubbles_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/exceptions.h>
22 #include <deal.II/base/tensor.h>
23 #include <deal.II/base/point.h>
24 #include <deal.II/base/polynomial.h>
25 #include <deal.II/base/tensor_product_polynomials.h>
26 #include <deal.II/base/utilities.h>
27 
28 #include <vector>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 
46 template <int dim>
48 {
49 public:
55  template <class Pol>
56  TensorProductPolynomialsBubbles (const std::vector<Pol> &pols);
57 
70  void compute (const Point<dim> &unit_point,
71  std::vector<double> &values,
72  std::vector<Tensor<1,dim> > &grads,
73  std::vector<Tensor<2,dim> > &grad_grads,
74  std::vector<Tensor<3,dim> > &third_derivatives,
75  std::vector<Tensor<4,dim> > &fourth_derivatives) const;
76 
89  double compute_value (const unsigned int i,
90  const Point<dim> &p) const;
91 
104  template <int order>
105  Tensor<order,dim> compute_derivative (const unsigned int i,
106  const Point<dim> &p) const;
107 
120  Tensor<1,dim> compute_grad (const unsigned int i,
121  const Point<dim> &p) const;
122 
135  Tensor<2,dim> compute_grad_grad (const unsigned int i,
136  const Point<dim> &p) const;
137 
144  unsigned int n () const;
145 };
146 
150 /* ---------------- template and inline functions ---------- */
151 
152 #ifndef DOXYGEN
153 
154 template <int dim>
155 template <class Pol>
156 inline
158 TensorProductPolynomialsBubbles(const std::vector<Pol> &pols)
159  :
160  TensorProductPolynomials<dim>(pols)
161 {
162  const unsigned int q_degree = this->polynomials.size()-1;
163  const unsigned int n_bubbles = ((q_degree<=1)?1:dim);
164  // append index for renumbering
165  for (unsigned int i=0; i<n_bubbles; ++i)
166  {
167  this->index_map.push_back(i+this->n_tensor_pols);
168  this->index_map_inverse.push_back(i+this->n_tensor_pols);
169  }
170 }
171 
172 
173 
174 template <int dim>
175 inline
176 unsigned int
178 {
179  return this->n_tensor_pols+dim;
180 }
181 
182 
183 
184 template <>
185 inline
186 unsigned int
188 {
190 }
191 
192 template <int dim>
193 template <int order>
196  const Point<dim> &p) const
197 {
198  const unsigned int q_degree = this->polynomials.size()-1;
199  const unsigned int max_q_indices = this->n_tensor_pols;
200  const unsigned int n_bubbles = ((q_degree<=1)?1:dim);
201  (void)n_bubbles;
202  Assert (i<max_q_indices+n_bubbles, ExcInternalError());
203 
204  // treat the regular basis functions
205  if (i<max_q_indices)
206  return this->TensorProductPolynomials<dim>::template compute_derivative<order>(i,p);
207 
208  const unsigned int comp = i - this->n_tensor_pols;
209 
210  Tensor<order,dim> derivative;
211  switch (order)
212  {
213  case 1:
214  {
215  Tensor<1,dim> &derivative_1 = *reinterpret_cast<Tensor<1,dim>*>(&derivative);
216 
217  for (unsigned int d=0; d<dim ; ++d)
218  {
219  derivative_1[d] = 1.;
220  //compute grad(4*\prod_{i=1}^d (x_i(1-x_i)))(p)
221  for (unsigned j=0; j<dim; ++j)
222  derivative_1[d] *= (d==j ? 4*(1-2*p(j)) : 4*p(j)*(1-p(j)));
223  // and multiply with (2*x_i-1)^{r-1}
224  for (unsigned int i=0; i<q_degree-1; ++i)
225  derivative_1[d]*=2*p(comp)-1;
226  }
227 
228  if (q_degree>=2)
229  {
230  //add \prod_{i=1}^d 4*(x_i(1-x_i))(p)
231  double value=1.;
232  for (unsigned int j=0; j < dim; ++j)
233  value*=4*p(j)*(1-p(j));
234  //and multiply with grad(2*x_i-1)^{r-1}
235  double tmp=value*2*(q_degree-1);
236  for (unsigned int i=0; i<q_degree-2; ++i)
237  tmp*=2*p(comp)-1;
238  derivative_1[comp]+=tmp;
239  }
240 
241  return derivative;
242  }
243  case 2:
244  {
245  Tensor<2,dim> &derivative_2 = *reinterpret_cast<Tensor<2,dim>*>(&derivative);
246 
247  double v [dim+1][3];
248  {
249  for (unsigned int c=0; c<dim; ++c)
250  {
251  v[c][0] = 4*p(c)*(1-p(c));
252  v[c][1] = 4*(1-2*p(c));
253  v[c][2] = -8;
254  }
255 
256  double tmp=1.;
257  for (unsigned int i=0; i<q_degree-1; ++i)
258  tmp *= 2*p(comp)-1;
259  v[dim][0] = tmp;
260 
261  if (q_degree>=2)
262  {
263  double tmp = 2*(q_degree-1);
264  for (unsigned int i=0; i<q_degree-2; ++i)
265  tmp *= 2*p(comp)-1;
266  v[dim][1] = tmp;
267  }
268  else
269  v[dim][1] = 0.;
270 
271  if (q_degree>=3)
272  {
273  double tmp=4*(q_degree-2)*(q_degree-1);
274  for (unsigned int i=0; i<q_degree-3; ++i)
275  tmp *= 2*p(comp)-1;
276  v[dim][2] = tmp;
277  }
278  else
279  v[dim][2] = 0.;
280  }
281 
282  //calculate (\partial_j \partial_k \psi) * monomial
283  Tensor<2,dim> grad_grad_1;
284  for (unsigned int d1=0; d1<dim; ++d1)
285  for (unsigned int d2=0; d2<dim; ++d2)
286  {
287  grad_grad_1[d1][d2] = v[dim][0];
288  for (unsigned int x=0; x<dim; ++x)
289  {
290  unsigned int derivative=0;
291  if (d1==x || d2==x)
292  {
293  if (d1==d2)
294  derivative=2;
295  else
296  derivative=1;
297  }
298  grad_grad_1[d1][d2] *= v[x][derivative];
299  }
300  }
301 
302  //calculate (\partial_j \psi) *(\partial_k monomial)
303  // and (\partial_k \psi) *(\partial_j monomial)
304  Tensor<2,dim> grad_grad_2;
305  Tensor<2,dim> grad_grad_3;
306  for (unsigned int d=0; d<dim; ++d)
307  {
308  grad_grad_2[d][comp] = v[dim][1];
309  grad_grad_3[comp][d] = v[dim][1];
310  for (unsigned int x=0; x<dim; ++x)
311  {
312  grad_grad_2[d][comp] *= v[x][d==x];
313  grad_grad_3[comp][d] *= v[x][d==x];
314  }
315  }
316 
317  //calculate \psi *(\partial j \partial_k monomial) and sum
318  double psi_value = 1.;
319  for (unsigned int x=0; x<dim; ++x)
320  psi_value *= v[x][0];
321 
322  for (unsigned int d1=0; d1<dim; ++d1)
323  for (unsigned int d2=0; d2<dim; ++d2)
324  derivative_2[d1][d2] = grad_grad_1[d1][d2]
325  +grad_grad_2[d1][d2]
326  +grad_grad_3[d1][d2];
327  derivative_2[comp][comp]+=psi_value*v[dim][2];
328 
329  return derivative;
330  }
331  default:
332  {
333  Assert (false, ExcNotImplemented());
334  return derivative;
335  }
336  }
337 }
338 
339 
340 #endif // DOXYGEN
341 DEAL_II_NAMESPACE_CLOSE
342 
343 #endif
static const unsigned int invalid_unsigned_int
Definition: types.h:170
std::vector< Polynomials::Polynomial< double > > polynomials
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const
#define Assert(cond, exc)
Definition: exceptions.h:313
std::vector< unsigned int > index_map_inverse
Tensor< order, dim > compute_derivative(const unsigned int i, const Point< dim > &p) const
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
TensorProductPolynomialsBubbles(const std::vector< Pol > &pols)
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
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const
double compute_value(const unsigned int i, const Point< dim > &p) const
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcInternalError()