Reference documentation for deal.II version 9.4.1
|
#include <deal.II/base/polynomial.h>
Public Member Functions | |
Hierarchical (const unsigned int p) | |
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 p) |
static const std::vector< double > & | get_coefficients (const unsigned int p) |
Static Private Attributes | |
static std::vector< std::unique_ptr< const std::vector< double > > > | recursive_coefficients |
Subscriptor functionality | |
Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class. | |
static std::mutex | mutex |
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 ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (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 int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
void | check_no_subscribers () const noexcept |
Hierarchical polynomials of arbitrary degree on [0,1]
.
When Constructing a Hierarchical polynomial of degree p
, the coefficients will be computed by a recursion formula. The coefficients are stored in a static data vector to be available when needed next time.
These hierarchical polynomials are based on those of Demkowicz, Oden, Rachowicz, and Hardy (CMAME 77 (1989) 79-112, Sec. 4). The first two polynomials are the standard linear shape functions given by \(\phi_{0}(x) = 1 - x\) and \(\phi_{1}(x) = x\). For \(l \geq 2\) we use the definitions \(\phi_{l}(x) = (2x-1)^l - 1, l = 2,4,6,...\) and \(\phi_{l}(x) = (2x-1)^l - (2x-1), l = 3,5,7,...\). These satisfy the recursion relations \(\phi_{l}(x) = (2x-1)\phi_{l-1}, l=3,5,7,...\) and \(\phi_{l}(x) = (2x-1)\phi_{l-1} + \phi_{2}, l=4,6,8,...\).
The degrees of freedom are the values at the vertices and the derivatives at the midpoint. Currently, we do not scale the polynomials in any way, although better conditioning of the element stiffness matrix could possibly be achieved with scaling.
Calling the constructor with a given index p
will generate the following: if p==0
, then the resulting polynomial is the linear function associated with the left vertex, if p==1
the one associated with the right vertex. For higher values of p
, you get the polynomial of degree p
that is orthogonal to all previous ones. Note that for p==0
you therefore do not get a polynomial of degree zero, but one of degree one. This is to allow generating a complete basis for polynomial spaces, by just iterating over the indices given to the constructor.
On the other hand, the function generate_complete_basis() creates a complete basis of given degree. In order to be consistent with the concept of a polynomial degree, if the given argument is zero, it does not return the linear polynomial described above, but rather a constant polynomial.
Definition at line 528 of file polynomial.h.
Polynomials::Hierarchical::Hierarchical | ( | const unsigned int | p | ) |
Constructor for polynomial of degree p
. There is an exception for p==0
, see the general documentation.
Definition at line 872 of file polynomial.cc.
|
static |
Return a vector of Hierarchical polynomial objects of degrees zero through degree
, which then spans the full space of polynomials up to the given degree. Note that there is an exception if the given degree
equals zero, see the general documentation of this class.
This function may be used to initialize the TensorProductPolynomials, AnisotropicPolynomials, and PolynomialSpace classes.
Definition at line 1029 of file polynomial.cc.
|
staticprivate |
Compute coefficients recursively.
Definition at line 879 of file polynomial.cc.
|
staticprivate |
Get coefficients for constructor. This way, it can use the non- standard constructor of Polynomial.
Definition at line 1013 of file polynomial.cc.
|
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.
|
inherited |
Return the values and the derivatives of the Polynomial at point x
. values[i], i=0,...,values.size()-1
includes the i
th 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.
|
inlineinherited |
Return the values and the derivatives of the Polynomial at point x
. values[i], i=0,...,n_derivatives
includes the i
th 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.
|
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.
|
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.
|
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.
|
staticprotectedinherited |
This function performs the actual scaling.
Definition at line 270 of file polynomial.cc.
|
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.
|
staticprotectedinherited |
This function performs the actual shift
Definition at line 277 of file polynomial.cc.
|
inherited |
Compute the derivative of a polynomial.
Definition at line 205 of file polynomial.cc.
|
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.
|
inherited |
Multiply with a scalar.
Definition at line 218 of file polynomial.cc.
|
inherited |
Multiply with another polynomial.
Definition at line 224 of file polynomial.cc.
|
inherited |
Add a second polynomial.
Definition at line 230 of file polynomial.cc.
|
inherited |
Subtract a second polynomial.
Definition at line 236 of file polynomial.cc.
|
inherited |
Test for equality of two polynomials.
Definition at line 242 of file polynomial.cc.
|
inherited |
Print coefficients.
Definition at line 248 of file polynomial.cc.
|
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.
|
virtualinherited |
Return an estimate (in bytes) for the memory consumption of this object.
Definition at line 263 of file polynomial.cc.
|
staticprotectedinherited |
Multiply polynomial by a factor.
Definition at line 283 of file polynomial.cc.
|
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.
|
staticprivate |
Vector with already computed coefficients. For each degree of the polynomial, we keep one pointer to the list of coefficients; we do so rather than keeping a vector of vectors in order to simplify programming multithread-safe. In order to avoid memory leak, we use a unique_ptr in order to correctly free the memory of the vectors when the global destructor is called.
Definition at line 573 of file polynomial.h.
|
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.
|
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.
|
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.
|
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.