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
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 <basix/finite-element.h>
10 #include <dolfinx/common/types.h>
11 #include <dolfinx/mesh/cell_types.h>
12 #include <functional>
13 #include <memory>
14 #include <numeric>
15 #include <vector>
16 #include <xtl/xspan.hpp>
17 
18 struct ufc_finite_element;
19 
20 namespace dolfinx::fem
21 {
25 {
26 public:
29  explicit FiniteElement(const ufc_finite_element& ufc_element);
30 
32  FiniteElement(const FiniteElement& element) = delete;
33 
35  FiniteElement(FiniteElement&& element) = default;
36 
38  virtual ~FiniteElement() = default;
39 
41  FiniteElement& operator=(const FiniteElement& element) = delete;
42 
44  FiniteElement& operator=(FiniteElement&& element) = default;
45 
48  std::string signature() const noexcept;
49 
52  mesh::CellType cell_shape() const noexcept;
53 
56  int space_dimension() const noexcept;
57 
62  int block_size() const noexcept;
63 
66  int value_size() const noexcept;
67 
71  int reference_value_size() const noexcept;
72 
75  int value_rank() const noexcept;
76 
78  int value_dimension(int i) const;
79 
82  std::string family() const noexcept;
83 
93  void tabulate(xt::xtensor<double, 4>& values, const xt::xtensor<double, 2>& X,
94  int order) const;
95 
97  void transform_reference_basis(xt::xtensor<double, 3>& values,
98  const xt::xtensor<double, 3>& reference_values,
99  const xt::xtensor<double, 3>& J,
100  const xtl::span<const double>& detJ,
101  const xt::xtensor<double, 3>& K) const;
102 
105  int num_sub_elements() const noexcept;
106 
108  const std::vector<std::shared_ptr<const FiniteElement>>&
109  sub_elements() const noexcept;
110 
112  std::size_t hash() const noexcept;
113 
115  std::shared_ptr<const FiniteElement>
116  extract_sub_element(const std::vector<int>& component) const;
117 
124  bool interpolation_ident() const noexcept;
125 
132  const xt::xtensor<double, 2>& interpolation_points() const;
133 
151  template <typename T>
152  constexpr void interpolate(const xt::xtensor<T, 2>& values,
153  xtl::span<T> dofs) const
154  {
155  if (!_element)
156  {
157  throw std::runtime_error("No underlying element for interpolation. "
158  "Cannot interpolate mixed elements directly.");
159  }
160 
161  _element->interpolate(tcb::make_span(dofs), tcb::make_span(values), _bs);
162  }
163 
177  bool needs_dof_transformations() const noexcept;
178 
192  bool needs_dof_permutations() const noexcept;
193 
208  template <typename T>
209  std::function<void(const xtl::span<T>&, const xtl::span<const std::uint32_t>&,
210  std::int32_t, int)>
211  get_dof_transformation_function(bool inverse = false, bool transpose = false,
212  bool scalar_element = false) const
213  {
215  {
216  // If no permutation needed, return function that does nothing
217  return [](const xtl::span<T>&, const xtl::span<const std::uint32_t>&,
218  std::int32_t, int)
219  {
220  // Do nothing
221  };
222  }
223 
224  if (_sub_elements.size() != 0)
225  {
226  if (_bs == 1)
227  {
228  // Mixed element
229  std::vector<std::function<void(const xtl::span<T>&,
230  const xtl::span<const std::uint32_t>&,
231  std::int32_t, int)>>
232  sub_element_functions;
233  std::vector<int> dims;
234  for (std::size_t i = 0; i < _sub_elements.size(); ++i)
235  {
236  sub_element_functions.push_back(
237  _sub_elements[i]->get_dof_transformation_function<T>(inverse,
238  transpose));
239  dims.push_back(_sub_elements[i]->space_dimension());
240  }
241 
242  return [dims, sub_element_functions](
243  const xtl::span<T>& data,
244  const xtl::span<const std::uint32_t>& cell_info,
245  std::int32_t cell, int block_size)
246  {
247  std::size_t start = 0;
248  for (std::size_t e = 0; e < sub_element_functions.size(); ++e)
249  {
250  const std::size_t width = dims[e] * block_size;
251  sub_element_functions[e](data.subspan(start, width), cell_info,
252  cell, block_size);
253  start += width;
254  }
255  };
256  }
257  else if (!scalar_element)
258  {
259  // Vector element
260  const std::function<void(const xtl::span<T>&,
261  const xtl::span<const std::uint32_t>&,
262  std::int32_t, int)>
263  sub_function = _sub_elements[0]->get_dof_transformation_function<T>(
264  inverse, transpose);
265  const int ebs = _bs;
266  return
267  [ebs, sub_function](const xtl::span<T>& data,
268  const xtl::span<const std::uint32_t>& cell_info,
269  std::int32_t cell, int data_block_size)
270  { sub_function(data, cell_info, cell, ebs * data_block_size); };
271  }
272  }
273  if (transpose)
274  {
275  if (inverse)
276  {
277  return [this](const xtl::span<T>& data,
278  const xtl::span<const std::uint32_t>& cell_info,
279  std::int32_t cell, int block_size)
280  {
281  apply_inverse_transpose_dof_transformation(data, cell_info[cell],
282  block_size);
283  };
284  }
285  else
286  {
287  return [this](const xtl::span<T>& data,
288  const xtl::span<const std::uint32_t>& cell_info,
289  std::int32_t cell, int block_size) {
290  apply_transpose_dof_transformation(data, cell_info[cell], block_size);
291  };
292  }
293  }
294  else
295  {
296  if (inverse)
297  {
298  return [this](const xtl::span<T>& data,
299  const xtl::span<const std::uint32_t>& cell_info,
300  std::int32_t cell, int block_size) {
301  apply_inverse_dof_transformation(data, cell_info[cell], block_size);
302  };
303  }
304  else
305  {
306  return [this](const xtl::span<T>& data,
307  const xtl::span<const std::uint32_t>& cell_info,
308  std::int32_t cell, int block_size)
309  { apply_dof_transformation(data, cell_info[cell], block_size); };
310  }
311  }
312  }
313 
329  template <typename T>
330  std::function<void(const xtl::span<T>&, const xtl::span<const std::uint32_t>&,
331  std::int32_t, int)>
333  bool transpose = false,
334  bool scalar_element
335  = false) const
336  {
338  {
339  // If no permutation needed, return function that does nothing
340  return [](const xtl::span<T>&, const xtl::span<const std::uint32_t>&,
341  std::int32_t, int)
342  {
343  // Do nothing
344  };
345  }
346  else if (_sub_elements.size() != 0)
347  {
348  if (_bs == 1)
349  {
350  // Mixed element
351  std::vector<std::function<void(const xtl::span<T>&,
352  const xtl::span<const std::uint32_t>&,
353  std::int32_t, int)>>
354  sub_element_functions;
355  for (std::size_t i = 0; i < _sub_elements.size(); ++i)
356  {
357  sub_element_functions.push_back(
358  _sub_elements[i]->get_dof_transformation_to_transpose_function<T>(
359  inverse, transpose));
360  }
361 
362  return [this, sub_element_functions](
363  const xtl::span<T>& data,
364  const xtl::span<const std::uint32_t>& cell_info,
365  std::int32_t cell, int block_size)
366  {
367  std::size_t start = 0;
368  for (std::size_t e = 0; e < sub_element_functions.size(); ++e)
369  {
370  const std::size_t width
371  = _sub_elements[e]->space_dimension() * block_size;
372  sub_element_functions[e](data.subspan(start, width), cell_info,
373  cell, block_size);
374  start += width;
375  }
376  };
377  }
378  else if (!scalar_element)
379  {
380  // Vector element
381  const std::function<void(const xtl::span<T>&,
382  const xtl::span<const std::uint32_t>&,
383  std::int32_t, int)>
384  sub_function = _sub_elements[0]->get_dof_transformation_function<T>(
385  inverse, transpose);
386  return [this,
387  sub_function](const xtl::span<T>& data,
388  const xtl::span<const std::uint32_t>& cell_info,
389  std::int32_t cell, int data_block_size)
390  {
391  const int ebs = block_size();
392  const std::size_t dof_count = data.size() / data_block_size;
393  for (int block = 0; block < data_block_size; ++block)
394  {
395  sub_function(data.subspan(block * dof_count, dof_count), cell_info,
396  cell, ebs);
397  }
398  };
399  }
400  }
401 
402  if (transpose)
403  {
404  if (inverse)
405  {
406  return [this](const xtl::span<T>& data,
407  const xtl::span<const std::uint32_t>& cell_info,
408  std::int32_t cell, int block_size)
409  {
411  data, cell_info[cell], block_size);
412  };
413  }
414  else
415  {
416  return [this](const xtl::span<T>& data,
417  const xtl::span<const std::uint32_t>& cell_info,
418  std::int32_t cell, int block_size)
419  {
421  block_size);
422  };
423  }
424  }
425  else
426  {
427  if (inverse)
428  {
429  return [this](const xtl::span<T>& data,
430  const xtl::span<const std::uint32_t>& cell_info,
431  std::int32_t cell, int block_size)
432  {
433  apply_inverse_dof_transformation_to_transpose(data, cell_info[cell],
434  block_size);
435  };
436  }
437  else
438  {
439  return [this](const xtl::span<T>& data,
440  const xtl::span<const std::uint32_t>& cell_info,
441  std::int32_t cell, int block_size) {
442  apply_dof_transformation_to_transpose(data, cell_info[cell],
443  block_size);
444  };
445  }
446  }
447  }
448 
454  template <typename T>
455  void apply_dof_transformation(const xtl::span<T>& data,
456  std::uint32_t cell_permutation,
457  int block_size) const
458  {
459  assert(_element);
460  _element->apply_dof_transformation(data, block_size, cell_permutation);
461  }
462 
470  template <typename T>
471  void
473  std::uint32_t cell_permutation,
474  int block_size) const
475  {
476  assert(_element);
477  _element->apply_inverse_transpose_dof_transformation(data, block_size,
478  cell_permutation);
479  }
480 
487  template <typename T>
488  void apply_transpose_dof_transformation(const xtl::span<T>& data,
489  std::uint32_t cell_permutation,
490  int block_size) const
491  {
492  assert(_element);
493  _element->apply_transpose_dof_transformation(data, block_size,
494  cell_permutation);
495  }
496 
503  template <typename T>
504  void apply_inverse_dof_transformation(const xtl::span<T>& data,
505  std::uint32_t cell_permutation,
506  int block_size) const
507  {
508  assert(_element);
509  _element->apply_inverse_dof_transformation(data, block_size,
510  cell_permutation);
511  }
512 
518  template <typename T>
519  void apply_dof_transformation_to_transpose(const xtl::span<T>& data,
520  std::uint32_t cell_permutation,
521  int block_size) const
522  {
523  assert(_element);
524  _element->apply_dof_transformation_to_transpose(data, block_size,
525  cell_permutation);
526  }
527 
533  template <typename T>
534  void
536  std::uint32_t cell_permutation,
537  int block_size) const
538  {
539  assert(_element);
540  _element->apply_inverse_dof_transformation_to_transpose(data, block_size,
541  cell_permutation);
542  }
543 
549  template <typename T>
551  const xtl::span<T>& data, std::uint32_t cell_permutation,
552  int block_size) const
553  {
554  assert(_element);
555  _element->apply_transpose_dof_transformation_to_transpose(data, block_size,
556  cell_permutation);
557  }
558 
564  template <typename T>
566  const xtl::span<T>& data, std::uint32_t cell_permutation,
567  int block_size) const
568  {
569  assert(_element);
570  _element->apply_inverse_transpose_dof_transformation_to_transpose(
571  data, block_size, cell_permutation);
572  }
573 
592  template <typename O, typename P, typename Q, typename T, typename S>
593  void map_pull_back(const O& u, const P& J, const Q& detJ, const T& K,
594  S&& U) const
595  {
596  assert(_element);
597  _element->map_pull_back_m(u, J, detJ, K, U);
598  }
599 
604  void permute_dofs(const xtl::span<std::int32_t>& doflist,
605  std::uint32_t cell_permutation) const;
606 
611  void unpermute_dofs(const xtl::span<std::int32_t>& doflist,
612  std::uint32_t cell_permutation) const;
613 
625  std::function<void(const xtl::span<std::int32_t>&, std::uint32_t)>
626  get_dof_permutation_function(bool inverse = false,
627  bool scalar_element = false) const;
628 
629 private:
630  std::string _signature, _family;
631 
632  mesh::CellType _cell_shape;
633 
634  int _tdim, _space_dim, _value_size, _reference_value_size;
635 
636  // List of sub-elements (if any)
637  std::vector<std::shared_ptr<const FiniteElement>> _sub_elements;
638 
639  // Simple hash of the signature string
640  std::size_t _hash;
641 
642  // Dimension of each value space
643  std::vector<int> _value_dimension;
644 
645  // Block size for VectorElements and TensorElements. This gives the
646  // number of DOFs colocated at each point.
647  int _bs;
648 
649  // Indicate whether the element needs permutations or transformations
650  bool _needs_dof_permutations;
651  bool _needs_dof_transformations;
652 
653  // Basix Element (nullptr for mixed elements)
654  std::unique_ptr<basix::FiniteElement> _element;
655 };
656 } // namespace dolfinx::fem
Finite Element, containing the dof layout on a reference element, and various methods for evaluating ...
Definition: FiniteElement.h:25
bool needs_dof_transformations() const noexcept
Check if DOF transformations are needed for this element.
Definition: FiniteElement.cpp:265
std::function< void(const xtl::span< T > &, const xtl::span< const std::uint32_t > &, std::int32_t, int)> get_dof_transformation_to_transpose_function(bool inverse=false, bool transpose=false, bool scalar_element=false) const
Return a function that applies DOF transformation to some transposed data.
Definition: FiniteElement.h:332
void apply_transpose_dof_transformation(const xtl::span< T > &data, std::uint32_t cell_permutation, int block_size) const
Apply transpose transformation to some data. For VectorElements, this applies the transformations for...
Definition: FiniteElement.h:488
int space_dimension() const noexcept
Dimension of the finite element function space.
Definition: FiniteElement.cpp:179
void apply_inverse_transpose_dof_transformation(const xtl::span< T > &data, std::uint32_t cell_permutation, int block_size) const
Apply inverse transpose transformation to some data. For VectorElements, this applies the transformat...
Definition: FiniteElement.h:472
std::function< void(const xtl::span< std::int32_t > &, std::uint32_t)> get_dof_permutation_function(bool inverse=false, bool scalar_element=false) const
Return a function that applies DOF transformation to some data.
Definition: FiniteElement.cpp:288
std::function< void(const xtl::span< T > &, const xtl::span< const std::uint32_t > &, std::int32_t, int)> get_dof_transformation_function(bool inverse=false, bool transpose=false, bool scalar_element=false) const
Return a function that applies DOF transformation to some data.
Definition: FiniteElement.h:211
FiniteElement & operator=(FiniteElement &&element)=default
Move assignment.
void apply_transpose_dof_transformation_to_transpose(const xtl::span< T > &data, std::uint32_t cell_permutation, int block_size) const
Apply transpose of transformation to some transposed data.
Definition: FiniteElement.h:550
const std::vector< std::shared_ptr< const FiniteElement > > & sub_elements() const noexcept
Subelements (if any)
Definition: FiniteElement.cpp:228
int reference_value_size() const noexcept
The value size, e.g. 1 for a scalar function, 2 for a 2D vector for the reference element.
Definition: FiniteElement.cpp:183
std::string family() const noexcept
The finite element family.
Definition: FiniteElement.cpp:203
void apply_inverse_dof_transformation_to_transpose(const xtl::span< T > &data, std::uint32_t cell_permutation, int block_size) const
Apply inverse of DOF transformation to some transposed data.
Definition: FiniteElement.h:535
std::size_t hash() const noexcept
Return simple hash of the signature string.
Definition: FiniteElement.cpp:234
bool needs_dof_permutations() const noexcept
Check if DOF permutations are needed for this element.
Definition: FiniteElement.cpp:270
int value_dimension(int i) const
Return the dimension of the value space for axis i.
Definition: FiniteElement.cpp:195
const xt::xtensor< double, 2 > & interpolation_points() const
Points on the reference cell at which an expression need to be evaluated in order to interpolate the ...
Definition: FiniteElement.cpp:253
void permute_dofs(const xtl::span< std::int32_t > &doflist, std::uint32_t cell_permutation) const
Permute the DOFs of the element.
Definition: FiniteElement.cpp:275
void apply_dof_transformation_to_transpose(const xtl::span< T > &data, std::uint32_t cell_permutation, int block_size) const
Apply DOF transformation to some tranposed data.
Definition: FiniteElement.h:519
bool interpolation_ident() const noexcept
Check if interpolation into the finite element space is an identity operation given the evaluation on...
Definition: FiniteElement.cpp:247
std::string signature() const noexcept
String identifying the finite element.
Definition: FiniteElement.cpp:172
void tabulate(xt::xtensor< double, 4 > &values, const xt::xtensor< double, 2 > &X, int order) const
Evaluate all derivatives of the basis functions up to given order at given points in reference cell.
Definition: FiniteElement.cpp:205
std::shared_ptr< const FiniteElement > extract_sub_element(const std::vector< int > &component) const
Extract sub finite element for component.
Definition: FiniteElement.cpp:237
void apply_inverse_dof_transformation(const xtl::span< T > &data, std::uint32_t cell_permutation, int block_size) const
Apply inverse transformation to some data. For VectorElements, this applies the transformations for t...
Definition: FiniteElement.h:504
FiniteElement(const ufc_finite_element &ufc_element)
Create finite element from UFC finite element.
Definition: FiniteElement.cpp:80
void unpermute_dofs(const xtl::span< std::int32_t > &doflist, std::uint32_t cell_permutation) const
Unpermute the DOFs of the element.
Definition: FiniteElement.cpp:281
int value_size() const noexcept
The value size, e.g. 1 for a scalar function, 2 for a 2D vector.
Definition: FiniteElement.cpp:181
int block_size() const noexcept
Block size of the finite element function space. For VectorElements and TensorElements,...
Definition: FiniteElement.cpp:193
int value_rank() const noexcept
Rank of the value space.
Definition: FiniteElement.cpp:188
FiniteElement & operator=(const FiniteElement &element)=delete
Copy assignment.
void map_pull_back(const O &u, const P &J, const Q &detJ, const T &K, S &&U) const
Pull back data from the physical element to the reference element. It can process batches of points t...
Definition: FiniteElement.h:593
void apply_inverse_transpose_dof_transformation_to_transpose(const xtl::span< T > &data, std::uint32_t cell_permutation, int block_size) const
Apply inverse transpose transformation to some transposed data.
Definition: FiniteElement.h:565
void apply_dof_transformation(const xtl::span< T > &data, std::uint32_t cell_permutation, int block_size) const
Apply DOF transformation to some data.
Definition: FiniteElement.h:455
mesh::CellType cell_shape() const noexcept
Cell shape.
Definition: FiniteElement.cpp:174
virtual ~FiniteElement()=default
Destructor.
int num_sub_elements() const noexcept
Get the number of sub elements (for a mixed element)
Definition: FiniteElement.cpp:222
constexpr void interpolate(const xt::xtensor< T, 2 > &values, xtl::span< T > dofs) const
Definition: FiniteElement.h:152
FiniteElement(FiniteElement &&element)=default
Move constructor.
FiniteElement(const FiniteElement &element)=delete
Copy constructor.
void transform_reference_basis(xt::xtensor< double, 3 > &values, const xt::xtensor< double, 3 > &reference_values, const xt::xtensor< double, 3 > &J, const xtl::span< const double > &detJ, const xt::xtensor< double, 3 > &K) const
Push basis functions forward to physical element.
Definition: FiniteElement.cpp:212
Finite element method functionality.
Definition: assemble_matrix_impl.h:23
CellType
Cell type identifier.
Definition: cell_types.h:22