Skip to content

Conversation

@ax3l
Copy link
Member

@ax3l ax3l commented Jan 29, 2021

This ports a prior empirical algorithm from libSplash to determine an optimal (large) chunk size for an HDF5 dataset based on its datatype and global extent.

Original implementation by Felix Schmitt @f-schmitt (ZIH, TU Dresden) in libSplash.

Original source:

Close #406
Related to #898 (improve HDF5 baseline performance)
Required for #510: basis to extend resizable data sets (#829) to HDF5

To Do

  • pass tests
  • measure performance (h5chunks_cori.sbatch.txt)
    • Cori seems to be still slow, also when playing further with OPENPMD_HDF5_INDEPENDENT="OFF" and OPENPMD_HDF5_ALIGNMENT="1048576" for our 8_benchmark_parallel -w benchmark
  • measure filesize overhead

Sample & bin directory du -hs bin samples with:

OPENPMD_HDF5_CHUNKS size
"auto"/unset 13M + 1.6G
"none" 11M + 1.6G

on my laptop (4KiB blocksize).

With "auto" the MPI.8_benchmark_parallel test is significantly slower. Changing the 4D test to 3D brings the difference down to about 20% slowdown #1010. (Maybe the 4th, 10-element dimension is sub-ideal for chunking?)

Follow-Ups

  • allow definition on a per-dataset basis
  • allow passing user-defined per-dataset Byte arrays for chunks

@ax3l ax3l requested a review from guj January 29, 2021 23:21
@ax3l ax3l force-pushed the topic-hdf5OptChunkSize branch from 287d7ee to fdbbf4f Compare January 29, 2021 23:23
@ax3l
Copy link
Member Author

ax3l commented Jan 31, 2021

At the moment, a naïve port of the logic causes a deadlock in the MPI tests.

Stuck only in hipace_like_write with

dims: 300 20
optimal chunks found: 256 16

After some tests, it looks like adding chunking makes the dataset declaration H5Dcreate collective in that case, urgh... (It's actually declared as always collective.)

Notes:

Chunk size cannot exceed the size of a fixed-size dataset. For example, a dataset consisting of a 5x4 fixed-size array cannot be defined with 10x10 chunks. 

//for( auto const& val : parameters.chunkSize )
// chunk_dims.push_back(static_cast< hsize_t >(val));

herr_t status = H5Pset_chunk(datasetCreationProperty, chunk_dims.size(), chunk_dims.data());
Copy link
Member Author

@ax3l ax3l Jan 31, 2021

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Fun fact: H5Pset_chunk_opts (HDF5 1.10.0+):

H5Pset_chunk_opts is used to specify storage options for chunks on the edge of a dataset’s dataspace. This capability allows the user to tune performance in cases where the dataset size may not be a multiple of the chunk size and the handling of partial edge chunks can impact performance. 

@ax3l ax3l force-pushed the topic-hdf5OptChunkSize branch 2 times, most recently from a625114 to 83695d8 Compare February 1, 2021 07:43
@ax3l ax3l requested a review from franzpoeschel February 1, 2021 07:43
@ax3l ax3l mentioned this pull request Feb 2, 2021
@ax3l ax3l force-pushed the topic-hdf5OptChunkSize branch from 83695d8 to 04c0164 Compare February 2, 2021 20:15
@ax3l ax3l force-pushed the topic-hdf5OptChunkSize branch from 04c0164 to 44f72a4 Compare February 16, 2021 19:09
@franzpoeschel
Copy link
Contributor

I've added some scaffolding for JSON options in HDF5

@ax3l ax3l force-pushed the topic-hdf5OptChunkSize branch from 0c2a55a to 34cde08 Compare April 6, 2021 06:03
@ax3l ax3l force-pushed the topic-hdf5OptChunkSize branch from 34cde08 to 22f70f0 Compare June 10, 2021 05:03
@ax3l ax3l changed the title [WIP] HDF5: Empiric for Optimal Chunk Size HDF5: Empiric for Optimal Chunk Size Jun 10, 2021
@ax3l
Copy link
Member Author

ax3l commented Jun 10, 2021

Finished the global options (JSON & env) to disable chucking when needed (mainly for HiPACE's legacy pipeline + potential regressions).

@ax3l ax3l force-pushed the topic-hdf5OptChunkSize branch 4 times, most recently from 368c2c4 to aea7e96 Compare June 11, 2021 19:35
@ax3l ax3l force-pushed the topic-hdf5OptChunkSize branch from aea7e96 to 2d34cf6 Compare June 12, 2021 20:39
@ax3l
Copy link
Member Author

ax3l commented Jun 14, 2021

It's a bit concerning that the clang sanitizer run parallel benchmark (8) runs into a time-out with the new patch:

15/32 Test #15: MPI.8_benchmark_parallel ...............***Timeout 1500.17 sec

But I don't see an immediately relation. 'll turn the chunking off for this one.

@ax3l ax3l force-pushed the topic-hdf5OptChunkSize branch from d834f9f to ff1d13d Compare June 14, 2021 23:06
@franzpoeschel
Copy link
Contributor

I tried out whether this PR alone already enables extensible datasets in HDF5, apparently not:

HDF5-DIAG: Error detected in HDF5 (1.10.4) thread 0:
  #000: ../../../src/H5D.c line 906 in H5Dset_extent(): unable to set extend dataset
    major: Dataset
    minor: Unable to initialize object 
  #001: ../../../src/H5Dint.c line 2615 in H5D__set_extent(): unable to modify size of dataspace
    major: Dataset
    minor: Unable to initialize object
  #002: ../../../src/H5S.c line 2005 in H5S_set_extent(): dimension cannot exceed the existing maximal size (new: 10 max: 5)
    major: Dataspace
    minor: Bad value

According to the documentation, only certain kinds of datasets can have their extents extended:

  1. A chunked dataset with unlimited dimensions
  2. A chunked dataset with fixed dimensions if the new dimension sizes are less than the maximum sizes set with maxdims (see H5Screate_simple)

Our datasets seem to fall into the second category. I guess, for extensible datasets, we would have to create unlimited datasets from the beginning?

@ax3l
Copy link
Member Author

ax3l commented Jun 21, 2021

In order to make a dataset resizable, we need to pass H5F_UNLIMITED into the dimensions of H5Screate_simple that should be resizable. Yes, we need chunking (this PR) to be able to use this (ref / PR).

Setting H5F_UNLIMITED on dimensions can have a negative performance implication, so users need to declare that intent upfront and then we just make all dims of a data set resizable.

This ports a prior empirical algorithm from libSplash to determine
an optimal (large) chunk size for an HDF5 dataset based on its
datatype and global extent.

Original implementation by Felix Schmitt @f-schmitt (ZIH, TU Dresden)
in
[libSplash](https://github.com/ComputationalRadiationPhysics/libSplash).

Original source:
- https://github.com/ComputationalRadiationPhysics/libSplash/blob/v1.7.0/src/DCDataSet.cpp
- https://github.com/ComputationalRadiationPhysics/libSplash/blob/v1.7.0/src/include/splash/core/DCHelper.hpp

Co-authored-by: Felix Schmitt <[email protected]>
franzpoeschel and others added 4 commits June 23, 2021 14:19
The parallel, independent I/O pattern here is corner-case for what
HDF5 can support, due to non-collective declarations of data sets.
Testing shows that it does not work with chunking.
Runs into timeout for unclear reasons with this patch:
```
15/32 Test openPMD#15: MPI.8_benchmark_parallel ...............***Timeout 1500.17 sec
```
@ax3l ax3l force-pushed the topic-hdf5OptChunkSize branch from ff1d13d to 16c14c7 Compare June 23, 2021 21:20
@ax3l ax3l merged commit 8c2d9ce into openPMD:dev Jun 24, 2021
@ax3l ax3l deleted the topic-hdf5OptChunkSize branch June 24, 2021 02:31
@ax3l ax3l mentioned this pull request Jun 25, 2021

A full configuration of the HDF5 backend:

.. literalinclude:: hdf5.json
Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@franzpoeschel I just realized I forgot to add a hdf5.json file here 🤪

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Ah, I accidentally named it json.json

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

-> #1169

All keys found under ``hdf5.dataset`` are applicable globally (future: as well as per dataset).
Explanation of the single keys:

* ``adios2.dataset.chunks``: This key contains options for data chunking via `H5Pset_chunk <https://support.hdfgroup.org/HDF5/doc/RM/H5P/H5Pset_chunk.htm>`__.
Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Ouch, that should read hdf5.dataset.chunks...

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

-> #1169

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Projects

None yet

Development

Successfully merging this pull request may close these issues.

HDF5: Chunking

2 participants