12#include <basix/mdspan.hpp>
15#include <dolfinx/common/IndexMap.h>
17#include <dolfinx/fem/FiniteElement.h>
18#include <dolfinx/fem/FunctionSpace.h>
19#include <dolfinx/mesh/Geometry.h>
20#include <dolfinx/mesh/Mesh.h>
21#include <dolfinx/mesh/Topology.h>
30template <std::
floating_po
int T>
37template <std::
floating_po
int T>
57std::tuple<std::vector<T>, std::array<std::size_t, 2>,
58 std::vector<std::int64_t>, std::vector<std::uint8_t>>
59tabulate_lagrange_dof_coordinates(
const fem::FunctionSpace<T>& V)
63 auto topology = mesh->topology();
65 const std::size_t gdim = mesh->geometry().dim();
66 const int tdim = topology->dim();
69 auto dofmap = V.dofmap();
71 auto map_dofs = dofmap->index_map;
73 const int index_map_bs = dofmap->index_map_bs();
74 const int dofmap_bs = dofmap->bs();
77 auto element = V.element();
79 const int e_bs = element->block_size();
80 const std::size_t scalar_dofs = element->space_dimension() / e_bs;
81 const std::size_t num_nodes
82 = index_map_bs * (map_dofs->size_local() + map_dofs->num_ghosts())
87 const auto [X, Xshape] = element->interpolation_points();
88 const fem::CoordinateElement<T>& cmap = mesh->geometry().cmap();
91 auto dofmap_x = mesh->geometry().dofmap();
92 std::span<const T> x_g = mesh->geometry().x();
93 const std::size_t num_dofs_g = cmap.dim();
95 std::span<const std::uint32_t> cell_info;
96 if (element->needs_dof_transformations())
98 mesh->topology_mutable()->create_entity_permutations();
99 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
104 auto apply_dof_transformation
107 using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
108 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
109 using cmdspan4_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
110 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 4>>;
113 const std::array<std::size_t, 4> phi_shape
114 = cmap.tabulate_shape(0, Xshape[0]);
115 std::vector<T> phi_b(
116 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
117 cmdspan4_t phi_full(phi_b.data(), phi_shape);
118 cmap.tabulate(0, X, Xshape, phi_b);
119 auto phi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
120 phi_full, 0, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
121 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
124 auto map = topology->index_map(tdim);
126 const std::int32_t num_cells = map->size_local() + map->num_ghosts();
127 std::vector<T> x_b(scalar_dofs * gdim);
128 mdspan2_t x(x_b.data(), scalar_dofs, gdim);
129 std::vector<T> coordinate_dofs_b(num_dofs_g * gdim);
130 mdspan2_t coordinate_dofs(coordinate_dofs_b.data(), num_dofs_g, gdim);
132 std::vector<T> coords(num_nodes * 3, 0.0);
133 std::array<std::size_t, 2> cshape = {num_nodes, 3};
134 for (std::int32_t c = 0; c < num_cells; ++c)
137 for (std::size_t i = 0; i < dofmap_x.extent(1); ++i)
138 for (std::size_t j = 0; j < gdim; ++j)
139 coordinate_dofs(i, j) = x_g[3 * dofmap_x(c, i) + j];
142 cmap.push_forward(x, coordinate_dofs, phi);
143 apply_dof_transformation(x_b, std::span(cell_info.data(), cell_info.size()),
147 auto dofs = dofmap->cell_dofs(c);
148 for (std::size_t i = 0; i < dofs.size(); ++i)
150 std::int32_t dof = dofs[i];
151 for (std::size_t j = 0; j < gdim; ++j)
152 coords[3 * dof + j] = x(i, j);
157 std::vector<std::int64_t> x_id(num_nodes);
158 std::array<std::int64_t, 2> range = map_dofs->local_range();
159 std::int32_t size_local = range[1] - range[0];
160 std::iota(x_id.begin(), std::next(x_id.begin(), size_local), range[0]);
161 std::span ghosts = map_dofs->ghosts();
162 std::copy(ghosts.begin(), ghosts.end(), std::next(x_id.begin(), size_local));
165 std::vector<std::uint8_t> id_ghost(num_nodes, 0);
166 std::fill(std::next(id_ghost.begin(), size_local), id_ghost.end(), 1);
168 return {std::move(coords), cshape, std::move(x_id), std::move(id_ghost)};
188std::tuple<std::vector<T>, std::array<std::size_t, 2>,
189 std::vector<std::int64_t>, std::vector<std::uint8_t>,
190 std::vector<std::int64_t>, std::array<std::size_t, 2>>
193 auto mesh = V.mesh();
195 auto topology = mesh->topology();
197 const int tdim = topology->dim();
200 if (V.element()->is_mixed())
201 throw std::runtime_error(
"Can't create VTK mesh from a mixed element");
203 const auto [x, xshape, x_id, x_ghost]
204 = impl::tabulate_lagrange_dof_coordinates(V);
205 auto map = topology->index_map(tdim);
206 const std::size_t num_cells = map->size_local() + map->num_ghosts();
209 auto dofmap = V.dofmap();
211 const int element_block_size = V.element()->block_size();
212 const std::uint32_t num_nodes
213 = V.element()->space_dimension() / element_block_size;
219 std::array<std::size_t, 2> shape = {num_cells, num_nodes};
220 std::vector<std::int64_t> vtk_topology(shape[0] * shape[1]);
221 for (std::size_t c = 0; c < shape[0]; ++c)
223 auto dofs = dofmap->cell_dofs(c);
224 for (std::size_t i = 0; i < dofs.size(); ++i)
225 vtk_topology[c * shape[1] + i] = dofs[vtkmap[i]];
228 return {std::move(x),
232 std::move(vtk_topology),
250std::pair<std::vector<std::int64_t>, std::array<std::size_t, 2>>
252 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
254 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
Degree-of-freedom map representations and tools.
This class represents a finite element function space defined by a mesh, a finite element,...
Definition FunctionSpace.h:34
std::vector< std::uint16_t > perm_vtk(mesh::CellType type, int num_nodes)
Definition cells.cpp:577
std::vector< std::uint16_t > transpose(std::span< const std::uint16_t > map)
Definition cells.cpp:654
std::tuple< std::vector< T >, std::array< std::size_t, 2 >, std::vector< std::int64_t >, std::vector< std::uint8_t >, std::vector< std::int64_t >, std::array< std::size_t, 2 > > vtk_mesh_from_space(const fem::FunctionSpace< T > &V)
Given a FunctionSpace, create a topology and geometry based on the dof coordinates.
Definition vtk_utils.h:191
std::pair< std::vector< std::int64_t >, std::array< std::size_t, 2 > > extract_vtk_connectivity(MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan< const std::int32_t, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents< std::size_t, 2 > > dofmap_x, mesh::CellType cell_type)
Extract the cell topology (connectivity) in VTK ordering for all cells the mesh. The 'topology' inclu...
Definition vtk_utils.cpp:23
CellType
Cell type identifier.
Definition cell_types.h:22
Top-level namespace.
Definition defines.h:12