Reference documentation for deal.II version 8.4.2
tensor_product_polynomials_bubbles.h
1 // ---------------------------------------------------------------------
2 // @f$Id@f$
3 //
4 // Copyright (C) 2012 - 2016 by the deal.II authors
5 //
6 // This file is part of the deal.II library.
7 //
8 // The deal.II library is free software; you can use it, redistribute
9 // it, and/or modify it under the terms of the GNU Lesser General
10 // Public License as published by the Free Software Foundation; either
11 // version 2.1 of the License, or (at your option) any later version.
12 // The full text of the license can be found in the file LICENSE at
13 // the top level of the deal.II distribution.
14 //
15 // ---------------------------------------------------------------------
16 
17 #ifndef dealii__tensor_product_polynomials_bubbles_h
18 #define dealii__tensor_product_polynomials_bubbles_h
19 
20 
21 #include <deal.II/base/config.h>
22 #include <deal.II/base/exceptions.h>
23 #include <deal.II/base/tensor.h>
24 #include <deal.II/base/point.h>
25 #include <deal.II/base/polynomial.h>
26 #include <deal.II/base/tensor_product_polynomials.h>
27 #include <deal.II/base/utilities.h>
28 
29 #include <vector>
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 
47 template <int dim>
49 {
50 public:
55  static const unsigned int dimension = dim;
56 
62  template <class Pol>
63  TensorProductPolynomialsBubbles (const std::vector<Pol> &pols);
64 
77  void compute (const Point<dim> &unit_point,
78  std::vector<double> &values,
79  std::vector<Tensor<1,dim> > &grads,
80  std::vector<Tensor<2,dim> > &grad_grads,
81  std::vector<Tensor<3,dim> > &third_derivatives,
82  std::vector<Tensor<4,dim> > &fourth_derivatives) const;
83 
96  double compute_value (const unsigned int i,
97  const Point<dim> &p) const;
98 
111  template <int order>
112  Tensor<order,dim> compute_derivative (const unsigned int i,
113  const Point<dim> &p) const;
114 
127  Tensor<1,dim> compute_grad (const unsigned int i,
128  const Point<dim> &p) const;
129 
142  Tensor<2,dim> compute_grad_grad (const unsigned int i,
143  const Point<dim> &p) const;
144 
151  unsigned int n () const;
152 };
153 
157 /* ---------------- template and inline functions ---------- */
158 
159 #ifndef DOXYGEN
160 
161 template <int dim>
162 template <class Pol>
163 inline
165 TensorProductPolynomialsBubbles(const std::vector<Pol> &pols)
166  :
167  TensorProductPolynomials<dim>(pols)
168 {
169  const unsigned int q_degree = this->polynomials.size()-1;
170  const unsigned int n_bubbles = ((q_degree<=1)?1:dim);
171  // append index for renumbering
172  for (unsigned int i=0; i<n_bubbles; ++i)
173  {
174  this->index_map.push_back(i+this->n_tensor_pols);
175  this->index_map_inverse.push_back(i+this->n_tensor_pols);
176  }
177 }
178 
179 
180 
181 template <int dim>
182 inline
183 unsigned int
185 {
186  return this->n_tensor_pols+dim;
187 }
188 
189 
190 
191 template <>
192 inline
193 unsigned int
195 {
197 }
198 
199 template <int dim>
200 template <int order>
203  const Point<dim> &p) const
204 {
205  const unsigned int q_degree = this->polynomials.size()-1;
206  const unsigned int max_q_indices = this->n_tensor_pols;
207  const unsigned int n_bubbles = ((q_degree<=1)?1:dim);
208  (void)n_bubbles;
209  Assert (i<max_q_indices+n_bubbles, ExcInternalError());
210 
211  // treat the regular basis functions
212  if (i<max_q_indices)
213  return this->TensorProductPolynomials<dim>::template compute_derivative<order>(i,p);
214 
215  const unsigned int comp = i - this->n_tensor_pols;
216 
217  Tensor<order,dim> derivative;
218  switch (order)
219  {
220  case 1:
221  {
222  Tensor<1,dim> &derivative_1 = *reinterpret_cast<Tensor<1,dim>*>(&derivative);
223 
224  for (unsigned int d=0; d<dim ; ++d)
225  {
226  derivative_1[d] = 1.;
227  //compute grad(4*\prod_{i=1}^d (x_i(1-x_i)))(p)
228  for (unsigned j=0; j<dim; ++j)
229  derivative_1[d] *= (d==j ? 4*(1-2*p(j)) : 4*p(j)*(1-p(j)));
230  // and multiply with (2*x_i-1)^{r-1}
231  for (unsigned int i=0; i<q_degree-1; ++i)
232  derivative_1[d]*=2*p(comp)-1;
233  }
234 
235  if (q_degree>=2)
236  {
237  //add \prod_{i=1}^d 4*(x_i(1-x_i))(p)
238  double value=1.;
239  for (unsigned int j=0; j < dim; ++j)
240  value*=4*p(j)*(1-p(j));
241  //and multiply with grad(2*x_i-1)^{r-1}
242  double tmp=value*2*(q_degree-1);
243  for (unsigned int i=0; i<q_degree-2; ++i)
244  tmp*=2*p(comp)-1;
245  derivative_1[comp]+=tmp;
246  }
247 
248  return derivative;
249  }
250  case 2:
251  {
252  Tensor<2,dim> &derivative_2 = *reinterpret_cast<Tensor<2,dim>*>(&derivative);
253 
254  double v [dim+1][3];
255  {
256  for (unsigned int c=0; c<dim; ++c)
257  {
258  v[c][0] = 4*p(c)*(1-p(c));
259  v[c][1] = 4*(1-2*p(c));
260  v[c][2] = -8;
261  }
262 
263  double tmp=1.;
264  for (unsigned int i=0; i<q_degree-1; ++i)
265  tmp *= 2*p(comp)-1;
266  v[dim][0] = tmp;
267 
268  if (q_degree>=2)
269  {
270  double tmp = 2*(q_degree-1);
271  for (unsigned int i=0; i<q_degree-2; ++i)
272  tmp *= 2*p(comp)-1;
273  v[dim][1] = tmp;
274  }
275  else
276  v[dim][1] = 0.;
277 
278  if (q_degree>=3)
279  {
280  double tmp=4*(q_degree-2)*(q_degree-1);
281  for (unsigned int i=0; i<q_degree-3; ++i)
282  tmp *= 2*p(comp)-1;
283  v[dim][2] = tmp;
284  }
285  else
286  v[dim][2] = 0.;
287  }
288 
289  //calculate (\partial_j \partial_k \psi) * monomial
290  Tensor<2,dim> grad_grad_1;
291  for (unsigned int d1=0; d1<dim; ++d1)
292  for (unsigned int d2=0; d2<dim; ++d2)
293  {
294  grad_grad_1[d1][d2] = v[dim][0];
295  for (unsigned int x=0; x<dim; ++x)
296  {
297  unsigned int derivative=0;
298  if (d1==x || d2==x)
299  {
300  if (d1==d2)
301  derivative=2;
302  else
303  derivative=1;
304  }
305  grad_grad_1[d1][d2] *= v[x][derivative];
306  }
307  }
308 
309  //calculate (\partial_j \psi) *(\partial_k monomial)
310  // and (\partial_k \psi) *(\partial_j monomial)
311  Tensor<2,dim> grad_grad_2;
312  Tensor<2,dim> grad_grad_3;
313  for (unsigned int d=0; d<dim; ++d)
314  {
315  grad_grad_2[d][comp] = v[dim][1];
316  grad_grad_3[comp][d] = v[dim][1];
317  for (unsigned int x=0; x<dim; ++x)
318  {
319  grad_grad_2[d][comp] *= v[x][d==x];
320  grad_grad_3[comp][d] *= v[x][d==x];
321  }
322  }
323 
324  //calculate \psi *(\partial j \partial_k monomial) and sum
325  double psi_value = 1.;
326  for (unsigned int x=0; x<dim; ++x)
327  psi_value *= v[x][0];
328 
329  for (unsigned int d1=0; d1<dim; ++d1)
330  for (unsigned int d2=0; d2<dim; ++d2)
331  derivative_2[d1][d2] = grad_grad_1[d1][d2]
332  +grad_grad_2[d1][d2]
333  +grad_grad_3[d1][d2];
334  derivative_2[comp][comp]+=psi_value*v[dim][2];
335 
336  return derivative;
337  }
338  default:
339  {
340  Assert (false, ExcNotImplemented());
341  return derivative;
342  }
343  }
344 }
345 
346 
347 #endif // DOXYGEN
348 DEAL_II_NAMESPACE_CLOSE
349 
350 #endif
static const unsigned int invalid_unsigned_int
Definition: types.h:164
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:294
std::vector< unsigned int > index_map_inverse
Tensor< order, dim > compute_derivative(const unsigned int i, const Point< dim > &p) const
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