DOLFINx
0.4.1
DOLFINx C++ interface
|
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 Function< T >>> &coefficients, const std::vector< std::shared_ptr< const Constant< T >>> &constants, const xt::xtensor< double, 2 > &X, const std::function< void(T *, const T *, const T *, const double *, const int *, const uint8_t *)> fn, const std::vector< int > &value_shape, const std::shared_ptr< const mesh::Mesh > &mesh=nullptr, const std::shared_ptr< const FunctionSpace > argument_function_space=nullptr) | |
Create an Expression. More... | |
Expression (Expression &&form)=default | |
Move constructor. | |
virtual | ~Expression ()=default |
Destructor. | |
std::shared_ptr< const fem::FunctionSpace > | argument_function_space () const |
Get argument function space. More... | |
const std::vector< std::shared_ptr< const fem::Function< T > > > & | coefficients () const |
Get coefficients. More... | |
const std::vector< std::shared_ptr< const fem::Constant< T > > > & | constants () const |
Get constants. More... | |
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. | |
template<typename U > | |
void | eval (const xtl::span< const std::int32_t > &cells, U &values) const |
Evaluate the expression on cells. More... | |
const std::function< void(T *, const T *, const T *, const double *, const int *, const uint8_t *)> & | get_tabulate_expression () const |
Get function for tabulate_expression. More... | |
std::shared_ptr< const mesh::Mesh > | mesh () const |
Get mesh. More... | |
int | value_size () const |
Get value size. More... | |
const std::vector< int > & | value_shape () const |
Get value shape. More... | |
const xt::xtensor< double, 2 > & | X () const |
Get evaluation points on reference cell. More... | |
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.
|
inline |
Create an Expression.
Users should prefer the create_expression factory functions
[in] | coefficients | Coefficients in the Expression |
[in] | constants | Constants in the Expression |
[in] | mesh | |
[in] | X | points on reference cell, number of points rows and tdim cols |
[in] | fn | function for tabulating expression |
[in] | value_shape | shape of expression evaluated at single point |
[in] | argument_function_space | Function space for Argument |
|
inline |
Get argument function space.
|
inline |
Get coefficients.
|
inline |
Get constants.
|
inline |
Evaluate the expression on cells.
[in] | cells | Cells on which to evaluate the Expression |
[out] | values | A 2D array to store the result. Caller is responsible for correct sizing which should be (num_cells, num_points * value_size * num_all_argument_dofs columns). |
|
inline |
Get function for tabulate_expression.
|
inline |
Get mesh.
|
inline |
Get value shape.
|
inline |
Get value size.
|
inline |
Get evaluation points on reference cell.