43 std::mutex coefficients_lock;
53 template <
typename number>
56 , in_lagrange_product_form(false)
62 template <
typename number>
64 : coefficients(n + 1, 0.)
65 , in_lagrange_product_form(false)
71 template <
typename number>
74 : in_lagrange_product_form(true)
80 number tmp_lagrange_weight = 1.;
81 for (
unsigned int i = 0; i < supp.size(); ++i)
85 tmp_lagrange_weight *= supp[
center](0) - supp[i](0);
89 Assert(std::fabs(tmp_lagrange_weight) > std::numeric_limits<number>::min(),
90 ExcMessage(
"Underflow in computation of Lagrange denominator."));
91 Assert(std::fabs(tmp_lagrange_weight) < std::numeric_limits<number>::max(),
92 ExcMessage(
"Overflow in computation of Lagrange denominator."));
99 template <
typename number>
105 value(x, values.size() - 1, values.data());
110 template <
typename number>
119 coefficients.resize(lagrange_support_points.size() + 1);
120 if (lagrange_support_points.size() == 0)
121 coefficients[0] = 1.;
124 coefficients[0] = -lagrange_support_points[0];
125 coefficients[1] = 1.;
126 for (
unsigned int i = 1; i < lagrange_support_points.size(); ++i)
128 coefficients[i + 1] = 1.;
129 for (
unsigned int j = i; j > 0; --j)
130 coefficients[j] = (-lagrange_support_points[i] * coefficients[j] +
131 coefficients[j - 1]);
132 coefficients[0] *= -lagrange_support_points[i];
135 for (
unsigned int i = 0; i < lagrange_support_points.size() + 1; ++i)
136 coefficients[i] *= lagrange_weight;
139 std::vector<number> new_points;
140 lagrange_support_points.swap(new_points);
141 in_lagrange_product_form =
false;
142 lagrange_weight = 1.;
147 template <
typename number>
153 for (
typename std::vector<number>::iterator c = coefficients.begin();
154 c != coefficients.end();
164 template <
typename number>
171 if (in_lagrange_product_form ==
true)
173 number inv_fact = number(1.) / factor;
174 number accumulated_fact = 1.;
175 for (
unsigned int i = 0; i < lagrange_support_points.size(); ++i)
177 lagrange_support_points[i] *= inv_fact;
178 accumulated_fact *= factor;
180 lagrange_weight *= accumulated_fact;
184 scale(coefficients, factor);
189 template <
typename number>
194 for (
typename std::vector<number>::iterator c = coefficients.begin();
195 c != coefficients.end();
202 template <
typename number>
206 if (in_lagrange_product_form ==
true)
207 lagrange_weight *= s;
210 for (
typename std::vector<number>::iterator c = coefficients.begin();
211 c != coefficients.end();
220 template <
typename number>
229 lagrange_support_points.insert(lagrange_support_points.end(),
235 else if (in_lagrange_product_form ==
true)
236 transform_into_standard_form();
241 std::unique_ptr<Polynomial<number>> q_data;
245 q_data = std::make_unique<Polynomial<number>>(p);
253 unsigned int new_degree = this->degree() + q->
degree();
255 std::vector<number> new_coefficients(new_degree + 1, 0.);
257 for (
unsigned int i = 0; i < q->
coefficients.size(); ++i)
258 for (
unsigned int j = 0; j < this->coefficients.size(); ++j)
259 new_coefficients[i + j] += this->coefficients[j] * q->
coefficients[i];
260 this->coefficients = new_coefficients;
267 template <
typename number>
279 if (in_lagrange_product_form ==
true)
280 transform_into_standard_form();
285 std::unique_ptr<Polynomial<number>> q_data;
289 q_data = std::make_unique<Polynomial<number>>(p);
301 for (
unsigned int i = 0; i < q->
coefficients.size(); ++i)
309 template <
typename number>
315 if (in_lagrange_product_form ==
true)
316 transform_into_standard_form();
321 std::unique_ptr<Polynomial<number>> q_data;
325 q_data = std::make_unique<Polynomial<number>>(p);
337 for (
unsigned int i = 0; i < q->
coefficients.size(); ++i)
345 template <
typename number>
357 else if (in_lagrange_product_form ==
true)
375 template <
typename number>
376 template <
typename number2>
379 const number2 offset)
387 std::vector<number2> new_coefficients(coefficients.begin(),
393 for (
unsigned int d = 1; d < new_coefficients.size(); ++d)
395 const unsigned int n = d;
400 unsigned int binomial_coefficient = 1;
405 number2 offset_power = offset;
412 for (
unsigned int k = 0; k < d; ++k)
418 binomial_coefficient = (binomial_coefficient * (n - k)) / (k + 1);
420 new_coefficients[d - k - 1] +=
421 new_coefficients[d] * binomial_coefficient * offset_power;
422 offset_power *= offset;
432 coefficients.assign(new_coefficients.begin(), new_coefficients.end());
437 template <
typename number>
438 template <
typename number2>
445 if (in_lagrange_product_form ==
true)
447 for (
unsigned int i = 0; i < lagrange_support_points.size(); ++i)
448 lagrange_support_points[i] -= offset;
452 shift(coefficients, offset);
457 template <
typename number>
466 std::unique_ptr<Polynomial<number>> q_data;
468 if (in_lagrange_product_form ==
true)
470 q_data = std::make_unique<Polynomial<number>>(*this);
477 std::vector<number> newcoefficients(q->
coefficients.size() - 1);
478 for (
unsigned int i = 1; i < q->
coefficients.size(); ++i)
479 newcoefficients[i - 1] = number(i) * q->
coefficients[i];
486 template <
typename number>
492 std::unique_ptr<Polynomial<number>> q_data;
494 if (in_lagrange_product_form ==
true)
496 q_data = std::make_unique<Polynomial<number>>(*this);
503 std::vector<number> newcoefficients(q->
coefficients.size() + 1);
504 newcoefficients[0] = 0.;
505 for (
unsigned int i = 0; i < q->
coefficients.size(); ++i)
506 newcoefficients[i + 1] = q->
coefficients[i] / number(i + 1.);
513 template <
typename number>
517 if (in_lagrange_product_form ==
true)
519 out << lagrange_weight;
520 for (
unsigned int i = 0; i < lagrange_support_points.size(); ++i)
521 out <<
" (x-" << lagrange_support_points[i] <<
")";
525 for (
int i = degree(); i >= 0; --i)
527 out << coefficients[i] <<
" x^" << i << std::endl;
532 template <
typename number>
546 template <
typename number>
550 std::vector<number> result(n + 1, 0.);
551 result[n] = coefficient;
557 template <
typename number>
559 :
Polynomial<number>(make_vector(n, coefficient))
564 template <
typename number>
565 std::vector<Polynomial<number>>
568 std::vector<Polynomial<number>> v;
569 v.reserve(degree + 1);
570 for (
unsigned int i = 0; i <= degree; ++i)
581 namespace LagrangeEquidistantImplementation
583 std::vector<Point<1>>
586 std::vector<Point<1>> points(n + 1);
587 const double one_over_n = 1. / n;
588 for (
unsigned int k = 0; k <= n; ++k)
589 points[k](0) =
static_cast<double>(k) * one_over_n;
598 const unsigned int support_point)
600 generate_equidistant_unit_points(n),
611 std::vector<double> new_support_points;
622 const unsigned int support_point,
623 std::vector<double> &a)
627 unsigned int n_functions = n + 1;
629 double const *x =
nullptr;
635 static const double x1[4] = {1.0, -1.0, 0.0, 1.0};
641 static const double x2[9] = {
642 1.0, -3.0, 2.0, 0.0, 4.0, -4.0, 0.0, -1.0, 2.0};
648 static const double x3[16] = {1.0,
672 for (
unsigned int i = 0; i < n_functions; ++i)
673 a[i] = x[support_point * n_functions + i];
678 std::vector<Polynomial<double>>
683 return std::vector<Polynomial<double>>(
689 std::vector<Polynomial<double>> v;
690 for (
unsigned int i = 0; i <=
degree; ++i)
701 std::vector<Polynomial<double>>
704 std::vector<Polynomial<double>> p;
705 p.reserve(points.size());
707 for (
unsigned int i = 0; i < points.size(); ++i)
708 p.emplace_back(points, i);
730 for (
unsigned int i = 0; i < k; ++i)
737 for (
unsigned int i = 0; i < k; ++i)
744 std::vector<Polynomial<double>>
747 std::vector<Polynomial<double>> v;
749 for (
unsigned int i = 0; i <=
degree; ++i)
799 std::vector<double> legendre_coefficients_tmp1(p);
800 std::vector<double> legendre_coefficients_tmp2(p - 1);
804 legendre_coefficients_tmp1[0] = 1.0;
806 for (
unsigned int i = 2; i < p; ++i)
808 for (
unsigned int j = 0; j < i - 1; ++j)
809 legendre_coefficients_tmp2[j] = legendre_coefficients_tmp1[j];
811 for (
unsigned int j = 0; j < i; ++j)
816 ((1.0 - 2 * i) * legendre_coefficients_tmp1[0] /
818 (1.0 - i) * legendre_coefficients_tmp2[0] /
822 for (
unsigned int j = 1; j < i - 1; ++j)
826 (2.0 * legendre_coefficients_tmp1[j - 1] -
827 legendre_coefficients_tmp1[j]) +
828 (1.0 - i) * legendre_coefficients_tmp2[j] /
833 (2.0 * legendre_coefficients_tmp1[i - 2] -
834 legendre_coefficients_tmp1[i - 1]) /
837 legendre_coefficients_tmp1[i - 1] / i;
840 for (
int i = p; i > 0; --i)
849 std::vector<Polynomial<double>>
852 std::vector<Polynomial<double>> basis(p + 1);
854 for (
unsigned int i = 0; i <= p; ++i)
867 std::vector<std::unique_ptr<const std::vector<double>>>
888 std::lock_guard<std::mutex> lock(coefficients_lock);
931 std::vector<double> c0(2);
935 std::vector<double> c1(2);
942 std::make_unique<const std::vector<double>>(std::move(c0));
944 std::make_unique<const std::vector<double>>(std::move(c1));
948 coefficients_lock.unlock();
950 coefficients_lock.lock();
952 std::vector<double> c2(3);
961 std::make_unique<const std::vector<double>>(std::move(c2));
973 coefficients_lock.unlock();
975 coefficients_lock.lock();
977 std::vector<double> ck(k + 1);
983 for (
unsigned int i = 1; i <= k - 1; ++i)
1005 std::make_unique<const std::vector<double>>(std::move(ck));
1012 const std::vector<double> &
1022 std::lock_guard<std::mutex> lock(coefficients_lock);
1028 std::vector<Polynomial<double>>
1041 return std::vector<Polynomial<double>>(
1045 std::vector<Polynomial<double>> v;
1047 for (
unsigned int i = 0; i <=
degree; ++i)
1102 (*this) *= legendre;
1108 std::vector<Polynomial<double>>
1113 "degrees less than three"));
1114 std::vector<Polynomial<double>> basis(n + 1);
1116 for (
unsigned int i = 0; i <= n; ++i)
1132 find_support_point_x_star(
const std::vector<double> &jacobi_roots)
1138 double guess_left = 0;
1139 double guess_right = 0.5;
1140 const unsigned int degree = jacobi_roots.size() + 3;
1152 double integral_left = 0, integral_right = 0;
1153 for (
unsigned int q = 0; q < gauss.size(); ++q)
1155 const double x = gauss.point(q)[0];
1156 double poly_val_common = x;
1157 for (
unsigned int j = 0; j < degree - 3; ++j)
1158 poly_val_common *= Utilities::fixed_power<2>(x - jacobi_roots[j]);
1159 poly_val_common *= Utilities::fixed_power<4>(x - 1.);
1161 gauss.weight(q) * (poly_val_common * (x - guess_left));
1163 gauss.weight(q) * (poly_val_common * (x - guess_right));
1168 return guess_right - (guess_right - guess_left) /
1169 (integral_right - integral_left) * integral_right;
1176 const unsigned int index)
1188 else if (degree == 1)
1201 else if (degree == 2)
1209 else if (index == 1)
1222 else if (degree == 3)
1239 else if (index == 1)
1249 else if (index == 2)
1256 else if (index == 3)
1287 std::vector<double> jacobi_roots =
1288 jacobi_polynomial_roots<double>(
degree - 3, 4, 4);
1297 const double auxiliary_zero =
1298 find_support_point_x_star(jacobi_roots);
1300 for (
unsigned int m = 0; m <
degree - 3; ++m)
1308 else if (index == 1)
1311 for (
unsigned int m = 0; m <
degree - 3; ++m)
1324 std::vector<Point<1>> points(
degree);
1326 for (
unsigned int i = 0; i <
degree; ++i)
1333 std::vector<double> value_and_grad(2);
1334 helper.
value(0., value_and_grad);
1338 const double auxiliary_zero =
1339 find_support_point_x_star(jacobi_roots);
1341 (1. / auxiliary_zero - value_and_grad[1] / value_and_grad[0]) /
1344 else if (index >= 2 && index <
degree - 1)
1348 for (
unsigned int m = 0, c = 2; m <
degree - 3; ++m)
1358 else if (index ==
degree - 1)
1362 for (
unsigned int m = 0; m <
degree - 3; ++m)
1366 std::vector<Point<1>> points(
degree);
1368 for (
unsigned int i = 0; i <
degree; ++i)
1375 std::vector<double> value_and_grad(2);
1376 helper.
value(1., value_and_grad);
1380 const double auxiliary_zero =
1381 find_support_point_x_star(jacobi_roots);
1383 (-1. / auxiliary_zero - value_and_grad[1] / value_and_grad[0]) /
1386 else if (index ==
degree)
1388 const double auxiliary_zero =
1389 find_support_point_x_star(jacobi_roots);
1392 for (
unsigned int m = 0; m <
degree - 3; ++m)
1404 std::vector<Polynomial<double>>
1407 std::vector<Polynomial<double>> basis(
degree + 1);
1409 for (
unsigned int i = 0; i <=
degree; ++i)
1422 template class Polynomial<float>;
1423 template class Polynomial<double>;
1424 template class Polynomial<long double>;
1439 template class Monomial<float>;
1440 template class Monomial<double>;
1441 template class Monomial<long double>;
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
static void multiply(std::vector< number > &coefficients, const number factor)
Polynomial< number > & operator*=(const double s)
virtual std::size_t memory_consumption() const
unsigned int degree() const
const std::vector< Point< dim > > & get_points() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcZero()
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcEmptyObject()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
std::enable_if_t< std::is_fundamental< T >::value, std::size_t > memory_consumption(const T &t)
std::vector< Point< 1 > > generate_equidistant_unit_points(const unsigned int n)
std::vector< Polynomial< double > > generate_complete_Lagrange_basis(const std::vector< Point< 1 > > &points)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)