9#include "ElementDofLayout.h" 
   12#include <basix/element-families.h> 
   13#include <basix/mdspan.hpp> 
   17#include <dolfinx/common/math.h> 
   18#include <dolfinx/mesh/cell_types.h> 
   25template <std::
floating_po
int T>
 
   36template <std::
floating_po
int T>
 
   51                    basix::element::lagrange_variant type
 
   52                    = basix::element::lagrange_variant::unset);
 
   74  basix::element::lagrange_variant 
variant() 
const;
 
   83                                            std::size_t num_points) 
const;
 
   93  void tabulate(
int nd, std::span<const T> X, std::array<std::size_t, 2> shape,
 
   94                std::span<T> basis) 
const;
 
  104  template <
typename U, 
typename V, 
typename W>
 
  107    math::dot(cell_geometry, dphi, J, 
true);
 
 
  113  template <
typename U, 
typename V>
 
  116    int gdim = J.extent(0);
 
  117    int tdim = K.extent(0);
 
 
  129  template <
typename U>
 
  133    static_assert(U::rank() == 2, 
"Must be rank 2");
 
  134    if (J.extent(0) == J.extent(1))
 
  138      assert(w.size() >= 2 * J.extent(0) * J.extent(1));
 
  139      using X = 
typename U::element_type;
 
  140      using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
 
  141          X, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
 
  142      mdspan2_t B(w.data(), J.extent(1), J.extent(0));
 
  143      mdspan2_t BA(w.data() + J.extent(0) * J.extent(1), B.extent(0),
 
  145      for (std::size_t i = 0; i < B.extent(0); ++i)
 
  146        for (std::size_t j = 0; j < B.extent(1); ++j)
 
  150      std::fill_n(BA.data_handle(), BA.size(), 0);
 
  152      return std::sqrt(math::det(BA));
 
 
  166  template <
typename U, 
typename V, 
typename W>
 
  169    for (std::size_t i = 0; i < x.extent(0); ++i)
 
  170      for (std::size_t j = 0; j < x.extent(1); ++j)
 
  174    math::dot(phi, cell_geometry, x);
 
 
  187  template <
typename U, 
typename V, 
typename W>
 
  191    assert(X.extent(0) == x.extent(0));
 
  192    assert(X.extent(1) == K.extent(0));
 
  193    assert(x.extent(1) == K.extent(1));
 
  194    for (std::size_t i = 0; i < X.extent(0); ++i)
 
  195      for (std::size_t j = 0; j < X.extent(1); ++j)
 
  199    for (std::size_t p = 0; p < x.extent(0); ++p)
 
  200      for (std::size_t i = 0; i < K.extent(0); ++i)
 
  201        for (std::size_t j = 0; j < K.extent(1); ++j)
 
  202          X(p, i) += K(i, j) * (x(p, j) - x0[j]);
 
 
  206  template <
typename X>
 
  207  using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
 
  208      X, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
 
  223                           double tol = 1.0e-6, 
int maxit = 15) 
const;
 
  226  void permute(std::span<std::int32_t> dofs, std::uint32_t cell_perm) 
const;
 
  229  void permute_inv(std::span<std::int32_t> dofs, std::uint32_t cell_perm) 
const;
 
  247  std::shared_ptr<const basix::FiniteElement<T>> _element;
 
 
Definition CoordinateElement.h:26
 
Definition CoordinateElement.h:38
 
ElementDofLayout create_dof_layout() const
Compute and return the dof layout.
Definition CoordinateElement.cpp:63
 
static void compute_jacobian(const U &dphi, const V &cell_geometry, W &&J)
Definition CoordinateElement.h:105
 
static void pull_back_affine(U &&X, const V &K, const std::array< T, 3 > &x0, const W &x)
Compute reference coordinates X for physical coordinates x for an affine map. For the affine case,...
Definition CoordinateElement.h:188
 
void tabulate(int nd, std::span< const T > X, std::array< std::size_t, 2 > shape, std::span< T > basis) const
Evaluate basis values and derivatives at set of points.
Definition CoordinateElement.cpp:54
 
CoordinateElement(std::shared_ptr< const basix::FiniteElement< T > > element)
Create a coordinate element from a Basix element.
Definition CoordinateElement.cpp:18
 
void permute_inv(std::span< std::int32_t > dofs, std::uint32_t cell_perm) const
Reverses a DOF permutation.
Definition CoordinateElement.cpp:171
 
basix::element::lagrange_variant variant() const
Variant of the element.
Definition CoordinateElement.cpp:201
 
mesh::CellType cell_shape() const
Cell shape.
Definition CoordinateElement.cpp:40
 
static void compute_jacobian_inverse(const U &J, V &&K)
Compute the inverse of the Jacobian.
Definition CoordinateElement.h:114
 
MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan< X, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents< std::size_t, 2 > > mdspan2_t
mdspan typedef
Definition CoordinateElement.h:207
 
std::array< std::size_t, 4 > tabulate_shape(std::size_t nd, std::size_t num_points) const
Shape of array to fill when calling tabulate.
Definition CoordinateElement.cpp:47
 
int degree() const
The polynomial degree of the element.
Definition CoordinateElement.cpp:187
 
void permute(std::span< std::int32_t > dofs, std::uint32_t cell_perm) const
Permute a list of DOF numbers on a cell.
Definition CoordinateElement.cpp:163
 
int dim() const
The dimension of the coordinate element space.
Definition CoordinateElement.cpp:194
 
void pull_back_nonaffine(mdspan2_t< T > X, mdspan2_t< const T > x, mdspan2_t< const T > cell_geometry, double tol=1.0e-6, int maxit=15) const
Compute reference coordinates X for physical coordinates x for a non-affine map.
Definition CoordinateElement.cpp:71
 
virtual ~CoordinateElement()=default
Destructor.
 
bool is_affine() const noexcept
Definition CoordinateElement.h:240
 
static void push_forward(U &&x, const V &cell_geometry, const W &phi)
Compute physical coordinates x for points X in the reference configuration.
Definition CoordinateElement.h:167
 
static double compute_jacobian_determinant(const U &J, std::span< typename U::value_type > w)
Compute the determinant of the Jacobian.
Definition CoordinateElement.h:131
 
bool needs_dof_permutations() const
Indicates whether the geometry DOF numbers on each cell need permuting.
Definition CoordinateElement.cpp:179
 
Definition ElementDofLayout.h:30
 
Finite element method functionality.
Definition assemble_matrix_impl.h:26
 
CellType
Cell type identifier.
Definition cell_types.h:22