DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
Functions
dolfinx::io::cells Namespace Reference

Functions for the re-ordering of input mesh topology to the DOLFINx ordering, and transpose orderings for file output. More...

Functions

std::vector< std::uint16_t > perm_vtk (mesh::CellType type, int num_nodes)
 
std::vector< std::uint16_t > perm_gmsh (mesh::CellType type, int num_nodes)
 
std::vector< std::uint16_t > transpose (std::span< const std::uint16_t > map)
 
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)
 
std::int8_t get_vtk_cell_type (mesh::CellType cell, int dim)
 

Detailed Description

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.

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

Function Documentation

◆ apply_permutation()

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 )

Permute cell topology by applying a permutation array for each cell

Parameters
[in]cellsArray of cell topologies, with each row representing a cell (row-major storage)
[in]shapeThe shape of the cells array
[in]pThe permutation array that maps a_p[i] = a[p[i]], where a_p is the permuted array
Returns
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.

◆ get_vtk_cell_type()

std::int8_t get_vtk_cell_type ( mesh::CellType cell,
int dim )

Get VTK cell identifier

Parameters
[in]cellThe cell type
[in]dimThe topological dimension of the cell
Returns
The VTK cell identifier

◆ perm_gmsh()

std::vector< std::uint16_t > perm_gmsh ( mesh::CellType type,
int num_nodes )

Permutation array to map from Gmsh to DOLFINx node ordering

Parameters
[in]typeThe cell shape
[in]num_nodes
Returns
Permutation array for permuting from Gmsh ordering to DOLFINx ordering, i.e. a_dolfin[i] = a_gmsh[p[i]] @details Ifp = [0, 2, 1, 3]anda = [10, 3, 4, 7], thena_p =[a[p[0]], a[p[1]], a[p[2]], a[p[3]]] = [10, 4, 3, 7]`

◆ perm_vtk()

std::vector< std::uint16_t > perm_vtk ( mesh::CellType type,
int num_nodes )

Permutation array to map from VTK to DOLFINx node ordering

Parameters
[in]typeThe cell shape
[in]num_nodesThe number of cell 'nodes'
Returns
Permutation array for permuting from VTK ordering to DOLFINx ordering, i.e. a_dolfin[i] = a_vtk[p[i]] @details Ifp = [0, 2, 1, 3]anda = [10, 3, 4, 7], thena_p =[a[p[0]], a[p[1]], a[p[2]], a[p[3]]] = [10, 4, 3, 7]`

◆ transpose()

std::vector< std::uint16_t > transpose ( std::span< const std::uint16_t > map)

Compute the transpose of a re-ordering map

Parameters
[in]mapA re-ordering map
Returns
Transpose of the map. E.g., is map = {1, 2, 3, 0}, the transpose will be {3 , 0, 1, 2 }.