ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
h5md_core.cpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2010-2022 The ESPResSo project
3 * Copyright (C) 2002,2003,2004,2005,2006,2007,2008,2009,2010
4 * Max-Planck-Institute for Polymer Research, Theory Group
5 *
6 * This file is part of ESPResSo.
7 *
8 * ESPResSo is free software: you can redistribute it and/or modify
9 * it under the terms of the GNU General Public License as published by
10 * the Free Software Foundation, either version 3 of the License, or
11 * (at your option) any later version.
12 *
13 * ESPResSo is distributed in the hope that it will be useful,
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 * GNU General Public License for more details.
17 *
18 * You should have received a copy of the GNU General Public License
19 * along with this program. If not, see <http://www.gnu.org/licenses/>.
20 */
21
22#include "hdf5_patches.hpp" // must appear first
23
24#include "h5md_core.hpp"
25#include "h5md_dataset.hpp"
27
28#include "BoxGeometry.hpp"
29#include "Particle.hpp"
31
32#include <config/version.hpp>
33
34#include <utils/Vector.hpp>
35
36#include <boost/array.hpp>
37#include <boost/mpi/collectives.hpp>
38#include <boost/multi_array.hpp>
39
40#include <highfive/boost.hpp>
41#include <highfive/highfive.hpp>
42
43#include <mpi.h>
44
45#include <algorithm>
46#include <cstddef>
47#include <filesystem>
48#include <fstream>
49#include <functional>
50#include <iterator>
51#include <memory>
52#include <ranges>
53#include <stdexcept>
54#include <string>
55#include <utility>
56#include <vector>
57
58namespace Writer {
59namespace H5md {
60
61using MultiArray3i = boost::multi_array<int, 3>;
68
69static std::unordered_map<std::string, H5MDOutputFields> const fields_map = {
70 {"all", H5MD_OUT_ALL},
71 {"particle.type", H5MD_OUT_TYPE},
72 {"particle.position", H5MD_OUT_POS},
73 {"particle.image", H5MD_OUT_IMG},
74 {"particle.velocity", H5MD_OUT_VEL},
75 {"particle.force", H5MD_OUT_FORCE},
76 {"particle.bonds", H5MD_OUT_BONDS},
77 {"particle.charge", H5MD_OUT_CHARGE},
78 {"particle.mass", H5MD_OUT_MASS},
79 {"box.length", H5MD_OUT_BOX_L},
80 {"lees_edwards.offset", H5MD_OUT_LE_OFF},
81 {"lees_edwards.direction", H5MD_OUT_LE_DIR},
82 {"lees_edwards.normal", H5MD_OUT_LE_NORMAL},
83};
84
85static auto fields_list_to_bitfield(std::vector<std::string> const &fields) {
86 unsigned int bitfield = H5MD_OUT_NONE;
87 for (auto const &field_name : fields) {
88 if (not fields_map.contains(field_name)) {
89 throw std::invalid_argument("Unknown field '" + field_name + "'");
90 }
91 bitfield |= fields_map.at(field_name);
92 }
93 return bitfield;
94}
95
96static void backup_file(std::filesystem::path const &from,
97 std::filesystem::path const &to) {
98 /*
99 * If the file itself *and* a backup file exists, something must
100 * have gone wrong.
101 */
102 auto constexpr option_fail_if_exists = std::filesystem::copy_options::none;
103 try {
104 std::filesystem::copy_file(from, to, option_fail_if_exists);
105 } catch (std::filesystem::filesystem_error const &) {
106 throw left_backupfile();
107 }
108}
109
110template <typename extent_type>
111static void extend_dataset(HighFive::DataSet &dataset,
112 extent_type const &change_extent) {
113 auto extents = dataset.getSpace().getDimensions();
114 auto const rank = extents.size();
115 /* Extend the dataset for another timestep */
116 for (auto i = 0u; i < rank; i++) {
117 extents[i] += change_extent[i];
118 }
119 dataset.resize(extents);
120}
121
122template <typename value_type, typename extent_type>
123static void write_dataset(value_type const &data, HighFive::DataSet &dataset,
124 extent_type const &offset, extent_type const &count) {
125 auto xfer_props = HighFive::DataTransferProps{};
126 xfer_props.add(HighFive::UseCollectiveIO{});
127 /* write the data to the dataset. */
128 dataset.select(offset, count).write(data, xfer_props);
129}
130
131template <typename value_type, typename extent_type>
132static void write_dataset(value_type const &data, HighFive::DataSet &dataset,
133 extent_type const &change_extent,
134 extent_type const &offset, extent_type const &count) {
135 extend_dataset(dataset, change_extent);
136 /* write the data to the dataset. */
137 write_dataset(data, dataset, offset, count);
138}
139
140static void write_script(HighFive::File &h5md_file,
141 std::filesystem::path const &script_path) {
142 if (!script_path.empty()) {
143 std::ifstream scriptfile(script_path);
144 std::string buffer(std::istreambuf_iterator<char>(scriptfile),
145 std::istreambuf_iterator<char>{});
146 h5md_file.createGroup("/parameters");
147 auto group = h5md_file.createGroup("/parameters/files");
148 group.createAttribute("script", buffer);
149 }
150}
151
152/* Initialize the file-related variables after parameters have been set. */
153void File::init_file() {
154 auto const file_exists = std::filesystem::exists(m_file_path);
155 auto const backup_file_exists = std::filesystem::exists(m_backup_path);
156 /* Perform a barrier synchronization. Otherwise one process might already
157 * create the file while another still checks for its existence. */
158 m_comm.barrier();
159 if (file_exists) {
160 if (m_h5md_specification.is_compliant(m_file_path)) {
161 /*
162 * If the file exists and has a valid H5MD structure, let's create a
163 * backup of it. This has the advantage, that the new file can
164 * just be deleted if the simulation crashes at some point and we
165 * still have a valid trajectory backed up, from which we can restart.
166 */
167 if (m_comm.rank() == 0)
168 backup_file(m_file_path, m_backup_path);
169 load_file();
170 } else {
171 throw incompatible_h5mdfile();
172 }
173 } else {
174 if (backup_file_exists)
175 throw left_backupfile();
176 create_file();
177 }
178}
179
180void File::load_datasets() {
181 auto &datasets = *m_datasets;
182 for (auto const &ds : m_h5md_specification.get_datasets()) {
183 if (ds.is_link)
184 continue;
185 auto path = ds.path();
186 datasets[path] = m_h5md_file->getDataSet(path);
187 }
188}
189
190void File::create_groups() {
191 for (auto const &ds : m_h5md_specification.get_datasets()) {
192 std::stringstream ss(ds.group);
193 std::string segment;
194 std::string current_path = "/";
195 while (std::getline(ss, segment, '/')) {
196 if (segment.empty())
197 continue;
198 current_path += "/" + segment;
199 if (!m_h5md_file->exist(current_path)) {
200 m_h5md_file->createGroup(current_path);
201 }
202 }
203 }
204}
205
206static std::vector<std::size_t> create_dims(hsize_t rank, hsize_t data_dim) {
207 if (rank == 3ul) {
208 return {0ul, 0ul, data_dim};
209 }
210 if (rank == 2ul) {
211 return {0ul, data_dim};
212 }
213 assert(rank == 1ul);
214 return {data_dim};
215}
216
217static std::vector<std::size_t> create_maxdims(hsize_t rank, hsize_t data_dim,
218 hsize_t max_dim) {
219 if (rank == 3ul) {
220 return {max_dim, max_dim, data_dim};
221 }
222 if (rank == 2ul) {
223 return {max_dim, max_dim};
224 }
225 assert(rank == 1ul);
226 return {max_dim};
227}
228
229static std::vector<hsize_t> create_chunk_dims(hsize_t rank, hsize_t data_dim,
230 hsize_t size) {
231 auto const chunk_size = (rank > 1ul) ? size : hsize_t{1ul};
232 if (rank == 3ul) {
233 return {1ul, chunk_size, data_dim};
234 }
235 if (rank == 2ul) {
236 return {1ul, chunk_size};
237 }
238 assert(rank == 1ul);
239 return {chunk_size};
240}
241
242void File::create_datasets() {
243 auto &datasets = *m_datasets;
244 for (auto const &ds : m_h5md_specification.get_datasets()) {
245 if (ds.is_link)
246 continue;
247 auto dims = create_dims(ds.rank, ds.data_dim);
248 auto maxdims = create_maxdims(ds.rank, ds.data_dim, H5S_UNLIMITED);
249 auto dataspace = HighFive::DataSpace(dims, maxdims);
250 auto const chunk_size = static_cast<hsize_t>(m_chunk_size);
251 auto const chunk = create_chunk_dims(ds.rank, ds.data_dim, chunk_size);
252 HighFive::DataSetCreateProps props;
253 props.add(HighFive::Chunking(chunk));
254 auto path = ds.path();
255 if (ds.type == H5T_NATIVE_INT) {
256 datasets.emplace(path,
257 m_h5md_file->createDataSet<int>(path, dataspace, props));
258 } else if (ds.type == H5T_NATIVE_DOUBLE) {
259 datasets.emplace(
260 path, m_h5md_file->createDataSet<double>(path, dataspace, props));
261 }
262 }
263}
264
265void File::load_file() {
266 HighFive::FileAccessProps fapl;
267 fapl.add(HighFive::MPIOFileAccess{m_comm, MPI_INFO_NULL});
268 fapl.add(HighFive::MPIOCollectiveMetadata{});
269 m_h5md_file = std::make_unique<HighFive::File>(
270 m_file_path.string(), HighFive::File::ReadWrite, fapl);
271 load_datasets();
272}
273
274static void write_attributes(HighFive::File &h5md_file) {
275 auto h5md_group = h5md_file.createGroup("h5md");
276 auto att = h5md_group.createAttribute<std::size_t>(
277 "version",
278 HighFive::DataSpace::From(std::array<std::size_t, 2>{{1ul, 1ul}}));
279 att.write(std::array<std::size_t, 2>{{1ul, 1ul}});
280 auto h5md_creator_group = h5md_group.createGroup("creator");
281 h5md_creator_group.createAttribute("name", "ESPResSo");
282 h5md_creator_group.createAttribute("version", ESPRESSO_VERSION);
283 auto h5md_author_group = h5md_group.createGroup("author");
284 h5md_author_group.createAttribute("name", "N/A");
285 auto box_path = "/particles/atoms/box";
286 if (h5md_file.exist(box_path)) {
287 auto group = h5md_file.getGroup(box_path);
288 group.createAttribute("dimension", 3);
289 group.createAttribute("boundary", "periodic");
290 }
291}
292
293void File::write_units() {
294 auto &datasets = *m_datasets;
295 if (!mass_unit().empty() and (m_fields & H5MD_OUT_MASS)) {
296 datasets.at("/particles/atoms/mass/value")
297 .createAttribute("unit", mass_unit());
298 }
299 if (!charge_unit().empty() and (m_fields & H5MD_OUT_CHARGE)) {
300 datasets.at("/particles/atoms/charge/value")
301 .createAttribute("unit", charge_unit());
302 }
303 if (!length_unit().empty() and (m_fields & H5MD_OUT_BOX_L)) {
304 datasets.at("/particles/atoms/position/value")
305 .createAttribute("unit", length_unit());
306 datasets.at("/particles/atoms/box/edges/value")
307 .createAttribute("unit", length_unit());
308 }
309 if (!length_unit().empty() and (m_fields & H5MD_OUT_LE_OFF)) {
310 datasets.at("/particles/atoms/lees_edwards/offset/value")
311 .createAttribute("unit", length_unit());
312 }
313 if (!velocity_unit().empty() and (m_fields & H5MD_OUT_VEL)) {
314 datasets.at("/particles/atoms/velocity/value")
315 .createAttribute("unit", velocity_unit());
316 }
317 if (!force_unit().empty() and (m_fields & H5MD_OUT_FORCE)) {
318 datasets.at("/particles/atoms/force/value")
319 .createAttribute("unit", force_unit());
320 }
321 if (!time_unit().empty()) {
322 datasets.at("/particles/atoms/id/time")
323 .createAttribute("unit", time_unit());
324 }
325}
326
327void File::create_hard_links() {
328 std::string path_step = "/particles/atoms/id/step";
329 std::string path_time = "/particles/atoms/id/time";
330 for (auto &ds : m_h5md_specification.get_datasets()) {
331 if (ds.is_link) {
332 char const *from = nullptr;
333 if (ds.name == "step") {
334 from = path_step.c_str();
335 } else if (ds.name == "time") {
336 from = path_time.c_str();
337 }
338 assert(from != nullptr);
339 if (H5Lcreate_hard(m_h5md_file->getId(), from, m_h5md_file->getId(),
340 ds.path().c_str(), H5P_DEFAULT, H5P_DEFAULT) < 0) {
341 throw std::runtime_error("Error creating hard link for " + ds.path());
342 }
343 }
344 }
345}
346
347void File::create_file() {
348 HighFive::FileAccessProps fapl;
349 fapl.add(HighFive::MPIOFileAccess{m_comm, MPI_INFO_NULL});
350 fapl.add(HighFive::MPIOCollectiveMetadata{});
351 m_h5md_file = std::make_unique<HighFive::File>(m_file_path.string(),
352 HighFive::File::Create, fapl);
353 write_script(*m_h5md_file, m_absolute_script_path);
354 create_groups();
355 create_datasets();
356 write_attributes(*m_h5md_file);
357 write_units();
358 create_hard_links();
359}
360
362 if (m_comm.rank() == 0) {
363 std::filesystem::remove(m_backup_path);
364 }
365}
366
367namespace detail {
368
369template <std::size_t rank> struct slice_info {};
370
371template <> struct slice_info<3> {
372 static auto extent(hsize_t n_part_diff) {
373 return Vector3s{1ul, n_part_diff, 0ul};
374 }
375 static constexpr auto count(std::size_t local_n_part) {
376 return Vector3s{1ul, local_n_part, 3ul};
377 }
378 static auto offset(hsize_t n_time_steps, hsize_t prefix) {
379 return Vector3s{n_time_steps, prefix, 0ul};
380 }
381 template <typename T>
382 static boost::multi_array<T, 3> reshape(std::vector<T> const &v1d,
383 Vector3s const &count) {
384 if (v1d.empty()) {
385 boost::multi_array<T, 3> data(boost::extents[0][0][0]);
386 return data;
387 }
388 auto const rows = count[1];
389 auto const cols = count[2];
390
391 boost::multi_array<T, 3> data(
392 boost::extents[1][static_cast<long>(rows)][static_cast<long>(cols)]);
393
394 for (std::size_t i = 0; i < rows; ++i) {
395 for (std::size_t j = 0; j < cols; ++j) {
396 data[0][i][j] = v1d[cols * i + j];
397 }
398 }
399
400 return data;
401 }
402};
403
404template <> struct slice_info<2> {
405 static auto extent(hsize_t n_part_diff) { return Vector2s{1ul, n_part_diff}; }
406 static constexpr auto count(std::size_t local_n) {
407 return Vector2s{1ul, local_n};
408 }
409 static auto offset(hsize_t n_time_steps, hsize_t prefix) {
410 return Vector2s{n_time_steps, prefix};
411 }
412 template <typename T>
413 static boost::multi_array<T, 2> reshape(std::vector<T> const &v1d,
414 Vector2s const &count) {
415 if (v1d.empty()) {
416 boost::multi_array<T, 2> data(boost::extents[0][0]);
417 return data;
418 }
419 auto const cols = count[1];
420
421 boost::multi_array<T, 2> data(boost::extents[1][static_cast<long>(cols)]);
422
423 for (std::size_t i = 0; i < cols; ++i) {
424 data[0][i] = v1d[i];
425 }
426
427 return data;
428 }
429};
430
431template <typename T> struct get_buffer_traits {};
432
433template <typename T>
434 requires std::is_arithmetic_v<T>
435struct get_buffer_traits<T> {
436 using type = T;
437 constexpr static std::size_t dim = 1ul;
438};
439
440template <typename T, std::size_t N>
441 requires std::is_arithmetic_v<T>
442struct get_buffer_traits<Utils::Vector<T, N>> {
443 using type = T;
444 constexpr static std::size_t dim = N;
445};
446
447template <typename Functor> class ParticleDataSerializer {
448 using RetVal = std::decay_t<std::invoke_result_t<Functor, Particle const &>>;
449 Functor m_getter;
450
451 template <typename T>
452 requires std::is_arithmetic_v<T>
453 void serialize(auto &buffer, T const &value) const {
454 buffer.emplace_back(value);
455 }
456
457 template <typename T, std::size_t N>
458 void serialize(auto &buffer, Utils::Vector<T, N> const &value) const {
459 buffer.insert(buffer.end(), value.cbegin(), value.cend());
460 }
461
462public:
463 explicit ParticleDataSerializer(Functor lambda) : m_getter{lambda} {}
464
465 auto operator()(ParticleRange const &particles) const {
466 auto constexpr value_dim = get_buffer_traits<RetVal>::dim;
467 std::vector<typename get_buffer_traits<RetVal>::type> buffer{};
468 buffer.reserve(particles.size() * value_dim);
469 for (auto const &p : particles) {
470 serialize(buffer, m_getter(p));
471 }
472 return buffer;
473 }
474};
475
476template <typename Functor> auto make_serializer(Functor lambda) {
477 return ParticleDataSerializer<Functor>{lambda};
478}
479template <typename RetVal>
480auto make_serializer(RetVal (Particle::*getter)() const) {
481 auto kernel = [getter](Particle const &p) -> RetVal { return (p.*getter)(); };
482 return ParticleDataSerializer<decltype(kernel)>{std::move(kernel)};
483}
484
485} // namespace detail
486
487template <std::size_t dim, typename Serializer>
488void write_td_particle_property(hsize_t prefix, hsize_t n_part_global,
489 ParticleRange const &particles,
490 HighFive::DataSet &dataset,
491 Serializer serializer) {
492 auto const n_part_local = static_cast<hsize_t>(particles.size());
493 auto const old_extents = dataset.getSpace().getDimensions();
494 auto const extent_n_part =
495 std::max(n_part_global, static_cast<hsize_t>(old_extents[1])) -
496 old_extents[1];
497 extend_dataset(dataset, detail::slice_info<dim>::extent(extent_n_part));
498 auto const count = detail::slice_info<dim>::count(n_part_local);
499 auto const offset = detail::slice_info<dim>::offset(old_extents[0], prefix);
500 HighFive::DataType dtype = dataset.getDataType();
501 auto buffer = serializer(particles);
502 write_dataset(detail::slice_info<dim>::reshape(buffer, count), dataset,
503 offset, count);
504}
505
506static void write_box(BoxGeometry const &box_geo, HighFive::DataSet &dataset) {
507 auto const extents = dataset.getSpace().getDimensions();
508 extend_dataset(dataset, Vector2hs{1ul, 0ul});
509 Vector2s const offset{extents[0], 0ul};
510 Vector2s const count{1ul, 3ul};
511 auto const data = box_geo.length().as_vector();
512 write_dataset(detail::slice_info<2>::reshape(data, count), dataset, offset,
513 count);
514}
515
516static void write_le_off(LeesEdwardsBC const &lebc,
517 HighFive::DataSet &dataset) {
518 auto const extents = dataset.getSpace().getDimensions();
519 extend_dataset(dataset, Vector2hs{1ul, 0ul});
520 Vector2s const offset{extents[0], 0ul};
521 Vector2s const count{1ul, 1ul};
522 auto const data = std::vector<double>{lebc.pos_offset};
523 write_dataset(detail::slice_info<2>::reshape(data, count), dataset, offset,
524 count);
525}
526
527static void write_le_dir(LeesEdwardsBC const &lebc,
528 HighFive::DataSet &dataset) {
529 auto const shear_direction = static_cast<int>(lebc.shear_direction);
530 auto const extents = dataset.getSpace().getDimensions();
531 extend_dataset(dataset, Vector2hs{1ul, 0ul});
532 Vector2s const offset{extents[0], 0ul};
533 Vector2s const count{1ul, 1ul};
534 auto const data = std::vector<int>{shear_direction};
535 write_dataset(detail::slice_info<2>::reshape(data, count), dataset, offset,
536 count);
537}
538
539static void write_le_normal(LeesEdwardsBC const &lebc,
540 HighFive::DataSet &dataset) {
541 auto const shear_plane_normal = static_cast<int>(lebc.shear_plane_normal);
542 auto const extents = dataset.getSpace().getDimensions();
543 extend_dataset(dataset, Vector2hs{1ul, 0ul});
544 Vector2s const offset{extents[0], 0ul};
545 Vector2s const count{1ul, 1ul};
546 auto const data = std::vector<int>{shear_plane_normal};
547 write_dataset(detail::slice_info<2>::reshape(data, count), dataset, offset,
548 count);
549}
550
551void File::write(const ParticleRange &particles, double time, int step,
552 BoxGeometry const &box_geo) {
553 auto &datasets = *m_datasets;
554 if (m_fields & H5MD_OUT_BOX_L) {
555 write_box(box_geo, datasets.at("/particles/atoms/box/edges/value"));
556 }
557 auto const &lebc = box_geo.lees_edwards_bc();
558 if (m_fields & H5MD_OUT_LE_OFF) {
559 write_le_off(lebc,
560 datasets.at("/particles/atoms/lees_edwards/offset/value"));
561 }
562 if (m_fields & H5MD_OUT_LE_DIR) {
563 write_le_dir(lebc,
564 datasets.at("/particles/atoms/lees_edwards/direction/value"));
565 }
566 if (m_fields & H5MD_OUT_LE_NORMAL) {
567 write_le_normal(lebc,
568 datasets.at("/particles/atoms/lees_edwards/normal/value"));
569 }
570
571 // calculate particle count and offset
572 static_assert(sizeof(hsize_t) == 8ul);
573 auto const n_part_local = static_cast<hsize_t>(particles.size());
574 hsize_t prefix{0ul};
575 BOOST_MPI_CHECK_RESULT(
576 MPI_Exscan, (&n_part_local, &prefix, 1, MPI_UINT64_T, MPI_SUM, m_comm));
577 auto const n_part_global =
578 boost::mpi::all_reduce(m_comm, n_part_local, std::plus<hsize_t>());
579
580 write_td_particle_property<2>(prefix, n_part_global, particles,
581 datasets.at("/particles/atoms/id/value"),
582 detail::make_serializer(&Particle::id));
583
584 {
585 HighFive::DataSet &dataset = datasets.at("/particles/atoms/id/value");
586 auto const extents = dataset.getSpace().getDimensions();
587 write_dataset(std::vector<double>{time},
588 datasets.at("/particles/atoms/id/time"), Vector1s{1},
589 Vector1s{extents[0]}, Vector1s{1});
590 write_dataset(std::vector<int>{step},
591 datasets.at("/particles/atoms/id/step"), Vector1s{1},
592 Vector1s{extents[0]}, Vector1s{1});
593 }
594
595 if (m_fields & H5MD_OUT_TYPE) {
596 write_td_particle_property<2>(prefix, n_part_global, particles,
597 datasets.at("/particles/atoms/species/value"),
598 detail::make_serializer(&Particle::type));
599 }
600 if (m_fields & H5MD_OUT_MASS) {
601 write_td_particle_property<2>(prefix, n_part_global, particles,
602 datasets.at("/particles/atoms/mass/value"),
603 detail::make_serializer(&Particle::mass));
604 }
605 if (m_fields & H5MD_OUT_POS) {
606 write_td_particle_property<3>(
607 prefix, n_part_global, particles,
608 datasets.at("/particles/atoms/position/value"),
609 detail::make_serializer([&](Particle const &p) {
610 return box_geo.folded_position(p.pos());
611 }));
612 }
613 if (m_fields & H5MD_OUT_IMG) {
614 write_td_particle_property<3>(
615 prefix, n_part_global, particles,
616 datasets.at("/particles/atoms/image/value"),
617 detail::make_serializer([&](Particle const &p) {
618 return box_geo.folded_image_box(p.pos(), p.image_box());
619 }));
620 }
621 if (m_fields & H5MD_OUT_VEL) {
622 write_td_particle_property<3>(
623 prefix, n_part_global, particles,
624 datasets.at("/particles/atoms/velocity/value"),
625 detail::make_serializer(&Particle::v));
626 }
627 if (m_fields & H5MD_OUT_FORCE) {
628 write_td_particle_property<3>(prefix, n_part_global, particles,
629 datasets.at("/particles/atoms/force/value"),
630 detail::make_serializer(&Particle::force));
631 }
632 if (m_fields & H5MD_OUT_CHARGE) {
633 write_td_particle_property<2>(prefix, n_part_global, particles,
634 datasets.at("/particles/atoms/charge/value"),
635 detail::make_serializer(&Particle::q));
636 }
637 if (m_fields & H5MD_OUT_BONDS) {
638 write_connectivity(particles);
639 }
640}
641
642void File::write_connectivity(const ParticleRange &particles) {
643 MultiArray3i bond(boost::extents[0][0][0]);
644 for (auto const &p : particles) {
645 auto nbonds_local = static_cast<decltype(bond)::index>(bond.shape()[1]);
646 for (auto const b : p.bonds()) {
647 auto const &partner_ids = b.partner_ids();
648 if (partner_ids.size() == 1u) {
649 bond.resize(boost::extents[1][nbonds_local + 1][2]);
650 bond[0][nbonds_local][0] = p.id();
651 bond[0][nbonds_local][1] = partner_ids[0];
652 nbonds_local++;
653 }
654 }
655 }
656
657 auto const n_bonds_local = static_cast<int>(bond.shape()[1]);
658 auto &datasets = *m_datasets;
659 int prefix_bonds = 0;
660 BOOST_MPI_CHECK_RESULT(
661 MPI_Exscan, (&n_bonds_local, &prefix_bonds, 1, MPI_INT, MPI_SUM, m_comm));
662 auto const n_bonds_total =
663 boost::mpi::all_reduce(m_comm, n_bonds_local, std::plus<int>());
664 auto const extents =
665 datasets.at("/connectivity/atoms/value").getSpace().getDimensions();
666 Vector3s offset_bonds = {extents[0], static_cast<std::size_t>(prefix_bonds),
667 0};
668 Vector3s count_bonds = {1, static_cast<std::size_t>(n_bonds_local), 2};
669 auto const n_bond_diff = std::max(static_cast<hsize_t>(n_bonds_total),
670 static_cast<hsize_t>(extents[1])) -
671 extents[1];
672 Vector3s change_extent_bonds = {1, static_cast<std::size_t>(n_bond_diff), 0};
673 write_dataset(bond, datasets.at("/connectivity/atoms/value"),
674 change_extent_bonds, offset_bonds, count_bonds);
675}
676
677void File::flush() { m_h5md_file->flush(); }
678
679std::vector<std::string> File::valid_fields() const {
680 auto const view = std::views::elements<0>(fields_map);
681 return {view.begin(), view.end()};
682}
683
684File::File(std::filesystem::path file_path, std::filesystem::path script_path,
685 std::vector<std::string> const &output_fields, std::string mass_unit,
686 std::string length_unit, std::string time_unit,
687 std::string force_unit, std::string velocity_unit,
688 std::string charge_unit, int chunk_size)
689 : m_file_path(std::move(file_path)), m_backup_path(m_file_path),
690 m_script_path(std::move(script_path)),
691 m_absolute_script_path(
692 m_script_path.empty()
693 ? std::filesystem::path()
694 : std::filesystem::weakly_canonical(m_script_path)),
695 m_mass_unit(std::move(mass_unit)), m_length_unit(std::move(length_unit)),
696 m_time_unit(std::move(time_unit)), m_force_unit(std::move(force_unit)),
697 m_velocity_unit(std::move(velocity_unit)),
698 m_charge_unit(std::move(charge_unit)), m_chunk_size(chunk_size),
699 m_comm(boost::mpi::communicator()),
700 m_fields(fields_list_to_bitfield(output_fields)),
701 m_datasets(std::make_unique<decltype(m_datasets)::element_type>()),
702 m_h5md_specification(m_fields) {
703 if (chunk_size <= 0) {
704 throw std::domain_error("Parameter 'chunk_size' must be > 0");
705 }
706 m_backup_path += ".bak";
707 init_file();
708}
709
710File::~File() = default;
711
712} /* namespace H5md */
713} /* namespace Writer */
Vector implementation and trait types for boost qvm interoperability.
Utils::Vector3d const & length() const
Box length.
LeesEdwardsBC const & lees_edwards_bc() const
A range of particles.
base_type::size_type size() const
DEVICE_QUALIFIER constexpr const_iterator cbegin() const noexcept
Definition Array.hpp:148
std::vector< T > as_vector() const
Definition Vector.hpp:140
DEVICE_QUALIFIER constexpr const_iterator cend() const noexcept
Definition Array.hpp:160
void write(const ParticleRange &particles, double time, int step, BoxGeometry const &geometry)
Write data to the hdf5 file.
auto const & chunk_size() const
Retrieve the set chunk size.
auto const & length_unit() const
Retrieve the set length unit.
auto const & time_unit() const
Retrieve the set time unit.
void close()
Method to perform the renaming of the temporary file from "filename" + ".bak" to "filename".
auto const & force_unit() const
Retrieve the set force unit.
auto const & mass_unit() const
Retrieve the set mass unit.
auto const & charge_unit() const
Retrieve the set charge unit.
auto const & velocity_unit() const
Retrieve the set velocity unit.
std::vector< std::string > valid_fields() const
Build the list of valid output fields.
File(std::filesystem::path file_path, std::filesystem::path script_path, std::vector< std::string > const &output_fields, std::string mass_unit, std::string length_unit, std::string time_unit, std::string force_unit, std::string velocity_unit, std::string charge_unit, int chunk_size)
Constructor.
void flush()
Method to enforce flushing the buffer to disk.
Communicator communicator
ParticleRange particles(std::span< Cell *const > cells)
boost::multi_array< int, 3 > MultiArray3i
Definition h5md_core.cpp:61
static void write_le_normal(LeesEdwardsBC const &lebc, HighFive::DataSet &dataset)
static std::vector< std::size_t > create_maxdims(hsize_t rank, hsize_t data_dim, hsize_t max_dim)
static auto fields_list_to_bitfield(std::vector< std::string > const &fields)
Definition h5md_core.cpp:85
static void write_box(BoxGeometry const &box_geo, HighFive::DataSet &dataset)
static void write_script(HighFive::File &h5md_file, std::filesystem::path const &script_path)
static std::unordered_map< std::string, H5MDOutputFields > const fields_map
Definition h5md_core.cpp:69
static void write_dataset(value_type const &data, HighFive::DataSet &dataset, extent_type const &offset, extent_type const &count)
static std::vector< hsize_t > create_chunk_dims(hsize_t rank, hsize_t data_dim, hsize_t size)
Utils::Vector< std::size_t, 2 > Vector2s
Definition h5md_core.cpp:66
static void extend_dataset(HighFive::DataSet &dataset, extent_type const &change_extent)
static void write_le_off(LeesEdwardsBC const &lebc, HighFive::DataSet &dataset)
static void backup_file(std::filesystem::path const &from, std::filesystem::path const &to)
Definition h5md_core.cpp:96
static void write_le_dir(LeesEdwardsBC const &lebc, HighFive::DataSet &dataset)
static void write_attributes(HighFive::File &h5md_file)
Utils::Vector< std::size_t, 3 > Vector3s
Definition h5md_core.cpp:67
void write_td_particle_property(hsize_t prefix, hsize_t n_part_global, ParticleRange const &particles, HighFive::DataSet &dataset, Serializer serializer)
static std::vector< std::size_t > create_dims(hsize_t rank, hsize_t data_dim)
void serialize(Archive &ar, std::tuple< T... > &pack, unsigned int const)
Serialize std::tuple.
STL namespace.
unsigned int shear_plane_normal
unsigned int shear_direction
Struct holding all information for one particle.
Definition Particle.hpp:450
auto const & mass() const
Definition Particle.hpp:507
auto const & q() const
Definition Particle.hpp:593
auto const & v() const
Definition Particle.hpp:488
auto const & type() const
Definition Particle.hpp:473
auto const & force() const
Definition Particle.hpp:490
auto const & id() const
Definition Particle.hpp:469
bool is_compliant(std::filesystem::path const &file) const