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>
109 template <
typename number>
118 coefficients.resize(lagrange_support_points.size() + 1);
119 if (lagrange_support_points.size() == 0)
120 coefficients[0] = 1.;
123 coefficients[0] = -lagrange_support_points[0];
124 coefficients[1] = 1.;
125 for (
unsigned int i = 1; i < lagrange_support_points.size(); ++i)
127 coefficients[i + 1] = 1.;
128 for (
unsigned int j = i; j > 0; --j)
129 coefficients[j] = (-lagrange_support_points[i] * coefficients[j] +
130 coefficients[j - 1]);
131 coefficients[0] *= -lagrange_support_points[i];
134 for (
unsigned int i = 0; i < lagrange_support_points.size() + 1; ++i)
135 coefficients[i] *= lagrange_weight;
138 std::vector<number> new_points;
139 lagrange_support_points.swap(new_points);
140 in_lagrange_product_form =
false;
141 lagrange_weight = 1.;
146 template <
typename number>
152 for (
typename std::vector<number>::iterator c = coefficients.begin();
153 c != coefficients.end();
163 template <
typename number>
170 if (in_lagrange_product_form ==
true)
172 number inv_fact = number(1.) / factor;
173 number accumulated_fact = 1.;
174 for (
unsigned int i = 0; i < lagrange_support_points.size(); ++i)
176 lagrange_support_points[i] *= inv_fact;
177 accumulated_fact *= factor;
179 lagrange_weight *= accumulated_fact;
183 scale(coefficients, factor);
188 template <
typename number>
193 for (
typename std::vector<number>::iterator c = coefficients.begin();
194 c != coefficients.end();
201 template <
typename number>
205 if (in_lagrange_product_form ==
true)
206 lagrange_weight *= s;
209 for (
typename std::vector<number>::iterator c = coefficients.begin();
210 c != coefficients.end();
219 template <
typename number>
228 lagrange_support_points.insert(lagrange_support_points.end(),
234 else if (in_lagrange_product_form ==
true)
235 transform_into_standard_form();
240 std::unique_ptr<Polynomial<number>> q_data;
244 q_data = std::make_unique<Polynomial<number>>(p);
252 unsigned int new_degree = this->degree() + q->
degree();
254 std::vector<number> new_coefficients(new_degree + 1, 0.);
256 for (
unsigned int i = 0; i < q->
coefficients.size(); ++i)
257 for (
unsigned int j = 0; j < this->coefficients.size(); ++j)
259 this->coefficients = new_coefficients;
266 template <
typename number>
278 if (in_lagrange_product_form ==
true)
279 transform_into_standard_form();
284 std::unique_ptr<Polynomial<number>> q_data;
288 q_data = std::make_unique<Polynomial<number>>(p);
300 for (
unsigned int i = 0; i < q->
coefficients.size(); ++i)
308 template <
typename number>
314 if (in_lagrange_product_form ==
true)
315 transform_into_standard_form();
320 std::unique_ptr<Polynomial<number>> q_data;
324 q_data = std::make_unique<Polynomial<number>>(p);
336 for (
unsigned int i = 0; i < q->
coefficients.size(); ++i)
344 template <
typename number>
356 else if (in_lagrange_product_form ==
true)
374 template <
typename number>
375 template <
typename number2>
378 const number2 offset)
386 std::vector<number2> new_coefficients(coefficients.begin(),
392 for (
unsigned int d = 1;
d < new_coefficients.size(); ++
d)
394 const unsigned int n =
d;
399 unsigned int binomial_coefficient = 1;
404 number2 offset_power = offset;
411 for (
unsigned int k = 0; k <
d; ++k)
417 binomial_coefficient = (binomial_coefficient * (n - k)) / (k + 1);
419 new_coefficients[
d - k - 1] +=
420 new_coefficients[
d] * binomial_coefficient * offset_power;
421 offset_power *= offset;
431 coefficients.assign(new_coefficients.begin(), new_coefficients.end());
436 template <
typename number>
437 template <
typename number2>
444 if (in_lagrange_product_form ==
true)
446 for (
unsigned int i = 0; i < lagrange_support_points.size(); ++i)
447 lagrange_support_points[i] -= offset;
451 shift(coefficients, offset);
456 template <
typename number>
465 std::unique_ptr<Polynomial<number>> q_data;
467 if (in_lagrange_product_form ==
true)
469 q_data = std::make_unique<Polynomial<number>>(*this);
476 std::vector<number> newcoefficients(q->
coefficients.size() - 1);
477 for (
unsigned int i = 1; i < q->
coefficients.size(); ++i)
478 newcoefficients[i - 1] = number(i) * q->
coefficients[i];
485 template <
typename number>
491 std::unique_ptr<Polynomial<number>> q_data;
493 if (in_lagrange_product_form ==
true)
495 q_data = std::make_unique<Polynomial<number>>(*this);
502 std::vector<number> newcoefficients(q->
coefficients.size() + 1);
503 newcoefficients[0] = 0.;
504 for (
unsigned int i = 0; i < q->
coefficients.size(); ++i)
505 newcoefficients[i + 1] = q->
coefficients[i] / number(i + 1.);
512 template <
typename number>
516 if (in_lagrange_product_form ==
true)
518 out << lagrange_weight;
519 for (
unsigned int i = 0; i < lagrange_support_points.size(); ++i)
520 out <<
" (x-" << lagrange_support_points[i] <<
")";
524 for (
int i = degree(); i >= 0; --i)
526 out << coefficients[i] <<
" x^" << i << std::endl;
531 template <
typename number>
545 template <
typename number>
549 std::vector<number> result(n + 1, 0.);
550 result[n] = coefficient;
556 template <
typename number>
558 :
Polynomial<number>(make_vector(n, coefficient))
563 template <
typename number>
564 std::vector<Polynomial<number>>
567 std::vector<Polynomial<number>> v;
568 v.reserve(degree + 1);
569 for (
unsigned int i = 0; i <= degree; ++i)
580 namespace LagrangeEquidistantImplementation
582 std::vector<Point<1>>
585 std::vector<Point<1>> points(n + 1);
586 const double one_over_n = 1. / n;
587 for (
unsigned int k = 0; k <= n; ++k)
588 points[k](0) =
static_cast<double>(k) * one_over_n;
597 const unsigned int support_point)
610 std::vector<double> new_support_points;
612 this->coefficients.resize(n + 1);
621 const unsigned int support_point,
622 std::vector<double> &a)
626 unsigned int n_functions = n + 1;
628 double const *x =
nullptr;
634 static const double x1[4] = {1.0, -1.0, 0.0, 1.0};
640 static const double x2[9] = {
641 1.0, -3.0, 2.0, 0.0, 4.0, -4.0, 0.0, -1.0, 2.0};
647 static const double x3[16] = {1.0,
671 for (
unsigned int i = 0; i < n_functions; ++i)
672 a[i] = x[support_point * n_functions + i];
677 std::vector<Polynomial<double>>
682 return std::vector<Polynomial<double>>(
688 std::vector<Polynomial<double>> v;
689 for (
unsigned int i = 0; i <=
degree; ++i)
700 std::vector<Polynomial<double>>
703 std::vector<Polynomial<double>> p;
704 p.reserve(points.size());
706 for (
unsigned int i = 0; i < points.size(); ++i)
707 p.emplace_back(points, i);
720 this->coefficients.clear();
729 for (
unsigned int i = 0; i < k; ++i)
736 for (
unsigned int i = 0; i < k; ++i)
743 std::vector<Polynomial<double>>
746 std::vector<Polynomial<double>> v;
748 for (
unsigned int i = 0; i <=
degree; ++i)
798 std::vector<double> legendre_coefficients_tmp1(p);
799 std::vector<double> legendre_coefficients_tmp2(p - 1);
803 legendre_coefficients_tmp1[0] = 1.0;
805 for (
unsigned int i = 2; i < p; ++i)
807 for (
unsigned int j = 0; j < i - 1; ++j)
808 legendre_coefficients_tmp2[j] = legendre_coefficients_tmp1[j];
810 for (
unsigned int j = 0; j < i; ++j)
815 ((1.0 - 2 * i) * legendre_coefficients_tmp1[0] /
817 (1.0 - i) * legendre_coefficients_tmp2[0] /
821 for (
unsigned int j = 1; j < i - 1; ++j)
825 (2.0 * legendre_coefficients_tmp1[j - 1] -
826 legendre_coefficients_tmp1[j]) +
827 (1.0 - i) * legendre_coefficients_tmp2[j] /
832 (2.0 * legendre_coefficients_tmp1[i - 2] -
833 legendre_coefficients_tmp1[i - 1]) /
836 legendre_coefficients_tmp1[i - 1] / i;
839 for (
int i = p; i > 0; --i)
848 std::vector<Polynomial<double>>
851 std::vector<Polynomial<double>> basis(p + 1);
853 for (
unsigned int i = 0; i <= p; ++i)
866 std::vector<std::unique_ptr<const std::vector<double>>>
887 std::lock_guard<std::mutex> lock(coefficients_lock);
930 std::vector<double> c0(2);
934 std::vector<double> c1(2);
941 std::make_unique<const std::vector<double>>(std::move(c0));
943 std::make_unique<const std::vector<double>>(std::move(c1));
947 coefficients_lock.unlock();
949 coefficients_lock.lock();
951 std::vector<double> c2(3);
960 std::make_unique<const std::vector<double>>(std::move(c2));
972 coefficients_lock.unlock();
974 coefficients_lock.lock();
976 std::vector<double> ck(k + 1);
982 for (
unsigned int i = 1; i <= k - 1; ++i)
1004 std::make_unique<const std::vector<double>>(std::move(ck));
1011 const std::vector<double> &
1021 std::lock_guard<std::mutex> lock(coefficients_lock);
1027 std::vector<Polynomial<double>>
1040 return std::vector<Polynomial<double>>(
1044 std::vector<Polynomial<double>> v;
1046 for (
unsigned int i = 0; i <=
degree; ++i)
1059 this->coefficients.clear();
1107 std::vector<Polynomial<double>>
1112 "degrees less than three"));
1113 std::vector<Polynomial<double>> basis(n + 1);
1115 for (
unsigned int i = 0; i <= n; ++i)
1131 find_support_point_x_star(
const std::vector<double> &jacobi_roots)
1137 double guess_left = 0;
1138 double guess_right = 0.5;
1139 const unsigned int degree = jacobi_roots.size() + 3;
1151 double integral_left = 0, integral_right = 0;
1152 for (
unsigned int q = 0; q < gauss.size(); ++q)
1154 const double x = gauss.point(q)[0];
1155 double poly_val_common = x;
1156 for (
unsigned int j = 0; j < degree - 3; ++j)
1157 poly_val_common *= Utilities::fixed_power<2>(x - jacobi_roots[j]);
1158 poly_val_common *= Utilities::fixed_power<4>(x - 1.);
1160 gauss.weight(q) * (poly_val_common * (x - guess_left));
1162 gauss.weight(q) * (poly_val_common * (x - guess_right));
1167 return guess_right - (guess_right - guess_left) /
1168 (integral_right - integral_left) * integral_right;
1175 const unsigned int index)
1180 this->coefficients.clear();
1187 else if (degree == 1)
1200 else if (degree == 2)
1208 else if (index == 1)
1221 else if (degree == 3)
1238 else if (index == 1)
1248 else if (index == 2)
1255 else if (index == 3)
1286 std::vector<double> jacobi_roots =
1287 jacobi_polynomial_roots<double>(
degree - 3, 4, 4);
1296 const double auxiliary_zero =
1297 find_support_point_x_star(jacobi_roots);
1299 for (
unsigned int m = 0; m <
degree - 3; m++)
1307 else if (index == 1)
1310 for (
unsigned int m = 0; m <
degree - 3; m++)
1323 std::vector<Point<1>> points(
degree);
1325 for (
unsigned int i = 0; i <
degree; ++i)
1332 std::vector<double> value_and_grad(2);
1333 helper.
value(0., value_and_grad);
1337 const double auxiliary_zero =
1338 find_support_point_x_star(jacobi_roots);
1340 (1. / auxiliary_zero - value_and_grad[1] / value_and_grad[0]) /
1343 else if (index >= 2 && index <
degree - 1)
1347 for (
unsigned int m = 0, c = 2; m <
degree - 3; m++)
1357 else if (index ==
degree - 1)
1361 for (
unsigned int m = 0; m <
degree - 3; m++)
1365 std::vector<Point<1>> points(
degree);
1367 for (
unsigned int i = 0; i <
degree; ++i)
1374 std::vector<double> value_and_grad(2);
1375 helper.
value(1., value_and_grad);
1379 const double auxiliary_zero =
1380 find_support_point_x_star(jacobi_roots);
1382 (-1. / auxiliary_zero - value_and_grad[1] / value_and_grad[0]) /
1385 else if (index ==
degree)
1387 const double auxiliary_zero =
1388 find_support_point_x_star(jacobi_roots);
1391 for (
unsigned int m = 0; m <
degree - 3; m++)
1403 std::vector<Polynomial<double>>
1406 std::vector<Polynomial<double>> basis(
degree + 1);
1408 for (
unsigned int i = 0; i <=
degree; ++i)
1421 template class Polynomial<float>;
1422 template class Polynomial<double>;
1423 template class Polynomial<long double>;
1438 template class Monomial<float>;
1439 template class Monomial<double>;
1440 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)
Expression fabs(const Expression &x)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
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)
double legendre(unsigned int l, double x)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)