16 #include <deal.II/base/polynomial.h> 17 #include <deal.II/base/point.h> 18 #include <deal.II/base/exceptions.h> 19 #include <deal.II/base/thread_management.h> 20 #include <deal.II/base/quadrature_lib.h> 26 DEAL_II_NAMESPACE_OPEN
51 template <
typename number>
55 in_lagrange_product_form (false),
61 template <
typename number>
64 coefficients (n+1, 0.),
65 in_lagrange_product_form (false),
71 template <
typename number>
73 const unsigned int center)
75 in_lagrange_product_form (true)
80 lagrange_support_points.reserve (supp.size()-1);
81 number tmp_lagrange_weight = 1.;
82 for (
unsigned int i=0; i<supp.size(); ++i)
85 lagrange_support_points.push_back(supp[i](0));
86 tmp_lagrange_weight *= supp[center](0) - supp[i](0);
90 Assert (std::fabs(tmp_lagrange_weight) > std::numeric_limits<number>::min(),
91 ExcMessage (
"Underflow in computation of Lagrange denominator."));
92 Assert (std::fabs(tmp_lagrange_weight) < std::numeric_limits<number>::max(),
93 ExcMessage (
"Overflow in computation of Lagrange denominator."));
95 lagrange_weight =
static_cast<number
>(1.)/tmp_lagrange_weight;
100 template <
typename number>
103 std::vector<number> &values)
const 107 value(x,values.size()-1,&values[0]);
112 template <
typename number>
115 const unsigned int n_derivatives,
116 number *values)
const 119 if (in_lagrange_product_form ==
true)
124 const unsigned int n_supp = lagrange_support_points.size();
125 switch (n_derivatives)
129 for (
unsigned int d=1; d<=n_derivatives; ++d)
131 for (
unsigned int i=0; i<n_supp; ++i)
133 const number v = x-lagrange_support_points[i];
141 for (
unsigned int k=n_derivatives; k>0; --k)
142 values[k] = (values[k] * v + values[k-1]);
154 number k_faculty = 1;
155 for (
unsigned int k=0; k<=n_derivatives; ++k)
157 values[k] *= k_faculty * lagrange_weight;
158 k_faculty *=
static_cast<number
>(k+1);
168 for (
unsigned int i=0; i<n_supp; ++i)
170 const number v = x-lagrange_support_points[i];
173 values[0] *= lagrange_weight;
179 for (
unsigned int i=0; i<n_supp; ++i)
181 const number v = x-lagrange_support_points[i];
182 values[1] = values[1] * v + values[0];
185 values[0] *= lagrange_weight;
186 values[1] *= lagrange_weight;
193 for (
unsigned int i=0; i<n_supp; ++i)
195 const number v = x-lagrange_support_points[i];
196 values[2] = values[2] * v + values[1];
197 values[1] = values[1] * v + values[0];
200 values[0] *= lagrange_weight;
201 values[1] *= lagrange_weight;
202 values[2] *=
static_cast<number
>(2) * lagrange_weight;
213 if (n_derivatives == 0)
215 values[0] = value(x);
221 const unsigned int m=coefficients.size();
222 std::vector<number> a(coefficients);
223 unsigned int j_faculty=1;
228 const unsigned int min_valuessize_m=std::min(n_derivatives+1, m);
229 for (
unsigned int j=0; j<min_valuessize_m; ++j)
231 for (
int k=m-2; k>=
static_cast<int>(j); --k)
233 values[j]=
static_cast<number
>(j_faculty)*a[j];
239 for (
unsigned int j=min_valuessize_m; j<=n_derivatives; ++j)
245 template <
typename number>
254 coefficients.resize (lagrange_support_points.size()+1);
255 if (lagrange_support_points.size() == 0)
256 coefficients[0] = 1.;
259 coefficients[0] = -lagrange_support_points[0];
260 coefficients[1] = 1.;
261 for (
unsigned int i=1; i<lagrange_support_points.size(); ++i)
263 coefficients[i+1] = 1.;
264 for (
unsigned int j=i; j>0; --j)
265 coefficients[j] = (-lagrange_support_points[i]*coefficients[j] +
267 coefficients[0] *= -lagrange_support_points[i];
270 for (
unsigned int i=0; i<lagrange_support_points.size()+1; ++i)
271 coefficients[i] *= lagrange_weight;
274 std::vector<number> new_points;
275 lagrange_support_points.swap(new_points);
276 in_lagrange_product_form =
false;
277 lagrange_weight = 1.;
282 template <
typename number>
288 for (
typename std::vector<number>::iterator c = coefficients.begin();
289 c != coefficients.end(); ++c)
298 template <
typename number>
305 if (in_lagrange_product_form ==
true)
307 number inv_fact = number(1.)/factor;
308 number accumulated_fact = 1.;
309 for (
unsigned int i=0; i<lagrange_support_points.size(); ++i)
311 lagrange_support_points[i] *= inv_fact;
312 accumulated_fact *= factor;
314 lagrange_weight *= accumulated_fact;
318 scale (coefficients, factor);
323 template <
typename number>
328 for (
typename std::vector<number>::iterator c = coefficients.begin();
329 c != coefficients.end(); ++c)
335 template <
typename number>
339 if (in_lagrange_product_form ==
true)
340 lagrange_weight *= s;
343 for (
typename std::vector<number>::iterator c = coefficients.begin();
344 c != coefficients.end(); ++c)
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_cxx11::shared_ptr<Polynomial<number> > q_data;
385 unsigned int new_degree = this->degree() + q->
degree();
387 std::vector<number> new_coefficients(new_degree+1, 0.);
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_cxx11::shared_ptr<Polynomial<number> > q_data;
441 template <
typename number>
447 if (in_lagrange_product_form ==
true)
448 transform_into_standard_form();
453 std_cxx11::shared_ptr<Polynomial<number> > q_data;
477 template <
typename number >
486 if (in_lagrange_product_form ==
true &&
490 else if (in_lagrange_product_form ==
true)
508 template <
typename number>
509 template <
typename number2>
512 const number2 offset)
520 std::vector<number2> new_coefficients(coefficients.begin(),
526 for (
unsigned int d=1; d<new_coefficients.size(); ++d)
528 const unsigned int n = d;
533 unsigned int binomial_coefficient = 1;
538 number2 offset_power = offset;
545 for (
unsigned int k=0; k<d; ++k)
551 binomial_coefficient = (binomial_coefficient*(n-k))/(k+1);
553 new_coefficients[d-k-1] += new_coefficients[d]
554 * binomial_coefficient
556 offset_power *= offset;
566 coefficients.assign(new_coefficients.begin(), new_coefficients.end());
571 template <
typename number>
572 template <
typename number2>
579 if (in_lagrange_product_form ==
true)
581 for (
unsigned int i=0; i<lagrange_support_points.size(); ++i)
582 lagrange_support_points[i] -= offset;
586 shift (coefficients, offset);
591 template <
typename number>
600 std_cxx11::shared_ptr<Polynomial<number> > q_data;
602 if (in_lagrange_product_form ==
true)
611 std::vector<number> newcoefficients (q->
coefficients.size()-1);
620 template <
typename number>
626 std_cxx11::shared_ptr<Polynomial<number> > q_data;
628 if (in_lagrange_product_form ==
true)
637 std::vector<number> newcoefficients (q->
coefficients.size()+1);
638 newcoefficients[0] = 0.;
647 template <
typename number>
651 if (in_lagrange_product_form ==
true)
653 out << lagrange_weight;
654 for (
unsigned int i=0; i<lagrange_support_points.size(); ++i)
655 out <<
" (x-" << lagrange_support_points[i] <<
")";
659 for (
int i=degree(); i>=0; --i)
661 out << coefficients[i] <<
" x^" << i << std::endl;
669 template <
typename number>
674 std::vector<number> result(n+1, 0.);
675 result[n] = coefficient;
681 template <
typename number>
684 :
Polynomial<number>(make_vector(n, coefficient))
689 template <
typename number>
690 std::vector<Polynomial<number> >
693 std::vector<Polynomial<number> > v;
695 for (
unsigned int i=0; i<=degree; ++i)
708 std::vector<Point<1> >
709 generate_equidistant_unit_points (
const unsigned int n)
711 std::vector<Point<1> > points (n+1);
712 const double one_over_n = 1./n;
713 for (
unsigned int k=0; k<=n; ++k)
714 points[k](0) =
static_cast<double>(k)*one_over_n;
723 const unsigned int support_point)
726 generate_equidistant_unit_points (n),
737 std::vector<double> new_support_points;
748 const unsigned int support_point,
749 std::vector<double> &a)
753 unsigned int n_functions=n+1;
754 Assert(support_point<n_functions,
762 static const double x1[4]=
772 static const double x2[9]=
783 static const double x3[16]=
785 1.0, -11.0/2.0, 9.0, -9.0/2.0,
786 0.0, 9.0, -45.0/2.0, 27.0/2.0,
787 0.0, -9.0/2.0, 18.0, -27.0/2.0,
788 0.0, 1.0, -9.0/2.0, 9.0/2.0
798 for (
unsigned int i=0; i<n_functions; ++i)
799 a[i]=x[support_point*n_functions+i];
804 std::vector<Polynomial<double> >
810 return std::vector<Polynomial<double> >
816 std::vector<Polynomial<double> > v;
817 for (
unsigned int i=0; i<=
degree; ++i)
828 std::vector<Polynomial<double> >
831 std::vector<Polynomial<double> > p;
832 p.reserve (points.size());
834 for (
unsigned int i=0; i<points.size(); ++i)
858 for (
unsigned int i=0; i<k; ++i)
865 for (
unsigned int i=0; i<k; ++i)
872 std::vector<Polynomial<double> >
875 std::vector<Polynomial<double> > v;
877 for (
unsigned int i=0; i<=
degree; ++i)
926 std::vector<double> legendre_coefficients_tmp1 (p);
927 std::vector<double> legendre_coefficients_tmp2 (p - 1);
931 legendre_coefficients_tmp1[0] = 1.0;
933 for (
unsigned int i = 2; i < p; ++i)
935 for (
unsigned int j = 0; j < i - 1; ++j)
936 legendre_coefficients_tmp2[j] = legendre_coefficients_tmp1[j];
938 for (
unsigned int j = 0; j < i; ++j)
941 coefficients[0] = std::sqrt (2 * i + 1.) * ((1.0 - 2 * i) * legendre_coefficients_tmp1[0] / std::sqrt (2 * i - 1.) + (1.0 - i) * legendre_coefficients_tmp2[0] / std::sqrt (2 * i - 3.)) / i;
943 for (
unsigned int j = 1; j < i - 1; ++j)
944 coefficients[j] = std::sqrt (2 * i + 1.) * (std::sqrt (2 * i - 1.) * (2.0 * legendre_coefficients_tmp1[j - 1] - legendre_coefficients_tmp1[j]) + (1.0 - i) * legendre_coefficients_tmp2[j] / std::sqrt (2 * i - 3.)) / i;
946 coefficients[i - 1] = std::sqrt (4 * i * i - 1.) * (2.0 * legendre_coefficients_tmp1[i - 2] - legendre_coefficients_tmp1[i - 1]) / i;
947 coefficients[i] = 2.0 * std::sqrt (4 * i * i - 1.) * legendre_coefficients_tmp1[i - 1] / i;
950 for (
int i = p; i > 0; --i)
961 std::vector<Polynomial<double> > basis (p + 1);
963 for (
unsigned int i = 0; i <= p; ++i)
976 std::vector<std_cxx11::shared_ptr<const std::vector<double> > >
1042 std::vector<double> *c0 =
new std::vector<double>(2);
1046 std::vector<double> *c1 =
new std::vector<double>(2);
1053 std_cxx11::shared_ptr<const std::vector<double> >(c0);
1055 std_cxx11::shared_ptr<const std::vector<double> >(c1);
1059 coefficients_lock.release ();
1061 coefficients_lock.acquire ();
1063 std::vector<double> *c2 =
new std::vector<double>(3);
1065 const double a = 1.;
1072 std_cxx11::shared_ptr<const std::vector<double> >(c2);
1084 coefficients_lock.release ();
1086 coefficients_lock.acquire ();
1088 std::vector<double> *ck =
new std::vector<double>(k+1);
1090 const double a = 1.;
1094 for (
unsigned int i=1; i<=k-1; ++i)
1116 std_cxx11::shared_ptr<const std::vector<double> >(ck);
1123 const std::vector<double> &
1139 std::vector<Polynomial<double> >
1152 return std::vector<Polynomial<double> >
1156 std::vector<Polynomial<double> > v;
1158 for (
unsigned int i=0; i<=
degree; ++i)
1214 (*this) *= legendre;
1220 std::vector<Polynomial<double> >
1223 std::vector<Polynomial<double> > basis (n + 1);
1225 for (
unsigned int i = 0; i <= n; ++i)
1237 template class Polynomial<float>;
1238 template class Polynomial<double>;
1239 template class Polynomial<long double>;
1248 template class Monomial<float>;
1249 template class Monomial<double>;
1250 template class Monomial<long double>;
1253 DEAL_II_NAMESPACE_CLOSE
void transform_into_standard_form()
LagrangeEquidistant(const unsigned int n, const unsigned int support_point)
void shift(const number2 offset)
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int p)
const std::vector< Point< dim > > & get_points() const
unsigned int degree() const
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
Legendre(const unsigned int p)
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
#define AssertIndexRange(index, range)
std::vector< Polynomial< double > > generate_complete_Lagrange_basis(const std::vector< Point< 1 > > &points)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
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)
bool in_lagrange_product_form
std::vector< double > compute_coefficients(const unsigned int p)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
HermiteInterpolation(const unsigned int p)
Lobatto(const unsigned int p=0)
static std::vector< Polynomial< number > > generate_complete_basis(const unsigned int degree)
Hierarchical(const unsigned int p)
std::vector< number > coefficients
static void compute_coefficients(const unsigned int p)
static const std::vector< double > & get_coefficients(const unsigned int p)
std::vector< number > lagrange_support_points
static ::ExceptionBase & ExcEmptyObject()
static ::ExceptionBase & ExcNotImplemented()
static std::vector< std_cxx11::shared_ptr< const std::vector< double > > > recursive_coefficients
static ::ExceptionBase & ExcZero()
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int p)
static ::ExceptionBase & ExcInternalError()