Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/v0.1.0/v0.9.0/cpp
DOLFINx  0.1.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/common/array2d.h>
12 #include <dolfinx/mesh/cell_types.h>
13 #include <functional>
14 #include <memory>
15 #include <string>
16 #include <xtensor/xtensor.hpp>
17 #include <xtl/xspan.hpp>
18 
19 namespace basix
20 {
21 class FiniteElement;
22 }
23 
24 namespace dolfinx::fem
25 {
26 
27 // FIXME: A dof layout on a reference cell needs to be defined.
29 
31 {
32 public:
35  explicit CoordinateElement(std::shared_ptr<basix::FiniteElement> element);
36 
40  CoordinateElement(mesh::CellType celltype, int degree);
41 
43  virtual ~CoordinateElement() = default;
44 
47  mesh::CellType cell_shape() const;
48 
50  int topological_dimension() const;
51 
59  xt::xtensor<double, 4> tabulate(int n, const xt::xtensor<double, 2>& X) const;
60 
70  void compute_jacobian(const xt::xtensor<double, 4>& dphi,
71  const xt::xtensor<double, 2>& cell_geometry,
72  xt::xtensor<double, 3>& J) const;
73 
82  void compute_jacobian_inverse(const xt::xtensor<double, 3>& J,
83  xt::xtensor<double, 3>& K) const;
84 
92  void compute_jacobian_determinant(const xt::xtensor<double, 3>& J,
93  xt::xtensor<double, 1>& detJ) const;
94 
97 
103  static void push_forward(xt::xtensor<double, 2>& x,
104  const xt::xtensor<double, 2>& cell_geometry,
105  const xt::xtensor<double, 2>& phi);
106 
109  void pull_back(xt::xtensor<double, 2>& X, xt::xtensor<double, 3>& J,
110  xt::xtensor<double, 1>& detJ, xt::xtensor<double, 3>& K,
111  const xt::xtensor<double, 2>& x,
112  const xt::xtensor<double, 2>& cell_geometry) const;
113 
115  void permute_dofs(xtl::span<std::int32_t> dofs,
116  const std::uint32_t cell_perm) const;
117 
119  void unpermute_dofs(xtl::span<std::int32_t> dofs,
120  const std::uint32_t cell_perm) const;
121 
124  bool needs_permutation_data() const;
125 
127  double non_affine_atol = 1.0e-8;
128 
131 
132 private:
133  // Flag denoting affine map
134  bool _is_affine;
135 
136  // Basix Element
137  std::shared_ptr<basix::FiniteElement> _element;
138 };
139 } // namespace dolfinx::fem
dolfinx::fem::CoordinateElement::dof_layout
ElementDofLayout dof_layout() const
Return the dof layout.
Definition: CoordinateElement.cpp:180
dolfinx::fem::CoordinateElement::compute_jacobian
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:81
dolfinx::mesh::CellType
CellType
Cell type identifier.
Definition: cell_types.h:21
dolfinx::fem::CoordinateElement::compute_jacobian_determinant
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:159
dolfinx::fem::CoordinateElement::compute_jacobian_inverse
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:122
dolfinx::fem::CoordinateElement::push_forward
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:197
dolfinx::fem::CoordinateElement::non_affine_max_its
int non_affine_max_its
Maximum number of iterations for non-affine Newton solver.
Definition: CoordinateElement.h:130
dolfinx::fem::CoordinateElement::non_affine_atol
double non_affine_atol
Absolute increment stopping criterium for non-affine Newton solver.
Definition: CoordinateElement.h:127
dolfinx::fem::CoordinateElement::~CoordinateElement
virtual ~CoordinateElement()=default
Destructor.
dolfinx::fem::CoordinateElement::CoordinateElement
CoordinateElement(std::shared_ptr< basix::FiniteElement > element)
Create a coordinate element from a Basix element.
Definition: CoordinateElement.cpp:36
dolfinx::fem::CoordinateElement::pull_back
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:212
dolfinx::fem::CoordinateElement::tabulate
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:76
dolfinx::fem::ElementDofLayout
The class represents the degree-of-freedom (dofs) for an element. Dofs are associated with a mesh ent...
Definition: ElementDofLayout.h:35
dolfinx::fem::CoordinateElement::unpermute_dofs
void unpermute_dofs(xtl::span< std::int32_t > dofs, const std::uint32_t cell_perm) const
Reverses a DOF permutation.
Definition: CoordinateElement.cpp:306
dolfinx::fem
Finite element method functionality.
Definition: assemble_matrix_impl.h:22
dolfinx::fem::CoordinateElement::needs_permutation_data
bool needs_permutation_data() const
Indicates whether the coordinate map needs permutation data passing in (for higher order geometries)
Definition: CoordinateElement.cpp:313
dolfinx::fem::CoordinateElement::topological_dimension
int topological_dimension() const
Return the topological dimension of the cell shape.
Definition: CoordinateElement.cpp:70
dolfinx::fem::CoordinateElement::permute_dofs
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:299
dolfinx::fem::CoordinateElement::cell_shape
mesh::CellType cell_shape() const
Cell shape.
Definition: CoordinateElement.cpp:52
dolfinx::fem::CoordinateElement
This class manages coordinate mappings for isoparametric cells.
Definition: CoordinateElement.h:30