16#ifndef dealii_polynomial_h
17#define dealii_polynomial_h
62 template <
typename number>
89 const unsigned int evaluation_point);
119 value(
const number x, std::vector<number> &
values)
const;
139 template <
typename Number2>
142 const unsigned int n_derivatives,
161 scale(
const number factor);
178 template <
typename number2>
180 shift(
const number2 offset);
229 print(std::ostream &out)
const;
236 template <
class Archive>
256 template <
typename number2>
308 template <
typename number>
316 Monomial(
const unsigned int n,
const double coefficient = 1.);
324 static std::vector<Polynomial<number>>
331 static std::vector<number>
332 make_vector(
unsigned int n,
const double coefficient);
369 static std::vector<Polynomial<double>>
379 const unsigned int support_point,
380 std::vector<double> &a);
391 std::vector<Polynomial<double>>
423 static std::vector<Polynomial<double>>
453 Lobatto(
const unsigned int p = 0);
459 static std::vector<Polynomial<double>>
528 static std::vector<Polynomial<double>>
542 static const std::vector<double> &
553 static std::vector<std::unique_ptr<const std::vector<double>>>
600 static std::vector<Polynomial<double>>
715 const unsigned int index);
721 static std::vector<Polynomial<double>>
736 template <
typename Number>
756 template <
typename Number>
770 template <
typename number>
772 : in_lagrange_product_form(false)
773 , lagrange_weight(1.)
778 template <
typename number>
782 if (in_lagrange_product_form ==
true)
784 return lagrange_support_points.size();
789 return coefficients.size() - 1;
795 template <
typename number>
799 if (in_lagrange_product_form ==
false)
804 const unsigned int m = coefficients.size();
805 number value = coefficients.back();
806 for (
int k = m - 2; k >= 0; --k)
807 value = value * x + coefficients[k];
813 const unsigned int m = lagrange_support_points.size();
815 for (
unsigned int j = 0; j < m; ++j)
816 value *= x - lagrange_support_points[j];
817 value *= lagrange_weight;
824 template <
typename number>
825 template <
typename Number2>
828 const unsigned int n_derivatives,
832 if (in_lagrange_product_form ==
true)
837 const unsigned int n_supp = lagrange_support_points.size();
838 const number weight = lagrange_weight;
839 switch (n_derivatives)
843 for (
unsigned int d = 1;
d <= n_derivatives; ++
d)
845 for (
unsigned int i = 0; i < n_supp; ++i)
847 const Number2 v = x - lagrange_support_points[i];
855 for (
unsigned int k = n_derivatives; k > 0; --k)
868 number k_factorial = 1;
869 for (
unsigned int k = 0; k <= n_derivatives; ++k)
871 values[k] *= k_factorial * weight;
872 k_factorial *=
static_cast<number
>(k + 1);
885 for (
unsigned int i = 0; i < n_supp; ++i)
887 const Number2 v = x - lagrange_support_points[i];
890 values[0] = weight * value;
897 Number2 derivative = 0.;
898 for (
unsigned int i = 0; i < n_supp; ++i)
900 const Number2 v = x - lagrange_support_points[i];
901 derivative = derivative * v + value;
904 values[0] = weight * value;
905 values[1] = weight * derivative;
912 Number2 derivative = 0.;
914 for (
unsigned int i = 0; i < n_supp; ++i)
916 const Number2 v = x - lagrange_support_points[i];
918 derivative = derivative * v + value;
921 values[0] = weight * value;
922 values[1] = weight * derivative;
934 const unsigned int m = coefficients.size();
935 std::vector<Number2> a(coefficients.size());
936 std::copy(coefficients.begin(), coefficients.end(), a.begin());
937 unsigned int j_factorial = 1;
942 const unsigned int min_valuessize_m =
std::min(n_derivatives + 1, m);
943 for (
unsigned int j = 0; j < min_valuessize_m; ++j)
945 for (
int k = m - 2; k >=
static_cast<int>(j); --k)
946 a[k] += x * a[k + 1];
947 values[j] =
static_cast<number
>(j_factorial) * a[j];
949 j_factorial *= j + 1;
953 for (
unsigned int j = min_valuessize_m; j <= n_derivatives; ++j)
959 template <
typename number>
960 template <
class Archive>
967 ar &in_lagrange_product_form;
968 ar &lagrange_support_points;
974 template <
typename Number>
981 Assert(alpha >= 0 && beta >= 0,
988 const Number xeval = Number(-1) + 2. * x;
994 p1 = ((alpha + beta + 2) * xeval + (alpha - beta)) / 2;
998 for (
unsigned int i = 1; i < degree; ++i)
1000 const Number v = 2 * i + (alpha + beta);
1001 const Number a1 = 2 * (i + 1) * (i + (alpha + beta + 1)) * v;
1002 const Number a2 = (v + 1) * (alpha * alpha - beta * beta);
1003 const Number a3 = v * (v + 1) * (v + 2);
1004 const Number a4 = 2 * (i + alpha) * (i + beta) * (v + 2);
1006 const Number pn = ((a2 + a3 * xeval) * p1 - a4 * p0) / a1;
1015 template <
typename Number>
1021 std::vector<Number> x(degree, 0.5);
1032 const Number tolerance =
1042 const unsigned int n_points = (alpha == beta ? degree / 2 : degree);
1043 for (
unsigned int k = 0; k < n_points; ++k)
1047 Number r = 0.5 - 0.5 *
std::cos(
static_cast<Number
>(2 * k + 1) /
1050 r = (r + x[k - 1]) / 2;
1053 for (
unsigned int it = 1; it < 1000; ++it)
1056 for (
unsigned int i = 0; i < k; ++i)
1057 s += 1. / (r - x[i]);
1061 (alpha + beta + degree + 1) *
1066 const Number delta = f / (f * s - J_x);
1074 if (it == converged + 1)
1079 ExcMessage(
"Newton iteration for zero of Jacobi polynomial "
1080 "did not converge."));
1086 for (
unsigned int k = n_points; k < degree; ++k)
1087 x[k] = 1.0 - x[degree - k - 1];
HermiteInterpolation(const unsigned int p)
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int p)
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
HermiteLikeInterpolation(const unsigned int degree, const unsigned int index)
static std::vector< std::unique_ptr< const std::vector< double > > > recursive_coefficients
Hierarchical(const unsigned int p)
static void compute_coefficients(const unsigned int p)
static const std::vector< double > & get_coefficients(const unsigned int p)
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
static void compute_coefficients(const unsigned int n, const unsigned int support_point, std::vector< double > &a)
LagrangeEquidistant(const unsigned int n, const unsigned int support_point)
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
Legendre(const unsigned int p)
std::vector< double > compute_coefficients(const unsigned int p)
Lobatto(const unsigned int p=0)
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int p)
static std::vector< Polynomial< number > > generate_complete_basis(const unsigned int degree)
Monomial(const unsigned int n, const double coefficient=1.)
static std::vector< number > make_vector(unsigned int n, const double coefficient)
number value(const number x) const
bool operator==(const Polynomial< number > &p) const
std::vector< number > coefficients
Polynomial< number > primitive() const
Polynomial< number > & operator+=(const Polynomial< number > &p)
Polynomial< number > derivative() const
void transform_into_standard_form()
void scale(const number factor)
Polynomial< number > & operator-=(const Polynomial< number > &p)
std::vector< number > lagrange_support_points
void shift(const number2 offset)
void print(std::ostream &out) const
bool in_lagrange_product_form
void serialize(Archive &ar, const unsigned int version)
static void multiply(std::vector< number > &coefficients, const number factor)
void value(const Number2 x, const unsigned int n_derivatives, Number2 *values) const
Polynomial< number > & operator*=(const double s)
virtual std::size_t memory_consumption() const
unsigned int degree() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcEmptyObject()
#define Assert(cond, exc)
static ::ExceptionBase & ExcMessage(std::string arg1)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > epsilon(const Tensor< 2, dim, Number > &Grad_u)
Number jacobi_polynomial_value(const unsigned int degree, const int alpha, const int beta, const Number x)
std::vector< Polynomial< double > > generate_complete_Lagrange_basis(const std::vector< Point< 1 > > &points)
std::vector< Number > jacobi_polynomial_roots(const unsigned int degree, const int alpha, const int beta)
void copy(const T *begin, const T *end, U *dest)
static constexpr double PI
static const unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > cos(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)