]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add support for all the scalar and complex types
authorDaniel Garcia-Sanchez <daniel.garcia-sanchez@insp.upmc.fr>
Thu, 2 Aug 2018 18:16:56 +0000 (20:16 +0200)
committerDaniel Garcia-Sanchez <daniel.garcia-sanchez@insp.upmc.fr>
Wed, 9 Jan 2019 17:09:15 +0000 (18:09 +0100)
Add FullMatrix support in write_dataset

Modify the namespace

Place get_hdf5_datatype in the internal namespace

include/deal.II/base/hdf5.h
source/base/hdf5.cc

index 08e9c981fba21f126d277a7a25a6b7a3a6e7a2d4..4d34cd66464d576677796b55966856e6c430d8d6 100644 (file)
@@ -32,7 +32,9 @@ DEAL_II_NAMESPACE_OPEN
  * This set of classes can be used to store the data in the hdf5 file format.
  *
  * The classes can be used to improve the interaction with python scripts.
- *
+ */
+// clang-format off
+/**
  * This python script writes the parameters for a deal.ii simulation:
  * @code
  * h5_file = h5py.File('simulation.hdf5','w')
@@ -45,8 +47,7 @@ DEAL_II_NAMESPACE_OPEN
  *
  * C++ deal.ii simulation with MPI HDF5:
  * @code
- * hdf5::File data_file("simulation.hdf5", MPI_COMM_WORLD,
- *                      hdf5::File::Mode::open);
+ * hdf5::File data_file("simulation.hdf5", MPI_COMM_WORLD, hdf5::File::Mode::open);
  * hdf5::Group data = data_file.group("data");
  *
  * // Read some parameters
@@ -56,15 +57,13 @@ DEAL_II_NAMESPACE_OPEN
  * const std::string simulation_type =
  * data.attr<std::string>("simulation_type");
  *
- * // Create some distributed datasets
- * hdf5::DataSet<double> frequency_dataset(data.template
- *                          create_dataset<double>("frequency",
- *                          std::vector<hsize_t>{nb_frequency_points}));
+ * // Create some distributed datasets.
+ * // Note that the keyword template is used between data and create because
+ * // the function create_dataset returns a dependent template.
+ * auto frequency_dataset = data.template create_dataset<double>("frequency", std::vector<hsize_t>{nb_frequency_points});
  *
- * hdf5::DataSet<std::complex<double>> displacement(parameters.data.template
- * create_dataset<std::complex<double>>("displacement",
- *                                      std::vector<hsize_t>{nb_slice_points,
- *                                      nb_frequency_points}));
+ * auto displacement = data.template create_dataset<std::complex<double>>("displacement",
+ *                                                           std::vector<hsize_t>{nb_slice_points, nb_frequency_points}));
  *
  * // ... do some calculations ...
  * // The displacement data is distributed across the MPI processes.
@@ -100,12 +99,26 @@ DEAL_II_NAMESPACE_OPEN
  * displacement = data['displacement'] # complex128 dtype
  * print('Degrees of freedom: ' + repr(data.attrs['degrees_of_freedom']))
  * @endcode
+ */
+// clang-format on
+/**
  *
  * @author Daniel Garcia-Sanchez, 2018
  */
-namespace hdf5
+namespace HDF5
 
 {
+  namespace internal
+  {
+    // This function gives the HDF5 datatype corresponding to the C++ type. In
+    // the case of std::complex types the HDF5 handlers are automatically freed
+    // using the destructor of std::shared_ptr.
+    template <typename T>
+    std::shared_ptr<hid_t>
+    get_hdf5_datatype();
+  } // namespace internal
+
+
   /**
    * General class for the HDF5 objects.
    *
@@ -204,7 +217,7 @@ namespace hdf5
      * User's Guide.
      */
     void
-    write_data(const dealii::FullMatrix<T> &data) const;
+    write_data(const FullMatrix<T> &data) const;
 
     /**
      * Writes data to a subset of the dataset. T can be double, int, unsigned
@@ -262,9 +275,9 @@ namespace hdf5
      * User's Guide.
      */
     void
-    write_data_hyperslab(const dealii::FullMatrix<T> &data,
-                         const std::vector<hsize_t>   offset,
-                         const std::vector<hsize_t>   count) const;
+    write_data_hyperslab(const FullMatrix<T> &      data,
+                         const std::vector<hsize_t> offset,
+                         const std::vector<hsize_t> count) const;
     void
                                write_data_none() const;
     const unsigned int         rank;
@@ -300,6 +313,25 @@ namespace hdf5
     Group
     create_group(const std::string name);
 
+    /**
+     * Creates a dataset. T can be double, int, unsigned int, bool or
+     * std::complex<double>.
+     *
+     * The keyword template has to be used between create_dataset and the name
+     * of the object because the function create_dataset returns a dependent
+     * template. See the example in HDF5.
+     *
+     * Datatype conversion takes place at the time of a read or write and is
+     * automatic. See the <a
+     * href="https://support.hdfgroup.org/HDF5/doc/UG/HDF5_Users_Guide-Responsive%20HTML5/index.html#t=HDF5_Users_Guide%2FDatatypes%2FHDF5_Datatypes.htm%23TOC_6_10_Data_Transferbc-26&rhtocid=6.5_2">"Data
+     * Transfer: Datatype Conversion and Selection"</a>  section in the HDF5
+     * User's Guide.
+     */
+    template <typename T>
+    DataSet<T>
+    create_dataset(const std::string          name,
+                   const std::vector<hsize_t> dimensions) const;
+
     /**
      * Creates and writes data to a dataset. T can be double, int, unsigned int,
      * bool or std::complex<double>.
@@ -312,12 +344,11 @@ namespace hdf5
      */
     template <typename T>
     void
-    write_dataset(const std::string name, const T &data) const;
-    template <typename T>
+    write_dataset(const std::string name, const std::vector<T> &data) const;
 
     /**
-     * Creates a dataset. T can be double, int, unsigned int, bool or
-     * std::complex<double>.
+     * Creates and writes data to a dataset. T can be double, int, unsigned int,
+     * bool or std::complex<double>.
      *
      * Datatype conversion takes place at the time of a read or write and is
      * automatic. See the <a
@@ -325,9 +356,9 @@ namespace hdf5
      * Transfer: Datatype Conversion and Selection"</a>  section in the HDF5
      * User's Guide.
      */
-    DataSet<T>
-    create_dataset(const std::string          name,
-                   const std::vector<hsize_t> dimensions) const;
+    template <typename T>
+    void
+    write_dataset(const std::string name, const FullMatrix<T> &data) const;
   };
 
   /**
@@ -356,7 +387,7 @@ namespace hdf5
      */
     File(const std::string name, const Mode mode = Mode::create);
   };
-} // namespace hdf5
+} // namespace HDF5
 
 DEAL_II_NAMESPACE_CLOSE
 
index 73384b4c788658def5d11794304de4b109dfdcc7..80c28624d7d4b6dcba6a06b8089e9bc9b3b226d8 100644 (file)
 
 DEAL_II_NAMESPACE_OPEN
 
-namespace hdf5
+namespace HDF5
 
 {
+  namespace internal
+  {
+    template <typename T>
+    std::shared_ptr<hid_t>
+    get_hdf5_datatype()
+    {
+      std::shared_ptr<hid_t> t_type;
+      if (std::is_same<T, float>::value)
+        {
+          t_type  = std::shared_ptr<hid_t>(new hid_t);
+          *t_type = H5T_NATIVE_FLOAT;
+        }
+      else if (std::is_same<T, double>::value)
+        {
+          t_type  = std::shared_ptr<hid_t>(new hid_t);
+          *t_type = H5T_NATIVE_DOUBLE;
+        }
+      else if (std::is_same<T, long double>::value)
+        {
+          t_type  = std::shared_ptr<hid_t>(new hid_t);
+          *t_type = H5T_NATIVE_LDOUBLE;
+        }
+      else if (std::is_same<T, int>::value)
+        {
+          t_type  = std::shared_ptr<hid_t>(new hid_t);
+          *t_type = H5T_NATIVE_INT;
+        }
+      else if (std::is_same<T, unsigned int>::value)
+        {
+          t_type  = std::shared_ptr<hid_t>(new hid_t);
+          *t_type = H5T_NATIVE_UINT;
+        }
+      else if (std::is_same<T, std::complex<float>>::value)
+        {
+          t_type  = std::shared_ptr<hid_t>(new hid_t, [](auto pointer) {
+            // Relase the HDF5 resource
+            H5Tclose(*pointer);
+            delete pointer;
+          });
+          *t_type = H5Tcreate(H5T_COMPOUND, sizeof(std::complex<float>));
+          //  The C++ standards committee agreed to mandate that the storage
+          //  format used for the std::complex type be binary-compatible with
+          //  the C99 type, i.e. an array T[2] with consecutive real [0] and
+          //  imaginary [1] parts.
+          H5Tinsert(*t_type, "r", 0, H5T_NATIVE_FLOAT);
+          H5Tinsert(*t_type, "i", sizeof(float), H5T_NATIVE_FLOAT);
+        }
+      else if (std::is_same<T, std::complex<double>>::value)
+        {
+          t_type  = std::shared_ptr<hid_t>(new hid_t, [](auto pointer) {
+            // Relase the HDF5 resource
+            H5Tclose(*pointer);
+            delete pointer;
+          });
+          *t_type = H5Tcreate(H5T_COMPOUND, sizeof(std::complex<double>));
+          //  The C++ standards committee agreed to mandate that the storage
+          //  format used for the std::complex type be binary-compatible with
+          //  the C99 type, i.e. an array T[2] with consecutive real [0] and
+          //  imaginary [1] parts.
+          H5Tinsert(*t_type, "r", 0, H5T_NATIVE_DOUBLE);
+          H5Tinsert(*t_type, "i", sizeof(double), H5T_NATIVE_DOUBLE);
+        }
+      else if (std::is_same<T, std::complex<long double>>::value)
+        {
+          t_type  = std::shared_ptr<hid_t>(new hid_t, [](auto pointer) {
+            // Relase the HDF5 resource
+            H5Tclose(*pointer);
+            delete pointer;
+          });
+          *t_type = H5Tcreate(H5T_COMPOUND, sizeof(std::complex<long double>));
+          //  The C++ standards committee agreed to mandate that the storage
+          //  format used for the std::complex type be binary-compatible with
+          //  the C99 type, i.e. an array T[2] with consecutive real [0] and
+          //  imaginary [1] parts.
+          H5Tinsert(*t_type, "r", 0, H5T_NATIVE_LDOUBLE);
+          H5Tinsert(*t_type, "i", sizeof(long double), H5T_NATIVE_LDOUBLE);
+        }
+      else
+        {
+          Assert(false, ExcInternalError());
+        }
+      return t_type;
+    }
+  } // namespace internal
+
+
   HDF5Object::HDF5Object(const std::string name, bool mpi)
     : name(name)
     , mpi(mpi)
@@ -40,43 +126,10 @@ namespace hdf5
   T
   HDF5Object::attr(const std::string attr_name) const
   {
-    std::shared_ptr<hid_t> t_type;
+    std::shared_ptr<hid_t> t_type = internal::get_hdf5_datatype<T>();
     T                      value;
     hid_t                  attr;
-    if (std::is_same<T, double>::value)
-      {
-        t_type  = std::shared_ptr<hid_t>(new hid_t);
-        *t_type = H5T_NATIVE_DOUBLE;
-      }
-    else if (std::is_same<T, int>::value)
-      {
-        t_type  = std::shared_ptr<hid_t>(new hid_t);
-        *t_type = H5T_NATIVE_INT;
-      }
-    else if (std::is_same<T, unsigned int>::value)
-      {
-        t_type  = std::shared_ptr<hid_t>(new hid_t);
-        *t_type = H5T_NATIVE_UINT;
-      }
-    else if (std::is_same<T, std::complex<double>>::value)
-      {
-        t_type  = std::shared_ptr<hid_t>(new hid_t, [](auto pointer) {
-          // Relase the HDF5 resource
-          H5Tclose(*pointer);
-          delete pointer;
-        });
-        *t_type = H5Tcreate(H5T_COMPOUND, sizeof(std::complex<double>));
-        //  The C++ standards committee agreed to mandate that the storage
-        //  format used for the std::complex type be binary-compatible with the
-        //  C99 type, i.e. an array T[2] with consecutive real [0] and imaginary
-        //  [1] parts.
-        H5Tinsert(*t_type, "r", 0, H5T_NATIVE_DOUBLE);
-        H5Tinsert(*t_type, "i", sizeof(double), H5T_NATIVE_DOUBLE);
-      }
-    else
-      {
-        Assert(false, ExcInternalError());
-      }
+
 
     attr = H5Aopen(*hdf5_reference, attr_name.data(), H5P_DEFAULT);
     H5Aread(attr, *t_type, &value);
@@ -151,7 +204,7 @@ namespace hdf5
   }
 
   template <>
-  std::vector<std::vector<double>>
+  FullMatrix<double>
   HDF5Object::attr(const std::string attr_name) const
   {
     hid_t attr;
@@ -162,8 +215,7 @@ namespace hdf5
     hsize_t dims[2];
     H5Sget_simple_extent_dims(attr_space, dims, NULL);
 
-    std::vector<std::vector<double>> matrix_value(dims[0],
-                                                  std::vector<double>(dims[1]));
+    FullMatrix<double> matrix_value(dims[0], dims[1]);
 
     double *hdf5_data = (double *)malloc(dims[0] * dims[1] * sizeof(double));
     H5Aread(attr, H5T_NATIVE_DOUBLE, hdf5_data);
@@ -190,42 +242,8 @@ namespace hdf5
   {
     hid_t                  attr;
     hid_t                  aid;
-    std::shared_ptr<hid_t> t_type;
+    std::shared_ptr<hid_t> t_type = internal::get_hdf5_datatype<T>();
 
-    if (std::is_same<T, double>::value)
-      {
-        t_type  = std::shared_ptr<hid_t>(new hid_t);
-        *t_type = H5T_NATIVE_DOUBLE;
-      }
-    else if (std::is_same<T, int>::value)
-      {
-        t_type  = std::shared_ptr<hid_t>(new hid_t);
-        *t_type = H5T_NATIVE_INT;
-      }
-    else if (std::is_same<T, unsigned int>::value)
-      {
-        t_type  = std::shared_ptr<hid_t>(new hid_t);
-        *t_type = H5T_NATIVE_UINT;
-      }
-    else if (std::is_same<T, std::complex<double>>::value)
-      {
-        t_type  = std::shared_ptr<hid_t>(new hid_t, [](auto pointer) {
-          // Relase the HDF5 resource
-          H5Tclose(*pointer);
-          delete pointer;
-        });
-        *t_type = H5Tcreate(H5T_COMPOUND, sizeof(std::complex<double>));
-        //  The C++ standards committee agreed to mandate that the storage
-        //  format used for the std::complex type be binary-compatible with the
-        //  C99 type, i.e. an array T[2] with consecutive real [0] and imaginary
-        //  [1] parts.
-        H5Tinsert(*t_type, "r", 0, H5T_NATIVE_DOUBLE);
-        H5Tinsert(*t_type, "i", sizeof(double), H5T_NATIVE_DOUBLE);
-      }
-    else
-      {
-        Assert(false, ExcInternalError());
-      }
 
     /*
      * Create scalar attribute.
@@ -256,6 +274,7 @@ namespace hdf5
     : HDF5Object(name, mpi)
     , rank(dimensions.size())
     , dimensions(dimensions)
+    , t_type(internal::get_hdf5_datatype<T>())
   {
     hdf5_reference = std::shared_ptr<hid_t>(new hid_t, [](auto pointer) {
       // Relase the HDF5 resource
@@ -268,40 +287,6 @@ namespace hdf5
       delete pointer;
     });
 
-    if (std::is_same<T, double>::value)
-      {
-        t_type  = std::shared_ptr<hid_t>(new hid_t);
-        *t_type = H5T_NATIVE_DOUBLE;
-      }
-    else if (std::is_same<T, int>::value)
-      {
-        t_type  = std::shared_ptr<hid_t>(new hid_t);
-        *t_type = H5T_NATIVE_INT;
-      }
-    else if (std::is_same<T, unsigned int>::value)
-      {
-        t_type  = std::shared_ptr<hid_t>(new hid_t);
-        *t_type = H5T_NATIVE_UINT;
-      }
-    else if (std::is_same<T, std::complex<double>>::value)
-      {
-        t_type  = std::shared_ptr<hid_t>(new hid_t, [](auto pointer) {
-          // Relase the HDF5 resource
-          H5Tclose(*pointer);
-          delete pointer;
-        });
-        *t_type = H5Tcreate(H5T_COMPOUND, sizeof(std::complex<double>));
-        //  The C++ standards committee agreed to mandate that the storage
-        //  format used for the std::complex type be binary-compatible with the
-        //  C99 type, i.e. an array T[2] with consecutive real [0] and imaginary
-        //  [1] parts.
-        H5Tinsert(*t_type, "r", 0, H5T_NATIVE_DOUBLE);
-        H5Tinsert(*t_type, "i", sizeof(double), H5T_NATIVE_DOUBLE);
-      }
-    else
-      {
-        Assert(false, ExcInternalError());
-      }
 
     total_size = 1;
     for (auto &&dimension : dimensions)
@@ -368,7 +353,7 @@ namespace hdf5
 
   template <typename T>
   void
-  DataSet<T>::write_data(const dealii::FullMatrix<T> &data) const
+  DataSet<T>::write_data(const FullMatrix<T> &data) const
   {
     AssertDimension(total_size, data.m() * data.n());
     hid_t plist;
@@ -477,9 +462,9 @@ namespace hdf5
 
   template <typename T>
   void
-  DataSet<T>::write_data_hyperslab(const dealii::FullMatrix<T> &data,
-                                   const std::vector<hsize_t>   offset,
-                                   const std::vector<hsize_t>   count) const
+  DataSet<T>::write_data_hyperslab(const FullMatrix<T> &      data,
+                                   const std::vector<hsize_t> offset,
+                                   const std::vector<hsize_t> count) const
   {
     AssertDimension(std::accumulate(count.begin(),
                                     count.end(),
@@ -597,34 +582,33 @@ namespace hdf5
     return Group(name, *this, mpi, Mode::create);
   }
 
-  template <>
+  template <typename T>
+  DataSet<T>
+  Group::create_dataset(const std::string          name,
+                        const std::vector<hsize_t> dimensions) const
+  {
+    return DataSet<T>(
+      name, *hdf5_reference, dimensions, mpi, DataSet<T>::Mode::create);
+  }
+
+  template <typename T>
   void
-  Group::write_dataset(const std::string          name,
-                       const std::vector<double> &data) const
+  Group::write_dataset(const std::string name, const std::vector<T> &data) const
   {
     std::vector<hsize_t> dimensions = {data.size()};
-    auto                 dataset    = create_dataset<double>(name, dimensions);
+    auto                 dataset    = create_dataset<T>(name, dimensions);
     dataset.write_data(data);
   }
 
-  template <>
+  template <typename T>
   void
-  Group::write_dataset(const std::string                 name,
-                       const dealii::FullMatrix<double> &data) const
+  Group::write_dataset(const std::string name, const FullMatrix<T> &data) const
   {
     std::vector<hsize_t> dimensions = {data.m(), data.n()};
-    auto                 dataset    = create_dataset<double>(name, dimensions);
+    auto                 dataset    = create_dataset<T>(name, dimensions);
     dataset.write_data(data);
   }
 
-  template <typename T>
-  DataSet<T>
-  Group::create_dataset(const std::string          name,
-                        const std::vector<hsize_t> dimensions) const
-  {
-    return DataSet<T>(
-      name, *hdf5_reference, dimensions, mpi, DataSet<T>::Mode::create);
-  }
 
   File::File(const std::string name,
              const bool        mpi,
@@ -684,36 +668,84 @@ namespace hdf5
 
 
   // explicit instantiations of classes
+  template class DataSet<float>;
   template class DataSet<double>;
+  template class DataSet<long double>;
+  template class DataSet<std::complex<float>>;
   template class DataSet<std::complex<double>>;
+  template class DataSet<std::complex<long double>>;
+  template class DataSet<int>;
+  template class DataSet<unsigned int>;
+
 
   // explicit instantiations of functions
+  template float
+  HDF5Object::attr<float>(const std::string attr_name) const;
   template double
   HDF5Object::attr<double>(const std::string attr_name) const;
+  template long double
+  HDF5Object::attr<long double>(const std::string attr_name) const;
+  template std::complex<float>
+  HDF5Object::attr<std::complex<float>>(const std::string attr_name) const;
+  template std::complex<double>
+  HDF5Object::attr<std::complex<double>>(const std::string attr_name) const;
+  template std::complex<long double>
+  HDF5Object::attr<std::complex<long double>>(
+    const std::string attr_name) const;
   template int
   HDF5Object::attr<int>(const std::string attr_name) const;
   template unsigned int
   HDF5Object::attr<unsigned int>(const std::string attr_name) const;
-  template std::complex<double>
-  HDF5Object::attr<std::complex<double>>(const std::string attr_name) const;
   // The specialization of HDF5Object::attr<std::string> has been defined above
 
+  template void
+  HDF5Object::write_attr<float>(const std::string attr_name, float value) const;
   template void
   HDF5Object::write_attr<double>(const std::string attr_name,
                                  double            value) const;
   template void
-  HDF5Object::write_attr<int>(const std::string attr_name, int value) const;
+  HDF5Object::write_attr<long double>(const std::string attr_name,
+                                      long double       value) const;
   template void
-  HDF5Object::write_attr<unsigned int>(const std::string attr_name,
-                                       unsigned int      value) const;
+  HDF5Object::write_attr<std::complex<float>>(const std::string   attr_name,
+                                              std::complex<float> value) const;
   template void
   HDF5Object::write_attr<std::complex<double>>(
     const std::string    attr_name,
     std::complex<double> value) const;
+  template void
+  HDF5Object::write_attr<std::complex<long double>>(
+    const std::string         attr_name,
+    std::complex<long double> value) const;
+  template void
+  HDF5Object::write_attr<int>(const std::string attr_name, int value) const;
+  template void
+  HDF5Object::write_attr<unsigned int>(const std::string attr_name,
+                                       unsigned int      value) const;
+
 
+  template DataSet<float>
+  Group::create_dataset<float>(const std::string          name,
+                               const std::vector<hsize_t> dimensions) const;
   template DataSet<double>
   Group::create_dataset<double>(const std::string          name,
                                 const std::vector<hsize_t> dimensions) const;
+  template DataSet<long double>
+  Group::create_dataset<long double>(
+    const std::string          name,
+    const std::vector<hsize_t> dimensions) const;
+  template DataSet<std::complex<float>>
+  Group::create_dataset<std::complex<float>>(
+    const std::string          name,
+    const std::vector<hsize_t> dimensions) const;
+  template DataSet<std::complex<double>>
+  Group::create_dataset<std::complex<double>>(
+    const std::string          name,
+    const std::vector<hsize_t> dimensions) const;
+  template DataSet<std::complex<long double>>
+  Group::create_dataset<std::complex<long double>>(
+    const std::string          name,
+    const std::vector<hsize_t> dimensions) const;
   template DataSet<int>
   Group::create_dataset<int>(const std::string          name,
                              const std::vector<hsize_t> dimensions) const;
@@ -721,12 +753,51 @@ namespace hdf5
   Group::create_dataset<unsigned int>(
     const std::string          name,
     const std::vector<hsize_t> dimensions) const;
-  template DataSet<std::complex<double>>
-  Group::create_dataset<std::complex<double>>(
-    const std::string          name,
-    const std::vector<hsize_t> dimensions) const;
 
+  template void
+  Group::write_dataset(const std::string         name,
+                       const std::vector<float> &data) const;
+  template void
+  Group::write_dataset(const std::string          name,
+                       const std::vector<double> &data) const;
+  template void
+  Group::write_dataset(const std::string               name,
+                       const std::vector<long double> &data) const;
+  template void
+  Group::write_dataset(const std::string                       name,
+                       const std::vector<std::complex<float>> &data) const;
+  template void
+  Group::write_dataset(const std::string                        name,
+                       const std::vector<std::complex<double>> &data) const;
+  template void
+  Group::write_dataset(
+    const std::string                             name,
+    const std::vector<std::complex<long double>> &data) const;
+  template void
+  Group::write_dataset(const std::string       name,
+                       const std::vector<int> &data) const;
+  template void
+  Group::write_dataset(const std::string                name,
+                       const std::vector<unsigned int> &data) const;
 
-} // namespace hdf5
+  template void
+  Group::write_dataset(const std::string        name,
+                       const FullMatrix<float> &data) const;
+  template void
+  Group::write_dataset(const std::string         name,
+                       const FullMatrix<double> &data) const;
+  template void
+  Group::write_dataset(const std::string              name,
+                       const FullMatrix<long double> &data) const;
+  template void
+  Group::write_dataset(const std::string                      name,
+                       const FullMatrix<std::complex<float>> &data) const;
+  template void
+  Group::write_dataset(const std::string                       name,
+                       const FullMatrix<std::complex<double>> &data) const;
+  template void
+  Group::write_dataset(const std::string                            name,
+                       const FullMatrix<std::complex<long double>> &data) const;
+} // namespace HDF5
 
 DEAL_II_NAMESPACE_CLOSE

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.