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().cmaps()[0];
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());
101 auto apply_dof_transformation
102 = element->template get_dof_transformation_function<T>();
104 using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
105 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
106 using cmdspan4_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
107 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 4>>;
110 const std::array<std::size_t, 4> phi_shape
111 = cmap.tabulate_shape(0, Xshape[0]);
112 std::vector<T> phi_b(
113 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
114 cmdspan4_t phi_full(phi_b.data(), phi_shape);
115 cmap.tabulate(0, X, Xshape, phi_b);
116 auto phi = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
117 submdspan(phi_full, 0, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
118 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
121 auto map = topology->index_map(tdim);
123 const std::int32_t num_cells = map->size_local() + map->num_ghosts();
124 std::vector<T> x_b(scalar_dofs * gdim);
125 mdspan2_t x(x_b.data(), scalar_dofs, gdim);
126 std::vector<T> coordinate_dofs_b(num_dofs_g * gdim);
127 mdspan2_t coordinate_dofs(coordinate_dofs_b.data(), num_dofs_g, gdim);
129 std::vector<T> coords(num_nodes * 3, 0.0);
130 std::array<std::size_t, 2> cshape = {num_nodes, 3};
131 for (std::int32_t c = 0; c < num_cells; ++c)
134 for (std::size_t i = 0; i < dofmap_x.extent(1); ++i)
135 for (std::size_t j = 0; j < gdim; ++j)
136 coordinate_dofs(i, j) = x_g[3 * dofmap_x(c, i) + j];
139 cmap.push_forward(x, coordinate_dofs, phi);
140 apply_dof_transformation(x_b, std::span(cell_info.data(), cell_info.size()),
144 auto dofs = dofmap->cell_dofs(c);
145 for (std::size_t i = 0; i < dofs.size(); ++i)
147 std::int32_t dof = dofs[i];
148 for (std::size_t j = 0; j < gdim; ++j)
149 coords[3 * dof + j] = x(i, j);
154 std::vector<std::int64_t> x_id(num_nodes);
155 std::array<std::int64_t, 2> range = map_dofs->local_range();
156 std::int32_t size_local = range[1] - range[0];
157 std::iota(x_id.begin(), std::next(x_id.begin(), size_local), range[0]);
158 const std::vector<std::int64_t>& ghosts = map_dofs->ghosts();
159 std::copy(ghosts.begin(), ghosts.end(), std::next(x_id.begin(), size_local));
162 std::vector<std::uint8_t> id_ghost(num_nodes, 0);
163 std::fill(std::next(id_ghost.begin(), size_local), id_ghost.end(), 1);
165 return {std::move(coords), cshape, std::move(x_id), std::move(id_ghost)};
185std::tuple<std::vector<T>, std::array<std::size_t, 2>,
186 std::vector<std::int64_t>, std::vector<std::uint8_t>,
187 std::vector<std::int64_t>, std::array<std::size_t, 2>>
190 auto mesh = V.
mesh();
192 auto topology = mesh->topology();
194 const int tdim = topology->dim();
198 throw std::runtime_error(
"Can't create VTK mesh from a mixed element");
200 const auto [x, xshape, x_id, x_ghost]
201 = impl::tabulate_lagrange_dof_coordinates(V);
202 auto map = topology->index_map(tdim);
203 const std::size_t num_cells = map->size_local() + map->num_ghosts();
208 const int element_block_size = V.
element()->block_size();
209 const std::uint32_t num_nodes
210 = V.
element()->space_dimension() / element_block_size;
216 std::array<std::size_t, 2> shape = {num_cells, num_nodes};
217 std::vector<std::int64_t> vtk_topology(shape[0] * shape[1]);
218 for (std::size_t c = 0; c < shape[0]; ++c)
220 auto dofs = dofmap->cell_dofs(c);
221 for (std::size_t i = 0; i < dofs.size(); ++i)
222 vtk_topology[c * shape[1] + i] = dofs[vtkmap[i]];
225 return {std::move(x),
229 std::move(vtk_topology),
247std::pair<std::vector<std::int64_t>, std::array<std::size_t, 2>>
249 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
251 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:35
std::shared_ptr< const DofMap > dofmap() const
The dofmap.
Definition FunctionSpace.h:309
std::shared_ptr< const FiniteElement< T > > element() const
The finite element.
Definition FunctionSpace.h:306
std::shared_ptr< const mesh::Mesh< T > > mesh() const
The mesh.
Definition FunctionSpace.h:303
FunctionSpace(U mesh, V element, W dofmap) -> FunctionSpace< typename std::remove_cvref< typename U::element_type >::type::geometry_type::value_type >
Type deduction.
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::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:188
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