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;
271 std::vector<double> &values,
344 unsigned int n ()
const;
350 const std::vector<std::vector<Polynomials::Polynomial<double> > >
polynomials;
365 unsigned int (&indices)[dim])
const;
383 template <
int dim,
typename PolynomialType>
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)
399 index_map_inverse[i]=i;
405 template <
int dim,
typename PolynomialType>
413 return n_tensor_pols;
418 template <
int dim,
typename PolynomialType>
420 const std::vector<unsigned int> &
427 template <
int dim,
typename PolynomialType>
429 const std::vector<unsigned int> &
432 return index_map_inverse;
435 template <
int dim,
typename PolynomialType>
439 (
const unsigned int i,
442 unsigned int indices[dim];
443 compute_index (i, indices);
447 std::vector<double> tmp (5);
448 for (
unsigned int d=0;
d<dim; ++
d)
450 polynomials[indices[
d]].value (p(d), tmp);
465 for (
unsigned int d=0;
d<dim; ++
d)
467 derivative_1[
d] = 1.;
468 for (
unsigned int x=0; x<dim; ++x)
470 unsigned int x_order=0;
473 derivative_1[
d] *= v[x][x_order];
482 for (
unsigned int d1=0; d1<dim; ++d1)
483 for (
unsigned int d2=0; d2<dim; ++d2)
485 derivative_2[d1][d2] = 1.;
486 for (
unsigned int x=0; x<dim; ++x)
488 unsigned int x_order=0;
489 if (d1==x) ++x_order;
490 if (d2==x) ++x_order;
492 derivative_2[d1][d2] *= v[x][x_order];
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)
505 derivative_3[d1][d2][d3] = 1.;
506 for (
unsigned int x=0; x<dim; ++x)
508 unsigned int x_order=0;
509 if (d1==x) ++x_order;
510 if (d2==x) ++x_order;
511 if (d3==x) ++x_order;
513 derivative_3[d1][d2][d3] *= v[x][x_order];
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)
527 derivative_4[d1][d2][d3][d4] = 1.;
528 for (
unsigned int x=0; x<dim; ++x)
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;
536 derivative_4[d1][d2][d3][d4] *= v[x][x_order];
556 unsigned int indices[dim];
557 compute_index (i, indices);
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]);
569 for (
unsigned int d=0;
d<dim; ++
d)
571 derivative_1[
d] = 1.;
572 for (
unsigned int x=0; x<dim; ++x)
574 unsigned int x_order=0;
577 derivative_1[
d] *= v[x][x_order];
586 for (
unsigned int d1=0; d1<dim; ++d1)
587 for (
unsigned int d2=0; d2<dim; ++d2)
589 derivative_2[d1][d2] = 1.;
590 for (
unsigned int x=0; x<dim; ++x)
592 unsigned int x_order=0;
593 if (d1==x) ++x_order;
594 if (d2==x) ++x_order;
596 derivative_2[d1][d2] *= v[x][x_order];
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)
609 derivative_3[d1][d2][d3] = 1.;
610 for (
unsigned int x=0; x<dim; ++x)
612 unsigned int x_order=0;
613 if (d1==x) ++x_order;
614 if (d2==x) ++x_order;
615 if (d3==x) ++x_order;
617 derivative_3[d1][d2][d3] *= v[x][x_order];
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)
631 derivative_4[d1][d2][d3][d4] = 1.;
632 for (
unsigned int x=0; x<dim; ++x)
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;
640 derivative_4[d1][d2][d3][d4] *= v[x][x_order];
657 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
AnisotropicPolynomials(const std::vector< std::vector< Polynomials::Polynomial< double > > > &pols)
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
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