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 <sstream>
45#include <string>
46#include <type_traits>
47#include <vector>
48
49namespace Cabana
50{
51namespace Experimental
52{
53namespace HDF5ParticleOutput
54{
55
56namespace Impl
57{
58// XDMF file creation routines
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 )
63{
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";
80 xdmf_file.close();
81}
82
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 )
88{
89 std::string AttributeType = "\"Scalar\"";
90 if ( dims2 != 0 )
91 AttributeType = "\"Tensor\"";
92 else if ( dims1 != 0 )
93 AttributeType = "\"Vector\"";
94
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=\""
100 << dims0;
101 if ( dims1 != 0 )
102 xdmf_file << " " << dims1;
103 if ( dims2 != 0 )
104 xdmf_file << " " << dims2;
105 xdmf_file << "\" DataType=\"" << dtype << "\" Precision=\"" << precision
106 << "\"";
107 xdmf_file << " Format=\"HDF\"> " << h5_file_name << ":/" << dataitem;
108 xdmf_file << " </DataItem>\n";
109 xdmf_file << " </Attribute>\n";
110 xdmf_file.close();
111}
112
113inline void writeXdmfFooter( const char* xml_file_name )
114{
115 std::ofstream xdmf_file( xml_file_name, std::ios::app );
116 xdmf_file << " </Grid>\n";
117 xdmf_file << " </Domain>\n</Xdmf>\n";
118 xdmf_file.close();
119}
121} // namespace Impl
122
135{
137 bool collective = false;
138
140 bool align = false;
142 unsigned long threshold = 0;
144 unsigned long alignment = 16777216;
145
147 bool meta_collective = true;
148
150 bool evict_on_close = false;
151
152#ifdef H5_HAVE_SUBFILING_VFD
153
155 bool subfiling = false;
156
157 // Optional subfiling file driver configuration parameters
158
160 int64_t subfiling_stripe_size = H5FD_SUBFILING_DEFAULT_STRIPE_SIZE;
161
163 int32_t subfiling_stripe_count = H5FD_SUBFILING_DEFAULT_STRIPE_COUNT;
164
166 int subfiling_ioc_selection = SELECT_IOC_ONE_PER_NODE;
167
169 int32_t subfiling_thread_pool_size = H5FD_IOC_DEFAULT_THREAD_POOL_SIZE;
170#endif
171};
172
174// Format traits for both HDF5 and XDMF.
175template <typename T>
176struct HDF5Traits;
177
178template <>
179struct HDF5Traits<int>
180{
181 static hid_t type( std::string* dtype, uint* precision )
182 {
183 *dtype = "Int";
184 *precision = sizeof( int );
185 return H5T_NATIVE_INT;
186 }
187};
188
189template <>
190struct HDF5Traits<unsigned int>
191{
192 static hid_t type( std::string* dtype, uint* precision )
193 {
194 *dtype = "UInt";
195 *precision = sizeof( unsigned int );
196 return H5T_NATIVE_UINT;
197 }
198};
199
200template <>
201struct HDF5Traits<long>
202{
203 static hid_t type( std::string* dtype, uint* precision )
204 {
205 *dtype = "Int";
206 *precision = sizeof( long );
207 return H5T_NATIVE_LONG;
208 }
209};
210
211template <>
212struct HDF5Traits<unsigned long>
213{
214 static hid_t type( std::string* dtype, uint* precision )
215 {
216 *dtype = "UInt";
217 *precision = sizeof( unsigned long );
218 return H5T_NATIVE_ULONG;
219 }
220};
221
222template <>
223struct HDF5Traits<float>
224{
225 static hid_t type( std::string* dtype, uint* precision )
226 {
227 *dtype = "Float";
228 *precision = sizeof( float );
229 return H5T_NATIVE_FLOAT;
230 }
231};
232
233template <>
234struct HDF5Traits<double>
235{
236 static hid_t type( std::string* dtype, uint* precision )
237 {
238 *dtype = "Float";
239 *precision = sizeof( double );
240 return H5T_NATIVE_DOUBLE;
241 }
242};
243
244namespace Impl
245{
246//---------------------------------------------------------------------------//
247// HDF5 (XDMF) Particle Field Output.
248//---------------------------------------------------------------------------//
249
250// Rank-0 field
251template <class SliceType>
252void writeFields(
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 )
259{
260 hid_t plist_id;
261 hid_t dset_id;
262 hid_t dcpl_id;
263 hid_t filespace_id;
264 hid_t memspace_id;
265
266 // HDF5 hyperslab parameters
267 hsize_t offset[1];
268 hsize_t dimsf[1];
269 hsize_t count[1];
270
271 offset[0] = n_offset;
272 count[0] = n_local;
273 dimsf[0] = n_global;
274
275 // Reorder in a contiguous blocked format.
276 Kokkos::View<typename SliceType::value_type*,
277 typename SliceType::memory_space>
278 view( Kokkos::ViewAllocateWithoutInitializing( "field" ), n_local );
279 copySliceToView( view, slice, 0, n_local );
280
281 // Mirror the field to the host.
282 auto host_view =
283 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), view );
284
285 std::string dtype;
286 uint precision = 0;
287 hid_t type_id =
288 HDF5Traits<typename SliceType::value_type>::type( &dtype, &precision );
289
290 filespace_id = H5Screate_simple( 1, dimsf, NULL );
291
292 dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
293 H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );
294
295 dset_id = H5Dcreate( file_id, slice.label().c_str(), type_id, filespace_id,
296 H5P_DEFAULT, dcpl_id, H5P_DEFAULT );
297
298 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, NULL, count,
299 NULL );
300
301 memspace_id = H5Screate_simple( 1, count, NULL );
302
303 plist_id = H5Pcreate( H5P_DATASET_XFER );
304 // Default IO in HDF5 is independent
305 if ( h5_config.collective )
306 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
307
308 H5Dwrite( dset_id, type_id, memspace_id, filespace_id, plist_id,
309 host_view.data() );
310
311 H5Pclose( plist_id );
312 H5Pclose( dcpl_id );
313 H5Sclose( memspace_id );
314 H5Dclose( dset_id );
315 H5Sclose( filespace_id );
316
317 if ( 0 == comm_rank )
318 {
319 hsize_t zero = 0;
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() );
323 }
324}
325
326// Rank-1 field
327template <class SliceType>
328void writeFields(
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 )
335{
336 hid_t plist_id;
337 hid_t dset_id;
338 hid_t dcpl_id;
339 hid_t filespace_id;
340 hid_t memspace_id;
341
342 // Reorder in a contiguous blocked format.
343 Kokkos::View<typename SliceType::value_type**, Kokkos::LayoutRight,
344 typename SliceType::memory_space>
345 view( Kokkos::ViewAllocateWithoutInitializing( "field" ), n_local,
346 slice.extent( 2 ) );
347 copySliceToView( view, slice, 0, n_local );
348
349 // Mirror the field to the host.
350 auto host_view =
351 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), view );
352
353 hsize_t offset[2];
354 hsize_t dimsf[2];
355 hsize_t dimsm[2];
356 hsize_t count[2];
357
358 dimsf[0] = n_global;
359 dimsf[1] = host_view.extent( 1 );
360 dimsm[0] = n_local;
361 dimsm[1] = host_view.extent( 1 );
362
363 offset[0] = n_offset;
364 offset[1] = 0;
365
366 count[0] = dimsm[0];
367 count[1] = dimsm[1];
368
369 std::string dtype;
370 uint precision;
371 hid_t type_id =
372 HDF5Traits<typename SliceType::value_type>::type( &dtype, &precision );
373
374 filespace_id = H5Screate_simple( 2, dimsf, NULL );
375
376 dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
377 H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );
378
379 dset_id = H5Dcreate( file_id, slice.label().c_str(), type_id, filespace_id,
380 H5P_DEFAULT, dcpl_id, H5P_DEFAULT );
381
382 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, NULL, count,
383 NULL );
384
385 memspace_id = H5Screate_simple( 2, dimsm, NULL );
386 plist_id = H5Pcreate( H5P_DATASET_XFER );
387 // Default IO in HDF5 is independent
388 if ( h5_config.collective )
389 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
390
391 H5Dwrite( dset_id, type_id, memspace_id, filespace_id, plist_id,
392 host_view.data() );
393
394 H5Pclose( plist_id );
395 H5Pclose( dcpl_id );
396 H5Sclose( memspace_id );
397 H5Dclose( dset_id );
398 H5Sclose( filespace_id );
399
400 if ( 0 == comm_rank )
401 {
402 hsize_t zero = 0;
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() );
406 }
407}
408
409// Rank-2 field
410template <class SliceType>
411void writeFields(
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 )
418{
419 hid_t plist_id;
420 hid_t dset_id;
421 hid_t dcpl_id;
422 hid_t filespace_id;
423 hid_t memspace_id;
424
425 // Reorder in a contiguous blocked format.
426 Kokkos::View<typename SliceType::value_type***, Kokkos::LayoutRight,
427 typename SliceType::memory_space>
428 view( Kokkos::ViewAllocateWithoutInitializing( "field" ), n_local,
429 slice.extent( 2 ), slice.extent( 3 ) );
430 copySliceToView( view, slice, 0, n_local );
431
432 // Mirror the field to the host.
433 auto host_view =
434 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), view );
435
436 hsize_t offset[3];
437 hsize_t dimsf[3];
438 hsize_t dimsm[3];
439 hsize_t count[3];
440
441 dimsf[0] = n_global;
442 dimsf[1] = host_view.extent( 1 );
443 dimsf[2] = host_view.extent( 2 );
444 dimsm[0] = n_local;
445 dimsm[1] = host_view.extent( 1 );
446 dimsm[2] = host_view.extent( 2 );
447
448 offset[0] = n_offset;
449 offset[1] = 0;
450 offset[2] = 0;
451
452 count[0] = dimsm[0];
453 count[1] = dimsm[1];
454 count[2] = dimsm[2];
455
456 std::string dtype;
457 uint precision;
458 hid_t type_id =
459 HDF5Traits<typename SliceType::value_type>::type( &dtype, &precision );
460
461 filespace_id = H5Screate_simple( 3, dimsf, NULL );
462
463 dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
464 H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );
465
466 dset_id = H5Dcreate( file_id, slice.label().c_str(), type_id, filespace_id,
467 H5P_DEFAULT, dcpl_id, H5P_DEFAULT );
468
469 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, NULL, count,
470 NULL );
471
472 memspace_id = H5Screate_simple( 3, dimsm, NULL );
473 plist_id = H5Pcreate( H5P_DATASET_XFER );
474 // Default IO in HDF5 is independent
475 if ( h5_config.collective )
476 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
477
478 H5Dwrite( dset_id, type_id, memspace_id, filespace_id, plist_id,
479 host_view.data() );
480
481 H5Pclose( plist_id );
482 H5Pclose( dcpl_id );
483 H5Sclose( memspace_id );
484 H5Dclose( dset_id );
485 H5Sclose( filespace_id );
486
487 if ( 0 == comm_rank )
488 {
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() );
492 }
493}
495} // namespace Impl
496
498inline void writeFields( HDF5Config, hid_t, std::size_t, std::size_t, hsize_t,
499 int, const char*, const char* )
500{
501}
502
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 )
509{
510 Impl::writeFields( h5_config, file_id, n_local, n_global, n_offset,
511 comm_rank, filename_hdf5, filename_xdmf, slice );
512}
513
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 )
520{
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... );
525}
526
527//---------------------------------------------------------------------------//
539template <class CoordSliceType, class... FieldSliceTypes>
540void writeTimeStep( HDF5Config h5_config, const std::string& prefix,
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 )
545{
546 Kokkos::Profiling::ScopedRegion region( "Cabana::HDF5ParticleOutput" );
547
548 hid_t plist_id;
549 hid_t dset_id;
550 hid_t dcpl_id;
551 hid_t file_id;
552 hid_t filespace_id;
553 hid_t memspace_id;
554
555 // HDF5 hyperslab parameters
556 hsize_t offset[2];
557 hsize_t dimsf[2];
558 hsize_t count[2];
559
560 int comm_rank;
561 MPI_Comm_rank( comm, &comm_rank );
562 int comm_size;
563 MPI_Comm_size( comm, &comm_size );
564
565 // Compose a data file name.
566 std::stringstream filename_hdf5;
567 filename_hdf5 << prefix << "_" << time_step_index << ".h5";
568
569 std::stringstream filename_xdmf;
570 filename_xdmf << prefix << "_" << time_step_index << ".xmf";
571
572 plist_id = H5Pcreate( H5P_FILE_ACCESS );
573 H5Pset_libver_bounds( plist_id, H5F_LIBVER_LATEST, H5F_LIBVER_LATEST );
574
575#if H5_VERSION_GE( 1, 10, 1 )
576 if ( h5_config.evict_on_close )
577 {
578 H5Pset_evict_on_close( plist_id, (hbool_t)1 );
579 }
580#endif
581
582#if H5_VERSION_GE( 1, 10, 0 )
583 if ( h5_config.collective )
584 {
585 H5Pset_all_coll_metadata_ops( plist_id, 1 );
586 H5Pset_coll_metadata_write( plist_id, 1 );
587 }
588#endif
589
590 if ( h5_config.align )
591 H5Pset_alignment( plist_id, h5_config.threshold, h5_config.alignment );
592
593#ifdef H5_HAVE_SUBFILING_VFD
594 if ( h5_config.subfiling )
595 {
596 H5FD_subfiling_config_t subfiling_config;
597 H5FD_ioc_config_t ioc_config;
598
599 H5FD_subfiling_config_t* subfiling_ptr = NULL;
600 H5FD_ioc_config_t* ioc_ptr = NULL;
601
602 // Get the default subfiling configuration parameters
603 hid_t fapl_id = H5I_INVALID_HID;
604
605 fapl_id = H5Pcreate( H5P_FILE_ACCESS );
606 H5Pget_fapl_subfiling( fapl_id, &subfiling_config );
607
608 if ( h5_config.subfiling_stripe_size !=
609 subfiling_config.shared_cfg.stripe_size )
610 {
611 subfiling_config.shared_cfg.stripe_size =
612 h5_config.subfiling_stripe_size;
613 if ( subfiling_ptr == NULL )
614 subfiling_ptr = &subfiling_config;
615 }
616 if ( h5_config.subfiling_stripe_count !=
617 subfiling_config.shared_cfg.stripe_count )
618 {
619 subfiling_config.shared_cfg.stripe_count =
620 h5_config.subfiling_stripe_count;
621 if ( subfiling_ptr == NULL )
622 subfiling_ptr = &subfiling_config;
623 }
624 if ( h5_config.subfiling_ioc_selection !=
625 (int)subfiling_config.shared_cfg.ioc_selection )
626 {
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;
631 }
632 if ( h5_config.subfiling_thread_pool_size !=
633 H5FD_IOC_DEFAULT_THREAD_POOL_SIZE )
634 {
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;
639 }
640 H5Pclose( fapl_id );
641
642 H5Pset_mpi_params( plist_id, comm, MPI_INFO_NULL );
643
644 if ( ioc_ptr != NULL )
645 H5Pset_fapl_ioc( subfiling_config.ioc_fapl_id, ioc_ptr );
646
647 H5Pset_fapl_subfiling( plist_id, subfiling_ptr );
648 }
649 else
650#endif
651 {
652 H5Pset_fapl_mpio( plist_id, comm, MPI_INFO_NULL );
653 }
654
655 file_id = H5Fcreate( filename_hdf5.str().c_str(), H5F_ACC_TRUNC,
656 H5P_DEFAULT, plist_id );
657 H5Pclose( plist_id );
658
659 // Write current simulation time
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 );
664 H5Aclose( attr_id );
665 H5Sclose( fspace );
666
667 // Reorder the coordinates in a blocked format.
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 );
679 } );
680
681 // Mirror the coordinates to the host.
682 auto host_coords =
683 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), coords_view );
684
685 std::vector<int> all_offsets( comm_size );
686 all_offsets[comm_rank] = n_local;
687
688 MPI_Allreduce( MPI_IN_PLACE, all_offsets.data(), comm_size, MPI_INT,
689 MPI_SUM, comm );
690
691 offset[0] = 0;
692 offset[1] = 0;
693
694 size_t n_global = 0;
695 for ( int i = 0; i < comm_size; i++ )
696 {
697 if ( i < comm_rank )
698 {
699 offset[0] += static_cast<hsize_t>( all_offsets[i] );
700 }
701 n_global += (size_t)all_offsets[i];
702 }
703 std::vector<int>().swap( all_offsets );
704
705 dimsf[0] = n_global;
706 dimsf[1] = 3;
707
708 filespace_id = H5Screate_simple( 2, dimsf, NULL );
709
710 count[0] = n_local;
711 count[1] = 3;
712
713 memspace_id = H5Screate_simple( 2, count, NULL );
714
715 plist_id = H5Pcreate( H5P_DATASET_XFER );
716
717 // Default IO in HDF5 is independent
718 if ( h5_config.collective )
719 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
720
721 std::string dtype;
722 uint precision;
723 hid_t type_id = HDF5Traits<typename CoordSliceType::value_type>::type(
724 &dtype, &precision );
725
726 dcpl_id = H5Pcreate( H5P_DATASET_CREATE );
727 H5Pset_fill_time( dcpl_id, H5D_FILL_TIME_NEVER );
728
729 dset_id = H5Dcreate( file_id, coords_slice.label().c_str(), type_id,
730 filespace_id, H5P_DEFAULT, dcpl_id, H5P_DEFAULT );
731
732 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, NULL, count,
733 NULL );
734
735 H5Dwrite( dset_id, type_id, memspace_id, filespace_id, plist_id,
736 host_coords.data() );
737 H5Dclose( dset_id );
738
739 H5Pclose( plist_id );
740 H5Pclose( dcpl_id );
741 H5Sclose( filespace_id );
742 H5Sclose( memspace_id );
743
744 if ( 0 == comm_rank )
745 {
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() );
750 }
751
752 // Add variables.
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(),
756 fields... );
757
758 H5Fclose( file_id );
759
760 if ( 0 == comm_rank )
761 Impl::writeXdmfFooter( filename_xdmf.str().c_str() );
762}
763
764//---------------------------------------------------------------------------//
765// HDF5 (XDMF) Particle Field Input.
766//---------------------------------------------------------------------------//
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 )
774{
775 // Read the field into a View.
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,
779 host_view.data() );
780
781 // Mirror the field and copy.
782 auto view = Kokkos::create_mirror_view_and_copy(
783 typename SliceType::memory_space(), host_view );
784 copyViewToSlice( slice, view, 0, n_local );
785}
786
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 )
794{
795 // Read the field into a View.
796 Kokkos::View<typename SliceType::value_type**, Kokkos::LayoutRight,
797 Kokkos::HostSpace>
798 host_view( Kokkos::ViewAllocateWithoutInitializing( "field" ), n_local,
799 slice.extent( 2 ) );
800 H5Dread( dset_id, dtype_id, memspace_id, filespace_id, plist_id,
801 host_view.data() );
802
803 // Mirror the field and copy.
804 auto view = Kokkos::create_mirror_view_and_copy(
805 typename SliceType::memory_space(), host_view );
806 copyViewToSlice( slice, view, 0, n_local );
807}
808
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 )
816{
817 // Read the field into a View.
818 Kokkos::View<typename SliceType::value_type***, Kokkos::LayoutRight,
819 Kokkos::HostSpace>
820 host_view( Kokkos::ViewAllocateWithoutInitializing( "field" ), n_local,
821 slice.extent( 2 ), slice.extent( 3 ) );
822 H5Dread( dset_id, dtype_id, memspace_id, filespace_id, plist_id,
823 host_view.data() );
824
825 // Mirror the field and copy.
826 auto view = Kokkos::create_mirror_view_and_copy(
827 typename SliceType::memory_space(), host_view );
828 copyViewToSlice( slice, view, 0, n_local );
829}
830
831//---------------------------------------------------------------------------//
843template <class FieldSliceType>
844void readTimeStep( HDF5Config h5_config, const std::string& prefix,
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 )
848{
849 Kokkos::Profiling::ScopedRegion region( "Cabana::HDF5ParticleInput" );
850
851 hid_t plist_id;
852 hid_t dset_id;
853 hid_t file_id;
854 hid_t dtype_id;
855 hid_t filespace_id;
856 hid_t memspace_id;
857 int ndims;
858
859 // HDF5 hyperslab parameters
860 hsize_t offset[3] = { 0, 0, 0 };
861 hsize_t dimsf[3] = { 0, 0, 0 };
862 hsize_t count[3] = { 0, 0, 0 };
863
864 int comm_rank;
865 MPI_Comm_rank( comm, &comm_rank );
866 int comm_size;
867 MPI_Comm_size( comm, &comm_size );
868
869 // Retrieve data file name.
870 std::stringstream filename_hdf5;
871 filename_hdf5 << prefix << "_" << time_step_index << ".h5";
872
873 plist_id = H5Pcreate( H5P_FILE_ACCESS );
874 H5Pset_fapl_mpio( plist_id, comm, MPI_INFO_NULL );
875
876#if H5_VERSION_GE( 1, 10, 0 )
877 if ( h5_config.collective )
878 {
879 H5Pset_all_coll_metadata_ops( plist_id, 1 );
880 }
881#endif
882
883 // Open the HDF5 file.
884 file_id = H5Fopen( filename_hdf5.str().c_str(), H5F_ACC_RDONLY, plist_id );
885 H5Pclose( plist_id );
886
887 // Get current simulation time associated with time_step_index
888 hid_t attr_id = H5Aopen( file_id, "Time", H5P_DEFAULT );
889 H5Aread( attr_id, H5T_NATIVE_DOUBLE, &time );
890 H5Aclose( attr_id );
891
892 // Open the dataset.
893 dset_id = H5Dopen( file_id, dataset_name.c_str(), H5P_DEFAULT );
894
895 // Get the datatype of the dataset.
896 dtype_id = H5Dget_type( dset_id );
897
898 // Get the dataspace of the dataset.
899 filespace_id = H5Dget_space( dset_id );
900
901 // Get the rank of the dataspace.
902 ndims = H5Sget_simple_extent_ndims( filespace_id );
903
904 // Get the extents of the file dataspace.
905 H5Sget_simple_extent_dims( filespace_id, dimsf, NULL );
906
907 std::vector<int> all_offsets( comm_size );
908 all_offsets[comm_rank] = n_local;
909
910 MPI_Allreduce( MPI_IN_PLACE, all_offsets.data(), comm_size, MPI_INT,
911 MPI_SUM, comm );
912
913 for ( int i = 0; i < comm_size; i++ )
914 {
915 if ( i < comm_rank )
916 {
917 offset[0] += static_cast<hsize_t>( all_offsets[i] );
918 }
919 }
920 std::vector<int>().swap( all_offsets );
921
922 count[0] = n_local;
923 count[1] = dimsf[1];
924 count[2] = dimsf[2];
925
926 memspace_id = H5Screate_simple( ndims, count, NULL );
927
928 plist_id = H5Pcreate( H5P_DATASET_XFER );
929
930 // Default IO in HDF5 is independent
931 if ( h5_config.collective )
932 H5Pset_dxpl_mpio( plist_id, H5FD_MPIO_COLLECTIVE );
933
934 H5Sselect_hyperslab( filespace_id, H5S_SELECT_SET, offset, NULL, count,
935 NULL );
936
937 readField( dset_id, dtype_id, memspace_id, filespace_id, plist_id, n_local,
938 field );
939
940 H5Pclose( plist_id );
941 H5Sclose( memspace_id );
942 H5Sclose( filespace_id );
943 H5Dclose( dset_id );
944 H5Fclose( file_id );
945}
946
947//---------------------------------------------------------------------------//
948
949} // namespace HDF5ParticleOutput
950} // namespace Experimental
951} // end namespace Cabana
952
953#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: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