Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/df/d34/ADIOS2Writers_8h_source.html
DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
ADIOS2Writers.h
Go to the documentation of this file.
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 <basix/mdspan.hpp>
14#include <cassert>
15#include <complex>
16#include <concepts>
17#include <dolfinx/common/IndexMap.h>
18#include <dolfinx/fem/DofMap.h>
19#include <dolfinx/fem/FiniteElement.h>
20#include <dolfinx/fem/Function.h>
21#include <dolfinx/mesh/Geometry.h>
22#include <dolfinx/mesh/Mesh.h>
23#include <filesystem>
24#include <memory>
25#include <mpi.h>
26#include <string>
27#include <type_traits>
28#include <variant>
29#include <vector>
30
33
34namespace dolfinx::fem
35{
36template <dolfinx::scalar T, std::floating_point U>
37class Function;
38}
39
40namespace dolfinx::io
41{
42namespace adios2_writer
43{
45template <std::floating_point T>
46using U = std::vector<std::variant<
47 std::shared_ptr<const fem::Function<float, T>>,
48 std::shared_ptr<const fem::Function<double, T>>,
49 std::shared_ptr<const fem::Function<std::complex<float>, T>>,
50 std::shared_ptr<const fem::Function<std::complex<double>, T>>>>;
51} // namespace adios2_writer
52
55{
56protected:
63 ADIOS2Writer(MPI_Comm comm, const std::filesystem::path& filename,
64 std::string tag, std::string engine);
65
67 ADIOS2Writer(ADIOS2Writer&& writer) = default;
68
70 ADIOS2Writer(const ADIOS2Writer&) = delete;
71
74
76 ADIOS2Writer& operator=(ADIOS2Writer&& writer) = default;
77
78 // Copy assignment
79 ADIOS2Writer& operator=(const ADIOS2Writer&) = delete;
80
81public:
83 void close();
84
85protected:
86 std::unique_ptr<adios2::ADIOS> _adios;
87 std::unique_ptr<adios2::IO> _io;
88 std::unique_ptr<adios2::Engine> _engine;
89};
90
92namespace impl_adios2
93{
96constexpr std::array field_ext = {"_real", "_imag"};
97
100template <class T>
101adios2::Attribute<T> define_attribute(adios2::IO& io, std::string name,
102 const T& value, std::string var_name = "",
103 std::string separator = "/")
104{
105 if (adios2::Attribute<T> attr = io.InquireAttribute<T>(name); attr)
106 return attr;
107 else
108 return io.DefineAttribute<T>(name, value, var_name, separator);
109}
110
113template <class T>
114adios2::Variable<T> define_variable(adios2::IO& io, std::string name,
115 const adios2::Dims& shape = adios2::Dims(),
116 const adios2::Dims& start = adios2::Dims(),
117 const adios2::Dims& count = adios2::Dims())
118{
119 if (adios2::Variable v = io.InquireVariable<T>(name); v)
120 {
121 if (v.Count() != count and v.ShapeID() == adios2::ShapeID::LocalArray)
122 v.SetSelection({start, count});
123 return v;
124 }
125 else
126 return io.DefineVariable<T>(name, shape, start, count);
127}
128
130template <std::floating_point T>
131std::shared_ptr<const mesh::Mesh<T>>
132extract_common_mesh(const typename adios2_writer::U<T>& u)
133{
134 // Extract mesh from first function
135 assert(!u.empty());
136 auto mesh = std::visit([](auto&& u) { return u->function_space()->mesh(); },
137 u.front());
138 assert(mesh);
139
140 // Check that all functions share the same mesh
141 for (auto& v : u)
142 {
143 std::visit(
144 [&mesh](auto&& u)
145 {
146 if (mesh != u->function_space()->mesh())
147 {
148 throw std::runtime_error(
149 "ADIOS2Writer only supports functions sharing the same mesh");
150 }
151 },
152 v);
153 }
154
155 return mesh;
156}
157
158} // namespace impl_adios2
159
161namespace impl_fides
162{
167void initialize_mesh_attributes(adios2::IO& io, mesh::CellType type);
168
173template <std::floating_point T>
175 const typename adios2_writer::U<T>& u)
176{
177 // Array of function (name, cell association types) for each function
178 // added to the file
179 std::vector<std::array<std::string, 2>> u_data;
180 for (auto& v : u)
181 {
182 std::visit(
183 [&u_data](auto&& u)
184 {
185 using U = std::decay_t<decltype(u)>;
186 using X = typename U::element_type;
187 if constexpr (std::is_floating_point_v<typename X::value_type>)
188 u_data.push_back({u->name, "points"});
189 else
190 {
191 u_data.push_back({u->name + impl_adios2::field_ext[0], "points"});
192 u_data.push_back({u->name + impl_adios2::field_ext[1], "points"});
193 }
194 },
195 v);
196 }
197
198 // Write field associations to file
199 if (adios2::Attribute<std::string> assc
200 = io.InquireAttribute<std::string>("Fides_Variable_Associations");
201 !assc)
202 {
203 std::vector<std::string> u_type;
204 std::transform(u_data.cbegin(), u_data.cend(), std::back_inserter(u_type),
205 [](auto& f) { return f[1]; });
206 io.DefineAttribute<std::string>("Fides_Variable_Associations",
207 u_type.data(), u_type.size());
208 }
209
210 // Write field pointers to file
211 if (adios2::Attribute<std::string> fields
212 = io.InquireAttribute<std::string>("Fides_Variable_List");
213 !fields)
214 {
215 std::vector<std::string> names;
216 std::transform(u_data.cbegin(), u_data.cend(), std::back_inserter(names),
217 [](auto& f) { return f[0]; });
218 io.DefineAttribute<std::string>("Fides_Variable_List", names.data(),
219 names.size());
220 }
221}
222
225template <typename T, std::floating_point U>
227{
228 auto V = u.function_space();
229 assert(V);
230 auto dofmap = V->dofmap();
231 assert(dofmap);
232 auto mesh = V->mesh();
233 assert(mesh);
234 const mesh::Geometry<U>& geometry = mesh->geometry();
235 auto topology = mesh->topology();
236 assert(topology);
237
238 // The Function and the mesh must have identical element_dof_layouts
239 // (up to the block size)
240 assert(dofmap->element_dof_layout() == geometry.cmap().create_dof_layout());
241
242 int tdim = topology->dim();
243 auto cell_map = topology->index_map(tdim);
244 assert(cell_map);
245 std::int32_t num_cells = cell_map->size_local() + cell_map->num_ghosts();
246
247 auto vertex_map = topology->index_map(0);
248 assert(vertex_map);
249 std::uint32_t num_vertices
250 = vertex_map->size_local() + vertex_map->num_ghosts();
251
252 int rank = V->value_shape().size();
253 std::uint32_t num_components = std::pow(3, rank);
254
255 // Get dof array and pack into array (padded where appropriate)
256 auto dofmap_x = geometry.dofmap();
257 const int bs = dofmap->bs();
258 std::span<const T> u_data = u.x()->array();
259 std::vector<T> data(num_vertices * num_components, 0);
260 for (std::int32_t c = 0; c < num_cells; ++c)
261 {
262 auto dofs = dofmap->cell_dofs(c);
263 auto dofs_x = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
264 dofmap_x, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
265 assert(dofs.size() == dofs_x.size());
266 for (std::size_t i = 0; i < dofs.size(); ++i)
267 for (int j = 0; j < bs; ++j)
268 data[num_components * dofs_x[i] + j] = u_data[bs * dofs[i] + j];
269 }
270
271 return data;
272}
273
280template <typename T, std::floating_point U>
281void write_data(adios2::IO& io, adios2::Engine& engine,
282 const fem::Function<T, U>& u)
283{
284 // FIXME: There is an implicit assumptions that u and the mesh have
285 // the same ElementDoflayout
286 auto V = u.function_space();
287 assert(V);
288 auto dofmap = V->dofmap();
289 assert(dofmap);
290 auto mesh = V->mesh();
291 assert(mesh);
292 const int gdim = mesh->geometry().dim();
293
294 // Vectors and tensor need padding in gdim < 3
295 int rank = V->value_shape().size();
296 bool need_padding = rank > 0 and gdim != 3 ? true : false;
297
298 // Get vertex data. If the mesh and function dofmaps are the same we
299 // can work directly with the dof array.
300 std::span<const T> data;
301 std::vector<T> _data;
302 auto eq_check = [](auto x, auto y) -> bool
303 {
304 return x.extents() == y.extents()
305 and std::equal(x.data_handle(), x.data_handle() + x.size(),
306 y.data_handle());
307 };
308
309 if (!need_padding and eq_check(mesh->geometry().dofmap(), dofmap->map()))
310 data = u.x()->array();
311 else
312 {
313 _data = impl_fides::pack_function_data(u);
314 data = std::span<const T>(_data);
315 }
316
317 auto vertex_map = mesh->topology()->index_map(0);
318 assert(vertex_map);
319 std::uint32_t num_vertices
320 = vertex_map->size_local() + vertex_map->num_ghosts();
321
322 // Write each real and imaginary part of the function
323 std::uint32_t num_components = std::pow(3, rank);
324 assert(data.size() % num_components == 0);
325 if constexpr (std::is_scalar_v<T>)
326 {
327 // ---- Real
328 adios2::Variable local_output = impl_adios2::define_variable<T>(
329 io, u.name, {}, {}, {num_vertices, num_components});
330
331 // To reuse out_data, we use sync mode here
332 engine.Put(local_output, data.data());
333 engine.PerformPuts();
334 }
335 else
336 {
337 // ---- Complex
338 using X = typename T::value_type;
339
340 std::vector<X> data_real(data.size()), data_imag(data.size());
341
342 adios2::Variable local_output_r = impl_adios2::define_variable<X>(
343 io, u.name + impl_adios2::field_ext[0], {}, {},
344 {num_vertices, num_components});
345 std::transform(data.begin(), data.end(), data_real.begin(),
346 [](auto x) -> X { return std::real(x); });
347 engine.Put(local_output_r, data_real.data());
348
349 adios2::Variable local_output_c = impl_adios2::define_variable<X>(
350 io, u.name + impl_adios2::field_ext[1], {}, {},
351 {num_vertices, num_components});
352 std::transform(data.begin(), data.end(), data_imag.begin(),
353 [](auto x) -> X { return std::imag(x); });
354 engine.Put(local_output_c, data_imag.data());
355 engine.PerformPuts();
356 }
357}
358
363template <std::floating_point T>
364void write_mesh(adios2::IO& io, adios2::Engine& engine,
365 const mesh::Mesh<T>& mesh)
366{
367 const mesh::Geometry<T>& geometry = mesh.geometry();
368 auto topology = mesh.topology();
369 assert(topology);
370
371 // "Put" geometry data
372 auto x_map = geometry.index_map();
373 std::uint32_t num_vertices = x_map->size_local() + x_map->num_ghosts();
374 adios2::Variable local_geometry = impl_adios2::define_variable<T>(
375 io, "points", {}, {}, {num_vertices, 3});
376 engine.Put(local_geometry, geometry.x().data());
377
378 // TODO: The DOLFINx and VTK topology are the same for some cell types
379 // - no need to repack via extract_vtk_connectivity in these cases
380
381 // Get topological dimension, number of cells and number of 'nodes' per
382 // cell, and compute 'VTK' connectivity
383 int tdim = topology->dim();
384 std::int32_t num_cells = topology->index_map(tdim)->size_local();
385 int num_nodes = geometry.cmap().dim();
386 auto [cells, shape] = io::extract_vtk_connectivity(mesh.geometry().dofmap(),
387 topology->cell_type());
388
389 // "Put" topology data in the result in the ADIOS2 file
390 adios2::Variable local_topology = impl_adios2::define_variable<std::int64_t>(
391 io, "connectivity", {}, {}, {std::size_t(num_cells * num_nodes)});
392 engine.Put(local_topology, cells.data());
393 engine.PerformPuts();
394}
395
396} // namespace impl_fides
397
399enum class FidesMeshPolicy
400{
401 update,
402 reuse
404};
405
409template <std::floating_point T>
410class FidesWriter : public ADIOS2Writer
411{
412public:
421 FidesWriter(MPI_Comm comm, const std::filesystem::path& filename,
422 std::shared_ptr<const mesh::Mesh<T>> mesh,
423 std::string engine = "BPFile")
424 : ADIOS2Writer(comm, filename, "Fides mesh writer", engine),
425 _mesh_reuse_policy(FidesMeshPolicy::update), _mesh(mesh)
426 {
427 assert(_io);
428 assert(mesh);
429 auto topology = mesh->topology();
430 assert(topology);
431 mesh::CellType type = topology->cell_type();
432 if (mesh->geometry().cmap().dim() != mesh::cell_num_entities(type, 0))
433 throw std::runtime_error("Fides only supports lowest-order meshes.");
434 impl_fides::initialize_mesh_attributes(*_io, type);
435 }
436
450 FidesWriter(MPI_Comm comm, const std::filesystem::path& filename,
451 const typename adios2_writer::U<T>& u, std::string engine,
452 const FidesMeshPolicy mesh_policy = FidesMeshPolicy::update)
453 : ADIOS2Writer(comm, filename, "Fides function writer", engine),
454 _mesh_reuse_policy(mesh_policy),
455 _mesh(impl_adios2::extract_common_mesh<T>(u)), _u(u)
456 {
457 if (u.empty())
458 throw std::runtime_error("FidesWriter fem::Function list is empty");
459
460 // Extract space and mesh from first function
461 auto V0 = std::visit([](auto& u) { return u->function_space().get(); },
462 u.front());
463 assert(V0);
464 auto mesh = V0->mesh().get();
465 assert(mesh);
466
467 // Extract element from first function
468 auto element0 = V0->element().get();
469 assert(element0);
470
471 // Check if function is mixed
472 if (element0->is_mixed())
473 {
474 throw std::runtime_error(
475 "Mixed functions are not supported by FidesWriter");
476 }
477
478 // FIXME: is the below check adequate for detecting a
479 // Lagrange element? Check that element is Lagrange
480 if (!element0->interpolation_ident())
481 {
482 throw std::runtime_error("Only Lagrange functions are supported. "
483 "Interpolate Functions before output.");
484 }
485
486 // Raise exception if element is DG 0
487 if (element0->space_dimension() / element0->block_size() == 1)
488 {
489 throw std::runtime_error(
490 "Piecewise constants are not (yet) supported by FidesWriter");
491 }
492
493 // Check that all functions are first order Lagrange
494 mesh::CellType cell_type = mesh->topology()->cell_type();
495 int num_vertices_per_cell = mesh::cell_num_entities(cell_type, 0);
496 for (auto& v : _u)
497 {
498 std::visit(
499 [V0, num_vertices_per_cell, element0](auto&& u)
500 {
501 auto element = u->function_space()->element();
502 assert(element);
503 if (*element != *element0)
504 {
505 throw std::runtime_error("All functions in FidesWriter must have "
506 "the same element type");
507 }
508 auto dof_layout
509 = u->function_space()->dofmap()->element_dof_layout();
510 int num_vertex_dofs = dof_layout.num_entity_dofs(0);
511 int num_dofs = element->space_dimension() / element->block_size();
512 if (num_dofs != num_vertices_per_cell or num_vertex_dofs != 1)
513 {
514 throw std::runtime_error("Only first order Lagrange spaces are "
515 "supported by FidesWriter");
516 }
517
518#ifndef NDEBUG
519 auto dmap0 = V0->dofmap()->map();
520 auto dmap = u->function_space()->dofmap()->map();
521 if (dmap0.size() != dmap.size()
522 or !std::equal(dmap0.data_handle(),
523 dmap0.data_handle() + dmap0.size(),
524 dmap.data_handle()))
525 {
526 throw std::runtime_error(
527 "All functions must have the same dofmap for FideWriter.");
528 }
529#endif
530 },
531 v);
532 }
533
534 auto topology = mesh->topology();
535 assert(topology);
536 mesh::CellType type = topology->cell_type();
537 if (mesh->geometry().cmap().dim() != mesh::cell_num_entities(type, 0))
538 throw std::runtime_error("Fides only supports lowest-order meshes.");
539 impl_fides::initialize_mesh_attributes(*_io, type);
540 impl_fides::initialize_function_attributes<T>(*_io, u);
541 }
542
553 FidesWriter(MPI_Comm comm, const std::filesystem::path& filename,
554 const typename adios2_writer::U<T>& u,
555 const FidesMeshPolicy mesh_policy = FidesMeshPolicy::update)
556 : FidesWriter(comm, filename, u, "BPFile", mesh_policy)
557 {
558 }
559
560 // Copy constructor
561 FidesWriter(const FidesWriter&) = delete;
562
564 FidesWriter(FidesWriter&& file) = default;
565
567 ~FidesWriter() = default;
568
570 FidesWriter& operator=(FidesWriter&&) = default;
571
572 // Copy assignment
573 FidesWriter& operator=(const FidesWriter&) = delete;
574
577 void write(double t)
578 {
579 assert(_io);
580 assert(_engine);
581 _engine->BeginStep();
582 adios2::Variable var_step
583 = impl_adios2::define_variable<double>(*_io, "step");
584 _engine->template Put<double>(var_step, t);
585 if (auto v = _io->template InquireVariable<std::int64_t>("connectivity");
586 !v or _mesh_reuse_policy == FidesMeshPolicy::update)
587 {
588 impl_fides::write_mesh(*_io, *_engine, *_mesh);
589 }
590
591 for (auto& v : _u)
592 {
593 std::visit([this](auto&& u)
594 { impl_fides::write_data(*_io, *_engine, *u); }, v);
595 }
596
597 _engine->EndStep();
598 }
599
600private:
601 // Control whether the mesh is written to file once or at every time
602 // step
603 FidesMeshPolicy _mesh_reuse_policy;
604
605 std::shared_ptr<const mesh::Mesh<T>> _mesh;
606 adios2_writer::U<T> _u;
607};
608
610namespace impl_vtx
611{
614std::stringstream create_vtk_schema(const std::vector<std::string>& point_data,
615 const std::vector<std::string>& cell_data);
616
618template <std::floating_point T>
619std::vector<std::string>
620extract_function_names(const typename adios2_writer::U<T>& u)
621{
622 std::vector<std::string> names;
623 for (auto& v : u)
624 {
625 std::visit(
626 [&names](auto&& u)
627 {
628 using U = std::decay_t<decltype(u)>;
629 using X = typename U::element_type;
630 if constexpr (std::is_floating_point_v<typename X::value_type>)
631 names.push_back(u->name);
632 else
633 {
634 names.push_back(u->name + impl_adios2::field_ext[0]);
635 names.push_back(u->name + impl_adios2::field_ext[1]);
636 }
637 },
638 v);
639 }
640
641 return names;
642}
643
653template <typename T, std::floating_point X>
654void vtx_write_data(adios2::IO& io, adios2::Engine& engine,
655 const fem::Function<T, X>& u)
656{
657 // Get function data array and information about layout
658 assert(u.x());
659 std::span<const T> u_vector = u.x()->array();
660 int rank = u.function_space()->value_shape().size();
661 std::uint32_t num_comp = std::pow(3, rank);
662 std::shared_ptr<const fem::DofMap> dofmap = u.function_space()->dofmap();
663 assert(dofmap);
664 std::shared_ptr<const common::IndexMap> index_map = dofmap->index_map;
665 assert(index_map);
666 int index_map_bs = dofmap->index_map_bs();
667 int dofmap_bs = dofmap->bs();
668 std::uint32_t num_dofs = index_map_bs
669 * (index_map->size_local() + index_map->num_ghosts())
670 / dofmap_bs;
671 if constexpr (std::is_scalar_v<T>)
672 {
673 // ---- Real
674 std::vector<T> data(num_dofs * num_comp, 0);
675 for (std::size_t i = 0; i < num_dofs; ++i)
676 for (int j = 0; j < index_map_bs; ++j)
677 data[i * num_comp + j] = u_vector[i * index_map_bs + j];
678
679 adios2::Variable output = impl_adios2::define_variable<T>(
680 io, u.name, {}, {}, {num_dofs, num_comp});
681 engine.Put(output, data.data(), adios2::Mode::Sync);
682 }
683 else
684 {
685 // ---- Complex
686 using U = typename T::value_type;
687
688 std::vector<U> data(num_dofs * num_comp, 0);
689 for (std::size_t i = 0; i < num_dofs; ++i)
690 for (int j = 0; j < index_map_bs; ++j)
691 data[i * num_comp + j] = std::real(u_vector[i * index_map_bs + j]);
692
693 adios2::Variable output_real = impl_adios2::define_variable<U>(
694 io, u.name + impl_adios2::field_ext[0], {}, {}, {num_dofs, num_comp});
695 engine.Put(output_real, data.data(), adios2::Mode::Sync);
696
697 std::fill(data.begin(), data.end(), 0);
698 for (std::size_t i = 0; i < num_dofs; ++i)
699 for (int j = 0; j < index_map_bs; ++j)
700 data[i * num_comp + j] = std::imag(u_vector[i * index_map_bs + j]);
701 adios2::Variable output_imag = impl_adios2::define_variable<U>(
702 io, u.name + impl_adios2::field_ext[1], {}, {}, {num_dofs, num_comp});
703 engine.Put(output_imag, data.data(), adios2::Mode::Sync);
704 }
705}
706
711template <std::floating_point T>
712void vtx_write_mesh(adios2::IO& io, adios2::Engine& engine,
713 const mesh::Mesh<T>& mesh)
714{
715 const mesh::Geometry<T>& geometry = mesh.geometry();
716 auto topology = mesh.topology();
717 assert(topology);
718
719 // "Put" geometry
720 std::shared_ptr<const common::IndexMap> x_map = geometry.index_map();
721 std::uint32_t num_vertices = x_map->size_local() + x_map->num_ghosts();
722 adios2::Variable local_geometry = impl_adios2::define_variable<T>(
723 io, "geometry", {}, {}, {num_vertices, 3});
724 engine.Put(local_geometry, geometry.x().data());
725
726 // Put number of nodes. The mesh data is written with local indices,
727 // therefore we need the ghost vertices.
728 adios2::Variable vertices = impl_adios2::define_variable<std::uint32_t>(
729 io, "NumberOfNodes", {adios2::LocalValueDim});
730 engine.Put<std::uint32_t>(vertices, num_vertices);
731
732 auto [vtkcells, shape]
733 = io::extract_vtk_connectivity(geometry.dofmap(), topology->cell_type());
734
735 // Add cell metadata
736 int tdim = topology->dim();
737 adios2::Variable cell_var = impl_adios2::define_variable<std::uint32_t>(
738 io, "NumberOfCells", {adios2::LocalValueDim});
739 engine.Put<std::uint32_t>(cell_var, shape[0]);
740 adios2::Variable celltype_var
741 = impl_adios2::define_variable<std::uint32_t>(io, "types");
742 engine.Put<std::uint32_t>(
743 celltype_var, cells::get_vtk_cell_type(topology->cell_type(), tdim));
744
745 // Pack mesh 'nodes'. Output is written as [N0, v0_0,...., v0_N0, N1,
746 // v1_0,...., v1_N1,....], where N is the number of cell nodes and v0,
747 // etc, is the node index
748 std::vector<std::int64_t> cells(shape[0] * (shape[1] + 1), shape[1]);
749 for (std::size_t c = 0; c < shape[0]; ++c)
750 {
751 std::span vtkcell(vtkcells.data() + c * shape[1], shape[1]);
752 std::span cell(cells.data() + c * (shape[1] + 1), shape[1] + 1);
753 std::copy(vtkcell.begin(), vtkcell.end(), std::next(cell.begin()));
754 }
755
756 // Put topology (nodes)
757 adios2::Variable local_topology = impl_adios2::define_variable<std::int64_t>(
758 io, "connectivity", {}, {}, {shape[0], shape[1] + 1});
759 engine.Put(local_topology, cells.data());
760
761 // Vertex global ids and ghost markers
762 adios2::Variable orig_id = impl_adios2::define_variable<std::int64_t>(
763 io, "vtkOriginalPointIds", {}, {}, {num_vertices});
764 engine.Put(orig_id, geometry.input_global_indices().data());
765
766 std::vector<std::uint8_t> x_ghost(num_vertices, 0);
767 std::fill(std::next(x_ghost.begin(), x_map->size_local()), x_ghost.end(), 1);
768 adios2::Variable ghost = impl_adios2::define_variable<std::uint8_t>(
769 io, "vtkGhostType", {}, {}, {x_ghost.size()});
770 engine.Put(ghost, x_ghost.data());
771 engine.PerformPuts();
772}
773
781template <std::floating_point T>
782std::pair<std::vector<std::int64_t>, std::vector<std::uint8_t>>
783vtx_write_mesh_from_space(adios2::IO& io, adios2::Engine& engine,
784 const fem::FunctionSpace<T>& V)
785{
786 auto mesh = V.mesh();
787 assert(mesh);
788 auto topology = mesh->topology();
789 assert(topology);
790 int tdim = topology->dim();
791
792 // Get a VTK mesh with points at the 'nodes'
793 auto [x, xshape, x_id, x_ghost, vtk, vtkshape] = io::vtk_mesh_from_space(V);
794
795 std::uint32_t num_dofs = xshape[0];
796
797 // -- Pack mesh 'nodes'. Output is written as [N0, v0_0,...., v0_N0, N1,
798 // v1_0,...., v1_N1,....], where N is the number of cell nodes and v0,
799 // etc, is the node index.
800
801 // Create vector, setting all entries to nodes per cell (vtk.shape(1))
802 std::vector<std::int64_t> cells(vtkshape[0] * (vtkshape[1] + 1), vtkshape[1]);
803
804 // Set the [v0_0,...., v0_N0, v1_0,...., v1_N1,....] data
805 for (std::size_t c = 0; c < vtkshape[0]; ++c)
806 {
807 std::span vtkcell(vtk.data() + c * vtkshape[1], vtkshape[1]);
808 std::span cell(cells.data() + c * (vtkshape[1] + 1), vtkshape[1] + 1);
809 std::copy(vtkcell.begin(), vtkcell.end(), std::next(cell.begin()));
810 }
811
812 // Define ADIOS2 variables for geometry, topology, celltypes and
813 // corresponding VTK data
814 adios2::Variable local_geometry
815 = impl_adios2::define_variable<T>(io, "geometry", {}, {}, {num_dofs, 3});
816 adios2::Variable local_topology = impl_adios2::define_variable<std::int64_t>(
817 io, "connectivity", {}, {}, {vtkshape[0], vtkshape[1] + 1});
818 adios2::Variable cell_type
819 = impl_adios2::define_variable<std::uint32_t>(io, "types");
820 adios2::Variable vertices = impl_adios2::define_variable<std::uint32_t>(
821 io, "NumberOfNodes", {adios2::LocalValueDim});
822 adios2::Variable elements = impl_adios2::define_variable<std::uint32_t>(
823 io, "NumberOfEntities", {adios2::LocalValueDim});
824
825 // Write mesh information to file
826 engine.Put<std::uint32_t>(vertices, num_dofs);
827 engine.Put<std::uint32_t>(elements, vtkshape[0]);
828 engine.Put<std::uint32_t>(
829 cell_type, cells::get_vtk_cell_type(topology->cell_type(), tdim));
830 engine.Put(local_geometry, x.data());
831 engine.Put(local_topology, cells.data());
832
833 // Node global ids
834 adios2::Variable orig_id = impl_adios2::define_variable<std::int64_t>(
835 io, "vtkOriginalPointIds", {}, {}, {x_id.size()});
836 engine.Put(orig_id, x_id.data());
837 adios2::Variable ghost = impl_adios2::define_variable<std::uint8_t>(
838 io, "vtkGhostType", {}, {}, {x_ghost.size()});
839 engine.Put(ghost, x_ghost.data());
840
841 engine.PerformPuts();
842 return {std::move(x_id), std::move(x_ghost)};
843}
844} // namespace impl_vtx
845
847enum class VTXMeshPolicy
848{
849 update,
850 reuse
852};
853
859template <std::floating_point T>
860class VTXWriter : public ADIOS2Writer
861{
862public:
874 VTXWriter(MPI_Comm comm, const std::filesystem::path& filename,
875 std::shared_ptr<const mesh::Mesh<T>> mesh,
876 std::string engine = "BPFile")
877 : ADIOS2Writer(comm, filename, "VTX mesh writer", engine), _mesh(mesh),
878 _mesh_reuse_policy(VTXMeshPolicy::update), _is_piecewise_constant(false)
879 {
880 // Define VTK scheme attribute for mesh
881 std::string vtk_scheme = impl_vtx::create_vtk_schema({}, {}).str();
882 impl_adios2::define_attribute<std::string>(*_io, "vtk.xml", vtk_scheme);
883 }
884
900 VTXWriter(MPI_Comm comm, const std::filesystem::path& filename,
901 const typename adios2_writer::U<T>& u, std::string engine,
902 VTXMeshPolicy mesh_policy = VTXMeshPolicy::update)
903 : ADIOS2Writer(comm, filename, "VTX function writer", engine),
904 _mesh(impl_adios2::extract_common_mesh<T>(u)),
905 _mesh_reuse_policy(mesh_policy), _u(u), _is_piecewise_constant(false)
906 {
907 if (u.empty())
908 throw std::runtime_error("VTXWriter fem::Function list is empty.");
909
910 // Extract space from first function
911 auto V0 = std::visit([](auto& u) { return u->function_space().get(); },
912 u.front());
913 assert(V0);
914 auto element0 = V0->element().get();
915 assert(element0);
916
917 // Check if function is mixed
918 if (element0->is_mixed())
919 {
920 throw std::runtime_error(
921 "Mixed functions are not supported by VTXWriter.");
922 }
923
924 // FIXME: is the below check adequate for detecting a Lagrange
925 // element?
926 // Check that element is Lagrange
927 if (!element0->interpolation_ident())
928 {
929 throw std::runtime_error(
930 "Only (discontinuous) Lagrange functions are "
931 "supported. Interpolate Functions before output.");
932 }
933
934 // Check if function is DG 0
935 if (element0->space_dimension() / element0->block_size() == 1)
936 _is_piecewise_constant = true;
937
938 // Check that all functions come from same element type
939 for (auto& v : _u)
940 {
941 std::visit(
942 [V0](auto& u)
943 {
944 auto element = u->function_space()->element();
945 assert(element);
946 if (*element != *V0->element().get())
947 {
948 throw std::runtime_error("All functions in VTXWriter must have "
949 "the same element type.");
950 }
951#ifndef NDEBUG
952 auto dmap0 = V0->dofmap()->map();
953 auto dmap = u->function_space()->dofmap()->map();
954 if (dmap0.size() != dmap.size()
955 or !std::equal(dmap0.data_handle(),
956 dmap0.data_handle() + dmap0.size(),
957 dmap.data_handle()))
958 {
959 throw std::runtime_error(
960 "All functions must have the same dofmap for VTXWriter.");
961 }
962#endif
963 },
964 v);
965 }
966
967 // Define VTK scheme attribute for set of functions
968 std::vector<std::string> names = impl_vtx::extract_function_names<T>(u);
969 std::string vtk_scheme;
970 if (_is_piecewise_constant)
971 vtk_scheme = impl_vtx::create_vtk_schema({}, names).str();
972 else
973 vtk_scheme = impl_vtx::create_vtk_schema(names, {}).str();
974
975 impl_adios2::define_attribute<std::string>(*_io, "vtk.xml", vtk_scheme);
976 }
977
992 VTXWriter(MPI_Comm comm, const std::filesystem::path& filename,
993 const typename adios2_writer::U<T>& u,
994 VTXMeshPolicy mesh_policy = VTXMeshPolicy::update)
995 : VTXWriter(comm, filename, u, "BPFile", mesh_policy)
996 {
997 }
998
999 // Copy constructor
1000 VTXWriter(const VTXWriter&) = delete;
1001
1003 VTXWriter(VTXWriter&& file) = default;
1004
1006 ~VTXWriter() = default;
1007
1009 VTXWriter& operator=(VTXWriter&&) = default;
1010
1011 // Copy assignment
1012 VTXWriter& operator=(const VTXWriter&) = delete;
1013
1016 void write(double t)
1017 {
1018 assert(_io);
1019 adios2::Variable var_step
1020 = impl_adios2::define_variable<double>(*_io, "step");
1021
1022 assert(_engine);
1023 _engine->BeginStep();
1024 _engine->template Put<double>(var_step, t);
1025
1026 // If we have no functions or DG functions write the mesh to file
1027 if (_is_piecewise_constant or _u.empty())
1028 {
1029 impl_vtx::vtx_write_mesh(*_io, *_engine, *_mesh);
1030 if (_is_piecewise_constant)
1031 {
1032 for (auto& v : _u)
1033 {
1034 std::visit([&](auto& u)
1035 { impl_vtx::vtx_write_data(*_io, *_engine, *u); }, v);
1036 }
1037 }
1038 }
1039 else
1040 {
1041 if (_mesh_reuse_policy == VTXMeshPolicy::update
1042 or !(_io->template InquireVariable<std::int64_t>("connectivity")))
1043 {
1044 // Write a single mesh for functions as they share finite
1045 // element
1046 std::tie(_x_id, _x_ghost) = std::visit(
1047 [&](auto& u)
1048 {
1049 return impl_vtx::vtx_write_mesh_from_space(*_io, *_engine,
1050 *u->function_space());
1051 },
1052 _u[0]);
1053 }
1054 else
1055 {
1056 // Node global ids
1057 adios2::Variable orig_id = impl_adios2::define_variable<std::int64_t>(
1058 *_io, "vtkOriginalPointIds", {}, {}, {_x_id.size()});
1059 _engine->Put(orig_id, _x_id.data());
1060 adios2::Variable ghost = impl_adios2::define_variable<std::uint8_t>(
1061 *_io, "vtkGhostType", {}, {}, {_x_ghost.size()});
1062 _engine->Put(ghost, _x_ghost.data());
1063 _engine->PerformPuts();
1064 }
1065
1066 // Write function data for each function to file
1067 for (auto& v : _u)
1068 {
1069 std::visit([&](auto& u)
1070 { impl_vtx::vtx_write_data(*_io, *_engine, *u); }, v);
1071 }
1072 }
1073
1074 _engine->EndStep();
1075 }
1076
1077private:
1078 std::shared_ptr<const mesh::Mesh<T>> _mesh;
1079 adios2_writer::U<T> _u;
1080
1081 // Control whether the mesh is written to file once or at every time
1082 // step
1083 VTXMeshPolicy _mesh_reuse_policy;
1084 std::vector<std::int64_t> _x_id;
1085 std::vector<std::uint8_t> _x_ghost;
1086
1087 // Special handling of piecewise constant functions
1088 bool _is_piecewise_constant;
1089};
1090
1092template <typename U, typename T>
1093VTXWriter(MPI_Comm comm, U filename, T mesh)
1094 -> VTXWriter<typename std::remove_cvref<
1095 typename T::element_type>::type::geometry_type::value_type>;
1096
1098template <typename U, typename T>
1099FidesWriter(MPI_Comm comm, U filename, T mesh)
1100 -> FidesWriter<typename std::remove_cvref<
1101 typename T::element_type>::type::geometry_type::value_type>;
1102
1103} // namespace dolfinx::io
1104
1105#endif
std::shared_ptr< const mesh::Mesh< T > > extract_common_mesh(const typename adios2_writer::U< T > &u)
Extract common mesh from list of Functions.
Definition ADIOS2Writers.h:132
std::pair< std::vector< std::int64_t >, std::vector< std::uint8_t > > vtx_write_mesh_from_space(adios2::IO &io, adios2::Engine &engine, const fem::FunctionSpace< T > &V)
Given a FunctionSpace, create a topology and geometry based on the function space dof coordinates....
Definition ADIOS2Writers.h:783
adios2::Variable< T > define_variable(adios2::IO &io, std::string name, const adios2::Dims &shape=adios2::Dims(), const adios2::Dims &start=adios2::Dims(), const adios2::Dims &count=adios2::Dims())
Definition ADIOS2Writers.h:114
void vtx_write_mesh(adios2::IO &io, adios2::Engine &engine, const mesh::Mesh< T > &mesh)
Definition ADIOS2Writers.h:712
adios2::Attribute< T > define_attribute(adios2::IO &io, std::string name, const T &value, std::string var_name="", std::string separator="/")
Definition ADIOS2Writers.h:101
std::vector< std::string > extract_function_names(const typename adios2_writer::U< T > &u)
Extract name of functions and split into real and imaginary component.
Definition ADIOS2Writers.h:620
constexpr std::array field_ext
Definition ADIOS2Writers.h:96
void write_data(adios2::IO &io, adios2::Engine &engine, const fem::Function< T, U > &u)
Definition ADIOS2Writers.h:281
void write_mesh(adios2::IO &io, adios2::Engine &engine, const mesh::Mesh< T > &mesh)
Write mesh geometry and connectivity (topology) for Fides.
Definition ADIOS2Writers.h:364
void vtx_write_data(adios2::IO &io, adios2::Engine &engine, const fem::Function< T, X > &u)
Definition ADIOS2Writers.h:654
std::vector< T > pack_function_data(const fem::Function< T, U > &u)
Definition ADIOS2Writers.h:226
std::stringstream create_vtk_schema(const std::vector< std::string > &point_data, const std::vector< std::string > &cell_data)
Definition ADIOS2Writers.cpp:94
void initialize_function_attributes(adios2::IO &io, const typename adios2_writer::U< T > &u)
Definition ADIOS2Writers.h:174
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
Definition Function.h:45
Base class for ADIOS2-based writers.
Definition ADIOS2Writers.h:55
ADIOS2Writer(MPI_Comm comm, const std::filesystem::path &filename, std::string tag, std::string engine)
Create an ADIOS2-based writer.
Definition ADIOS2Writers.cpp:52
ADIOS2Writer & operator=(ADIOS2Writer &&writer)=default
Move assignment.
void close()
Close the file.
Definition ADIOS2Writers.cpp:64
ADIOS2Writer(const ADIOS2Writer &)=delete
Copy constructor.
ADIOS2Writer(ADIOS2Writer &&writer)=default
Move constructor.
~ADIOS2Writer()
Destructor.
Definition ADIOS2Writers.cpp:62
Geometry stores the geometry imposed on a mesh.
Definition Geometry.h:33
std::shared_ptr< const common::IndexMap > index_map() const
Index map.
Definition Geometry.h:159
const fem::CoordinateElement< value_type > & cmap() const
The element that describes the geometry map.
Definition Geometry.h:180
MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan< const std::int32_t, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents< std::size_t, 2 > > dofmap() const
DofMap for the geometry.
Definition Geometry.h:123
const std::vector< std::int64_t > & input_global_indices() const
Global user indices.
Definition Geometry.h:201
std::span< const value_type > x() const
Access geometry degrees-of-freedom data (const version).
Definition Geometry.h:168
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition Mesh.h:23
Finite element method functionality.
Definition assemble_matrix_impl.h:26
Support for file IO.
Definition ADIOS2Writers.h:41
CellType
Cell type identifier.
Definition cell_types.h:22