File: cells.h

package info (click to toggle)
fenics-dolfinx 1%3A0.9.0-11
  • links: PTS, VCS
  • area: main
  • in suites: forky
  • size: 5,376 kB
  • sloc: cpp: 33,701; python: 22,338; makefile: 230; sh: 171; xml: 55
file content (167 lines) | stat: -rw-r--r-- 6,707 bytes parent folder | download | duplicates (4)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
// Copyright (C) 2019 Jorgen S. Dokken
//
// This file is part of DOLFINx (https://www.fenicsproject.org)
//
// SPDX-License-Identifier:    LGPL-3.0-or-later

#pragma once

#include <array>
#include <cstdint>
#include <dolfinx/mesh/cell_types.h>
#include <span>
#include <vector>

/// @brief Functions for the re-ordering of input mesh topology to the
/// DOLFINx ordering, and transpose orderings for file output.
///
/// Basix ordering is used for the geometry nodes, and is shown below
/// for a range of cell types.
/*!
@verbatim
Triangle:               Triangle6:          Triangle10:
v
^
|
2                       2                    2
|`\                     |`\                  | \
|  `\                   |  `\                6   4
|    `\                 4    `3              |     \
|      `\               |      `\            5   9   3
|        `\             |        `\          |         \
0----------1 --> u      0-----5----1         0---7---8---1


Quadrilateral:         Quadrilateral9:         Quadrilateral16:
v
^
|
2-----------3          2-----7-----3           2--10--11---3
|           |          |           |           |           |
|           |          |           |           7  14  15   9
|           |          5     8     6           |           |
|           |          |           |           6  12  13   8
|           |          |           |           |           |
0-----------1 --> u    0-----4-----1           0---4---5---1


Tetrahedron:                    Tetrahedron10:               Tetrahedron20
            v
            /
          2                            2                         2
        ,/|`\                        ,/|`\                     ,/|`\
      ,/  |  `\                    ,/  |  `\                 13  |  `9
      ,/    '.   `\                ,8    '.   `6           ,/     4   `\
    ,/       |     `\            ,/       4     `\       12    19 |     `8
  ,/         |       `\        ,/         |       `\   ,/         |       `\
0-----------'.--------1 -> u 0--------9--'.--------1 0-----14----'.--15----1
  `\.         |      ,/        `\.         |      ,/   `\.  17     |  16 ,/
    `\.      |    ,/             `\.      |    ,5       10.   18 5    ,6
        `\.   '. ,/                  `7.   '. ,/            `\.   '.  7
          `\. |/                       `\. |/                 11. |/
              `3                            `3                    `3
                `\.
                    w


Hexahedron:          Hexahedron27:
        w
    6----------7               6----19----7
    /|   ^   v /|              /|         /|
  / |   |  / / |            17 |  25    18|
  /  |   | / /  |            / 14    24 /  15
4----------5   |           4----16----5   |
|   |   +--|---|--> u      |22 |  26  | 23|
|   2------+---3           |   2----13+---3
|  /       |  /           10  / 21   12  /
| /        | /             | 9    20  | 11
|/         |/              |/         |/
0----------1               0-----8----1


Prism:                      Prism15:
            w
            ^
            |
            3                       3
          ,/|`\                   ,/|`\
        ,/  |  `\               12  |  13
      ,/    |    `\           ,/    |    `\
    4------+------5         4------14-----5
    |      |      |         |      8      |
    |    ,/|`\    |         |      |      |
    |  ,/  |  `\  |         |      |      |
    |,/    0    `\|         10     0      11
  ./|    ,/ `\    |\        |    ,/ `\    |
  /  |  ,/     `\  | `\      |  ,6     `7  |
u    |,/         `\|   v     |,/         `\|
    1-------------2         1------9------2


Pyramid:                     Pyramid13:
                4                            4
              ,/|\                         ,/|\
            ,/ .'|\                      ,/ .'|\
  v      ,/   | | \                   ,/   | | \
    \.  ,/    .' | `.                ,/    .' | `.
      \.      |  '.  \             11      |  12  \
    ,/  \.  .'  w |   \          ,/       .'   |   \
  ,/      \. |  ^ |    \       ,/         7    |    9
2----------\'--|-3    `.     2-------10-.'----3     `.
  `\       |  \.|  `\    \      `\        |      `\    \
    `\     .'   +----`\ - \ -> u  `6     .'         8   \
      `\   |           `\  \        `\   |           `\  \
        `\.'               `\         `\.'               `\
          0-----------------1           0--------5--------1
@endverbatim
*/
namespace dolfinx::io::cells
{

/// @brief Permutation array to map from VTK to DOLFINx node ordering.
///
/// @param[in] type The cell shape
/// @param[in] num_nodes The number of cell 'nodes'
/// @return Permutation array @p for permuting from VTK ordering to
/// DOLFINx ordering, i.e. `a_dolfin[i] = a_vtk[p[i]]`.
/// @details If `p = [0, 2, 1, 3]` and `a = [10, 3, 4, 7]`, then `a_p
/// =[a[p[0]], a[p[1]], a[p[2]], a[p[3]]] = [10, 4, 3, 7]`
std::vector<std::uint16_t> perm_vtk(mesh::CellType type, int num_nodes);

/// @brief Permutation array to map from Gmsh to DOLFINx node ordering.
///
/// @param[in] type The cell shape
/// @param[in] num_nodes
/// @return Permutation array @p for permuting from Gmsh ordering to
/// DOLFINx ordering, i.e. `a_dolfin[i] = a_gmsh[p[i]]`.
/// @details If `p = [0, 2, 1, 3]` and `a = [10, 3, 4, 7]`, then `a_p
/// =[a[p[0]], a[p[1]], a[p[2]], a[p[3]]] = [10, 4, 3, 7]`
std::vector<std::uint16_t> perm_gmsh(mesh::CellType type, int num_nodes);

/// @brief Compute the transpose of a re-ordering map.
///
/// @param[in] map A re-ordering map
/// @return Transpose of the `map`. E.g., is `map = {1, 2, 3, 0}`, the
/// transpose will be `{3 , 0, 1, 2 }`.
std::vector<std::uint16_t> transpose(std::span<const std::uint16_t> map);

/// Permute cell topology by applying a permutation array for each cell
/// @param[in] cells Array of cell topologies, with each row
/// representing a cell (row-major storage)
/// @param[in] shape The shape of the `cells` array
/// @param[in] p The permutation array that maps `a_p[i] = a[p[i]]`,
/// where `a_p` is the permuted array
/// @return Permuted cell topology, where for a cell `v_new[i] =
/// v_old[map[i]]`. The storage is row-major and the shape is the same
/// as `cells`.
std::vector<std::int64_t> apply_permutation(std::span<const std::int64_t> cells,
                                            std::array<std::size_t, 2> shape,
                                            std::span<const std::uint16_t> p);

/// Get VTK cell identifier
/// @param[in] cell The cell type
/// @param[in] dim The topological dimension of the cell
/// @return The VTK cell identifier
std::int8_t get_vtk_cell_type(mesh::CellType cell, int dim);

} // namespace dolfinx::io::cells