16#ifndef dealii_tensor_product_polynomials_h
17#define dealii_tensor_product_polynomials_h
75template <
int dim,
typename PolynomialType = Polynomials::Polynomial<
double>>
110 const std::vector<unsigned int> &
116 const std::vector<unsigned int> &
133 std::vector<double> & values,
239 virtual std::unique_ptr<ScalarPolynomialsBase<dim>>
240 clone()
const override;
252 std::vector<PolynomialType>
279 std::array<unsigned int, dim> &indices)
const;
359 std::vector<double> & values,
465 virtual std::unique_ptr<ScalarPolynomialsBase<dim>>
466 clone()
const override;
472 const std::vector<std::vector<Polynomials::Polynomial<double>>>
polynomials;
482 std::array<unsigned int, dim> &indices)
const;
500template <
int dim,
typename PolynomialType>
503 const std::vector<Pol> &pols)
505 , polynomials(pols.begin(), pols.end())
506 , index_map(this->n())
507 , index_map_inverse(this->n())
511 for (
unsigned int i = 0; i < this->
n(); ++i)
514 index_map_inverse[i] = i;
519template <
int dim,
typename PolynomialType>
520inline const std::vector<unsigned int> &
527template <
int dim,
typename PolynomialType>
528inline const std::vector<unsigned int> &
531 return index_map_inverse;
535template <
int dim,
typename PolynomialType>
539 return "TensorProductPolynomials";
543template <
int dim,
typename PolynomialType>
547 const unsigned int i,
550 std::array<unsigned int, dim> indices;
551 compute_index(i, indices);
555 std::vector<double> tmp(5);
556 for (
unsigned int d = 0;
d < dim; ++
d)
558 polynomials[indices[
d]].value(p(d), tmp);
574 for (
unsigned int d = 0;
d < dim; ++
d)
576 derivative_1[
d] = 1.;
577 for (
unsigned int x = 0; x < dim; ++x)
579 unsigned int x_order = 0;
583 derivative_1[
d] *= v[x][x_order];
593 for (
unsigned int d1 = 0; d1 < dim; ++d1)
594 for (
unsigned int d2 = 0; d2 < dim; ++d2)
596 derivative_2[d1][d2] = 1.;
597 for (
unsigned int x = 0; x < dim; ++x)
599 unsigned int x_order = 0;
605 derivative_2[d1][d2] *= v[x][x_order];
615 for (
unsigned int d1 = 0; d1 < dim; ++d1)
616 for (
unsigned int d2 = 0; d2 < dim; ++d2)
617 for (
unsigned int d3 = 0; d3 < dim; ++d3)
619 derivative_3[d1][d2][d3] = 1.;
620 for (
unsigned int x = 0; x < dim; ++x)
622 unsigned int x_order = 0;
630 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;
658 derivative_4[d1][d2][d3][d4] *= v[x][x_order];
678 compute_derivative(
const unsigned int,
const Point<0> &)
const
687template <
int dim,
typename PolynomialType>
690 const unsigned int i,
693 return compute_derivative<1>(i, p);
698template <
int dim,
typename PolynomialType>
701 const unsigned int i,
704 return compute_derivative<2>(i, p);
709template <
int dim,
typename PolynomialType>
712 const unsigned int i,
715 return compute_derivative<3>(i, p);
720template <
int dim,
typename PolynomialType>
723 const unsigned int i,
726 return compute_derivative<4>(i, p);
737 std::array<unsigned int, dim> indices;
738 compute_index(i, indices);
740 std::vector<std::vector<double>> v(dim, std::vector<double>(order + 1));
741 for (
unsigned int d = 0;
d < dim; ++
d)
742 polynomials[d][indices[d]].
value(p(d), v[d]);
751 for (
unsigned int d = 0;
d < dim; ++
d)
753 derivative_1[
d] = 1.;
754 for (
unsigned int x = 0; x < dim; ++x)
756 unsigned int x_order = 0;
760 derivative_1[
d] *= v[x][x_order];
770 for (
unsigned int d1 = 0; d1 < dim; ++d1)
771 for (
unsigned int d2 = 0; d2 < dim; ++d2)
773 derivative_2[d1][d2] = 1.;
774 for (
unsigned int x = 0; x < dim; ++x)
776 unsigned int x_order = 0;
782 derivative_2[d1][d2] *= v[x][x_order];
792 for (
unsigned int d1 = 0; d1 < dim; ++d1)
793 for (
unsigned int d2 = 0; d2 < dim; ++d2)
794 for (
unsigned int d3 = 0; d3 < dim; ++d3)
796 derivative_3[d1][d2][d3] = 1.;
797 for (
unsigned int x = 0; x < dim; ++x)
799 unsigned int x_order = 0;
807 derivative_3[d1][d2][d3] *= v[x][x_order];
817 for (
unsigned int d1 = 0; d1 < dim; ++d1)
818 for (
unsigned int d2 = 0; d2 < dim; ++d2)
819 for (
unsigned int d3 = 0; d3 < dim; ++d3)
820 for (
unsigned int d4 = 0; d4 < dim; ++d4)
822 derivative_4[d1][d2][d3][d4] = 1.;
823 for (
unsigned int x = 0; x < dim; ++x)
825 unsigned int x_order = 0;
835 derivative_4[d1][d2][d3][d4] *= v[x][x_order];
869 return compute_derivative<1>(i, p);
879 return compute_derivative<2>(i, p);
889 return compute_derivative<3>(i, p);
899 return compute_derivative<4>(i, p);
908 return "AnisotropicPolynomials";
void evaluate(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 override
double compute_value(const unsigned int i, const Point< dim > &p) const override
const std::vector< std::vector< Polynomials::Polynomial< double > > > polynomials
virtual Tensor< 1, dim > compute_1st_derivative(const unsigned int i, const Point< dim > &p) const override
static unsigned int get_n_tensor_pols(const std::vector< std::vector< Polynomials::Polynomial< double > > > &pols)
virtual Tensor< 3, dim > compute_3rd_derivative(const unsigned int i, const Point< dim > &p) const override
void compute_index(const unsigned int i, std::array< unsigned int, dim > &indices) const
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const override
virtual Tensor< 2, dim > compute_2nd_derivative(const unsigned int i, const Point< dim > &p) const override
Tensor< order, dim > compute_derivative(const unsigned int i, const Point< dim > &p) const
std::string name() const override
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const override
virtual Tensor< 4, dim > compute_4th_derivative(const unsigned int i, const Point< dim > &p) const override
virtual std::unique_ptr< ScalarPolynomialsBase< dim > > clone() const override
void evaluate(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 override
void output_indices(std::ostream &out) const
void compute_index(const unsigned int i, std::array< unsigned int, dim > &indices) const
double compute_value(const unsigned int i, const Point< dim > &p) const override
virtual std::size_t memory_consumption() const override
virtual Tensor< 2, dim > compute_2nd_derivative(const unsigned int i, const Point< dim > &p) const override
const std::vector< unsigned int > & get_numbering() const
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const override
Tensor< order, dim > compute_derivative(const unsigned int i, const Point< dim > &p) const
virtual std::unique_ptr< ScalarPolynomialsBase< dim > > clone() const override
virtual Tensor< 1, dim > compute_1st_derivative(const unsigned int i, const Point< dim > &p) const override
std::vector< unsigned int > index_map
std::vector< unsigned int > index_map_inverse
std::vector< PolynomialType > get_underlying_polynomials() const
virtual Tensor< 3, dim > compute_3rd_derivative(const unsigned int i, const Point< dim > &p) const override
std::vector< PolynomialType > polynomials
void set_numbering(const std::vector< unsigned int > &renumber)
const std::vector< unsigned int > & get_numbering_inverse() const
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const override
std::string name() const override
static constexpr unsigned int dimension
virtual Tensor< 4, dim > compute_4th_derivative(const unsigned int i, const Point< dim > &p) const override
TensorProductPolynomials(const std::vector< Pol > &pols)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertThrow(cond, exc)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
typename internal::ndarray::HelperArray< T, Ns... >::type ndarray