Reference documentation for deal.II version 8.5.1
polynomial.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
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>
21 
22 #include <cmath>
23 #include <algorithm>
24 #include <limits>
25 
26 DEAL_II_NAMESPACE_OPEN
27 
28 
29 
30 // have a lock that guarantees that at most one thread is changing and
31 // accessing the @p{coefficients} arrays of classes implementing
32 // polynomials with tables. make this lock local to this file.
33 //
34 // having only one lock for all of these classes is probably not going
35 // to be a problem since we only need it on very rare occasions. if
36 // someone finds this is a bottleneck, feel free to replace it by a
37 // more fine-grained solution
38 namespace
39 {
40  Threads::Mutex coefficients_lock;
41 }
42 
43 
44 
45 namespace Polynomials
46 {
47 
48 // -------------------- class Polynomial ---------------- //
49 
50 
51  template <typename number>
52  Polynomial<number>::Polynomial (const std::vector<number> &a)
53  :
54  coefficients (a),
55  in_lagrange_product_form (false),
56  lagrange_weight (1.)
57  {}
58 
59 
60 
61  template <typename number>
62  Polynomial<number>::Polynomial (const unsigned int n)
63  :
64  coefficients (n+1, 0.),
65  in_lagrange_product_form (false),
66  lagrange_weight (1.)
67  {}
68 
69 
70 
71  template <typename number>
72  Polynomial<number>::Polynomial (const std::vector<Point<1> > &supp,
73  const unsigned int center)
74  :
75  in_lagrange_product_form (true)
76  {
77  Assert (supp.size()>0, ExcEmptyObject());
78  AssertIndexRange (center, supp.size());
79 
80  lagrange_support_points.reserve (supp.size()-1);
81  number tmp_lagrange_weight = 1.;
82  for (unsigned int i=0; i<supp.size(); ++i)
83  if (i!=center)
84  {
85  lagrange_support_points.push_back(supp[i](0));
86  tmp_lagrange_weight *= supp[center](0) - supp[i](0);
87  }
88 
89  // check for underflow and overflow
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."));
94 
95  lagrange_weight = static_cast<number>(1.)/tmp_lagrange_weight;
96  }
97 
98 
99 
100  template <typename number>
101  void
102  Polynomial<number>::value (const number x,
103  std::vector<number> &values) const
104  {
105  Assert (values.size() > 0, ExcZero());
106 
107  value(x,values.size()-1,&values[0]);
108  }
109 
110 
111 
112  template <typename number>
113  void
114  Polynomial<number>::value (const number x,
115  const unsigned int n_derivatives,
116  number *values) const
117  {
118  // evaluate Lagrange polynomial and derivatives
119  if (in_lagrange_product_form == true)
120  {
121  // to compute the value and all derivatives of a polynomial of the
122  // form (x-x_1)*(x-x_2)*...*(x-x_n), expand the derivatives like
123  // automatic differentiation does.
124  const unsigned int n_supp = lagrange_support_points.size();
125  switch (n_derivatives)
126  {
127  default:
128  values[0] = 1;
129  for (unsigned int d=1; d<=n_derivatives; ++d)
130  values[d] = 0;
131  for (unsigned int i=0; i<n_supp; ++i)
132  {
133  const number v = x-lagrange_support_points[i];
134 
135  // multiply by (x-x_i) and compute action on all derivatives,
136  // too (inspired from automatic differentiation: implement the
137  // product rule for the old value and the new variable 'v',
138  // i.e., expand value v and derivative one). since we reuse a
139  // value from the next lower derivative from the steps before,
140  // need to start from the highest derivative
141  for (unsigned int k=n_derivatives; k>0; --k)
142  values[k] = (values[k] * v + values[k-1]);
143  values[0] *= v;
144  }
145  // finally, multiply by the weight in the Lagrange
146  // denominator. Could be done instead of setting values[0] = 1
147  // above, but that gives different accumulation of round-off
148  // errors (multiplication is not associative) compared to when we
149  // computed the weight, and hence a basis function might not be
150  // exactly one at the center point, which is nice to have. We also
151  // multiply derivatives by k! to transform the product p_n =
152  // p^(n)(x)/k! into the actual form of the derivative
153  {
154  number k_faculty = 1;
155  for (unsigned int k=0; k<=n_derivatives; ++k)
156  {
157  values[k] *= k_faculty * lagrange_weight;
158  k_faculty *= static_cast<number>(k+1);
159  }
160  }
161  break;
162 
163  // manually implement case 0 (values only), case 1 (value + first
164  // derivative), and case 2 (up to second derivative) since they
165  // might be called often. then, we can unroll the loop.
166  case 0:
167  values[0] = 1;
168  for (unsigned int i=0; i<n_supp; ++i)
169  {
170  const number v = x-lagrange_support_points[i];
171  values[0] *= v;
172  }
173  values[0] *= lagrange_weight;
174  break;
175 
176  case 1:
177  values[0] = 1;
178  values[1] = 0;
179  for (unsigned int i=0; i<n_supp; ++i)
180  {
181  const number v = x-lagrange_support_points[i];
182  values[1] = values[1] * v + values[0];
183  values[0] *= v;
184  }
185  values[0] *= lagrange_weight;
186  values[1] *= lagrange_weight;
187  break;
188 
189  case 2:
190  values[0] = 1;
191  values[1] = 0;
192  values[2] = 0;
193  for (unsigned int i=0; i<n_supp; ++i)
194  {
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];
198  values[0] *= v;
199  }
200  values[0] *= lagrange_weight;
201  values[1] *= lagrange_weight;
202  values[2] *= static_cast<number>(2) * lagrange_weight;
203  break;
204  }
205  return;
206  }
207 
208  Assert (coefficients.size() > 0, ExcEmptyObject());
209 
210  // if we only need the value, then call the other function since that is
211  // significantly faster (there is no need to allocate and free memory,
212  // which is really expensive compared to all the other operations!)
213  if (n_derivatives == 0)
214  {
215  values[0] = value(x);
216  return;
217  };
218 
219  // if there are derivatives needed, then do it properly by the full Horner
220  // scheme
221  const unsigned int m=coefficients.size();
222  std::vector<number> a(coefficients);
223  unsigned int j_faculty=1;
224 
225  // loop over all requested derivatives. note that derivatives @p{j>m} are
226  // necessarily zero, as they differentiate the polynomial more often than
227  // the highest power is
228  const unsigned int min_valuessize_m=std::min(n_derivatives+1, m);
229  for (unsigned int j=0; j<min_valuessize_m; ++j)
230  {
231  for (int k=m-2; k>=static_cast<int>(j); --k)
232  a[k]+=x*a[k+1];
233  values[j]=static_cast<number>(j_faculty)*a[j];
234 
235  j_faculty*=j+1;
236  }
237 
238  // fill higher derivatives by zero
239  for (unsigned int j=min_valuessize_m; j<=n_derivatives; ++j)
240  values[j] = 0;
241  }
242 
243 
244 
245  template <typename number>
246  void
248  {
249  // should only be called when the product form is active
250  Assert (in_lagrange_product_form == true, ExcInternalError());
251  Assert (coefficients.size() == 0, ExcInternalError());
252 
253  // compute coefficients by expanding the product (x-x_i) term by term
254  coefficients.resize (lagrange_support_points.size()+1);
255  if (lagrange_support_points.size() == 0)
256  coefficients[0] = 1.;
257  else
258  {
259  coefficients[0] = -lagrange_support_points[0];
260  coefficients[1] = 1.;
261  for (unsigned int i=1; i<lagrange_support_points.size(); ++i)
262  {
263  coefficients[i+1] = 1.;
264  for (unsigned int j=i; j>0; --j)
265  coefficients[j] = (-lagrange_support_points[i]*coefficients[j] +
266  coefficients[j-1]);
267  coefficients[0] *= -lagrange_support_points[i];
268  }
269  }
270  for (unsigned int i=0; i<lagrange_support_points.size()+1; ++i)
271  coefficients[i] *= lagrange_weight;
272 
273  // delete the product form data
274  std::vector<number> new_points;
275  lagrange_support_points.swap(new_points);
276  in_lagrange_product_form = false;
277  lagrange_weight = 1.;
278  }
279 
280 
281 
282  template <typename number>
283  void
284  Polynomial<number>::scale (std::vector<number> &coefficients,
285  const number factor)
286  {
287  number f = 1.;
288  for (typename std::vector<number>::iterator c = coefficients.begin();
289  c != coefficients.end(); ++c)
290  {
291  *c *= f;
292  f *= factor;
293  }
294  }
295 
296 
297 
298  template <typename number>
299  void
300  Polynomial<number>::scale (const number factor)
301  {
302  // to scale (x-x_0)*(x-x_1)*...*(x-x_n), scale
303  // support points by 1./factor and the weight
304  // likewise
305  if (in_lagrange_product_form == true)
306  {
307  number inv_fact = number(1.)/factor;
308  number accumulated_fact = 1.;
309  for (unsigned int i=0; i<lagrange_support_points.size(); ++i)
310  {
311  lagrange_support_points[i] *= inv_fact;
312  accumulated_fact *= factor;
313  }
314  lagrange_weight *= accumulated_fact;
315  }
316  // otherwise, use the function above
317  else
318  scale (coefficients, factor);
319  }
320 
321 
322 
323  template <typename number>
324  void
325  Polynomial<number>::multiply (std::vector<number> &coefficients,
326  const number factor)
327  {
328  for (typename std::vector<number>::iterator c = coefficients.begin();
329  c != coefficients.end(); ++c)
330  *c *= factor;
331  }
332 
333 
334 
335  template <typename number>
338  {
339  if (in_lagrange_product_form == true)
340  lagrange_weight *= s;
341  else
342  {
343  for (typename std::vector<number>::iterator c = coefficients.begin();
344  c != coefficients.end(); ++c)
345  *c *= s;
346  }
347  return *this;
348  }
349 
350 
351 
352  template <typename number>
355  {
356  // if we are in Lagrange form, just append the
357  // new points
358  if (in_lagrange_product_form == true && p.in_lagrange_product_form == true)
359  {
360  lagrange_weight *= p.lagrange_weight;
361  lagrange_support_points.insert (lagrange_support_points.end(),
362  p.lagrange_support_points.begin(),
363  p.lagrange_support_points.end());
364  }
365 
366  // cannot retain product form, recompute...
367  else if (in_lagrange_product_form == true)
368  transform_into_standard_form();
369 
370  // need to transform p into standard form as
371  // well if necessary. copy the polynomial to
372  // do this
373  std_cxx11::shared_ptr<Polynomial<number> > q_data;
374  const Polynomial<number> *q = 0;
375  if (p.in_lagrange_product_form == true)
376  {
377  q_data.reset (new Polynomial<number>(p));
379  q = q_data.get();
380  }
381  else
382  q = &p;
383 
384  // Degree of the product
385  unsigned int new_degree = this->degree() + q->degree();
386 
387  std::vector<number> new_coefficients(new_degree+1, 0.);
388 
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;
393 
394  return *this;
395  }
396 
397 
398 
399  template <typename number>
402  {
403  // Lagrange product form cannot reasonably be
404  // retained after polynomial addition. we
405  // could in theory check if either this
406  // polynomial or the other is a zero
407  // polynomial and retain it, but we actually
408  // currently (r23974) assume that the addition
409  // of a zero polynomial changes the state and
410  // tests equivalence.
411  if (in_lagrange_product_form == true)
412  transform_into_standard_form();
413 
414  // need to transform p into standard form as
415  // well if necessary. copy the polynomial to
416  // do this
417  std_cxx11::shared_ptr<Polynomial<number> > q_data;
418  const Polynomial<number> *q = 0;
419  if (p.in_lagrange_product_form == true)
420  {
421  q_data.reset (new Polynomial<number>(p));
423  q = q_data.get();
424  }
425  else
426  q = &p;
427 
428  // if necessary expand the number
429  // of coefficients we store
430  if (q->coefficients.size() > coefficients.size())
431  coefficients.resize (q->coefficients.size(), 0.);
432 
433  for (unsigned int i=0; i<q->coefficients.size(); ++i)
434  coefficients[i] += q->coefficients[i];
435 
436  return *this;
437  }
438 
439 
440 
441  template <typename number>
444  {
445  // Lagrange product form cannot reasonably be
446  // retained after polynomial addition
447  if (in_lagrange_product_form == true)
448  transform_into_standard_form();
449 
450  // need to transform p into standard form as
451  // well if necessary. copy the polynomial to
452  // do this
453  std_cxx11::shared_ptr<Polynomial<number> > q_data;
454  const Polynomial<number> *q = 0;
455  if (p.in_lagrange_product_form == true)
456  {
457  q_data.reset (new Polynomial<number>(p));
459  q = q_data.get();
460  }
461  else
462  q = &p;
463 
464  // if necessary expand the number
465  // of coefficients we store
466  if (q->coefficients.size() > coefficients.size())
467  coefficients.resize (q->coefficients.size(), 0.);
468 
469  for (unsigned int i=0; i<q->coefficients.size(); ++i)
470  coefficients[i] -= q->coefficients[i];
471 
472  return *this;
473  }
474 
475 
476 
477  template <typename number >
478  bool
480  {
481  // need to distinguish a few cases based on
482  // whether we are in product form or not. two
483  // polynomials can still be the same when they
484  // are on different forms, but the expansion
485  // is the same
486  if (in_lagrange_product_form == true &&
487  p.in_lagrange_product_form == true)
488  return ((lagrange_weight == p.lagrange_weight) &&
489  (lagrange_support_points == p.lagrange_support_points));
490  else if (in_lagrange_product_form == true)
491  {
492  Polynomial<number> q = *this;
494  return (q.coefficients == p.coefficients);
495  }
496  else if (p.in_lagrange_product_form == true)
497  {
498  Polynomial<number> q = p;
500  return (q.coefficients == coefficients);
501  }
502  else
503  return (p.coefficients == coefficients);
504  }
505 
506 
507 
508  template <typename number>
509  template <typename number2>
510  void
511  Polynomial<number>::shift(std::vector<number> &coefficients,
512  const number2 offset)
513  {
514  // too many coefficients cause overflow in
515  // the binomial coefficient used below
516  Assert (coefficients.size() < 31, ExcNotImplemented());
517 
518  // Copy coefficients to a vector of
519  // accuracy given by the argument
520  std::vector<number2> new_coefficients(coefficients.begin(),
521  coefficients.end());
522 
523  // Traverse all coefficients from
524  // c_1. c_0 will be modified by
525  // higher degrees, only.
526  for (unsigned int d=1; d<new_coefficients.size(); ++d)
527  {
528  const unsigned int n = d;
529  // Binomial coefficients are
530  // needed for the
531  // computation. The rightmost
532  // value is unity.
533  unsigned int binomial_coefficient = 1;
534 
535  // Powers of the offset will be
536  // needed and computed
537  // successively.
538  number2 offset_power = offset;
539 
540  // Compute (x+offset)^d
541  // and modify all values c_k
542  // with k<d.
543  // The coefficient in front of
544  // x^d is not modified in this step.
545  for (unsigned int k=0; k<d; ++k)
546  {
547  // Recursion from Bronstein
548  // Make sure no remainders
549  // occur in integer
550  // division.
551  binomial_coefficient = (binomial_coefficient*(n-k))/(k+1);
552 
553  new_coefficients[d-k-1] += new_coefficients[d]
554  * binomial_coefficient
555  * offset_power;
556  offset_power *= offset;
557  }
558  // The binomial coefficient
559  // should have gone through a
560  // whole row of Pascal's
561  // triangle.
562  Assert (binomial_coefficient == 1, ExcInternalError());
563  }
564 
565  // copy new elements to old vector
566  coefficients.assign(new_coefficients.begin(), new_coefficients.end());
567  }
568 
569 
570 
571  template <typename number>
572  template <typename number2>
573  void
574  Polynomial<number>::shift (const number2 offset)
575  {
576  // shift is simple for a polynomial in product
577  // form, (x-x_0)*(x-x_1)*...*(x-x_n). just add
578  // offset to all shifts
579  if (in_lagrange_product_form == true)
580  {
581  for (unsigned int i=0; i<lagrange_support_points.size(); ++i)
582  lagrange_support_points[i] -= offset;
583  }
584  else
585  // do the shift in any case
586  shift (coefficients, offset);
587  }
588 
589 
590 
591  template <typename number>
594  {
595  // no simple form possible for Lagrange
596  // polynomial on product form
597  if (degree() == 0)
598  return Monomial<number>(0, 0.);
599 
600  std_cxx11::shared_ptr<Polynomial<number> > q_data;
601  const Polynomial<number> *q = 0;
602  if (in_lagrange_product_form == true)
603  {
604  q_data.reset (new Polynomial<number>(*this));
606  q = q_data.get();
607  }
608  else
609  q = this;
610 
611  std::vector<number> newcoefficients (q->coefficients.size()-1);
612  for (unsigned int i=1 ; i<q->coefficients.size() ; ++i)
613  newcoefficients[i-1] = number(i) * q->coefficients[i];
614 
615  return Polynomial<number> (newcoefficients);
616  }
617 
618 
619 
620  template <typename number>
623  {
624  // no simple form possible for Lagrange
625  // polynomial on product form
626  std_cxx11::shared_ptr<Polynomial<number> > q_data;
627  const Polynomial<number> *q = 0;
628  if (in_lagrange_product_form == true)
629  {
630  q_data.reset (new Polynomial<number>(*this));
632  q = q_data.get();
633  }
634  else
635  q = this;
636 
637  std::vector<number> newcoefficients (q->coefficients.size()+1);
638  newcoefficients[0] = 0.;
639  for (unsigned int i=0 ; i<q->coefficients.size() ; ++i)
640  newcoefficients[i+1] = q->coefficients[i]/number(i+1.);
641 
642  return Polynomial<number> (newcoefficients);
643  }
644 
645 
646 
647  template <typename number>
648  void
649  Polynomial<number>::print (std::ostream &out) const
650  {
651  if (in_lagrange_product_form == true)
652  {
653  out << lagrange_weight;
654  for (unsigned int i=0; i<lagrange_support_points.size(); ++i)
655  out << " (x-" << lagrange_support_points[i] << ")";
656  out << std::endl;
657  }
658  else
659  for (int i=degree(); i>=0; --i)
660  {
661  out << coefficients[i] << " x^" << i << std::endl;
662  }
663  }
664 
665 
666 
667 // ------------------ class Monomial -------------------------- //
668 
669  template <typename number>
670  std::vector<number>
672  double coefficient)
673  {
674  std::vector<number> result(n+1, 0.);
675  result[n] = coefficient;
676  return result;
677  }
678 
679 
680 
681  template <typename number>
682  Monomial<number>::Monomial (unsigned int n,
683  double coefficient)
684  : Polynomial<number>(make_vector(n, coefficient))
685  {}
686 
687 
688 
689  template <typename number>
690  std::vector<Polynomial<number> >
691  Monomial<number>::generate_complete_basis (const unsigned int degree)
692  {
693  std::vector<Polynomial<number> > v;
694  v.reserve(degree+1);
695  for (unsigned int i=0; i<=degree; ++i)
696  v.push_back (Monomial<number>(i));
697  return v;
698  }
699 
700 
701 
702 // ------------------ class LagrangeEquidistant --------------- //
703 
704  namespace internal
705  {
706  namespace LagrangeEquidistant
707  {
708  std::vector<Point<1> >
709  generate_equidistant_unit_points (const unsigned int n)
710  {
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;
715  return points;
716  }
717  }
718  }
719 
720 
721 
723  const unsigned int support_point)
724  :
726  generate_equidistant_unit_points (n),
727  support_point)
728  {
729  Assert (coefficients.size() == 0, ExcInternalError());
730 
731  // For polynomial order up to 3, we have precomputed weights. Use these
732  // weights instead of the product form
733  if (n <= 3)
734  {
735  this->in_lagrange_product_form = false;
736  this->lagrange_weight = 1.;
737  std::vector<double> new_support_points;
738  this->lagrange_support_points.swap(new_support_points);
739  this->coefficients.resize (n+1);
740  compute_coefficients(n, support_point, this->coefficients);
741  }
742  }
743 
744 
745 
746  void
748  const unsigned int support_point,
749  std::vector<double> &a)
750  {
751  Assert(support_point<n+1, ExcIndexRange(support_point, 0, n+1));
752 
753  unsigned int n_functions=n+1;
754  Assert(support_point<n_functions,
755  ExcIndexRange(support_point, 0, n_functions));
756  double const *x=0;
757 
758  switch (n)
759  {
760  case 1:
761  {
762  static const double x1[4]=
763  {
764  1.0, -1.0,
765  0.0, 1.0
766  };
767  x=&x1[0];
768  break;
769  }
770  case 2:
771  {
772  static const double x2[9]=
773  {
774  1.0, -3.0, 2.0,
775  0.0, 4.0, -4.0,
776  0.0, -1.0, 2.0
777  };
778  x=&x2[0];
779  break;
780  }
781  case 3:
782  {
783  static const double x3[16]=
784  {
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
789  };
790  x=&x3[0];
791  break;
792  }
793  default:
794  Assert(false, ExcInternalError())
795  }
796 
797  Assert(x!=0, ExcInternalError());
798  for (unsigned int i=0; i<n_functions; ++i)
799  a[i]=x[support_point*n_functions+i];
800  }
801 
802 
803 
804  std::vector<Polynomial<double> >
806  generate_complete_basis (const unsigned int degree)
807  {
808  if (degree==0)
809  // create constant polynomial
810  return std::vector<Polynomial<double> >
811  (1, Polynomial<double> (std::vector<double> (1,1.)));
812  else
813  {
814  // create array of Lagrange
815  // polynomials
816  std::vector<Polynomial<double> > v;
817  for (unsigned int i=0; i<=degree; ++i)
818  v.push_back(LagrangeEquidistant(degree,i));
819  return v;
820  }
821  }
822 
823 
824 
825 //----------------------------------------------------------------------//
826 
827 
828  std::vector<Polynomial<double> >
829  generate_complete_Lagrange_basis (const std::vector<Point<1> > &points)
830  {
831  std::vector<Polynomial<double> > p;
832  p.reserve (points.size());
833 
834  for (unsigned int i=0; i<points.size(); ++i)
835  p.push_back (Polynomial<double> (points, i));
836  return p;
837  }
838 
839 
840 
841 // ------------------ class Legendre --------------- //
842 
843 
844 
845  Legendre::Legendre (const unsigned int k)
846  :
847  Polynomial<double> (0)
848  {
849  this->coefficients.clear();
850  this->in_lagrange_product_form = true;
851  this->lagrange_support_points.resize(k);
852 
853  // the roots of a Legendre polynomial are exactly the points in the
854  // Gauss-Legendre quadrature formula
855  if (k > 0)
856  {
857  QGauss<1> gauss(k);
858  for (unsigned int i=0; i<k; ++i)
859  this->lagrange_support_points[i] = gauss.get_points()[i][0];
860  }
861 
862  // compute the abscissa in zero of the product of monomials. The exact
863  // value should be sqrt(2*k+1), so set the weight to that value.
864  double prod = 1.;
865  for (unsigned int i=0; i<k; ++i)
866  prod *= this->lagrange_support_points[i];
867  this->lagrange_weight = std::sqrt(double(2*k+1)) / prod;
868  }
869 
870 
871 
872  std::vector<Polynomial<double> >
873  Legendre::generate_complete_basis (const unsigned int degree)
874  {
875  std::vector<Polynomial<double> > v;
876  v.reserve(degree+1);
877  for (unsigned int i=0; i<=degree; ++i)
878  v.push_back (Legendre(i));
879  return v;
880  }
881 
882 
883 
884 // ------------------ class Lobatto -------------------- //
885 
886 
887  Lobatto::Lobatto (const unsigned int p) : Polynomial<double> (compute_coefficients (p))
888  {
889  }
890 
891  std::vector<double> Lobatto::compute_coefficients (const unsigned int p)
892  {
893  switch (p)
894  {
895  case 0:
896  {
897  std::vector<double> coefficients (2);
898 
899  coefficients[0] = 1.0;
900  coefficients[1] = -1.0;
901  return coefficients;
902  }
903 
904  case 1:
905  {
906  std::vector<double> coefficients (2);
907 
908  coefficients[0] = 0.0;
909  coefficients[1] = 1.0;
910  return coefficients;
911  }
912 
913  case 2:
914  {
915  std::vector<double> coefficients (3);
916 
917  coefficients[0] = 0.0;
918  coefficients[1] = -1.0 * std::sqrt (3.);
919  coefficients[2] = std::sqrt (3.);
920  return coefficients;
921  }
922 
923  default:
924  {
925  std::vector<double> coefficients (p + 1);
926  std::vector<double> legendre_coefficients_tmp1 (p);
927  std::vector<double> legendre_coefficients_tmp2 (p - 1);
928 
929  coefficients[0] = -1.0 * std::sqrt (3.);
930  coefficients[1] = 2.0 * std::sqrt (3.);
931  legendre_coefficients_tmp1[0] = 1.0;
932 
933  for (unsigned int i = 2; i < p; ++i)
934  {
935  for (unsigned int j = 0; j < i - 1; ++j)
936  legendre_coefficients_tmp2[j] = legendre_coefficients_tmp1[j];
937 
938  for (unsigned int j = 0; j < i; ++j)
939  legendre_coefficients_tmp1[j] = coefficients[j];
940 
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;
942 
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;
945 
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;
948  }
949 
950  for (int i = p; i > 0; --i)
951  coefficients[i] = coefficients[i - 1] / i;
952 
953  coefficients[0] = 0.0;
954  return coefficients;
955  }
956  }
957  }
958 
959  std::vector<Polynomial<double> > Lobatto::generate_complete_basis (const unsigned int p)
960  {
961  std::vector<Polynomial<double> > basis (p + 1);
962 
963  for (unsigned int i = 0; i <= p; ++i)
964  basis[i] = Lobatto (i);
965 
966  return basis;
967  }
968 
969 
970 
971 // ------------------ class Hierarchical --------------- //
972 
973 
974 // Reserve space for polynomials up to degree 19. Should be sufficient
975 // for the start.
976  std::vector<std_cxx11::shared_ptr<const std::vector<double> > >
978 
979 
980 
981  Hierarchical::Hierarchical (const unsigned int k)
982  :
983  Polynomial<double> (get_coefficients(k))
984  {}
985 
986 
987 
988  void
989  Hierarchical::compute_coefficients (const unsigned int k_)
990  {
991  unsigned int k = k_;
992 
993  // first make sure that no other
994  // thread intercepts the operation
995  // of this function
996  // for this, acquire the lock
997  // until we quit this function
998  Threads::Mutex::ScopedLock lock(coefficients_lock);
999 
1000  // The first 2 coefficients
1001  // are hard-coded
1002  if (k==0)
1003  k=1;
1004  // check: does the information
1005  // already exist?
1006  if ( (recursive_coefficients.size() < k+1) ||
1007  ((recursive_coefficients.size() >= k+1) &&
1008  (recursive_coefficients[k].get() == 0)) )
1009  // no, then generate the
1010  // respective coefficients
1011  {
1012  // make sure that there is enough
1013  // space in the array for the
1014  // coefficients, so we have to resize
1015  // it to size k+1
1016 
1017  // but it's more complicated than
1018  // that: we call this function
1019  // recursively, so if we simply
1020  // resize it to k+1 here, then
1021  // compute the coefficients for
1022  // degree k-1 by calling this
1023  // function recursively, then it will
1024  // reset the size to k -- not enough
1025  // for what we want to do below. the
1026  // solution therefore is to only
1027  // resize the size if we are going to
1028  // *increase* it
1029  if (recursive_coefficients.size() < k+1)
1030  recursive_coefficients.resize (k+1);
1031 
1032  if (k<=1)
1033  {
1034  // create coefficients
1035  // vectors for k=0 and k=1
1036  //
1037  // allocate the respective
1038  // amount of memory and
1039  // later assign it to the
1040  // coefficients array to
1041  // make it const
1042  std::vector<double> *c0 = new std::vector<double>(2);
1043  (*c0)[0] = 1.;
1044  (*c0)[1] = -1.;
1045 
1046  std::vector<double> *c1 = new std::vector<double>(2);
1047  (*c1)[0] = 0.;
1048  (*c1)[1] = 1.;
1049 
1050  // now make these arrays
1051  // const
1053  std_cxx11::shared_ptr<const std::vector<double> >(c0);
1055  std_cxx11::shared_ptr<const std::vector<double> >(c1);
1056  }
1057  else if (k==2)
1058  {
1059  coefficients_lock.release ();
1061  coefficients_lock.acquire ();
1062 
1063  std::vector<double> *c2 = new std::vector<double>(3);
1064 
1065  const double a = 1.; //1./8.;
1066 
1067  (*c2)[0] = 0.*a;
1068  (*c2)[1] = -4.*a;
1069  (*c2)[2] = 4.*a;
1070 
1072  std_cxx11::shared_ptr<const std::vector<double> >(c2);
1073  }
1074  else
1075  {
1076  // for larger numbers,
1077  // compute the coefficients
1078  // recursively. to do so,
1079  // we have to release the
1080  // lock temporarily to
1081  // allow the called
1082  // function to acquire it
1083  // itself
1084  coefficients_lock.release ();
1085  compute_coefficients(k-1);
1086  coefficients_lock.acquire ();
1087 
1088  std::vector<double> *ck = new std::vector<double>(k+1);
1089 
1090  const double a = 1.; //1./(2.*k);
1091 
1092  (*ck)[0] = - a*(*recursive_coefficients[k-1])[0];
1093 
1094  for (unsigned int i=1; i<=k-1; ++i)
1095  (*ck)[i] = a*( 2.*(*recursive_coefficients[k-1])[i-1]
1096  - (*recursive_coefficients[k-1])[i] );
1097 
1098  (*ck)[k] = a*2.*(*recursive_coefficients[k-1])[k-1];
1099  // for even degrees, we need
1100  // to add a multiple of
1101  // basis fcn phi_2
1102  if ( (k%2) == 0 )
1103  {
1104  double b = 1.; //8.;
1105  //for (unsigned int i=1; i<=k; i++)
1106  // b /= 2.*i;
1107 
1108  (*ck)[1] += b*(*recursive_coefficients[2])[1];
1109  (*ck)[2] += b*(*recursive_coefficients[2])[2];
1110  }
1111  // finally assign the newly
1112  // created vector to the
1113  // const pointer in the
1114  // coefficients array
1116  std_cxx11::shared_ptr<const std::vector<double> >(ck);
1117  };
1118  };
1119  }
1120 
1121 
1122 
1123  const std::vector<double> &
1124  Hierarchical::get_coefficients (const unsigned int k)
1125  {
1126  // first make sure the coefficients
1127  // get computed if so necessary
1129 
1130  // then get a pointer to the array
1131  // of coefficients. do that in a MT
1132  // safe way
1133  Threads::Mutex::ScopedLock lock (coefficients_lock);
1134  return *recursive_coefficients[k];
1135  }
1136 
1137 
1138 
1139  std::vector<Polynomial<double> >
1140  Hierarchical::generate_complete_basis (const unsigned int degree)
1141  {
1142  if (degree==0)
1143  // create constant
1144  // polynomial. note that we
1145  // can't use the other branch
1146  // of the if-statement, since
1147  // calling the constructor of
1148  // this class with argument
1149  // zero does _not_ create the
1150  // constant polynomial, but
1151  // rather 1-x
1152  return std::vector<Polynomial<double> >
1153  (1, Polynomial<double> (std::vector<double> (1,1.)));
1154  else
1155  {
1156  std::vector<Polynomial<double> > v;
1157  v.reserve(degree+1);
1158  for (unsigned int i=0; i<=degree; ++i)
1159  v.push_back (Hierarchical(i));
1160  return v;
1161  }
1162  }
1163 
1164 
1165 
1166 // ------------------ HermiteInterpolation --------------- //
1167 
1169  :
1170  Polynomial<double>(0)
1171  {
1172  this->coefficients.clear();
1173  this->in_lagrange_product_form = true;
1174 
1175  this->lagrange_support_points.resize(3);
1176  if (p==0)
1177  {
1178  this->lagrange_support_points[0] = -0.5;
1179  this->lagrange_support_points[1] = 1.;
1180  this->lagrange_support_points[2] = 1.;
1181  this->lagrange_weight = 2.;
1182  }
1183  else if (p==1)
1184  {
1185  this->lagrange_support_points[0] = 0.;
1186  this->lagrange_support_points[1] = 0.;
1187  this->lagrange_support_points[2] = 1.5;
1188  this->lagrange_weight = -2.;
1189  }
1190  else if (p==2)
1191  {
1192  this->lagrange_support_points[0] = 0.;
1193  this->lagrange_support_points[1] = 1.;
1194  this->lagrange_support_points[2] = 1.;
1195  }
1196  else if (p==3)
1197  {
1198  this->lagrange_support_points[0] = 0.;
1199  this->lagrange_support_points[1] = 0.;
1200  this->lagrange_support_points[2] = 1.;
1201  }
1202  else
1203  {
1204  this->lagrange_support_points.resize(4);
1205  this->lagrange_support_points[0] = 0.;
1206  this->lagrange_support_points[1] = 0.;
1207  this->lagrange_support_points[2] = 1.;
1208  this->lagrange_support_points[3] = 1.;
1209  this->lagrange_weight = 16.;
1210 
1211  if (p>4)
1212  {
1213  Legendre legendre(p-4);
1214  (*this) *= legendre;
1215  }
1216  }
1217  }
1218 
1219 
1220  std::vector<Polynomial<double> >
1222  {
1223  std::vector<Polynomial<double> > basis (n + 1);
1224 
1225  for (unsigned int i = 0; i <= n; ++i)
1226  basis[i] = HermiteInterpolation (i);
1227 
1228  return basis;
1229  }
1230 
1231 }
1232 
1233 // ------------------ explicit instantiations --------------- //
1234 
1235 namespace Polynomials
1236 {
1237  template class Polynomial<float>;
1238  template class Polynomial<double>;
1239  template class Polynomial<long double>;
1240 
1241  template void Polynomial<float>::shift(const float offset);
1242  template void Polynomial<float>::shift(const double offset);
1243  template void Polynomial<double>::shift(const double offset);
1244  template void Polynomial<long double>::shift(const long double offset);
1245  template void Polynomial<float>::shift(const long double offset);
1246  template void Polynomial<double>::shift(const long double offset);
1247 
1248  template class Monomial<float>;
1249  template class Monomial<double>;
1250  template class Monomial<long double>;
1251 }
1252 
1253 DEAL_II_NAMESPACE_CLOSE
void transform_into_standard_form()
Definition: polynomial.cc:247
LagrangeEquidistant(const unsigned int n, const unsigned int support_point)
Definition: polynomial.cc:722
void shift(const number2 offset)
Definition: polynomial.cc:574
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int p)
Definition: polynomial.cc:959
const std::vector< Point< dim > > & get_points() const
unsigned int degree() const
Definition: polynomial.h:610
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
Definition: polynomial.cc:1140
Legendre(const unsigned int p)
Definition: polynomial.cc:845
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
Definition: polynomial.cc:873
#define AssertIndexRange(index, range)
Definition: exceptions.h:1170
std::vector< Polynomial< double > > generate_complete_Lagrange_basis(const std::vector< Point< 1 > > &points)
Definition: polynomial.cc:829
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)
Definition: polynomial.cc:747
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
Definition: polynomial.cc:806
Definition: point.h:89
std::vector< double > compute_coefficients(const unsigned int p)
Definition: polynomial.cc:891
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:313
HermiteInterpolation(const unsigned int p)
Definition: polynomial.cc:1168
Lobatto(const unsigned int p=0)
Definition: polynomial.cc:887
static std::vector< Polynomial< number > > generate_complete_basis(const unsigned int degree)
Definition: polynomial.cc:691
Hierarchical(const unsigned int p)
Definition: polynomial.cc:981
std::vector< number > coefficients
Definition: polynomial.h:250
static void compute_coefficients(const unsigned int p)
Definition: polynomial.cc:989
static const std::vector< double > & get_coefficients(const unsigned int p)
Definition: polynomial.cc:1124
std::vector< number > lagrange_support_points
Definition: polynomial.h:262
static ::ExceptionBase & ExcEmptyObject()
static ::ExceptionBase & ExcNotImplemented()
static std::vector< std_cxx11::shared_ptr< const std::vector< double > > > recursive_coefficients
Definition: polynomial.h:537
static ::ExceptionBase & ExcZero()
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int p)
Definition: polynomial.cc:1221
static ::ExceptionBase & ExcInternalError()