Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
hdf5.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2018 - 2019 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/config.h>
17 
18 #ifdef DEAL_II_WITH_HDF5
19 
20 # include <deal.II/base/array_view.h>
21 # include <deal.II/base/hdf5.h>
22 
23 # include <deal.II/lac/full_matrix.h>
24 # include <deal.II/lac/vector.h>
25 
26 # include <hdf5.h>
27 
28 # include <memory>
29 # include <numeric>
30 # include <vector>
31 
33 
34 namespace HDF5
35 {
36  namespace internal
37  {
38  // This function gives the HDF5 datatype corresponding to the C++ type. In
39  // the case of std::complex types the HDF5 handlers are automatically freed
40  // using the destructor of std::shared_ptr.
41  // std::shared_ptr is used instead of std::unique_ptr because the destructor
42  // of std::shared_ptr doesn't have to be defined in the template argument.
43  // In the other hand, the destructor of std::unique has to be defined in the
44  // template argument. Native types such as H5T_NATIVE_DOUBLE do not
45  // require a destructor, but compound types such as std::complex<double>
46  // require a destructor to free the HDF5 resources.
47  template <typename number>
48  std::shared_ptr<hid_t>
50  {
51  std::shared_ptr<hid_t> t_type;
53  {
54  return std::make_shared<hid_t>(H5T_NATIVE_FLOAT);
55  }
57  {
58  return std::make_shared<hid_t>(H5T_NATIVE_DOUBLE);
59  }
61  {
62  return std::make_shared<hid_t>(H5T_NATIVE_INT);
63  }
65  {
66  return std::make_shared<hid_t>(H5T_NATIVE_UINT);
67  }
68  else if (std::is_same<number, std::complex<float>>::value)
69  {
70  t_type = std::shared_ptr<hid_t>(new hid_t, [](hid_t *pointer) {
71  // Release the HDF5 resource
72  const herr_t ret = H5Tclose(*pointer);
73  AssertNothrow(ret >= 0, ExcInternalError());
74  (void)ret;
75  delete pointer;
76  });
77  *t_type = H5Tcreate(H5T_COMPOUND, sizeof(std::complex<float>));
78  // The C++ standards committee agreed to mandate that the storage
79  // format used for the std::complex type be binary-compatible with
80  // the C99 type, i.e. an array T[2] with consecutive real [0] and
81  // imaginary [1] parts.
82  herr_t ret = H5Tinsert(*t_type, "r", 0, H5T_NATIVE_FLOAT);
83  Assert(ret >= 0, ExcInternalError());
84  ret = H5Tinsert(*t_type, "i", sizeof(float), H5T_NATIVE_FLOAT);
85  Assert(ret >= 0, ExcInternalError());
86  (void)ret;
87  return t_type;
88  }
89  else if (std::is_same<number, std::complex<double>>::value)
90  {
91  t_type = std::shared_ptr<hid_t>(new hid_t, [](hid_t *pointer) {
92  // Release the HDF5 resource
93  const herr_t ret = H5Tclose(*pointer);
94  AssertNothrow(ret >= 0, ExcInternalError());
95  (void)ret;
96  delete pointer;
97  });
98  *t_type = H5Tcreate(H5T_COMPOUND, sizeof(std::complex<double>));
99  // The C++ standards committee agreed to mandate that the storage
100  // format used for the std::complex type be binary-compatible with
101  // the C99 type, i.e. an array T[2] with consecutive real [0] and
102  // imaginary [1] parts.
103  herr_t ret = H5Tinsert(*t_type, "r", 0, H5T_NATIVE_DOUBLE);
104  Assert(ret >= 0, ExcInternalError());
105  ret = H5Tinsert(*t_type, "i", sizeof(double), H5T_NATIVE_DOUBLE);
106  Assert(ret >= 0, ExcInternalError());
107  (void)ret;
108  return t_type;
109  }
110 
111  // The function should not reach this point
112  Assert(false, ExcInternalError());
113  return t_type;
114  }
115 
116 
117 
118  // Several HDF5 functions such as H5Screate_simple() require a
119  // one-dimensional array that specifies the size of each dimension of the
120  // container, see:
121  // https://support.hdfgroup.org/HDF5/doc1.8/RM/RM_H5S.html#Dataspace-CreateSimple
122  // The function get_container_dimensions returns a vector with the
123  // dimensions of the container.
124  // For a std::vector the function returns std::vector<hsize_t>{vector_size}
125  // For a Vector the function returns std::vector<hsize_t>{vector_size}
126  // For a FullMatrix the function returns std::vector<hsize_t>{rows, columns}
127  template <typename number>
128  std::vector<hsize_t>
129  get_container_dimensions(const std::vector<number> &data)
130  {
131  std::vector<hsize_t> dimensions = {data.size()};
132  return dimensions;
133  }
134 
135 
136 
137  template <typename number>
138  std::vector<hsize_t>
140  {
141  std::vector<hsize_t> dimensions = {data.size()};
142  return dimensions;
143  }
144 
145 
146 
147  template <typename number>
148  std::vector<hsize_t>
150  {
151  std::vector<hsize_t> dimensions = {data.m(), data.n()};
152  return dimensions;
153  }
154 
155 
156 
157  // This function returns the total size of the container.
158  // For a std::vector the function returns int(vector_size)
159  // For a Vector the function returns int(vector_size)
160  // For a FullMatrix the function returns int(rows*columns)
161  template <typename number>
162  unsigned int
163  get_container_size(const std::vector<number> &data)
164  {
165  return static_cast<unsigned int>(data.size());
166  }
167 
168 
169 
170  template <typename number>
171  unsigned int
173  {
174  return static_cast<unsigned int>(data.size());
175  }
176 
177 
178 
179  template <typename number>
180  unsigned int
182  {
183  return static_cast<unsigned int>(data.m() * data.n());
184  }
185 
186 
187 
188  // This function initializes and returns a container of type std::vector,
189  // Vector or FullMatrix. The function does not set the values of the
190  // elements of the container. The container can store data of a HDF5 dataset
191  // or a HDF5 selection. The dimensions parameter holds the dimensions of the
192  // HDF5 dataset or selection.
193  //
194  // In the case of a std::vector, the size of the vector will be the total
195  // size given by dimensions. For example in the case of a dataset of rank 3,
196  // the dimensions are std::vector<hsize_t>{dim_0,dim_1,dim_2}. The size of
197  // the returned std::vector will be dim_0*dim_1*dim_2
198  //
199  // In the case of a ::Vector, the size of the returned ::Vector
200  // will be as well dim_0*dim_1*dim_2
201  //
202  // A FullMatrix can store only data of HDF5 datasets with rank 2. The size
203  // of the FullMatrix will be FullMatrix(dim_0,dim_2)
204  template <typename Container>
205  typename std::enable_if<
206  std::is_same<Container,
207  std::vector<typename Container::value_type>>::value,
208  Container>::type
209  initialize_container(const std::vector<hsize_t> &dimensions)
210  {
211  return Container(std::accumulate(
212  dimensions.begin(), dimensions.end(), 1, std::multiplies<int>()));
213  }
214 
215 
216 
217  template <typename Container>
218  typename std::enable_if<
219  std::is_same<Container, Vector<typename Container::value_type>>::value,
220  Container>::type
221  initialize_container(const std::vector<hsize_t> &dimensions)
222  {
223  return Container(std::accumulate(
224  dimensions.begin(), dimensions.end(), 1, std::multiplies<int>()));
225  }
226 
227 
228 
229  template <typename Container>
230  typename std::enable_if<
231  std::is_same<Container,
233  Container>::type
234  initialize_container(const std::vector<hsize_t> &dimensions)
235  {
236  // If the rank is higher than 2, then remove single-dimensional entries
237  // from the shape defined by dimensions. This is equivalent to the squeeze
238  // function of python/numpy. For example the following code would convert
239  // the vector {1,3,1,2} to {3,2}
240  std::vector<hsize_t> squeezed_dimensions;
241 
242  if (dimensions.size() > 2)
243  {
244  for (const auto &dimension : dimensions)
245  {
246  if (dimension > 1)
247  squeezed_dimensions.push_back(dimension);
248  }
249  }
250  else
251  {
252  squeezed_dimensions = dimensions;
253  }
254 
255  AssertDimension(squeezed_dimensions.size(), 2);
256  return Container(squeezed_dimensions[0], squeezed_dimensions[1]);
257  }
258 
259 
260  // This helper function sets the property list of the read and write
261  // operations of DataSet. A property list has to be created for the MPI
262  // driver. For the serial driver the default H5P_DEFAULT can be used. In
263  // addition H5Pset_dxpl_mpio is used to set the MPI mode to collective.
264  void
265  set_plist(hid_t &plist, const bool mpi)
266  {
267  if (mpi)
268  {
269 # ifdef DEAL_II_WITH_MPI
270  plist = H5Pcreate(H5P_DATASET_XFER);
271  Assert(plist >= 0, ExcInternalError());
272  const herr_t ret = H5Pset_dxpl_mpio(plist, H5FD_MPIO_COLLECTIVE);
273  (void)ret;
274  Assert(ret >= 0, ExcInternalError());
275 # else
276  AssertThrow(false, ExcNotImplemented());
277 # endif
278  }
279  else
280  {
281  plist = H5P_DEFAULT;
282  }
283 
284  (void)plist;
285  (void)mpi;
286  }
287 
288 
289  // This helper function release the property list handler of the read and
290  // write operations of DataSet. For the serial there is no need to release
291  // the property list handler because H5P_DEFAULT has been used. If
292  // query_io_mode is True then H5Pget_mpio_actual_io_mode and
293  // H5Pget_mpio_no_collective_cause are used to check if the operation has
294  // been collective.
295  void
296  release_plist(hid_t & plist,
297  H5D_mpio_actual_io_mode_t &io_mode,
298  uint32_t & local_no_collective_cause,
299  uint32_t & global_no_collective_cause,
300  const bool mpi,
301  const bool query_io_mode)
302  {
303  if (mpi)
304  {
305 # ifdef DEAL_II_WITH_MPI
306  herr_t ret;
307  (void)ret;
308  if (query_io_mode)
309  {
310  ret = H5Pget_mpio_actual_io_mode(plist, &io_mode);
311  Assert(ret >= 0, ExcInternalError());
312  ret =
313  H5Pget_mpio_no_collective_cause(plist,
314  &local_no_collective_cause,
315  &global_no_collective_cause);
316  Assert(ret >= 0, ExcInternalError());
317  }
318  ret = H5Pclose(plist);
319  Assert(ret >= 0, ExcInternalError());
320 # else
321  AssertThrow(false, ExcNotImplemented());
322 # endif
323  }
324 
325  (void)plist;
326  (void)io_mode;
327  (void)local_no_collective_cause;
328  (void)global_no_collective_cause;
329  (void)mpi;
330  (void)query_io_mode;
331  }
332 
333 
334  // Convert a HDF5 no_collective_cause code to a human readable string
335  std::string
336  no_collective_cause_to_string(const uint32_t no_collective_cause)
337  {
338  std::string message;
339 
340  auto append_to_message = [&message](const char *p) {
341  if (message.length() > 0)
342  message += ", ";
343  message += p;
344  };
345 
346  // The first is not a bitmask comparison, the rest are bitmask
347  // comparisons.
348  // https://support.hdfgroup.org/HDF5/doc/RM/RM_H5P.html#Property-GetMpioNoCollectiveCause
349  // See H5Ppublic.h
350  // Hex codes are used because the HDF5 Group can deprecate some of the
351  // enum codes. For example the enum code H5D_MPIO_FILTERS is not defined
352  // in 1.10.2 because it is possible to use compressed datasets with the
353  // MPI/IO driver.
354 
355  // H5D_MPIO_COLLECTIVE
356  if (no_collective_cause == 0x00)
357  {
358  append_to_message("H5D_MPIO_COLLECTIVE");
359  }
360  // H5D_MPIO_SET_INDEPENDENT
361  if ((no_collective_cause & 0x01) == 0x01)
362  {
363  append_to_message("H5D_MPIO_SET_INDEPENDENT");
364  }
365  // H5D_MPIO_DATATYPE_CONVERSION
366  if ((no_collective_cause & 0x02) == 0x02)
367  {
368  append_to_message("H5D_MPIO_DATATYPE_CONVERSION");
369  }
370  // H5D_MPIO_DATA_TRANSFORMS
371  if ((no_collective_cause & 0x04) == 0x04)
372  {
373  append_to_message("H5D_MPIO_DATA_TRANSFORMS");
374  }
375  // H5D_MPIO_NOT_SIMPLE_OR_SCALAR_DATASPACES
376  if ((no_collective_cause & 0x10) == 0x10)
377  {
378  append_to_message("H5D_MPIO_NOT_SIMPLE_OR_SCALAR_DATASPACES");
379  }
380  // H5D_MPIO_NOT_CONTIGUOUS_OR_CHUNKED_DATASET
381  if ((no_collective_cause & 0x20) == 0x20)
382  {
383  append_to_message("H5D_MPIO_NOT_CONTIGUOUS_OR_CHUNKED_DATASET");
384  }
385  // H5D_MPIO_FILTERS
386  if ((no_collective_cause & 0x40) == 0x40)
387  {
388  append_to_message("H5D_MPIO_FILTERS");
389  }
390  return message;
391  }
392  } // namespace internal
393 
394 
395 
396  HDF5Object::HDF5Object(const std::string &name, const bool mpi)
397  : name(name)
398  , mpi(mpi)
399  {}
400 
401 
402 
403  template <typename T>
404  T
405  HDF5Object::get_attribute(const std::string &attr_name) const
406  {
407  const std::shared_ptr<hid_t> t_type = internal::get_hdf5_datatype<T>();
408  T value;
409  hid_t attr;
410  herr_t ret;
411 
412  attr = H5Aopen(*hdf5_reference, attr_name.data(), H5P_DEFAULT);
413  Assert(attr >= 0, ExcMessage("Error at H5Aopen"));
414  (void)ret;
415  ret = H5Aread(attr, *t_type, &value);
416  Assert(ret >= 0, ExcMessage("Error at H5Aread"));
417  (void)ret;
418  ret = H5Aclose(attr);
419  Assert(ret >= 0, ExcMessage("Error at H5Aclose"));
420 
421  return value;
422  }
423 
424 
425 
426  template <>
427  bool
428  HDF5Object::get_attribute(const std::string &attr_name) const
429  {
430  // The enum field generated by h5py can be casted to int
431  int int_value;
432  hid_t attr;
433  herr_t ret;
434 
435  attr = H5Aopen(*hdf5_reference, attr_name.data(), H5P_DEFAULT);
436  Assert(attr >= 0, ExcMessage("Error at H5Aopen"));
437  ret = H5Aread(attr, H5T_NATIVE_INT, &int_value);
438  (void)ret;
439  Assert(ret >= 0, ExcMessage("Error at H5Aread"));
440  ret = H5Aclose(attr);
441  (void)ret;
442  Assert(ret >= 0, ExcMessage("Error at H5Aclose"));
443 
444  return (int_value != 0);
445  }
446 
447 
448 
449  template <>
450  std::string
451  HDF5Object::get_attribute(const std::string &attr_name) const
452  {
453  // Reads a UTF8 variable string
454  //
455  // code inspired from
456  // https://support.hdfgroup.org/ftp/HDF5/examples/misc-examples/vlstratt.c
457  //
458  // In the case of a variable length string the user does not have to reserve
459  // memory for string_out. H5Aread will reserve the memory and the
460  // user has to free the memory.
461  //
462  // Todo:
463  // - Use H5Dvlen_reclaim instead of free
464 
465  char * string_out;
466  hid_t attr;
467  hid_t type;
468  herr_t ret;
469 
470  /* Create a datatype to refer to. */
471  type = H5Tcopy(H5T_C_S1);
472  Assert(type >= 0, ExcInternalError());
473 
474  // Python strings are encoded in UTF8
475  ret = H5Tset_cset(type, H5T_CSET_UTF8);
476  Assert(type >= 0, ExcInternalError());
477 
478  ret = H5Tset_size(type, H5T_VARIABLE);
479  Assert(ret >= 0, ExcInternalError());
480 
481  attr = H5Aopen(*hdf5_reference, attr_name.data(), H5P_DEFAULT);
482  Assert(attr >= 0, ExcInternalError());
483 
484  ret = H5Aread(attr, type, &string_out);
485  Assert(ret >= 0, ExcInternalError());
486 
487  std::string string_value(string_out);
488  // The memory of the variable length string has to be freed.
489  // H5Dvlen_reclaim could be also used
490  free(string_out);
491  ret = H5Tclose(type);
492  Assert(ret >= 0, ExcInternalError());
493 
494  ret = H5Aclose(attr);
495  Assert(ret >= 0, ExcInternalError());
496 
497  (void)ret;
498  return string_value;
499  }
500 
501 
502 
503  template <typename T>
504  void
505  HDF5Object::set_attribute(const std::string &attr_name, const T value)
506  {
507  hid_t attr;
508  hid_t aid;
509  herr_t ret;
510 
511  const std::shared_ptr<hid_t> t_type = internal::get_hdf5_datatype<T>();
512 
513 
514  /*
515  * Create scalar attribute.
516  */
517  aid = H5Screate(H5S_SCALAR);
518  Assert(aid >= 0, ExcMessage("Error at H5Screate"));
519  attr = H5Acreate2(*hdf5_reference,
520  attr_name.data(),
521  *t_type,
522  aid,
523  H5P_DEFAULT,
524  H5P_DEFAULT);
525  Assert(attr >= 0, ExcMessage("Error at H5Acreate2"));
526 
527  /*
528  * Write scalar attribute.
529  */
530  ret = H5Awrite(attr, *t_type, &value);
531  Assert(ret >= 0, ExcMessage("Error at H5Awrite"));
532 
533  ret = H5Sclose(aid);
534  Assert(ret >= 0, ExcMessage("Error at H5Sclose"));
535  ret = H5Aclose(attr);
536  Assert(ret >= 0, ExcMessage("Error at H5Aclose"));
537 
538  (void)ret;
539  }
540 
541 
542 
543  template <>
544  void
545  HDF5Object::set_attribute(const std::string &attr_name,
546  const std::string value) // NOLINT
547  {
548  // Writes a UTF8 variable string
549  //
550  // code inspired from
551  // https://support.hdfgroup.org/ftp/HDF5/examples/misc-examples/vlstratt.c
552 
553  hid_t attr;
554  hid_t aid;
555  hid_t t_type;
556  herr_t ret;
557 
558  /* Create a datatype to refer to. */
559  t_type = H5Tcopy(H5T_C_S1);
560  Assert(t_type >= 0, ExcInternalError());
561 
562  // Python strings are encoded in UTF8
563  ret = H5Tset_cset(t_type, H5T_CSET_UTF8);
564  Assert(t_type >= 0, ExcInternalError());
565 
566  ret = H5Tset_size(t_type, H5T_VARIABLE);
567  Assert(ret >= 0, ExcInternalError());
568 
569  /*
570  * Create scalar attribute.
571  */
572  aid = H5Screate(H5S_SCALAR);
573  Assert(aid >= 0, ExcMessage("Error at H5Screate"));
574  attr = H5Acreate2(
575  *hdf5_reference, attr_name.data(), t_type, aid, H5P_DEFAULT, H5P_DEFAULT);
576  Assert(attr >= 0, ExcMessage("Error at H5Acreate2"));
577 
578  /*
579  * Write scalar attribute.
580  * In most of the cases H5Awrite and H5Dwrite take a pointer to the data.
581  * But in the particular case of a variable length string, H5Awrite takes
582  * the address of the pointer of the string.
583  */
584  const char *c_string_value = value.c_str();
585  ret = H5Awrite(attr, t_type, &c_string_value);
586  Assert(ret >= 0, ExcInternalError());
587 
588  ret = H5Sclose(aid);
589  Assert(ret >= 0, ExcMessage("Error at H5Sclose"));
590  ret = H5Aclose(attr);
591  Assert(ret >= 0, ExcMessage("Error at H5Aclose"));
592 
593  (void)ret;
594  }
595 
596 
597 
598  std::string
600  {
601  return name;
602  }
603 
604 
605 
606  DataSet::DataSet(const std::string &name,
607  const hid_t & parent_group_id,
608  const bool mpi)
609  : HDF5Object(name, mpi)
610  , query_io_mode(false)
611  , io_mode(H5D_MPIO_NO_COLLECTIVE)
612  , local_no_collective_cause(H5D_MPIO_SET_INDEPENDENT)
613  , global_no_collective_cause(H5D_MPIO_SET_INDEPENDENT)
614  {
615  hdf5_reference = std::shared_ptr<hid_t>(new hid_t, [](hid_t *pointer) {
616  // Release the HDF5 resource
617  const herr_t ret = H5Dclose(*pointer);
618  AssertNothrow(ret >= 0, ExcInternalError());
619  (void)ret;
620  delete pointer;
621  });
622  dataspace = std::shared_ptr<hid_t>(new hid_t, [](hid_t *pointer) {
623  // Release the HDF5 resource
624  const herr_t ret = H5Sclose(*pointer);
625  AssertNothrow(ret >= 0, ExcInternalError());
626  (void)ret;
627  delete pointer;
628  });
629 
630  // rank_ret can take a negative value if the functions
631  // H5Sget_simple_extent_ndims and H5Sget_simple_extent_dims fail. rank is
632  // unsigned int therefore it can not be used to store the return value of
633  // H5Sget_simple_extent_ndims and H5Sget_simple_extent_dims
634  int rank_ret;
635 
636  *hdf5_reference = H5Dopen2(parent_group_id, name.data(), H5P_DEFAULT);
637  Assert(*hdf5_reference >= 0, ExcMessage("Error at H5Dopen2"));
638  *dataspace = H5Dget_space(*hdf5_reference);
639  Assert(*dataspace >= 0, ExcMessage("Error at H5Dget_space"));
640  rank_ret = H5Sget_simple_extent_ndims(*dataspace);
641  Assert(rank_ret >= 0, ExcInternalError());
642  rank = rank_ret;
643  dimensions.resize(rank);
644  rank_ret =
645  H5Sget_simple_extent_dims(*dataspace, dimensions.data(), nullptr);
646  AssertDimension(rank_ret, static_cast<int>(rank));
647 
648  size = 1;
649  for (const auto &dimension : dimensions)
650  {
651  size *= dimension;
652  }
653  }
654 
655 
656 
657  DataSet::DataSet(const std::string & name,
658  const hid_t & parent_group_id,
659  const std::vector<hsize_t> & dimensions,
660  const std::shared_ptr<hid_t> &t_type,
661  const bool mpi)
662  : HDF5Object(name, mpi)
663  , rank(dimensions.size())
664  , dimensions(dimensions)
665  , query_io_mode(false)
666  , io_mode(H5D_MPIO_NO_COLLECTIVE)
667  , local_no_collective_cause(H5D_MPIO_SET_INDEPENDENT)
668  , global_no_collective_cause(H5D_MPIO_SET_INDEPENDENT)
669  {
670  hdf5_reference = std::shared_ptr<hid_t>(new hid_t, [](hid_t *pointer) {
671  // Release the HDF5 resource
672  const herr_t ret = H5Dclose(*pointer);
673  AssertNothrow(ret >= 0, ExcInternalError());
674  (void)ret;
675  delete pointer;
676  });
677  dataspace = std::shared_ptr<hid_t>(new hid_t, [](hid_t *pointer) {
678  // Release the HDF5 resource
679  const herr_t ret = H5Sclose(*pointer);
680  AssertNothrow(ret >= 0, ExcInternalError());
681  (void)ret;
682  delete pointer;
683  });
684 
685  *dataspace = H5Screate_simple(rank, dimensions.data(), nullptr);
686  Assert(*dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
687 
688  *hdf5_reference = H5Dcreate2(parent_group_id,
689  name.data(),
690  *t_type,
691  *dataspace,
692  H5P_DEFAULT,
693  H5P_DEFAULT,
694  H5P_DEFAULT);
695  Assert(*hdf5_reference >= 0, ExcMessage("Error at H5Dcreate2"));
696 
697  size = 1;
698  for (const auto &dimension : dimensions)
699  {
700  size *= dimension;
701  }
702  }
703 
704 
705 
706  template <typename Container>
707  Container
709  {
710  const std::shared_ptr<hid_t> t_type =
711  internal::get_hdf5_datatype<typename Container::value_type>();
712  hid_t plist;
713  herr_t ret;
714 
715  Container data = internal::initialize_container<Container>(dimensions);
716 
717  internal::set_plist(plist, mpi);
718 
719  ret = H5Dread(*hdf5_reference,
720  *t_type,
721  H5S_ALL,
722  H5S_ALL,
723  plist,
724  make_array_view(data).data());
725  Assert(ret >= 0, ExcInternalError());
726 
727 
729  io_mode,
732  mpi,
733  query_io_mode);
734 
735  (void)ret;
736  return data;
737  }
738 
739 
740 
741  template <typename Container>
742  Container
743  DataSet::read_selection(const std::vector<hsize_t> &coordinates)
744  {
745  Assert(coordinates.size() % rank == 0,
746  ExcMessage(
747  "The dimension of coordinates has to be divisible by the rank"));
748  const std::shared_ptr<hid_t> t_type =
749  internal::get_hdf5_datatype<typename Container::value_type>();
750  hid_t plist;
751  hid_t memory_dataspace;
752  herr_t ret;
753 
754  std::vector<hsize_t> data_dimensions{
755  static_cast<hsize_t>(coordinates.size() / rank)};
756 
757  Container data = internal::initialize_container<Container>(data_dimensions);
758 
759  memory_dataspace = H5Screate_simple(1, data_dimensions.data(), nullptr);
760  Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
761  ret = H5Sselect_elements(*dataspace,
762  H5S_SELECT_SET,
763  data.size(),
764  coordinates.data());
765  Assert(ret >= 0, ExcMessage("Error at H5Sselect_elements"));
766 
767  internal::set_plist(plist, mpi);
768 
769  ret = H5Dread(*hdf5_reference,
770  *t_type,
771  memory_dataspace,
772  *dataspace,
773  plist,
774  make_array_view(data).data());
775  Assert(ret >= 0, ExcMessage("Error at H5Dread"));
776 
778  io_mode,
781  mpi,
782  query_io_mode);
783 
784  ret = H5Sclose(memory_dataspace);
785  Assert(ret >= 0, ExcMessage("Error at H5SClose"));
786 
787  (void)ret;
788  return data;
789  }
790 
791 
792 
793  template <typename Container>
794  Container
795  DataSet::read_hyperslab(const std::vector<hsize_t> &offset,
796  const std::vector<hsize_t> &count)
797  {
798  const std::shared_ptr<hid_t> t_type =
799  internal::get_hdf5_datatype<typename Container::value_type>();
800  hid_t plist;
801  hid_t memory_dataspace;
802  herr_t ret;
803 
804  // In this particular overload of read_hyperslab the data_dimensions are
805  // the same as count
806  std::vector<hsize_t> data_dimensions = count;
807 
808  Container data = internal::initialize_container<Container>(data_dimensions);
809 
810  memory_dataspace =
811  H5Screate_simple(data_dimensions.size(), data_dimensions.data(), nullptr);
812  Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
813  ret = H5Sselect_hyperslab(*dataspace,
814  H5S_SELECT_SET,
815  offset.data(),
816  nullptr,
817  count.data(),
818  nullptr);
819  Assert(ret >= 0, ExcMessage("Error at H5Sselect_hyperslab"));
820 
821  internal::set_plist(plist, mpi);
822 
823  ret = H5Dread(*hdf5_reference,
824  *t_type,
825  memory_dataspace,
826  *dataspace,
827  plist,
828  make_array_view(data).data());
829  Assert(ret >= 0, ExcMessage("Error at H5Dread"));
830 
832  io_mode,
835  mpi,
836  query_io_mode);
837 
838  ret = H5Sclose(memory_dataspace);
839  Assert(ret >= 0, ExcMessage("Error at H5SClose"));
840 
841  (void)ret;
842  return data;
843  }
844 
845 
846 
847  template <typename Container>
848  Container
849  DataSet::read_hyperslab(const std::vector<hsize_t> &data_dimensions,
850  const std::vector<hsize_t> &offset,
851  const std::vector<hsize_t> &stride,
852  const std::vector<hsize_t> &count,
853  const std::vector<hsize_t> &block)
854  {
855  const std::shared_ptr<hid_t> t_type =
856  internal::get_hdf5_datatype<typename Container::value_type>();
857  hid_t plist;
858  hid_t memory_dataspace;
859  herr_t ret;
860 
861  Container data = internal::initialize_container<Container>(data_dimensions);
862 
863  memory_dataspace =
864  H5Screate_simple(data_dimensions.size(), data_dimensions.data(), nullptr);
865  Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
866  ret = H5Sselect_hyperslab(*dataspace,
867  H5S_SELECT_SET,
868  offset.data(),
869  stride.data(),
870  count.data(),
871  block.data());
872  Assert(ret >= 0, ExcMessage("Error at H5Sselect_hyperslab"));
873 
874  internal::set_plist(plist, mpi);
875 
876  ret = H5Dread(*hdf5_reference,
877  *t_type,
878  memory_dataspace,
879  *dataspace,
880  plist,
881  make_array_view(data).data());
882  Assert(ret >= 0, ExcMessage("Error at H5Dread"));
883 
885  io_mode,
888  mpi,
889  query_io_mode);
890 
891  ret = H5Sclose(memory_dataspace);
892  Assert(ret >= 0, ExcMessage("Error at H5SClose"));
893 
894  (void)ret;
895  return data;
896  }
897 
898 
899 
900  template <typename number>
901  void
903  {
904  const std::shared_ptr<hid_t> t_type = internal::get_hdf5_datatype<number>();
905  const std::vector<hsize_t> data_dimensions = {0};
906 
907  hid_t memory_dataspace;
908  hid_t plist;
909  herr_t ret;
910 
911  memory_dataspace = H5Screate_simple(1, data_dimensions.data(), nullptr);
912  Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
913  ret = H5Sselect_none(*dataspace);
914  Assert(ret >= 0, ExcMessage("H5Sselect_none"));
915 
916  internal::set_plist(plist, mpi);
917 
918  // The pointer of data can safely be nullptr, see the discussion at the HDF5
919  // forum:
920  // https://forum.hdfgroup.org/t/parallel-i-o-does-not-support-filters-yet/884/17
921  ret = H5Dread(
922  *hdf5_reference, *t_type, memory_dataspace, *dataspace, plist, nullptr);
923  Assert(ret >= 0, ExcMessage("Error at H5Dread"));
924 
926  io_mode,
929  mpi,
930  query_io_mode);
931 
932  ret = H5Sclose(memory_dataspace);
933  Assert(ret >= 0, ExcMessage("Error at H5SClose"));
934 
935  (void)ret;
936  }
937 
938 
939 
940  template <typename Container>
941  void
942  DataSet::write(const Container &data)
943  {
945  const std::shared_ptr<hid_t> t_type =
946  internal::get_hdf5_datatype<typename Container::value_type>();
947  hid_t plist;
948  herr_t ret;
949 
950  internal::set_plist(plist, mpi);
951 
952  ret = H5Dwrite(*hdf5_reference,
953  *t_type,
954  H5S_ALL,
955  H5S_ALL,
956  plist,
957  make_array_view(data).data());
958  Assert(ret >= 0, ExcMessage("Error at H5Dwrite"));
959 
961  io_mode,
964  mpi,
965  query_io_mode);
966 
967  (void)ret;
968  }
969 
970 
971 
972  template <typename Container>
973  void
974  DataSet::write_selection(const Container & data,
975  const std::vector<hsize_t> &coordinates)
976  {
977  AssertDimension(coordinates.size(), data.size() * rank);
978  const std::shared_ptr<hid_t> t_type =
979  internal::get_hdf5_datatype<typename Container::value_type>();
980  const std::vector<hsize_t> data_dimensions =
982 
983  hid_t memory_dataspace;
984  hid_t plist;
985  herr_t ret;
986 
987 
988  memory_dataspace = H5Screate_simple(1, data_dimensions.data(), nullptr);
989  Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
990  ret = H5Sselect_elements(*dataspace,
991  H5S_SELECT_SET,
992  data.size(),
993  coordinates.data());
994  Assert(ret >= 0, ExcMessage("Error at H5Sselect_elements"));
995 
996  internal::set_plist(plist, mpi);
997 
998  ret = H5Dwrite(*hdf5_reference,
999  *t_type,
1000  memory_dataspace,
1001  *dataspace,
1002  plist,
1003  make_array_view(data).data());
1004  Assert(ret >= 0, ExcMessage("Error at H5Dwrite"));
1005 
1007  io_mode,
1010  mpi,
1011  query_io_mode);
1012 
1013  ret = H5Sclose(memory_dataspace);
1014  Assert(ret >= 0, ExcMessage("Error at H5SClose"));
1015 
1016  (void)ret;
1017  }
1018 
1019 
1020 
1021  template <typename Container>
1022  void
1023  DataSet::write_hyperslab(const Container & data,
1024  const std::vector<hsize_t> &offset,
1025  const std::vector<hsize_t> &count)
1026  {
1027  AssertDimension(std::accumulate(count.begin(),
1028  count.end(),
1029  1,
1030  std::multiplies<unsigned int>()),
1032  const std::shared_ptr<hid_t> t_type =
1033  internal::get_hdf5_datatype<typename Container::value_type>();
1034  // In this particular overload of write_hyperslab the data_dimensions are
1035  // the same as count
1036  const std::vector<hsize_t> &data_dimensions = count;
1037 
1038  hid_t memory_dataspace;
1039  hid_t plist;
1040  herr_t ret;
1041 
1042  memory_dataspace =
1043  H5Screate_simple(data_dimensions.size(), data_dimensions.data(), nullptr);
1044  Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
1045  ret = H5Sselect_hyperslab(*dataspace,
1046  H5S_SELECT_SET,
1047  offset.data(),
1048  nullptr,
1049  count.data(),
1050  nullptr);
1051  Assert(ret >= 0, ExcMessage("Error at H5Sselect_hyperslab"));
1052 
1053  internal::set_plist(plist, mpi);
1054 
1055  ret = H5Dwrite(*hdf5_reference,
1056  *t_type,
1057  memory_dataspace,
1058  *dataspace,
1059  plist,
1060  make_array_view(data).data());
1061  Assert(ret >= 0, ExcMessage("Error at H5Dwrite"));
1062 
1064  io_mode,
1067  mpi,
1068  query_io_mode);
1069 
1070  ret = H5Sclose(memory_dataspace);
1071  Assert(ret >= 0, ExcMessage("Error at H5Sclose"));
1072 
1073  (void)ret;
1074  }
1075 
1076 
1077 
1078  template <typename Container>
1079  void
1080  DataSet::write_hyperslab(const Container & data,
1081  const std::vector<hsize_t> &data_dimensions,
1082  const std::vector<hsize_t> &offset,
1083  const std::vector<hsize_t> &stride,
1084  const std::vector<hsize_t> &count,
1085  const std::vector<hsize_t> &block)
1086  {
1087  const std::shared_ptr<hid_t> t_type =
1088  internal::get_hdf5_datatype<typename Container::value_type>();
1089 
1090  hid_t memory_dataspace;
1091  hid_t plist;
1092  herr_t ret;
1093 
1094  memory_dataspace =
1095  H5Screate_simple(data_dimensions.size(), data_dimensions.data(), nullptr);
1096  Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
1097  ret = H5Sselect_hyperslab(*dataspace,
1098  H5S_SELECT_SET,
1099  offset.data(),
1100  stride.data(),
1101  count.data(),
1102  block.data());
1103  Assert(ret >= 0, ExcMessage("Error at H5Sselect_hyperslab"));
1104 
1105  internal::set_plist(plist, mpi);
1106 
1107  ret = H5Dwrite(*hdf5_reference,
1108  *t_type,
1109  memory_dataspace,
1110  *dataspace,
1111  plist,
1112  make_array_view(data).data());
1113  Assert(ret >= 0, ExcMessage("Error at H5Dwrite"));
1114 
1116  io_mode,
1119  mpi,
1120  query_io_mode);
1121 
1122  ret = H5Sclose(memory_dataspace);
1123  Assert(ret >= 0, ExcMessage("Error at H5Sclose"));
1124 
1125  (void)ret;
1126  }
1127 
1128 
1129 
1130  template <typename number>
1131  void
1133  {
1134  std::shared_ptr<hid_t> t_type = internal::get_hdf5_datatype<number>();
1135  std::vector<hsize_t> data_dimensions = {0};
1136 
1137  hid_t memory_dataspace;
1138  hid_t plist;
1139  herr_t ret;
1140 
1141  memory_dataspace = H5Screate_simple(1, data_dimensions.data(), nullptr);
1142  Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
1143  ret = H5Sselect_none(*dataspace);
1144  Assert(ret >= 0, ExcMessage("Error at H5PSselect_none"));
1145 
1146  internal::set_plist(plist, mpi);
1147 
1148  // The pointer of data can safely be nullptr, see the discussion at the HDF5
1149  // forum:
1150  // https://forum.hdfgroup.org/t/parallel-i-o-does-not-support-filters-yet/884/17
1151  ret = H5Dwrite(
1152  *hdf5_reference, *t_type, memory_dataspace, *dataspace, plist, nullptr);
1153  Assert(ret >= 0, ExcMessage("Error at H5Dwrite"));
1154 
1156  io_mode,
1159  mpi,
1160  query_io_mode);
1161 
1162  ret = H5Sclose(memory_dataspace);
1163  Assert(ret >= 0, ExcMessage("Error at H5Sclose"));
1164 
1165  (void)ret;
1166  }
1167 
1168 
1169 
1170  void
1171  DataSet::set_query_io_mode(const bool new_query_io_mode)
1172  {
1173  query_io_mode = new_query_io_mode;
1174  }
1175 
1176 
1177 
1178  std::vector<hsize_t>
1180  {
1181  return dimensions;
1182  }
1183 
1184 
1185 
1186  std::string
1188  {
1190  ExcMessage("query_io_mode should be true in order to use io_mode"));
1191  switch (io_mode)
1192  {
1193  case (H5D_MPIO_NO_COLLECTIVE):
1194  return std::string("H5D_MPIO_NO_COLLECTIVE");
1195  break;
1196  case (H5D_MPIO_CHUNK_INDEPENDENT):
1197  return std::string("H5D_MPIO_CHUNK_INDEPENDENT");
1198  break;
1199  case (H5D_MPIO_CHUNK_COLLECTIVE):
1200  return std::string("H5D_MPIO_CHUNK_COLLECTIVE");
1201  break;
1202  case (H5D_MPIO_CHUNK_MIXED):
1203  return std::string("H5D_MPIO_CHUNK_MIXED");
1204  break;
1205  case (H5D_MPIO_CONTIGUOUS_COLLECTIVE):
1206  return std::string("H5D_MPIO_CONTIGUOUS_COLLECTIVE");
1207  break;
1208  default:
1209  Assert(false, ExcInternalError());
1210  return std::string("Internal error");
1211  break;
1212  }
1213  // The function should not reach this line.
1214  Assert(false, ExcInternalError());
1215  return std::string("Internal error");
1216  }
1217 
1218 
1219 
1220  H5D_mpio_actual_io_mode_t
1222  {
1224  ExcMessage(
1225  "query_io_mode should be true in order to use get_io_mode"));
1226  return io_mode;
1227  }
1228 
1229 
1230 
1231  std::string
1233  {
1234  Assert(
1235  query_io_mode,
1236  ExcMessage(
1237  "query_io_mode should be true in order to use get_local_no_collective_cause()"));
1239  }
1240 
1241 
1242 
1243  uint32_t
1245  {
1246  Assert(
1247  query_io_mode,
1248  ExcMessage(
1249  "query_io_mode should be true in order to use get_local_no_collective_cause()"));
1251  }
1252 
1253 
1254 
1255  std::string
1257  {
1258  Assert(
1259  query_io_mode,
1260  ExcMessage(
1261  "query_io_mode should be true in order to use get_global_no_collective_cause()"));
1263  }
1264 
1265 
1266 
1267  uint32_t
1269  {
1270  Assert(
1271  query_io_mode,
1272  ExcMessage(
1273  "query_io_mode should be true in order to use get_global_no_collective_cause()"));
1275  }
1276 
1277 
1278  bool
1280  {
1281  return query_io_mode;
1282  }
1283 
1284 
1285 
1286  unsigned int
1288  {
1289  return size;
1290  }
1291 
1292 
1293 
1294  unsigned int
1296  {
1297  return rank;
1298  }
1299 
1300 
1301 
1302  Group::Group(const std::string & name,
1303  const Group & parentGroup,
1304  const bool mpi,
1305  const GroupAccessMode mode)
1306  : HDF5Object(name, mpi)
1307  {
1308  hdf5_reference = std::shared_ptr<hid_t>(new hid_t, [](hid_t *pointer) {
1309  // Release the HDF5 resource
1310  const herr_t ret = H5Gclose(*pointer);
1311  AssertNothrow(ret >= 0, ExcInternalError());
1312  (void)ret;
1313  delete pointer;
1314  });
1315  switch (mode)
1316  {
1317  case (GroupAccessMode::open):
1318  *hdf5_reference =
1319  H5Gopen2(*(parentGroup.hdf5_reference), name.data(), H5P_DEFAULT);
1320  Assert(*hdf5_reference >= 0, ExcMessage("Error at H5Gopen2"));
1321  break;
1322  case (GroupAccessMode::create):
1323  *hdf5_reference = H5Gcreate2(*(parentGroup.hdf5_reference),
1324  name.data(),
1325  H5P_DEFAULT,
1326  H5P_DEFAULT,
1327  H5P_DEFAULT);
1328  Assert(*hdf5_reference >= 0, ExcMessage("Error at H5Gcreate2"));
1329  break;
1330  default:
1331  Assert(false, ExcInternalError());
1332  break;
1333  }
1334  }
1335 
1336 
1337 
1338  Group::Group(const std::string &name, const bool mpi)
1339  : HDF5Object(name, mpi)
1340  {}
1341 
1342 
1343 
1344  Group
1345  Group::open_group(const std::string &name) const
1346  {
1347  return {name, *this, mpi, GroupAccessMode::open};
1348  }
1349 
1350 
1351 
1352  Group
1353  Group::create_group(const std::string &name) const
1354  {
1355  return {name, *this, mpi, GroupAccessMode::create};
1356  }
1357 
1358 
1359 
1360  DataSet
1361  Group::open_dataset(const std::string &name) const
1362  {
1363  return {name, *hdf5_reference, mpi};
1364  }
1365 
1366 
1367 
1368  template <typename number>
1369  DataSet
1370  Group::create_dataset(const std::string & name,
1371  const std::vector<hsize_t> &dimensions) const
1372  {
1373  std::shared_ptr<hid_t> t_type = internal::get_hdf5_datatype<number>();
1374  return {name, *hdf5_reference, dimensions, t_type, mpi};
1375  }
1376 
1377 
1378 
1379  template <typename Container>
1380  void
1381  Group::write_dataset(const std::string &name, const Container &data) const
1382  {
1383  std::vector<hsize_t> dimensions = internal::get_container_dimensions(data);
1384  auto dataset =
1385  create_dataset<typename Container::value_type>(name, dimensions);
1386  dataset.write(data);
1387  }
1388 
1389 
1390 
1391  File::File(const std::string &name, const FileAccessMode mode)
1392  : File(name,
1393  mode,
1394  false,
1395 # ifdef DEAL_II_WITH_MPI
1396  MPI_COMM_NULL
1397 # else
1398  0
1399 # endif
1400  )
1401  {}
1402 
1403 
1404 
1405  File::File(const std::string & name,
1406  const FileAccessMode mode,
1407  const MPI_Comm mpi_communicator)
1408  : File(name, mode, true, mpi_communicator)
1409  {}
1410 
1411 
1412 
1413  File::File(const std::string & name,
1414  const FileAccessMode mode,
1415  const bool mpi,
1416  const MPI_Comm mpi_communicator)
1417  : Group(name, mpi)
1418  {
1419  hdf5_reference = std::shared_ptr<hid_t>(new hid_t, [](hid_t *pointer) {
1420  // Release the HDF5 resource
1421  const herr_t err = H5Fclose(*pointer);
1422  AssertNothrow(err >= 0, ExcInternalError());
1423  (void)err;
1424  delete pointer;
1425  });
1426 
1427  hid_t plist;
1428  herr_t ret;
1429 
1430  if (mpi)
1431  {
1432 # ifdef DEAL_II_WITH_MPI
1433 # ifdef H5_HAVE_PARALLEL
1434  const MPI_Info info = MPI_INFO_NULL;
1435 
1436  plist = H5Pcreate(H5P_FILE_ACCESS);
1437  Assert(plist >= 0, ExcMessage("Error at H5Pcreate"));
1438  ret = H5Pset_fapl_mpio(plist, mpi_communicator, info);
1439  Assert(ret >= 0, ExcMessage("Error at H5Pset_fapl_mpio"));
1440 # else
1441  AssertThrow(false, ExcMessage("HDF5 parallel support is disabled."));
1442 # endif // H5_HAVE_PARALLEL
1443 # else
1444  AssertThrow(false, ExcMessage("MPI support is disabled."));
1445 # endif // DEAL_II_WITH_MPI
1446  }
1447  else
1448  {
1449  plist = H5P_DEFAULT;
1450  }
1451 
1452  switch (mode)
1453  {
1454  case (FileAccessMode::open):
1455  *hdf5_reference = H5Fopen(name.data(), H5F_ACC_RDWR, plist);
1456  Assert(*hdf5_reference >= 0, ExcMessage("Error at H5Fopen"));
1457  break;
1458  case (FileAccessMode::create):
1459  *hdf5_reference =
1460  H5Fcreate(name.data(), H5F_ACC_TRUNC, H5P_DEFAULT, plist);
1461  Assert(*hdf5_reference >= 0, ExcMessage("Error at H5Fcreate"));
1462  break;
1463  default:
1464  Assert(false, ExcInternalError());
1465  break;
1466  }
1467 
1468  if (mpi)
1469  {
1470  // Release the HDF5 resource
1471  ret = H5Pclose(plist);
1472  Assert(ret >= 0, ExcMessage("Error at H5Pclose"));
1473  }
1474 
1475  (void)ret;
1476  (void)mpi_communicator;
1477  }
1478 
1479 
1480 
1481 # ifndef DOXYGEN
1482 
1483  // instantiations of functions
1484 
1485 # include "hdf5.inst"
1486 
1487  // Instantiations for int and unsigned int can be found below. Only
1488  // std::vector is provided with instantiations with int and unsigned int.
1489  // Instantiations of FullMatrix and Vector are provided with float, double,
1490  // std::complex<float> and std::complex<double> in hdf5.inst
1491 
1492  template int
1493  HDF5Object::get_attribute<int>(const std::string &attr_name) const;
1494  template unsigned int
1495  HDF5Object::get_attribute<unsigned int>(const std::string &attr_name) const;
1496  // The specializations of HDF5Object::get_attribute<std::string>
1497  // and HDF5Object::get_attribute<bool> have been defined above
1498 
1499  template void
1500  HDF5Object::set_attribute<int>(const std::string &attr_name, int value);
1501  template void
1502  HDF5Object::set_attribute<unsigned int>(const std::string &attr_name,
1503  unsigned int value);
1504 
1505  template std::vector<int>
1506  DataSet::read<std::vector<int>>();
1507  template std::vector<unsigned int>
1508  DataSet::read<std::vector<unsigned int>>();
1509 
1510  template std::vector<int>
1511  DataSet::read_selection<std::vector<int>>(
1512  const std::vector<hsize_t> &coordinates);
1513  template std::vector<unsigned int>
1514  DataSet::read_selection<std::vector<unsigned int>>(
1515  const std::vector<hsize_t> &coordinates);
1516 
1517  template std::vector<int>
1518  DataSet::read_hyperslab<std::vector<int>>(const std::vector<hsize_t> &offset,
1519  const std::vector<hsize_t> &count);
1520  template std::vector<unsigned int>
1521  DataSet::read_hyperslab<std::vector<unsigned int>>(
1522  const std::vector<hsize_t> &offset,
1523  const std::vector<hsize_t> &count);
1524 
1525  template std::vector<int>
1526  DataSet::read_hyperslab<std::vector<int>>(
1527  const std::vector<hsize_t> &data_dimensions,
1528  const std::vector<hsize_t> &offset,
1529  const std::vector<hsize_t> &stride,
1530  const std::vector<hsize_t> &count,
1531  const std::vector<hsize_t> &block);
1532  template std::vector<unsigned int>
1533  DataSet::read_hyperslab<std::vector<unsigned int>>(
1534  const std::vector<hsize_t> &data_dimensions,
1535  const std::vector<hsize_t> &offset,
1536  const std::vector<hsize_t> &stride,
1537  const std::vector<hsize_t> &count,
1538  const std::vector<hsize_t> &block);
1539 
1540  template void
1541  DataSet::read_none<int>();
1542  template void
1543  DataSet::read_none<unsigned int>();
1544 
1545  template void
1546  DataSet::write<std::vector<int>>(const std::vector<int> &data);
1547  template void
1548  DataSet::write<std::vector<unsigned int>>(
1549  const std::vector<unsigned int> &data);
1550 
1551  template void
1552  DataSet::write_selection<std::vector<int>>(
1553  const std::vector<int> & data,
1554  const std::vector<hsize_t> &coordinates);
1555  template void
1556  DataSet::write_selection<std::vector<unsigned int>>(
1557  const std::vector<unsigned int> &data,
1558  const std::vector<hsize_t> & coordinates);
1559 
1560  template void
1561  DataSet::write_hyperslab<std::vector<int>>(const std::vector<int> & data,
1562  const std::vector<hsize_t> &offset,
1563  const std::vector<hsize_t> &count);
1564  template void
1565  DataSet::write_hyperslab<std::vector<unsigned int>>(
1566  const std::vector<unsigned int> &data,
1567  const std::vector<hsize_t> & offset,
1568  const std::vector<hsize_t> & count);
1569 
1570  template void
1571  DataSet::write_hyperslab<std::vector<int>>(
1572  const std::vector<int> & data,
1573  const std::vector<hsize_t> &data_dimensions,
1574  const std::vector<hsize_t> &offset,
1575  const std::vector<hsize_t> &stride,
1576  const std::vector<hsize_t> &count,
1577  const std::vector<hsize_t> &block);
1578  template void
1579  DataSet::write_hyperslab<std::vector<unsigned int>>(
1580  const std::vector<unsigned int> &data,
1581  const std::vector<hsize_t> & data_dimensions,
1582  const std::vector<hsize_t> & offset,
1583  const std::vector<hsize_t> & stride,
1584  const std::vector<hsize_t> & count,
1585  const std::vector<hsize_t> & block);
1586 
1587  template void
1588  DataSet::write_none<int>();
1589  template void
1590  DataSet::write_none<unsigned int>();
1591 
1592  template DataSet
1593  Group::create_dataset<int>(const std::string & name,
1594  const std::vector<hsize_t> &dimensions) const;
1595  template DataSet
1596  Group::create_dataset<unsigned int>(
1597  const std::string & name,
1598  const std::vector<hsize_t> &dimensions) const;
1599 
1600  template void
1601  Group::write_dataset<std::vector<int>>(const std::string & name,
1602  const std::vector<int> &data) const;
1603  template void
1604  Group::write_dataset<std::vector<unsigned int>>(
1605  const std::string & name,
1606  const std::vector<unsigned int> &data) const;
1607 
1608 # endif // DOXYGEN
1609 
1610 } // namespace HDF5
1611 
1613 
1614 #endif // DEAL_II_WITH_HDF5
HDF5::HDF5Object::get_attribute
T get_attribute(const std::string &attr_name) const
Definition: hdf5.cc:405
array_view.h
HDF5::File::FileAccessMode::open
@ open
HDF5::File::FileAccessMode::create
@ create
HDF5::DataSet::DataSet
DataSet(const std::string &name, const hid_t &parent_group_id, bool mpi)
Definition: hdf5.cc:606
HDF5::DataSet::read
Container read()
Definition: hdf5.cc:708
HDF5::File
Definition: hdf5.h:1067
HDF5::HDF5Object::name
const std::string name
Definition: hdf5.h:394
HDF5::HDF5Object::set_attribute
void set_attribute(const std::string &attr_name, const T value)
Definition: hdf5.cc:505
HDF5::Group::open_group
Group open_group(const std::string &name) const
Definition: hdf5.cc:1345
HDF5::internal::get_container_dimensions
std::vector< hsize_t > get_container_dimensions(const std::vector< number > &data)
Definition: hdf5.cc:129
HDF5::internal::initialize_container
std::enable_if< std::is_same< Container, std::vector< typename Container::value_type > >::value, Container >::type initialize_container(const std::vector< hsize_t > &dimensions)
Definition: hdf5.cc:209
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
HDF5::Group::create_dataset
DataSet create_dataset(const std::string &name, const std::vector< hsize_t > &dimensions) const
Definition: hdf5.cc:1370
FullMatrix::m
size_type m() const
HDF5::DataSet::read_selection
Container read_selection(const std::vector< hsize_t > &coordinates)
Definition: hdf5.cc:743
HDF5
Definition: hdf5.h:331
HDF5::DataSet
Definition: hdf5.h:416
MPI_Comm
make_array_view
ArrayView< typename std::remove_reference< typename std::iterator_traits< Iterator >::reference >::type, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
Definition: array_view.h:607
FullMatrix::n
size_type n() const
Vector::size
size_type size() const
HDF5::internal::set_plist
void set_plist(hid_t &plist, const bool mpi)
Definition: hdf5.cc:265
HDF5::DataSet::set_query_io_mode
void set_query_io_mode(const bool new_query_io_mode)
Definition: hdf5.cc:1171
HDF5::internal::get_hdf5_datatype
std::shared_ptr< hid_t > get_hdf5_datatype()
Definition: hdf5.cc:49
HDF5::HDF5Object
Definition: hdf5.h:338
HDF5::HDF5Object::mpi
const bool mpi
Definition: hdf5.h:408
HDF5::DataSet::local_no_collective_cause
uint32_t local_no_collective_cause
Definition: hdf5.h:953
HDF5::DataSet::get_dimensions
std::vector< hsize_t > get_dimensions() const
Definition: hdf5.cc:1179
HDF5::DataSet::get_global_no_collective_cause
std::string get_global_no_collective_cause()
Definition: hdf5.cc:1256
HDF5::Group::GroupAccessMode
GroupAccessMode
Definition: hdf5.h:974
HDF5::DataSet::read_hyperslab
Container read_hyperslab(const std::vector< hsize_t > &offset, const std::vector< hsize_t > &count)
Definition: hdf5.cc:795
Utilities::CUDA::free
void free(T *&pointer)
Definition: cuda.h:96
LAPACKSupport::T
static const char T
Definition: lapack_support.h:163
HDF5::DataSet::get_rank
unsigned int get_rank() const
Definition: hdf5.cc:1295
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
HDF5::DataSet::write
void write(const Container &data)
Definition: hdf5.cc:942
HDF5::DataSet::get_size
unsigned int get_size() const
Definition: hdf5.cc:1287
HDF5::DataSet::write_hyperslab
void write_hyperslab(const Container &data, const std::vector< hsize_t > &offset, const std::vector< hsize_t > &count)
Definition: hdf5.cc:1023
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
HDF5::DataSet::get_local_no_collective_cause
std::string get_local_no_collective_cause()
Definition: hdf5.cc:1232
HDF5::DataSet::get_query_io_mode
bool get_query_io_mode() const
Definition: hdf5.cc:1279
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
HDF5::Group::open_dataset
DataSet open_dataset(const std::string &name) const
Definition: hdf5.cc:1361
HDF5::Group::GroupAccessMode::create
@ create
HDF5::DataSet::get_global_no_collective_cause_as_hdf5_type
uint32_t get_global_no_collective_cause_as_hdf5_type()
Definition: hdf5.cc:1268
HDF5::DataSet::read_none
void read_none()
Definition: hdf5.cc:902
HDF5::DataSet::global_no_collective_cause
uint32_t global_no_collective_cause
Definition: hdf5.h:960
HDF5::HDF5Object::get_name
std::string get_name() const
Definition: hdf5.cc:599
HDF5::DataSet::get_local_no_collective_cause_as_hdf5_type
uint32_t get_local_no_collective_cause_as_hdf5_type()
Definition: hdf5.cc:1244
HDF5::DataSet::dataspace
std::shared_ptr< hid_t > dataspace
Definition: hdf5.h:925
HDF5::Group::write_dataset
void write_dataset(const std::string &name, const Container &data) const
Definition: hdf5.cc:1381
value
static const bool value
Definition: dof_tools_constraints.cc:433
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
HDF5::DataSet::rank
unsigned int rank
Definition: hdf5.h:914
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
HDF5::File::FileAccessMode
FileAccessMode
Definition: hdf5.h:1073
HDF5::Group::GroupAccessMode::open
@ open
vector.h
HDF5::internal::get_container_size
unsigned int get_container_size(const std::vector< number > &data)
Definition: hdf5.cc:163
HDF5::Group
Definition: hdf5.h:968
HDF5::HDF5Object::HDF5Object
HDF5Object(const std::string &name, const bool mpi)
Definition: hdf5.cc:396
HDF5::internal::release_plist
void release_plist(hid_t &plist, H5D_mpio_actual_io_mode_t &io_mode, uint32_t &local_no_collective_cause, uint32_t &global_no_collective_cause, const bool mpi, const bool query_io_mode)
Definition: hdf5.cc:296
HDF5::Group::Group
Group(const std::string &name, const Group &parent_group, const bool mpi, const GroupAccessMode mode)
Definition: hdf5.cc:1302
HDF5::File::File
File(const std::string &name, const FileAccessMode mode)
Definition: hdf5.cc:1391
config.h
AssertNothrow
#define AssertNothrow(cond, exc)
Definition: exceptions.h:1483
hdf5.h
FullMatrix
Definition: full_matrix.h:71
HDF5::DataSet::write_none
void write_none()
Definition: hdf5.cc:1132
HDF5::DataSet::get_io_mode
std::string get_io_mode()
Definition: hdf5.cc:1187
internal
Definition: aligned_vector.h:369
HDF5::DataSet::write_selection
void write_selection(const Container &data, const std::vector< hsize_t > &coordinates)
Definition: hdf5.cc:974
HDF5::HDF5Object::hdf5_reference
std::shared_ptr< hid_t > hdf5_reference
Definition: hdf5.h:403
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
HDF5::DataSet::dimensions
std::vector< hsize_t > dimensions
Definition: hdf5.h:920
HDF5::DataSet::query_io_mode
bool query_io_mode
Definition: hdf5.h:941
HDF5::DataSet::io_mode
H5D_mpio_actual_io_mode_t io_mode
Definition: hdf5.h:946
Vector< number >
HDF5::DataSet::size
unsigned int size
Definition: hdf5.h:930
AssertThrow
#define AssertThrow(cond, exc)
Definition: exceptions.h:1531
HDF5::DataSet::get_io_mode_as_hdf5_type
H5D_mpio_actual_io_mode_t get_io_mode_as_hdf5_type()
Definition: hdf5.cc:1221
full_matrix.h
HDF5::Group::create_group
Group create_group(const std::string &name) const
Definition: hdf5.cc:1353
HDF5::internal::no_collective_cause_to_string
std::string no_collective_cause_to_string(const uint32_t no_collective_cause)
Definition: hdf5.cc:336