Skip to content

Commit

Permalink
Merge pull request #644 from brtnfld/hdf5.opt
Browse files Browse the repository at this point in the history
Adds HDF5 option to use subfiling
  • Loading branch information
streeve authored Jan 12, 2024
2 parents 6fd0c34 + 07df756 commit 9e099a6
Show file tree
Hide file tree
Showing 6 changed files with 474 additions and 7 deletions.
33 changes: 33 additions & 0 deletions .github/workflows/CI.yml
Original file line number Diff line number Diff line change
Expand Up @@ -345,6 +345,39 @@ jobs:
${hypre_cmake_opts[@]}
cmake --build build --parallel 2
cmake --install build
- name: Checkout HDF5
if: ${{ matrix.hdf5 == 'HDF5' && matrix.distro == 'ubuntu:latest' }}
uses: actions/checkout@v3
with:
repository: HDFGroup/hdf5
ref: hdf5-1_14_3
path: hdf5-1_14_3
- name: Build HDF5
if: ${{ matrix.hdf5 == 'HDF5' && matrix.distro == 'ubuntu:latest' }}
working-directory: hdf5-1_14_3
run: |
mkdir hdf5
cd hdf5
export HDF5_ROOT=$HOME/hdf5
echo "HDF5_ROOT=$HDF5_ROOT" >> $GITHUB_ENV
cmake -C ../config/cmake/cacheinit.cmake -G "Unix Makefiles" \
-DCMAKE_C_COMPILER="mpicc" \
-DCMAKE_INSTALL_PREFIX=$HDF5_ROOT \
-DHDF5_BUILD_FORTRAN:BOOL=OFF \
-DHDF5_BUILD_JAVA:BOOL=OFF \
-DHDF5_BUILD_CPP_LIB:BOOL=OFF \
-DHDF5_BUILD_HL_LIB:BOOL=OFF \
-DHDF5_ENABLE_SZIP_SUPPORT:BOOL=OFF \
-DHDF5_ENABLE_Z_LIB_SUPPORT:BOOL=OFF \
-DBUILD_TESTING:BOOL=OFF \
-DHDF5_BUILD_EXAMPLES:BOOL=OFF \
-DHDF5_BUILD_HL_LIB:BOOL=OFF \
-DHDF5_BUILD_TOOLS:BOOL=OFF \
-DHDF5_ENABLE_PARALLEL:BOOL=ON \
-DHDF5_ENABLE_SUBFILING_VFD:BOOL=ON\
-DBUILD_SHARED_LIBS:BOOL=ON \
..
make -j 2 install
- name: Checkout ALL
if: ${{ matrix.liball == 'libALL' }}
run: |
Expand Down
119 changes: 114 additions & 5 deletions core/src/Cabana_HDF5ParticleOutput.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -31,6 +31,11 @@
#include <Kokkos_Core.hpp>

#include <hdf5.h>
#ifdef H5_HAVE_SUBFILING_VFD
#include "H5FDioc.h" /* Private header for the IOC VFD */
#include "H5FDsubfiling.h" /* Private header for the subfiling VFD */
#endif

#include <mpi.h>

#include <fstream>
Expand Down Expand Up @@ -142,6 +147,26 @@ struct HDF5Config

//! Cause all metadata for an object to be evicted from the cache
bool evict_on_close = false;

#ifdef H5_HAVE_SUBFILING_VFD

//! Use the subfiling file driver
bool subfiling = false;

// Optional subfiling file driver configuration parameters

//! Size (in bytes) of data stripes in subfiles
int64_t subfiling_stripe_size = H5FD_SUBFILING_DEFAULT_STRIPE_SIZE;

//! Target number of subfiles to use
int32_t subfiling_stripe_count = H5FD_SUBFILING_DEFAULT_STRIPE_COUNT;

//! The method to use for selecting MPI ranks to be I/O concentrators.
int subfiling_ioc_selection = SELECT_IOC_ONE_PER_NODE;

//! Number of I/O concentrator worker threads to use
int32_t subfiling_thread_pool_size = H5FD_IOC_DEFAULT_THREAD_POOL_SIZE;
#endif
};

//! \cond Impl
Expand Down Expand Up @@ -233,6 +258,7 @@ void writeFields(
{
hid_t plist_id;
hid_t dset_id;
hid_t dcpl_id;
hid_t filespace_id;
hid_t memspace_id;

Expand Down Expand Up @@ -261,8 +287,12 @@ void writeFields(
HDF5Traits<typename SliceType::value_type>::type( &dtype, &precision );

filespace_id = H5Screate_simple( 1, dimsf, NULL );

dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );

dset_id = H5Dcreate( file_id, slice.label().c_str(), type_id, filespace_id,
H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT );
H5P_DEFAULT, dcpl_id, H5P_DEFAULT );

H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, NULL, count,
NULL );
Expand All @@ -278,6 +308,7 @@ void writeFields(
host_view.data() );

H5Pclose( plist_id );
H5Pclose( dcpl_id );
H5Sclose( memspace_id );
H5Dclose( dset_id );
H5Sclose( filespace_id );
Expand All @@ -303,6 +334,7 @@ void writeFields(
{
hid_t plist_id;
hid_t dset_id;
hid_t dcpl_id;
hid_t filespace_id;
hid_t memspace_id;

Expand Down Expand Up @@ -339,8 +371,12 @@ void writeFields(
HDF5Traits<typename SliceType::value_type>::type( &dtype, &precision );

filespace_id = H5Screate_simple( 2, dimsf, NULL );

dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );

dset_id = H5Dcreate( file_id, slice.label().c_str(), type_id, filespace_id,
H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT );
H5P_DEFAULT, dcpl_id, H5P_DEFAULT );

H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, NULL, count,
NULL );
Expand All @@ -355,6 +391,7 @@ void writeFields(
host_view.data() );

H5Pclose( plist_id );
H5Pclose( dcpl_id );
H5Sclose( memspace_id );
H5Dclose( dset_id );
H5Sclose( filespace_id );
Expand All @@ -380,6 +417,7 @@ void writeFields(
{
hid_t plist_id;
hid_t dset_id;
hid_t dcpl_id;
hid_t filespace_id;
hid_t memspace_id;

Expand Down Expand Up @@ -420,8 +458,12 @@ void writeFields(
HDF5Traits<typename SliceType::value_type>::type( &dtype, &precision );

filespace_id = H5Screate_simple( 3, dimsf, NULL );

dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );

dset_id = H5Dcreate( file_id, slice.label().c_str(), type_id, filespace_id,
H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT );
H5P_DEFAULT, dcpl_id, H5P_DEFAULT );

H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, NULL, count,
NULL );
Expand All @@ -436,6 +478,7 @@ void writeFields(
host_view.data() );

H5Pclose( plist_id );
H5Pclose( dcpl_id );
H5Sclose( memspace_id );
H5Dclose( dset_id );
H5Sclose( filespace_id );
Expand Down Expand Up @@ -503,6 +546,7 @@ void writeTimeStep( HDF5Config h5_config, const std::string& prefix,

hid_t plist_id;
hid_t dset_id;
hid_t dcpl_id;
hid_t file_id;
hid_t filespace_id;
hid_t memspace_id;
Expand All @@ -525,7 +569,6 @@ void writeTimeStep( HDF5Config h5_config, const std::string& prefix,
filename_xdmf << prefix << "_" << time_step_index << ".xmf";

plist_id = H5Pcreate( H5P_FILE_ACCESS );
H5Pset_fapl_mpio( plist_id, comm, MPI_INFO_NULL );
H5Pset_libver_bounds( plist_id, H5F_LIBVER_LATEST, H5F_LIBVER_LATEST );

#if H5_VERSION_GE( 1, 10, 1 )
Expand All @@ -546,6 +589,68 @@ void writeTimeStep( HDF5Config h5_config, const std::string& prefix,
if ( h5_config.align )
H5Pset_alignment( plist_id, h5_config.threshold, h5_config.alignment );

#ifdef H5_HAVE_SUBFILING_VFD
if ( h5_config.subfiling )
{
H5FD_subfiling_config_t subfiling_config;
H5FD_ioc_config_t ioc_config;

H5FD_subfiling_config_t* subfiling_ptr = NULL;
H5FD_ioc_config_t* ioc_ptr = NULL;

// Get the default subfiling configuration parameters
hid_t fapl_id = H5I_INVALID_HID;

fapl_id = H5Pcreate( H5P_FILE_ACCESS );
H5Pget_fapl_subfiling( fapl_id, &subfiling_config );

if ( h5_config.subfiling_stripe_size !=
subfiling_config.shared_cfg.stripe_size )
{
subfiling_config.shared_cfg.stripe_size =
h5_config.subfiling_stripe_size;
if ( subfiling_ptr == NULL )
subfiling_ptr = &subfiling_config;
}
if ( h5_config.subfiling_stripe_count !=
subfiling_config.shared_cfg.stripe_count )
{
subfiling_config.shared_cfg.stripe_count =
h5_config.subfiling_stripe_count;
if ( subfiling_ptr == NULL )
subfiling_ptr = &subfiling_config;
}
if ( h5_config.subfiling_ioc_selection !=
(int)subfiling_config.shared_cfg.ioc_selection )
{
subfiling_config.shared_cfg.ioc_selection =
(H5FD_subfiling_ioc_select_t)h5_config.subfiling_ioc_selection;
if ( subfiling_ptr == NULL )
subfiling_ptr = &subfiling_config;
}
if ( h5_config.subfiling_thread_pool_size !=
H5FD_IOC_DEFAULT_THREAD_POOL_SIZE )
{
H5Pget_fapl_ioc( fapl_id, &ioc_config );
ioc_config.thread_pool_size = h5_config.subfiling_thread_pool_size;
if ( ioc_ptr == NULL )
ioc_ptr = &ioc_config;
}
H5Pclose( fapl_id );

H5Pset_mpi_params( plist_id, comm, MPI_INFO_NULL );

if ( ioc_ptr != NULL )
H5Pset_fapl_ioc( subfiling_config.ioc_fapl_id, ioc_ptr );

H5Pset_fapl_subfiling( plist_id, subfiling_ptr );
}
else
#endif
{
H5Pset_fapl_mpio( plist_id, comm, MPI_INFO_NULL );
}

file_id = H5Fcreate( filename_hdf5.str().c_str(), H5F_ACC_TRUNC,
H5P_DEFAULT, plist_id );
H5Pclose( plist_id );
Expand Down Expand Up @@ -617,8 +722,11 @@ void writeTimeStep( HDF5Config h5_config, const std::string& prefix,
hid_t type_id = HDF5Traits<typename CoordSliceType::value_type>::type(
&dtype, &precision );

dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );

dset_id = H5Dcreate( file_id, coords_slice.label().c_str(), type_id,
filespace_id, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT );
filespace_id, H5P_DEFAULT, dcpl_id, H5P_DEFAULT );

H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, NULL, count,
NULL );
Expand All @@ -628,6 +736,7 @@ void writeTimeStep( HDF5Config h5_config, const std::string& prefix,
H5Dclose( dset_id );

H5Pclose( plist_id );
H5Pclose( dcpl_id );
H5Sclose( filespace_id );
H5Sclose( memspace_id );

Expand Down
4 changes: 2 additions & 2 deletions example/core_tutorial/13_hdf5_output/hdf5_output_example.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -16,7 +16,7 @@
//---------------------------------------------------------------------------//
// HDF5 output example
//---------------------------------------------------------------------------//
void hdf5_output()
void hdf5Output()
{
/*
HDF5 is a parallel file format for large datasets. In this example, we
Expand Down Expand Up @@ -157,7 +157,7 @@ int main( int argc, char* argv[] )
MPI_Init( &argc, &argv );
Kokkos::initialize( argc, argv );

hdf5_output();
hdf5Output();

Kokkos::finalize();

Expand Down
22 changes: 22 additions & 0 deletions example/core_tutorial/13_hdf5_output_advanced/CMakeLists.txt
Original file line number Diff line number Diff line change
@@ -0,0 +1,22 @@
############################################################################
# Copyright (c) 2018-2022 by the Cabana authors #
# All rights reserved. #
# #
# This file is part of the Cabana library. Cabana is distributed under a #
# BSD 3-clause license. For the licensing terms see the LICENSE file in #
# the top-level directory. #
# #
# SPDX-License-Identifier: BSD-3-Clause #
############################################################################

add_executable(AdvancedHDF5Output advanced_hdf5_output_example.cpp)
target_link_libraries(AdvancedHDF5Output Cabana::Core)
set(HDF5_BIN_DIR ${HDF5_ROOT}/bin)
if(EXISTS ${HDF5_BIN_DIR}/h5fuse.sh)
file(COPY ${HDF5_BIN_DIR}/h5fuse.sh DESTINATION ${CMAKE_CURRENT_BINARY_DIR})
elseif(EXISTS ${HDF5_BIN_DIR}/h5fuse)
file(COPY ${HDF5_BIN_DIR}/h5fuse DESTINATION ${CMAKE_CURRENT_BINARY_DIR})
endif()
add_test(NAME Core_tutorial_13_subfiling COMMAND ${MPIEXEC} ${MPIEXEC_NUMPROC_FLAG}
${MPIEXEC_MAX_NUMPROCS} ${MPIEXEC_PREFLAGS} $<TARGET_FILE:AdvancedHDF5Output> ${MPIEXEC_POSTFLAGS})
set_tests_properties(Core_tutorial_13_subfiling PROPERTIES PROCESSORS ${MPIEXEC_MAX_NUMPROCS})
Loading

0 comments on commit 9e099a6

Please sign in to comment.