Represents a mathematical expression evaluated at a pre-defined set of points on the reference cell.
More...
|
| Expression (const std::vector< std::shared_ptr< const Function< scalar_type, geometry_type > > > &coefficients, const std::vector< std::shared_ptr< const Constant< scalar_type > > > &constants, std::span< const geometry_type > X, std::array< std::size_t, 2 > Xshape, std::function< void(scalar_type *, const scalar_type *, const scalar_type *, const geometry_type *, const int *, const uint8_t *)> fn, const std::vector< int > &value_shape, std::shared_ptr< const FunctionSpace< geometry_type > > argument_function_space=nullptr) |
| Create an Expression.
|
|
| Expression (Expression &&e)=default |
| Move constructor.
|
|
virtual | ~Expression ()=default |
| Destructor.
|
|
std::shared_ptr< const FunctionSpace< geometry_type > > | argument_function_space () const |
| Get argument function space.
|
|
const std::vector< std::shared_ptr< const Function< scalar_type, geometry_type > > > & | coefficients () const |
| Get coefficients.
|
|
const std::vector< std::shared_ptr< const Constant< scalar_type > > > & | constants () const |
| Get constants.
|
|
std::vector< int > | coefficient_offsets () const |
| Offset for each coefficient expansion array on a cell.
|
|
void | eval (const mesh::Mesh< geometry_type > &mesh, std::span< const std::int32_t > cells, std::span< scalar_type > values, std::array< std::size_t, 2 > vshape) const |
| Evaluate Expression on cells.
|
|
const std::function< void(scalar_type *, const scalar_type *, const scalar_type *, const geometry_type *, const int *, const uint8_t *)> & | get_tabulate_expression () const |
| Get function for tabulate_expression.
|
|
int | value_size () const |
| Get value size.
|
|
const std::vector< int > & | value_shape () const |
| Get value shape.
|
|
std::pair< std::vector< geometry_type >, std::array< std::size_t, 2 > > | X () const |
| Evaluation points on the reference cell.
|
|
template<
dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_type_t<T>>
class dolfinx::fem::Expression< T, U >
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.
The 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.
- Template Parameters
-
T | The scalar type |
U | The mesh geometry scalar type |