| 
    DOLFINx
    0.1.0
    
   DOLFINx C++ interface 
   | 
 
This class manages coordinate mappings for isoparametric cells. More...
#include <CoordinateElement.h>
Public Member Functions | |
| CoordinateElement (std::shared_ptr< basix::FiniteElement > element) | |
| Create a coordinate element from a Basix element.  More... | |
| CoordinateElement (mesh::CellType celltype, int degree) | |
| Create a Lagrage coordinate element.  More... | |
| virtual | ~CoordinateElement ()=default | 
| Destructor.  | |
| mesh::CellType | cell_shape () const | 
| Cell shape.  More... | |
| int | topological_dimension () const | 
| Return the topological dimension of the cell shape.  | |
| xt::xtensor< double, 4 > | tabulate (int n, const xt::xtensor< double, 2 > &X) const | 
| Compute basis values and derivatives at set of points.  More... | |
| void | compute_jacobian (const xt::xtensor< double, 4 > &dphi, const xt::xtensor< double, 2 > &cell_geometry, xt::xtensor< double, 3 > &J) const | 
| Compute Jacobian for a cell with given geometry using the basis functions and first order derivatives.  More... | |
| void | compute_jacobian_inverse (const xt::xtensor< double, 3 > &J, xt::xtensor< double, 3 > &K) const | 
| Compute the inverse of the Jacobian. If the coordinate element is affine, it computes the inverse at only one point.  More... | |
| void | compute_jacobian_determinant (const xt::xtensor< double, 3 > &J, xt::xtensor< double, 1 > &detJ) const | 
| Compute the determinant of the Jacobian. If the coordinate element is affine, it computes the determinant at only one point.  More... | |
| ElementDofLayout | dof_layout () const | 
| Return the dof layout.  | |
| void | pull_back (xt::xtensor< double, 2 > &X, xt::xtensor< double, 3 > &J, xt::xtensor< double, 1 > &detJ, xt::xtensor< double, 3 > &K, const xt::xtensor< double, 2 > &x, const xt::xtensor< double, 2 > &cell_geometry) const | 
| Compute reference coordinates X, and J, detJ and K for physical coordinates x.  | |
| void | permute_dofs (xtl::span< std::int32_t > dofs, const std::uint32_t cell_perm) const | 
| Permutes a list of DOF numbers on a cell.  | |
| void | unpermute_dofs (xtl::span< std::int32_t > dofs, const std::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)  | |
Static Public Member Functions | |
| static void | push_forward (xt::xtensor< double, 2 > &x, const xt::xtensor< double, 2 > &cell_geometry, const xt::xtensor< double, 2 > &phi) | 
| Compute physical coordinates x for points X in the reference configuration.  More... | |
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.
      
  | 
  explicit | 
Create a coordinate element from a Basix element.
| [in] | element | Element from Basix | 
| CoordinateElement::CoordinateElement | ( | mesh::CellType | celltype, | 
| int | degree | ||
| ) | 
Create a Lagrage coordinate element.
| [in] | celltype | The cell shape | 
| [in] | degree | Polynomial degree of the map | 
| mesh::CellType CoordinateElement::cell_shape | ( | ) | const | 
Cell shape.
| void CoordinateElement::compute_jacobian | ( | const xt::xtensor< double, 4 > & | dphi, | 
| const xt::xtensor< double, 2 > & | cell_geometry, | ||
| xt::xtensor< double, 3 > & | J | ||
| ) | const | 
Compute Jacobian for a cell with given geometry using the basis functions and first order derivatives.
| [in] | dphi | Pre-computed first order derivatives of basis functions at set of points. The shape of dphi is (tdim, number of points, number of basis fn, 1). | 
| [in] | cell_geometry | Coordinates/geometry | 
| [in,out] | J | The Jacobian The shape of J is (number of points, geometric dimension, topological dimenson). | 
| void CoordinateElement::compute_jacobian_determinant | ( | const xt::xtensor< double, 3 > & | J, | 
| xt::xtensor< double, 1 > & | detJ | ||
| ) | const | 
Compute the determinant of the Jacobian. If the coordinate element is affine, it computes the determinant at only one point.
| [in] | J | Polynomial degree of the map The shape of J is (number of points, geometric dimension, topological dimenson). | 
| [in,out] | detJ | Jacobian The shape of detJ is (number of points) | 
| void CoordinateElement::compute_jacobian_inverse | ( | const xt::xtensor< double, 3 > & | J, | 
| xt::xtensor< double, 3 > & | K | ||
| ) | const | 
Compute the inverse of the Jacobian. If the coordinate element is affine, it computes the inverse at only one point.
| [in] | J | The Jacobian The shape of J is (number of points, geometric dimension, topological dimenson). | 
| [in,out] | K | The inverse of the Jacobian The shape of J is (number of points, tpological dimension, geometrical dimenson). | 
      
  | 
  static | 
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 | 
| xt::xtensor< double, 4 > CoordinateElement::tabulate | ( | int | n, | 
| const xt::xtensor< double, 2 > & | X | ||
| ) | const | 
Compute basis values and derivatives at set of points.
| [in] | n | The order of derivatives, up to and including, to compute. Use 0 for the basis functions only. | 
| [in] | X | The points at which to compute the basis functions. The shape of x is (number of points, geometric dimension). | 
 1.8.17