DOLFIN
DOLFIN C++ interface
Public Member Functions | List of all members
dolfin::FiniteElement Class Reference

This is a wrapper for a UFC finite element (ufc::finite_element). More...

#include <FiniteElement.h>

Public Member Functions

 FiniteElement (std::shared_ptr< const ufc::finite_element > element)
 
virtual ~FiniteElement ()
 Destructor.
 
std::string signature () const
 
ufc::shape cell_shape () const
 
std::size_t topological_dimension () const
 
virtual unsigned int geometric_dimension () const
 
std::size_t space_dimension () const
 
std::size_t value_rank () const
 Return the rank of the value space.
 
std::size_t value_dimension (std::size_t i) const
 Return the dimension of the value space for axis i.
 
void evaluate_basis (std::size_t i, double *values, const double *x, const double *coordinate_dofs, int cell_orientation) const
 Evaluate basis function i at given point in cell.
 
void evaluate_basis_all (double *values, const double *x, const double *coordinate_dofs, int cell_orientation) const
 Evaluate all basis functions at given point in cell.
 
void evaluate_basis_derivatives (unsigned int i, unsigned int n, double *values, const double *x, const double *coordinate_dofs, int cell_orientation) const
 Evaluate order n derivatives of basis function i at given point in cell.
 
void evaluate_basis_derivatives_all (unsigned int n, double *values, const double *x, const double *coordinate_dofs, int cell_orientation) const
 
double evaluate_dof (std::size_t i, const ufc::function &function, const double *coordinate_dofs, int cell_orientation, const ufc::cell &c) const
 Evaluate linear functional for dof i on the function f.
 
void evaluate_dofs (double *values, const ufc::function &f, const double *coordinate_dofs, int cell_orientation, const ufc::cell &c) const
 Evaluate linear functionals for all dofs on the function f.
 
void interpolate_vertex_values (double *vertex_values, double *coefficients, const double *coordinate_dofs, int cell_orientation) const
 
void tabulate_dof_coordinates (boost::multi_array< double, 2 > &coordinates, const std::vector< double > &coordinate_dofs, const Cell &cell) const
 
std::size_t num_sub_elements () const
 
std::size_t hash () const
 Return simple hash of the signature string.
 
std::shared_ptr< const FiniteElementcreate_sub_element (std::size_t i) const
 Create a new finite element for sub element i (for a mixed element)
 
std::shared_ptr< const FiniteElementcreate () const
 Create a new class instance.
 
std::shared_ptr< const FiniteElementextract_sub_element (const std::vector< std::size_t > &component) const
 Extract sub finite element for component.
 
std::shared_ptr< const ufc::finite_element > ufc_element () const
 

Detailed Description

This is a wrapper for a UFC finite element (ufc::finite_element).

Constructor & Destructor Documentation

FiniteElement::FiniteElement ( std::shared_ptr< const ufc::finite_element >  element)

Create finite element from UFC finite element (data may be shared)

Parameters
element(ufc::finite_element) UFC finite element

Member Function Documentation

ufc::shape dolfin::FiniteElement::cell_shape ( ) const
inline

Return the cell shape

Returns
ufc::shape
void dolfin::FiniteElement::evaluate_basis_derivatives_all ( unsigned int  n,
double *  values,
const double *  x,
const double *  coordinate_dofs,
int  cell_orientation 
) const
inline

Evaluate order n derivatives of all basis functions at given point in cell

virtual unsigned int dolfin::FiniteElement::geometric_dimension ( ) const
inlinevirtual

Return the geometric dimension of the cell shape

Returns
unsigned int
void dolfin::FiniteElement::interpolate_vertex_values ( double *  vertex_values,
double *  coefficients,
const double *  coordinate_dofs,
int  cell_orientation 
) const
inline

Interpolate vertex values from dof values

Parameters
vertex_values(double*)
coefficients(double*)
coordinate_dofs(const double*)
cell_orientation(int)
std::size_t dolfin::FiniteElement::num_sub_elements ( ) const
inline

Return the number of sub elements (for a mixed element)

Returns
std::size_t number of sub-elements
std::string dolfin::FiniteElement::signature ( ) const
inline

Return a string identifying the finite element

Returns
std::string
std::size_t dolfin::FiniteElement::space_dimension ( ) const
inline

Return the dimension of the finite element function space

Returns
std::size_t
void FiniteElement::tabulate_dof_coordinates ( boost::multi_array< double, 2 > &  coordinates,
const std::vector< double > &  coordinate_dofs,
const Cell cell 
) const

Tabulate the coordinates of all dofs on an element

Parameters
[in,out]coordinates(boost::multi_array<double, 2>) The coordinates of all dofs on a cell.
[in]coordinate_dofs(std::vector<double>) The cell coordinates
[in]cell(Cell) The cell.
std::size_t dolfin::FiniteElement::topological_dimension ( ) const
inline

Return the topological dimension of the cell shape

Returns
std::size_t
std::shared_ptr<const ufc::finite_element> dolfin::FiniteElement::ufc_element ( ) const
inline

Return underlying UFC element. Intended for libray usage only and may change.


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