Cabana 0.8.0-dev
 
Loading...
Searching...
No Matches
Cabana_HDF5ParticleOutput.hpp
Go to the documentation of this file.
1/****************************************************************************
2 * Copyright (c) 2018-2023 by the Cabana authors *
3 * All rights reserved. *
4 * *
5 * This file is part of the Cabana library. Cabana is distributed under a *
6 * BSD 3-clause license. For the licensing terms see the LICENSE file in *
7 * the top-level directory. *
8 * *
9 * SPDX-License-Identifier: BSD-3-Clause *
10 ****************************************************************************/
11
12/****************************************************************************
13 * Copyright (c) 2022 by the Picasso authors *
14 * All rights reserved. *
15 * *
16 * This file is part of the Picasso library. Picasso is distributed under a *
17 * BSD 3-clause license. For the licensing terms see the LICENSE file in *
18 * the top-level directory. *
19 * *
20 * SPDX-License-Identifier: BSD-3-Clause *
21 ****************************************************************************/
22
27
28#ifndef CABANA_HDF5PARTICLEOUTPUT_HPP
29#define CABANA_HDF5PARTICLEOUTPUT_HPP
30
31#include <Kokkos_Core.hpp>
32#include <Kokkos_Profiling_ScopedRegion.hpp>
33
34#include <hdf5.h>
35#ifdef H5_HAVE_SUBFILING_VFD
36#include "H5FDioc.h" /* Private header for the IOC VFD */
37#include "H5FDsubfiling.h" /* Private header for the subfiling VFD */
38#endif
39
40#include <mpi.h>
41
42#include <fstream>
43#include <iostream>
44#include <limits>
45#include <sstream>
46#include <string>
47#include <type_traits>
48#include <vector>
49
50namespace Cabana
51{
52namespace Experimental
53{
54namespace HDF5ParticleOutput
55{
56
57namespace Impl
58{
59// XDMF file creation routines
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 )
65{
66 std::ofstream xdmf_file( xml_file_name, std::ios::trunc );
67 // Set precision to guarantee that conversion to text and back is exact.
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" )
79 << "\">\n";
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";
86 xdmf_file.close();
87}
88
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 )
94{
95 std::string AttributeType = "\"Scalar\"";
96 if ( dims2 != 0 )
97 AttributeType = "\"Tensor\"";
98 else if ( dims1 != 0 )
99 AttributeType = "\"Vector\"";
100
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=\""
106 << dims0;
107 if ( dims1 != 0 )
108 xdmf_file << " " << dims1;
109 if ( dims2 != 0 )
110 xdmf_file << " " << dims2;
111 xdmf_file << "\" DataType=\"" << dtype << "\" Precision=\"" << precision
112 << "\"";
113 xdmf_file << " Format=\"HDF\"> " << h5_file_name << ":/" << dataitem;
114 xdmf_file << " </DataItem>\n";
115 xdmf_file << " </Attribute>\n";
116 xdmf_file.close();
117}
118
119inline void writeXdmfFooter( const char* xml_file_name )
120{
121 std::ofstream xdmf_file( xml_file_name, std::ios::app );
122 xdmf_file << " </Grid>\n";
123 xdmf_file << " </Domain>\n</Xdmf>\n";
124 xdmf_file.close();
125}
127} // namespace Impl
128
141{
143 bool collective = false;
144
146 bool align = false;
148 unsigned long threshold = 0;
150 unsigned long alignment = 16777216;
151
153 bool meta_collective = true;
154
156 bool evict_on_close = false;
157
158#ifdef H5_HAVE_SUBFILING_VFD
159
161 bool subfiling = false;
162
163 // Optional subfiling file driver configuration parameters
164
166 int64_t subfiling_stripe_size = H5FD_SUBFILING_DEFAULT_STRIPE_SIZE;
167
169 int32_t subfiling_stripe_count = H5FD_SUBFILING_DEFAULT_STRIPE_COUNT;
170
172 int subfiling_ioc_selection = SELECT_IOC_ONE_PER_NODE;
173
175 int32_t subfiling_thread_pool_size = H5FD_IOC_DEFAULT_THREAD_POOL_SIZE;
176#endif
177};
178
180// Format traits for both HDF5 and XDMF.
181template <typename T>
182struct HDF5Traits;
183
184template <>
185struct HDF5Traits<int>
186{
187 static hid_t type( std::string* dtype, uint* precision )
188 {
189 *dtype = "Int";
190 *precision = sizeof( int );
191 return H5T_NATIVE_INT;
192 }
193};
194
195template <>
196struct HDF5Traits<unsigned int>
197{
198 static hid_t type( std::string* dtype, uint* precision )
199 {
200 *dtype = "UInt";
201 *precision = sizeof( unsigned int );
202 return H5T_NATIVE_UINT;
203 }
204};
205
206template <>
207struct HDF5Traits<long>
208{
209 static hid_t type( std::string* dtype, uint* precision )
210 {
211 *dtype = "Int";
212 *precision = sizeof( long );
213 return H5T_NATIVE_LONG;
214 }
215};
216
217template <>
218struct HDF5Traits<unsigned long>
219{
220 static hid_t type( std::string* dtype, uint* precision )
221 {
222 *dtype = "UInt";
223 *precision = sizeof( unsigned long );
224 return H5T_NATIVE_ULONG;
225 }
226};
227
228template <>
229struct HDF5Traits<float>
230{
231 static hid_t type( std::string* dtype, uint* precision )
232 {
233 *dtype = "Float";
234 *precision = sizeof( float );
235 return H5T_NATIVE_FLOAT;
236 }
237};
238
239template <>
240struct HDF5Traits<double>
241{
242 static hid_t type( std::string* dtype, uint* precision )
243 {
244 *dtype = "Float";
245 *precision = sizeof( double );
246 return H5T_NATIVE_DOUBLE;
247 }
248};
249
250namespace Impl
251{
252//---------------------------------------------------------------------------//
253// HDF5 (XDMF) Particle Field Output.
254//---------------------------------------------------------------------------//
255
256// Rank-0 field
257template <class SliceType>
258void writeFields(
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 )
265{
266 hid_t plist_id;
267 hid_t dset_id;
268 hid_t dcpl_id;
269 hid_t filespace_id;
270 hid_t memspace_id;
271
272 // HDF5 hyperslab parameters
273 hsize_t offset[1];
274 hsize_t dimsf[1];
275 hsize_t count[1];
276
277 offset[0] = n_offset;
278 count[0] = n_local;
279 dimsf[0] = n_global;
280
281 // Reorder in a contiguous blocked format.
282 Kokkos::View<typename SliceType::value_type*,
283 typename SliceType::memory_space>
284 view( Kokkos::ViewAllocateWithoutInitializing( "field" ), n_local );
285 copySliceToView( view, slice, 0, n_local );
286
287 // Mirror the field to the host.
288 auto host_view =
289 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), view );
290
291 std::string dtype;
292 uint precision = 0;
293 hid_t type_id =
294 HDF5Traits<typename SliceType::value_type>::type( &dtype, &precision );
295
296 filespace_id = H5Screate_simple( 1, dimsf, nullptr );
297
298 dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
299 H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );
300
301 dset_id = H5Dcreate( file_id, slice.label().c_str(), type_id, filespace_id,
302 H5P_DEFAULT, dcpl_id, H5P_DEFAULT );
303
304 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, nullptr, count,
305 nullptr );
306
307 memspace_id = H5Screate_simple( 1, count, nullptr );
308
309 plist_id = H5Pcreate( H5P_DATASET_XFER );
310 // Default IO in HDF5 is independent
311 if ( h5_config.collective )
312 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
313
314 H5Dwrite( dset_id, type_id, memspace_id, filespace_id, plist_id,
315 host_view.data() );
316
317 H5Pclose( plist_id );
318 H5Pclose( dcpl_id );
319 H5Sclose( memspace_id );
320 H5Dclose( dset_id );
321 H5Sclose( filespace_id );
322
323 if ( 0 == comm_rank )
324 {
325 hsize_t zero = 0;
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() );
329 }
330}
331
332// Rank-1 field
333template <class SliceType>
334void writeFields(
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 )
341{
342 hid_t plist_id;
343 hid_t dset_id;
344 hid_t dcpl_id;
345 hid_t filespace_id;
346 hid_t memspace_id;
347
348 // Reorder in a contiguous blocked format.
349 Kokkos::View<typename SliceType::value_type**, Kokkos::LayoutRight,
350 typename SliceType::memory_space>
351 view( Kokkos::ViewAllocateWithoutInitializing( "field" ), n_local,
352 slice.extent( 2 ) );
353 copySliceToView( view, slice, 0, n_local );
354
355 // Mirror the field to the host.
356 auto host_view =
357 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), view );
358
359 hsize_t offset[2];
360 hsize_t dimsf[2];
361 hsize_t dimsm[2];
362 hsize_t count[2];
363
364 dimsf[0] = n_global;
365 dimsf[1] = host_view.extent( 1 );
366 dimsm[0] = n_local;
367 dimsm[1] = host_view.extent( 1 );
368
369 offset[0] = n_offset;
370 offset[1] = 0;
371
372 count[0] = dimsm[0];
373 count[1] = dimsm[1];
374
375 std::string dtype;
376 uint precision;
377 hid_t type_id =
378 HDF5Traits<typename SliceType::value_type>::type( &dtype, &precision );
379
380 filespace_id = H5Screate_simple( 2, dimsf, nullptr );
381
382 dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
383 H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );
384
385 dset_id = H5Dcreate( file_id, slice.label().c_str(), type_id, filespace_id,
386 H5P_DEFAULT, dcpl_id, H5P_DEFAULT );
387
388 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, nullptr, count,
389 nullptr );
390
391 memspace_id = H5Screate_simple( 2, dimsm, nullptr );
392 plist_id = H5Pcreate( H5P_DATASET_XFER );
393 // Default IO in HDF5 is independent
394 if ( h5_config.collective )
395 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
396
397 H5Dwrite( dset_id, type_id, memspace_id, filespace_id, plist_id,
398 host_view.data() );
399
400 H5Pclose( plist_id );
401 H5Pclose( dcpl_id );
402 H5Sclose( memspace_id );
403 H5Dclose( dset_id );
404 H5Sclose( filespace_id );
405
406 if ( 0 == comm_rank )
407 {
408 hsize_t zero = 0;
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() );
412 }
413}
414
415// Rank-2 field
416template <class SliceType>
417void writeFields(
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 )
424{
425 hid_t plist_id;
426 hid_t dset_id;
427 hid_t dcpl_id;
428 hid_t filespace_id;
429 hid_t memspace_id;
430
431 // Reorder in a contiguous blocked format.
432 Kokkos::View<typename SliceType::value_type***, Kokkos::LayoutRight,
433 typename SliceType::memory_space>
434 view( Kokkos::ViewAllocateWithoutInitializing( "field" ), n_local,
435 slice.extent( 2 ), slice.extent( 3 ) );
436 copySliceToView( view, slice, 0, n_local );
437
438 // Mirror the field to the host.
439 auto host_view =
440 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), view );
441
442 hsize_t offset[3];
443 hsize_t dimsf[3];
444 hsize_t dimsm[3];
445 hsize_t count[3];
446
447 dimsf[0] = n_global;
448 dimsf[1] = host_view.extent( 1 );
449 dimsf[2] = host_view.extent( 2 );
450 dimsm[0] = n_local;
451 dimsm[1] = host_view.extent( 1 );
452 dimsm[2] = host_view.extent( 2 );
453
454 offset[0] = n_offset;
455 offset[1] = 0;
456 offset[2] = 0;
457
458 count[0] = dimsm[0];
459 count[1] = dimsm[1];
460 count[2] = dimsm[2];
461
462 std::string dtype;
463 uint precision;
464 hid_t type_id =
465 HDF5Traits<typename SliceType::value_type>::type( &dtype, &precision );
466
467 filespace_id = H5Screate_simple( 3, dimsf, nullptr );
468
469 dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
470 H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );
471
472 dset_id = H5Dcreate( file_id, slice.label().c_str(), type_id, filespace_id,
473 H5P_DEFAULT, dcpl_id, H5P_DEFAULT );
474
475 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, nullptr, count,
476 nullptr );
477
478 memspace_id = H5Screate_simple( 3, dimsm, nullptr );
479 plist_id = H5Pcreate( H5P_DATASET_XFER );
480 // Default IO in HDF5 is independent
481 if ( h5_config.collective )
482 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
483
484 H5Dwrite( dset_id, type_id, memspace_id, filespace_id, plist_id,
485 host_view.data() );
486
487 H5Pclose( plist_id );
488 H5Pclose( dcpl_id );
489 H5Sclose( memspace_id );
490 H5Dclose( dset_id );
491 H5Sclose( filespace_id );
492
493 if ( 0 == comm_rank )
494 {
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() );
498 }
499}
501} // namespace Impl
502
504inline void writeFields( HDF5Config, hid_t, std::size_t, std::size_t, hsize_t,
505 int, const char*, const char* )
506{
507}
508
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 )
515{
516 Impl::writeFields( h5_config, file_id, n_local, n_global, n_offset,
517 comm_rank, filename_hdf5, filename_xdmf, slice );
518}
519
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 )
526{
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... );
531}
532
533//---------------------------------------------------------------------------//
545template <class CoordSliceType, class... FieldSliceTypes>
546void writeTimeStep( HDF5Config h5_config, const std::string& prefix,
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 )
551{
552 Kokkos::Profiling::ScopedRegion region( "Cabana::HDF5ParticleOutput" );
553
554 hid_t plist_id;
555 hid_t dset_id;
556 hid_t dcpl_id;
557 hid_t file_id;
558 hid_t filespace_id;
559 hid_t memspace_id;
560
561 // HDF5 hyperslab parameters
562 hsize_t offset[2];
563 hsize_t dimsf[2];
564 hsize_t count[2];
565
566 int comm_rank;
567 MPI_Comm_rank( comm, &comm_rank );
568 int comm_size;
569 MPI_Comm_size( comm, &comm_size );
570
571 // Compose a data file name.
572 std::stringstream filename_hdf5;
573 filename_hdf5 << prefix << "_" << time_step_index << ".h5";
574
575 std::stringstream filename_xdmf;
576 filename_xdmf << prefix << "_" << time_step_index << ".xmf";
577
578 plist_id = H5Pcreate( H5P_FILE_ACCESS );
579 H5Pset_libver_bounds( plist_id, H5F_LIBVER_LATEST, H5F_LIBVER_LATEST );
580
581#if H5_VERSION_GE( 1, 10, 1 )
582 if ( h5_config.evict_on_close )
583 {
584 H5Pset_evict_on_close( plist_id, true );
585 }
586#endif
587
588#if H5_VERSION_GE( 1, 10, 0 )
589 if ( h5_config.collective )
590 {
591 H5Pset_all_coll_metadata_ops( plist_id, true );
592 H5Pset_coll_metadata_write( plist_id, true );
593 }
594#endif
595
596 if ( h5_config.align )
597 H5Pset_alignment( plist_id, h5_config.threshold, h5_config.alignment );
598
599#ifdef H5_HAVE_SUBFILING_VFD
600 if ( h5_config.subfiling )
601 {
602 H5FD_subfiling_config_t subfiling_config;
603 H5FD_ioc_config_t ioc_config;
604
605 H5FD_subfiling_config_t* subfiling_ptr = nullptr;
606 H5FD_ioc_config_t* ioc_ptr = nullptr;
607
608 // Get the default subfiling configuration parameters
609 hid_t fapl_id = H5I_INVALID_HID;
610
611 fapl_id = H5Pcreate( H5P_FILE_ACCESS );
612 H5Pget_fapl_subfiling( fapl_id, &subfiling_config );
613
614 if ( h5_config.subfiling_stripe_size !=
615 subfiling_config.shared_cfg.stripe_size )
616 {
617 subfiling_config.shared_cfg.stripe_size =
618 h5_config.subfiling_stripe_size;
619 if ( subfiling_ptr == nullptr )
620 subfiling_ptr = &subfiling_config;
621 }
622 if ( h5_config.subfiling_stripe_count !=
623 subfiling_config.shared_cfg.stripe_count )
624 {
625 subfiling_config.shared_cfg.stripe_count =
626 h5_config.subfiling_stripe_count;
627 if ( subfiling_ptr == nullptr )
628 subfiling_ptr = &subfiling_config;
629 }
630 if ( h5_config.subfiling_ioc_selection !=
631 (int)subfiling_config.shared_cfg.ioc_selection )
632 {
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;
637 }
638 if ( h5_config.subfiling_thread_pool_size !=
639 H5FD_IOC_DEFAULT_THREAD_POOL_SIZE )
640 {
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;
645 }
646 H5Pclose( fapl_id );
647
648 H5Pset_mpi_params( plist_id, comm, MPI_INFO_NULL );
649
650 if ( ioc_ptr != nullptr )
651 H5Pset_fapl_ioc( subfiling_config.ioc_fapl_id, ioc_ptr );
652
653 H5Pset_fapl_subfiling( plist_id, subfiling_ptr );
654 }
655 else
656#endif
657 {
658 H5Pset_fapl_mpio( plist_id, comm, MPI_INFO_NULL );
659 }
660
661 file_id = H5Fcreate( filename_hdf5.str().c_str(), H5F_ACC_TRUNC,
662 H5P_DEFAULT, plist_id );
663 H5Pclose( plist_id );
664
665 // Write current simulation time
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 );
670 H5Aclose( attr_id );
671 H5Sclose( fspace );
672
673 // Reorder the coordinates in a blocked format.
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 );
685 } );
686
687 // Mirror the coordinates to the host.
688 auto host_coords =
689 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), coords_view );
690
691 std::vector<int> all_offsets( comm_size );
692 all_offsets[comm_rank] = n_local;
693
694 MPI_Allreduce( MPI_IN_PLACE, all_offsets.data(), comm_size, MPI_INT,
695 MPI_SUM, comm );
696
697 offset[0] = 0;
698 offset[1] = 0;
699
700 size_t n_global = 0;
701 for ( int i = 0; i < comm_size; i++ )
702 {
703 if ( i < comm_rank )
704 {
705 offset[0] += static_cast<hsize_t>( all_offsets[i] );
706 }
707 n_global += (size_t)all_offsets[i];
708 }
709 std::vector<int>().swap( all_offsets );
710
711 dimsf[0] = n_global;
712 dimsf[1] = coords_slice.extent( 2 );
713
714 filespace_id = H5Screate_simple( 2, dimsf, nullptr );
715
716 count[0] = n_local;
717 count[1] = coords_slice.extent( 2 );
718
719 memspace_id = H5Screate_simple( 2, count, nullptr );
720
721 plist_id = H5Pcreate( H5P_DATASET_XFER );
722
723 // Default IO in HDF5 is independent
724 if ( h5_config.collective )
725 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
726
727 std::string dtype;
728 uint precision;
729 hid_t type_id = HDF5Traits<typename CoordSliceType::value_type>::type(
730 &dtype, &precision );
731
732 dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
733 H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );
734
735 dset_id = H5Dcreate( file_id, coords_slice.label().c_str(), type_id,
736 filespace_id, H5P_DEFAULT, dcpl_id, H5P_DEFAULT );
737
738 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, nullptr, count,
739 nullptr );
740
741 H5Dwrite( dset_id, type_id, memspace_id, filespace_id, plist_id,
742 host_coords.data() );
743 H5Dclose( dset_id );
744
745 H5Pclose( plist_id );
746 H5Pclose( dcpl_id );
747 H5Sclose( filespace_id );
748 H5Sclose( memspace_id );
749
750 if ( 0 == comm_rank )
751 {
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() );
756 }
757
758 // Add variables.
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(),
762 fields... );
763
764 H5Fclose( file_id );
765
766 if ( 0 == comm_rank )
767 Impl::writeXdmfFooter( filename_xdmf.str().c_str() );
768}
769
770//---------------------------------------------------------------------------//
771// HDF5 (XDMF) Particle Field Input.
772//---------------------------------------------------------------------------//
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 )
780{
781 // Read the field into a View.
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,
785 host_view.data() );
786
787 // Mirror the field and copy.
788 auto view = Kokkos::create_mirror_view_and_copy(
789 typename SliceType::memory_space(), host_view );
790 copyViewToSlice( slice, view, 0, n_local );
791}
792
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 )
800{
801 // Read the field into a View.
802 Kokkos::View<typename SliceType::value_type**, Kokkos::LayoutRight,
803 Kokkos::HostSpace>
804 host_view( Kokkos::ViewAllocateWithoutInitializing( "field" ), n_local,
805 slice.extent( 2 ) );
806 H5Dread( dset_id, dtype_id, memspace_id, filespace_id, plist_id,
807 host_view.data() );
808
809 // Mirror the field and copy.
810 auto view = Kokkos::create_mirror_view_and_copy(
811 typename SliceType::memory_space(), host_view );
812 copyViewToSlice( slice, view, 0, n_local );
813}
814
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 )
822{
823 // Read the field into a View.
824 Kokkos::View<typename SliceType::value_type***, Kokkos::LayoutRight,
825 Kokkos::HostSpace>
826 host_view( Kokkos::ViewAllocateWithoutInitializing( "field" ), n_local,
827 slice.extent( 2 ), slice.extent( 3 ) );
828 H5Dread( dset_id, dtype_id, memspace_id, filespace_id, plist_id,
829 host_view.data() );
830
831 // Mirror the field and copy.
832 auto view = Kokkos::create_mirror_view_and_copy(
833 typename SliceType::memory_space(), host_view );
834 copyViewToSlice( slice, view, 0, n_local );
835}
836
837//---------------------------------------------------------------------------//
849template <class FieldSliceType>
850void readTimeStep( HDF5Config h5_config, const std::string& prefix,
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 )
854{
855 Kokkos::Profiling::ScopedRegion region( "Cabana::HDF5ParticleInput" );
856
857 hid_t plist_id;
858 hid_t dset_id;
859 hid_t file_id;
860 hid_t dtype_id;
861 hid_t filespace_id;
862 hid_t memspace_id;
863 int ndims;
864
865 // HDF5 hyperslab parameters
866 hsize_t offset[3] = { 0, 0, 0 };
867 hsize_t dimsf[3] = { 0, 0, 0 };
868 hsize_t count[3] = { 0, 0, 0 };
869
870 int comm_rank;
871 MPI_Comm_rank( comm, &comm_rank );
872 int comm_size;
873 MPI_Comm_size( comm, &comm_size );
874
875 // Retrieve data file name.
876 std::stringstream filename_hdf5;
877 filename_hdf5 << prefix << "_" << time_step_index << ".h5";
878
879 plist_id = H5Pcreate( H5P_FILE_ACCESS );
880 H5Pset_fapl_mpio( plist_id, comm, MPI_INFO_NULL );
881
882#if H5_VERSION_GE( 1, 10, 0 )
883 if ( h5_config.collective )
884 {
885 H5Pset_all_coll_metadata_ops( plist_id, true );
886 }
887#endif
888
889 // Open the HDF5 file.
890 file_id = H5Fopen( filename_hdf5.str().c_str(), H5F_ACC_RDONLY, plist_id );
891 H5Pclose( plist_id );
892
893 // Get current simulation time associated with time_step_index
894 hid_t attr_id = H5Aopen( file_id, "Time", H5P_DEFAULT );
895 H5Aread( attr_id, H5T_NATIVE_DOUBLE, &time );
896 H5Aclose( attr_id );
897
898 // Open the dataset.
899 dset_id = H5Dopen( file_id, dataset_name.c_str(), H5P_DEFAULT );
900
901 // Get the datatype of the dataset.
902 dtype_id = H5Dget_type( dset_id );
903
904 // Get the dataspace of the dataset.
905 filespace_id = H5Dget_space( dset_id );
906
907 // Get the rank of the dataspace.
908 ndims = H5Sget_simple_extent_ndims( filespace_id );
909
910 // Get the extents of the file dataspace.
911 H5Sget_simple_extent_dims( filespace_id, dimsf, nullptr );
912
913 std::vector<int> all_offsets( comm_size );
914 all_offsets[comm_rank] = n_local;
915
916 MPI_Allreduce( MPI_IN_PLACE, all_offsets.data(), comm_size, MPI_INT,
917 MPI_SUM, comm );
918
919 for ( int i = 0; i < comm_size; i++ )
920 {
921 if ( i < comm_rank )
922 {
923 offset[0] += static_cast<hsize_t>( all_offsets[i] );
924 }
925 }
926 std::vector<int>().swap( all_offsets );
927
928 count[0] = n_local;
929 count[1] = dimsf[1];
930 count[2] = dimsf[2];
931
932 memspace_id = H5Screate_simple( ndims, count, nullptr );
933
934 plist_id = H5Pcreate( H5P_DATASET_XFER );
935
936 // Default IO in HDF5 is independent
937 if ( h5_config.collective )
938 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
939
940 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, nullptr, count,
941 nullptr );
942
943 readField( dset_id, dtype_id, memspace_id, filespace_id, plist_id, n_local,
944 field );
945
946 H5Pclose( plist_id );
947 H5Sclose( memspace_id );
948 H5Sclose( filespace_id );
949 H5Dclose( dset_id );
950 H5Fclose( file_id );
951}
952
953//---------------------------------------------------------------------------//
954
955} // namespace HDF5ParticleOutput
956} // namespace Experimental
957} // end namespace Cabana
958
959#endif // CABANA_HDF5PARTICLEOUTPUT_HPP
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