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> 25 DEAL_II_NAMESPACE_OPEN
50 template <
typename number>
54 in_lagrange_product_form (false),
60 template <
typename number>
63 coefficients (n+1, 0.),
64 in_lagrange_product_form (false),
70 template <
typename number>
72 const unsigned int center)
74 in_lagrange_product_form (true)
76 Assert (supp.size()>0, ExcEmptyObject());
79 lagrange_support_points.reserve (supp.size()-1);
80 number tmp_lagrange_weight = 1.;
81 for (
unsigned int i=0; i<supp.size(); ++i)
84 lagrange_support_points.push_back(supp[i](0));
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."));
94 lagrange_weight =
static_cast<number
>(1.)/tmp_lagrange_weight;
99 template <
typename number>
102 std::vector<number> &values)
const 104 Assert (values.size() > 0, ExcZero());
105 const unsigned int values_size=values.size();
108 if (in_lagrange_product_form ==
true)
113 const unsigned int n_supp = lagrange_support_points.size();
118 for (
unsigned int d=1; d<values_size; ++d)
120 for (
unsigned int i=0; i<n_supp; ++i)
122 const number v = x-lagrange_support_points[i];
130 for (
unsigned int k=values_size-1; k>0; --k)
131 values[k] = (values[k] * v + values[k-1]);
143 number k_faculty = 1;
144 for (
unsigned int k=0; k<values_size; ++k)
146 values[k] *= k_faculty * lagrange_weight;
147 k_faculty *=
static_cast<number
>(k+1);
157 for (
unsigned int i=0; i<n_supp; ++i)
159 const number v = x-lagrange_support_points[i];
162 values[0] *= lagrange_weight;
168 for (
unsigned int i=0; i<n_supp; ++i)
170 const number v = x-lagrange_support_points[i];
171 values[1] = values[1] * v + values[0];
174 values[0] *= lagrange_weight;
175 values[1] *= lagrange_weight;
182 for (
unsigned int i=0; i<n_supp; ++i)
184 const number v = x-lagrange_support_points[i];
185 values[2] = values[2] * v + values[1];
186 values[1] = values[1] * v + values[0];
189 values[0] *= lagrange_weight;
190 values[1] *= lagrange_weight;
191 values[2] *=
static_cast<number
>(2) * lagrange_weight;
197 Assert (coefficients.size() > 0, ExcEmptyObject());
202 if (values_size == 1)
204 values[0] = value(x);
210 const unsigned int m=coefficients.size();
211 std::vector<number> a(coefficients);
212 unsigned int j_faculty=1;
217 const unsigned int min_valuessize_m=std::min(values_size, m);
218 for (
unsigned int j=0; j<min_valuessize_m; ++j)
220 for (
int k=m-2; k>=
static_cast<int>(j); --k)
222 values[j]=
static_cast<number
>(j_faculty)*a[j];
228 for (
unsigned int j=min_valuessize_m; j<values_size; ++j)
234 template <
typename number>
239 Assert (in_lagrange_product_form ==
true, ExcInternalError());
240 Assert (coefficients.size() == 0, ExcInternalError());
243 coefficients.resize (lagrange_support_points.size()+1);
244 if (lagrange_support_points.size() == 0)
245 coefficients[0] = 1.;
248 coefficients[0] = -lagrange_support_points[0];
249 coefficients[1] = 1.;
250 for (
unsigned int i=1; i<lagrange_support_points.size(); ++i)
252 coefficients[i+1] = 1.;
253 for (
unsigned int j=i; j>0; --j)
254 coefficients[j] = (-lagrange_support_points[i]*coefficients[j] +
256 coefficients[0] *= -lagrange_support_points[i];
259 for (
unsigned int i=0; i<lagrange_support_points.size()+1; ++i)
260 coefficients[i] *= lagrange_weight;
263 std::vector<number> new_points;
264 lagrange_support_points.swap(new_points);
265 in_lagrange_product_form =
false;
266 lagrange_weight = 1.;
271 template <
typename number>
277 for (
typename std::vector<number>::iterator c = coefficients.begin();
278 c != coefficients.end(); ++c)
287 template <
typename number>
294 if (in_lagrange_product_form ==
true)
296 number inv_fact = number(1.)/factor;
297 number accumulated_fact = 1.;
298 for (
unsigned int i=0; i<lagrange_support_points.size(); ++i)
300 lagrange_support_points[i] *= inv_fact;
301 accumulated_fact *= factor;
303 lagrange_weight *= accumulated_fact;
307 scale (coefficients, factor);
312 template <
typename number>
317 for (
typename std::vector<number>::iterator c = coefficients.begin();
318 c != coefficients.end(); ++c)
324 template <
typename number>
328 if (in_lagrange_product_form ==
true)
329 lagrange_weight *= s;
332 for (
typename std::vector<number>::iterator c = coefficients.begin();
333 c != coefficients.end(); ++c)
341 template <
typename number>
350 lagrange_support_points.insert (lagrange_support_points.end(),
356 else if (in_lagrange_product_form ==
true)
357 transform_into_standard_form();
362 std_cxx11::shared_ptr<Polynomial<number> > q_data;
374 unsigned int new_degree = this->degree() + q->
degree();
376 std::vector<number> new_coefficients(new_degree+1, 0.);
379 for (
unsigned int j=0; j<this->coefficients.size(); ++j)
380 new_coefficients[i+j] += this->coefficients[j]*q->
coefficients[i];
381 this->coefficients = new_coefficients;
388 template <
typename number>
400 if (in_lagrange_product_form ==
true)
401 transform_into_standard_form();
406 std_cxx11::shared_ptr<Polynomial<number> > q_data;
430 template <
typename number>
436 if (in_lagrange_product_form ==
true)
437 transform_into_standard_form();
442 std_cxx11::shared_ptr<Polynomial<number> > q_data;
466 template <
typename number >
475 if (in_lagrange_product_form ==
true &&
479 else if (in_lagrange_product_form ==
true)
497 template <
typename number>
498 template <
typename number2>
501 const number2 offset)
505 Assert (coefficients.size() < 31, ExcNotImplemented());
509 std::vector<number2> new_coefficients(coefficients.begin(),
515 for (
unsigned int d=1; d<new_coefficients.size(); ++d)
517 const unsigned int n = d;
522 unsigned int binomial_coefficient = 1;
527 number2 offset_power = offset;
534 for (
unsigned int k=0; k<d; ++k)
540 binomial_coefficient = (binomial_coefficient*(n-k))/(k+1);
542 new_coefficients[d-k-1] += new_coefficients[d]
543 * binomial_coefficient
545 offset_power *= offset;
551 Assert (binomial_coefficient == 1, ExcInternalError());
555 coefficients.assign(new_coefficients.begin(), new_coefficients.end());
560 template <
typename number>
561 template <
typename number2>
568 if (in_lagrange_product_form ==
true)
570 for (
unsigned int i=0; i<lagrange_support_points.size(); ++i)
571 lagrange_support_points[i] -= offset;
575 shift (coefficients, offset);
580 template <
typename number>
589 std_cxx11::shared_ptr<Polynomial<number> > q_data;
591 if (in_lagrange_product_form ==
true)
600 std::vector<number> newcoefficients (q->
coefficients.size()-1);
609 template <
typename number>
615 std_cxx11::shared_ptr<Polynomial<number> > q_data;
617 if (in_lagrange_product_form ==
true)
626 std::vector<number> newcoefficients (q->
coefficients.size()+1);
627 newcoefficients[0] = 0.;
636 template <
typename number>
640 if (in_lagrange_product_form ==
true)
642 out << lagrange_weight;
643 for (
unsigned int i=0; i<lagrange_support_points.size(); ++i)
644 out <<
" (x-" << lagrange_support_points[i] <<
")";
648 for (
int i=degree(); i>=0; --i)
650 out << coefficients[i] <<
" x^" << i << std::endl;
658 template <
typename number>
663 std::vector<number> result(n+1, 0.);
664 result[n] = coefficient;
670 template <
typename number>
673 :
Polynomial<number>(make_vector(n, coefficient))
678 template <
typename number>
679 std::vector<Polynomial<number> >
682 std::vector<Polynomial<number> > v;
684 for (
unsigned int i=0; i<=
degree; ++i)
697 std::vector<Point<1> >
698 generate_equidistant_unit_points (
const unsigned int n)
700 std::vector<Point<1> > points (n+1);
701 const double one_over_n = 1./n;
702 for (
unsigned int k=0; k<=n; ++k)
703 points[k](0) =
static_cast<double>(k)*one_over_n;
712 const unsigned int support_point)
715 generate_equidistant_unit_points (n),
726 std::vector<double> new_support_points;
737 const unsigned int support_point,
738 std::vector<double> &a)
740 Assert(support_point<n+1, ExcIndexRange(support_point, 0, n+1));
742 unsigned int n_functions=n+1;
743 Assert(support_point<n_functions,
744 ExcIndexRange(support_point, 0, n_functions));
751 static const double x1[4]=
761 static const double x2[9]=
772 static const double x3[16]=
774 1.0, -11.0/2.0, 9.0, -9.0/2.0,
775 0.0, 9.0, -45.0/2.0, 27.0/2.0,
776 0.0, -9.0/2.0, 18.0, -27.0/2.0,
777 0.0, 1.0, -9.0/2.0, 9.0/2.0
783 Assert(
false, ExcInternalError())
786 Assert(x!=0, ExcInternalError());
787 for (
unsigned int i=0; i<n_functions; ++i)
788 a[i]=x[support_point*n_functions+i];
793 std::vector<Polynomial<double> >
799 return std::vector<Polynomial<double> >
805 std::vector<Polynomial<double> > v;
806 for (
unsigned int i=0; i<=
degree; ++i)
817 std::vector<Polynomial<double> >
820 std::vector<Polynomial<double> > p;
821 p.reserve (points.size());
823 for (
unsigned int i=0; i<points.size(); ++i)
835 std::vector<std_cxx11::shared_ptr<const std::vector<double> > >
837 std::vector<std_cxx11::shared_ptr<const std::vector<double> > >
856 #ifdef DEAL_II_LONG_DOUBLE_LOOP_BUG 857 typedef double SHIFT_TYPE;
859 typedef long double SHIFT_TYPE;
879 std_cxx11::shared_ptr<
const std::vector<double> >())))
913 std::vector<double> *c0 =
new std::vector<double>(1);
916 std::vector<double> *c1 =
new std::vector<double>(2);
927 std_cxx11::shared_ptr<const std::vector<double> >(c0);
929 std_cxx11::shared_ptr<const std::vector<double> >(c1);
933 c0 =
new std::vector<double>(*c0);
934 c1 =
new std::vector<double>(*c1);
954 coefficients_lock.release ();
956 coefficients_lock.acquire ();
958 std::vector<double> *ck =
new std::vector<double>(k+1);
960 const double a = 1./(k);
961 const double b = a*(2*k-1);
962 const double c = a*(k-1);
966 for (
unsigned int i=1 ; i<= k-2 ; ++i)
977 std_cxx11::shared_ptr<const std::vector<double> >(ck);
980 ck =
new std::vector<double>(*ck);
985 std_cxx11::shared_ptr<const std::vector<double> >(ck);
992 const std::vector<double> &
1008 std::vector<Polynomial<double> >
1011 std::vector<Polynomial<double> > v;
1012 v.reserve(degree+1);
1013 for (
unsigned int i=0; i<=
degree; ++i)
1035 coefficients[0] = 1.0;
1036 coefficients[1] = -1.0;
1044 coefficients[0] = 0.0;
1045 coefficients[1] = 1.0;
1053 coefficients[0] = 0.0;
1054 coefficients[1] = -1.0 * std::sqrt (3.);
1055 coefficients[2] = std::sqrt (3.);
1062 std::vector<double> legendre_coefficients_tmp1 (p);
1063 std::vector<double> legendre_coefficients_tmp2 (p - 1);
1065 coefficients[0] = -1.0 * std::sqrt (3.);
1066 coefficients[1] = 2.0 * std::sqrt (3.);
1067 legendre_coefficients_tmp1[0] = 1.0;
1069 for (
unsigned int i = 2; i < p; ++i)
1071 for (
unsigned int j = 0; j < i - 1; ++j)
1072 legendre_coefficients_tmp2[j] = legendre_coefficients_tmp1[j];
1074 for (
unsigned int j = 0; j < i; ++j)
1075 legendre_coefficients_tmp1[j] = coefficients[j];
1077 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;
1079 for (
unsigned int j = 1; j < i - 1; ++j)
1080 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;
1082 coefficients[i - 1] = std::sqrt (4 * i * i - 1.) * (2.0 * legendre_coefficients_tmp1[i - 2] - legendre_coefficients_tmp1[i - 1]) / i;
1083 coefficients[i] = 2.0 * std::sqrt (4 * i * i - 1.) * legendre_coefficients_tmp1[i - 1] / i;
1086 for (
int i = p; i > 0; --i)
1087 coefficients[i] = coefficients[i - 1] / i;
1089 coefficients[0] = 0.0;
1097 std::vector<Polynomial<double> > basis (p + 1);
1099 for (
unsigned int i = 0; i <= p; ++i)
1112 std::vector<std_cxx11::shared_ptr<const std::vector<double> > >
1127 unsigned int k = k_;
1178 std::vector<double> *c0 =
new std::vector<double>(2);
1182 std::vector<double> *c1 =
new std::vector<double>(2);
1189 std_cxx11::shared_ptr<const std::vector<double> >(c0);
1191 std_cxx11::shared_ptr<const std::vector<double> >(c1);
1195 coefficients_lock.release ();
1197 coefficients_lock.acquire ();
1199 std::vector<double> *c2 =
new std::vector<double>(3);
1201 const double a = 1.;
1208 std_cxx11::shared_ptr<const std::vector<double> >(c2);
1220 coefficients_lock.release ();
1222 coefficients_lock.acquire ();
1224 std::vector<double> *ck =
new std::vector<double>(k+1);
1226 const double a = 1.;
1230 for (
unsigned int i=1; i<=k-1; ++i)
1252 std_cxx11::shared_ptr<const std::vector<double> >(ck);
1259 const std::vector<double> &
1275 std::vector<Polynomial<double> >
1288 return std::vector<Polynomial<double> >
1292 std::vector<Polynomial<double> > v;
1293 v.reserve(degree+1);
1294 for (
unsigned int i=0; i<=
degree; ++i)
1337 (*this) *= legendre;
1343 std::vector<Polynomial<double> >
1346 std::vector<Polynomial<double> > basis (n + 1);
1348 for (
unsigned int i = 0; i <= n; ++i)
1375 DEAL_II_NAMESPACE_CLOSE
void scale(const number factor)
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)
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)
::ExceptionBase & ExcMessage(std::string arg1)
#define AssertIndexRange(index, range)
std::vector< Polynomial< double > > generate_complete_Lagrange_basis(const std::vector< Point< 1 > > &points)
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
static std::vector< std_cxx11::shared_ptr< const std::vector< double > > > shifted_coefficients
std::vector< double > compute_coefficients(const unsigned int p)
static void multiply(std::vector< number > &coefficients, const number factor)
#define Assert(cond, exc)
HermiteInterpolation(const unsigned int p)
Lobatto(const unsigned int p=0)
static std::vector< std_cxx11::shared_ptr< const std::vector< double > > > recursive_coefficients
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 k)
static const std::vector< double > & get_coefficients(const unsigned int p)
std::vector< number > lagrange_support_points
static void compute_coefficients(const unsigned int p)
static std::vector< std_cxx11::shared_ptr< const std::vector< double > > > recursive_coefficients
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int p)