libROM  v1.0
Data-driven physical simulation library
HDFDatabaseMPIO.cpp
1 
11 // Description: The concrete database implementation using HDF5.
12 
13 #include "HDFDatabaseMPIO.h"
14 #include "utils/mpi_utils.h"
15 #include <iostream>
16 
17 namespace CAROM {
18 
20 #if HDF5_IS_PARALLEL
21  HDFDatabase(),
22  d_rank(-1),
23  d_comm(MPI_COMM_NULL)
24 #else
25  HDFDatabase()
26 #endif
27 {}
28 
29 #if HDF5_IS_PARALLEL
30 
31 bool
33  const std::string& file_name,
34  const MPI_Comm comm)
35 {
36  int mpi_init;
37  MPI_Initialized(&mpi_init);
38  if (mpi_init) {
39  CAROM_VERIFY(comm != MPI_COMM_NULL);
40  d_comm = comm;
41  MPI_Comm_rank(d_comm, &d_rank);
42  }
43  else
44  d_rank = 0;
45 
46  /* Create a single file with a name compatible to HDFDatabase. */
47  std::string file_name_ext(file_name + ".000000");
48  if (d_rank == 0)
49  Database::create(file_name_ext);
50  CAROM_VERIFY(!file_name.empty());
51 
52  hid_t plist_id;
53  /*
54  * Set up file access property list with parallel I/O access
55  */
56  plist_id = H5Pcreate(H5P_FILE_ACCESS);
57  H5Pset_fapl_mpio(plist_id, d_comm, MPI_INFO_NULL);
58  /*
59  * OPTIONAL: It is generally recommended to set collective
60  * metadata reads/writes on FAPL to perform metadata reads
61  * collectively, which usually allows datasets
62  * to perform better at scale, although it is not
63  * strictly necessary.
64  */
65  H5Pset_all_coll_metadata_ops(plist_id, true);
66  H5Pset_coll_metadata_write(plist_id, true);
67 
68  hid_t file_id = H5Fcreate(file_name_ext.c_str(),
69  H5F_ACC_TRUNC,
70  H5P_DEFAULT,
71  plist_id);
72  bool result = file_id >= 0;
73  CAROM_VERIFY(result);
74  d_is_file = true;
75  d_file_id = file_id;
76  d_group_id = file_id;
77 
78  H5Pclose(plist_id);
79 
80  return result;
81 }
82 
83 bool
85  const std::string& file_name,
86  const std::string& type,
87  const MPI_Comm comm)
88 {
89  int mpi_init;
90  MPI_Initialized(&mpi_init);
91  if (mpi_init) {
92  CAROM_VERIFY(comm != MPI_COMM_NULL);
93  d_comm = comm;
94  MPI_Comm_rank(d_comm, &d_rank);
95  }
96  else
97  d_rank = 0;
98 
99  std::string file_name_ext(file_name + ".000000");
100  if (d_rank == 0)
101  Database::open(file_name_ext, type);
102  CAROM_VERIFY(!file_name.empty());
103  CAROM_VERIFY(type == "r" || type == "wr");
104 
105  hid_t plist_id;
106  /*
107  * Set up file access property list with parallel I/O access
108  */
109  plist_id = H5Pcreate(H5P_FILE_ACCESS);
110  H5Pset_fapl_mpio(plist_id, d_comm, MPI_INFO_NULL);
111  /*
112  * OPTIONAL: It is generally recommended to set collective
113  * metadata reads/writes on FAPL to perform metadata reads
114  * collectively, which usually allows datasets
115  * to perform better at scale, although it is not
116  * strictly necessary.
117  */
118  H5Pset_all_coll_metadata_ops(plist_id, true);
119  H5Pset_coll_metadata_write(plist_id, true);
120 
121  hid_t file_id;
122  if (type == "r")
123  {
124  file_id = H5Fopen(file_name_ext.c_str(),
125  H5F_ACC_RDONLY,
126  plist_id);
127  }
128  else if (type == "wr")
129  {
130  file_id = H5Fopen(file_name_ext.c_str(),
131  H5F_ACC_RDWR,
132  plist_id);
133  }
134  bool result = file_id >= 0;
135  CAROM_VERIFY(result);
136  d_is_file = true;
137  d_file_id = file_id;
138  d_group_id = file_id;
139 
140  H5Pclose(plist_id);
141 
142  return result;
143 }
144 
145 void
146 HDFDatabaseMPIO::putIntegerArray_parallel(
147  const std::string& key,
148  const int* const data,
149  int nelem_local)
150 {
151  CAROM_VERIFY(!key.empty());
152  CAROM_VERIFY(data != nullptr);
153  CAROM_VERIFY(nelem_local >= 0);
154 
155  /* determine global nelements and offsets */
156  std::vector<int> offsets;
157  int nelements = CAROM::get_global_offsets(nelem_local, offsets, d_comm);
158  CAROM_VERIFY(nelements > 0);
159 
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);
164 
165 #if (H5_VERS_MAJOR > 1) || ((H5_VERS_MAJOR == 1) && (H5_VERS_MINOR > 6))
166  hid_t dataset = H5Dcreate(d_group_id,
167  key.c_str(),
168  H5T_NATIVE_INT,
169  filespace,
170  H5P_DEFAULT,
171  H5P_DEFAULT,
172  H5P_DEFAULT);
173 #else
174  CAROM_ERROR("HDFDatabaseMPIO is not compatible with current HDF5 version!\n");
175 #endif
176  CAROM_VERIFY(dataset >= 0);
177  H5Sclose(filespace);
178 
179  /*
180  * Each process defines dataset in memory and writes it to the hyperslab
181  * in the file.
182  */
183  /* hyperslab selection parameters */
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);
189 
190  /*
191  * Select hyperslab in the file.
192  */
193  filespace = H5Dget_space(dataset);
194  if (nelem_local > 0)
195  H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);
196  else
197  H5Sselect_none(filespace);
198 
199  /*
200  * Create property list for collective dataset write.
201  */
202  hid_t plist_id = H5Pcreate(H5P_DATASET_XFER);
203  H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);
204 
205  herr_t errf = H5Dwrite(dataset,
206  H5T_NATIVE_INT,
207  memspace,
208  filespace,
209  plist_id,
210  data);
211  CAROM_VERIFY(errf >= 0);
212 
213  // Write attribute so we know what kind of data this is.
214  writeAttribute(KEY_INT_ARRAY, dataset);
215 
216  errf = H5Sclose(filespace);
217  CAROM_VERIFY(errf >= 0);
218 
219  errf = H5Sclose(memspace);
220  CAROM_VERIFY(errf >= 0);
221 
222  errf = H5Pclose(plist_id);
223  CAROM_VERIFY(errf >= 0);
224 
225  errf = H5Dclose(dataset);
226  CAROM_VERIFY(errf >= 0);
227 #ifndef DEBUG_CHECK_ASSERTIONS
228  CAROM_NULL_USE(errf);
229 #endif
230 }
231 
232 void
233 HDFDatabaseMPIO::putDoubleArray_parallel(
234  const std::string& key,
235  const double* const data,
236  int nelem_local)
237 {
238  CAROM_VERIFY(!key.empty());
239  CAROM_VERIFY(data != nullptr);
240  CAROM_VERIFY(nelem_local >= 0);
241 
242  /* determine global nelements and offsets */
243  std::vector<int> offsets;
244  int nelements = CAROM::get_global_offsets(nelem_local, offsets, d_comm);
245  CAROM_VERIFY(nelements > 0);
246 
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);
251 
252 #if (H5_VERS_MAJOR > 1) || ((H5_VERS_MAJOR == 1) && (H5_VERS_MINOR > 6))
253  hid_t dataset = H5Dcreate(d_group_id,
254  key.c_str(),
255  H5T_NATIVE_DOUBLE,
256  filespace,
257  H5P_DEFAULT,
258  H5P_DEFAULT,
259  H5P_DEFAULT);
260 #else
261  CAROM_ERROR("HDFDatabaseMPIO is not compatible with current HDF5 version!\n");
262 #endif
263  CAROM_VERIFY(dataset >= 0);
264  H5Sclose(filespace);
265 
266  /*
267  * Each process defines dataset in memory and writes it to the hyperslab
268  * in the file.
269  */
270  /* hyperslab selection parameters */
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);
276 
277  /*
278  * Select hyperslab in the file.
279  */
280  filespace = H5Dget_space(dataset);
281  if (nelem_local > 0)
282  H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);
283  else
284  H5Sselect_none(filespace);
285 
286  /*
287  * Create property list for collective dataset write.
288  */
289  hid_t plist_id = H5Pcreate(H5P_DATASET_XFER);
290  H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);
291 
292  herr_t errf = H5Dwrite(dataset,
293  H5T_NATIVE_DOUBLE,
294  memspace,
295  filespace,
296  plist_id,
297  data);
298  CAROM_VERIFY(errf >= 0);
299 
300  // Write attribute so we know what kind of data this is.
302 
303  errf = H5Sclose(filespace);
304  CAROM_VERIFY(errf >= 0);
305 
306  errf = H5Sclose(memspace);
307  CAROM_VERIFY(errf >= 0);
308 
309  errf = H5Pclose(plist_id);
310  CAROM_VERIFY(errf >= 0);
311 
312  errf = H5Dclose(dataset);
313  CAROM_VERIFY(errf >= 0);
314 #ifndef DEBUG_CHECK_ASSERTIONS
315  CAROM_NULL_USE(errf);
316 #endif
317 }
318 
319 void
320 HDFDatabaseMPIO::getIntegerArray_parallel(
321  const std::string& key,
322  int* data,
323  int nelem_local)
324 {
325  CAROM_VERIFY(nelem_local >= 0);
326  /* determine global nelements and offsets */
327  std::vector<int> offsets;
328  int nelements = CAROM::get_global_offsets(nelem_local, offsets, d_comm);
329  if (nelements == 0) return;
330 
331  CAROM_VERIFY(!key.empty());
332 #ifndef DEBUG_CHECK_ASSERTIONS
333  CAROM_NULL_USE(nelem_local);
334 #endif
335 
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);
338 #else
339  CAROM_ERROR("HDFDatabaseMPIO is not compatible with current HDF5 version!\n");
340 #endif
341  CAROM_VERIFY(dset >= 0);
342 
343  hid_t filespace = H5Dget_space(dset);
344  CAROM_VERIFY(filespace >= 0);
345 
346  hsize_t nsel = H5Sget_select_npoints(filespace);
347  CAROM_VERIFY(static_cast<int>(nsel) == nelements);
348 
349  const int dim_rank = 1;
350  H5Sclose(filespace);
351 
352  /*
353  * Each process defines dataset in memory and writes it to the hyperslab
354  * in the file.
355  */
356  /* hyperslab selection parameters */
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);
362 
363  /*
364  * Select hyperslab in the file.
365  */
366  filespace = H5Dget_space(dset);
367  if (nelem_local > 0)
368  H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);
369  else
370  H5Sselect_none(filespace);
371 
372  /*
373  * Create property list for collective dataset write.
374  */
375  hid_t plist_id = H5Pcreate(H5P_DATASET_XFER);
376  H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);
377 
378  herr_t errf;
379  errf = H5Dread(dset, H5T_NATIVE_INT, memspace, filespace, plist_id, data);
380  CAROM_VERIFY(errf >= 0);
381 
382  errf = H5Sclose(filespace);
383  CAROM_VERIFY(errf >= 0);
384 
385  errf = H5Sclose(memspace);
386  CAROM_VERIFY(errf >= 0);
387 
388  errf = H5Pclose(plist_id);
389  CAROM_VERIFY(errf >= 0);
390 
391  errf = H5Dclose(dset);
392  CAROM_VERIFY(errf >= 0);
393 #ifndef DEBUG_CHECK_ASSERTIONS
394  CAROM_NULL_USE(errf);
395 #endif
396 }
397 
398 void
399 HDFDatabaseMPIO::getDoubleArray_parallel(
400  const std::string& key,
401  double* data,
402  int nelem_local)
403 {
404  CAROM_VERIFY(nelem_local >= 0);
405  /* determine global nelements and offsets */
406  std::vector<int> offsets;
407  int nelements = CAROM::get_global_offsets(nelem_local, offsets, d_comm);
408  if (nelements == 0) return;
409 
410  CAROM_VERIFY(!key.empty());
411 #ifndef DEBUG_CHECK_ASSERTIONS
412  CAROM_NULL_USE(nelements);
413 #endif
414 
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);
417 #else
418  hid_t dset = H5Dopen(d_group_id, key.c_str());
419 #endif
420  CAROM_VERIFY(dset >= 0);
421 
422  hid_t filespace = H5Dget_space(dset);
423  CAROM_VERIFY(filespace >= 0);
424 
425  hsize_t nsel = H5Sget_select_npoints(filespace);
426  CAROM_VERIFY(static_cast<int>(nsel) == nelements);
427 
428  const int dim_rank = 1;
429  H5Sclose(filespace);
430 
431  /*
432  * Each process defines dataset in memory and writes it to the hyperslab
433  * in the file.
434  */
435  /* hyperslab selection parameters */
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);
441 
442  /*
443  * Select hyperslab in the file.
444  */
445  filespace = H5Dget_space(dset);
446  if (nelem_local > 0)
447  H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);
448  else
449  H5Sselect_none(filespace);
450 
451  /*
452  * Create property list for collective dataset write.
453  */
454  hid_t plist_id = H5Pcreate(H5P_DATASET_XFER);
455  H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);
456 
457  herr_t errf;
458  errf = H5Dread(dset, H5T_NATIVE_DOUBLE, memspace, filespace, plist_id, data);
459  CAROM_VERIFY(errf >= 0);
460 
461  errf = H5Sclose(filespace);
462  CAROM_VERIFY(errf >= 0);
463 
464  errf = H5Sclose(memspace);
465  CAROM_VERIFY(errf >= 0);
466 
467  errf = H5Pclose(plist_id);
468  CAROM_VERIFY(errf >= 0);
469 
470  errf = H5Dclose(dset);
471  CAROM_VERIFY(errf >= 0);
472 #ifndef DEBUG_CHECK_ASSERTIONS
473  CAROM_NULL_USE(errf);
474 #endif
475 }
476 
477 void
478 HDFDatabaseMPIO::getDoubleArray_parallel(
479  const std::string& key,
480  double* data,
481  int nelem_local,
482  const std::vector<int>& idx_local)
483 {
484  std::vector<int> idx_tmp = idx_local;
485  if (idx_local.size() == 0)
486  {
487  for (int i = 0; i < nelem_local; i++)
488  idx_tmp.push_back(i);
489  }
490 
491  std::vector<double> alldata(nelem_local);
492  getDoubleArray_parallel(key, alldata.data(), nelem_local);
493  int k = 0;
494  for (int i = 0; i < nelem_local; ++i)
495  {
496  if (idx_tmp[k] == i)
497  {
498  data[k++] = alldata[i];
499  }
500  if (k == idx_tmp.size())
501  {
502  break;
503  }
504  }
505  CAROM_VERIFY(k == idx_tmp.size());
506 }
507 
508 void
509 HDFDatabaseMPIO::getDoubleArray_parallel(
510  const std::string& key,
511  double* data,
512  int nelem_local,
513  int block_offset_global,
514  int block_size_global,
515  int stride_global)
516 {
517  CAROM_VERIFY(nelem_local >= 0);
518  CAROM_VERIFY(CAROM::is_same(stride_global, d_comm));
519  CAROM_VERIFY(CAROM::is_same(block_size_global, d_comm));
520  CAROM_VERIFY(CAROM::is_same(block_offset_global, d_comm));
521  /* We assume stride sets the maximum block size */
522  CAROM_VERIFY(block_offset_global + block_size_global <= stride_global);
523  /* determine global nelements and offsets */
524  hsize_t num_local_blocks = nelem_local / block_size_global;
525  std::vector<int> global_offsets;
526  int nelements = CAROM::get_global_offsets(num_local_blocks * stride_global,
527  global_offsets, d_comm);
528 
529  CAROM_VERIFY(!key.empty());
530 #ifndef DEBUG_CHECK_ASSERTIONS
531  CAROM_NULL_USE(nelements);
532 #endif
533 
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);
536 #else
537  hid_t dset = H5Dopen(d_group_id, key.c_str());
538 #endif
539  CAROM_VERIFY(dset >= 0);
540 
541  hid_t filespace = H5Dget_space(dset);
542  CAROM_VERIFY(filespace >= 0);
543 
544  hsize_t nsel = H5Sget_select_npoints(filespace);
545  CAROM_VERIFY((nsel == 0) || (nsel == nelements));
546 
547  const int dim_rank = 1;
548  H5Sclose(filespace);
549 
550  herr_t errf;
551  if (nsel > 0) {
552  /*
553  * Each process defines dataset in memory and writes it to the hyperslab
554  * in the file.
555  */
556  /* hyperslab selection parameters */
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);
563 
564  /*
565  * Select hyperslab in the file.
566  */
567  filespace = H5Dget_space(dset);
568  if (nelem_local > 0)
569  H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset,
570  strides, num_blocks, block_sizes);
571  else
572  H5Sselect_none(filespace);
573 
574  /*
575  * Create property list for collective dataset write.
576  */
577  hid_t plist_id = H5Pcreate(H5P_DATASET_XFER);
578  H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);
579 
580  errf = H5Dread(dset, H5T_NATIVE_DOUBLE, memspace, filespace, plist_id, data);
581  CAROM_VERIFY(errf >= 0);
582 
583  errf = H5Sclose(memspace);
584  CAROM_VERIFY(errf >= 0);
585 
586  errf = H5Pclose(plist_id);
587  CAROM_VERIFY(errf >= 0);
588  }
589 
590  errf = H5Sclose(filespace);
591  CAROM_VERIFY(errf >= 0);
592 
593  errf = H5Dclose(dset);
594  CAROM_VERIFY(errf >= 0);
595 #ifndef DEBUG_CHECK_ASSERTIONS
596  CAROM_NULL_USE(errf);
597 #endif
598 }
599 
600 void
602  int type_key,
603  hid_t dataset_id)
604 {
605  CAROM_VERIFY(CAROM::is_same(type_key, d_comm));
606 
607  hid_t attr_id = H5Screate(H5S_SCALAR);
608  CAROM_VERIFY(attr_id >= 0);
609 
610 #if (H5_VERS_MAJOR > 1) || ((H5_VERS_MAJOR == 1) && (H5_VERS_MINOR > 6))
611  hid_t attr = H5Acreate(dataset_id,
612  "Type",
613  H5T_NATIVE_INT,
614  attr_id,
615  H5P_DEFAULT,
616  H5P_DEFAULT);
617 #else
618  CAROM_ERROR("HDFDatabaseMPIO is not compatible with current HDF5 version!\n");
619 #endif
620  CAROM_VERIFY(attr >= 0);
621 
622  /* only the last process writes the attribute */
623  herr_t errf = H5Awrite(attr, H5T_NATIVE_INT, &type_key);
624  CAROM_VERIFY(errf >= 0);
625 
626  errf = H5Aclose(attr);
627  CAROM_VERIFY(errf >= 0);
628 
629  errf = H5Sclose(attr_id);
630  CAROM_VERIFY(errf >= 0);
631 #ifndef DEBUG_CHECK_ASSERTIONS
632  CAROM_NULL_USE(errf);
633 #endif
634 }
635 
636 #endif
637 
638 }
virtual bool create(const std::string &file_name, const MPI_Comm comm=MPI_COMM_NULL)
Creates a new database file with the supplied name.
Definition: Database.cpp:29
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.
Definition: Database.cpp:38
bool d_is_file
True if the HDF5 database is mounted to a file.
Definition: HDFDatabase.h:323
static const int KEY_DOUBLE_ARRAY
The key representing a double array.
Definition: HDFDatabase.h:338
int d_rank
MPI process rank.
Definition: HDFDatabase.h:348
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.
Definition: HDFDatabase.cpp:45
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.
Definition: HDFDatabase.cpp:76
hid_t d_group_id
ID of the group attached to the database.
Definition: HDFDatabase.h:333
static const int KEY_INT_ARRAY
The key representing an integer array.
Definition: HDFDatabase.h:343
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.
Definition: HDFDatabase.h:328
HDFDatabaseMPIO()
Default constructor.
bool is_same(int x, const MPI_Comm &comm)
Check if an integer is equal over all MPI processes.
Definition: mpi_utils.cpp:72
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...
Definition: mpi_utils.cpp:41