Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/v0.3.0/v0.9.0/cpp
DOLFINx  0.3.0
DOLFINx C++ interface
CoordinateElement.h
1 // Copyright (C) 2018-2020 Garth N. Wells and Chris N. Richardson
2 //
3 // This file is part of DOLFINx (https://www.fenicsproject.org)
4 //
5 // SPDX-License-Identifier: LGPL-3.0-or-later
6 
7 #pragma once
8 
9 #include "ElementDofLayout.h"
10 #include <cstdint>
11 #include <dolfinx/mesh/cell_types.h>
12 #include <functional>
13 #include <memory>
14 #include <xtensor/xtensor.hpp>
15 #include <xtl/xspan.hpp>
16 
17 namespace basix
18 {
19 class FiniteElement;
20 }
21 
22 namespace dolfinx::fem
23 {
24 
25 // FIXME: A dof layout on a reference cell needs to be defined.
27 
29 {
30 public:
33  explicit CoordinateElement(std::shared_ptr<basix::FiniteElement> element);
34 
38  CoordinateElement(mesh::CellType celltype, int degree);
39 
41  virtual ~CoordinateElement() = default;
42 
45  mesh::CellType cell_shape() const;
46 
48  int topological_dimension() const;
49 
57  xt::xtensor<double, 4> tabulate(int n, const xt::xtensor<double, 2>& X) const;
58 
68  void compute_jacobian(const xt::xtensor<double, 4>& dphi,
69  const xt::xtensor<double, 2>& cell_geometry,
70  xt::xtensor<double, 3>& J) const;
71 
79  void compute_jacobian_inverse(const xt::xtensor<double, 3>& J,
80  xt::xtensor<double, 3>& K) const;
81 
88  void compute_jacobian_determinant(const xt::xtensor<double, 3>& J,
89  xt::xtensor<double, 1>& detJ) const;
90 
93 
99  static void push_forward(xt::xtensor<double, 2>& x,
100  const xt::xtensor<double, 2>& cell_geometry,
101  const xt::xtensor<double, 2>& phi);
102 
105  void pull_back(xt::xtensor<double, 2>& X, xt::xtensor<double, 3>& J,
106  xt::xtensor<double, 1>& detJ, xt::xtensor<double, 3>& K,
107  const xt::xtensor<double, 2>& x,
108  const xt::xtensor<double, 2>& cell_geometry) const;
109 
111  void permute_dofs(xtl::span<std::int32_t> dofs,
112  const std::uint32_t cell_perm) const;
113 
115  void unpermute_dofs(xtl::span<std::int32_t> dofs,
116  const std::uint32_t cell_perm) const;
117 
123  bool needs_dof_permutations() const;
124 
126  double non_affine_atol = 1.0e-8;
127 
130 
131 private:
132  // Flag denoting affine map
133  bool _is_affine;
134 
135  // Basix Element
136  std::shared_ptr<basix::FiniteElement> _element;
137 };
138 } // namespace dolfinx::fem
This class manages coordinate mappings for isoparametric cells.
Definition: CoordinateElement.h:29
void unpermute_dofs(xtl::span< std::int32_t > dofs, const std::uint32_t cell_perm) const
Reverses a DOF permutation.
Definition: CoordinateElement.cpp:302
virtual ~CoordinateElement()=default
Destructor.
int topological_dimension() const
Return the topological dimension of the cell shape.
Definition: CoordinateElement.cpp:73
double non_affine_atol
Absolute increment stopping criterium for non-affine Newton solver.
Definition: CoordinateElement.h:126
bool needs_dof_permutations() const
Indicates whether the geometry DOF numbers on each cell need permuting.
Definition: CoordinateElement.cpp:309
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.
Definition: CoordinateElement.cpp:208
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...
Definition: CoordinateElement.cpp:84
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.
Definition: CoordinateElement.cpp:295
int non_affine_max_its
Maximum number of iterations for non-affine Newton solver.
Definition: CoordinateElement.h:129
ElementDofLayout dof_layout() const
Return the dof layout.
Definition: CoordinateElement.cpp:183
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.
Definition: CoordinateElement.cpp:193
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 ...
Definition: CoordinateElement.cpp:125
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 determi...
Definition: CoordinateElement.cpp:162
xt::xtensor< double, 4 > tabulate(int n, const xt::xtensor< double, 2 > &X) const
Compute basis values and derivatives at set of points.
Definition: CoordinateElement.cpp:79
CoordinateElement(std::shared_ptr< basix::FiniteElement > element)
Create a coordinate element from a Basix element.
Definition: CoordinateElement.cpp:36
mesh::CellType cell_shape() const
Cell shape.
Definition: CoordinateElement.cpp:55
The class represents the degree-of-freedom (dofs) for an element. Dofs are associated with a mesh ent...
Definition: ElementDofLayout.h:34
Finite element method functionality.
Definition: assemble_matrix_impl.h:23
CellType
Cell type identifier.
Definition: cell_types.h:22