Reference documentation for deal.II version 9.0.0
tensor_product_polynomials.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2018 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 
253 template <int dim>
255 {
256 public:
269  AnisotropicPolynomials (const std::vector<std::vector<Polynomials::Polynomial<double> > > &base_polynomials);
270 
284  void compute (const Point<dim> &unit_point,
285  std::vector<double> &values,
286  std::vector<Tensor<1,dim> > &grads,
287  std::vector<Tensor<2,dim> > &grad_grads,
288  std::vector<Tensor<3,dim> > &third_derivatives,
289  std::vector<Tensor<4,dim> > &fourth_derivatives) const;
290 
303  double compute_value (const unsigned int i,
304  const Point<dim> &p) const;
305 
320  template <int order>
321  Tensor<order,dim> compute_derivative (const unsigned int i,
322  const Point<dim> &p) const;
323 
336  Tensor<1,dim> compute_grad (const unsigned int i,
337  const Point<dim> &p) const;
338 
351  Tensor<2,dim> compute_grad_grad (const unsigned int i,
352  const Point<dim> &p) const;
353 
358  unsigned int n () const;
359 
360 private:
364  const std::vector<std::vector<Polynomials::Polynomial<double> > > polynomials;
365 
370  const unsigned int n_tensor_pols;
371 
378  void compute_index (const unsigned int i,
379  unsigned int (&indices)[dim]) const;
380 
384  static
385  unsigned int
386  get_n_tensor_pols (const std::vector<std::vector<Polynomials::Polynomial<double> > > &pols);
387 };
388 
391 #ifndef DOXYGEN
392 
393 
394 /* ---------------- template and inline functions ---------- */
395 
396 
397 template <int dim, typename PolynomialType>
398 template <class Pol>
399 inline
401 TensorProductPolynomials(const std::vector<Pol> &pols)
402  :
403  polynomials (pols.begin(), pols.end()),
404  n_tensor_pols(Utilities::fixed_power<dim>(pols.size())),
405  index_map(n_tensor_pols),
406  index_map_inverse(n_tensor_pols)
407 {
408  // per default set this index map to identity. This map can be changed by
409  // the user through the set_numbering() function
410  for (unsigned int i=0; i<n_tensor_pols; ++i)
411  {
412  index_map[i]=i;
413  index_map_inverse[i]=i;
414  }
415 }
416 
417 
418 
419 template <int dim, typename PolynomialType>
420 inline
421 unsigned int
423 {
424  if (dim == 0)
426  else
427  return n_tensor_pols;
428 }
429 
430 
431 
432 template <int dim, typename PolynomialType>
433 inline
434 const std::vector<unsigned int> &
436 {
437  return index_map;
438 }
439 
440 
441 template <int dim, typename PolynomialType>
442 inline
443 const std::vector<unsigned int> &
445 {
446  return index_map_inverse;
447 }
448 
449 template <int dim, typename PolynomialType>
450 template <int order>
453 (const unsigned int i,
454  const Point<dim> &p) const
455 {
456  unsigned int indices[dim];
457  compute_index (i, indices);
458 
459  double v [dim][5];
460  {
461  std::vector<double> tmp (5);
462  for (unsigned int d=0; d<dim; ++d)
463  {
464  polynomials[indices[d]].value (p(d), tmp);
465  v[d][0] = tmp[0];
466  v[d][1] = tmp[1];
467  v[d][2] = tmp[2];
468  v[d][3] = tmp[3];
469  v[d][4] = tmp[4];
470  }
471  }
472 
473  Tensor<order,dim> derivative;
474  switch (order)
475  {
476  case 1:
477  {
478  Tensor<1,dim> &derivative_1 = *reinterpret_cast<Tensor<1,dim>*>(&derivative);
479  for (unsigned int d=0; d<dim; ++d)
480  {
481  derivative_1[d] = 1.;
482  for (unsigned int x=0; x<dim; ++x)
483  {
484  unsigned int x_order=0;
485  if (d==x) ++x_order;
486 
487  derivative_1[d] *= v[x][x_order];
488  }
489  }
490 
491  return derivative;
492  }
493  case 2:
494  {
495  Tensor<2,dim> &derivative_2 = *reinterpret_cast<Tensor<2,dim>*>(&derivative);
496  for (unsigned int d1=0; d1<dim; ++d1)
497  for (unsigned int d2=0; d2<dim; ++d2)
498  {
499  derivative_2[d1][d2] = 1.;
500  for (unsigned int x=0; x<dim; ++x)
501  {
502  unsigned int x_order=0;
503  if (d1==x) ++x_order;
504  if (d2==x) ++x_order;
505 
506  derivative_2[d1][d2] *= v[x][x_order];
507  }
508  }
509 
510  return derivative;
511  }
512  case 3:
513  {
514  Tensor<3,dim> &derivative_3 = *reinterpret_cast<Tensor<3,dim>*>(&derivative);
515  for (unsigned int d1=0; d1<dim; ++d1)
516  for (unsigned int d2=0; d2<dim; ++d2)
517  for (unsigned int d3=0; d3<dim; ++d3)
518  {
519  derivative_3[d1][d2][d3] = 1.;
520  for (unsigned int x=0; x<dim; ++x)
521  {
522  unsigned int x_order=0;
523  if (d1==x) ++x_order;
524  if (d2==x) ++x_order;
525  if (d3==x) ++x_order;
526 
527  derivative_3[d1][d2][d3] *= v[x][x_order];
528  }
529  }
530 
531  return derivative;
532  }
533  case 4:
534  {
535  Tensor<4,dim> &derivative_4 = *reinterpret_cast<Tensor<4,dim>*>(&derivative);
536  for (unsigned int d1=0; d1<dim; ++d1)
537  for (unsigned int d2=0; d2<dim; ++d2)
538  for (unsigned int d3=0; d3<dim; ++d3)
539  for (unsigned int d4=0; d4<dim; ++d4)
540  {
541  derivative_4[d1][d2][d3][d4] = 1.;
542  for (unsigned int x=0; x<dim; ++x)
543  {
544  unsigned int x_order=0;
545  if (d1==x) ++x_order;
546  if (d2==x) ++x_order;
547  if (d3==x) ++x_order;
548  if (d4==x) ++x_order;
549 
550  derivative_4[d1][d2][d3][d4] *= v[x][x_order];
551  }
552  }
553 
554  return derivative;
555  }
556  default:
557  {
558  Assert (false, ExcNotImplemented());
559  return derivative;
560  }
561  }
562 }
563 
564 template <int dim>
565 template <int order>
568  const Point<dim> &p) const
569 {
570  unsigned int indices[dim];
571  compute_index (i, indices);
572 
573  std::vector<std::vector<double> > v(dim, std::vector<double> (order+1));
574  for (unsigned int d=0; d<dim; ++d)
575  polynomials[d][indices[d]].value(p(d), v[d]);
576 
577  Tensor<order,dim> derivative;
578  switch (order)
579  {
580  case 1:
581  {
582  Tensor<1,dim> &derivative_1 = *reinterpret_cast<Tensor<1,dim>*>(&derivative);
583  for (unsigned int d=0; d<dim; ++d)
584  {
585  derivative_1[d] = 1.;
586  for (unsigned int x=0; x<dim; ++x)
587  {
588  unsigned int x_order=0;
589  if (d==x) ++x_order;
590 
591  derivative_1[d] *= v[x][x_order];
592  }
593  }
594 
595  return derivative;
596  }
597  case 2:
598  {
599  Tensor<2,dim> &derivative_2 = *reinterpret_cast<Tensor<2,dim>*>(&derivative);
600  for (unsigned int d1=0; d1<dim; ++d1)
601  for (unsigned int d2=0; d2<dim; ++d2)
602  {
603  derivative_2[d1][d2] = 1.;
604  for (unsigned int x=0; x<dim; ++x)
605  {
606  unsigned int x_order=0;
607  if (d1==x) ++x_order;
608  if (d2==x) ++x_order;
609 
610  derivative_2[d1][d2] *= v[x][x_order];
611  }
612  }
613 
614  return derivative;
615  }
616  case 3:
617  {
618  Tensor<3,dim> &derivative_3 = *reinterpret_cast<Tensor<3,dim>*>(&derivative);
619  for (unsigned int d1=0; d1<dim; ++d1)
620  for (unsigned int d2=0; d2<dim; ++d2)
621  for (unsigned int d3=0; d3<dim; ++d3)
622  {
623  derivative_3[d1][d2][d3] = 1.;
624  for (unsigned int x=0; x<dim; ++x)
625  {
626  unsigned int x_order=0;
627  if (d1==x) ++x_order;
628  if (d2==x) ++x_order;
629  if (d3==x) ++x_order;
630 
631  derivative_3[d1][d2][d3] *= v[x][x_order];
632  }
633  }
634 
635  return derivative;
636  }
637  case 4:
638  {
639  Tensor<4,dim> &derivative_4 = *reinterpret_cast<Tensor<4,dim>*>(&derivative);
640  for (unsigned int d1=0; d1<dim; ++d1)
641  for (unsigned int d2=0; d2<dim; ++d2)
642  for (unsigned int d3=0; d3<dim; ++d3)
643  for (unsigned int d4=0; d4<dim; ++d4)
644  {
645  derivative_4[d1][d2][d3][d4] = 1.;
646  for (unsigned int x=0; x<dim; ++x)
647  {
648  unsigned int x_order=0;
649  if (d1==x) ++x_order;
650  if (d2==x) ++x_order;
651  if (d3==x) ++x_order;
652  if (d4==x) ++x_order;
653 
654  derivative_4[d1][d2][d3][d4] *= v[x][x_order];
655  }
656  }
657 
658  return derivative;
659  }
660  default:
661  {
662  Assert (false, ExcNotImplemented());
663  return derivative;
664  }
665  }
666 }
667 
668 
669 
670 #endif // DOXYGEN
671 DEAL_II_NAMESPACE_CLOSE
672 
673 #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
Tensor< order, dim > compute_derivative(const unsigned int i, const Point< dim > &p) const
static const unsigned int invalid_unsigned_int
Definition: types.h:173
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:1142
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
AnisotropicPolynomials(const std::vector< std::vector< Polynomials::Polynomial< double > > > &base_polynomials)
Definition: cuda.h:31
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