Reference documentation for deal.II version 9.6.0
|
#include <deal.II/base/polynomials_abf.h>
Public Member Functions | |
PolynomialsABF (const unsigned int k) | |
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 |
unsigned int | n () const |
unsigned int | degree () const |
Static Public Member Functions | |
static unsigned int | n_polynomials (const unsigned int degree) |
Private Attributes | |
const AnisotropicPolynomials< dim > | polynomial_space |
Threads::Mutex | mutex |
std::vector< double > | p_values |
std::vector< Tensor< 1, dim > > | p_grads |
std::vector< Tensor< 2, dim > > | p_grad_grads |
std::vector< Tensor< 3, dim > > | p_third_derivatives |
std::vector< Tensor< 4, dim > > | p_fourth_derivatives |
const unsigned int | polynomial_degree |
const unsigned int | n_pols |
This class implements the Hdiv-conforming, vector-valued Arnold-Boffi-Falk polynomials as described in the article by Arnold, Boffi, and Falk: Quadrilateral H(div) finite elements, SIAM J. Numer. Anal., Vol. 42, No. 6, pp.2429-2451
The ABF polynomials are constructed such that the divergence is in the tensor product polynomial space Qk. Therefore, the polynomial order of each component must be two orders higher in the corresponding direction, yielding the polynomial spaces (Qk+2,k, Qk,k+2) and (Qk+2,k,k, Qk,k+2,k, Qk,k,k+2) in 2d and 3d, resp.
Definition at line 52 of file polynomials_abf.h.
PolynomialsABF< dim >::PolynomialsABF | ( | const unsigned int | k | ) |
Constructor. Creates all basis functions for Raviart-Thomas polynomials of given degree.
Definition at line 49 of file polynomials_abf.cc.
|
overridevirtual |
Compute the value and the first and second 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.
Implements TensorPolynomialsBase< dim >.
Definition at line 63 of file polynomials_abf.cc.
|
inlineoverridevirtual |
Return the name of the space, which is ABF
.
Implements TensorPolynomialsBase< dim >.
Definition at line 142 of file polynomials_abf.h.
|
static |
Return the number of polynomials in the space RT(degree)
without requiring to build an object of PolynomialsABF. This is required by the FiniteElement classes.
Definition at line 156 of file polynomials_abf.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 184 of file polynomials_abf.cc.
|
inlineinherited |
Return the number of polynomials.
Definition at line 151 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 160 of file tensor_polynomials_base.h.
|
private |
An object representing the polynomial space for a single component. We can re-use it for the other vector components by rotating the coordinates of the evaluation point.
Definition at line 106 of file polynomials_abf.h.
|
mutableprivate |
A mutex that guards the following scratch arrays.
Definition at line 111 of file polynomials_abf.h.
|
mutableprivate |
Auxiliary memory.
Definition at line 116 of file polynomials_abf.h.
|
mutableprivate |
Auxiliary memory.
Definition at line 121 of file polynomials_abf.h.
|
mutableprivate |
Auxiliary memory.
Definition at line 126 of file polynomials_abf.h.
|
mutableprivate |
Auxiliary memory.
Definition at line 131 of file polynomials_abf.h.
|
mutableprivate |
Auxiliary memory.
Definition at line 136 of file polynomials_abf.h.
|
privateinherited |
The highest polynomial degree of this functions represented by this object.
Definition at line 139 of file tensor_polynomials_base.h.
|
privateinherited |
The number of polynomials represented by this object.
Definition at line 144 of file tensor_polynomials_base.h.