Skip to content

Commit c9a6719

Browse files
committed
Add functions to write CG1 and DG0 functions to VTKHDF
1 parent 100c21c commit c9a6719

File tree

3 files changed

+152
-1
lines changed

3 files changed

+152
-1
lines changed

cpp/dolfinx/io/VTKHDF.h

Lines changed: 129 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -6,9 +6,13 @@
66

77
#include "HDF5Interface.h"
88
#include <algorithm>
9+
#include <basix/element-families.h>
10+
#include <basix/finite-element.h>
911
#include <concepts>
1012
#include <dolfinx/common/IndexMap.h>
13+
#include <dolfinx/fem/Function.h>
1114
#include <dolfinx/io/cells.h>
15+
#include <dolfinx/io/utils.h>
1216
#include <dolfinx/mesh/Mesh.h>
1317
#include <dolfinx/mesh/Topology.h>
1418
#include <dolfinx/mesh/utils.h>
@@ -298,6 +302,131 @@ void write_data(std::string point_or_cell, std::string filename,
298302
hdf5::close_file(h5file);
299303
}
300304

305+
/// @brief Write a function to VTKHDF.
306+
///
307+
/// Adds a function to an existing VTKHDF file, which already contains a mesh.
308+
///
309+
/// @tparam U Scalar type.
310+
/// @param[in] filename File for output.
311+
/// @param[in] mesh Mesh, which must be the same as the original mesh
312+
/// used in the file.
313+
/// @param[in] u Function to write to file.
314+
/// @param[in] time Timestamp.
315+
///
316+
/// @note Mesh must be written to file first using `VTKHDF::write_mesh`.
317+
/// @note Only one dataset "u" can be written per file at present, with
318+
/// multiple timesteps.
319+
/// @note Limited support for floating point types at present (no
320+
/// complex number support). This function only supports DG0 and CG1 functions.
321+
template <std::floating_point U>
322+
void write_function(std::string filename, const mesh::Mesh<U>& mesh,
323+
const fem::Function<U>& u, double time)
324+
{
325+
auto dofmap = u.function_space()->dofmap();
326+
assert(dofmap);
327+
const int bs = dofmap->bs();
328+
329+
auto map_c = mesh.topology()->index_map(mesh.topology()->dim());
330+
assert(map_c);
331+
332+
std::shared_ptr<const fem::FiniteElement<U>> element
333+
= u.function_space()->element();
334+
assert(element);
335+
336+
std::span<const std::size_t> value_shape = element->value_shape();
337+
int rank = value_shape.size();
338+
std::int32_t num_components = std::reduce(
339+
value_shape.begin(), value_shape.end(), 1, std::multiplies{});
340+
341+
std::span<const U> x = u.x()->array();
342+
343+
// Check that it is a Lagrange family element
344+
if (element->basix_element().family() != basix::element::family::P)
345+
{
346+
throw std::runtime_error("Unsupported function space. Only DG0 and CG1 are "
347+
"supported at the moment.");
348+
}
349+
350+
// DG0
351+
if (element->basix_element().degree() == 0)
352+
{
353+
const std::int32_t num_local_cells = map_c->size_local();
354+
std::vector<U> data(num_local_cells * num_components);
355+
356+
for (std::int32_t c = 0; c < num_local_cells; ++c)
357+
{
358+
auto dofs = dofmap->cell_dofs(c);
359+
assert(dofs.size() == 1);
360+
for (std::size_t i = 0; i < dofs.size(); ++i)
361+
{
362+
std::copy_n(std::cbegin(x) + bs * dofs[i], bs,
363+
std::begin(data) + num_components * c);
364+
}
365+
}
366+
367+
io::VTKHDF::write_data<U>("Cell", filename, mesh, data, time);
368+
}
369+
// CG1
370+
else if (element->basix_element().discontinuous() == false
371+
and element->basix_element().degree() == 1)
372+
{
373+
auto map_x = mesh.geometry().index_map();
374+
assert(map_x);
375+
376+
auto& geometry = mesh.geometry();
377+
auto& cmap = geometry.cmap();
378+
int cmap_dim = cmap.dim();
379+
int cell_dim = element->space_dimension() / element->block_size();
380+
if (cmap_dim != cell_dim)
381+
{
382+
throw std::runtime_error("Degree of output Function must be the same as "
383+
"mesh degree. Maybe the "
384+
"Function needs to be interpolated?");
385+
}
386+
387+
// Check that dofmap layouts are equal and check Lagrange variants
388+
if (dofmap->element_dof_layout() != cmap.create_dof_layout())
389+
{
390+
throw std::runtime_error("Function and Mesh dof layouts do not match. "
391+
"Maybe the Function needs to be interpolated?");
392+
}
393+
if (cmap.degree() > 2
394+
and element->basix_element().lagrange_variant() != cmap.variant())
395+
{
396+
throw std::runtime_error("Mismatch in Lagrange family. Maybe the "
397+
"Function needs to be interpolated?");
398+
}
399+
400+
std::int32_t num_cells = map_c->size_local() + map_c->num_ghosts();
401+
std::int32_t num_local_points = map_x->size_local();
402+
403+
// Get dof array and pack into array (padded where appropriate)
404+
auto dofmap_x = geometry.dofmap();
405+
std::vector<U> data(num_local_points * num_components);
406+
for (std::int32_t c = 0; c < num_cells; ++c)
407+
{
408+
auto dofs = dofmap->cell_dofs(c);
409+
auto dofs_x = md::submdspan(dofmap_x, c, md::full_extent);
410+
assert(dofs.size() == dofs_x.size());
411+
for (std::size_t i = 0; i < dofs.size(); ++i)
412+
{
413+
if (dofs_x[i] < num_local_points)
414+
{
415+
std::copy_n(std::cbegin(x) + bs * dofs[i], bs,
416+
std::begin(data) + num_components * dofs_x[i]);
417+
}
418+
}
419+
}
420+
421+
io::VTKHDF::write_data<U>("Point", filename, mesh, data, time);
422+
}
423+
else
424+
{
425+
throw std::runtime_error("Unsupported function space. Only DG0 and CG1 are "
426+
"supported at the moment.");
427+
}
428+
}
429+
301430
/// @brief Read a mesh from a VTKHDF format file.
302431
///
303432
/// @tparam U Scalar type of mesh

python/dolfinx/io/vtkhdf.py

Lines changed: 21 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -17,11 +17,19 @@
1717
read_vtkhdf_mesh_float32,
1818
read_vtkhdf_mesh_float64,
1919
write_vtkhdf_data,
20+
write_vtkhdf_function,
2021
write_vtkhdf_mesh,
2122
)
23+
from dolfinx.fem import Function
2224
from dolfinx.mesh import Mesh
2325

24-
__all__ = ["read_mesh", "write_cell_data", "write_mesh", "write_point_data"]
26+
__all__ = [
27+
"read_mesh",
28+
"write_cell_data",
29+
"write_function",
30+
"write_mesh",
31+
"write_point_data",
32+
]
2533

2634

2735
def read_mesh(
@@ -100,3 +108,15 @@ def write_cell_data(filename: str | Path, mesh: Mesh, data: npt.NDArray, time: f
100108
time: Timestamp.
101109
"""
102110
write_vtkhdf_data("Cell", filename, mesh._cpp_object, data, time)
111+
112+
113+
def write_function(filename: str | Path, mesh: Mesh, u: Function, time: float):
114+
"""Write a function to file. Only CG1 and DG0 functions are supported.
115+
116+
Args:
117+
filename: File to write to.
118+
mesh: Mesh.
119+
u: Function to write.
120+
time: Timestamp.
121+
"""
122+
write_vtkhdf_function(filename, mesh._cpp_object, u._cpp_object, time)

python/dolfinx/wrappers/io.cpp

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -225,6 +225,8 @@ void io(nb::module_& m)
225225
.def("write_vtkhdf_mesh", &dolfinx::io::VTKHDF::write_mesh<float>);
226226
m.def("write_vtkhdf_data", &dolfinx::io::VTKHDF::write_data<double>);
227227
m.def("write_vtkhdf_data", &dolfinx::io::VTKHDF::write_data<float>);
228+
m.def("write_vtkhdf_function", &dolfinx::io::VTKHDF::write_function<float>);
229+
m.def("write_vtkhdf_function", &dolfinx::io::VTKHDF::write_function<double>);
228230
m.def("read_vtkhdf_mesh_float64",
229231
[](MPICommWrapper comm, const std::string& filename, std::size_t gdim,
230232
std::optional<std::int32_t> max_facet_to_cell_links)

0 commit comments

Comments
 (0)