Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/d9/de1/cells_8h_source.html
DOLFINx 0.7.3
DOLFINx C++ interface
Loading...
Searching...
No Matches
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 <array>
10#include <cstdint>
11#include <dolfinx/mesh/cell_types.h>
12#include <span>
13#include <vector>
14
18{
19/*
20 The basix ordering is used for the geometry nodes, and is shown below
21 for a range of cell types.
22
23 Triangle: Triangle6: Triangle10:
24 v
25 ^
26 |
27 2 2 2
28 |`\ |`\ | \
29 | `\ | `\ 6 4
30 | `\ 4 `3 | \
31 | `\ | `\ 5 9 3
32 | `\ | `\ | \
33 0----------1 --> u 0-----5----1 0---7---8---1
34
35 Quadrilateral: Quadrilateral9: Quadrilateral16:
36 v
37 ^
38 |
39 2-----------3 2-----7-----3 2--10--11---3
40 | | | | | |
41 | | | | 7 14 15 9
42 | | 5 8 6 | |
43 | | | | 6 12 13 8
44 | | | | | |
45 0-----------1 --> u 0-----4-----1 0---4---5---1
46
47 Tetrahedron: Tetrahedron10: Tetrahedron20
48 v
49 /
50 2 2 2
51 ,/|`\ ,/|`\ ,/|`\
52 ,/ | `\ ,/ | `\ 13 | `9
53 ,/ '. `\ ,8 '. `6 ,/ 4 `\
54 ,/ | `\ ,/ 4 `\ 12 19 | `8
55 ,/ | `\ ,/ | `\ ,/ | `\
56 0-----------'.--------1 -> u 0--------9--'.--------1 0-----14----'.--15----1
57 `\. | ,/ `\. | ,/ `\. 17 | 16 ,/
58 `\. | ,/ `\. | ,5 10. 18 5 ,6
59 `\. '. ,/ `7. '. ,/ `\. '. 7
60 `\. |/ `\. |/ 11. |/
61 `3 `3 `3
62 `\.
63 w
64
65 Hexahedron: Hexahedron27:
66 w
67 6----------7 6----19----7
68 /| ^ v /| /| /|
69 / | | / / | 17 | 25 18|
70 / | | / / | / 14 24 / 15
71 4----------5 | 4----16----5 |
72 | | +--|---|--> u |22 | 26 | 23|
73 | 2------+---3 | 2----13+---3
74 | / | / 10 / 21 12 /
75 | / | / | 9 20 | 11
76 |/ |/ |/ |/
77 0----------1 0-----8----1
78
79 Prism: Prism15:
80 w
81 ^
82 |
83 3 3
84 ,/|`\ ,/|`\
85 ,/ | `\ 12 | 13
86 ,/ | `\ ,/ | `\
87 4------+------5 4------14-----5
88 | | | | 8 |
89 | ,/|`\ | | | |
90 | ,/ | `\ | | | |
91 |,/ 0 `\| 10 0 11
92 ./| ,/ `\ |\ | ,/ `\ |
93 / | ,/ `\ | `\ | ,6 `7 |
94 u |,/ `\| v |,/ `\|
95 1-------------2 1------9------2
96
97 Pyramid: Pyramid13:
98 4 4
99 ,/|\ ,/|\
100 ,/ .'|\ ,/ .'|\
101 v ,/ | | \ ,/ | | \
102 \. ,/ .' | `. ,/ .' | `.
103 \. | '. \ 11 | 12 \
104 ,/ \. .' w | \ ,/ .' | \
105 ,/ \. | ^ | \ ,/ 7 | 9
106 2----------\'--|-3 `. 2-------10-.'----3 `.
107 `\ | \.| `\ \ `\ | `\ \
108 `\ .' +----`\ - \ -> u `6 .' 8 \
109 `\ | `\ \ `\ | `\ \
110 `\.' `\ `\.' `\
111 0-----------------1 0--------5--------1
112
113*/
114
123std::vector<std::uint8_t> perm_vtk(mesh::CellType type, int num_nodes);
124
133std::vector<std::uint8_t> perm_gmsh(mesh::CellType type, int num_nodes);
134
140std::vector<std::uint8_t> transpose(const std::vector<std::uint8_t>& map);
141
151std::vector<std::int64_t>
152apply_permutation(const std::span<const std::int64_t>& cells,
153 std::array<std::size_t, 2> shape,
154 const std::span<const std::uint8_t>& p);
155
160std::int8_t get_vtk_cell_type(mesh::CellType cell, int dim);
161
162} // namespace dolfinx::io::cells
Functions for the re-ordering of input mesh topology to the DOLFINx ordering, and transpose orderings...
Definition cells.h:18
std::vector< std::uint8_t > transpose(const std::vector< std::uint8_t > &map)
Compute the transpose of a re-ordering map.
Definition cells.cpp:427
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:350
std::vector< std::int64_t > apply_permutation(const std::span< const std::int64_t > &cells, std::array< std::size_t, 2 > shape, const std::span< const std::uint8_t > &p)
Permute cell topology by applying a permutation array for each cell.
Definition cells.cpp:436
std::int8_t get_vtk_cell_type(mesh::CellType cell, int dim)
Get VTK cell identifier.
Definition cells.cpp:455
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:388
CellType
Cell type identifier.
Definition cell_types.h:22