Reference documentation for deal.II version 9.5.0
|
#include <deal.II/base/polynomials_raviart_thomas.h>
Public Member Functions | |
PolynomialsRaviartThomas (const unsigned int degree_normal, const unsigned int degree_tangential) | |
PolynomialsRaviartThomas (const unsigned int k) | |
PolynomialsRaviartThomas (const PolynomialsRaviartThomas &other)=default | |
void | evaluate (const Point< dim > &unit_point, std::vector< Tensor< 1, dim > > &values, std::vector< Tensor< 2, dim > > &grads, std::vector< Tensor< 3, dim > > &grad_grads, std::vector< Tensor< 4, dim > > &third_derivatives, std::vector< Tensor< 5, dim > > &fourth_derivatives) const override |
std::string | name () const override |
virtual std::unique_ptr< TensorPolynomialsBase< dim > > | clone () const override |
std::vector< Point< dim > > | get_polynomial_support_points () const |
unsigned int | n () const |
unsigned int | degree () const |
Static Public Member Functions | |
static unsigned int | n_polynomials (const unsigned int normal_degree, const unsigned int tangential_degree) |
static unsigned int | n_polynomials (const unsigned int degree) |
static std::vector< unsigned int > | get_lexicographic_numbering (const unsigned int normal_degree, const unsigned int tangential_degree) |
Private Attributes | |
const unsigned int | normal_degree |
const unsigned int | tangential_degree |
const AnisotropicPolynomials< dim > | polynomial_space |
std::vector< unsigned int > | lexicographic_to_hierarchic |
std::vector< unsigned int > | hierarchic_to_lexicographic |
std::array< std::vector< unsigned int >, dim > | renumber_aniso |
const unsigned int | polynomial_degree |
const unsigned int | n_pols |
This class implements the Hdiv-conforming, vector-valued Raviart-Thomas polynomials as described in the book by Brezzi and Fortin.
The Raviart-Thomas polynomials are constructed such that the divergence is in the tensor product polynomial space Qk. Therefore, the polynomial order of each component must be one order higher in the corresponding direction, yielding the polynomial spaces (Qk+1,k, Qk,k+1) and (Qk+1,k,k, Qk,k+1,k, Qk,k,k+1) in 2d and 3d, resp.
Definition at line 50 of file polynomials_raviart_thomas.h.
PolynomialsRaviartThomas< dim >::PolynomialsRaviartThomas | ( | const unsigned int | degree_normal, |
const unsigned int | degree_tangential | ||
) |
Constructor. Creates all basis functions for Raviart-Thomas polynomials of given degree in normal and tangential directions. The usual FE_RaviartThomas and FE_RaviartThomasNodal classes will use degree + 1
and degree
in the two directions, respectively.
Definition at line 64 of file polynomials_raviart_thomas.cc.
PolynomialsRaviartThomas< dim >::PolynomialsRaviartThomas | ( | const unsigned int | k | ) |
Constructor, using the common Raviart-Thomas space of degree k + 1
in normal direction and k
in the tangential directions.
Definition at line 126 of file polynomials_raviart_thomas.cc.
|
default |
Copy constructor.
|
overridevirtual |
Compute the value and derivatives of each Raviart-Thomas polynomial at unit_point
.
The size of the vectors must either be zero or equal n()
. In the first case, the function will not compute these values.
If you need values or derivatives of all tensor product polynomials then use this function, rather than using any of the compute_value
, compute_grad
or compute_grad_grad
functions, see below, in a loop over all tensor product polynomials.
Implements TensorPolynomialsBase< dim >.
Definition at line 134 of file polynomials_raviart_thomas.cc.
|
overridevirtual |
Return the name of the space, which is PolynomialsRaviartThomas
.
Implements TensorPolynomialsBase< dim >.
Definition at line 228 of file polynomials_raviart_thomas.cc.
|
static |
Return the number of polynomials in the space without requiring to build an object of PolynomialsRaviartThomas. This is required by the FiniteElement classes.
Definition at line 246 of file polynomials_raviart_thomas.cc.
|
static |
Variant of the n_polynomials() function taking only a single argument degree
, assuming degree + 1
in the normal direction and degree
in the tangential directions.
Definition at line 237 of file polynomials_raviart_thomas.cc.
|
static |
Compute the lexicographic to hierarchic numbering underlying this class, computed as a free function.
Definition at line 258 of file polynomials_raviart_thomas.cc.
|
overridevirtual |
A sort of virtual copy constructor, this function returns a copy of the polynomial space object. Derived classes need to override the function here in this base class and return an object of the same type as the derived class.
Some places in the library, for example the constructors of FE_PolyTensor, need to make copies of polynomial spaces without knowing their exact type. They do so through this function.
Implements TensorPolynomialsBase< dim >.
Definition at line 323 of file polynomials_raviart_thomas.cc.
std::vector< Point< dim > > PolynomialsRaviartThomas< dim >::get_polynomial_support_points |
Compute the generalized support points in the ordering used by the polynomial shape functions. Note that these points are not support points in the classical sense as the Lagrange polynomials of the different components have different points, which need to be combined in terms of Piola transforms.
Definition at line 332 of file polynomials_raviart_thomas.cc.
|
inlineinherited |
Return the number of polynomials.
Definition at line 157 of file tensor_polynomials_base.h.
|
inlineinherited |
Return the highest polynomial degree of polynomials represented by this class. A derived class may override this if its value is different from my_degree
.
Definition at line 166 of file tensor_polynomials_base.h.
|
private |
The given degree in the normal direction.
Definition at line 148 of file polynomials_raviart_thomas.h.
|
private |
The given degree in the tangential direction.
Definition at line 153 of file polynomials_raviart_thomas.h.
|
private |
An object representing the polynomial space for a single component. We can re-use it by rotating the coordinates of the evaluation point.
Definition at line 159 of file polynomials_raviart_thomas.h.
|
private |
Renumbering from lexicographic to hierarchic order.
Definition at line 164 of file polynomials_raviart_thomas.h.
|
private |
Renumbering from hierarchic to lexicographic order. Inverse of lexicographic_to_hierarchic.
Definition at line 170 of file polynomials_raviart_thomas.h.
|
private |
Renumbering from shifted polynomial spaces to lexicographic one.
Definition at line 175 of file polynomials_raviart_thomas.h.
|
privateinherited |
The highest polynomial degree of this functions represented by this object.
Definition at line 145 of file tensor_polynomials_base.h.
|
privateinherited |
The number of polynomials represented by this object.
Definition at line 150 of file tensor_polynomials_base.h.