
Dimension independent data for finite elements. See the derived class FiniteElement class for information on its use. All its data are available to the implementation in a concrete finite element class.
Definition at line 181 of file fe_base.h.
| 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 | ( | ) |
Default constructor. Constructs an element with no dofs. Checking n_dofs_per_cell() is therefore a good way to check if something went wrong.
| 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 unsigned int | n_blocks = numbers::invalid_unsigned_int |
||
| ) |
Constructor, computing all necessary values from the distribution of dofs to geometrcal objects.
| dofs_per_object | Number of dofs on geometrical objects for each dimension. In this vector, entry 0 refers to dofs on vertices, entry 1 on lines and so on. Its length must be dim+1. |
| n_components | Number of vector components of the element. |
| degree | Maximal polynomial degree in a single direction. |
| conformity | The finite element space has continuity of this Sobolev space. |
| n_blocks | obsolete and ignored. |
| 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.
Referenced by MeshWorker::IntegrationInfo< dim, spacedim >::initialize().
| unsigned int FiniteElementData< dim >::n_blocks | ( | ) | const |
Number of blocks.
| const BlockIndices& FiniteElementData< dim >::block_indices | ( | ) | const |
Detailed infromation on block sizes.
| bool FiniteElementData< dim >::is_primitive | ( | ) | const |
Return whether the entire finite element is primitive, in the sense that all its shape functions are primitive. If the finite element is scalar, then this is always the case.
Since this is an extremely common operation, the result is cached in the cached_primitivity variable which is computed in the constructor.
| 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.
Referenced by MeshWorker::IntegrationInfoBox< dim, spacedim >::initialize().
| bool FiniteElementData< dim >::conforms | ( | const | Conformity ) | const |
Test whether a finite element space conforms to a certain Sobolev space.
| unsigned int FiniteElementData< dim >::face_to_cell_index | ( | const unsigned int | face_index, |
| const unsigned int | face, | ||
| const bool | face_orientation = true, |
||
| const bool | face_flip = false, |
||
| const bool | face_rotation = false |
||
| ) | const |
Given an index in the natural ordering of indices on a face, return the index of the same degree of freedom on the cell.
Referenced by BlockList::initialize_vertex_patches_mg().
| unsigned int FiniteElementData< dim >::face_to_equivalent_cell_index | ( | const unsigned int | index ) | const |
Given an index in the natural ordering of indices on a face, return the index of an equivalent degree of freedom on the cell.
To explain the concept, consider the case where we would like to know whether a degree of freedom on a face is primitive. Unfortunately, the is_primitive() function in the FiniteElement class takes a cell index, so we would need to find the cell index of the shape function that corresponds to the present face index. This function does that.
Code implementing this would then look like this:
for (i=0; i<dofs_per_face; ++i) if (fe.is_primitive(fe.face_to_equivalent_cell_index(i))) ... do whatever
| bool FiniteElementData< dim >::operator== | ( | const FiniteElementData< dim > & | ) | const |
Comparison operator.
| void FiniteElementData< dim >::set_primitivity | ( | const bool | value ) | [protected] |
Set the primitivity of the element. This is usually done by the constructor of a derived class. See primitive for details.
const unsigned int FiniteElementData< dim >::dimension = dim [static] |
| const unsigned int FiniteElementData< dim >::dofs_per_vertex |
Number of degrees of freedom on a vertex.
Definition at line 302 of file fe_base.h.
Referenced by BlockList::initialize_vertex_patches_mg().
| 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 |
Number of degrees of freedom on a face. This is the accumulated number of degrees of freedom on all the objects of dimension up to dim-1 constituting a face.
Definition at line 363 of file fe_base.h.
Referenced by BlockList::initialize_vertex_patches_mg().
| 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 |
Maximal polynomial degree of a shape function in a single coordinate direction.
Reimplemented in FE_DGPNonparametric< dim, spacedim >.
| const Conformity FiniteElementData< dim >::conforming_space |
| 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.
bool FiniteElementData< dim >::cached_primitivity [private] |
documentation generated on Fri Feb 3 2012 06:04:08 by
doxygen
1.7.2