Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Public Member Functions | Private Member Functions | Static Private Member Functions | List of all members
FE_P1NC Class Reference

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

Inheritance diagram for FE_P1NC:
[legend]

Public Member Functions

 FE_P1NC ()
 
virtual std::string get_name () const override
 
virtual UpdateFlags requires_update_flags (const UpdateFlags flags) const override
 
virtual std::unique_ptr< FiniteElement< 2, 2 > > clone () const override
 
virtual ~FE_P1NC () override=default
 
- Public Member Functions inherited from FiniteElement< 2, 2 >
 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 intoperator^ (const unsigned int multiplicity) const
 
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone () const=0
 
virtual std::string get_name () const=0
 
const FiniteElement< dim, spacedim > & operator[] (const unsigned int fe_index) const
 
virtual double shape_value (const unsigned int i, const Point< dim > &p) const
 
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 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 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 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 final
 
virtual FiniteElementDomination::Domination compare_for_domination (const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const
 
virtual bool operator== (const FiniteElement< dim, spacedim > &fe) const
 
bool operator!= (const FiniteElement< dim, spacedim > &) const
 
std::pair< unsigned int, unsigned intsystem_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 intface_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 intsystem_to_base_index (const unsigned int index) const
 
std::pair< std::pair< unsigned int, unsigned int >, unsigned intface_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 intcomponent_to_base_index (const unsigned int component) const
 
std::pair< unsigned int, unsigned intblock_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
 
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
 
GeometryPrimitive get_associated_geometry_primitive (const unsigned int cell_dof_index) const
 
virtual void convert_generalized_support_point_values_to_dof_values (const std::vector< Vector< double >> &support_point_values, std::vector< double > &nodal_values) const
 
virtual std::size_t memory_consumption () const
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") 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

virtual std::unique_ptr< FiniteElement< 2, 2 >::InternalDataBaseget_data (const UpdateFlags update_flags, const Mapping< 2, 2 > &, const Quadrature< 2 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &output_data) const override
 
virtual std::unique_ptr< FiniteElement< 2, 2 >::InternalDataBaseget_face_data (const UpdateFlags update_flags, const Mapping< 2, 2 > &, const Quadrature< 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &output_data) const override
 
virtual std::unique_ptr< FiniteElement< 2, 2 >::InternalDataBaseget_subface_data (const UpdateFlags update_flags, const Mapping< 2, 2 > &, const Quadrature< 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &output_data) const override
 
virtual void fill_fe_values (const Triangulation< 2, 2 >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< 2 > &quadrature, const Mapping< 2, 2 > &mapping, const Mapping< 2, 2 >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< 2, 2 > &mapping_data, const FiniteElement< 2, 2 >::InternalDataBase &fe_internal, internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &output_data) const override
 
virtual void fill_fe_face_values (const Triangulation< 2, 2 >::cell_iterator &cell, const unsigned int face_no, const Quadrature< 1 > &quadrature, const Mapping< 2, 2 > &mapping, const Mapping< 2, 2 >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< 2, 2 > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &output_data) const override
 
virtual void fill_fe_subface_values (const Triangulation< 2, 2 >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< 1 > &quadrature, const Mapping< 2, 2 > &mapping, const Mapping< 2, 2 >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< 2, 2 > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &output_data) const override
 
void initialize_constraints ()
 

Static Private Member Functions

static std::vector< unsigned intget_dpo_vector ()
 
static std::array< std::array< double, 3 >, 4 > get_linear_shape_coefficients (const Triangulation< 2, 2 >::cell_iterator &cell)
 

Additional Inherited Members

- Public Types inherited from FiniteElementData< dim >
enum  Conformity {
  unknown = 0x00, L2 = 0x01, Hcurl = 0x02, Hdiv = 0x04,
  H1 = Hcurl | Hdiv, H2 = 0x0e
}
 
- Static Public Member Functions inherited from FiniteElement< 2, 2 >
static ::ExceptionBaseExcShapeFunctionNotPrimitive (int arg1)
 
static ::ExceptionBaseExcFENotPrimitive ()
 
static ::ExceptionBaseExcUnitShapeValuesDoNotExist ()
 
static ::ExceptionBaseExcFEHasNoSupportPoints ()
 
static ::ExceptionBaseExcEmbeddingVoid ()
 
static ::ExceptionBaseExcProjectionVoid ()
 
static ::ExceptionBaseExcWrongInterfaceMatrixSize (int arg1, int arg2)
 
static ::ExceptionBaseExcInterpolationNotImplemented ()
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (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< 2, 2 >
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< 2, 2 >
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< InternalDataBaseget_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< InternalDataBaseget_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< InternalDataBaseget_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< 2, 2 >
static std::vector< unsigned intcompute_n_nonzero_components (const std::vector< ComponentMask > &nonzero_components)
 
- Protected Attributes inherited from FiniteElement< 2, 2 >
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
 
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< boolrestriction_is_additive_flags
 
const std::vector< ComponentMasknonzero_components
 
const std::vector< unsigned intn_nonzero_components_table
 
const bool cached_primitivity
 

Detailed Description

Implementation of the scalar version of the P1 nonconforming finite element, a piecewise linear element on quadrilaterals in 2D. This implementation is only for 2D cells in a 2D space (i.e., codimension 0).

Unlike the usual continuous, \(H^1\) conforming finite elements, the P1 nonconforming element does not enforce continuity across edges. However, it requires the continuity in an integral sense: any function in the space should have the same integral values on two sides of the common edge shared by two adjacent elements.

Thus, each function in the nonconforming element space can be discontinuous, and consequently not included in \(H^1_0\), just like the basis functions in Discontinuous Galerkin (DG) finite element spaces. On the other hand, basis functions in DG spaces are completely discontinuous across edges without any relation between the values from both sides. This is a reason why usual weak formulations for DG schemes contain additional penalty terms for jump across edges to control discontinuity. However, nonconforming elements usually do not need additional terms in their weak formulations because their integrals along edges are the same from both sides, i.e., there is some level of continuity.

Dice Rule

Since any function in the P1 nonconforming space is piecewise linear on each element, the function value at the midpoint of each edge is same as the mean value on the edge. Thus the continuity of the integral value across each edge is equivalent to the continuity of the midpoint value of each edge in this case.

Thus for the P1 nonconforming element, the function values at midpoints on edges of a cell are important. The first attempt to define (local) degrees of freedom (DoFs) on a quadrilateral is by using midpoint values of a function.

However, these 4 functionals are not linearly independent because a linear function on 2D is uniquely determined by only 3 independent values. A simple observation reads that any linear function on a quadrilateral should satisfy the 'dice rule': the sum of two function values at the midpoints of the edge pair on opposite sides of a cell is equal to the sum of those at the midpoints of the other edge pair. This is called the 'dice rule' because the number of points on opposite sides of a dice always adds up to the same number as well (in the case of dice, to seven).

In formulas, the dice rule is written as \(\phi(m_0) + \phi(m_1) = \phi(m_2) + \phi(m_3)\) for all \(\phi\) in the function space where \(m_j\) is the midpoint of the edge \(e_j\). Here, we assume the standard numbering convention for edges used in deal.II and described in class GeometryInfo.

Conversely if 4 values at midpoints satisfying the dice rule are given, then there always exists the unique linear function which coincides with 4 midpoints values.

Due to the dice rule, three values at any three midpoints can determine the last value at the last midpoint. It means that the number of independent local functionals on a cell is 3, and this is also the dimension of the linear polynomial space on a cell in 2D.

Shape functions

Before introducing the degrees of freedom, we present 4 local shape functions on a cell. Due to the dice rule, we need a special construction for shape functions. Although the following 4 shape functions are not linearly independent within a cell, they are helpful to define the global basis functions which are linearly independent on the whole domain. Again, we assume the standard numbering for vertices used in deal.II.

*  2---------|---------3
*  |                   |
*  |                   |
*  |                   |
*  |                   |
*  -                   -
*  |                   |
*  |                   |
*  |                   |
*  |                   |
*  0---------|---------1
*  

For each vertex \(v_j\) of given cell, there are two edges of which \(v_j\) is one of end points. Consider a linear function such that it has value 0.5 at the midpoints of two adjacent edges, and 0.0 at the two midpoints of the other edges. Note that the set of these values satisfies the dice rule which is described above. We denote such a function associated with vertex \(v_j\) by \(\phi_j\). Then the set of 4 shape functions is a partition of unity on a cell: \(\sum_{j=0}^{3} \phi_j = 1\). (This is easy to see: at each edge midpoint, the sum of the four function adds up to one because two functions have value 0.5 and the other value 0.0. Because the function is globally linear, the only function that can have value 1 at four points must also be globally equal to one.)

The following figures represent \(\phi_j\) for \(j=0,\cdots,3\) with their midpoint values:

The local DoFs are defined by the coefficients of the shape functions associated with vertices, respectively. Although these 4 local DoFs are not linearly independent within a single cell, this definition is a good start point for the definition of the global DoFs.

We want to emphasize that the shape functions are constructed on each cell, not on the reference cell only. Usual finite elements are defined based on a 'parametric' concept. It means that a function space for a finite element is defined on one reference cell, and it is transformed into each cell via a mapping from the reference cell. However the P1 nonconforming element does not follow such concept. It defines a function space with linear shape functions on each cell without any help of a function space on the reference cell. In other words, the element is defined in real space, not via a mapping from a reference cell. In this, it is similar to the FE_DGPNonparametric element.

Thus this implementation does not have to compute shape values on the reference cell. Rather, the shape values are computed by construction of the shape functions on each cell independently.

Degrees of freedom

We next have to consider the global basis functions for the element because the system of equations which we ultimately have to solve is for a global system, not local. The global basis functions associated with a node are defined by a cell-wise composition of local shape functions associated with the node on each element.

There is a theoretical result about the linear independence of the global basis functions depending on the type of the boundary condition we consider.

When homogeneous Dirichlet boundary conditions are given, the global basis functions associated with interior nodes are linearly independent. Then, the number of DoFs is equal to the number of interior nodes, and consequently the same as the number of DoFs for the standard bilinear \(Q_1\) finite element.

When Neumann boundary conditions are given, the global basis functions associated with all nodes (including boundary nodes) are actually not linearly independent. There exists one redundancy. Thus in this case, the number of DoFs is equal to the number of all nodes minus 1. This is, again as for the regular \(Q_1\) element.

Unit support points

For a smooth function, we construct a piecewise linear function which belongs to the element space by using its nodal values as DoF values.

Note that for the P1 nonconforming element, two nodal values of a smooth function and its interpolant do not coincide in general, in contrast to ordinary Lagrange finite elements. Of course, it is meaningless to refer 'nodal value' because the element space has nonconformity. But it is also true even though the single global basis function associated with a node is considered the unique 'nodal value' at the node. For instance, consider the basis function associated with a node. Consider two lines representing the level sets for value 0.5 and 0, respectively, by connecting two midpoints. Then we cut the quad into two sub-triangles by the diagonal which is placed along those two lines. It gives another level set for value 0.25 which coincides with the cutting diagonal. Therefore these three level sets are all parallel in the quad and it gives the value 0.75 at the base node, not value 1. This is true whether the quadrilateral is a rectangle, parallelogram, or any other shape.

Reference

The original paper for the P1 nonconforming element by Park and Sheen is accessible at https://doi.org/10.1137/S0036142902404923 , see [54] .

Author
Jaeryun Yim, 2015, 2016.

Definition at line 259 of file fe_p1nc.h.

Constructor & Destructor Documentation

◆ FE_P1NC()

FE_P1NC::FE_P1NC ( )

Constructor for the P1 nonconforming element. It is only for 2D and codimension = 0.

Definition at line 24 of file fe_p1nc.cc.

◆ ~FE_P1NC()

virtual FE_P1NC::~FE_P1NC ( )
overridevirtualdefault

Destructor.

Member Function Documentation

◆ get_name()

std::string FE_P1NC::get_name ( ) const
overridevirtual

Definition at line 42 of file fe_p1nc.cc.

◆ requires_update_flags()

UpdateFlags FE_P1NC::requires_update_flags ( const UpdateFlags  update_flags) const
overridevirtual

Given a set of update flags, compute which other quantities also need to be computed in order to satisfy the request by the given flags. Then return the combination of the original set of flags and those just computed.

As an example, if update_flags contains update_gradients a finite element class will typically require the computation of the inverse of the Jacobian matrix in order to rotate the gradient of shape functions on the reference cell to the real cell. It would then return not just update_gradients, but also update_covariant_transformation, the flag that makes the mapping class produce the inverse of the Jacobian matrix.

An extensive discussion of the interaction between this function and FEValues can be found in the How Mapping, FiniteElement, and FEValues work together documentation module.

See also
UpdateFlags

Implements FiniteElement< 2, 2 >.

Definition at line 50 of file fe_p1nc.cc.

◆ clone()

std::unique_ptr< FiniteElement< 2, 2 > > FE_P1NC::clone ( ) const
overridevirtual

Definition at line 69 of file fe_p1nc.cc.

◆ get_dpo_vector()

std::vector< unsigned int > FE_P1NC::get_dpo_vector ( )
staticprivate

Return the vector consists of the numbers of degrees of freedom per objects.

Definition at line 77 of file fe_p1nc.cc.

◆ get_linear_shape_coefficients()

std::array< std::array< double, 3 >, 4 > FE_P1NC::get_linear_shape_coefficients ( const Triangulation< 2, 2 >::cell_iterator &  cell)
staticprivate

Return the coefficients of 4 local linear shape functions \(\phi_j(x,y) = a x + b y + c\) on given cell. For each local shape function, the array consists of three coefficients is in order of a,b and c.

Definition at line 89 of file fe_p1nc.cc.

◆ get_data()

std::unique_ptr< FiniteElement< 2, 2 >::InternalDataBase > FE_P1NC::get_data ( const UpdateFlags  update_flags,
const Mapping< 2, 2 > &  ,
const Quadrature< 2 > &  quadrature,
::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &  output_data 
) const
overrideprivatevirtual

Do the work which is needed before cellwise data computation. Since the shape functions are constructed independently on each cell, the data on the reference cell is not necessary. It returns an empty variable type of @ InternalDataBase and updates @ update_flags, and computes trivially zero Hessian for each cell if it is needed.

Definition at line 136 of file fe_p1nc.cc.

◆ get_face_data()

std::unique_ptr< FiniteElement< 2, 2 >::InternalDataBase > FE_P1NC::get_face_data ( const UpdateFlags  update_flags,
const Mapping< 2, 2 > &  ,
const Quadrature< 1 > &  quadrature,
::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &  output_data 
) const
overrideprivatevirtual

Definition at line 161 of file fe_p1nc.cc.

◆ get_subface_data()

std::unique_ptr< FiniteElement< 2, 2 >::InternalDataBase > FE_P1NC::get_subface_data ( const UpdateFlags  update_flags,
const Mapping< 2, 2 > &  ,
const Quadrature< 1 > &  quadrature,
::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &  output_data 
) const
overrideprivatevirtual

Definition at line 186 of file fe_p1nc.cc.

◆ fill_fe_values()

void FE_P1NC::fill_fe_values ( const Triangulation< 2, 2 >::cell_iterator &  cell,
const CellSimilarity::Similarity  cell_similarity,
const Quadrature< 2 > &  quadrature,
const Mapping< 2, 2 > &  mapping,
const Mapping< 2, 2 >::InternalDataBase mapping_internal,
const internal::FEValuesImplementation::MappingRelatedData< 2, 2 > &  mapping_data,
const FiniteElement< 2, 2 >::InternalDataBase fe_internal,
internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &  output_data 
) const
overrideprivatevirtual

Compute the data on the current cell.

Definition at line 211 of file fe_p1nc.cc.

◆ fill_fe_face_values()

void FE_P1NC::fill_fe_face_values ( const Triangulation< 2, 2 >::cell_iterator &  cell,
const unsigned int  face_no,
const Quadrature< 1 > &  quadrature,
const Mapping< 2, 2 > &  mapping,
const Mapping< 2, 2 >::InternalDataBase mapping_internal,
const ::internal::FEValuesImplementation::MappingRelatedData< 2, 2 > &  mapping_data,
const InternalDataBase fe_internal,
::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &  output_data 
) const
overrideprivatevirtual

Compute the data on the face of the current cell.

Definition at line 249 of file fe_p1nc.cc.

◆ fill_fe_subface_values()

void FE_P1NC::fill_fe_subface_values ( const Triangulation< 2, 2 >::cell_iterator &  cell,
const unsigned int  face_no,
const unsigned int  sub_no,
const Quadrature< 1 > &  quadrature,
const Mapping< 2, 2 > &  mapping,
const Mapping< 2, 2 >::InternalDataBase mapping_internal,
const ::internal::FEValuesImplementation::MappingRelatedData< 2, 2 > &  mapping_data,
const InternalDataBase fe_internal,
::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &  output_data 
) const
overrideprivatevirtual

Compute the data on the subface of the current cell.

Definition at line 293 of file fe_p1nc.cc.

◆ initialize_constraints()

void FE_P1NC::initialize_constraints ( )
private

Create the constraints matrix for hanging edges.

Definition at line 340 of file fe_p1nc.cc.


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