DOLFIN-X
DOLFIN-X C++ interface
|
This class manages coordinate mappings for isoparametric cells. More...
#include <CoordinateElement.h>
Public Member Functions | |
CoordinateElement (std::shared_ptr< basix::FiniteElement > element, int geometric_dimension, const std::string &signature, const ElementDofLayout &dof_layout) | |
Create a coordinate element. More... | |
virtual | ~CoordinateElement ()=default |
Destructor. | |
std::string | signature () const |
String identifying the finite element. More... | |
mesh::CellType | cell_shape () const |
Cell shape. More... | |
int | topological_dimension () const |
Return the topological dimension of the cell shape. | |
int | geometric_dimension () const |
Return the geometric dimension of the cell shape. | |
xt::xtensor< double, 4 > | tabulate (int n, const array2d< double > &X) const |
Tabulate shape functions up to n-th order derivative at points X in the reference geometry Note: Dynamic allocation. | |
void | compute_jacobian_data (const xt::xtensor< double, 4 > &tabulated_data, const array2d< double > &X, const array2d< double > &cell_geometry, std::vector< double > &J, tcb::span< double > detJ, std::vector< double > &K) const |
Compute J, K and detJ for a cell with given geometry, and the basis functions and first order derivatives at points X. | |
const ElementDofLayout & | dof_layout () const |
Return the dof layout. | |
void | push_forward (array2d< double > &x, const array2d< double > &cell_geometry, const xt::xtensor< double, 2 > &phi) const |
Compute physical coordinates x for points X in the reference configuration. More... | |
void | compute_reference_geometry (array2d< double > &X, std::vector< double > &J, tcb::span< double > detJ, std::vector< double > &K, const array2d< double > &x, const array2d< double > &cell_geometry) const |
Compute reference coordinates X, and J, detJ and K for physical coordinates x. | |
void | permute_dofs (std::int32_t *dofs, const uint32_t cell_perm) const |
Permutes a list of DOF numbers on a cell. | |
void | unpermute_dofs (std::int32_t *dofs, const uint32_t cell_perm) const |
Reverses a DOF permutation. | |
bool | needs_permutation_data () const |
Indicates whether the coordinate map needs permutation data passing in (for higher order geometries) | |
Public Attributes | |
double | non_affine_atol = 1.0e-8 |
Absolute increment stopping criterium for non-affine Newton solver. | |
int | non_affine_max_its = 10 |
Maximum number of iterations for non-affine Newton solver. | |
This class manages coordinate mappings for isoparametric cells.
CoordinateElement::CoordinateElement | ( | std::shared_ptr< basix::FiniteElement > | element, |
int | geometric_dimension, | ||
const std::string & | signature, | ||
const ElementDofLayout & | dof_layout | ||
) |
Create a coordinate element.
[in] | element | Element from basix |
[in] | geometric_dimension | Geometric dimension |
[in] | signature | Signature string description of coordinate map |
[in] | dof_layout | Layout of the geometry degrees-of-freedom |
mesh::CellType CoordinateElement::cell_shape | ( | ) | const |
Cell shape.
void CoordinateElement::push_forward | ( | array2d< double > & | x, |
const array2d< double > & | cell_geometry, | ||
const xt::xtensor< double, 2 > & | phi | ||
) | const |
Compute physical coordinates x for points X in the reference configuration.
[in,out] | x | The physical coordinates of the reference points X |
[in] | cell_geometry | The cell node coordinates (physical) |
[in] | phi | Tabulated basis functions at reference points X |
std::string dolfinx::fem::CoordinateElement::signature | ( | ) | const |
String identifying the finite element.