Reference documentation for deal.II version 9.1.1
|
#include <deal.II/base/polynomial.h>
Public Member Functions | |
Hierarchical (const unsigned int p) | |
Public Member Functions inherited from Polynomials::Polynomial< double > | |
Polynomial (const std::vector< double > &coefficients) | |
Polynomial (const unsigned int n) | |
Polynomial (const std::vector< Point< 1 >> &lagrange_support_points, const unsigned int evaluation_point) | |
Polynomial () | |
double | value (const double x) const |
void | value (const double x, std::vector< double > &values) const |
void | value (const double x, const unsigned int n_derivatives, double *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) |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
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 |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Static Public Member Functions | |
static std::vector< Polynomial< double > > | generate_complete_basis (const unsigned int degree) |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
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 |
Additional Inherited Members | |
Protected Member Functions inherited from Polynomials::Polynomial< double > | |
void | transform_into_standard_form () |
Static Protected Member Functions inherited from Polynomials::Polynomial< double > | |
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 inherited from Polynomials::Polynomial< double > | |
std::vector< double > | coefficients |
bool | in_lagrange_product_form |
std::vector< double > | lagrange_support_points |
double | lagrange_weight |
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 501 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 993 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 1150 of file polynomial.cc.
|
staticprivate |
Compute coefficients recursively.
Definition at line 1000 of file polynomial.cc.
|
staticprivate |
Get coefficients for constructor. This way, it can use the non- standard constructor of Polynomial.
Definition at line 1134 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 546 of file polynomial.h.