16 #ifndef dealii_tensor_product_polynomials_h 17 #define dealii_tensor_product_polynomials_h 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> 29 DEAL_II_NAMESPACE_OPEN
63 template <
int dim,
typename PolynomialType=Polynomials::Polynomial<
double> >
91 void set_numbering(
const std::vector<unsigned int> &renumber);
116 std::vector<double> &values,
189 unsigned int n ()
const;
221 unsigned int (&indices)[(dim>0?dim:1)])
const;
285 std::vector<double> &values,
358 unsigned int n ()
const;
364 const std::vector<std::vector<Polynomials::Polynomial<double> > >
polynomials;
379 unsigned int (&indices)[dim])
const;
397 template <
int dim,
typename PolynomialType>
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)
413 index_map_inverse[i]=i;
419 template <
int dim,
typename PolynomialType>
427 return n_tensor_pols;
432 template <
int dim,
typename PolynomialType>
434 const std::vector<unsigned int> &
441 template <
int dim,
typename PolynomialType>
443 const std::vector<unsigned int> &
446 return index_map_inverse;
449 template <
int dim,
typename PolynomialType>
453 (
const unsigned int i,
456 unsigned int indices[dim];
457 compute_index (i, indices);
461 std::vector<double> tmp (5);
462 for (
unsigned int d=0;
d<dim; ++
d)
464 polynomials[indices[
d]].value (p(d), tmp);
479 for (
unsigned int d=0;
d<dim; ++
d)
481 derivative_1[
d] = 1.;
482 for (
unsigned int x=0; x<dim; ++x)
484 unsigned int x_order=0;
487 derivative_1[
d] *= v[x][x_order];
496 for (
unsigned int d1=0; d1<dim; ++d1)
497 for (
unsigned int d2=0; d2<dim; ++d2)
499 derivative_2[d1][d2] = 1.;
500 for (
unsigned int x=0; x<dim; ++x)
502 unsigned int x_order=0;
503 if (d1==x) ++x_order;
504 if (d2==x) ++x_order;
506 derivative_2[d1][d2] *= v[x][x_order];
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)
519 derivative_3[d1][d2][d3] = 1.;
520 for (
unsigned int x=0; x<dim; ++x)
522 unsigned int x_order=0;
523 if (d1==x) ++x_order;
524 if (d2==x) ++x_order;
525 if (d3==x) ++x_order;
527 derivative_3[d1][d2][d3] *= v[x][x_order];
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)
541 derivative_4[d1][d2][d3][d4] = 1.;
542 for (
unsigned int x=0; x<dim; ++x)
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;
550 derivative_4[d1][d2][d3][d4] *= v[x][x_order];
570 unsigned int indices[dim];
571 compute_index (i, indices);
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]);
583 for (
unsigned int d=0;
d<dim; ++
d)
585 derivative_1[
d] = 1.;
586 for (
unsigned int x=0; x<dim; ++x)
588 unsigned int x_order=0;
591 derivative_1[
d] *= v[x][x_order];
600 for (
unsigned int d1=0; d1<dim; ++d1)
601 for (
unsigned int d2=0; d2<dim; ++d2)
603 derivative_2[d1][d2] = 1.;
604 for (
unsigned int x=0; x<dim; ++x)
606 unsigned int x_order=0;
607 if (d1==x) ++x_order;
608 if (d2==x) ++x_order;
610 derivative_2[d1][d2] *= v[x][x_order];
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)
623 derivative_3[d1][d2][d3] = 1.;
624 for (
unsigned int x=0; x<dim; ++x)
626 unsigned int x_order=0;
627 if (d1==x) ++x_order;
628 if (d2==x) ++x_order;
629 if (d3==x) ++x_order;
631 derivative_3[d1][d2][d3] *= v[x][x_order];
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)
645 derivative_4[d1][d2][d3][d4] = 1.;
646 for (
unsigned int x=0; x<dim; ++x)
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;
654 derivative_4[d1][d2][d3][d4] *= v[x][x_order];
671 DEAL_II_NAMESPACE_CLOSE
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
Tensor< order, dim > compute_derivative(const unsigned int i, const Point< dim > &p) const
unsigned int n_tensor_pols
static const unsigned int invalid_unsigned_int
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)
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)
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()
const unsigned int n_tensor_pols
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