32#include <config/version.hpp>
36#include <boost/array.hpp>
37#include <boost/mpi/collectives.hpp>
38#include <boost/multi_array.hpp>
40#if defined(__GNUC__) or defined(__GNUG__)
42#pragma GCC diagnostic push
43#pragma GCC diagnostic ignored "-Wuninitialized"
45#if __has_include(<highfive/boost.hpp>)
46#include <highfive/boost.hpp>
48#include <highfive/highfive.hpp>
49#if defined(__GNUC__) or defined(__GNUG__)
50#pragma GCC diagnostic pop
79static std::unordered_map<std::string, H5MDOutputFields>
const fields_map = {
99 throw std::invalid_argument(
"Unknown field '" +
field_name +
"'");
107 std::filesystem::path
const &
to) {
115 }
catch (std::filesystem::filesystem_error
const &) {
120template <
typename extent_type>
124 auto const rank =
extents.size();
126 for (
auto i = 0
u; i < rank; i++) {
132template <
typename value_type,
typename extent_type>
135 auto xfer_props = HighFive::DataTransferProps{};
141template <
typename value_type,
typename extent_type>
151 std::filesystem::path
const &script_path) {
152 if (!script_path.empty()) {
154 std::string buffer(std::istreambuf_iterator<char>(
scriptfile),
155 std::istreambuf_iterator<char>{});
157 auto group =
h5md_file.createGroup(
"/parameters/files");
158 group.createAttribute(
"script", buffer);
163void File::init_file() {
164 auto const file_exists = std::filesystem::exists(m_file_path);
177 if (m_comm.rank() == 0)
181 throw incompatible_h5mdfile();
185 throw left_backupfile();
190void File::load_datasets() {
192 for (
auto const &
ds : m_h5md_specification.get_datasets()) {
195 auto path =
ds.path();
196 datasets[path] = m_h5md_file->getDataSet(path);
200void File::create_groups() {
201 for (
auto const &
ds : m_h5md_specification.get_datasets()) {
202 std::stringstream
ss(
ds.group);
218 return {0
ul, 0
ul, data_dim};
221 return {0
ul, data_dim};
241 auto const chunk_size = (rank > 1ul) ? size :
hsize_t{1ul};
243 return {1ul, chunk_size, data_dim};
246 return {1ul, chunk_size};
252void File::create_datasets() {
254 for (
auto const &
ds : m_h5md_specification.get_datasets()) {
262 HighFive::DataSetCreateProps
props;
264 auto path =
ds.path();
270 path, m_h5md_file->createDataSet<
double>(path,
dataspace,
props));
275void File::load_file() {
276 HighFive::FileAccessProps
fapl;
278 fapl.add(HighFive::MPIOCollectiveMetadata{});
279 m_h5md_file = std::make_unique<HighFive::File>(
280 m_file_path.string(), HighFive::File::ReadWrite,
fapl);
288 HighFive::DataSpace::From(std::array<std::size_t, 2>{{1ul, 1ul}}));
289 att.write(std::array<std::size_t, 2>{{1ul, 1ul}});
295 auto box_path =
"/particles/atoms/box";
298 group.createAttribute(
"dimension", 3);
299 group.createAttribute(
"boundary",
"periodic");
303void File::write_units() {
306 datasets.at(
"/particles/atoms/mass/value")
310 datasets.at(
"/particles/atoms/charge/value")
314 datasets.at(
"/particles/atoms/position/value")
316 datasets.at(
"/particles/atoms/box/edges/value")
320 datasets.at(
"/particles/atoms/lees_edwards/offset/value")
324 datasets.at(
"/particles/atoms/velocity/value")
328 datasets.at(
"/particles/atoms/force/value")
332 datasets.at(
"/particles/atoms/id/time")
337void File::create_hard_links() {
338 std::string
path_step =
"/particles/atoms/id/step";
339 std::string
path_time =
"/particles/atoms/id/time";
340 for (
auto &
ds : m_h5md_specification.get_datasets()) {
342 char const *
from =
nullptr;
343 if (
ds.name ==
"step") {
345 }
else if (
ds.name ==
"time") {
351 throw std::runtime_error(
"Error creating hard link for " +
ds.path());
357void File::create_file() {
358 HighFive::FileAccessProps
fapl;
360 fapl.add(HighFive::MPIOCollectiveMetadata{});
361 m_h5md_file = std::make_unique<HighFive::File>(m_file_path.string(),
362 HighFive::File::Create,
fapl);
372 if (m_comm.rank() == 0) {
373 std::filesystem::remove(m_backup_path);
379template <std::
size_t rank>
struct slice_info {};
381template <>
struct slice_info<3> {
391 template <
typename T>
392 static boost::multi_array<T, 3> reshape(std::vector<T>
const &
v1d,
395 boost::multi_array<T, 3> data(boost::extents[0][0][0]);
398 auto const rows = count[1];
399 auto const cols = count[2];
401 boost::multi_array<T, 3> data(
402 boost::extents[1][
static_cast<long>(rows)][
static_cast<long>(cols)]);
404 for (std::size_t i = 0; i < rows; ++i) {
405 for (std::size_t
j = 0;
j < cols; ++
j) {
406 data[0][i][
j] =
v1d[cols * i +
j];
414template <>
struct slice_info<2> {
416 static constexpr auto count(std::size_t
local_n) {
422 template <
typename T>
423 static boost::multi_array<T, 2> reshape(std::vector<T>
const &
v1d,
426 boost::multi_array<T, 2> data(boost::extents[0][0]);
429 auto const cols = count[1];
431 boost::multi_array<T, 2> data(boost::extents[1][
static_cast<long>(cols)]);
433 for (std::size_t i = 0; i < cols; ++i) {
441template <
typename T>
struct get_buffer_traits {};
444 requires std::is_arithmetic_v<T>
445struct get_buffer_traits<T> {
447 constexpr static std::size_t dim = 1ul;
450template <
typename T, std::
size_t N>
451 requires std::is_arithmetic_v<T>
452struct get_buffer_traits<
Utils::Vector<T, N>> {
454 constexpr static std::size_t dim =
N;
457template <
typename Functor>
class ParticleDataSerializer {
458 using RetVal = std::decay_t<std::invoke_result_t<Functor, Particle const &>>;
461 template <
typename T>
462 requires std::is_arithmetic_v<T>
463 void serialize(
auto &buffer, T
const &value)
const {
464 buffer.emplace_back(value);
467 template <
typename T, std::
size_t N>
469 buffer.insert(buffer.end(), value.
cbegin(), value.
cend());
473 explicit ParticleDataSerializer(
Functor lambda) : m_getter{lambda} {}
476 auto constexpr value_dim = get_buffer_traits<RetVal>::dim;
477 std::vector<typename get_buffer_traits<RetVal>::type> buffer{};
486template <
typename Functor>
auto make_serializer(
Functor lambda) {
489template <
typename RetVal>
492 return ParticleDataSerializer<
decltype(kernel)>{std::move(kernel)};
497template <std::
size_t dim,
typename Serializer>
508 auto const count = detail::slice_info<dim>::count(
n_part_local);
532 auto const data = std::vector<double>{
lebc.pos_offset};
539 auto const shear_direction =
static_cast<int>(
lebc.shear_direction);
544 auto const data = std::vector<int>{shear_direction};
551 auto const shear_plane_normal =
static_cast<int>(
lebc.shear_plane_normal);
556 auto const data = std::vector<int>{shear_plane_normal};
570 datasets.at(
"/particles/atoms/lees_edwards/offset/value"));
574 datasets.at(
"/particles/atoms/lees_edwards/direction/value"));
578 datasets.at(
"/particles/atoms/lees_edwards/normal/value"));
582 static_assert(
sizeof(
hsize_t) == 8ul);
588 boost::mpi::all_reduce(m_comm,
n_part_local, std::plus<hsize_t>());
591 datasets.at(
"/particles/atoms/id/value"),
607 datasets.at(
"/particles/atoms/species/value"),
612 datasets.at(
"/particles/atoms/mass/value"),
618 datasets.at(
"/particles/atoms/position/value"),
619 detail::make_serializer([&](
Particle const &p) {
620 return box_geo.folded_position(p.pos());
626 datasets.at(
"/particles/atoms/image/value"),
627 detail::make_serializer([&](
Particle const &p) {
628 return box_geo.folded_image_box(p.pos(), p.image_box());
634 datasets.at(
"/particles/atoms/velocity/value"),
639 datasets.at(
"/particles/atoms/force/value"),
644 datasets.at(
"/particles/atoms/charge/value"),
648 write_connectivity(particles);
652void File::write_connectivity(
const ParticleRange &particles) {
654 for (
auto const &p : particles) {
656 for (
auto const b : p.bonds()) {
657 auto const &partner_ids = b.partner_ids();
658 if (partner_ids.size() == 1u) {
673 boost::mpi::all_reduce(m_comm,
n_bonds_local, std::plus<int>());
675 datasets.at(
"/connectivity/atoms/value").getSpace().getDimensions();
694File::File(std::filesystem::path file_path, std::filesystem::path script_path,
695 std::vector<std::string>
const &
output_fields, std::string mass_unit,
696 std::string length_unit, std::string time_unit,
697 std::string force_unit, std::string velocity_unit,
698 std::string charge_unit,
int chunk_size)
699 : m_file_path(
std::
move(file_path)), m_backup_path(m_file_path),
700 m_script_path(
std::
move(script_path)),
701 m_absolute_script_path(
702 m_script_path.empty()
705 m_mass_unit(
std::
move(mass_unit)), m_length_unit(
std::
move(length_unit)),
706 m_time_unit(
std::
move(time_unit)), m_force_unit(
std::
move(force_unit)),
707 m_velocity_unit(
std::
move(velocity_unit)),
708 m_charge_unit(
std::
move(charge_unit)), m_chunk_size(chunk_size),
711 m_h5md_specification(m_fields) {
713 throw std::domain_error(
"Parameter 'chunk_size' must be > 0");
715 m_backup_path +=
".bak";
Vector implementation and trait types for boost qvm interoperability.
Utils::Vector3d const & length() const
Box length.
LeesEdwardsBC const & lees_edwards_bc() const
base_type::size_type size() const
DEVICE_QUALIFIER constexpr const_iterator cbegin() const noexcept
std::vector< T > as_vector() const
DEVICE_QUALIFIER constexpr const_iterator cend() const noexcept
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.
cudaStream_t stream[1]
CUDA streams for parallel computing on CPU and GPU.
Communicator communicator
ParticleRange particles(std::span< Cell *const > cells)
boost::multi_array< int, 3 > MultiArray3i
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)
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
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
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)
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
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.
Struct holding all information for one particle.
auto const & mass() const
auto const & type() const
auto const & force() const
bool is_compliant(std::filesystem::path const &file) const