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
cells.h
1 // Copyright (C) 2019 Jorgen S. Dokken
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 <cstdint>
10 #include <dolfinx/mesh/cell_types.h>
11 #include <vector>
12 #include <xtensor/xtensor.hpp>
13 
14 namespace dolfinx::mesh
15 {
16 class Mesh;
17 }
18 
22 {
23 /*
24  The basix ordering is used for the geometry nodes, and is shown below
25  for a range of cell types.
26 
27  Triangle: Triangle6: Triangle10:
28  v
29  ^
30  |
31  2 2 2
32  |`\ |`\ | \
33  | `\ | `\ 6 4
34  | `\ 4 `3 | \
35  | `\ | `\ 5 9 3
36  | `\ | `\ | \
37  0----------1 --> u 0-----5----1 0---7---8---1
38 
39  Quadrilateral: Quadrilateral9: Quadrilateral16:
40  v
41  ^
42  |
43  2-----------3 2-----7-----3 2--10--11---3
44  | | | | | |
45  | | | | 7 14 15 9
46  | | 5 8 6 | |
47  | | | | 6 12 13 8
48  | | | | | |
49  0-----------1 --> u 0-----4-----1 0---4---5---1
50 
51  Tetrahedron: Tetrahedron10: Tetrahedron20
52  v
53  /
54  2 2 2
55  ,/|`\ ,/|`\ ,/|`\
56  ,/ | `\ ,/ | `\ 13 | `9
57  ,/ '. `\ ,8 '. `6 ,/ 4 `\
58  ,/ | `\ ,/ 4 `\ 12 19 | `8
59  ,/ | `\ ,/ | `\ ,/ | `\
60  0-----------'.--------1 -> u 0--------9--'.--------1 0-----14----'.--15----1
61  `\. | ,/ `\. | ,/ `\. 17 | 16 ,/
62  `\. | ,/ `\. | ,5 10. 18 5 ,6
63  `\. '. ,/ `7. '. ,/ `\. '. 7
64  `\. |/ `\. |/ 11. |/
65  `3 `3 `3
66  `\.
67  w
68 
69  Hexahedron: Hexahedron27:
70  w
71  6----------7 6----19----7
72  /| ^ v /| /| /|
73  / | | / / | 17 | 25 18|
74  / | | / / | / 14 24 / 15
75  4----------5 | 4----16----5 |
76  | | +--|---|--> u |22 | 26 | 23|
77  | 2------+---3 | 2----13+---3
78  | / | / 10 / 21 12 /
79  | / | / | 9 20 | 11
80  |/ |/ |/ |/
81  0----------1 0-----8----1
82 
83 */
84 
93 std::vector<std::uint8_t> perm_vtk(mesh::CellType type, int num_nodes);
94 
103 std::vector<std::uint8_t> perm_gmsh(mesh::CellType type, int num_nodes);
104 
110 std::vector<std::uint8_t> transpose(const std::vector<std::uint8_t>& map);
111 
119 xt::xtensor<std::int64_t, 2>
120 compute_permutation(const xt::xtensor<std::int64_t, 2>& cells,
121  const std::vector<std::uint8_t>& p);
122 
127 std::int8_t get_vtk_cell_type(const dolfinx::mesh::Mesh& mesh, int dim);
128 
129 } // namespace dolfinx::io::cells
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition: Mesh.h:53
Functions for the re-ordering of input mesh topology to the DOLFINx ordering, and transpose orderings...
Definition: cells.h:22
std::vector< std::uint8_t > perm_gmsh(mesh::CellType type, int num_nodes)
Permutation array to map from Gmsh to DOLFINx node ordering.
Definition: cells.cpp:312
std::int8_t get_vtk_cell_type(const dolfinx::mesh::Mesh &mesh, int dim)
Get VTK cell identifier.
Definition: cells.cpp:368
std::vector< std::uint8_t > perm_vtk(mesh::CellType type, int num_nodes)
Permutation array to map from VTK to DOLFINx node ordering.
Definition: cells.cpp:280
xt::xtensor< std::int64_t, 2 > compute_permutation(const xt::xtensor< std::int64_t, 2 > &cells, const std::vector< std::uint8_t > &p)
Permute cell topology by applying a permutation array for each cell.
Definition: cells.cpp:354
std::vector< std::uint8_t > transpose(const std::vector< std::uint8_t > &map)
Compute the transpose of a re-ordering map.
Definition: cells.cpp:345
Mesh data structures and algorithms on meshes.
Definition: DirichletBC.h:20
CellType
Cell type identifier.
Definition: cell_types.h:22