28#ifndef CABANA_HDF5PARTICLEOUTPUT_HPP
29#define CABANA_HDF5PARTICLEOUTPUT_HPP
31#include <Kokkos_Core.hpp>
32#include <Kokkos_Profiling_ScopedRegion.hpp>
35#ifdef H5_HAVE_SUBFILING_VFD
37#include "H5FDsubfiling.h"
54namespace HDF5ParticleOutput
61inline void writeXdmfHeader(
const char* xml_file_name,
const double time,
62 hsize_t dims0, hsize_t dims1,
const char* dtype,
63 uint precision,
const char* h5_file_name,
64 const char* coords_name )
66 std::ofstream xdmf_file( xml_file_name, std::ios::trunc );
68 xdmf_file.precision( std::numeric_limits<double>::max_digits10 );
69 xdmf_file <<
"<?xml version=\"1.0\" ?>\n";
70 xdmf_file <<
"<!DOCTYPE Xdmf SYSTEM \"Xdmf.dtd\" []>\n";
71 xdmf_file <<
"<Xdmf Version=\"2.0\">\n";
72 xdmf_file <<
" <Domain>\n";
73 xdmf_file <<
" <Grid Name=\"points\" GridType=\"Uniform\">\n";
74 xdmf_file <<
" <Time Value=\"" << time <<
"\"/>\n";
75 xdmf_file <<
" <Topology TopologyType=\"Polyvertex\"";
76 xdmf_file <<
" Dimensions=\"" << dims0 <<
"\"";
77 xdmf_file <<
" NodesPerElement=\"1\"> </Topology>\n";
78 xdmf_file <<
" <Geometry Type=\"" << ( dims1 == 3 ?
"XYZ" :
"XY" )
80 xdmf_file <<
" <DataItem Dimensions=\"" << dims0 <<
" " << dims1;
81 xdmf_file <<
"\" NumberType=\"" << dtype;
82 xdmf_file <<
"\" Precision=\"" << precision;
83 xdmf_file <<
"\" Format=\"HDF\"> " << h5_file_name <<
":/" << coords_name;
84 xdmf_file <<
" </DataItem>\n";
85 xdmf_file <<
" </Geometry>\n";
89inline void writeXdmfAttribute(
const char* xml_file_name,
90 const char* field_name, hsize_t dims0,
91 hsize_t dims1, hsize_t dims2,
const char* dtype,
92 uint precision,
const char* h5_file_name,
93 const char* dataitem )
95 std::string AttributeType =
"\"Scalar\"";
97 AttributeType =
"\"Tensor\"";
98 else if ( dims1 != 0 )
99 AttributeType =
"\"Vector\"";
101 std::ofstream xdmf_file( xml_file_name, std::ios::app );
102 xdmf_file <<
" <Attribute AttributeType =" << AttributeType
103 <<
" Center=\"Node\"";
104 xdmf_file <<
" Name=\"" << field_name <<
"\">\n";
105 xdmf_file <<
" <DataItem ItemType=\"Uniform\" Dimensions=\""
108 xdmf_file <<
" " << dims1;
110 xdmf_file <<
" " << dims2;
111 xdmf_file <<
"\" DataType=\"" << dtype <<
"\" Precision=\"" << precision
113 xdmf_file <<
" Format=\"HDF\"> " << h5_file_name <<
":/" << dataitem;
114 xdmf_file <<
" </DataItem>\n";
115 xdmf_file <<
" </Attribute>\n";
119inline void writeXdmfFooter(
const char* xml_file_name )
121 std::ofstream xdmf_file( xml_file_name, std::ios::app );
122 xdmf_file <<
" </Grid>\n";
123 xdmf_file <<
" </Domain>\n</Xdmf>\n";
158#ifdef H5_HAVE_SUBFILING_VFD
161 bool subfiling =
false;
166 int64_t subfiling_stripe_size = H5FD_SUBFILING_DEFAULT_STRIPE_SIZE;
169 int32_t subfiling_stripe_count = H5FD_SUBFILING_DEFAULT_STRIPE_COUNT;
172 int subfiling_ioc_selection = SELECT_IOC_ONE_PER_NODE;
175 int32_t subfiling_thread_pool_size = H5FD_IOC_DEFAULT_THREAD_POOL_SIZE;
185struct HDF5Traits<int>
187 static hid_t type( std::string* dtype, uint* precision )
190 *precision =
sizeof( int );
191 return H5T_NATIVE_INT;
196struct HDF5Traits<unsigned int>
198 static hid_t type( std::string* dtype, uint* precision )
201 *precision =
sizeof(
unsigned int );
202 return H5T_NATIVE_UINT;
207struct HDF5Traits<long>
209 static hid_t type( std::string* dtype, uint* precision )
212 *precision =
sizeof( long );
213 return H5T_NATIVE_LONG;
218struct HDF5Traits<unsigned long>
220 static hid_t type( std::string* dtype, uint* precision )
223 *precision =
sizeof(
unsigned long );
224 return H5T_NATIVE_ULONG;
229struct HDF5Traits<float>
231 static hid_t type( std::string* dtype, uint* precision )
234 *precision =
sizeof( float );
235 return H5T_NATIVE_FLOAT;
240struct HDF5Traits<double>
242 static hid_t type( std::string* dtype, uint* precision )
245 *precision =
sizeof( double );
246 return H5T_NATIVE_DOUBLE;
257template <
class SliceType>
259 HDF5Config h5_config, hid_t file_id, std::size_t n_local,
260 std::size_t n_global, hsize_t n_offset,
int comm_rank,
261 const char* filename_hdf5,
const char* filename_xdmf,
262 const SliceType&
slice,
263 typename std::enable_if<
264 2 == SliceType::kokkos_view::traits::dimension::rank,
int*>::type = 0 )
277 offset[0] = n_offset;
282 Kokkos::View<
typename SliceType::value_type*,
283 typename SliceType::memory_space>
284 view( Kokkos::ViewAllocateWithoutInitializing(
"field" ), n_local );
289 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), view );
294 HDF5Traits<typename SliceType::value_type>::type( &dtype, &precision );
296 filespace_id = H5Screate_simple( 1, dimsf,
nullptr );
298 dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
299 H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );
301 dset_id = H5Dcreate( file_id,
slice.label().c_str(), type_id, filespace_id,
302 H5P_DEFAULT, dcpl_id, H5P_DEFAULT );
304 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset,
nullptr, count,
307 memspace_id = H5Screate_simple( 1, count,
nullptr );
309 plist_id = H5Pcreate( H5P_DATASET_XFER );
311 if ( h5_config.collective )
312 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
314 H5Dwrite( dset_id, type_id, memspace_id, filespace_id, plist_id,
317 H5Pclose( plist_id );
319 H5Sclose( memspace_id );
321 H5Sclose( filespace_id );
323 if ( 0 == comm_rank )
326 Impl::writeXdmfAttribute(
327 filename_xdmf,
slice.label().c_str(), dimsf[0], zero, zero,
328 dtype.c_str(), precision, filename_hdf5,
slice.label().c_str() );
333template <
class SliceType>
335 HDF5Config h5_config, hid_t file_id, std::size_t n_local,
336 std::size_t n_global, hsize_t n_offset,
int comm_rank,
337 const char* filename_hdf5,
const char* filename_xdmf,
338 const SliceType&
slice,
339 typename std::enable_if<
340 3 == SliceType::kokkos_view::traits::dimension::rank,
int*>::type = 0 )
349 Kokkos::View<
typename SliceType::value_type**, Kokkos::LayoutRight,
350 typename SliceType::memory_space>
351 view( Kokkos::ViewAllocateWithoutInitializing(
"field" ), n_local,
357 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), view );
365 dimsf[1] = host_view.extent( 1 );
367 dimsm[1] = host_view.extent( 1 );
369 offset[0] = n_offset;
378 HDF5Traits<typename SliceType::value_type>::type( &dtype, &precision );
380 filespace_id = H5Screate_simple( 2, dimsf,
nullptr );
382 dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
383 H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );
385 dset_id = H5Dcreate( file_id,
slice.label().c_str(), type_id, filespace_id,
386 H5P_DEFAULT, dcpl_id, H5P_DEFAULT );
388 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset,
nullptr, count,
391 memspace_id = H5Screate_simple( 2, dimsm,
nullptr );
392 plist_id = H5Pcreate( H5P_DATASET_XFER );
394 if ( h5_config.collective )
395 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
397 H5Dwrite( dset_id, type_id, memspace_id, filespace_id, plist_id,
400 H5Pclose( plist_id );
402 H5Sclose( memspace_id );
404 H5Sclose( filespace_id );
406 if ( 0 == comm_rank )
409 Impl::writeXdmfAttribute(
410 filename_xdmf,
slice.label().c_str(), dimsf[0], dimsf[1], zero,
411 dtype.c_str(), precision, filename_hdf5,
slice.label().c_str() );
416template <
class SliceType>
418 HDF5Config h5_config, hid_t file_id, std::size_t n_local,
419 std::size_t n_global, hsize_t n_offset,
int comm_rank,
420 const char* filename_hdf5,
const char* filename_xdmf,
421 const SliceType&
slice,
422 typename std::enable_if<
423 4 == SliceType::kokkos_view::traits::dimension::rank,
int*>::type = 0 )
432 Kokkos::View<
typename SliceType::value_type***, Kokkos::LayoutRight,
433 typename SliceType::memory_space>
434 view( Kokkos::ViewAllocateWithoutInitializing(
"field" ), n_local,
440 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), view );
448 dimsf[1] = host_view.extent( 1 );
449 dimsf[2] = host_view.extent( 2 );
451 dimsm[1] = host_view.extent( 1 );
452 dimsm[2] = host_view.extent( 2 );
454 offset[0] = n_offset;
465 HDF5Traits<typename SliceType::value_type>::type( &dtype, &precision );
467 filespace_id = H5Screate_simple( 3, dimsf,
nullptr );
469 dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
470 H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );
472 dset_id = H5Dcreate( file_id,
slice.label().c_str(), type_id, filespace_id,
473 H5P_DEFAULT, dcpl_id, H5P_DEFAULT );
475 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset,
nullptr, count,
478 memspace_id = H5Screate_simple( 3, dimsm,
nullptr );
479 plist_id = H5Pcreate( H5P_DATASET_XFER );
481 if ( h5_config.collective )
482 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
484 H5Dwrite( dset_id, type_id, memspace_id, filespace_id, plist_id,
487 H5Pclose( plist_id );
489 H5Sclose( memspace_id );
491 H5Sclose( filespace_id );
493 if ( 0 == comm_rank )
495 Impl::writeXdmfAttribute(
496 filename_xdmf,
slice.label().c_str(), dimsf[0], dimsf[1], dimsf[2],
497 dtype.c_str(), precision, filename_hdf5,
slice.label().c_str() );
504inline void writeFields( HDF5Config, hid_t, std::size_t, std::size_t, hsize_t,
505 int,
const char*,
const char* )
510template <
class SliceType>
511void writeFields( HDF5Config h5_config, hid_t file_id, std::size_t n_local,
512 std::size_t n_global, hsize_t n_offset,
int comm_rank,
513 const char* filename_hdf5,
const char* filename_xdmf,
514 const SliceType&
slice )
516 Impl::writeFields( h5_config, file_id, n_local, n_global, n_offset,
517 comm_rank, filename_hdf5, filename_xdmf,
slice );
521template <
class SliceType,
class... FieldSliceTypes>
522void writeFields( HDF5Config h5_config, hid_t file_id, std::size_t n_local,
523 std::size_t n_global, hsize_t n_offset,
int comm_rank,
524 const char* filename_hdf5,
const char* filename_xdmf,
525 const SliceType&
slice, FieldSliceTypes&&... fields )
527 Impl::writeFields( h5_config, file_id, n_local, n_global, n_offset,
528 comm_rank, filename_hdf5, filename_xdmf,
slice );
529 writeFields( h5_config, file_id, n_local, n_global, n_offset, comm_rank,
530 filename_hdf5, filename_xdmf, fields... );
545template <
class CoordSliceType,
class... FieldSliceTypes>
547 MPI_Comm comm,
const int time_step_index,
const double time,
548 const std::size_t n_local,
549 const CoordSliceType& coords_slice,
550 FieldSliceTypes&&... fields )
552 Kokkos::Profiling::ScopedRegion region(
"Cabana::HDF5ParticleOutput" );
567 MPI_Comm_rank( comm, &comm_rank );
569 MPI_Comm_size( comm, &comm_size );
572 std::stringstream filename_hdf5;
573 filename_hdf5 << prefix <<
"_" << time_step_index <<
".h5";
575 std::stringstream filename_xdmf;
576 filename_xdmf << prefix <<
"_" << time_step_index <<
".xmf";
578 plist_id = H5Pcreate( H5P_FILE_ACCESS );
579 H5Pset_libver_bounds( plist_id, H5F_LIBVER_LATEST, H5F_LIBVER_LATEST );
581#if H5_VERSION_GE( 1, 10, 1 )
582 if ( h5_config.evict_on_close )
584 H5Pset_evict_on_close( plist_id,
true );
588#if H5_VERSION_GE( 1, 10, 0 )
589 if ( h5_config.collective )
591 H5Pset_all_coll_metadata_ops( plist_id,
true );
592 H5Pset_coll_metadata_write( plist_id,
true );
596 if ( h5_config.align )
597 H5Pset_alignment( plist_id, h5_config.threshold, h5_config.alignment );
599#ifdef H5_HAVE_SUBFILING_VFD
600 if ( h5_config.subfiling )
602 H5FD_subfiling_config_t subfiling_config;
603 H5FD_ioc_config_t ioc_config;
605 H5FD_subfiling_config_t* subfiling_ptr =
nullptr;
606 H5FD_ioc_config_t* ioc_ptr =
nullptr;
609 hid_t fapl_id = H5I_INVALID_HID;
611 fapl_id = H5Pcreate( H5P_FILE_ACCESS );
612 H5Pget_fapl_subfiling( fapl_id, &subfiling_config );
614 if ( h5_config.subfiling_stripe_size !=
615 subfiling_config.shared_cfg.stripe_size )
617 subfiling_config.shared_cfg.stripe_size =
618 h5_config.subfiling_stripe_size;
619 if ( subfiling_ptr ==
nullptr )
620 subfiling_ptr = &subfiling_config;
622 if ( h5_config.subfiling_stripe_count !=
623 subfiling_config.shared_cfg.stripe_count )
625 subfiling_config.shared_cfg.stripe_count =
626 h5_config.subfiling_stripe_count;
627 if ( subfiling_ptr ==
nullptr )
628 subfiling_ptr = &subfiling_config;
630 if ( h5_config.subfiling_ioc_selection !=
631 (
int)subfiling_config.shared_cfg.ioc_selection )
633 subfiling_config.shared_cfg.ioc_selection =
634 (H5FD_subfiling_ioc_select_t)h5_config.subfiling_ioc_selection;
635 if ( subfiling_ptr ==
nullptr )
636 subfiling_ptr = &subfiling_config;
638 if ( h5_config.subfiling_thread_pool_size !=
639 H5FD_IOC_DEFAULT_THREAD_POOL_SIZE )
641 H5Pget_fapl_ioc( fapl_id, &ioc_config );
642 ioc_config.thread_pool_size = h5_config.subfiling_thread_pool_size;
643 if ( ioc_ptr ==
nullptr )
644 ioc_ptr = &ioc_config;
648 H5Pset_mpi_params( plist_id, comm, MPI_INFO_NULL );
650 if ( ioc_ptr !=
nullptr )
651 H5Pset_fapl_ioc( subfiling_config.ioc_fapl_id, ioc_ptr );
653 H5Pset_fapl_subfiling( plist_id, subfiling_ptr );
658 H5Pset_fapl_mpio( plist_id, comm, MPI_INFO_NULL );
661 file_id = H5Fcreate( filename_hdf5.str().c_str(), H5F_ACC_TRUNC,
662 H5P_DEFAULT, plist_id );
663 H5Pclose( plist_id );
666 hid_t fspace = H5Screate( H5S_SCALAR );
667 hid_t attr_id = H5Acreate( file_id,
"Time", H5T_NATIVE_DOUBLE, fspace,
668 H5P_DEFAULT, H5P_DEFAULT );
669 H5Awrite( attr_id, H5T_NATIVE_DOUBLE, &time );
674 Kokkos::View<
typename CoordSliceType::value_type**, Kokkos::LayoutRight,
675 typename CoordSliceType::memory_space>
676 coords_view( Kokkos::ViewAllocateWithoutInitializing(
"coords" ),
677 coords_slice.size(), coords_slice.extent( 2 ) );
678 Kokkos::parallel_for(
679 "Cabana::HDF5ParticleOutput::copyCoords",
680 Kokkos::RangePolicy<typename CoordSliceType::execution_space>(
681 0, coords_slice.size() ),
682 KOKKOS_LAMBDA(
const int i ) {
683 for ( std::size_t d0 = 0; d0 < coords_slice.extent( 2 ); ++d0 )
684 coords_view( i, d0 ) = coords_slice( i, d0 );
689 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), coords_view );
691 std::vector<int> all_offsets( comm_size );
692 all_offsets[comm_rank] = n_local;
694 MPI_Allreduce( MPI_IN_PLACE, all_offsets.data(), comm_size, MPI_INT,
701 for (
int i = 0; i < comm_size; i++ )
705 offset[0] +=
static_cast<hsize_t
>( all_offsets[i] );
707 n_global += (size_t)all_offsets[i];
709 std::vector<int>().swap( all_offsets );
712 dimsf[1] = coords_slice.extent( 2 );
714 filespace_id = H5Screate_simple( 2, dimsf,
nullptr );
717 count[1] = coords_slice.extent( 2 );
719 memspace_id = H5Screate_simple( 2, count,
nullptr );
721 plist_id = H5Pcreate( H5P_DATASET_XFER );
724 if ( h5_config.collective )
725 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
729 hid_t type_id = HDF5Traits<typename CoordSliceType::value_type>::type(
730 &dtype, &precision );
732 dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
733 H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );
735 dset_id = H5Dcreate( file_id, coords_slice.label().c_str(), type_id,
736 filespace_id, H5P_DEFAULT, dcpl_id, H5P_DEFAULT );
738 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset,
nullptr, count,
741 H5Dwrite( dset_id, type_id, memspace_id, filespace_id, plist_id,
742 host_coords.data() );
745 H5Pclose( plist_id );
747 H5Sclose( filespace_id );
748 H5Sclose( memspace_id );
750 if ( 0 == comm_rank )
752 Impl::writeXdmfHeader( filename_xdmf.str().c_str(), time, dimsf[0],
753 dimsf[1], dtype.c_str(), precision,
754 filename_hdf5.str().c_str(),
755 coords_slice.label().c_str() );
759 hsize_t n_offset = offset[0];
760 writeFields( h5_config, file_id, n_local, n_global, n_offset, comm_rank,
761 filename_hdf5.str().c_str(), filename_xdmf.str().c_str(),
766 if ( 0 == comm_rank )
767 Impl::writeXdmfFooter( filename_xdmf.str().c_str() );
774template <
class SliceType>
776 hid_t dset_id, hid_t dtype_id, hid_t memspace_id, hid_t filespace_id,
777 hid_t plist_id, std::size_t n_local,
const SliceType&
slice,
778 typename std::enable_if<
779 2 == SliceType::kokkos_view::traits::dimension::rank,
int*>::type = 0 )
782 Kokkos::View<typename SliceType::value_type*, Kokkos::HostSpace> host_view(
783 Kokkos::ViewAllocateWithoutInitializing(
"field" ), n_local );
784 H5Dread( dset_id, dtype_id, memspace_id, filespace_id, plist_id,
788 auto view = Kokkos::create_mirror_view_and_copy(
789 typename SliceType::memory_space(), host_view );
794template <
class SliceType>
796 hid_t dset_id, hid_t dtype_id, hid_t memspace_id, hid_t filespace_id,
797 hid_t plist_id, std::size_t n_local,
const SliceType&
slice,
798 typename std::enable_if<
799 3 == SliceType::kokkos_view::traits::dimension::rank,
int*>::type = 0 )
802 Kokkos::View<
typename SliceType::value_type**, Kokkos::LayoutRight,
804 host_view( Kokkos::ViewAllocateWithoutInitializing(
"field" ), n_local,
806 H5Dread( dset_id, dtype_id, memspace_id, filespace_id, plist_id,
810 auto view = Kokkos::create_mirror_view_and_copy(
811 typename SliceType::memory_space(), host_view );
816template <
class SliceType>
818 hid_t dset_id, hid_t dtype_id, hid_t memspace_id, hid_t filespace_id,
819 hid_t plist_id, std::size_t n_local,
const SliceType&
slice,
820 typename std::enable_if<
821 4 == SliceType::kokkos_view::traits::dimension::rank,
int*>::type = 0 )
824 Kokkos::View<
typename SliceType::value_type***, Kokkos::LayoutRight,
826 host_view( Kokkos::ViewAllocateWithoutInitializing(
"field" ), n_local,
828 H5Dread( dset_id, dtype_id, memspace_id, filespace_id, plist_id,
832 auto view = Kokkos::create_mirror_view_and_copy(
833 typename SliceType::memory_space(), host_view );
849template <
class FieldSliceType>
851 MPI_Comm comm,
const int time_step_index,
852 const std::size_t n_local,
const std::string& dataset_name,
853 double& time, FieldSliceType& field )
855 Kokkos::Profiling::ScopedRegion region(
"Cabana::HDF5ParticleInput" );
866 hsize_t offset[3] = { 0, 0, 0 };
867 hsize_t dimsf[3] = { 0, 0, 0 };
868 hsize_t count[3] = { 0, 0, 0 };
871 MPI_Comm_rank( comm, &comm_rank );
873 MPI_Comm_size( comm, &comm_size );
876 std::stringstream filename_hdf5;
877 filename_hdf5 << prefix <<
"_" << time_step_index <<
".h5";
879 plist_id = H5Pcreate( H5P_FILE_ACCESS );
880 H5Pset_fapl_mpio( plist_id, comm, MPI_INFO_NULL );
882#if H5_VERSION_GE( 1, 10, 0 )
883 if ( h5_config.collective )
885 H5Pset_all_coll_metadata_ops( plist_id,
true );
890 file_id = H5Fopen( filename_hdf5.str().c_str(), H5F_ACC_RDONLY, plist_id );
891 H5Pclose( plist_id );
894 hid_t attr_id = H5Aopen( file_id,
"Time", H5P_DEFAULT );
895 H5Aread( attr_id, H5T_NATIVE_DOUBLE, &time );
899 dset_id = H5Dopen( file_id, dataset_name.c_str(), H5P_DEFAULT );
902 dtype_id = H5Dget_type( dset_id );
905 filespace_id = H5Dget_space( dset_id );
908 ndims = H5Sget_simple_extent_ndims( filespace_id );
911 H5Sget_simple_extent_dims( filespace_id, dimsf,
nullptr );
913 std::vector<int> all_offsets( comm_size );
914 all_offsets[comm_rank] = n_local;
916 MPI_Allreduce( MPI_IN_PLACE, all_offsets.data(), comm_size, MPI_INT,
919 for (
int i = 0; i < comm_size; i++ )
923 offset[0] +=
static_cast<hsize_t
>( all_offsets[i] );
926 std::vector<int>().swap( all_offsets );
932 memspace_id = H5Screate_simple( ndims, count,
nullptr );
934 plist_id = H5Pcreate( H5P_DATASET_XFER );
937 if ( h5_config.collective )
938 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
940 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset,
nullptr, count,
943 readField( dset_id, dtype_id, memspace_id, filespace_id, plist_id, n_local,
946 H5Pclose( plist_id );
947 H5Sclose( memspace_id );
948 H5Sclose( filespace_id );
void writeFields(HDF5Config, hid_t, std::size_t, std::size_t, hsize_t, int, const char *, const char *)
Write particle data to HDF5 output. Empty overload if only writing coords.
Definition Cabana_HDF5ParticleOutput.hpp:504
void readField(hid_t dset_id, hid_t dtype_id, hid_t memspace_id, hid_t filespace_id, hid_t plist_id, std::size_t n_local, const SliceType &slice, typename std::enable_if< 2==SliceType::kokkos_view::traits::dimension::rank, int * >::type=0)
Read particle data from HDF5 output. Rank-0.
Definition Cabana_HDF5ParticleOutput.hpp:775
void readTimeStep(HDF5Config h5_config, const std::string &prefix, MPI_Comm comm, const int time_step_index, const std::size_t n_local, const std::string &dataset_name, double &time, FieldSliceType &field)
Read particle output from an HDF5 file.
Definition Cabana_HDF5ParticleOutput.hpp:850
void writeTimeStep(HDF5Config h5_config, const std::string &prefix, MPI_Comm comm, const int time_step_index, const double time, const std::size_t n_local, const CoordSliceType &coords_slice, FieldSliceTypes &&... fields)
Write particle output in HDF5 format.
Definition Cabana_HDF5ParticleOutput.hpp:546
Core: particle data structures and algorithms.
Definition Cabana_AoSoA.hpp:36
void copyViewToSlice(ExecutionSpace exec_space, SliceType &slice, const ViewType &view, const std::size_t begin, const std::size_t end, typename std::enable_if< 2==SliceType::kokkos_view::traits::dimension::rank, int * >::type=0)
Copy from View to slice. Rank-0.
Definition Cabana_Slice.hpp:939
void copySliceToView(ExecutionSpace exec_space, ViewType &view, const SliceType &slice, const std::size_t begin, const std::size_t end, typename std::enable_if< 2==SliceType::kokkos_view::traits::dimension::rank, int * >::type=0)
Copy from slice to View. Rank-0.
Definition Cabana_Slice.hpp:877
AoSoA_t::template member_slice_type< M > slice(const AoSoA_t &aosoa, const std::string &slice_label="")
Create a slice from an AoSoA.
Definition Cabana_AoSoA.hpp:77
HDF5 tuning settings.
Definition Cabana_HDF5ParticleOutput.hpp:141
bool meta_collective
Sets metadata I/O mode operations to collective or independent (default)
Definition Cabana_HDF5ParticleOutput.hpp:153
bool align
Set alignment on or off.
Definition Cabana_HDF5ParticleOutput.hpp:146
bool evict_on_close
Cause all metadata for an object to be evicted from the cache.
Definition Cabana_HDF5ParticleOutput.hpp:156
unsigned long threshold
Threshold for aligning file objects.
Definition Cabana_HDF5ParticleOutput.hpp:148
unsigned long alignment
Alignment value.
Definition Cabana_HDF5ParticleOutput.hpp:150
bool collective
I/O transfer mode to collective or independent (default)
Definition Cabana_HDF5ParticleOutput.hpp:143