42 std::mutex coefficients_lock;
52 template <
typename number>
55 , in_lagrange_product_form(false)
61 template <
typename number>
63 : coefficients(n + 1, 0.)
64 , in_lagrange_product_form(false)
70 template <
typename number>
73 : in_lagrange_product_form(true)
78 lagrange_support_points.reserve(supp.size() - 1);
79 number tmp_lagrange_weight = 1.;
80 for (
unsigned int i = 0; i < supp.size(); ++i)
83 lagrange_support_points.push_back(supp[i](0));
84 tmp_lagrange_weight *= supp[
center](0) - supp[i](0);
89 ExcMessage(
"Underflow in computation of Lagrange denominator."));
91 ExcMessage(
"Overflow in computation of Lagrange denominator."));
93 lagrange_weight =
static_cast<number
>(1.) / tmp_lagrange_weight;
98 template <
typename number>
104 value(x, values.size() - 1, values.data());
109 template <
typename number>
112 const unsigned int n_derivatives,
113 number * values)
const
116 if (in_lagrange_product_form ==
true)
121 const unsigned int n_supp = lagrange_support_points.size();
122 switch (n_derivatives)
126 for (
unsigned int d = 1;
d <= n_derivatives; ++
d)
128 for (
unsigned int i = 0; i < n_supp; ++i)
130 const number v = x - lagrange_support_points[i];
138 for (
unsigned int k = n_derivatives; k > 0; --k)
139 values[k] = (values[k] * v + values[k - 1]);
151 number k_faculty = 1;
152 for (
unsigned int k = 0; k <= n_derivatives; ++k)
154 values[k] *= k_faculty * lagrange_weight;
155 k_faculty *=
static_cast<number
>(k + 1);
165 for (
unsigned int i = 0; i < n_supp; ++i)
167 const number v = x - lagrange_support_points[i];
170 values[0] *= lagrange_weight;
176 for (
unsigned int i = 0; i < n_supp; ++i)
178 const number v = x - lagrange_support_points[i];
179 values[1] = values[1] * v + values[0];
182 values[0] *= lagrange_weight;
183 values[1] *= lagrange_weight;
190 for (
unsigned int i = 0; i < n_supp; ++i)
192 const number v = x - lagrange_support_points[i];
193 values[2] = values[2] * v + values[1];
194 values[1] = values[1] * v + values[0];
197 values[0] *= lagrange_weight;
198 values[1] *= lagrange_weight;
199 values[2] *=
static_cast<number
>(2) * lagrange_weight;
210 if (n_derivatives == 0)
212 values[0] =
value(x);
218 const unsigned int m = coefficients.size();
219 std::vector<number> a(coefficients);
220 unsigned int j_faculty = 1;
225 const unsigned int min_valuessize_m =
std::min(n_derivatives + 1, m);
226 for (
unsigned int j = 0; j < min_valuessize_m; ++j)
228 for (
int k = m - 2; k >=
static_cast<int>(j); --k)
229 a[k] += x * a[k + 1];
230 values[j] =
static_cast<number
>(j_faculty) * a[j];
236 for (
unsigned int j = min_valuessize_m; j <= n_derivatives; ++j)
242 template <
typename number>
251 coefficients.resize(lagrange_support_points.size() + 1);
252 if (lagrange_support_points.size() == 0)
253 coefficients[0] = 1.;
256 coefficients[0] = -lagrange_support_points[0];
257 coefficients[1] = 1.;
258 for (
unsigned int i = 1; i < lagrange_support_points.size(); ++i)
260 coefficients[i + 1] = 1.;
261 for (
unsigned int j = i; j > 0; --j)
262 coefficients[j] = (-lagrange_support_points[i] * coefficients[j] +
263 coefficients[j - 1]);
264 coefficients[0] *= -lagrange_support_points[i];
267 for (
unsigned int i = 0; i < lagrange_support_points.size() + 1; ++i)
268 coefficients[i] *= lagrange_weight;
271 std::vector<number> new_points;
272 lagrange_support_points.swap(new_points);
273 in_lagrange_product_form =
false;
274 lagrange_weight = 1.;
279 template <
typename number>
285 for (
typename std::vector<number>::iterator c = coefficients.begin();
286 c != coefficients.end();
296 template <
typename number>
303 if (in_lagrange_product_form ==
true)
305 number inv_fact = number(1.) / factor;
306 number accumulated_fact = 1.;
307 for (
unsigned int i = 0; i < lagrange_support_points.size(); ++i)
309 lagrange_support_points[i] *= inv_fact;
310 accumulated_fact *= factor;
312 lagrange_weight *= accumulated_fact;
316 scale(coefficients, factor);
321 template <
typename number>
326 for (
typename std::vector<number>::iterator c = coefficients.begin();
327 c != coefficients.end();
334 template <
typename number>
338 if (in_lagrange_product_form ==
true)
339 lagrange_weight *= s;
342 for (
typename std::vector<number>::iterator c = coefficients.begin();
343 c != coefficients.end();
352 template <
typename number>
361 lagrange_support_points.insert(lagrange_support_points.end(),
367 else if (in_lagrange_product_form ==
true)
368 transform_into_standard_form();
373 std::unique_ptr<Polynomial<number>> q_data;
377 q_data = std_cxx14::make_unique<Polynomial<number>>(p);
378 q_data->transform_into_standard_form();
385 unsigned int new_degree = this->degree() + q->
degree();
387 std::vector<number> new_coefficients(new_degree + 1, 0.);
389 for (
unsigned int i = 0; i < q->
coefficients.size(); ++i)
390 for (
unsigned int j = 0; j < this->coefficients.size(); ++j)
391 new_coefficients[i + j] += this->coefficients[j] * q->
coefficients[i];
392 this->coefficients = new_coefficients;
399 template <
typename number>
411 if (in_lagrange_product_form ==
true)
412 transform_into_standard_form();
417 std::unique_ptr<Polynomial<number>> q_data;
421 q_data = std_cxx14::make_unique<Polynomial<number>>(p);
422 q_data->transform_into_standard_form();
433 for (
unsigned int i = 0; i < q->
coefficients.size(); ++i)
441 template <
typename number>
447 if (in_lagrange_product_form ==
true)
448 transform_into_standard_form();
453 std::unique_ptr<Polynomial<number>> q_data;
457 q_data = std_cxx14::make_unique<Polynomial<number>>(p);
458 q_data->transform_into_standard_form();
469 for (
unsigned int i = 0; i < q->
coefficients.size(); ++i)
477 template <
typename number>
489 else if (in_lagrange_product_form ==
true)
507 template <
typename number>
508 template <
typename number2>
511 const number2 offset)
519 std::vector<number2> new_coefficients(coefficients.begin(),
525 for (
unsigned int d = 1;
d < new_coefficients.size(); ++
d)
527 const unsigned int n =
d;
532 unsigned int binomial_coefficient = 1;
537 number2 offset_power = offset;
544 for (
unsigned int k = 0; k <
d; ++k)
550 binomial_coefficient = (binomial_coefficient * (n - k)) / (k + 1);
552 new_coefficients[
d - k - 1] +=
553 new_coefficients[
d] * binomial_coefficient * offset_power;
554 offset_power *= offset;
564 coefficients.assign(new_coefficients.begin(), new_coefficients.end());
569 template <
typename number>
570 template <
typename number2>
577 if (in_lagrange_product_form ==
true)
579 for (
unsigned int i = 0; i < lagrange_support_points.size(); ++i)
580 lagrange_support_points[i] -= offset;
584 shift(coefficients, offset);
589 template <
typename number>
598 std::unique_ptr<Polynomial<number>> q_data;
600 if (in_lagrange_product_form ==
true)
602 q_data = std_cxx14::make_unique<Polynomial<number>>(*this);
603 q_data->transform_into_standard_form();
609 std::vector<number> newcoefficients(q->
coefficients.size() - 1);
610 for (
unsigned int i = 1; i < q->
coefficients.size(); ++i)
611 newcoefficients[i - 1] = number(i) * q->
coefficients[i];
618 template <
typename number>
624 std::unique_ptr<Polynomial<number>> q_data;
626 if (in_lagrange_product_form ==
true)
628 q_data = std_cxx14::make_unique<Polynomial<number>>(*this);
629 q_data->transform_into_standard_form();
635 std::vector<number> newcoefficients(q->
coefficients.size() + 1);
636 newcoefficients[0] = 0.;
637 for (
unsigned int i = 0; i < q->
coefficients.size(); ++i)
638 newcoefficients[i + 1] = q->
coefficients[i] / number(i + 1.);
645 template <
typename number>
649 if (in_lagrange_product_form ==
true)
651 out << lagrange_weight;
652 for (
unsigned int i = 0; i < lagrange_support_points.size(); ++i)
653 out <<
" (x-" << lagrange_support_points[i] <<
")";
657 for (
int i = degree(); i >= 0; --i)
659 out << coefficients[i] <<
" x^" << i << std::endl;
664 template <
typename number>
678 template <
typename number>
682 std::vector<number> result(n + 1, 0.);
683 result[n] = coefficient;
689 template <
typename number>
691 :
Polynomial<number>(make_vector(n, coefficient))
696 template <
typename number>
697 std::vector<Polynomial<number>>
700 std::vector<Polynomial<number>> v;
701 v.reserve(degree + 1);
702 for (
unsigned int i = 0; i <= degree; ++i)
713 namespace LagrangeEquidistantImplementation
715 std::vector<Point<1>>
718 std::vector<Point<1>> points(n + 1);
719 const double one_over_n = 1. / n;
720 for (
unsigned int k = 0; k <= n; ++k)
721 points[k](0) =
static_cast<double>(k) * one_over_n;
730 const unsigned int support_point)
743 std::vector<double> new_support_points;
754 const unsigned int support_point,
755 std::vector<double> &a)
759 unsigned int n_functions = n + 1;
761 double const *x =
nullptr;
767 static const double x1[4] = {1.0, -1.0, 0.0, 1.0};
773 static const double x2[9] = {
774 1.0, -3.0, 2.0, 0.0, 4.0, -4.0, 0.0, -1.0, 2.0};
780 static const double x3[16] = {1.0,
804 for (
unsigned int i = 0; i < n_functions; ++i)
805 a[i] = x[support_point * n_functions + i];
810 std::vector<Polynomial<double>>
815 return std::vector<Polynomial<double>>(
821 std::vector<Polynomial<double>> v;
822 for (
unsigned int i = 0; i <=
degree; ++i)
833 std::vector<Polynomial<double>>
836 std::vector<Polynomial<double>> p;
837 p.reserve(points.size());
839 for (
unsigned int i = 0; i < points.size(); ++i)
840 p.emplace_back(points, i);
862 for (
unsigned int i = 0; i < k; ++i)
869 for (
unsigned int i = 0; i < k; ++i)
876 std::vector<Polynomial<double>>
879 std::vector<Polynomial<double>> v;
881 for (
unsigned int i = 0; i <=
degree; ++i)
931 std::vector<double> legendre_coefficients_tmp1(p);
932 std::vector<double> legendre_coefficients_tmp2(p - 1);
936 legendre_coefficients_tmp1[0] = 1.0;
938 for (
unsigned int i = 2; i < p; ++i)
940 for (
unsigned int j = 0; j < i - 1; ++j)
941 legendre_coefficients_tmp2[j] = legendre_coefficients_tmp1[j];
943 for (
unsigned int j = 0; j < i; ++j)
947 std::sqrt(2 * i + 1.) *
948 ((1.0 - 2 * i) * legendre_coefficients_tmp1[0] /
949 std::sqrt(2 * i - 1.) +
950 (1.0 - i) * legendre_coefficients_tmp2[0] /
951 std::sqrt(2 * i - 3.)) /
954 for (
unsigned int j = 1; j < i - 1; ++j)
956 std::sqrt(2 * i + 1.) *
957 (std::sqrt(2 * i - 1.) *
958 (2.0 * legendre_coefficients_tmp1[j - 1] -
959 legendre_coefficients_tmp1[j]) +
960 (1.0 - i) * legendre_coefficients_tmp2[j] /
961 std::sqrt(2 * i - 3.)) /
965 (2.0 * legendre_coefficients_tmp1[i - 2] -
966 legendre_coefficients_tmp1[i - 1]) /
969 legendre_coefficients_tmp1[i - 1] / i;
972 for (
int i = p; i > 0; --i)
981 std::vector<Polynomial<double>>
984 std::vector<Polynomial<double>> basis(p + 1);
986 for (
unsigned int i = 0; i <= p; ++i)
999 std::vector<std::unique_ptr<const std::vector<double>>>
1013 unsigned int k = k_;
1020 std::lock_guard<std::mutex> lock(coefficients_lock);
1063 std::vector<double> c0(2);
1067 std::vector<double> c1(2);
1074 std_cxx14::make_unique<const std::vector<double>>(std::move(c0));
1076 std_cxx14::make_unique<const std::vector<double>>(std::move(c1));
1080 coefficients_lock.unlock();
1082 coefficients_lock.lock();
1084 std::vector<double> c2(3);
1086 const double a = 1.;
1093 std_cxx14::make_unique<const std::vector<double>>(std::move(c2));
1105 coefficients_lock.unlock();
1107 coefficients_lock.lock();
1109 std::vector<double> ck(k + 1);
1111 const double a = 1.;
1115 for (
unsigned int i = 1; i <= k - 1; ++i)
1137 std_cxx14::make_unique<const std::vector<double>>(std::move(ck));
1144 const std::vector<double> &
1154 std::lock_guard<std::mutex> lock(coefficients_lock);
1160 std::vector<Polynomial<double>>
1173 return std::vector<Polynomial<double>>(
1177 std::vector<Polynomial<double>> v;
1179 for (
unsigned int i = 0; i <=
degree; ++i)
1234 (*this) *= legendre;
1240 std::vector<Polynomial<double>>
1245 "degrees less than three"));
1246 std::vector<Polynomial<double>> basis(n + 1);
1248 for (
unsigned int i = 0; i <= n; ++i)
1264 find_support_point_x_star(
const std::vector<double> &jacobi_roots)
1270 double guess_left = 0;
1271 double guess_right = 0.5;
1272 const unsigned int degree = jacobi_roots.size() + 3;
1284 double integral_left = 0, integral_right = 0;
1285 for (
unsigned int q = 0; q < gauss.size(); ++q)
1287 const double x = gauss.point(q)[0];
1288 double poly_val_common = x;
1289 for (
unsigned int j = 0; j < degree - 3; ++j)
1290 poly_val_common *= Utilities::fixed_power<2>(x - jacobi_roots[j]);
1291 poly_val_common *= Utilities::fixed_power<4>(x - 1.);
1293 gauss.weight(q) * (poly_val_common * (x - guess_left));
1295 gauss.weight(q) * (poly_val_common * (x - guess_right));
1300 return guess_right - (guess_right - guess_left) /
1301 (integral_right - integral_left) * integral_right;
1308 const unsigned int index)
1320 else if (degree == 1)
1333 else if (degree == 2)
1341 else if (index == 1)
1354 else if (degree == 3)
1371 else if (index == 1)
1381 else if (index == 2)
1388 else if (index == 3)
1419 std::vector<double> jacobi_roots =
1420 jacobi_polynomial_roots<double>(
degree - 3, 4, 4);
1429 const double auxiliary_zero =
1430 find_support_point_x_star(jacobi_roots);
1432 for (
unsigned int m = 0; m <
degree - 3; m++)
1440 else if (index == 1)
1443 for (
unsigned int m = 0; m <
degree - 3; m++)
1456 std::vector<Point<1>> points(
degree);
1458 for (
unsigned int i = 0; i <
degree; ++i)
1465 std::vector<double> value_and_grad(2);
1466 helper.
value(0., value_and_grad);
1467 Assert(std::abs(value_and_grad[0]) > 1
e-10,
1470 const double auxiliary_zero =
1471 find_support_point_x_star(jacobi_roots);
1473 (1. / auxiliary_zero - value_and_grad[1] / value_and_grad[0]) /
1476 else if (index >= 2 && index <
degree - 1)
1480 for (
unsigned int m = 0, c = 2; m <
degree - 3; m++)
1490 else if (index ==
degree - 1)
1494 for (
unsigned int m = 0; m <
degree - 3; m++)
1498 std::vector<Point<1>> points(
degree);
1500 for (
unsigned int i = 0; i <
degree; ++i)
1507 std::vector<double> value_and_grad(2);
1508 helper.
value(1., value_and_grad);
1509 Assert(std::abs(value_and_grad[0]) > 1
e-10,
1512 const double auxiliary_zero =
1513 find_support_point_x_star(jacobi_roots);
1515 (-1. / auxiliary_zero - value_and_grad[1] / value_and_grad[0]) /
1518 else if (index ==
degree)
1520 const double auxiliary_zero =
1521 find_support_point_x_star(jacobi_roots);
1524 for (
unsigned int m = 0; m <
degree - 3; m++)
1536 std::vector<Polynomial<double>>
1539 std::vector<Polynomial<double>> basis(
degree + 1);
1541 for (
unsigned int i = 0; i <=
degree; ++i)
1554 template class Polynomial<float>;
1555 template class Polynomial<double>;
1556 template class Polynomial<long double>;
1571 template class Monomial<float>;
1572 template class Monomial<double>;
1573 template class Monomial<long double>;