|
Reference documentation for deal.II version 9.2.0
|
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\)
\(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\)
\(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\)
\(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Go to the documentation of this file.
16 #ifndef dealii_polynomial_h
17 #define dealii_polynomial_h
61 template <
typename number>
88 const unsigned int evaluation_point);
103 value(
const number x)
const;
116 value(
const number x, std::vector<number> &values)
const;
130 value(
const number x,
131 const unsigned int n_derivatives,
132 number * values)
const;
150 scale(
const number factor);
167 template <
typename number2>
169 shift(
const number2 offset);
218 print(std::ostream &out)
const;
224 template <
class Archive>
226 serialize(Archive &ar,
const unsigned int version);
244 template <
typename number2>
298 template <
typename number>
306 Monomial(
const unsigned int n,
const double coefficient = 1.);
314 static std::vector<Polynomial<number>>
321 static std::vector<number>
322 make_vector(
unsigned int n,
const double coefficient);
361 static std::vector<Polynomial<double>>
371 const unsigned int support_point,
372 std::vector<double> &a);
383 std::vector<Polynomial<double>>
417 static std::vector<Polynomial<double>>
449 Lobatto(
const unsigned int p = 0);
455 static std::vector<Polynomial<double>>
526 static std::vector<Polynomial<double>>
540 static const std::vector<double> &
551 static std::vector<std::unique_ptr<const std::vector<double>>>
601 static std::vector<Polynomial<double>>
719 const unsigned int index);
725 static std::vector<Polynomial<double>>
740 template <
typename Number>
760 template <
typename Number>
774 template <
typename number>
776 : in_lagrange_product_form(false)
777 , lagrange_weight(1.)
782 template <
typename number>
786 if (in_lagrange_product_form ==
true)
788 return lagrange_support_points.size();
793 return coefficients.size() - 1;
799 template <
typename number>
803 if (in_lagrange_product_form ==
false)
808 const unsigned int m = coefficients.size();
809 number
value = coefficients.back();
810 for (
int k = m - 2; k >= 0; --k)
817 const unsigned int m = lagrange_support_points.size();
819 for (
unsigned int j = 0; j < m; ++j)
820 value *= x - lagrange_support_points[j];
821 value *= lagrange_weight;
828 template <
typename number>
829 template <
class Archive>
836 ar &in_lagrange_product_form;
837 ar &lagrange_support_points;
843 template <
typename Number>
850 Assert(alpha >= 0 && beta >= 0,
857 const Number xeval = Number(-1) + 2. * x;
863 p1 = ((alpha + beta + 2) * xeval + (alpha - beta)) / 2;
867 for (
unsigned int i = 1; i < degree; ++i)
869 const Number v = 2 * i + (alpha + beta);
870 const Number a1 = 2 * (i + 1) * (i + (alpha + beta + 1)) * v;
871 const Number a2 = (v + 1) * (alpha * alpha - beta * beta);
872 const Number a3 = v * (v + 1) * (v + 2);
873 const Number a4 = 2 * (i + alpha) * (i + beta) * (v + 2);
875 const Number pn = ((a2 + a3 * xeval) * p1 - a4 * p0) / a1;
884 template <
typename Number>
890 std::vector<Number> x(degree, 0.5);
901 const Number tolerance =
911 const unsigned int n_points = (alpha == beta ? degree / 2 : degree);
912 for (
unsigned int k = 0; k < n_points; ++k)
916 Number r = 0.5 - 0.5 * std::cos(
static_cast<Number
>(2 * k + 1) /
919 r = (r + x[k - 1]) / 2;
922 for (
unsigned int it = 1; it < 1000; ++it)
925 for (
unsigned int i = 0; i < k; ++i)
926 s += 1. / (r - x[i]);
930 (alpha + beta + degree + 1) *
935 const Number delta = f / (f * s - J_x);
938 std::abs(delta) < tolerance)
943 if (it == converged + 1)
948 ExcMessage(
"Newton iteration for zero of Jacobi polynomial "
949 "did not converge."));
955 for (
unsigned int k = n_points; k < degree; ++k)
956 x[k] = 1.0 - x[degree - k - 1];
SymmetricTensor< 2, dim, Number > epsilon(const Tensor< 2, dim, Number > &Grad_u)
static const std::vector< double > & get_coefficients(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 p)
static ::ExceptionBase & ExcNotImplemented()
void serialize(Archive &ar, const unsigned int version)
static void compute_coefficients(const unsigned int p)
void transform_into_standard_form()
Legendre(const unsigned int p)
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
Monomial(const unsigned int n, const double coefficient=1.)
virtual std::size_t memory_consumption() const
std::vector< double > compute_coefficients(const unsigned int p)
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
std::vector< Polynomial< double > > generate_complete_Lagrange_basis(const std::vector< Point< 1 >> &points)
void scale(const number factor)
std::vector< Number > jacobi_polynomial_roots(const unsigned int degree, const int alpha, const int beta)
Number jacobi_polynomial_value(const unsigned int degree, const int alpha, const int beta, const Number x)
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
std::vector< number > lagrange_support_points
HermiteInterpolation(const unsigned int p)
static std::vector< Polynomial< number > > generate_complete_basis(const unsigned int degree)
static ::ExceptionBase & ExcMessage(std::string arg1)
static std::vector< number > make_vector(unsigned int n, const double coefficient)
static ::ExceptionBase & ExcEmptyObject()
#define DEAL_II_NAMESPACE_OPEN
LagrangeEquidistant(const unsigned int n, const unsigned int support_point)
Hierarchical(const unsigned int p)
Polynomial< number > & operator+=(const Polynomial< number > &p)
void print(std::ostream &out) const
Polynomial< number > & operator-=(const Polynomial< number > &p)
#define Assert(cond, exc)
Polynomial< number > & operator*=(const double s)
bool operator==(const Polynomial< number > &p) const
Lobatto(const unsigned int p=0)
void shift(const number2 offset)
static void multiply(std::vector< number > &coefficients, const number factor)
static const unsigned int invalid_unsigned_int
unsigned int degree() const
HermiteLikeInterpolation(const unsigned int degree, const unsigned int index)
static void compute_coefficients(const unsigned int n, const unsigned int support_point, std::vector< double > &a)
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
#define DEAL_II_NAMESPACE_CLOSE
static constexpr double PI
static std::vector< std::unique_ptr< const std::vector< double > > > recursive_coefficients
Polynomial< number > primitive() const
std::vector< number > coefficients
Polynomial< number > derivative() const
T max(const T &t, const MPI_Comm &mpi_communicator)
bool in_lagrange_product_form
number value(const number x) const