![]() |
Reference documentation for deal.II version 8.4.2
|
#include <deal.II/fe/fe_raviart_thomas.h>
Public Member Functions | |
| FE_RaviartThomasNodal (const unsigned int p) | |
| virtual std::string | get_name () const |
| virtual FiniteElement< dim > * | clone () 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 |
| virtual bool | hp_constraints_are_implemented () const |
Public Member Functions inherited from FE_PolyTensor< PolynomialsRaviartThomas< dim >, dim > | |
| FE_PolyTensor (const unsigned int degree, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components) | |
| virtual double | shape_value (const unsigned int i, const Point< dim > &p) const |
| virtual Tensor< 1, dim > | shape_grad (const unsigned int i, const Point< dim > &p) const |
| virtual Tensor< 2, dim > | shape_grad_grad (const unsigned int i, const Point< dim > &p) const |
Public Member Functions inherited from FiniteElement< dim, 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 (ExcInterpolationNotImplemented) | |
| DeclExceptionMsg (ExcUnitShapeValuesDoNotExist, "You are trying to access the values or derivatives of shape functions " "on the reference cell of an element that does not define its shape " "functions through mapping from the reference cell. Consequently, " "you cannot ask for shape function values or derivatives on the " "reference cell.") | |
| DeclExceptionMsg (ExcFEHasNoSupportPoints, "You are trying to access the support points of a finite " "element that either has no support points at all, or for " "which the corresponding tables have not been implemented.") | |
| DeclExceptionMsg (ExcEmbeddingVoid, "You are trying to access the matrices that describe how " "to embed a finite element function on one cell into the " "finite element space on one of its children (i.e., the " "'embedding' or 'prolongation' matrices). However, the " "current finite element can either not define this sort of " "operation, or it has not yet been implemented.") | |
| DeclExceptionMsg (ExcProjectionVoid, "You are trying to access the matrices that describe how " "to restrict a finite element function from the children " "of one cell to the finite element space defined on their " "parent (i.e., the 'restriction' or 'projection' matrices). " "However, the current finite element can either not define " "this sort of operation, or it has not yet been " "implemented.") | |
| DeclException2 (ExcWrongInterfaceMatrixSize, int, int,<< "The interface matrix has a size of "<< arg1<< "x"<< arg2<< ", which is not reasonable for the current element " "in the present dimension.") | |
| virtual Tensor< 3, dim > | shape_3rd_derivative (const unsigned int i, const Point< dim > &p) const |
| virtual Tensor< 3, dim > | shape_3rd_derivative_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const |
| virtual Tensor< 4, dim > | shape_4th_derivative (const unsigned int i, const Point< dim > &p) const |
| virtual Tensor< 4, dim > | shape_4th_derivative_component (const unsigned int i, const Point< dim > &p, const unsigned int component) 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 |
| 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 |
| 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_line_dof_index_for_line_orientation (const unsigned int index, const bool line_orientation) const |
| const ComponentMask & | get_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 |
| virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > | get_constant_modes () 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 |
| GeometryPrimitive | get_associated_geometry_primitive (const unsigned int cell_dof_index) const |
Public Member Functions inherited from Subscriptor | |
| Subscriptor () | |
| Subscriptor (const Subscriptor &) | |
| virtual | ~Subscriptor () |
| Subscriptor & | operator= (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."<< "\"<< "(Additional information: "<< arg3<< ")\"<< "See the entry in the Frequently Asked Questions of "<< "deal.II (linked to from http://www.dealii.org/) for "<< "a lot more information on what this error means and "<< "how to fix programs in which it happens.") | |
| DeclException2 (ExcNoSubscriber, char *, char *,<< "No subscriber with identifier <"<< arg2<< "> subscribes to this object of class "<< arg1<< ". Consequently, it cannot be unsubscribed.") | |
| template<class Archive > | |
| void | serialize (Archive &ar, const unsigned int version) |
Public Member Functions inherited from 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()) | |
| 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 |
| bool | is_primitive () const |
| unsigned int | tensor_degree () const |
| bool | conforms (const Conformity) const |
| bool | operator== (const FiniteElementData &) const |
Private Member Functions | |
| virtual bool | has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const |
| void | initialize_support_points (const unsigned int rt_degree) |
Static Private Member Functions | |
| static std::vector< unsigned int > | get_dpo_vector (const unsigned int degree) |
| static std::vector< bool > | get_ria_vector (const unsigned int degree) |
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 |
| const BlockIndices | block_indices_data |
Static Public Attributes inherited from FiniteElement< dim, dim > | |
| static const unsigned int | space_dimension |
Static Public Attributes inherited from FiniteElementData< dim > | |
| static const unsigned int | dimension = dim |
Protected Member Functions inherited from FiniteElement< dim, 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 |
| virtual InternalDataBase * | get_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const=0 |
| virtual InternalDataBase * | get_face_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const |
| virtual InternalDataBase * | get_subface_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const |
| virtual void | fill_fe_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValues::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const=0 |
| virtual void | fill_fe_face_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim-1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValues::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const=0 |
| virtual void | fill_fe_subface_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim-1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValues::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const=0 |
Protected Member Functions inherited from FiniteElementData< dim > | |
| void | set_primitivity (const bool value) |
Static Protected Member Functions inherited from FiniteElement< dim, dim > | |
| static std::vector< unsigned int > | compute_n_nonzero_components (const std::vector< ComponentMask > &nonzero_components) |
Protected Attributes inherited from FE_PolyTensor< PolynomialsRaviartThomas< dim >, dim > | |
| MappingType | mapping_type |
| PolynomialsRaviartThomas< dim > | poly_space |
| FullMatrix< double > | inverse_node_matrix |
| Point< dim > | cached_point |
| std::vector< Tensor< 1, dim > > | cached_values |
| std::vector< Tensor< 2, dim > > | cached_grads |
| std::vector< Tensor< 3, dim > > | cached_grad_grads |
Protected Attributes inherited from FiniteElement< dim, dim > | |
| std::vector< std::vector< FullMatrix< double > > > | restriction |
| std::vector< std::vector< FullMatrix< double > > > | prolongation |
| FullMatrix< double > | interface_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, int > | adjust_quad_dof_index_for_face_orientation_table |
| std::vector< int > | adjust_line_dof_index_for_line_orientation_table |
| std::vector< std::pair< unsigned int, unsigned int > > | system_to_component_table |
| std::vector< std::pair< unsigned int, unsigned int > > | face_system_to_component_table |
| std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > | system_to_base_table |
| std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > | face_system_to_base_table |
| BlockIndices | base_to_block_indices |
| std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > | component_to_base_table |
| const std::vector< bool > | restriction_is_additive_flags |
| const std::vector< ComponentMask > | nonzero_components |
| const std::vector< unsigned int > | n_nonzero_components_table |
The Raviart-Thomas elements with node functionals defined as point values in Gauss points.
For this Raviart-Thomas element, the node values are not cell and face moments with respect to certain polynomials, but the values in quadrature points. Following the general scheme for numbering degrees of freedom, the node values on edges are first, edge by edge, according to the natural ordering of the edges of a cell. The interior degrees of freedom are last.
For an RT-element of degree k, we choose (k+1)d-1 Gauss points on each face. These points are ordered lexicographically with respect to the orientation of the face. This way, the normal component which is in Qk is uniquely determined. Furthermore, since this Gauss-formula is exact on Q2k+1, these node values correspond to the exact integration of the moments of the RT-space.
In the interior of the cells, the moments are with respect to an anisotropic Qk space, where the test functions are one degree lower in the direction corresponding to the vector component under consideration. This is emulated by using an anisotropic Gauss formula for integration.
Definition at line 246 of file fe_raviart_thomas.h.
| FE_RaviartThomasNodal< dim >::FE_RaviartThomasNodal | ( | const unsigned int | p | ) |
Constructor for the Raviart-Thomas element of degree p.
Definition at line 36 of file fe_raviart_thomas_nodal.cc.
|
virtual |
Return a string that uniquely identifies a finite element. This class returns FE_RaviartThomasNodal<dim>(degree), with dim and degree replaced by appropriate values.
Implements FiniteElement< dim, dim >.
Definition at line 108 of file fe_raviart_thomas_nodal.cc.
|
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 copies of finite elements without knowing their exact type. They do so through this function.
Implements FiniteElement< dim, dim >.
Definition at line 130 of file fe_raviart_thomas_nodal.cc.
|
virtual |
Interpolate a set of scalar values, computed in the generalized support points.
Reimplemented from FiniteElement< dim, dim >.
Definition at line 292 of file fe_raviart_thomas_nodal.cc.
|
virtual |
Interpolate a set of vector values, computed in the generalized support points.
Since a finite element often only interpolates part of a vector, offset is used to determine the first component of the vector to be interpolated. Maybe consider changing your data structures to use the next function.
Reimplemented from FiniteElement< dim, dim >.
Definition at line 302 of file fe_raviart_thomas_nodal.cc.
|
virtual |
Interpolate a set of vector values, computed in the generalized support points.
Reimplemented from FiniteElement< dim, dim >.
Definition at line 350 of file fe_raviart_thomas_nodal.cc.
|
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, dim >.
Definition at line 397 of file fe_raviart_thomas_nodal.cc.
|
staticprivate |
Only for internal use. Its full name is get_dofs_per_object_vector function and it creates the dofs_per_object vector that is needed within the constructor to be passed to the constructor of FiniteElementData.
Definition at line 208 of file fe_raviart_thomas_nodal.cc.
|
staticprivate |
Compute the vector used for the restriction_is_additive field passed to the base class's constructor.
Definition at line 241 of file fe_raviart_thomas_nodal.cc.
|
privatevirtual |
This function returns true, if the shape function shape_index has non-zero function values somewhere on the face face_index.
Right now, this is only implemented for RT0 in 1D. Otherwise, returns always true.
Reimplemented from FiniteElement< dim, dim >.
Definition at line 263 of file fe_raviart_thomas_nodal.cc.
|
private |
Initialize the FiniteElement<dim>::generalized_support_points and FiniteElement<dim>::generalized_face_support_points fields. Called from the constructor.
See the glossary entry on generalized support points for more information.
Definition at line 144 of file fe_raviart_thomas_nodal.cc.
1.8.12