Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/v0.1.0/v0.9.0/cpp
DOLFINx  0.1.0
DOLFINx C++ interface
Public Types | Public Member Functions | List of all members
dolfinx::fem::Expression< T > Class Template Reference

Represents a mathematical expression evaluated at a pre-defined set of points on the reference cell. This class closely follows the concept of a UFC Expression. More...

#include <Expression.h>

Public Types

using scalar_type = T
 Scalar type (T).
 

Public Member Functions

 Expression (const std::vector< std::shared_ptr< const fem::Function< T >>> &coefficients, const std::vector< std::shared_ptr< const fem::Constant< T >>> &constants, const std::shared_ptr< const mesh::Mesh > &mesh, const array2d< double > &X, const std::function< void(T *, const T *, const T *, const double *)> fn, const std::size_t value_size)
 Create Expression. More...
 
 Expression (Expression &&form)=default
 Move constructor.
 
virtual ~Expression ()=default
 Destructor.
 
const std::vector< std::shared_ptr< const fem::Function< T > > > & coefficients () const
 Access coefficients.
 
std::vector< int > coefficient_offsets () const
 Offset for each coefficient expansion array on a cell. Used to pack data for multiple coefficients in a flat array. The last entry is the size required to store all coefficients.
 
void eval (const xtl::span< const std::int32_t > &active_cells, array2d< T > &values) const
 Evaluate the expression on cells. More...
 
const std::function< void(T *, const T *, const T *, const double *)> & get_tabulate_expression () const
 Get function for tabulate_expression. More...
 
const std::vector< std::shared_ptr< const fem::Constant< T > > > & constants () const
 Access constants. More...
 
std::shared_ptr< const mesh::Meshmesh () const
 Get mesh. More...
 
const array2d< double > & x () const
 Get evaluation points on reference cell. More...
 
const std::size_t value_size () const
 Get value size. More...
 
const std::size_t num_points () const
 Get number of points. More...
 

Detailed Description

template<typename T>
class dolfinx::fem::Expression< T >

Represents a mathematical expression evaluated at a pre-defined set of points on the reference cell. This class closely follows the concept of a UFC Expression.

This functionality can be used to evaluate a gradient of a Function at quadrature points in all cells. This evaluated gradient can then be used as input in to a non-FEniCS function that calculates a material constitutive model.

Constructor & Destructor Documentation

◆ Expression()

template<typename T >
dolfinx::fem::Expression< T >::Expression ( const std::vector< std::shared_ptr< const fem::Function< T >>> &  coefficients,
const std::vector< std::shared_ptr< const fem::Constant< T >>> &  constants,
const std::shared_ptr< const mesh::Mesh > &  mesh,
const array2d< double > &  X,
const std::function< void(T *, const T *, const T *, const double *)>  fn,
const std::size_t  value_size 
)
inline

Create Expression.

Parameters
[in]coefficientsCoefficients in the Expression
[in]constantsConstants in the Expression
[in]mesh
[in]Xpoints on reference cell, number of points rows and tdim cols
[in]fnfunction for tabulating expression
[in]value_sizesize of expression evaluated at single point

Member Function Documentation

◆ constants()

template<typename T >
const std::vector<std::shared_ptr<const fem::Constant<T> > >& dolfinx::fem::Expression< T >::constants ( ) const
inline

Access constants.

Returns
Vector of attached constants with their names. Names are used to set constants in user's c++ code. Index in the vector is the position of the constant in the original (nonsimplified) form.

◆ eval()

template<typename T >
void dolfinx::fem::Expression< T >::eval ( const xtl::span< const std::int32_t > &  active_cells,
array2d< T > &  values 
) const
inline

Evaluate the expression on cells.

Parameters
[in]active_cellsCells on which to evaluate the Expression
[out]valuesTo store the result. Caller responsible for correct sizing which should be num_cells rows by num_points*value_size columns.

◆ get_tabulate_expression()

template<typename T >
const std::function<void(T*, const T*, const T*, const double*)>& dolfinx::fem::Expression< T >::get_tabulate_expression ( ) const
inline

Get function for tabulate_expression.

Parameters
[out]fnFunction to tabulate expression.

◆ mesh()

template<typename T >
std::shared_ptr<const mesh::Mesh> dolfinx::fem::Expression< T >::mesh ( ) const
inline

Get mesh.

Returns
The mesh

◆ num_points()

template<typename T >
const std::size_t dolfinx::fem::Expression< T >::num_points ( ) const
inline

Get number of points.

Returns
number of points

◆ value_size()

template<typename T >
const std::size_t dolfinx::fem::Expression< T >::value_size ( ) const
inline

Get value size.

Returns
value_size

◆ x()

template<typename T >
const array2d<double>& dolfinx::fem::Expression< T >::x ( ) const
inline

Get evaluation points on reference cell.

Returns
Evaluation points

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