DOLFIN-X
DOLFIN-X C++ interface
FiniteElement.h
1 // Copyright (C) 2008-2020 Anders Logg and Garth N. Wells
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 <dolfinx/common/types.h>
10 #include <dolfinx/mesh/cell_types.h>
11 #include <functional>
12 #include <memory>
13 #include <petscsys.h>
14 #include <unsupported/Eigen/CXX11/Tensor>
15 #include <vector>
16 
17 struct ufc_coordinate_mapping;
18 struct ufc_finite_element;
19 
20 namespace dolfinx::fem
21 {
22 
26 {
27 public:
30  FiniteElement(const ufc_finite_element& ufc_element);
31 
33  virtual ~FiniteElement() = default;
34 
37  std::string signature() const;
38 
41  mesh::CellType cell_shape() const;
42 
45  int space_dimension() const;
46 
49  int value_size() const;
50 
54  int reference_value_size() const;
55 
58  int value_rank() const;
59 
61  int value_dimension(int i) const;
62 
63  // FIXME: Is this well-defined? What does it do on non-simplex
64  // elements?
66  int degree() const;
67 
70  std::string family() const;
71 
73  // reference_values[num_points][num_dofs][reference_value_size]
75  Eigen::Tensor<double, 3, Eigen::RowMajor>& reference_values,
76  const Eigen::Ref<const Eigen::Array<
77  double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>>& X) const;
78 
81  Eigen::Tensor<double, 3, Eigen::RowMajor>& values,
82  const Eigen::Tensor<double, 3, Eigen::RowMajor>& reference_values,
83  const Eigen::Ref<const Eigen::Array<double, Eigen::Dynamic,
84  Eigen::Dynamic, Eigen::RowMajor>>& X,
85  const Eigen::Tensor<double, 3, Eigen::RowMajor>& J,
86  const Eigen::Ref<const Eigen::Array<double, Eigen::Dynamic, 1>>& detJ,
87  const Eigen::Tensor<double, 3, Eigen::RowMajor>& K,
88  const std::uint32_t permutation_info) const;
89 
92  Eigen::Tensor<double, 4, Eigen::RowMajor>& values, std::size_t order,
93  const Eigen::Tensor<double, 4, Eigen::RowMajor>& reference_values,
94  const Eigen::Ref<const Eigen::Array<double, Eigen::Dynamic,
95  Eigen::Dynamic, Eigen::RowMajor>>& X,
96  const Eigen::Tensor<double, 3, Eigen::RowMajor>& J,
97  const Eigen::Ref<const Eigen::Array<double, Eigen::Dynamic, 1>>& detJ,
98  const Eigen::Tensor<double, 3, Eigen::RowMajor>& K,
99  const std::uint32_t permutation_info) const;
100 
103  bool has_dof_reference_coordinates() const noexcept;
104 
107  const Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>&
109 
112  void transform_values(
113  PetscScalar* reference_values,
114  const Eigen::Ref<const Eigen::Array<PetscScalar, Eigen::Dynamic,
115  Eigen::Dynamic, Eigen::RowMajor>>&
116  physical_values,
117  const Eigen::Ref<const Eigen::Array<double, Eigen::Dynamic,
118  Eigen::Dynamic, Eigen::RowMajor>>&
119  coordinate_dofs) const;
120 
122  int num_sub_elements() const;
123 
125  std::size_t hash() const;
126 
128  std::shared_ptr<const FiniteElement>
129  extract_sub_element(const std::vector<int>& component) const;
130 
131 private:
132  std::string _signature, _family;
133 
134  mesh::CellType _cell_shape;
135 
136  int _tdim, _space_dim, _value_size, _reference_value_size, _degree;
137 
138  // Dof coordinates on the reference element
139  bool _has_refX;
140  Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> _refX;
141 
142  // List of sub-elements (if any)
143  std::vector<std::shared_ptr<const FiniteElement>> _sub_elements;
144 
145  // Recursively extract sub finite element
146  static std::shared_ptr<const FiniteElement>
147  extract_sub_element(const FiniteElement& finite_element,
148  const std::vector<int>& component);
149 
150  // Simple hash of the signature string
151  std::size_t _hash;
152 
153  // Dimension of each value space
154  std::vector<int> _value_dimension;
155 
156  // Functions for basis and derivatives evaluation
157  std::function<int(double*, int, const double*)> _evaluate_reference_basis;
158 
159  std::function<int(double*, int, int, const double*)>
160  _evaluate_reference_basis_derivatives;
161 
162  std::function<int(double*, int, int, const double*, const double*,
163  const double*, const double*, const double*,
164  const std::uint32_t)>
165  _transform_reference_basis_derivatives;
166 
167  std::function<int(ufc_scalar_t*, const ufc_scalar_t*, const double*,
168  const ufc_coordinate_mapping*)>
169  _transform_values;
170 };
171 } // namespace dolfinx::fem
dolfinx::fem::FiniteElement::FiniteElement
FiniteElement(const ufc_finite_element &ufc_element)
Create finite element from UFC finite element.
Definition: FiniteElement.cpp:19
dolfinx::fem::FiniteElement::extract_sub_element
std::shared_ptr< const FiniteElement > extract_sub_element(const std::vector< int > &component) const
Extract sub finite element for component.
Definition: FiniteElement.cpp:207
dolfinx::fem::FiniteElement::value_dimension
int value_dimension(int i) const
Return the dimension of the value space for axis i.
Definition: FiniteElement.cpp:99
dolfinx::fem::FiniteElement::evaluate_reference_basis
void evaluate_reference_basis(Eigen::Tensor< double, 3, Eigen::RowMajor > &reference_values, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &X) const
Evaluate all basis functions at given point in reference cell.
Definition: FiniteElement.cpp:110
dolfinx::fem::FiniteElement::num_sub_elements
int num_sub_elements() const
Return the number of sub elements (for a mixed element)
Definition: FiniteElement.cpp:202
dolfinx::fem::FiniteElement::dof_reference_coordinates
const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > & dof_reference_coordinates() const
Tabulate the reference coordinates of all dofs on an element.
Definition: FiniteElement.cpp:176
dolfinx::mesh::CellType
CellType
Cell type identifier.
Definition: cell_types.h:22
dolfinx::fem::FiniteElement::cell_shape
mesh::CellType cell_shape() const
Cell shape.
Definition: FiniteElement.cpp:84
dolfinx::fem::FiniteElement::transform_values
void transform_values(PetscScalar *reference_values, const Eigen::Ref< const Eigen::Array< PetscScalar, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &physical_values, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &coordinate_dofs) const
Map values of field from physical to reference space which has been evaluated at points given by dof_...
Definition: FiniteElement.cpp:187
dolfinx::fem::FiniteElement
Finite Element, containing the dof layout on a reference element, and various methods for evaluating ...
Definition: FiniteElement.h:25
dolfinx::fem::FiniteElement::transform_reference_basis
void transform_reference_basis(Eigen::Tensor< double, 3, Eigen::RowMajor > &values, const Eigen::Tensor< double, 3, Eigen::RowMajor > &reference_values, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &X, const Eigen::Tensor< double, 3, Eigen::RowMajor > &J, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, 1 >> &detJ, const Eigen::Tensor< double, 3, Eigen::RowMajor > &K, const std::uint32_t permutation_info) const
Push basis functions forward to physical element.
Definition: FiniteElement.cpp:126
dolfinx::fem::FiniteElement::value_rank
int value_rank() const
Rank of the value space.
Definition: FiniteElement.cpp:97
dolfinx::fem::FiniteElement::hash
std::size_t hash() const
Return simple hash of the signature string.
Definition: FiniteElement.cpp:204
dolfinx::fem::FiniteElement::~FiniteElement
virtual ~FiniteElement()=default
Destructor.
dolfinx::fem::FiniteElement::transform_reference_basis_derivatives
void transform_reference_basis_derivatives(Eigen::Tensor< double, 4, Eigen::RowMajor > &values, std::size_t order, const Eigen::Tensor< double, 4, Eigen::RowMajor > &reference_values, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >> &X, const Eigen::Tensor< double, 3, Eigen::RowMajor > &J, const Eigen::Ref< const Eigen::Array< double, Eigen::Dynamic, 1 >> &detJ, const Eigen::Tensor< double, 3, Eigen::RowMajor > &K, const std::uint32_t permutation_info) const
Push basis function (derivatives) forward to physical element.
Definition: FiniteElement.cpp:148
dolfinx::fem::FiniteElement::space_dimension
int space_dimension() const
Dimension of the finite element function space.
Definition: FiniteElement.cpp:88
dolfinx::fem::FiniteElement::signature
std::string signature() const
String identifying the finite element.
Definition: FiniteElement.cpp:82
dolfinx::fem::FiniteElement::family
std::string family() const
The finite element family.
Definition: FiniteElement.cpp:108
dolfinx::fem
Finite element method functionality.
Definition: assemble_matrix_impl.h:34
dolfinx::fem::FiniteElement::value_size
int value_size() const
The value size, e.g. 1 for a scalar function, 2 for a 2D vector.
Definition: FiniteElement.cpp:90
dolfinx::fem::FiniteElement::degree
int degree() const
Return the maximum polynomial degree.
Definition: FiniteElement.cpp:106
dolfinx::fem::FiniteElement::reference_value_size
int reference_value_size() const
The value size, e.g. 1 for a scalar function, 2 for a 2D vector for the reference element.
Definition: FiniteElement.cpp:92
dolfinx::fem::FiniteElement::has_dof_reference_coordinates
bool has_dof_reference_coordinates() const noexcept
Check if reference coordinates for dofs are defined.
Definition: FiniteElement.cpp:170