Reference documentation for deal.II version 8.1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Member Functions | List of all members
FE_Nothing< dim > Class Template Reference

#include <fe_nothing.h>

Inheritance diagram for FE_Nothing< dim >:
[legend]

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
 
- Public Member Functions inherited from FiniteElement< dim >
 FiniteElement (const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
 
virtual ~FiniteElement ()
 
const FiniteElement< dim,
spacedim > & 
operator[] (const unsigned int fe_index) const
 
bool operator== (const FiniteElement< dim, spacedim > &) const
 
virtual std::size_t memory_consumption () const
 
 DeclException1 (ExcShapeFunctionNotPrimitive, int,<< "The shape function with index "<< arg1<< " is not primitive, i.e. it is vector-valued and "<< "has more than one non-zero vector component. This "<< "function cannot be called for these shape functions. "<< "Maybe you want to use the same function with the "<< "_component suffix?")
 
 DeclException0 (ExcFENotPrimitive)
 
 DeclException0 (ExcUnitShapeValuesDoNotExist)
 
 DeclException0 (ExcFEHasNoSupportPoints)
 
 DeclException0 (ExcEmbeddingVoid)
 
 DeclException0 (ExcProjectionVoid)
 
 DeclException0 (ExcConstraintsVoid)
 
 DeclException0 (ExcInterpolationNotImplemented)
 
 DeclException0 (ExcBoundaryFaceUsed)
 
 DeclException0 (ExcJacobiDeterminantHasWrongSign)
 
 DeclException2 (ExcWrongInterfaceMatrixSize, int, int,<< "The interface matrix has a size of "<< arg1<< "x"<< arg2<< ", which is not reasonable in the present dimension.")
 
 DeclException2 (ExcComponentIndexInvalid, int, int,<< "The component-index pair ("<< arg1<< ", "<< arg2<< ") is invalid, i.e. non-existent.")
 
virtual double shape_value_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const
 
virtual Tensor< 1, dim > shape_grad (const unsigned int i, const Point< dim > &p) const
 
virtual Tensor< 1, dim > shape_grad_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const
 
virtual Tensor< 2, dim > shape_grad_grad (const unsigned int i, const Point< dim > &p) const
 
virtual Tensor< 2, dim > shape_grad_grad_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const
 
virtual bool has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const
 
virtual const FullMatrix
< double > & 
get_restriction_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
 
virtual const FullMatrix
< double > & 
get_prolongation_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
 
bool prolongation_is_implemented () const
 
bool isotropic_prolongation_is_implemented () const
 
bool restriction_is_implemented () const
 
bool isotropic_restriction_is_implemented () const
 
bool restriction_is_additive (const unsigned int index) const
 
const FullMatrix< double > & constraints (const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
 
bool constraints_are_implemented (const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
 
virtual void get_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
 
virtual void get_face_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
 
virtual void get_subface_interpolation_matrix (const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
 
virtual std::vector< std::pair
< unsigned int, unsigned int > > 
hp_vertex_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const
 
virtual std::vector< std::pair
< unsigned int, unsigned int > > 
hp_line_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const
 
virtual std::vector< std::pair
< unsigned int, unsigned int > > 
hp_quad_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const
 
virtual
FiniteElementDomination::Domination 
compare_for_face_domination (const FiniteElement< dim, spacedim > &fe_other) const
 
std::pair< unsigned int,
unsigned int
system_to_component_index (const unsigned int index) const
 
unsigned int component_to_system_index (const unsigned int component, const unsigned int index) const
 
std::pair< unsigned int,
unsigned int
face_system_to_component_index (const unsigned int index) const
 
virtual unsigned int face_to_cell_index (const unsigned int face_dof_index, const unsigned int face, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false) const
 
unsigned int adjust_quad_dof_index_for_face_orientation (const unsigned int index, const bool face_orientation, const bool face_flip, const bool face_rotation) const
 
unsigned int adjust_line_dof_index_for_line_orientation (const unsigned int index, const bool line_orientation) const
 
const ComponentMaskget_nonzero_components (const unsigned int i) const
 
unsigned int n_nonzero_components (const unsigned int i) const
 
bool is_primitive (const unsigned int i) const
 
unsigned int n_base_elements () const
 
virtual const FiniteElement
< dim, spacedim > & 
base_element (const unsigned int index) const
 
unsigned int element_multiplicity (const unsigned int index) const
 
std::pair< std::pair< unsigned
int, unsigned int >, unsigned
int
system_to_base_index (const unsigned int index) const
 
std::pair< std::pair< unsigned
int, unsigned int >, unsigned
int
face_system_to_base_index (const unsigned int index) const
 
types::global_dof_index first_block_of_base (const unsigned int b) const
 
std::pair< unsigned int,
unsigned int
component_to_base_index (const unsigned int component) const
 
std::pair< unsigned int,
unsigned int
block_to_base_index (const unsigned int block) const
 
std::pair< unsigned int,
types::global_dof_index
system_to_block_index (const unsigned int component) const
 
unsigned int component_to_block_index (const unsigned int component) const
 
ComponentMask component_mask (const FEValuesExtractors::Scalar &scalar) const
 
ComponentMask component_mask (const FEValuesExtractors::Vector &vector) const
 
ComponentMask component_mask (const FEValuesExtractors::SymmetricTensor< 2 > &sym_tensor) const
 
ComponentMask component_mask (const BlockMask &block_mask) const
 
BlockMask block_mask (const FEValuesExtractors::Scalar &scalar) const
 
BlockMask block_mask (const FEValuesExtractors::Vector &vector) const
 
BlockMask block_mask (const FEValuesExtractors::SymmetricTensor< 2 > &sym_tensor) const
 
BlockMask block_mask (const ComponentMask &component_mask) const
 
const std::vector< Point< dim > > & get_unit_support_points () const
 
bool has_support_points () const
 
virtual Point< dim > unit_support_point (const unsigned int index) const
 
const std::vector< Point< dim-1 > > & get_unit_face_support_points () const
 
bool has_face_support_points () const
 
virtual Point< dim-1 > unit_face_support_point (const unsigned int index) const
 
const std::vector< Point< dim > > & get_generalized_support_points () const
 
bool has_generalized_support_points () const
 
const std::vector< Point< dim-1 > > & get_generalized_face_support_points () const
 
bool has_generalized_face_support_points () const
 
virtual void interpolate (std::vector< double > &local_dofs, const std::vector< double > &values) const
 
virtual void interpolate (std::vector< double > &local_dofs, const std::vector< Vector< double > > &values, unsigned int offset=0) const
 
virtual void interpolate (std::vector< double > &local_dofs, const VectorSlice< const std::vector< std::vector< double > > > &values) const
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
void subscribe (const char *identifier=0) const
 
void unsubscribe (const char *identifier=0) const
 
unsigned int n_subscriptions () const
 
void list_subscribers () const
 
 DeclException3 (ExcInUse, int, char *, std::string &,<< "Object of class "<< arg2<< " is still used by "<< arg1<< " other objects.\n"<< "(Additional information: "<< arg3<< ")\n"<< "Note the entry in the Frequently Asked Questions of "<< "deal.II (linked to from http://www.dealii.org/) for "<< "more information on what this error means.")
 
 DeclException2 (ExcNoSubscriber, char *, char *,<< "No subscriber with identifier \""<< arg2<< "\" did subscribe to this object of class "<< arg1)
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
- Public Member Functions inherited from FiniteElementData< dim >
 FiniteElementData ()
 
 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)
 
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 BlockIndicesblock_indices () const
 
bool is_primitive () const
 
unsigned int tensor_degree () const
 
bool conforms (const Conformity) const
 
bool operator== (const FiniteElementData &) const
 

Additional Inherited Members

- Public Types inherited from FiniteElementData< dim >
enum  Conformity {
  unknown = 0x00, L2 = 0x01, Hcurl = 0x02, Hdiv = 0x04,
  H1 = Hcurl | Hdiv, H2 = 0x0e
}
 
- Public Attributes inherited from FiniteElementData< dim >
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
 
BlockIndices block_indices_data
 
- Static Public Attributes inherited from FiniteElementData< dim >
static const unsigned int dimension = dim
 
- Protected Member Functions inherited from FiniteElement< dim >
void reinit_restriction_and_prolongation_matrices (const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
 
TableIndices< 2 > interface_constraints_size () const
 
void compute_2nd (const Mapping< dim, spacedim > &mapping, const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int offset, typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, InternalDataBase &fe_internal, FEValuesData< dim, spacedim > &data) const
 
- Protected Member Functions inherited from FiniteElementData< dim >
void set_primitivity (const bool value)
 
- Static Protected Member Functions inherited from FiniteElement< dim >
static std::vector< unsigned intcompute_n_nonzero_components (const std::vector< ComponentMask > &nonzero_components)
 
- Protected Attributes inherited from FiniteElement< dim >
std::vector< std::vector
< FullMatrix< double > > > 
restriction
 
std::vector< std::vector
< FullMatrix< double > > > 
prolongation
 
FullMatrix< doubleinterface_constraints
 
std::vector< Point< dim > > unit_support_points
 
std::vector< Point< dim-1 > > unit_face_support_points
 
std::vector< Point< dim > > generalized_support_points
 
std::vector< Point< dim-1 > > generalized_face_support_points
 
Table< 2, intadjust_quad_dof_index_for_face_orientation_table
 
std::vector< intadjust_line_dof_index_for_line_orientation_table
 

Detailed Description

template<int dim>
class FE_Nothing< dim >

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. An interesting application for this element is also presented in the paper A. Cangiani, J. Chapman, E. Georgoulis, M. Jensen: Implementation of the Continuous-Discontinuous Galerkin Finite Element Method, arXiv:1201.2878v1 [math.NA], 2012 (see http://arxiv.org/abs/1201.2878).

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.

Author
Joshua White, Wolfgang Bangerth

Definition at line 81 of file fe_nothing.h.

Constructor & Destructor Documentation

template<int dim>
FE_Nothing< dim >::FE_Nothing ( unsigned int  n_components = 1)

Constructor. Argument denotes the number of components to give this finite element (default = 1).

Member Function Documentation

template<int dim>
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 >.

template<int 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 >.

template<int 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 >.

template<int 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 >.

template<int 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 >.

template<int 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

Fill the fields of FEValues. This function performs all the operations needed to compute the data of an FEValues object.

In the current case, this function returns no meaningful information, since the element has no degrees of freedom.

template<int dim>
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.

template<int dim>
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.

template<int dim>
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.

template<int dim>
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.

template<int dim>
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 >.

template<int 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.

template<int dim>
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.


The documentation for this class was generated from the following file: