13 #include "HDFDatabaseMPIO.h"
14 #include "utils/mpi_utils.h"
33 const std::string& file_name,
37 MPI_Initialized(&mpi_init);
39 CAROM_VERIFY(comm != MPI_COMM_NULL);
41 MPI_Comm_rank(d_comm, &
d_rank);
47 std::string file_name_ext(file_name +
".000000");
50 CAROM_VERIFY(!file_name.empty());
56 plist_id = H5Pcreate(H5P_FILE_ACCESS);
57 H5Pset_fapl_mpio(plist_id, d_comm, MPI_INFO_NULL);
65 H5Pset_all_coll_metadata_ops(plist_id,
true);
66 H5Pset_coll_metadata_write(plist_id,
true);
68 hid_t file_id = H5Fcreate(file_name_ext.c_str(),
72 bool result = file_id >= 0;
85 const std::string& file_name,
86 const std::string& type,
90 MPI_Initialized(&mpi_init);
92 CAROM_VERIFY(comm != MPI_COMM_NULL);
94 MPI_Comm_rank(d_comm, &
d_rank);
99 std::string file_name_ext(file_name +
".000000");
102 CAROM_VERIFY(!file_name.empty());
103 CAROM_VERIFY(type ==
"r" || type ==
"wr");
109 plist_id = H5Pcreate(H5P_FILE_ACCESS);
110 H5Pset_fapl_mpio(plist_id, d_comm, MPI_INFO_NULL);
118 H5Pset_all_coll_metadata_ops(plist_id,
true);
119 H5Pset_coll_metadata_write(plist_id,
true);
124 file_id = H5Fopen(file_name_ext.c_str(),
128 else if (type ==
"wr")
130 file_id = H5Fopen(file_name_ext.c_str(),
134 bool result = file_id >= 0;
135 CAROM_VERIFY(result);
146 HDFDatabaseMPIO::putIntegerArray_parallel(
147 const std::string& key,
148 const int*
const data,
151 CAROM_VERIFY(!key.empty());
152 CAROM_VERIFY(data !=
nullptr);
153 CAROM_VERIFY(nelem_local >= 0);
156 std::vector<int> offsets;
158 CAROM_VERIFY(nelements > 0);
160 const int dim_rank = 1;
161 hsize_t dim[dim_rank] = {
static_cast<hsize_t
>(nelements) };
162 hid_t filespace = H5Screate_simple(dim_rank, dim, NULL);
163 CAROM_VERIFY(filespace >= 0);
165 #if (H5_VERS_MAJOR > 1) || ((H5_VERS_MAJOR == 1) && (H5_VERS_MINOR > 6))
174 CAROM_ERROR(
"HDFDatabaseMPIO is not compatible with current HDF5 version!\n");
176 CAROM_VERIFY(dataset >= 0);
184 hsize_t count[dim_rank];
185 hsize_t offset[dim_rank];
186 count[0] = nelem_local;
187 offset[0] = offsets[
d_rank];
188 hid_t memspace = H5Screate_simple(dim_rank, count, NULL);
193 filespace = H5Dget_space(dataset);
195 H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);
197 H5Sselect_none(filespace);
202 hid_t plist_id = H5Pcreate(H5P_DATASET_XFER);
203 H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);
205 herr_t errf = H5Dwrite(dataset,
211 CAROM_VERIFY(errf >= 0);
216 errf = H5Sclose(filespace);
217 CAROM_VERIFY(errf >= 0);
219 errf = H5Sclose(memspace);
220 CAROM_VERIFY(errf >= 0);
222 errf = H5Pclose(plist_id);
223 CAROM_VERIFY(errf >= 0);
225 errf = H5Dclose(dataset);
226 CAROM_VERIFY(errf >= 0);
227 #ifndef DEBUG_CHECK_ASSERTIONS
228 CAROM_NULL_USE(errf);
233 HDFDatabaseMPIO::putDoubleArray_parallel(
234 const std::string& key,
235 const double*
const data,
238 CAROM_VERIFY(!key.empty());
239 CAROM_VERIFY(data !=
nullptr);
240 CAROM_VERIFY(nelem_local >= 0);
243 std::vector<int> offsets;
245 CAROM_VERIFY(nelements > 0);
247 const int dim_rank = 1;
248 hsize_t dim[dim_rank] = {
static_cast<hsize_t
>(nelements) };
249 hid_t filespace = H5Screate_simple(dim_rank, dim, 0);
250 CAROM_VERIFY(filespace >= 0);
252 #if (H5_VERS_MAJOR > 1) || ((H5_VERS_MAJOR == 1) && (H5_VERS_MINOR > 6))
261 CAROM_ERROR(
"HDFDatabaseMPIO is not compatible with current HDF5 version!\n");
263 CAROM_VERIFY(dataset >= 0);
271 hsize_t count[dim_rank];
272 hsize_t offset[dim_rank];
273 count[0] = nelem_local;
274 offset[0] = offsets[
d_rank];
275 hid_t memspace = H5Screate_simple(dim_rank, count, NULL);
280 filespace = H5Dget_space(dataset);
282 H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);
284 H5Sselect_none(filespace);
289 hid_t plist_id = H5Pcreate(H5P_DATASET_XFER);
290 H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);
292 herr_t errf = H5Dwrite(dataset,
298 CAROM_VERIFY(errf >= 0);
303 errf = H5Sclose(filespace);
304 CAROM_VERIFY(errf >= 0);
306 errf = H5Sclose(memspace);
307 CAROM_VERIFY(errf >= 0);
309 errf = H5Pclose(plist_id);
310 CAROM_VERIFY(errf >= 0);
312 errf = H5Dclose(dataset);
313 CAROM_VERIFY(errf >= 0);
314 #ifndef DEBUG_CHECK_ASSERTIONS
315 CAROM_NULL_USE(errf);
320 HDFDatabaseMPIO::getIntegerArray_parallel(
321 const std::string& key,
325 CAROM_VERIFY(nelem_local >= 0);
327 std::vector<int> offsets;
329 if (nelements == 0)
return;
331 CAROM_VERIFY(!key.empty());
332 #ifndef DEBUG_CHECK_ASSERTIONS
333 CAROM_NULL_USE(nelem_local);
336 #if (H5_VERS_MAJOR > 1) || ((H5_VERS_MAJOR == 1) && (H5_VERS_MINOR > 6))
337 hid_t dset = H5Dopen(
d_group_id, key.c_str(), H5P_DEFAULT);
339 CAROM_ERROR(
"HDFDatabaseMPIO is not compatible with current HDF5 version!\n");
341 CAROM_VERIFY(dset >= 0);
343 hid_t filespace = H5Dget_space(dset);
344 CAROM_VERIFY(filespace >= 0);
346 hsize_t nsel = H5Sget_select_npoints(filespace);
347 CAROM_VERIFY(
static_cast<int>(nsel) == nelements);
349 const int dim_rank = 1;
357 hsize_t count[dim_rank];
358 hsize_t offset[dim_rank];
359 count[0] = nelem_local;
360 offset[0] = offsets[
d_rank];
361 hid_t memspace = H5Screate_simple(dim_rank, count, NULL);
366 filespace = H5Dget_space(dset);
368 H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);
370 H5Sselect_none(filespace);
375 hid_t plist_id = H5Pcreate(H5P_DATASET_XFER);
376 H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);
379 errf = H5Dread(dset, H5T_NATIVE_INT, memspace, filespace, plist_id, data);
380 CAROM_VERIFY(errf >= 0);
382 errf = H5Sclose(filespace);
383 CAROM_VERIFY(errf >= 0);
385 errf = H5Sclose(memspace);
386 CAROM_VERIFY(errf >= 0);
388 errf = H5Pclose(plist_id);
389 CAROM_VERIFY(errf >= 0);
391 errf = H5Dclose(dset);
392 CAROM_VERIFY(errf >= 0);
393 #ifndef DEBUG_CHECK_ASSERTIONS
394 CAROM_NULL_USE(errf);
399 HDFDatabaseMPIO::getDoubleArray_parallel(
400 const std::string& key,
404 CAROM_VERIFY(nelem_local >= 0);
406 std::vector<int> offsets;
408 if (nelements == 0)
return;
410 CAROM_VERIFY(!key.empty());
411 #ifndef DEBUG_CHECK_ASSERTIONS
412 CAROM_NULL_USE(nelements);
415 #if (H5_VERS_MAJOR > 1) || ((H5_VERS_MAJOR == 1) && (H5_VERS_MINOR > 6))
416 hid_t dset = H5Dopen(
d_group_id, key.c_str(), H5P_DEFAULT);
418 hid_t dset = H5Dopen(
d_group_id, key.c_str());
420 CAROM_VERIFY(dset >= 0);
422 hid_t filespace = H5Dget_space(dset);
423 CAROM_VERIFY(filespace >= 0);
425 hsize_t nsel = H5Sget_select_npoints(filespace);
426 CAROM_VERIFY(
static_cast<int>(nsel) == nelements);
428 const int dim_rank = 1;
436 hsize_t count[dim_rank];
437 hsize_t offset[dim_rank];
438 count[0] = nelem_local;
439 offset[0] = offsets[
d_rank];
440 hid_t memspace = H5Screate_simple(dim_rank, count, NULL);
445 filespace = H5Dget_space(dset);
447 H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);
449 H5Sselect_none(filespace);
454 hid_t plist_id = H5Pcreate(H5P_DATASET_XFER);
455 H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);
458 errf = H5Dread(dset, H5T_NATIVE_DOUBLE, memspace, filespace, plist_id, data);
459 CAROM_VERIFY(errf >= 0);
461 errf = H5Sclose(filespace);
462 CAROM_VERIFY(errf >= 0);
464 errf = H5Sclose(memspace);
465 CAROM_VERIFY(errf >= 0);
467 errf = H5Pclose(plist_id);
468 CAROM_VERIFY(errf >= 0);
470 errf = H5Dclose(dset);
471 CAROM_VERIFY(errf >= 0);
472 #ifndef DEBUG_CHECK_ASSERTIONS
473 CAROM_NULL_USE(errf);
478 HDFDatabaseMPIO::getDoubleArray_parallel(
479 const std::string& key,
482 const std::vector<int>& idx_local)
484 std::vector<int> idx_tmp = idx_local;
485 if (idx_local.size() == 0)
487 for (
int i = 0; i < nelem_local; i++)
488 idx_tmp.push_back(i);
491 std::vector<double> alldata(nelem_local);
492 getDoubleArray_parallel(key, alldata.data(), nelem_local);
494 for (
int i = 0; i < nelem_local; ++i)
498 data[k++] = alldata[i];
500 if (k == idx_tmp.size())
505 CAROM_VERIFY(k == idx_tmp.size());
509 HDFDatabaseMPIO::getDoubleArray_parallel(
510 const std::string& key,
513 int block_offset_global,
514 int block_size_global,
517 CAROM_VERIFY(nelem_local >= 0);
522 CAROM_VERIFY(block_offset_global + block_size_global <= stride_global);
524 hsize_t num_local_blocks = nelem_local / block_size_global;
525 std::vector<int> global_offsets;
527 global_offsets, d_comm);
529 CAROM_VERIFY(!key.empty());
530 #ifndef DEBUG_CHECK_ASSERTIONS
531 CAROM_NULL_USE(nelements);
534 #if (H5_VERS_MAJOR > 1) || ((H5_VERS_MAJOR == 1) && (H5_VERS_MINOR > 6))
535 hid_t dset = H5Dopen(
d_group_id, key.c_str(), H5P_DEFAULT);
537 hid_t dset = H5Dopen(
d_group_id, key.c_str());
539 CAROM_VERIFY(dset >= 0);
541 hid_t filespace = H5Dget_space(dset);
542 CAROM_VERIFY(filespace >= 0);
544 hsize_t nsel = H5Sget_select_npoints(filespace);
545 CAROM_VERIFY((nsel == 0) || (nsel == nelements));
547 const int dim_rank = 1;
557 hsize_t num_blocks[dim_rank] = {num_local_blocks};
558 hsize_t buffer_array_size[dim_rank] = {(hsize_t) nelem_local};
559 hsize_t offset[dim_rank] = {(hsize_t) global_offsets[
d_rank] + block_offset_global};
560 hsize_t strides[dim_rank] = {(hsize_t) stride_global};
561 hsize_t block_sizes[1] = {(hsize_t) block_size_global};
562 hid_t memspace = H5Screate_simple(dim_rank, buffer_array_size, NULL);
567 filespace = H5Dget_space(dset);
569 H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset,
570 strides, num_blocks, block_sizes);
572 H5Sselect_none(filespace);
577 hid_t plist_id = H5Pcreate(H5P_DATASET_XFER);
578 H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);
580 errf = H5Dread(dset, H5T_NATIVE_DOUBLE, memspace, filespace, plist_id, data);
581 CAROM_VERIFY(errf >= 0);
583 errf = H5Sclose(memspace);
584 CAROM_VERIFY(errf >= 0);
586 errf = H5Pclose(plist_id);
587 CAROM_VERIFY(errf >= 0);
590 errf = H5Sclose(filespace);
591 CAROM_VERIFY(errf >= 0);
593 errf = H5Dclose(dset);
594 CAROM_VERIFY(errf >= 0);
595 #ifndef DEBUG_CHECK_ASSERTIONS
596 CAROM_NULL_USE(errf);
607 hid_t attr_id = H5Screate(H5S_SCALAR);
608 CAROM_VERIFY(attr_id >= 0);
610 #if (H5_VERS_MAJOR > 1) || ((H5_VERS_MAJOR == 1) && (H5_VERS_MINOR > 6))
611 hid_t attr = H5Acreate(dataset_id,
618 CAROM_ERROR(
"HDFDatabaseMPIO is not compatible with current HDF5 version!\n");
620 CAROM_VERIFY(attr >= 0);
623 herr_t errf = H5Awrite(attr, H5T_NATIVE_INT, &type_key);
624 CAROM_VERIFY(errf >= 0);
626 errf = H5Aclose(attr);
627 CAROM_VERIFY(errf >= 0);
629 errf = H5Sclose(attr_id);
630 CAROM_VERIFY(errf >= 0);
631 #ifndef DEBUG_CHECK_ASSERTIONS
632 CAROM_NULL_USE(errf);
virtual bool create(const std::string &file_name, const MPI_Comm comm=MPI_COMM_NULL)
Creates a new database file with the supplied name.
virtual bool open(const std::string &file_name, const std::string &type, const MPI_Comm comm=MPI_COMM_NULL)
Opens an existing database file with the supplied name.
bool d_is_file
True if the HDF5 database is mounted to a file.
static const int KEY_DOUBLE_ARRAY
The key representing a double array.
int d_rank
MPI process rank.
virtual bool create(const std::string &file_name, const MPI_Comm comm=MPI_COMM_NULL) override
Creates a new HDF5 database file with the supplied name.
virtual bool open(const std::string &file_name, const std::string &type, const MPI_Comm comm=MPI_COMM_NULL) override
Opens an existing HDF5 database file with the supplied name.
hid_t d_group_id
ID of the group attached to the database.
static const int KEY_INT_ARRAY
The key representing an integer array.
virtual void writeAttribute(int type_key, hid_t dataset_id)
Write an attribute to the specified dataset.
hid_t d_file_id
ID of file attached to database or -1 if not mounted to a file.
HDFDatabaseMPIO()
Default constructor.
bool is_same(int x, const MPI_Comm &comm)
Check if an integer is equal over all MPI processes.
int get_global_offsets(const int local_dim, std::vector< int > &offsets, const MPI_Comm &comm)
Save integer offsets for each MPI rank under MPI communicator comm, where each rank as the size of lo...