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