Reference documentation for deal.II version 9.2.0
|
#include <deal.II/base/polynomials_bernardi_raugel.h>
Public Member Functions | |
PolynomialsBernardiRaugel (const unsigned int k) | |
std::string | name () const override |
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 |
virtual std::unique_ptr< TensorPolynomialsBase< dim > > | clone () const override |
Public Member Functions inherited from TensorPolynomialsBase< dim > | |
TensorPolynomialsBase (const unsigned int deg, const unsigned int n_polynomials) | |
TensorPolynomialsBase (TensorPolynomialsBase< dim > &&)=default | |
TensorPolynomialsBase (const TensorPolynomialsBase< dim > &)=default | |
virtual | ~TensorPolynomialsBase ()=default |
unsigned int | n () const |
unsigned int | degree () const |
Static Public Member Functions | |
static unsigned int | n_polynomials (const unsigned int k) |
Static Private Member Functions | |
static std::vector< std::vector< Polynomials::Polynomial< double > > > | create_polynomials_Q () |
static std::vector< std::vector< Polynomials::Polynomial< double > > > | create_polynomials_bubble () |
Private Attributes | |
const AnisotropicPolynomials< dim > | polynomial_space_Q |
const AnisotropicPolynomials< dim > | polynomial_space_bubble |
This class implements the Bernardi-Raugel polynomials similarly to the description in the Mathematics of Computation paper from 1985 by Christine Bernardi and Geneviève Raugel.
The Bernardi-Raugel polynomials are originally defined as an enrichment of the \((P_1)^d\) elements on simplicial meshes for Stokes problems by the addition of bubble functions, yielding a locking-free finite element which is a subset of \((P_2)^d\) elements. This implementation is an enrichment of \((Q_1)^d\) elements which is a subset of \((Q_2)^d\) elements for quadrilateral and hexahedral meshes.
The \(BR_1\) bubble functions are defined to have magnitude 1 at the center of face \(e_i\) and direction \(\mathbf{n}_i\) normal to face \(e_i\), and magnitude 0 on all other vertices and faces. Ordering is consistent with the face numbering in GeometryInfo. The vector \(\mathbf{n}_i\) points in the positive axis direction and not necessarily normal to the element for consistent orientation across edges.
\(x=0\) edge: \(\mathbf{p}_1 = \mathbf{n}_1 (1-x)(y)(1-y)\)
@f$x=1@f$ edge: @f$\mathbf{p}_2 = \mathbf{n}_2 (x)(y)(1-y)@f$ @f$y=0@f$ edge: @f$\mathbf{p}_3 = \mathbf{n}_3 (x)(1-x)(1-y)@f$ @f$y=1@f$ edge: @f$\mathbf{p}_4 = \mathbf{n}_4 (x)(1-x)(y)@f$
\(x=0\) edge: \(\mathbf{p}_1 = \mathbf{n}_1 (1-x)(y)(1-y)(z)(1-z)\)
@f$x=1@f$ edge: @f$\mathbf{p}_2 = \mathbf{n}_2 (x)(y)(1-y)(z)(1-z)@f$ @f$y=0@f$ edge: @f$\mathbf{p}_3 = \mathbf{n}_3 (x)(1-x)(1-y)(z)(1-z)@f$ @f$y=1@f$ edge: @f$\mathbf{p}_4 = \mathbf{n}_4 (x)(1-x)(y)(z)(1-z)@f$ @f$z=0@f$ edge: @f$\mathbf{p}_5 = \mathbf{n}_5 (x)(1-x)(y)(1-y)(1-z)@f$ @f$z=1@f$ edge: @f$\mathbf{p}_6 = \mathbf{n}_6 (x)(1-x)(y)(1-y)(z)@f$
Then the \(BR_1(E)\) polynomials are defined on quadrilaterals and hexahedra by \(BR_1(E) = Q_1(E) \oplus \mbox{span}\{\mathbf{p}_i, i=1,...,2d\}\).
Definition at line 87 of file polynomials_bernardi_raugel.h.
PolynomialsBernardiRaugel< dim >::PolynomialsBernardiRaugel | ( | const unsigned int | k | ) |
Constructor. Creates all basis functions for Bernardi-Raugel polynomials of given degree.
k=1
. Definition at line 24 of file polynomials_bernardi_raugel.cc.
|
inlineoverridevirtual |
Return the name of the space, which is BernardiRaugel
.
Implements TensorPolynomialsBase< dim >.
Definition at line 173 of file polynomials_bernardi_raugel.h.
|
overridevirtual |
Compute the value and derivatives of each Bernardi-Raugel 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 70 of file polynomials_bernardi_raugel.cc.
|
static |
Return the number of polynomials in the space BR(degree)
without requiring to build an object of PolynomialsBernardiRaugel. This is required by the FiniteElement classes.
Definition at line 239 of file polynomials_bernardi_raugel.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 255 of file polynomials_bernardi_raugel.cc.
|
staticprivate |
A static member function that creates the polynomial space we use to initialize the polynomial_space_Q member variable.
Definition at line 55 of file polynomials_bernardi_raugel.cc.
|
staticprivate |
A static member function that creates the polynomial space we use to initialize the polynomial_space_bubble member variable.
Definition at line 33 of file polynomials_bernardi_raugel.cc.
|
private |
An object representing the polynomial space of Q functions which forms the BR
polynomials through outer products of these with the corresponding unit ijk vectors.
Definition at line 146 of file polynomials_bernardi_raugel.h.
|
private |
An object representing the polynomial space of bubble functions which forms the BR
polynomials through outer products of these with the corresponding normals.
Definition at line 153 of file polynomials_bernardi_raugel.h.