
Public Member Functions | |
| FE_Nothing (unsigned int n_components=1) | |
| virtual FiniteElement< dim > * | clone () const |
| virtual std::string | get_name () const |
| virtual UpdateFlags | update_once (const UpdateFlags flags) const |
| virtual UpdateFlags | update_each (const UpdateFlags flags) const |
| virtual double | shape_value (const unsigned int i, const Point< dim > &p) const |
| virtual void | fill_fe_values (const Mapping< dim > &mapping, const typename Triangulation< dim >::cell_iterator &cell, const Quadrature< dim > &quadrature, typename Mapping< dim >::InternalDataBase &mapping_data, typename Mapping< dim >::InternalDataBase &fedata, FEValuesData< dim, dim > &data, CellSimilarity::Similarity &cell_similarity) const |
| virtual void | fill_fe_face_values (const Mapping< dim > &mapping, const typename Triangulation< dim >::cell_iterator &cell, const unsigned int face, const Quadrature< dim-1 > &quadrature, typename Mapping< dim >::InternalDataBase &mapping_data, typename Mapping< dim >::InternalDataBase &fedata, FEValuesData< dim, dim > &data) const |
| virtual void | fill_fe_subface_values (const Mapping< dim > &mapping, const typename Triangulation< dim >::cell_iterator &cell, const unsigned int face, const unsigned int subface, const Quadrature< dim-1 > &quadrature, typename Mapping< dim >::InternalDataBase &mapping_data, typename Mapping< dim >::InternalDataBase &fedata, FEValuesData< dim, dim > &data) const |
| virtual Mapping< dim > ::InternalDataBase * | get_data (const UpdateFlags update_flags, const Mapping< dim > &mapping, const Quadrature< dim > &quadrature) const |
| virtual FiniteElementDomination::Domination | compare_for_face_domination (const FiniteElement< dim > &fe_other) const |
| virtual std::vector< std::pair < unsigned int, unsigned int > > | hp_vertex_dof_identities (const FiniteElement< dim > &fe_other) const |
| virtual std::vector< std::pair < unsigned int, unsigned int > > | hp_line_dof_identities (const FiniteElement< dim > &fe_other) const |
| virtual std::vector< std::pair < unsigned int, unsigned int > > | hp_quad_dof_identities (const FiniteElement< dim > &fe_other) const |
| virtual bool | hp_constraints_are_implemented () const |
| virtual void | get_face_interpolation_matrix (const FiniteElement< dim > &source_fe, FullMatrix< double > &interpolation_matrix) const |
| virtual void | get_subface_interpolation_matrix (const FiniteElement< dim > &source_fe, const unsigned int index, FullMatrix< double > &interpolation_matrix) const |
Definition of a finite element with zero degrees of freedom. This class is useful (in the context of an hp method) to represent empty cells in the triangulation on which no degrees of freedom should be allocated, or to describe a field that is extended by zero to a part of the domain where we don't need it. Thus a triangulation may be divided into two regions: an active region where normal elements are used, and an inactive region where FE_Nothing elements are used. The hp::DoFHandler will therefore assign no degrees of freedom to the FE_Nothing cells, and this subregion is therefore implicitly deleted from the computation. step-46 shows a use case for this element.
Note that some care must be taken that the resulting mesh topology continues to make sense when FE_Nothing elements are introduced. This is particularly true when dealing with hanging node constraints, because the library makes some basic assumptions about the nature of those constraints. The following geometries are acceptable:
+---------+----+----+ | | 0 | | | 1 +----+----+ | | 0 | | +---------+----+----+
+---------+----+----+ | | 1 | | | 0 +----+----+ | | 1 | | +---------+----+----+
Here, 0 denotes an FE_Nothing cell, and 1 denotes some other element type. The library has no difficulty computing the necessary hanging node constraints in these cases (i.e. no constraint). However, the following geometry is NOT acceptable (at least in the current implementation):
+---------+----+----+ | | 0 | | | 1 +----+----+ | | 1 | | +---------+----+----+
The distinction lies in the mixed nature of the child faces, a case we have not implemented as of yet.
Definition at line 73 of file fe_nothing.h.
| FE_Nothing< dim >::FE_Nothing | ( | unsigned int | n_components = 1 ) |
Constructor. Argument denotes the number of components to give this finite element (default = 1).
| virtual FiniteElement<dim>* FE_Nothing< dim >::clone | ( | ) | const [virtual] |
A sort of virtual copy constructor. Some places in the library, for example the constructors of FESystem as well as the hp::FECollection class, need to make copied of finite elements without knowing their exact type. They do so through this function.
Implements FiniteElement< dim >.
| virtual std::string FE_Nothing< dim >::get_name | ( | ) | const [virtual] |
Return a string that uniquely identifies a finite element. In this case it is FE_Nothing<dim>.
Implements FiniteElement< dim >.
| virtual UpdateFlags FE_Nothing< dim >::update_once | ( | const UpdateFlags | flags ) | const [virtual] |
Determine the values a finite element should compute on initialization of data for FEValues.
Given a set of flags indicating what quantities are requested from a FEValues object, update_once() and update_each() compute which values must really be computed. Then, the fill_*_values functions are called with the result of these.
In this case, since the element has zero degrees of freedom and no information can be computed on it, this function simply returns the default (empty) set of update flags.
Implements FiniteElement< dim >.
| virtual UpdateFlags FE_Nothing< dim >::update_each | ( | const UpdateFlags | flags ) | const [virtual] |
Complementary function for update_once().
While update_once() returns the values to be computed on the unit cell for yielding the required data, this function determines the values that must be recomputed on each cell.
Refer to update_once() for more details.
Implements FiniteElement< dim >.
| virtual double FE_Nothing< dim >::shape_value | ( | const unsigned int | i, |
| const Point< dim > & | p | ||
| ) | const [virtual] |
Return the value of the ith shape function at the point p. p is a point on the reference element. Because the current element has no degrees of freedom, this function should obviously not be called in practice. All this function really does, therefore, is trigger an exception.
Reimplemented from FiniteElement< dim >.
| virtual void FE_Nothing< dim >::fill_fe_values | ( | const Mapping< dim > & | mapping, |
| const typename Triangulation< dim >::cell_iterator & | cell, | ||
| const Quadrature< dim > & | quadrature, | ||
| typename Mapping< dim >::InternalDataBase & | mapping_data, | ||
| typename Mapping< dim >::InternalDataBase & | fedata, | ||
| FEValuesData< dim, dim > & | data, | ||
| CellSimilarity::Similarity & | cell_similarity | ||
| ) | const [virtual] |
| virtual void FE_Nothing< dim >::fill_fe_face_values | ( | const Mapping< dim > & | mapping, |
| const typename Triangulation< dim >::cell_iterator & | cell, | ||
| const unsigned int | face, | ||
| const Quadrature< dim-1 > & | quadrature, | ||
| typename Mapping< dim >::InternalDataBase & | mapping_data, | ||
| typename Mapping< dim >::InternalDataBase & | fedata, | ||
| FEValuesData< dim, dim > & | data | ||
| ) | const [virtual] |
Fill the fields of FEFaceValues. This function performs all the operations needed to compute the data of an FEFaceValues object.
In the current case, this function returns no meaningful information, since the element has no degrees of freedom.
| virtual void FE_Nothing< dim >::fill_fe_subface_values | ( | const Mapping< dim > & | mapping, |
| const typename Triangulation< dim >::cell_iterator & | cell, | ||
| const unsigned int | face, | ||
| const unsigned int | subface, | ||
| const Quadrature< dim-1 > & | quadrature, | ||
| typename Mapping< dim >::InternalDataBase & | mapping_data, | ||
| typename Mapping< dim >::InternalDataBase & | fedata, | ||
| FEValuesData< dim, dim > & | data | ||
| ) | const [virtual] |
Fill the fields of FESubFaceValues. This function performs all the operations needed to compute the data of an FESubFaceValues object.
In the current case, this function returns no meaningful information, since the element has no degrees of freedom.
| virtual Mapping<dim>::InternalDataBase* FE_Nothing< dim >::get_data | ( | const UpdateFlags | update_flags, |
| const Mapping< dim > & | mapping, | ||
| const Quadrature< dim > & | quadrature | ||
| ) | const [virtual] |
Prepare internal data structures and fill in values independent of the cell. Returns a pointer to an object of which the caller of this function then has to assume ownership (which includes destruction when it is no more needed).
In the current case, this function just returns a default pointer, since no meaningful data exists for this element.
| virtual FiniteElementDomination::Domination FE_Nothing< dim >::compare_for_face_domination | ( | const FiniteElement< dim > & | fe_other ) | const [virtual] |
Return whether this element dominates the one given as argument when they meet at a common face, whether it is the other way around, whether neither dominates, or if either could dominate.
For a definition of domination, see FiniteElementBase::Domination and in particular the hp paper.
In the current case, this element is always assumed to dominate, unless it is also of type FE_Nothing(). In that situation, either element can dominate.
| virtual std::vector<std::pair<unsigned int, unsigned int> > FE_Nothing< dim >::hp_vertex_dof_identities | ( | const FiniteElement< dim > & | fe_other ) | const [virtual] |
| virtual std::vector<std::pair<unsigned int, unsigned int> > FE_Nothing< dim >::hp_line_dof_identities | ( | const FiniteElement< dim > & | fe_other ) | const [virtual] |
| virtual std::vector<std::pair<unsigned int, unsigned int> > FE_Nothing< dim >::hp_quad_dof_identities | ( | const FiniteElement< dim > & | fe_other ) | const [virtual] |
| virtual bool FE_Nothing< dim >::hp_constraints_are_implemented | ( | ) | const [virtual] |
Return whether this element implements its hanging node constraints in the new way, which has to be used to make elements "hp compatible". That means, the element properly implements the get_face_interpolation_matrix and get_subface_interpolation_matrix methods. Therefore the return value also indicates whether a call to the get_face_interpolation_matrix() method and the get_subface_interpolation_matrix() method will generate an error or not.
Currently the main purpose of this function is to allow the make_hanging_node_constraints method to decide whether the new procedures, which are supposed to work in the hp framework can be used, or if the old well verified but not hp capable functions should be used. Once the transition to the new scheme for computing the interface constraints is complete, this function will be superfluous and will probably go away.
Derived classes should implement this function accordingly. The default assumption is that a finite element does not provide hp capable face interpolation, and the default implementation therefore returns false.
Reimplemented from FiniteElement< dim >.
| virtual void FE_Nothing< dim >::get_face_interpolation_matrix | ( | const FiniteElement< dim > & | source_fe, |
| FullMatrix< double > & | interpolation_matrix | ||
| ) | const [virtual] |
Return the matrix interpolating from a face of of one element to the face of the neighboring element. The size of the matrix is then source.dofs_per_face times this->dofs_per_face.
Since the current finite element has no degrees of freedom, the interpolation matrix is necessarily empty.
| virtual void FE_Nothing< dim >::get_subface_interpolation_matrix | ( | const FiniteElement< dim > & | source_fe, |
| const unsigned int | index, | ||
| FullMatrix< double > & | interpolation_matrix | ||
| ) | const [virtual] |
Return the matrix interpolating from a face of of one element to the subface of the neighboring element. The size of the matrix is then source.dofs_per_face times this->dofs_per_face.
Since the current finite element has no degrees of freedom, the interpolation matrix is necessarily empty.
documentation generated on Fri Feb 3 2012 06:04:07 by
doxygen
1.7.2