DOLFINx 0.10.0.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
ADIOS2Writers.h
1// Copyright (C) 2021-2023 Jørgen S. Dokken and Garth N. Wells
2//
3// This file is part of DOLFINX (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#pragma once
8
9#ifdef HAS_ADIOS2
10
11#include "vtk_utils.h"
12#include <adios2.h>
13#include <algorithm>
14#include <basix/mdspan.hpp>
15#include <cassert>
16#include <complex>
17#include <concepts>
18#include <dolfinx/common/IndexMap.h>
19#include <dolfinx/fem/DofMap.h>
20#include <dolfinx/fem/FiniteElement.h>
21#include <dolfinx/fem/Function.h>
22#include <dolfinx/mesh/Geometry.h>
23#include <dolfinx/mesh/Mesh.h>
24#include <filesystem>
25#include <memory>
26#include <mpi.h>
27#include <string>
28#include <type_traits>
29#include <variant>
30#include <vector>
31
34
35namespace dolfinx::fem
36{
37template <dolfinx::scalar T, std::floating_point U>
38class Function;
39}
40
41namespace dolfinx::io
42{
43namespace adios2_writer
44{
46template <std::floating_point T>
47using U = std::vector<std::variant<
48 std::shared_ptr<const fem::Function<float, T>>,
49 std::shared_ptr<const fem::Function<double, T>>,
50 std::shared_ptr<const fem::Function<std::complex<float>, T>>,
51 std::shared_ptr<const fem::Function<std::complex<double>, T>>>>;
52} // namespace adios2_writer
53
55class ADIOS2Writer
56{
57protected:
64 ADIOS2Writer(MPI_Comm comm, const std::filesystem::path& filename,
65 std::string tag, std::string engine);
66
68 ADIOS2Writer(ADIOS2Writer&& writer) = default;
69
71 ADIOS2Writer(const ADIOS2Writer&) = delete;
72
74 ~ADIOS2Writer();
75
77 ADIOS2Writer& operator=(ADIOS2Writer&& writer) = default;
78
79 // Copy assignment
80 ADIOS2Writer& operator=(const ADIOS2Writer&) = delete;
81
82public:
84 void close();
85
86protected:
87 std::unique_ptr<adios2::ADIOS> _adios;
88 std::unique_ptr<adios2::IO> _io;
89 std::unique_ptr<adios2::Engine> _engine;
90};
91
93namespace impl_adios2
94{
97constexpr std::array field_ext = {"_real", "_imag"};
98
101template <class T>
102adios2::Attribute<T> define_attribute(adios2::IO& io, std::string name,
103 const T& value, std::string var_name = "",
104 std::string separator = "/")
105{
106 if (adios2::Attribute<T> attr = io.InquireAttribute<T>(name); attr)
107 return attr;
108 else
109 return io.DefineAttribute<T>(name, value, var_name, separator);
110}
111
114template <class T>
115adios2::Variable<T> define_variable(adios2::IO& io, std::string name,
116 const adios2::Dims& shape = adios2::Dims(),
117 const adios2::Dims& start = adios2::Dims(),
118 const adios2::Dims& count = adios2::Dims())
119{
120 if (adios2::Variable v = io.InquireVariable<T>(name); v)
121 {
122 if (v.Count() != count and v.ShapeID() == adios2::ShapeID::LocalArray)
123 v.SetSelection({start, count});
124 return v;
125 }
126 else
127 return io.DefineVariable<T>(name, shape, start, count);
128}
129
131template <std::floating_point T>
132std::shared_ptr<const mesh::Mesh<T>>
133extract_common_mesh(const typename adios2_writer::U<T>& u)
134{
135 // Extract mesh from first function
136 assert(!u.empty());
137 auto mesh = std::visit([](auto&& u) { return u->function_space()->mesh(); },
138 u.front());
139 assert(mesh);
140
141 // Check that all functions share the same mesh
142 for (auto& v : u)
143 {
144 std::visit(
145 [&mesh](auto&& u)
146 {
147 if (mesh != u->function_space()->mesh())
148 {
149 throw std::runtime_error(
150 "ADIOS2Writer only supports functions sharing the same mesh");
151 }
152 },
153 v);
154 }
155
156 return mesh;
157}
158
159} // namespace impl_adios2
160
162namespace impl_vtx
163{
166std::stringstream create_vtk_schema(const std::vector<std::string>& point_data,
167 const std::vector<std::string>& cell_data);
168
170template <std::floating_point T>
171std::vector<std::string>
172extract_function_names(const typename adios2_writer::U<T>& u)
173{
174 std::vector<std::string> names;
175 for (auto& v : u)
176 {
177 std::visit(
178 [&names](auto&& u)
179 {
180 using U = std::decay_t<decltype(u)>;
181 using X = typename U::element_type;
182 if constexpr (std::is_floating_point_v<typename X::value_type>)
183 names.push_back(u->name);
184 else
185 {
186 names.push_back(u->name + impl_adios2::field_ext[0]);
187 names.push_back(u->name + impl_adios2::field_ext[1]);
188 }
189 },
190 v);
191 }
192
193 return names;
194}
195
205template <typename T, std::floating_point X>
206void vtx_write_data(adios2::IO& io, adios2::Engine& engine,
207 const fem::Function<T, X>& u)
208{
209 // Get function data array and information about layout
210 assert(u.x());
211 std::span<const T> u_vector = u.x()->array();
212
213 // Pad to 3D if vector/tensor is product of dimensions is smaller than
214 // 3**rank to ensure that we can visualize them correctly in Paraview
215 std::span<const std::size_t> value_shape
216 = u.function_space()->element()->value_shape();
217 int rank = value_shape.size();
218 int num_comp = std::reduce(value_shape.begin(), value_shape.end(), 1,
219 std::multiplies{});
220 if (num_comp < std::pow(3, rank))
221 num_comp = std::pow(3, rank);
222
223 std::shared_ptr<const fem::DofMap> dofmap = u.function_space()->dofmap();
224 assert(dofmap);
225 std::shared_ptr<const common::IndexMap> index_map = dofmap->index_map;
226 assert(index_map);
227 int index_map_bs = dofmap->index_map_bs();
228 int dofmap_bs = dofmap->bs();
229 std::uint32_t num_dofs = index_map_bs
230 * (index_map->size_local() + index_map->num_ghosts())
231 / dofmap_bs;
232 if constexpr (std::is_scalar_v<T>)
233 {
234 // ---- Real
235 std::vector<T> data(num_dofs * num_comp, 0);
236 for (std::size_t i = 0; i < num_dofs; ++i)
237 for (int j = 0; j < index_map_bs; ++j)
238 data[i * num_comp + j] = u_vector[i * index_map_bs + j];
239
240 adios2::Variable output = impl_adios2::define_variable<T>(
241 io, u.name, {}, {}, {num_dofs, num_comp});
242 engine.Put(output, data.data(), adios2::Mode::Sync);
243 }
244 else
245 {
246 // ---- Complex
247 using U = typename T::value_type;
248
249 std::vector<U> data(num_dofs * num_comp, 0);
250 for (std::size_t i = 0; i < num_dofs; ++i)
251 for (int j = 0; j < index_map_bs; ++j)
252 data[i * num_comp + j] = std::real(u_vector[i * index_map_bs + j]);
253
254 adios2::Variable output_real = impl_adios2::define_variable<U>(
255 io, u.name + impl_adios2::field_ext[0], {}, {}, {num_dofs, num_comp});
256 engine.Put(output_real, data.data(), adios2::Mode::Sync);
257
258 std::ranges::fill(data, 0);
259 for (std::size_t i = 0; i < num_dofs; ++i)
260 for (int j = 0; j < index_map_bs; ++j)
261 data[i * num_comp + j] = std::imag(u_vector[i * index_map_bs + j]);
262 adios2::Variable output_imag = impl_adios2::define_variable<U>(
263 io, u.name + impl_adios2::field_ext[1], {}, {}, {num_dofs, num_comp});
264 engine.Put(output_imag, data.data(), adios2::Mode::Sync);
265 }
266}
267
272template <std::floating_point T>
273void vtx_write_mesh(adios2::IO& io, adios2::Engine& engine,
274 const mesh::Mesh<T>& mesh)
275{
276 const mesh::Geometry<T>& geometry = mesh.geometry();
277 auto topology = mesh.topology();
278 assert(topology);
279
280 // "Put" geometry
281 std::shared_ptr<const common::IndexMap> x_map = geometry.index_map();
282 std::uint32_t num_vertices = x_map->size_local() + x_map->num_ghosts();
283 adios2::Variable local_geometry = impl_adios2::define_variable<T>(
284 io, "geometry", {}, {}, {num_vertices, 3});
285 engine.Put(local_geometry, geometry.x().data());
286
287 // Put number of nodes. The mesh data is written with local indices,
288 // therefore we need the ghost vertices.
289 adios2::Variable vertices = impl_adios2::define_variable<std::uint32_t>(
290 io, "NumberOfNodes", {adios2::LocalValueDim});
291 engine.Put<std::uint32_t>(vertices, num_vertices);
292
293 auto [vtkcells, shape]
294 = io::extract_vtk_connectivity(geometry.dofmap(), topology->cell_type());
295
296 // Add cell metadata
297 int tdim = topology->dim();
298 adios2::Variable cell_var = impl_adios2::define_variable<std::uint32_t>(
299 io, "NumberOfCells", {adios2::LocalValueDim});
300 engine.Put<std::uint32_t>(cell_var, shape[0]);
301 adios2::Variable celltype_var
302 = impl_adios2::define_variable<std::uint32_t>(io, "types");
303 engine.Put<std::uint32_t>(
304 celltype_var, cells::get_vtk_cell_type(topology->cell_type(), tdim));
305
306 // Pack mesh 'nodes'. Output is written as [N0, v0_0,...., v0_N0, N1,
307 // v1_0,...., v1_N1,....], where N is the number of cell nodes and v0,
308 // etc, is the node index
309 std::vector<std::int64_t> cells(shape[0] * (shape[1] + 1), shape[1]);
310 for (std::size_t c = 0; c < shape[0]; ++c)
311 {
312 std::span vtkcell(vtkcells.data() + c * shape[1], shape[1]);
313 std::span cell(cells.data() + c * (shape[1] + 1), shape[1] + 1);
314 std::ranges::copy(vtkcell, std::next(cell.begin()));
315 }
316
317 // Put topology (nodes)
318 adios2::Variable local_topology = impl_adios2::define_variable<std::int64_t>(
319 io, "connectivity", {}, {}, {shape[0], shape[1] + 1});
320 engine.Put(local_topology, cells.data());
321
322 // Vertex global ids and ghost markers
323 adios2::Variable orig_id = impl_adios2::define_variable<std::int64_t>(
324 io, "vtkOriginalPointIds", {}, {}, {num_vertices});
325 engine.Put(orig_id, geometry.input_global_indices().data());
326
327 std::vector<std::uint8_t> x_ghost(num_vertices, 0);
328 std::fill(std::next(x_ghost.begin(), x_map->size_local()), x_ghost.end(), 1);
329 adios2::Variable ghost = impl_adios2::define_variable<std::uint8_t>(
330 io, "vtkGhostType", {}, {}, {x_ghost.size()});
331 engine.Put(ghost, x_ghost.data());
332 engine.PerformPuts();
333}
334
342template <std::floating_point T>
343std::pair<std::vector<std::int64_t>, std::vector<std::uint8_t>>
344vtx_write_mesh_from_space(adios2::IO& io, adios2::Engine& engine,
345 const fem::FunctionSpace<T>& V)
346{
347 auto mesh = V.mesh();
348 assert(mesh);
349 auto topology = mesh->topology();
350 assert(topology);
351 int tdim = topology->dim();
352
353 // Get a VTK mesh with points at the 'nodes'
354 auto [x, xshape, x_id, x_ghost, vtk, vtkshape] = io::vtk_mesh_from_space(V);
355
356 std::uint32_t num_dofs = xshape[0];
357
358 // -- Pack mesh 'nodes'. Output is written as [N0, v0_0,...., v0_N0, N1,
359 // v1_0,...., v1_N1,....], where N is the number of cell nodes and v0,
360 // etc, is the node index.
361
362 // Create vector, setting all entries to nodes per cell (vtk.shape(1))
363 std::vector<std::int64_t> cells(vtkshape[0] * (vtkshape[1] + 1), vtkshape[1]);
364
365 // Set the [v0_0,...., v0_N0, v1_0,...., v1_N1,....] data
366 for (std::size_t c = 0; c < vtkshape[0]; ++c)
367 {
368 std::span vtkcell(vtk.data() + c * vtkshape[1], vtkshape[1]);
369 std::span cell(cells.data() + c * (vtkshape[1] + 1), vtkshape[1] + 1);
370 std::ranges::copy(vtkcell, std::next(cell.begin()));
371 }
372
373 // Define ADIOS2 variables for geometry, topology, celltypes and
374 // corresponding VTK data
375 adios2::Variable local_geometry
376 = impl_adios2::define_variable<T>(io, "geometry", {}, {}, {num_dofs, 3});
377 adios2::Variable local_topology = impl_adios2::define_variable<std::int64_t>(
378 io, "connectivity", {}, {}, {vtkshape[0], vtkshape[1] + 1});
379 adios2::Variable cell_type
380 = impl_adios2::define_variable<std::uint32_t>(io, "types");
381 adios2::Variable vertices = impl_adios2::define_variable<std::uint32_t>(
382 io, "NumberOfNodes", {adios2::LocalValueDim});
383 adios2::Variable elements = impl_adios2::define_variable<std::uint32_t>(
384 io, "NumberOfEntities", {adios2::LocalValueDim});
385
386 // Write mesh information to file
387 engine.Put<std::uint32_t>(vertices, num_dofs);
388 engine.Put<std::uint32_t>(elements, vtkshape[0]);
389 engine.Put<std::uint32_t>(
390 cell_type, cells::get_vtk_cell_type(topology->cell_type(), tdim));
391 engine.Put(local_geometry, x.data());
392 engine.Put(local_topology, cells.data());
393
394 // Node global ids
395 adios2::Variable orig_id = impl_adios2::define_variable<std::int64_t>(
396 io, "vtkOriginalPointIds", {}, {}, {x_id.size()});
397 engine.Put(orig_id, x_id.data());
398 adios2::Variable ghost = impl_adios2::define_variable<std::uint8_t>(
399 io, "vtkGhostType", {}, {}, {x_ghost.size()});
400 engine.Put(ghost, x_ghost.data());
401
402 engine.PerformPuts();
403 return {std::move(x_id), std::move(x_ghost)};
404}
405} // namespace impl_vtx
406
408enum class VTXMeshPolicy
409{
410 update,
411 reuse
413};
414
420template <std::floating_point T>
421class VTXWriter : public ADIOS2Writer
422{
423public:
435 VTXWriter(MPI_Comm comm, const std::filesystem::path& filename,
436 std::shared_ptr<const mesh::Mesh<T>> mesh,
437 std::string engine = "BPFile")
438 : ADIOS2Writer(comm, filename, "VTX mesh writer", engine), _mesh(mesh),
439 _mesh_reuse_policy(VTXMeshPolicy::update), _is_piecewise_constant(false)
440 {
441 // Define VTK scheme attribute for mesh
442 std::string vtk_scheme = impl_vtx::create_vtk_schema({}, {}).str();
443 impl_adios2::define_attribute<std::string>(*_io, "vtk.xml", vtk_scheme);
444 }
445
461 VTXWriter(MPI_Comm comm, const std::filesystem::path& filename,
462 const typename adios2_writer::U<T>& u, std::string engine,
463 VTXMeshPolicy mesh_policy = VTXMeshPolicy::update)
464 : ADIOS2Writer(comm, filename, "VTX function writer", engine),
465 _mesh(impl_adios2::extract_common_mesh<T>(u)),
466 _mesh_reuse_policy(mesh_policy), _u(u), _is_piecewise_constant(false)
467 {
468 if (u.empty())
469 throw std::runtime_error("VTXWriter fem::Function list is empty.");
470
471 // Extract space from first function
472 auto V0 = std::visit([](auto& u) { return u->function_space().get(); },
473 u.front());
474 assert(V0);
475 auto element0 = V0->element().get();
476 assert(element0);
477
478 // Check if function is mixed
479 if (element0->is_mixed())
480 {
481 throw std::runtime_error(
482 "Mixed functions are not supported by VTXWriter.");
483 }
484
485 // FIXME: is the below check adequate for detecting a Lagrange
486 // element?
487 // Check that element is Lagrange
488 if (!element0->interpolation_ident())
489 {
490 throw std::runtime_error(
491 "Only (discontinuous) Lagrange functions are "
492 "supported. Interpolate Functions before output.");
493 }
494
495 // Check if function is DG 0
496 if (element0->space_dimension() / element0->block_size() == 1)
497 _is_piecewise_constant = true;
498
499 // Check that all functions come from same element type
500 for (auto& v : _u)
501 {
502 std::visit(
503 [V0](auto& u)
504 {
505 auto element = u->function_space()->element();
506 assert(element);
507 if (*element != *V0->element().get())
508 {
509 throw std::runtime_error("All functions in VTXWriter must have "
510 "the same element type.");
511 }
512#ifndef NDEBUG
513 auto dmap0 = V0->dofmap()->map();
514 auto dmap = u->function_space()->dofmap()->map();
515 if (dmap0.size() != dmap.size()
516 or !std::equal(dmap0.data_handle(),
517 dmap0.data_handle() + dmap0.size(),
518 dmap.data_handle()))
519 {
520 throw std::runtime_error(
521 "All functions must have the same dofmap for VTXWriter.");
522 }
523#endif
524 },
525 v);
526 }
527
528 // Define VTK scheme attribute for set of functions
529 std::vector<std::string> names = impl_vtx::extract_function_names<T>(u);
530 std::string vtk_scheme;
531 if (_is_piecewise_constant)
532 vtk_scheme = impl_vtx::create_vtk_schema({}, names).str();
533 else
534 vtk_scheme = impl_vtx::create_vtk_schema(names, {}).str();
535
536 impl_adios2::define_attribute<std::string>(*_io, "vtk.xml", vtk_scheme);
537 }
538
553 VTXWriter(MPI_Comm comm, const std::filesystem::path& filename,
554 const typename adios2_writer::U<T>& u,
555 VTXMeshPolicy mesh_policy = VTXMeshPolicy::update)
556 : VTXWriter(comm, filename, u, "BPFile", mesh_policy)
557 {
558 }
559
560 // Copy constructor
561 VTXWriter(const VTXWriter&) = delete;
562
564 VTXWriter(VTXWriter&& file) = default;
565
567 ~VTXWriter() = default;
568
570 VTXWriter& operator=(VTXWriter&&) = default;
571
572 // Copy assignment
573 VTXWriter& operator=(const VTXWriter&) = delete;
574
577 void write(double t)
578 {
579 assert(_io);
580 adios2::Variable var_step
581 = impl_adios2::define_variable<double>(*_io, "step");
582
583 assert(_engine);
584 _engine->BeginStep();
585 _engine->template Put<double>(var_step, t);
586
587 // If we have no functions or DG functions write the mesh to file
588 if (_is_piecewise_constant or _u.empty())
589 {
590 impl_vtx::vtx_write_mesh(*_io, *_engine, *_mesh);
591 if (_is_piecewise_constant)
592 {
593 for (auto& v : _u)
594 {
595 std::visit([&](auto& u)
596 { impl_vtx::vtx_write_data(*_io, *_engine, *u); }, v);
597 }
598 }
599 }
600 else
601 {
602 if (_mesh_reuse_policy == VTXMeshPolicy::update
603 or !(_io->template InquireVariable<std::int64_t>("connectivity")))
604 {
605 // Write a single mesh for functions as they share finite
606 // element
607 std::tie(_x_id, _x_ghost) = std::visit(
608 [&](auto& u)
609 {
610 return impl_vtx::vtx_write_mesh_from_space(*_io, *_engine,
611 *u->function_space());
612 },
613 _u[0]);
614 }
615 else
616 {
617 // Node global ids
618 adios2::Variable orig_id = impl_adios2::define_variable<std::int64_t>(
619 *_io, "vtkOriginalPointIds", {}, {}, {_x_id.size()});
620 _engine->Put(orig_id, _x_id.data());
621 adios2::Variable ghost = impl_adios2::define_variable<std::uint8_t>(
622 *_io, "vtkGhostType", {}, {}, {_x_ghost.size()});
623 _engine->Put(ghost, _x_ghost.data());
624 _engine->PerformPuts();
625 }
626
627 // Write function data for each function to file
628 for (auto& v : _u)
629 {
630 std::visit([&](auto& u)
631 { impl_vtx::vtx_write_data(*_io, *_engine, *u); }, v);
632 }
633 }
634
635 _engine->EndStep();
636 }
637
638private:
639 std::shared_ptr<const mesh::Mesh<T>> _mesh;
640 adios2_writer::U<T> _u;
641
642 // Control whether the mesh is written to file once or at every time
643 // step
644 VTXMeshPolicy _mesh_reuse_policy;
645 std::vector<std::int64_t> _x_id;
646 std::vector<std::uint8_t> _x_ghost;
647
648 // Special handling of piecewise constant functions
649 bool _is_piecewise_constant;
650};
651
653template <typename U, typename T>
654VTXWriter(MPI_Comm comm, U filename, T mesh)
655 -> VTXWriter<typename std::remove_cvref<
656 typename T::element_type>::type::geometry_type::value_type>;
657
658} // namespace dolfinx::io
659
660#endif
Degree-of-freedom map representations and tools.
int rank(MPI_Comm comm)
Return process rank for the communicator.
Definition MPI.cpp:64
void cells(la::SparsityPattern &pattern, std::array< std::span< const std::int32_t >, 2 > cells, std::array< std::reference_wrapper< const DofMap >, 2 > dofmaps)
Iterate over cells and insert entries into sparsity pattern.
Definition sparsitybuild.cpp:16
Finite element method functionality.
Definition assemble_matrix_impl.h:26
Support for file IO.
Definition cells.h:119