Reference documentation for deal.II version 9.1.0-pre
Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes | Friends | List of all members
FE_RaviartThomas< dim > Class Template Reference

#include <deal.II/fe/fe_raviart_thomas.h>

Inheritance diagram for FE_RaviartThomas< dim >:
[legend]

Public Member Functions

 FE_RaviartThomas (const unsigned int p)
 
virtual std::string get_name () const override
 
virtual std::unique_ptr< FiniteElement< dim, dim > > clone () const override
 
virtual bool has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const override
 
virtual void convert_generalized_support_point_values_to_dof_values (const std::vector< Vector< double >> &support_point_values, std::vector< double > &nodal_values) const override
 
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes () const override
 
virtual std::size_t memory_consumption () const override
 
- 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 override
 
virtual Tensor< 1, dim > shape_grad (const unsigned int i, const Point< dim > &p) const override
 
virtual Tensor< 2, dim > shape_grad_grad (const unsigned int i, const Point< dim > &p) const override
 
- 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)
 
 FiniteElement (FiniteElement< dim, spacedim > &&)=default
 
 FiniteElement (const FiniteElement< dim, spacedim > &)=default
 
virtual ~FiniteElement () override=default
 
std::pair< std::unique_ptr< FiniteElement< dim, spacedim > >, unsigned int > operator^ (const unsigned int multiplicity) const
 
const FiniteElement< dim, spacedim > & operator[] (const unsigned int fe_index) const
 
virtual bool operator== (const FiniteElement< dim, spacedim > &fe) const
 
bool operator!= (const FiniteElement< dim, spacedim > &) const
 
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 bool hp_constraints_are_implemented () 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 ComponentMaskget_nonzero_components (const unsigned int i) const
 
unsigned int n_nonzero_components (const unsigned int i) const
 
bool is_primitive () 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
 
const FiniteElement< dim, spacedim > & get_sub_fe (const ComponentMask &mask) const
 
virtual const FiniteElement< dim, spacedim > & get_sub_fe (const unsigned int first_component, const unsigned int n_selected_components) 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_indexsystem_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
 
GeometryPrimitive get_associated_geometry_primitive (const unsigned int cell_dof_index) const
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
void subscribe (const char *identifier=nullptr) const
 
void unsubscribe (const char *identifier=nullptr) const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
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 BlockIndicesblock_indices () const
 
unsigned int tensor_degree () const
 
bool conforms (const Conformity) const
 
bool operator== (const FiniteElementData &) const
 

Private Member Functions

void initialize_support_points (const unsigned int rt_degree)
 
void initialize_restriction ()
 

Static Private Member Functions

static std::vector< unsigned int > get_dpo_vector (const unsigned int degree)
 

Private Attributes

Table< 2, double > boundary_weights
 
Table< 3, double > interior_weights
 

Friends

template<int dim1>
class FE_RaviartThomas
 

Additional Inherited Members

- Public Types inherited from FiniteElementData< dim >
- Static Public Member Functions inherited from FiniteElement< dim, dim >
static::ExceptionBase & ExcShapeFunctionNotPrimitive (int arg1)
 
static::ExceptionBase & ExcFENotPrimitive ()
 
static::ExceptionBase & ExcUnitShapeValuesDoNotExist ()
 
static::ExceptionBase & ExcFEHasNoSupportPoints ()
 
static::ExceptionBase & ExcEmbeddingVoid ()
 
static::ExceptionBase & ExcProjectionVoid ()
 
static::ExceptionBase & ExcWrongInterfaceMatrixSize (int arg1, int arg2)
 
static::ExceptionBase & ExcInterpolationNotImplemented ()
 
- Static Public Member Functions inherited from Subscriptor
static::ExceptionBase & ExcInUse (int arg1, std::string arg2, std::string arg3)
 
static::ExceptionBase & ExcNoSubscriber (std::string arg1, std::string arg2)
 
- 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 std::unique_ptr< InternalDataBase > get_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature,::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
 
virtual std::unique_ptr< InternalDataBase > get_face_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature,::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
 
virtual std::unique_ptr< InternalDataBase > get_subface_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature,::internal::FEValuesImplementation::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::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal,::internal::FEValuesImplementation::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::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal,::internal::FEValuesImplementation::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::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal,::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
 
- 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
 
Threads::Mutex cache_mutex
 
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< ComponentMasknonzero_components
 
const std::vector< unsigned int > n_nonzero_components_table
 
const bool cached_primitivity
 

Detailed Description

template<int dim>
class FE_RaviartThomas< dim >

Implementation of Raviart-Thomas (RT) elements, conforming with the space Hdiv. These elements generate vector fields with normal components continuous between mesh cells.

We follow the usual definition of the degree of RT elements, which denotes the polynomial degree of the largest complete polynomial subspace contained in the RT space. Then, approximation order of the function itself is degree+1, as with usual polynomial spaces. The numbering so chosen implies the sequence

\[ Q_{k+1} \stackrel{\text{grad}}{\rightarrow} \text{Nedelec}_k \stackrel{\text{curl}}{\rightarrow} \text{RaviartThomas}_k \stackrel{\text{div}}{\rightarrow} DGQ_{k} \]

The lowest order element is consequently FE_RaviartThomas(0).

This class is not implemented for the codimension one case (spacedim != dim).

Todo:
Even if this element is implemented for two and three space dimensions, the definition of the node values relies on consistently oriented faces in 3D. Therefore, care should be taken on complicated meshes.

Interpolation

The interpolation operators associated with the RT element are constructed such that interpolation and computing the divergence are commuting operations. We require this from interpolating arbitrary functions as well as the restriction matrices. It can be achieved by two interpolation schemes, the simplified one in FE_RaviartThomasNodal and the original one here:

Node values on edges/faces

On edges or faces, the node values are the moments of the normal component of the interpolated function with respect to the traces of the RT polynomials. Since the normal trace of the RT space of degree k on an edge/face is the space Qk, the moments are taken with respect to this space.

Interior node values

Higher order RT spaces have interior nodes. These are moments taken with respect to the gradient of functions in Qk on the cell (this space is the matching space for RTk in a mixed formulation).

Generalized support points

The node values above rely on integrals, which will be computed by quadrature rules themselves. The generalized support points are a set of points such that this quadrature can be performed with sufficient accuracy. The points needed are those of QGaussk+1 on each face as well as QGaussk+1 in the interior of the cell (or none for RT0).

Author
Guido Kanschat, 2005, based on previous work by Wolfgang Bangerth.

Definition at line 105 of file fe_raviart_thomas.h.

Constructor & Destructor Documentation

template<int dim>
FE_RaviartThomas< dim >::FE_RaviartThomas ( const unsigned int  p)

Constructor for the Raviart-Thomas element of degree p.

Definition at line 47 of file fe_raviart_thomas.cc.

Member Function Documentation

template<int dim>
std::string FE_RaviartThomas< dim >::get_name ( ) const
overridevirtual

Return a string that uniquely identifies a finite element. This class returns FE_RaviartThomas<dim>(degree), with dim and degree replaced by appropriate values.

Implements FiniteElement< dim, dim >.

Definition at line 115 of file fe_raviart_thomas.cc.

template<int dim>
std::unique_ptr< FiniteElement< dim, dim > > FE_RaviartThomas< dim >::clone ( ) const
overridevirtual

A sort of virtual copy constructor, this function returns a copy of the finite element object. Derived classes need to override the function here in this base class and return an object of the same type as the derived class.

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 137 of file fe_raviart_thomas.cc.

template<int dim>
bool FE_RaviartThomas< dim >::has_support_on_face ( const unsigned int  shape_index,
const unsigned int  face_index 
) const
overridevirtual

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 442 of file fe_raviart_thomas.cc.

template<int dim>
void FE_RaviartThomas< dim >::convert_generalized_support_point_values_to_dof_values ( const std::vector< Vector< double >> &  support_point_values,
std::vector< double > &  nodal_values 
) const
overridevirtual

Given the values of a function \(f(\mathbf x)\) at the (generalized) support points of the reference cell, this function then computes what the nodal values of the element are, i.e., \(\Psi_i[f]\), where \(\Psi_i\) are the node functionals of the element (see also Node values or node functionals). The values \(\Psi_i[f]\) are then the expansion coefficients for the shape functions of the finite element function that interpolates the given function \(f(x)\), i.e., \( f_h(\mathbf x) = \sum_i \Psi_i[f] \varphi_i(\mathbf x) \) is the finite element interpolant of \(f\) with the current element. The operation described here is used, for example, in the FETools::compute_node_matrix() function.

In more detail, let us assume that the generalized support points (see this glossary entry ) of the current element are \(\hat{\mathbf x}_i\) and that the node functionals associated with the current element are \(\Psi_i[\cdot]\). Then, the fact that the element is based on generalized support points, implies that if we apply \(\Psi_i\) to a (possibly vector-valued) finite element function \(\varphi\), the result must have the form \(\Psi_i[\varphi] = f_i(\varphi(\hat{\mathbf x}_i))\) – in other words, the value of the node functional \(\Psi_i\) applied to \(\varphi\) only depends on the values of \(\varphi\) at \(\hat{\mathbf x}_i\) and not on values anywhere else, or integrals of \(\varphi\), or any other kind of information.

The exact form of \(f_i\) depends on the element. For example, for scalar Lagrange elements, we have that in fact \(\Psi_i[\varphi] = \varphi(\hat{\mathbf x}_i)\). If you combine multiple scalar Lagrange elements via an FESystem object, then \(\Psi_i[\varphi] = \varphi(\hat{\mathbf x}_i)_{c(i)}\) where \(c(i)\) is the result of the FiniteElement::system_to_component_index() function's return value's first component. In these two cases, \(f_i\) is therefore simply the identity (in the scalar case) or a function that selects a particular vector component of its argument. On the other hand, for Raviart-Thomas elements, one would have that \(f_i(\mathbf y) = \mathbf y \cdot \mathbf n_i\) where \(\mathbf n_i\) is the normal vector of the face at which the shape function is defined.

Given all of this, what this function does is the following: If you input a list of values of a function \(\varphi\) at all generalized support points (where each value is in fact a vector of values with as many components as the element has), then this function returns a vector of values obtained by applying the node functionals to these values. In other words, if you pass in \(\{\varphi(\hat{\mathbf x}_i)\}_{i=0}^{N-1}\) then you will get out a vector \(\{\Psi[\varphi]\}_{i=0}^{N-1}\) where \(N\) equals dofs_per_cell.

Parameters
[in]support_point_valuesAn array of size dofs_per_cell (which equals the number of points the get_generalized_support_points() function will return) where each element is a vector with as many entries as the element has vector components. This array should contain the values of a function at the generalized support points of the current element.
[out]nodal_valuesAn array of size dofs_per_cell that contains the node functionals of the element applied to the given function.
Note
It is safe to call this function for (transformed) values on the real cell only for elements with trivial MappingType. For all other elements (for example for H(curl), or H(div) conforming elements) vector values have to be transformed to the reference cell first.
Given what the function is supposed to do, the function clearly can only work for elements that actually implement (generalized) support points. Elements that do not have generalized support points – e.g., elements whose nodal functionals evaluate integrals or moments of functions (such as FE_Q_Hierarchical) – can in general not make sense of the operation that is required for this function. They consequently may not implement it.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 486 of file fe_raviart_thomas.cc.

template<int dim>
std::pair< Table< 2, bool >, std::vector< unsigned int > > FE_RaviartThomas< dim >::get_constant_modes ( ) const
overridevirtual

Return a list of constant modes of the element. This method is currently not correctly implemented because it returns ones for all components.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 420 of file fe_raviart_thomas.cc.

template<int dim>
std::size_t FE_RaviartThomas< dim >::memory_consumption ( ) const
overridevirtual

Determine an estimate for the memory consumption (in bytes) of this object.

This function is made virtual, since finite element objects are usually accessed through pointers to their base class, rather than the class itself.

Reimplemented from FiniteElement< dim, dim >.

Definition at line 529 of file fe_raviart_thomas.cc.

template<int dim>
std::vector< unsigned int > FE_RaviartThomas< dim >::get_dpo_vector ( const unsigned int  degree)
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 398 of file fe_raviart_thomas.cc.

template<int dim>
void FE_RaviartThomas< dim >::initialize_support_points ( const unsigned int  rt_degree)
private

Initialize the generalized_support_points field of the FiniteElement class and fill the tables with interpolation weights (boundary_weights and interior_weights). Called from the constructor.

Definition at line 150 of file fe_raviart_thomas.cc.

template<int dim>
void FE_RaviartThomas< dim >::initialize_restriction ( )
private

Initialize the interpolation from functions on refined mesh cells onto the father cell. According to the philosophy of the Raviart-Thomas element, this restriction operator preserves the divergence of a function weakly.

Definition at line 266 of file fe_raviart_thomas.cc.

Friends And Related Function Documentation

template<int dim>
template<int dim1>
friend class FE_RaviartThomas
friend

Allow access from other dimensions.

Definition at line 204 of file fe_raviart_thomas.h.

Member Data Documentation

template<int dim>
Table<2, double> FE_RaviartThomas< dim >::boundary_weights
private

These are the factors multiplied to a function in the generalized_face_support_points when computing the integration. They are organized such that there is one row for each generalized face support point and one column for each degree of freedom on the face.

See the glossary entry on generalized support points for more information.

Definition at line 190 of file fe_raviart_thomas.h.

template<int dim>
Table<3, double> FE_RaviartThomas< dim >::interior_weights
private

Precomputed factors for interpolation of interior degrees of freedom. The rationale for this Table is the same as for boundary_weights. Only, this table has a third coordinate for the space direction of the component evaluated.

Definition at line 198 of file fe_raviart_thomas.h.


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