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"
53namespace HDF5ParticleOutput
60inline void writeXdmfHeader(
const char* xml_file_name, hsize_t dims0,
61 hsize_t dims1,
const char* dtype, uint precision,
62 const char* h5_file_name,
const char* coords_name )
64 std::ofstream xdmf_file( xml_file_name, std::ios::trunc );
65 xdmf_file <<
"<?xml version=\"1.0\" ?>\n";
66 xdmf_file <<
"<!DOCTYPE Xdmf SYSTEM \"Xdmf.dtd\" []>\n";
67 xdmf_file <<
"<Xdmf Version=\"2.0\">\n";
68 xdmf_file <<
" <Domain>\n";
69 xdmf_file <<
" <Grid Name=\"points\" GridType=\"Uniform\">\n";
70 xdmf_file <<
" <Topology TopologyType=\"Polyvertex\"";
71 xdmf_file <<
" Dimensions=\"" << dims0 <<
"\"";
72 xdmf_file <<
" NodesPerElement=\"1\"> </Topology>\n";
73 xdmf_file <<
" <Geometry Type=\"XYZ\">\n";
74 xdmf_file <<
" <DataItem Dimensions=\"" << dims0 <<
" " << dims1;
75 xdmf_file <<
"\" NumberType=\"" << dtype;
76 xdmf_file <<
"\" Precision=\"" << precision;
77 xdmf_file <<
"\" Format=\"HDF\"> " << h5_file_name <<
":/" << coords_name;
78 xdmf_file <<
" </DataItem>\n";
79 xdmf_file <<
" </Geometry>\n";
83inline void writeXdmfAttribute(
const char* xml_file_name,
84 const char* field_name, hsize_t dims0,
85 hsize_t dims1, hsize_t dims2,
const char* dtype,
86 uint precision,
const char* h5_file_name,
87 const char* dataitem )
89 std::string AttributeType =
"\"Scalar\"";
91 AttributeType =
"\"Tensor\"";
92 else if ( dims1 != 0 )
93 AttributeType =
"\"Vector\"";
95 std::ofstream xdmf_file( xml_file_name, std::ios::app );
96 xdmf_file <<
" <Attribute AttributeType =" << AttributeType
97 <<
" Center=\"Node\"";
98 xdmf_file <<
" Name=\"" << field_name <<
"\">\n";
99 xdmf_file <<
" <DataItem ItemType=\"Uniform\" Dimensions=\""
102 xdmf_file <<
" " << dims1;
104 xdmf_file <<
" " << dims2;
105 xdmf_file <<
"\" DataType=\"" << dtype <<
"\" Precision=\"" << precision
107 xdmf_file <<
" Format=\"HDF\"> " << h5_file_name <<
":/" << dataitem;
108 xdmf_file <<
" </DataItem>\n";
109 xdmf_file <<
" </Attribute>\n";
113inline void writeXdmfFooter(
const char* xml_file_name )
115 std::ofstream xdmf_file( xml_file_name, std::ios::app );
116 xdmf_file <<
" </Grid>\n";
117 xdmf_file <<
" </Domain>\n</Xdmf>\n";
152#ifdef H5_HAVE_SUBFILING_VFD
155 bool subfiling =
false;
160 int64_t subfiling_stripe_size = H5FD_SUBFILING_DEFAULT_STRIPE_SIZE;
163 int32_t subfiling_stripe_count = H5FD_SUBFILING_DEFAULT_STRIPE_COUNT;
166 int subfiling_ioc_selection = SELECT_IOC_ONE_PER_NODE;
169 int32_t subfiling_thread_pool_size = H5FD_IOC_DEFAULT_THREAD_POOL_SIZE;
179struct HDF5Traits<int>
181 static hid_t type( std::string* dtype, uint* precision )
184 *precision =
sizeof( int );
185 return H5T_NATIVE_INT;
190struct HDF5Traits<unsigned int>
192 static hid_t type( std::string* dtype, uint* precision )
195 *precision =
sizeof(
unsigned int );
196 return H5T_NATIVE_UINT;
201struct HDF5Traits<long>
203 static hid_t type( std::string* dtype, uint* precision )
206 *precision =
sizeof( long );
207 return H5T_NATIVE_LONG;
212struct HDF5Traits<unsigned long>
214 static hid_t type( std::string* dtype, uint* precision )
217 *precision =
sizeof(
unsigned long );
218 return H5T_NATIVE_ULONG;
223struct HDF5Traits<float>
225 static hid_t type( std::string* dtype, uint* precision )
228 *precision =
sizeof( float );
229 return H5T_NATIVE_FLOAT;
234struct HDF5Traits<double>
236 static hid_t type( std::string* dtype, uint* precision )
239 *precision =
sizeof( double );
240 return H5T_NATIVE_DOUBLE;
251template <
class SliceType>
253 HDF5Config h5_config, hid_t file_id, std::size_t n_local,
254 std::size_t n_global, hsize_t n_offset,
int comm_rank,
255 const char* filename_hdf5,
const char* filename_xdmf,
256 const SliceType&
slice,
257 typename std::enable_if<
258 2 == SliceType::kokkos_view::traits::dimension::rank,
int*>::type = 0 )
271 offset[0] = n_offset;
276 Kokkos::View<
typename SliceType::value_type*,
277 typename SliceType::memory_space>
278 view( Kokkos::ViewAllocateWithoutInitializing(
"field" ), n_local );
283 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), view );
288 HDF5Traits<typename SliceType::value_type>::type( &dtype, &precision );
290 filespace_id = H5Screate_simple( 1, dimsf, NULL );
292 dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
293 H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );
295 dset_id = H5Dcreate( file_id,
slice.label().c_str(), type_id, filespace_id,
296 H5P_DEFAULT, dcpl_id, H5P_DEFAULT );
298 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, NULL, count,
301 memspace_id = H5Screate_simple( 1, count, NULL );
303 plist_id = H5Pcreate( H5P_DATASET_XFER );
305 if ( h5_config.collective )
306 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
308 H5Dwrite( dset_id, type_id, memspace_id, filespace_id, plist_id,
311 H5Pclose( plist_id );
313 H5Sclose( memspace_id );
315 H5Sclose( filespace_id );
317 if ( 0 == comm_rank )
320 Impl::writeXdmfAttribute(
321 filename_xdmf,
slice.label().c_str(), dimsf[0], zero, zero,
322 dtype.c_str(), precision, filename_hdf5,
slice.label().c_str() );
327template <
class SliceType>
329 HDF5Config h5_config, hid_t file_id, std::size_t n_local,
330 std::size_t n_global, hsize_t n_offset,
int comm_rank,
331 const char* filename_hdf5,
const char* filename_xdmf,
332 const SliceType&
slice,
333 typename std::enable_if<
334 3 == SliceType::kokkos_view::traits::dimension::rank,
int*>::type = 0 )
343 Kokkos::View<
typename SliceType::value_type**, Kokkos::LayoutRight,
344 typename SliceType::memory_space>
345 view( Kokkos::ViewAllocateWithoutInitializing(
"field" ), n_local,
351 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), view );
359 dimsf[1] = host_view.extent( 1 );
361 dimsm[1] = host_view.extent( 1 );
363 offset[0] = n_offset;
372 HDF5Traits<typename SliceType::value_type>::type( &dtype, &precision );
374 filespace_id = H5Screate_simple( 2, dimsf, NULL );
376 dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
377 H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );
379 dset_id = H5Dcreate( file_id,
slice.label().c_str(), type_id, filespace_id,
380 H5P_DEFAULT, dcpl_id, H5P_DEFAULT );
382 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, NULL, count,
385 memspace_id = H5Screate_simple( 2, dimsm, NULL );
386 plist_id = H5Pcreate( H5P_DATASET_XFER );
388 if ( h5_config.collective )
389 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
391 H5Dwrite( dset_id, type_id, memspace_id, filespace_id, plist_id,
394 H5Pclose( plist_id );
396 H5Sclose( memspace_id );
398 H5Sclose( filespace_id );
400 if ( 0 == comm_rank )
403 Impl::writeXdmfAttribute(
404 filename_xdmf,
slice.label().c_str(), dimsf[0], dimsf[1], zero,
405 dtype.c_str(), precision, filename_hdf5,
slice.label().c_str() );
410template <
class SliceType>
412 HDF5Config h5_config, hid_t file_id, std::size_t n_local,
413 std::size_t n_global, hsize_t n_offset,
int comm_rank,
414 const char* filename_hdf5,
const char* filename_xdmf,
415 const SliceType&
slice,
416 typename std::enable_if<
417 4 == SliceType::kokkos_view::traits::dimension::rank,
int*>::type = 0 )
426 Kokkos::View<
typename SliceType::value_type***, Kokkos::LayoutRight,
427 typename SliceType::memory_space>
428 view( Kokkos::ViewAllocateWithoutInitializing(
"field" ), n_local,
434 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), view );
442 dimsf[1] = host_view.extent( 1 );
443 dimsf[2] = host_view.extent( 2 );
445 dimsm[1] = host_view.extent( 1 );
446 dimsm[2] = host_view.extent( 2 );
448 offset[0] = n_offset;
459 HDF5Traits<typename SliceType::value_type>::type( &dtype, &precision );
461 filespace_id = H5Screate_simple( 3, dimsf, NULL );
463 dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
464 H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );
466 dset_id = H5Dcreate( file_id,
slice.label().c_str(), type_id, filespace_id,
467 H5P_DEFAULT, dcpl_id, H5P_DEFAULT );
469 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, NULL, count,
472 memspace_id = H5Screate_simple( 3, dimsm, NULL );
473 plist_id = H5Pcreate( H5P_DATASET_XFER );
475 if ( h5_config.collective )
476 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
478 H5Dwrite( dset_id, type_id, memspace_id, filespace_id, plist_id,
481 H5Pclose( plist_id );
483 H5Sclose( memspace_id );
485 H5Sclose( filespace_id );
487 if ( 0 == comm_rank )
489 Impl::writeXdmfAttribute(
490 filename_xdmf,
slice.label().c_str(), dimsf[0], dimsf[1], dimsf[2],
491 dtype.c_str(), precision, filename_hdf5,
slice.label().c_str() );
498inline void writeFields( HDF5Config, hid_t, std::size_t, std::size_t, hsize_t,
499 int,
const char*,
const char* )
504template <
class SliceType>
505void writeFields( HDF5Config h5_config, hid_t file_id, std::size_t n_local,
506 std::size_t n_global, hsize_t n_offset,
int comm_rank,
507 const char* filename_hdf5,
const char* filename_xdmf,
508 const SliceType&
slice )
510 Impl::writeFields( h5_config, file_id, n_local, n_global, n_offset,
511 comm_rank, filename_hdf5, filename_xdmf,
slice );
515template <
class SliceType,
class... FieldSliceTypes>
516void writeFields( HDF5Config h5_config, hid_t file_id, std::size_t n_local,
517 std::size_t n_global, hsize_t n_offset,
int comm_rank,
518 const char* filename_hdf5,
const char* filename_xdmf,
519 const SliceType&
slice, FieldSliceTypes&&... fields )
521 Impl::writeFields( h5_config, file_id, n_local, n_global, n_offset,
522 comm_rank, filename_hdf5, filename_xdmf,
slice );
523 writeFields( h5_config, file_id, n_local, n_global, n_offset, comm_rank,
524 filename_hdf5, filename_xdmf, fields... );
539template <
class CoordSliceType,
class... FieldSliceTypes>
541 MPI_Comm comm,
const int time_step_index,
const double time,
542 const std::size_t n_local,
543 const CoordSliceType& coords_slice,
544 FieldSliceTypes&&... fields )
546 Kokkos::Profiling::ScopedRegion region(
"Cabana::HDF5ParticleOutput" );
561 MPI_Comm_rank( comm, &comm_rank );
563 MPI_Comm_size( comm, &comm_size );
566 std::stringstream filename_hdf5;
567 filename_hdf5 << prefix <<
"_" << time_step_index <<
".h5";
569 std::stringstream filename_xdmf;
570 filename_xdmf << prefix <<
"_" << time_step_index <<
".xmf";
572 plist_id = H5Pcreate( H5P_FILE_ACCESS );
573 H5Pset_libver_bounds( plist_id, H5F_LIBVER_LATEST, H5F_LIBVER_LATEST );
575#if H5_VERSION_GE( 1, 10, 1 )
576 if ( h5_config.evict_on_close )
578 H5Pset_evict_on_close( plist_id, (hbool_t)1 );
582#if H5_VERSION_GE( 1, 10, 0 )
583 if ( h5_config.collective )
585 H5Pset_all_coll_metadata_ops( plist_id, 1 );
586 H5Pset_coll_metadata_write( plist_id, 1 );
590 if ( h5_config.align )
591 H5Pset_alignment( plist_id, h5_config.threshold, h5_config.alignment );
593#ifdef H5_HAVE_SUBFILING_VFD
594 if ( h5_config.subfiling )
596 H5FD_subfiling_config_t subfiling_config;
597 H5FD_ioc_config_t ioc_config;
599 H5FD_subfiling_config_t* subfiling_ptr = NULL;
600 H5FD_ioc_config_t* ioc_ptr = NULL;
603 hid_t fapl_id = H5I_INVALID_HID;
605 fapl_id = H5Pcreate( H5P_FILE_ACCESS );
606 H5Pget_fapl_subfiling( fapl_id, &subfiling_config );
608 if ( h5_config.subfiling_stripe_size !=
609 subfiling_config.shared_cfg.stripe_size )
611 subfiling_config.shared_cfg.stripe_size =
612 h5_config.subfiling_stripe_size;
613 if ( subfiling_ptr == NULL )
614 subfiling_ptr = &subfiling_config;
616 if ( h5_config.subfiling_stripe_count !=
617 subfiling_config.shared_cfg.stripe_count )
619 subfiling_config.shared_cfg.stripe_count =
620 h5_config.subfiling_stripe_count;
621 if ( subfiling_ptr == NULL )
622 subfiling_ptr = &subfiling_config;
624 if ( h5_config.subfiling_ioc_selection !=
625 (
int)subfiling_config.shared_cfg.ioc_selection )
627 subfiling_config.shared_cfg.ioc_selection =
628 (H5FD_subfiling_ioc_select_t)h5_config.subfiling_ioc_selection;
629 if ( subfiling_ptr == NULL )
630 subfiling_ptr = &subfiling_config;
632 if ( h5_config.subfiling_thread_pool_size !=
633 H5FD_IOC_DEFAULT_THREAD_POOL_SIZE )
635 H5Pget_fapl_ioc( fapl_id, &ioc_config );
636 ioc_config.thread_pool_size = h5_config.subfiling_thread_pool_size;
637 if ( ioc_ptr == NULL )
638 ioc_ptr = &ioc_config;
642 H5Pset_mpi_params( plist_id, comm, MPI_INFO_NULL );
644 if ( ioc_ptr != NULL )
645 H5Pset_fapl_ioc( subfiling_config.ioc_fapl_id, ioc_ptr );
647 H5Pset_fapl_subfiling( plist_id, subfiling_ptr );
652 H5Pset_fapl_mpio( plist_id, comm, MPI_INFO_NULL );
655 file_id = H5Fcreate( filename_hdf5.str().c_str(), H5F_ACC_TRUNC,
656 H5P_DEFAULT, plist_id );
657 H5Pclose( plist_id );
660 hid_t fspace = H5Screate( H5S_SCALAR );
661 hid_t attr_id = H5Acreate( file_id,
"Time", H5T_NATIVE_DOUBLE, fspace,
662 H5P_DEFAULT, H5P_DEFAULT );
663 H5Awrite( attr_id, H5T_NATIVE_DOUBLE, &time );
668 Kokkos::View<
typename CoordSliceType::value_type**, Kokkos::LayoutRight,
669 typename CoordSliceType::memory_space>
670 coords_view( Kokkos::ViewAllocateWithoutInitializing(
"coords" ),
671 coords_slice.size(), coords_slice.extent( 2 ) );
672 Kokkos::parallel_for(
673 "Cabana::HDF5ParticleOutput::writeCoords",
674 Kokkos::RangePolicy<typename CoordSliceType::execution_space>(
675 0, coords_slice.size() ),
676 KOKKOS_LAMBDA(
const int i ) {
677 for ( std::size_t d0 = 0; d0 < coords_slice.extent( 2 ); ++d0 )
678 coords_view( i, d0 ) = coords_slice( i, d0 );
683 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), coords_view );
685 std::vector<int> all_offsets( comm_size );
686 all_offsets[comm_rank] = n_local;
688 MPI_Allreduce( MPI_IN_PLACE, all_offsets.data(), comm_size, MPI_INT,
695 for (
int i = 0; i < comm_size; i++ )
699 offset[0] +=
static_cast<hsize_t
>( all_offsets[i] );
701 n_global += (size_t)all_offsets[i];
703 std::vector<int>().swap( all_offsets );
708 filespace_id = H5Screate_simple( 2, dimsf, NULL );
713 memspace_id = H5Screate_simple( 2, count, NULL );
715 plist_id = H5Pcreate( H5P_DATASET_XFER );
718 if ( h5_config.collective )
719 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
723 hid_t type_id = HDF5Traits<typename CoordSliceType::value_type>::type(
724 &dtype, &precision );
726 dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
727 H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );
729 dset_id = H5Dcreate( file_id, coords_slice.label().c_str(), type_id,
730 filespace_id, H5P_DEFAULT, dcpl_id, H5P_DEFAULT );
732 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, NULL, count,
735 H5Dwrite( dset_id, type_id, memspace_id, filespace_id, plist_id,
736 host_coords.data() );
739 H5Pclose( plist_id );
741 H5Sclose( filespace_id );
742 H5Sclose( memspace_id );
744 if ( 0 == comm_rank )
746 Impl::writeXdmfHeader( filename_xdmf.str().c_str(), dimsf[0], dimsf[1],
747 dtype.c_str(), precision,
748 filename_hdf5.str().c_str(),
749 coords_slice.label().c_str() );
753 hsize_t n_offset = offset[0];
754 writeFields( h5_config, file_id, n_local, n_global, n_offset, comm_rank,
755 filename_hdf5.str().c_str(), filename_xdmf.str().c_str(),
760 if ( 0 == comm_rank )
761 Impl::writeXdmfFooter( filename_xdmf.str().c_str() );
768template <
class SliceType>
770 hid_t dset_id, hid_t dtype_id, hid_t memspace_id, hid_t filespace_id,
771 hid_t plist_id, std::size_t n_local,
const SliceType&
slice,
772 typename std::enable_if<
773 2 == SliceType::kokkos_view::traits::dimension::rank,
int*>::type = 0 )
776 Kokkos::View<typename SliceType::value_type*, Kokkos::HostSpace> host_view(
777 Kokkos::ViewAllocateWithoutInitializing(
"field" ), n_local );
778 H5Dread( dset_id, dtype_id, memspace_id, filespace_id, plist_id,
782 auto view = Kokkos::create_mirror_view_and_copy(
783 typename SliceType::memory_space(), host_view );
788template <
class SliceType>
790 hid_t dset_id, hid_t dtype_id, hid_t memspace_id, hid_t filespace_id,
791 hid_t plist_id, std::size_t n_local,
const SliceType&
slice,
792 typename std::enable_if<
793 3 == SliceType::kokkos_view::traits::dimension::rank,
int*>::type = 0 )
796 Kokkos::View<
typename SliceType::value_type**, Kokkos::LayoutRight,
798 host_view( Kokkos::ViewAllocateWithoutInitializing(
"field" ), n_local,
800 H5Dread( dset_id, dtype_id, memspace_id, filespace_id, plist_id,
804 auto view = Kokkos::create_mirror_view_and_copy(
805 typename SliceType::memory_space(), host_view );
810template <
class SliceType>
812 hid_t dset_id, hid_t dtype_id, hid_t memspace_id, hid_t filespace_id,
813 hid_t plist_id, std::size_t n_local,
const SliceType&
slice,
814 typename std::enable_if<
815 4 == SliceType::kokkos_view::traits::dimension::rank,
int*>::type = 0 )
818 Kokkos::View<
typename SliceType::value_type***, Kokkos::LayoutRight,
820 host_view( Kokkos::ViewAllocateWithoutInitializing(
"field" ), n_local,
822 H5Dread( dset_id, dtype_id, memspace_id, filespace_id, plist_id,
826 auto view = Kokkos::create_mirror_view_and_copy(
827 typename SliceType::memory_space(), host_view );
843template <
class FieldSliceType>
845 MPI_Comm comm,
const int time_step_index,
846 const std::size_t n_local,
const std::string& dataset_name,
847 double& time, FieldSliceType& field )
849 Kokkos::Profiling::ScopedRegion region(
"Cabana::HDF5ParticleInput" );
860 hsize_t offset[3] = { 0, 0, 0 };
861 hsize_t dimsf[3] = { 0, 0, 0 };
862 hsize_t count[3] = { 0, 0, 0 };
865 MPI_Comm_rank( comm, &comm_rank );
867 MPI_Comm_size( comm, &comm_size );
870 std::stringstream filename_hdf5;
871 filename_hdf5 << prefix <<
"_" << time_step_index <<
".h5";
873 plist_id = H5Pcreate( H5P_FILE_ACCESS );
874 H5Pset_fapl_mpio( plist_id, comm, MPI_INFO_NULL );
876#if H5_VERSION_GE( 1, 10, 0 )
877 if ( h5_config.collective )
879 H5Pset_all_coll_metadata_ops( plist_id, 1 );
884 file_id = H5Fopen( filename_hdf5.str().c_str(), H5F_ACC_RDONLY, plist_id );
885 H5Pclose( plist_id );
888 hid_t attr_id = H5Aopen( file_id,
"Time", H5P_DEFAULT );
889 H5Aread( attr_id, H5T_NATIVE_DOUBLE, &time );
893 dset_id = H5Dopen( file_id, dataset_name.c_str(), H5P_DEFAULT );
896 dtype_id = H5Dget_type( dset_id );
899 filespace_id = H5Dget_space( dset_id );
902 ndims = H5Sget_simple_extent_ndims( filespace_id );
905 H5Sget_simple_extent_dims( filespace_id, dimsf, NULL );
907 std::vector<int> all_offsets( comm_size );
908 all_offsets[comm_rank] = n_local;
910 MPI_Allreduce( MPI_IN_PLACE, all_offsets.data(), comm_size, MPI_INT,
913 for (
int i = 0; i < comm_size; i++ )
917 offset[0] +=
static_cast<hsize_t
>( all_offsets[i] );
920 std::vector<int>().swap( all_offsets );
926 memspace_id = H5Screate_simple( ndims, count, NULL );
928 plist_id = H5Pcreate( H5P_DATASET_XFER );
931 if ( h5_config.collective )
932 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
934 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, NULL, count,
937 readField( dset_id, dtype_id, memspace_id, filespace_id, plist_id, n_local,
940 H5Pclose( plist_id );
941 H5Sclose( memspace_id );
942 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:498
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:769
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:844
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:540
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:932
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:870
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:135
bool meta_collective
Sets metadata I/O mode operations to collective or independent (default)
Definition Cabana_HDF5ParticleOutput.hpp:147
bool align
Set alignment on or off.
Definition Cabana_HDF5ParticleOutput.hpp:140
bool evict_on_close
Cause all metadata for an object to be evicted from the cache.
Definition Cabana_HDF5ParticleOutput.hpp:150
unsigned long threshold
Threshold for aligning file objects.
Definition Cabana_HDF5ParticleOutput.hpp:142
unsigned long alignment
Alignment value.
Definition Cabana_HDF5ParticleOutput.hpp:144
bool collective
I/O transfer mode to collective or independent (default)
Definition Cabana_HDF5ParticleOutput.hpp:137