XDMF Format
Using the XDMF format for visualization and time series data.
Introduction
The eXtensible Data Model and Format (XDMF) is a data format designed for scientific computing. It uses a split-file architecture:
- An XML file (
.xdmf) that describes the structure of the data - One or more HDF5 files (
.h5) that store the actual numerical data
This design combines human-readable metadata with efficient binary storage, making it ideal for large-scale simulations. XDMF files can be directly opened in ParaView for visualization.
The XDMF Class
The Rodin::
- Create an XDMF writer
- Associate a mesh with a named grid
- Add fields (scalar, vector, tensor) to the grid
- Write snapshots, optionally with time stamps
Basic Usage
#include <Rodin/Geometry.h> #include <Rodin/Variational.h> #include <Rodin/IO/XDMF.h> using namespace Rodin; using namespace Rodin::Geometry; using namespace Rodin::Variational; int main() { Mesh mesh; mesh = mesh.UniformGrid(Polytope::Type::Triangle, {16, 16}); P1 Vh(mesh); GridFunction u(Vh); u = [](const Geometry::Point& p) { return p.x() * p.y(); }; IO::XDMF xdmf("MyOutput"); xdmf.grid().setMesh(mesh).add("u", u); xdmf.write(); return 0; }
This creates:
MyOutput.xdmf— XML metadata file (open this in ParaView)MyOutput.Grid.mesh.h5— HDF5 file containing the mesh dataMyOutput.Grid.u.0.h5— HDF5 file containing the field data
Working with Grids
A grid in XDMF is a named collection of a mesh and associated fields. You obtain a grid handle from the XDMF writer:
IO::XDMF xdmf("Simulation"); auto grid = xdmf.grid(); // default grid name auto grid = xdmf.grid("Domain"); // named grid
Setting the Mesh
Each grid must have a mesh associated with it before writing:
auto grid = xdmf.grid(); grid.setMesh(mesh);
Adding Fields
Fields are added using the add() method. The field can be a named GridFunction or provided with an explicit name:
// Option 1: Name the grid function and add it GridFunction u(Vh); u.setName("velocity"); grid.add(u); // Option 2: Provide a name when adding grid.add("pressure", p);
Static Output
For a single snapshot (no time evolution), simply call write() without a time parameter:
IO::XDMF xdmf("StaticResult"); xdmf.grid().setMesh(mesh).add("solution", u); xdmf.write();
Time-Dependent Output
XDMF excels at storing time-dependent data. Call write(t) with a time value for each snapshot:
IO::XDMF xdmf("TimeEvolution"); auto grid = xdmf.grid(); grid.setMesh(mesh); grid.add(u); for (size_t k = 0; k < numSteps; ++k) { Real t = k * dt; // Update your solution... u = [t](const Geometry::Point& p) { return std::sin(Math::Constants::pi() * p.x()) * std::exp(-t); }; // Write snapshot at time t xdmf.write(t); }
When opened in ParaView, you can use the time controls to step through the evolution of the solution.
Mesh Policy
By default, the mesh is considered static (exported once). If your mesh changes over time (e.g., mesh adaptation), you can configure the writer to export the mesh at each time step:
IO::XDMF xdmf("AdaptiveSim"); auto grid = xdmf.grid(); grid.setMesh(mesh, IO::XDMF::MeshPolicy::Transient); grid.add(u);
Field Types
XDMF supports different types of fields:
Scalar Fields
Scalar fields have one component per node:
P1 Vh(mesh); // scalar P1 space GridFunction temperature(Vh); temperature = [](const Geometry::Point& p) { return p.x() + p.y(); }; grid.add("temperature", temperature);
Vector Fields
Vector fields have multiple components per node. Create the finite element space with the number of components:
P1 Vh(mesh, 2); // 2D vector P1 space GridFunction velocity(Vh); velocity = [](const Geometry::Point& p) { return Math::Vector<Real>{{ -p.y(), p.x() }}; }; grid.add("velocity", velocity);
File Patterns
The XDMF writer uses customizable file naming patterns. The default patterns are:
| File | Default Pattern | Example |
|---|---|---|
| XDMF XML | {stem}.xdmf | MyOutput.xdmf |
| Mesh HDF5 | {stem}.{grid}.mesh.h5 | MyOutput.Grid.mesh.h5 |
| Field HDF5 | {stem}.{grid}.{name}.{index}.h5 | MyOutput.Grid.u.0.h5 |
Available placeholders:
{stem}— The base name passed to the XDMF constructor{grid}— The name of the grid{name}— The name of the field{index}— The time step index{rank}— The MPI rank (for distributed output)
Opening Files in ParaView
To visualize XDMF output in ParaView:
- Open ParaView
- Use File > Open and select the
.xdmffile - Click Apply in the Properties panel
- Select a field from the dropdown menu to colorize the mesh
- For time-dependent data, use the time controls at the top
Complete Example
Here is a complete example that solves the Poisson equation and exports the result using XDMF:
#include <Rodin/Solver.h> #include <Rodin/Geometry.h> #include <Rodin/Variational.h> #include <Rodin/IO/XDMF.h> using namespace Rodin; using namespace Rodin::Geometry; using namespace Rodin::Variational; int main() { // Create mesh Mesh mesh; mesh = mesh.UniformGrid(Polytope::Type::Triangle, {32, 32}); mesh.getConnectivity().compute(1, 2); // Setup problem P1 Vh(mesh); TrialFunction u(Vh); TestFunction v(Vh); Problem poisson(u, v); poisson = Integral(Grad(u), Grad(v)) - Integral(RealFunction(1.0), v) + DirichletBC(u, Zero()); // Solve Solver::CG(poisson).solve(); // Export to XDMF IO::XDMF xdmf("PoissonResult"); xdmf.grid().setMesh(mesh).add("solution", u.getSolution()); xdmf.write(); return 0; }