Reference documentation for deal.II version 8.5.1
tensor_product_polynomials.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 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_h
17 #define dealii__tensor_product_polynomials_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/utilities.h>
26 
27 #include <vector>
28 
29 DEAL_II_NAMESPACE_OPEN
30 
63 template <int dim, typename PolynomialType=Polynomials::Polynomial<double> >
65 {
66 public:
71  static const unsigned int dimension = dim;
72 
79  template <class Pol>
80  TensorProductPolynomials (const std::vector<Pol> &pols);
81 
85  void output_indices(std::ostream &out) const;
86 
91  void set_numbering(const std::vector<unsigned int> &renumber);
92 
96  const std::vector<unsigned int> &get_numbering() const;
97 
101  const std::vector<unsigned int> &get_numbering_inverse() const;
102 
115  void compute (const Point<dim> &unit_point,
116  std::vector<double> &values,
117  std::vector<Tensor<1,dim> > &grads,
118  std::vector<Tensor<2,dim> > &grad_grads,
119  std::vector<Tensor<3,dim> > &third_derivatives,
120  std::vector<Tensor<4,dim> > &fourth_derivatives) const;
121 
134  double compute_value (const unsigned int i,
135  const Point<dim> &p) const;
136 
151  template <int order>
152  Tensor<order,dim> compute_derivative (const unsigned int i,
153  const Point<dim> &p) const;
154 
167  Tensor<1,dim> compute_grad (const unsigned int i,
168  const Point<dim> &p) const;
169 
182  Tensor<2,dim> compute_grad_grad (const unsigned int i,
183  const Point<dim> &p) const;
184 
189  unsigned int n () const;
190 
191 
192 protected:
196  std::vector<PolynomialType> polynomials;
197 
201  unsigned int n_tensor_pols;
202 
206  std::vector<unsigned int> index_map;
207 
211  std::vector<unsigned int> index_map_inverse;
212 
219  // fix to avoid compiler warnings about zero length arrays
220  void compute_index (const unsigned int i,
221  unsigned int (&indices)[(dim>0?dim:1)]) const;
222 };
223 
224 
225 
244 template <int dim>
246 {
247 public:
255  AnisotropicPolynomials (const std::vector<std::vector<Polynomials::Polynomial<double> > > &pols);
256 
270  void compute (const Point<dim> &unit_point,
271  std::vector<double> &values,
272  std::vector<Tensor<1,dim> > &grads,
273  std::vector<Tensor<2,dim> > &grad_grads,
274  std::vector<Tensor<3,dim> > &third_derivatives,
275  std::vector<Tensor<4,dim> > &fourth_derivatives) const;
276 
289  double compute_value (const unsigned int i,
290  const Point<dim> &p) const;
291 
306  template <int order>
307  Tensor<order,dim> compute_derivative (const unsigned int i,
308  const Point<dim> &p) const;
309 
322  Tensor<1,dim> compute_grad (const unsigned int i,
323  const Point<dim> &p) const;
324 
337  Tensor<2,dim> compute_grad_grad (const unsigned int i,
338  const Point<dim> &p) const;
339 
344  unsigned int n () const;
345 
346 private:
350  const std::vector<std::vector<Polynomials::Polynomial<double> > > polynomials;
351 
356  const unsigned int n_tensor_pols;
357 
364  void compute_index (const unsigned int i,
365  unsigned int (&indices)[dim]) const;
366 
370  static
371  unsigned int
372  get_n_tensor_pols (const std::vector<std::vector<Polynomials::Polynomial<double> > > &pols);
373 };
374 
377 #ifndef DOXYGEN
378 
379 
380 /* ---------------- template and inline functions ---------- */
381 
382 
383 template <int dim, typename PolynomialType>
384 template <class Pol>
385 inline
387 TensorProductPolynomials(const std::vector<Pol> &pols)
388  :
389  polynomials (pols.begin(), pols.end()),
390  n_tensor_pols(Utilities::fixed_power<dim>(pols.size())),
391  index_map(n_tensor_pols),
392  index_map_inverse(n_tensor_pols)
393 {
394  // per default set this index map to identity. This map can be changed by
395  // the user through the set_numbering() function
396  for (unsigned int i=0; i<n_tensor_pols; ++i)
397  {
398  index_map[i]=i;
399  index_map_inverse[i]=i;
400  }
401 }
402 
403 
404 
405 template <int dim, typename PolynomialType>
406 inline
407 unsigned int
409 {
410  if (dim == 0)
412  else
413  return n_tensor_pols;
414 }
415 
416 
417 
418 template <int dim, typename PolynomialType>
419 inline
420 const std::vector<unsigned int> &
422 {
423  return index_map;
424 }
425 
426 
427 template <int dim, typename PolynomialType>
428 inline
429 const std::vector<unsigned int> &
431 {
432  return index_map_inverse;
433 }
434 
435 template <int dim, typename PolynomialType>
436 template <int order>
439 (const unsigned int i,
440  const Point<dim> &p) const
441 {
442  unsigned int indices[dim];
443  compute_index (i, indices);
444 
445  double v [dim][5];
446  {
447  std::vector<double> tmp (5);
448  for (unsigned int d=0; d<dim; ++d)
449  {
450  polynomials[indices[d]].value (p(d), tmp);
451  v[d][0] = tmp[0];
452  v[d][1] = tmp[1];
453  v[d][2] = tmp[2];
454  v[d][3] = tmp[3];
455  v[d][4] = tmp[4];
456  }
457  }
458 
459  Tensor<order,dim> derivative;
460  switch (order)
461  {
462  case 1:
463  {
464  Tensor<1,dim> &derivative_1 = *reinterpret_cast<Tensor<1,dim>*>(&derivative);
465  for (unsigned int d=0; d<dim; ++d)
466  {
467  derivative_1[d] = 1.;
468  for (unsigned int x=0; x<dim; ++x)
469  {
470  unsigned int x_order=0;
471  if (d==x) ++x_order;
472 
473  derivative_1[d] *= v[x][x_order];
474  }
475  }
476 
477  return derivative;
478  }
479  case 2:
480  {
481  Tensor<2,dim> &derivative_2 = *reinterpret_cast<Tensor<2,dim>*>(&derivative);
482  for (unsigned int d1=0; d1<dim; ++d1)
483  for (unsigned int d2=0; d2<dim; ++d2)
484  {
485  derivative_2[d1][d2] = 1.;
486  for (unsigned int x=0; x<dim; ++x)
487  {
488  unsigned int x_order=0;
489  if (d1==x) ++x_order;
490  if (d2==x) ++x_order;
491 
492  derivative_2[d1][d2] *= v[x][x_order];
493  }
494  }
495 
496  return derivative;
497  }
498  case 3:
499  {
500  Tensor<3,dim> &derivative_3 = *reinterpret_cast<Tensor<3,dim>*>(&derivative);
501  for (unsigned int d1=0; d1<dim; ++d1)
502  for (unsigned int d2=0; d2<dim; ++d2)
503  for (unsigned int d3=0; d3<dim; ++d3)
504  {
505  derivative_3[d1][d2][d3] = 1.;
506  for (unsigned int x=0; x<dim; ++x)
507  {
508  unsigned int x_order=0;
509  if (d1==x) ++x_order;
510  if (d2==x) ++x_order;
511  if (d3==x) ++x_order;
512 
513  derivative_3[d1][d2][d3] *= v[x][x_order];
514  }
515  }
516 
517  return derivative;
518  }
519  case 4:
520  {
521  Tensor<4,dim> &derivative_4 = *reinterpret_cast<Tensor<4,dim>*>(&derivative);
522  for (unsigned int d1=0; d1<dim; ++d1)
523  for (unsigned int d2=0; d2<dim; ++d2)
524  for (unsigned int d3=0; d3<dim; ++d3)
525  for (unsigned int d4=0; d4<dim; ++d4)
526  {
527  derivative_4[d1][d2][d3][d4] = 1.;
528  for (unsigned int x=0; x<dim; ++x)
529  {
530  unsigned int x_order=0;
531  if (d1==x) ++x_order;
532  if (d2==x) ++x_order;
533  if (d3==x) ++x_order;
534  if (d4==x) ++x_order;
535 
536  derivative_4[d1][d2][d3][d4] *= v[x][x_order];
537  }
538  }
539 
540  return derivative;
541  }
542  default:
543  {
544  Assert (false, ExcNotImplemented());
545  return derivative;
546  }
547  }
548 }
549 
550 template <int dim>
551 template <int order>
554  const Point<dim> &p) const
555 {
556  unsigned int indices[dim];
557  compute_index (i, indices);
558 
559  std::vector<std::vector<double> > v(dim, std::vector<double> (order+1));
560  for (unsigned int d=0; d<dim; ++d)
561  polynomials[d][indices[d]].value(p(d), v[d]);
562 
563  Tensor<order,dim> derivative;
564  switch (order)
565  {
566  case 1:
567  {
568  Tensor<1,dim> &derivative_1 = *reinterpret_cast<Tensor<1,dim>*>(&derivative);
569  for (unsigned int d=0; d<dim; ++d)
570  {
571  derivative_1[d] = 1.;
572  for (unsigned int x=0; x<dim; ++x)
573  {
574  unsigned int x_order=0;
575  if (d==x) ++x_order;
576 
577  derivative_1[d] *= v[x][x_order];
578  }
579  }
580 
581  return derivative;
582  }
583  case 2:
584  {
585  Tensor<2,dim> &derivative_2 = *reinterpret_cast<Tensor<2,dim>*>(&derivative);
586  for (unsigned int d1=0; d1<dim; ++d1)
587  for (unsigned int d2=0; d2<dim; ++d2)
588  {
589  derivative_2[d1][d2] = 1.;
590  for (unsigned int x=0; x<dim; ++x)
591  {
592  unsigned int x_order=0;
593  if (d1==x) ++x_order;
594  if (d2==x) ++x_order;
595 
596  derivative_2[d1][d2] *= v[x][x_order];
597  }
598  }
599 
600  return derivative;
601  }
602  case 3:
603  {
604  Tensor<3,dim> &derivative_3 = *reinterpret_cast<Tensor<3,dim>*>(&derivative);
605  for (unsigned int d1=0; d1<dim; ++d1)
606  for (unsigned int d2=0; d2<dim; ++d2)
607  for (unsigned int d3=0; d3<dim; ++d3)
608  {
609  derivative_3[d1][d2][d3] = 1.;
610  for (unsigned int x=0; x<dim; ++x)
611  {
612  unsigned int x_order=0;
613  if (d1==x) ++x_order;
614  if (d2==x) ++x_order;
615  if (d3==x) ++x_order;
616 
617  derivative_3[d1][d2][d3] *= v[x][x_order];
618  }
619  }
620 
621  return derivative;
622  }
623  case 4:
624  {
625  Tensor<4,dim> &derivative_4 = *reinterpret_cast<Tensor<4,dim>*>(&derivative);
626  for (unsigned int d1=0; d1<dim; ++d1)
627  for (unsigned int d2=0; d2<dim; ++d2)
628  for (unsigned int d3=0; d3<dim; ++d3)
629  for (unsigned int d4=0; d4<dim; ++d4)
630  {
631  derivative_4[d1][d2][d3][d4] = 1.;
632  for (unsigned int x=0; x<dim; ++x)
633  {
634  unsigned int x_order=0;
635  if (d1==x) ++x_order;
636  if (d2==x) ++x_order;
637  if (d3==x) ++x_order;
638  if (d4==x) ++x_order;
639 
640  derivative_4[d1][d2][d3][d4] *= v[x][x_order];
641  }
642  }
643 
644  return derivative;
645  }
646  default:
647  {
648  Assert (false, ExcNotImplemented());
649  return derivative;
650  }
651  }
652 }
653 
654 
655 
656 #endif // DOXYGEN
657 DEAL_II_NAMESPACE_CLOSE
658 
659 #endif
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const
Tensor< order, dim > compute_derivative(const unsigned int i, const Point< dim > &p) const
unsigned int n() const
AnisotropicPolynomials(const std::vector< std::vector< Polynomials::Polynomial< double > > > &pols)
Tensor< order, dim > compute_derivative(const unsigned int i, const Point< dim > &p) const
static const unsigned int invalid_unsigned_int
Definition: types.h:170
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const
std::vector< PolynomialType > polynomials
void compute_index(const unsigned int i, unsigned int(&indices)[(dim >0?dim:1)]) const
void output_indices(std::ostream &out) const
void set_numbering(const std::vector< unsigned int > &renumber)
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
TensorProductPolynomials(const std::vector< Pol > &pols)
void compute_index(const unsigned int i, unsigned int(&indices)[dim]) const
static unsigned int get_n_tensor_pols(const std::vector< std::vector< Polynomials::Polynomial< double > > > &pols)
#define Assert(cond, exc)
Definition: exceptions.h:313
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
const std::vector< std::vector< Polynomials::Polynomial< double > > > polynomials
std::vector< unsigned int > index_map_inverse
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
const std::vector< unsigned int > & get_numbering() const
Definition: mpi.h:48
static const unsigned int dimension
std::vector< unsigned int > index_map
double compute_value(const unsigned int i, const Point< dim > &p) const
const std::vector< unsigned int > & get_numbering_inverse() const
static ::ExceptionBase & ExcNotImplemented()
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
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const