Reference documentation for deal.II version 9.1.1
|
#include <deal.II/base/polynomials_abf.h>
Public Member Functions | |
PolynomialsABF (const unsigned int k) | |
void | compute (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 |
unsigned int | n () const |
unsigned int | degree () const |
std::string | name () const |
Static Public Member Functions | |
static unsigned int | compute_n_pols (unsigned int degree) |
Private Attributes | |
const unsigned int | my_degree |
const AnisotropicPolynomials< dim > | polynomial_space |
unsigned int | n_pols |
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 |
This class implements the Hdiv-conforming, vector-valued Arnold-Boffi-Falk polynomials as described in the article by Arnold-Boffi- 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 54 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.
void PolynomialsABF< dim >::compute | ( | 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 |
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.
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.
Definition at line 64 of file polynomials_abf.cc.
|
inline |
Return the number of ABF polynomials.
Definition at line 166 of file polynomials_abf.h.
|
inline |
Return the degree of the ABF space, which is two less than the highest polynomial degree.
Definition at line 174 of file polynomials_abf.h.
|
inline |
Return the name of the space, which is ABF
.
Definition at line 182 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 157 of file polynomials_abf.cc.
|
private |
The degree of this object as given to the constructor.
Definition at line 118 of file polynomials_abf.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 125 of file polynomials_abf.h.
|
private |
Number of Raviart-Thomas polynomials.
Definition at line 130 of file polynomials_abf.h.
|
mutableprivate |
A mutex that guards the following scratch arrays.
Definition at line 135 of file polynomials_abf.h.
|
mutableprivate |
Auxiliary memory.
Definition at line 140 of file polynomials_abf.h.
|
mutableprivate |
Auxiliary memory.
Definition at line 145 of file polynomials_abf.h.
|
mutableprivate |
Auxiliary memory.
Definition at line 150 of file polynomials_abf.h.
|
mutableprivate |
Auxiliary memory.
Definition at line 155 of file polynomials_abf.h.
|
mutableprivate |
Auxiliary memory.
Definition at line 160 of file polynomials_abf.h.