DOLFIN-X
DOLFIN-X C++ interface
Public Member Functions | List of all members
dolfinx::fem::FiniteElement Class Reference

Finite Element, containing the dof layout on a reference element, and various methods for evaluating and transforming the basis. More...

#include <FiniteElement.h>

Public Member Functions

 FiniteElement (const ufc_finite_element &ufc_element)
 Create finite element from UFC finite element. More...
 
virtual ~FiniteElement ()=default
 Destructor.
 
std::string signature () const
 String identifying the finite element. More...
 
mesh::CellType cell_shape () const
 Cell shape. More...
 
int space_dimension () const
 Dimension of the finite element function space. More...
 
int value_size () const
 The value size, e.g. 1 for a scalar function, 2 for a 2D vector. More...
 
int reference_value_size () const
 The value size, e.g. 1 for a scalar function, 2 for a 2D vector for the reference element. More...
 
int value_rank () const
 Rank of the value space. More...
 
int value_dimension (int i) const
 Return the dimension of the value space for axis i.
 
int degree () const
 Return the maximum polynomial degree.
 
std::string family () const
 The finite element family. More...
 
void evaluate_reference_basis (Eigen::Tensor< double, 3, Eigen::RowMajor > &reference_values, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &X) const
 Evaluate all basis functions at given point in reference cell.
 
void transform_reference_basis (Eigen::Tensor< double, 3, Eigen::RowMajor > &values, const Eigen::Tensor< double, 3, Eigen::RowMajor > &reference_values, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &X, const Eigen::Tensor< double, 3, Eigen::RowMajor > &J, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, 1 >> &detJ, const Eigen::Tensor< double, 3, Eigen::RowMajor > &K, const std::uint32_t permutation_info) const
 Push basis functions forward to physical element.
 
void transform_reference_basis_derivatives (Eigen::Tensor< double, 4, Eigen::RowMajor > &values, std::size_t order, const Eigen::Tensor< double, 4, Eigen::RowMajor > &reference_values, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &X, const Eigen::Tensor< double, 3, Eigen::RowMajor > &J, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, 1 >> &detJ, const Eigen::Tensor< double, 3, Eigen::RowMajor > &K, const std::uint32_t permutation_info) const
 Push basis function (derivatives) forward to physical element.
 
bool has_dof_reference_coordinates () const noexcept
 Check if reference coordinates for dofs are defined. More...
 
const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > & dof_reference_coordinates () const
 Tabulate the reference coordinates of all dofs on an element. More...
 
void transform_values (PetscScalar *reference_values, const Eigen::Ref< const Eigen::Array< PetscScalar, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &physical_values, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &coordinate_dofs) const
 Map values of field from physical to reference space which has been evaluated at points given by dof_reference_coordinates()
 
int num_sub_elements () const
 Return the number of sub elements (for a mixed element)
 
std::size_t hash () const
 Return simple hash of the signature string.
 
std::shared_ptr< const FiniteElementextract_sub_element (const std::vector< int > &component) const
 Extract sub finite element for component.
 

Detailed Description

Finite Element, containing the dof layout on a reference element, and various methods for evaluating and transforming the basis.

Constructor & Destructor Documentation

◆ FiniteElement()

FiniteElement::FiniteElement ( const ufc_finite_element &  ufc_element)

Create finite element from UFC finite element.

Parameters
[in]ufc_elementUFC finite element

Member Function Documentation

◆ cell_shape()

mesh::CellType FiniteElement::cell_shape ( ) const

Cell shape.

Returns
Element cell shape

◆ dof_reference_coordinates()

const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > & FiniteElement::dof_reference_coordinates ( ) const

Tabulate the reference coordinates of all dofs on an element.

Returns
The coordinates of all dofs on the reference cell

◆ family()

std::string FiniteElement::family ( ) const

The finite element family.

Returns
The string of the finite element family

◆ has_dof_reference_coordinates()

bool FiniteElement::has_dof_reference_coordinates ( ) const
noexcept

Check if reference coordinates for dofs are defined.

Returns
True if the dof coordinates are available

◆ reference_value_size()

int FiniteElement::reference_value_size ( ) const

The value size, e.g. 1 for a scalar function, 2 for a 2D vector for the reference element.

Returns
The value size for the reference element

◆ signature()

std::string FiniteElement::signature ( ) const

String identifying the finite element.

Returns
Element signature

◆ space_dimension()

int FiniteElement::space_dimension ( ) const

Dimension of the finite element function space.

Returns
Dimension of the finite element space

◆ value_rank()

int FiniteElement::value_rank ( ) const

Rank of the value space.

Returns
The value rank

◆ value_size()

int FiniteElement::value_size ( ) const

The value size, e.g. 1 for a scalar function, 2 for a 2D vector.

Returns
The value size

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