Reference documentation for deal.II version 9.4.1
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Member Functions | Static Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Static Private Member Functions | List of all members
Polynomials::LagrangeEquidistant Class Reference

#include <deal.II/base/polynomial.h>

Inheritance diagram for Polynomials::LagrangeEquidistant:
[legend]

Public Member Functions

 LagrangeEquidistant (const unsigned int n, const unsigned int support_point)
 
double value (const double x) const
 
void value (const double x, std::vector< double > &values) const
 
void value (const Number2 x, const unsigned int n_derivatives, Number2 *values) const
 
void values_of_array (const std::array< Number2, n_entries > &points, const unsigned int n_derivatives, std::array< Number2, n_entries > *values) const
 
unsigned int degree () const
 
void scale (const double factor)
 
void shift (const number2 offset)
 
Polynomial< double > derivative () const
 
Polynomial< double > primitive () const
 
Polynomial< double > & operator*= (const double s)
 
Polynomial< double > & operator*= (const Polynomial< double > &p)
 
Polynomial< double > & operator+= (const Polynomial< double > &p)
 
Polynomial< double > & operator-= (const Polynomial< double > &p)
 
bool operator== (const Polynomial< double > &p) const
 
void print (std::ostream &out) const
 
void serialize (Archive &ar, const unsigned int version)
 
virtual std::size_t memory_consumption () const
 

Static Public Member Functions

static std::vector< Polynomial< double > > generate_complete_basis (const unsigned int degree)
 

Protected Member Functions

void transform_into_standard_form ()
 

Static Protected Member Functions

static void scale (std::vector< double > &coefficients, const double factor)
 
static void shift (std::vector< double > &coefficients, const number2 shift)
 
static void multiply (std::vector< double > &coefficients, const double factor)
 

Protected Attributes

std::vector< double > coefficients
 
bool in_lagrange_product_form
 
std::vector< double > lagrange_support_points
 
double lagrange_weight
 

Static Private Member Functions

static void compute_coefficients (const unsigned int n, const unsigned int support_point, std::vector< double > &a)
 

Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 
static std::mutex mutex
 
void check_no_subscribers () const noexcept
 

Detailed Description

Lagrange polynomials with equidistant interpolation points in [0,1]. The polynomial of degree n has got n+1 interpolation points. The interpolation points are sorted in ascending order. This order gives an index to each interpolation point. A Lagrangian polynomial equals to 1 at its ‘support point’, and 0 at all other interpolation points. For example, if the degree is 3, and the support point is 1, then the polynomial represented by this object is cubic and its value is 1 at the point x=1/3, and zero at the point x=0, x=2/3, and x=1. All the polynomials have polynomial degree equal to degree, but together they span the entire space of polynomials of degree less than or equal degree.

The Lagrange polynomials are implemented up to degree 10.

Definition at line 370 of file polynomial.h.

Constructor & Destructor Documentation

◆ LagrangeEquidistant()

Polynomials::LagrangeEquidistant::LagrangeEquidistant ( const unsigned int  n,
const unsigned int  support_point 
)

Constructor. Takes the degree n of the Lagrangian polynomial and the index support_point of the support point. Fills the coefficients of the base class Polynomial.

Definition at line 597 of file polynomial.cc.

Member Function Documentation

◆ generate_complete_basis()

std::vector< Polynomial< double > > Polynomials::LagrangeEquidistant::generate_complete_basis ( const unsigned int  degree)
static

Return a vector of polynomial objects of degree degree, which then spans the full space of polynomials up to the given degree. The polynomials are generated by calling the constructor of this class with the same degree but support point running from zero to degree. This function may be used to initialize the TensorProductPolynomials and PolynomialSpace classes.

Definition at line 679 of file polynomial.cc.

◆ compute_coefficients()

void Polynomials::LagrangeEquidistant::compute_coefficients ( const unsigned int  n,
const unsigned int  support_point,
std::vector< double > &  a 
)
staticprivate

Compute the coefficients of the base class Polynomial. This function is static to allow to be called in the constructor.

Definition at line 621 of file polynomial.cc.

◆ value() [1/3]

double Polynomials::Polynomial< double >::value ( const double  x) const
inlineinherited

Return the value of this polynomial at the given point.

This function uses the most numerically stable evaluation algorithm for the provided form of the polynomial. If the polynomial is in the product form of roots, the evaluation is based on products of the form (x - x_i), whereas the Horner scheme is used for polynomials in the coefficient form.

Definition at line 107 of file polynomial.h.

◆ value() [2/3]

void Polynomials::Polynomial< double >::value ( const double  x,
std::vector< double > &  values 
) const
inherited

Return the values and the derivatives of the Polynomial at point x. values[i], i=0,...,values.size()-1 includes the ith derivative. The number of derivatives to be computed is thus determined by the size of the array passed.

This function uses the Horner scheme for numerical stability of the evaluation for polynomials in the coefficient form or the product of terms involving the roots if that representation is used.

Definition at line 120 of file polynomial.cc.

◆ value() [3/3]

void Polynomials::Polynomial< double >::value ( const Number2  x,
const unsigned int  n_derivatives,
Number2 *  values 
) const
inlineinherited

Return the values and the derivatives of the Polynomial at point x. values[i], i=0,...,n_derivatives includes the ith derivative. The number of derivatives to be computed is determined by n_derivatives and values has to provide sufficient space for n_derivatives + 1 values.

This function uses the most numerically stable evaluation algorithm for the provided form of the polynomial. If the polynomial is in the product form of roots, the evaluation is based on products of the form (x - x_i), whereas the Horner scheme is used for polynomials in the coefficient form.

The template type Number2 must implement arithmetic operations such as additions or multiplication with the type number of the polynomial, and must be convertible from number by operator=.

Definition at line 142 of file polynomial.h.

◆ values_of_array()

void Polynomials::Polynomial< double >::values_of_array ( const std::array< Number2, n_entries > &  points,
const unsigned int  n_derivatives,
std::array< Number2, n_entries > *  values 
) const
inlineinherited

Similar to the function above, but evaluate the polynomials on several positions at once, as described by the array argument points. This function is can be faster than the other function when the same polynomial should be evaluated on several positions at once, e.g., the x,y,z coordinates of a point for tensor-product polynomials.

The template type Number2 must implement arithmetic operations such as additions or multiplication with the type number of the polynomial, and must be convertible from number by operator=.

Definition at line 160 of file polynomial.h.

◆ degree()

unsigned int Polynomials::Polynomial< double >::degree
inlineinherited

Degree of the polynomial. This is the degree reflected by the number of coefficients provided by the constructor. Leading non-zero coefficients are not treated separately.

Definition at line 170 of file polynomial.h.

◆ scale() [1/2]

void Polynomials::Polynomial< double >::scale ( const double  factor)
inherited

Scale the abscissa of the polynomial. Given the polynomial p(t) and the scaling t = ax, then the result of this operation is the polynomial q, such that q(x) = p(t).

The operation is performed in place.

Definition at line 180 of file polynomial.cc.

◆ scale() [2/2]

void Polynomials::Polynomial< double >::scale ( std::vector< double > &  coefficients,
const double  factor 
)
staticprotectedinherited

This function performs the actual scaling.

Definition at line 270 of file polynomial.cc.

◆ shift() [1/2]

void Polynomials::Polynomial< double >::shift ( const number2  offset)
inherited

Shift the abscissa oft the polynomial. Given the polynomial p(t) and the shift t = x + a, then the result of this operation is the polynomial q, such that q(x) = p(t).

The template parameter allows to compute the new coefficients with higher accuracy, since all computations are performed with type number2. This may be necessary, since this operation involves a big number of additions. On a Sun Sparc Ultra with Solaris 2.8, the difference between double and long double was not significant, though.

The operation is performed in place, i.e. the coefficients of the present object are changed.

Definition at line 199 of file polynomial.cc.

◆ shift() [2/2]

void Polynomials::Polynomial< double >::shift ( std::vector< double > &  coefficients,
const number2  shift 
)
staticprotectedinherited

This function performs the actual shift

Definition at line 277 of file polynomial.cc.

◆ derivative()

Polynomial< double > Polynomials::Polynomial< double >::derivative
inherited

Compute the derivative of a polynomial.

Definition at line 205 of file polynomial.cc.

◆ primitive()

Polynomial< double > Polynomials::Polynomial< double >::primitive
inherited

Compute the primitive of a polynomial. the coefficient of the zero order term of the polynomial is zero.

Definition at line 212 of file polynomial.cc.

◆ operator*=() [1/2]

Polynomial< double > & Polynomials::Polynomial< double >::operator*= ( const double  s)
inherited

Multiply with a scalar.

Definition at line 218 of file polynomial.cc.

◆ operator*=() [2/2]

Polynomial< double > & Polynomials::Polynomial< double >::operator*= ( const Polynomial< double > &  p)
inherited

Multiply with another polynomial.

Definition at line 224 of file polynomial.cc.

◆ operator+=()

Polynomial< double > & Polynomials::Polynomial< double >::operator+= ( const Polynomial< double > &  p)
inherited

Add a second polynomial.

Definition at line 230 of file polynomial.cc.

◆ operator-=()

Polynomial< double > & Polynomials::Polynomial< double >::operator-= ( const Polynomial< double > &  p)
inherited

Subtract a second polynomial.

Definition at line 236 of file polynomial.cc.

◆ operator==()

bool Polynomials::Polynomial< double >::operator== ( const Polynomial< double > &  p) const
inherited

Test for equality of two polynomials.

Definition at line 242 of file polynomial.cc.

◆ print()

void Polynomials::Polynomial< double >::print ( std::ostream &  out) const
inherited

Print coefficients.

Definition at line 248 of file polynomial.cc.

◆ serialize()

void Polynomials::Polynomial< double >::serialize ( Archive &  ar,
const unsigned int  version 
)
inlineinherited

Write or read the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.

Definition at line 257 of file polynomial.h.

◆ memory_consumption()

std::size_t Polynomials::Polynomial< double >::memory_consumption
virtualinherited

Return an estimate (in bytes) for the memory consumption of this object.

Definition at line 263 of file polynomial.cc.

◆ multiply()

void Polynomials::Polynomial< double >::multiply ( std::vector< double > &  coefficients,
const double  factor 
)
staticprotectedinherited

Multiply polynomial by a factor.

Definition at line 283 of file polynomial.cc.

◆ transform_into_standard_form()

void Polynomials::Polynomial< double >::transform_into_standard_form
protectedinherited

Transform polynomial form of product of linear factors into standard form, \(\sum_i a_i x^i\). Deletes all data structures related to the product form.

Definition at line 291 of file polynomial.cc.

Member Data Documentation

◆ coefficients

std::vector<double > Polynomials::Polynomial< double >::coefficients
protectedinherited

Coefficients of the polynomial \(\sum_i a_i x^i\). This vector is filled by the constructor of this class and may be passed down by derived classes.

This vector cannot be constant since we want to allow copying of polynomials.

Definition at line 301 of file polynomial.h.

◆ in_lagrange_product_form

bool Polynomials::Polynomial< double >::in_lagrange_product_form
protectedinherited

Stores whether the polynomial is in Lagrange product form, i.e., constructed as a product \((x-x_0) (x-x_1) \ldots (x-x_n)/c\), or not.

Definition at line 307 of file polynomial.h.

◆ lagrange_support_points

std::vector<double > Polynomials::Polynomial< double >::lagrange_support_points
protectedinherited

If the polynomial is in Lagrange product form, i.e., constructed as a product \((x-x_0) (x-x_1) \ldots (x-x_n)/c\), store the shifts \(x_i\).

Definition at line 313 of file polynomial.h.

◆ lagrange_weight

double Polynomials::Polynomial< double >::lagrange_weight
protectedinherited

If the polynomial is in Lagrange product form, i.e., constructed as a product \((x-x_0) (x-x_1) \ldots (x-x_n)/c\), store the weight c.

Definition at line 319 of file polynomial.h.


The documentation for this class was generated from the following files: