Reference documentation for deal.II version 9.0.0
|
#include <deal.II/fe/fe_base.h>
Public Types | |
enum | Conformity { unknown = 0x00, L2 = 0x01, Hcurl = 0x02, Hdiv = 0x04, H1 = Hcurl | Hdiv, H2 = 0x0e } |
Public Member Functions | |
FiniteElementData (const std::vector< unsigned int > &dofs_per_object, const unsigned int n_components, const unsigned int degree, const Conformity conformity=unknown, const BlockIndices &block_indices=BlockIndices()) | |
unsigned int | n_dofs_per_vertex () const |
unsigned int | n_dofs_per_line () const |
unsigned int | n_dofs_per_quad () const |
unsigned int | n_dofs_per_hex () const |
unsigned int | n_dofs_per_face () const |
unsigned int | n_dofs_per_cell () const |
template<int structdim> | |
unsigned int | n_dofs_per_object () const |
unsigned int | n_components () const |
unsigned int | n_blocks () const |
const BlockIndices & | block_indices () const |
unsigned int | tensor_degree () const |
bool | conforms (const Conformity) const |
bool | operator== (const FiniteElementData &) const |
Public Attributes | |
const unsigned int | dofs_per_vertex |
const unsigned int | dofs_per_line |
const unsigned int | dofs_per_quad |
const unsigned int | dofs_per_hex |
const unsigned int | first_line_index |
const unsigned int | first_quad_index |
const unsigned int | first_hex_index |
const unsigned int | first_face_line_index |
const unsigned int | first_face_quad_index |
const unsigned int | dofs_per_face |
const unsigned int | dofs_per_cell |
const unsigned int | components |
const unsigned int | degree |
const Conformity | conforming_space |
const BlockIndices | block_indices_data |
Static Public Attributes | |
static const unsigned int | dimension = dim |
A class that declares a number of scalar constant variables that describe basic properties of a finite element implementation. This includes, for example, the number of degrees of freedom per vertex, line, or cell; the number of vector components; etc.
The kind of information stored here is computed during initialization of a finite element object and is passed down to this class via its constructor. The data stored by this class is part of the public interface of the FiniteElement class (which derives from the current class). See there for more information.
enum FiniteElementData::Conformity |
Enumerator for the different types of continuity a finite element may have. Continuity is measured by the Sobolev space containing the constructed finite element space and is also called this way.
Note that certain continuities may imply others. For instance, a function in H1 is in Hcurl and Hdiv as well.
If you are interested in continuity in the classical sense, then the following relations hold:
H1 implies that the function is continuous over cell boundaries.
H2 implies that the function is continuously differentiable over cell boundaries.
In order to test if a finite element conforms to a certain space, use FiniteElementData<dim>::conforms().
FiniteElementData< dim >::FiniteElementData | ( | const std::vector< unsigned int > & | dofs_per_object, |
const unsigned int | n_components, | ||
const unsigned int | degree, | ||
const Conformity | conformity = unknown , |
||
const BlockIndices & | block_indices = BlockIndices() |
||
) |
Constructor, computing all necessary values from the distribution of dofs to geometrical objects.
[in] | dofs_per_object | A vector that describes the number of degrees of freedom on geometrical objects for each dimension. This vector must have size dim+1, and entry 0 describes the number of degrees of freedom per vertex, entry 1 the number of degrees of freedom per line, etc. As an example, for the common \(Q_1\) Lagrange element in 2d, this vector would have elements (1,0,0) . On the other hand, for a \(Q_3\) element in 3d, it would have entries (1,2,4,8) . |
[in] | n_components | Number of vector components of the element. |
[in] | degree | The maximal polynomial degree of any of the shape functions of this element in any variable on the reference element. For example, for the \(Q_1\) element (in any space dimension), this would be one; this is so despite the fact that the element has a shape function of the form \(\hat x\hat y\) (in 2d) and \(\hat x\hat y\hat z\) (in 3d), which, although quadratic and cubic polynomials, are still only linear in each reference variable separately. The information provided by this variable is typically used in determining what an appropriate quadrature formula is. |
[in] | conformity | A variable describing which Sobolev space this element conforms to. For example, the \(Q_p\) Lagrange elements (implemented by the FE_Q class) are \(H^1\) conforming, whereas the Raviart-Thomas element (implemented by the FE_RaviartThomas class) is \(H_\text{div}\) conforming; finally, completely discontinuous elements (implemented by the FE_DGQ class) are only \(L_2\) conforming. |
[in] | block_indices | An argument that describes how the base elements of a finite element are grouped. The default value constructs a single block that consists of all dofs_per_cell degrees of freedom. This is appropriate for all "atomic" elements (including non-primitive ones) and these can therefore omit this argument. On the other hand, composed elements such as FESystem will want to pass a different value here. |
Definition at line 23 of file fe_data.cc.
unsigned int FiniteElementData< dim >::n_dofs_per_vertex | ( | ) | const |
Number of dofs per vertex.
unsigned int FiniteElementData< dim >::n_dofs_per_line | ( | ) | const |
Number of dofs per line. Not including dofs on lower dimensional objects.
unsigned int FiniteElementData< dim >::n_dofs_per_quad | ( | ) | const |
Number of dofs per quad. Not including dofs on lower dimensional objects.
unsigned int FiniteElementData< dim >::n_dofs_per_hex | ( | ) | const |
Number of dofs per hex. Not including dofs on lower dimensional objects.
unsigned int FiniteElementData< dim >::n_dofs_per_face | ( | ) | const |
Number of dofs per face, accumulating degrees of freedom of all lower dimensional objects.
unsigned int FiniteElementData< dim >::n_dofs_per_cell | ( | ) | const |
Number of dofs per cell, accumulating degrees of freedom of all lower dimensional objects.
unsigned int FiniteElementData< dim >::n_dofs_per_object | ( | ) | const |
Return the number of degrees per structdim-dimensional object. For structdim==0, the function therefore returns dofs_per_vertex, for structdim==1 dofs_per_line, etc. This function is mostly used to allow some template trickery for functions that should work on all sorts of objects without wanting to use the different names (vertex, line, ...) associated with these objects.
unsigned int FiniteElementData< dim >::n_components | ( | ) | const |
Number of components. See the glossary for more information.
unsigned int FiniteElementData< dim >::n_blocks | ( | ) | const |
Number of blocks. See the glossary for more information.
const BlockIndices& FiniteElementData< dim >::block_indices | ( | ) | const |
Detailed information on block sizes.
unsigned int FiniteElementData< dim >::tensor_degree | ( | ) | const |
Maximal polynomial degree of a shape function in a single coordinate direction.
This function can be used to determine the optimal quadrature rule.
bool FiniteElementData< dim >::conforms | ( | const Conformity | ) | const |
Test whether a finite element space conforms to a certain Sobolev space.
bool FiniteElementData< dim >::operator== | ( | const FiniteElementData< dim > & | f | ) | const |
Comparison operator.
Definition at line 72 of file fe_data.cc.
|
static |
const unsigned int FiniteElementData< dim >::dofs_per_vertex |
const unsigned int FiniteElementData< dim >::dofs_per_line |
const unsigned int FiniteElementData< dim >::dofs_per_quad |
const unsigned int FiniteElementData< dim >::dofs_per_hex |
const unsigned int FiniteElementData< dim >::first_line_index |
const unsigned int FiniteElementData< dim >::first_quad_index |
const unsigned int FiniteElementData< dim >::first_hex_index |
const unsigned int FiniteElementData< dim >::first_face_line_index |
const unsigned int FiniteElementData< dim >::first_face_quad_index |
const unsigned int FiniteElementData< dim >::dofs_per_face |
const unsigned int FiniteElementData< dim >::dofs_per_cell |
const unsigned int FiniteElementData< dim >::components |
Number of vector components of this finite element, and dimension of the image space. For vector-valued finite elements (i.e. when this number is greater than one), the number of vector components is in many cases equal to the number of base elements glued together with the help of the FESystem class. However, for elements like the Nedelec element, the number is greater than one even though we only have one base element.
const unsigned int FiniteElementData< dim >::degree |
const Conformity FiniteElementData< dim >::conforming_space |
const BlockIndices FiniteElementData< dim >::block_indices_data |
Storage for an object describing the sizes of each block of a compound element. For an element which is not an FESystem, this contains only a single block with length dofs_per_cell.