isce3  0.1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
Serialization.h
Go to the documentation of this file.
1 //-*- C++ -*-
2 //-*- coding: utf-8 -*-
3 //
4 // Author: Bryan V. Riel
5 // Copyright 2017-2018
6 
11 #pragma once
12 
13 #include <iostream>
14 #include <sstream>
15 
16 // isce3::core
17 #include <isce3/core/DateTime.h>
18 
19 // isce3::io
20 #include <isce3/io/IH5.h>
21 
23 namespace isce3 {
25  namespace io {
26 
33  template <typename H5obj>
34  inline bool exists(H5obj & h5obj, const std::string & name,
35  const std::string start = ".", const std::string type = "BOTH") {
36  // Search for objects
37  std::vector<std::string> objs = h5obj.find(name, start, type, "FULL");
38  // Check size of vector
39  if (objs.size() > 0) {
40  return true;
41  } else {
42  return false;
43  }
44  }
45 
51  template <typename H5obj, typename T>
52  inline void loadFromH5(H5obj & h5obj, const std::string & datasetPath, T & v) {
53  // Open dataset
54  isce3::io::IDataSet dataset = h5obj.openDataSet(datasetPath);
55  // Read the scalar dataset
56  dataset.read(v);
57  }
58 
64  template <typename H5obj, typename T>
65  inline void loadFromH5(H5obj & h5obj, const std::string & datasetPath,
66  std::vector<T> & v) {
67  // Open dataset
68  isce3::io::IDataSet dataset = h5obj.openDataSet(datasetPath);
69  // Read the vector dataset
70  dataset.read(v);
71  }
72 
78  template <typename H5obj, typename T>
79  inline void loadFromH5(H5obj & h5obj, const std::string & datasetPath,
80  std::valarray<T> & v) {
81  // Open dataset
82  isce3::io::IDataSet dataset = h5obj.openDataSet(datasetPath);
83  // Read the valarray dataset
84  dataset.read(v);
85  }
86 
92  template <typename H5obj, typename T>
93  inline void loadFromH5(H5obj & h5obj, const std::string & datasetPath,
95  // Open dataset
96  isce3::io::IDataSet dataset = h5obj.openDataSet(datasetPath);
97  // Read the Matrix dataset using raw pointer interface
98  dataset.read(m.data());
99  }
100 
107  template<typename H5obj, typename T>
108  inline void saveToH5(H5obj& h5obj, const std::string &datasetPath,
109  const T &val, const std::string & units="") {
110  //Check for existence of dataset
111  if (exists(h5obj, datasetPath)) {
112  return;
113  }
114  //Create dataset
115  isce3::io::IDataSet dset = h5obj.createDataSet(datasetPath, val);
116  //Update units attribute if long enough
117  if (units.length() > 0) {
118  dset.createAttribute("units", units);
119  }
120  }
121 
128  template <typename H5obj, typename T>
129  inline void saveToH5(H5obj & h5obj, const std::string & datasetPath,
130  const std::vector<T> & v, const std::string & units = "") {
131  // Check for existence of dataset
132  if (exists(h5obj, datasetPath)) {
133  return;
134  }
135  // Create dataset
136  isce3::io::IDataSet dset = h5obj.createDataSet(datasetPath, v);
137  // Update units attribute if long enough
138  if (units.length() > 0) {
139  dset.createAttribute("units", units);
140  }
141  }
142 
149  template <typename H5obj, typename T>
150  inline void saveToH5(H5obj & h5obj, const std::string & datasetPath,
151  const std::valarray<T> & v, const std::string & units = "") {
152  // Check for existence of dataset
153  if (exists(h5obj, datasetPath)) {
154  return;
155  }
156  // Create dataset
157  isce3::io::IDataSet dset = h5obj.createDataSet(datasetPath, v);
158  // Update units attribute if long enough
159  if (units.length() > 0) {
160  dset.createAttribute("units", units);
161  }
162  }
163 
170  template <typename H5obj, typename T, size_t S>
171  inline void saveToH5(H5obj & h5obj, const std::string & datasetPath,
172  const std::vector<T> & v, std::array<size_t, S> dims,
173  const std::string & units = "") {
174  // Check for existence of dataset
175  if (exists(h5obj, datasetPath)) {
176  return;
177  }
178  // Create dataset
179  isce3::io::IDataSet dset = h5obj.createDataSet(datasetPath, v, dims);
180  // Update units attribute if long enough
181  if (units.length() > 0) {
182  dset.createAttribute("units", units);
183  }
184  }
185 
192  template <typename H5obj, typename T, size_t S>
193  inline void saveToH5(H5obj & h5obj, const std::string & datasetPath,
194  const std::valarray<T> & v, std::array<size_t, S> dims,
195  const std::string & units = "") {
196  // Check for existence of dataset
197  if (exists(h5obj, datasetPath)) {
198  return;
199  }
200  // Create dataset
201  isce3::io::IDataSet dset = h5obj.createDataSet(datasetPath, v, dims);
202  // Update units attribute if long enough
203  if (units.length() > 0) {
204  dset.createAttribute("units", units);
205  }
206  }
207 
214  template <typename H5obj, typename T>
215  inline void saveToH5(H5obj & h5obj, const std::string & datasetPath,
216  const isce3::core::Matrix<T> & mat,
217  const std::string & units = "") {
218  // Check for existence of dataset
219  if (exists(h5obj, datasetPath)) {
220  return;
221  }
222  // Create dataset
223  std::array<size_t, 2> dims{mat.length(), mat.width()};
224  isce3::io::IDataSet dset = h5obj.createDataSet(datasetPath, mat.data(), dims);
225  // Update units attribute if long enough
226  if (units.length() > 0) {
227  dset.createAttribute("units", units);
228  }
229  }
230 
235  template <typename H5obj>
236  inline std::vector<int> getImageDims(H5obj & h5obj, const std::string & datasetPath) {
237  // Open dataset
238  isce3::io::IDataSet dataset = h5obj.openDataSet(datasetPath);
239  // Get dimensions
240  return dataset.getDimensions();
241  }
242 
248  template <typename H5obj>
249  inline isce3::core::DateTime getRefEpoch(H5obj & h5obj, const std::string & datasetPath) {
250 
251  // Open the dataset
252  isce3::io::IDataSet dset = h5obj.openDataSet(datasetPath);
253 
254  // Get the attribute
255  std::string unitAttr;
256  std::string attribute("units");
257  dset.read(unitAttr, attribute);
258 
259  // Parse it
260  std::string dummy1, dummy2, datestr, timestr;
261  std::istringstream ss(unitAttr);
262  ss >> dummy1 >> dummy2 >> datestr >> timestr;
263  std::string inputstr = datestr + "T" + timestr + ".00";
264  isce3::core::DateTime epoch(inputstr);
265 
266  // Done
267  return epoch;
268  }
269 
275  template <typename H5obj>
276  inline void setRefEpoch(H5obj & h5obj, const std::string & datasetPath,
277  const isce3::core::DateTime & refEpoch) {
278 
279  // Open the dataset
280  isce3::io::IDataSet dset = h5obj.openDataSet(datasetPath);
281 
282  // Need to create string representation of DateTime manually
283  char buffer[40];
284  sprintf(buffer,
285  "seconds since %04d-%02d-%02d %02d:%02d:%02d",
286  refEpoch.year,
287  refEpoch.months,
288  refEpoch.days,
289  refEpoch.hours,
290  refEpoch.minutes,
291  refEpoch.seconds);
292  std::string unitsAttr{buffer};
293 
294  // Save buffer to attribute
295  dset.createAttribute("units", unitsAttr);
296  }
297 
298  }
299 }
std::vector< int > getDimensions(const std::string &v="")
Get the size of each dimension of the dataset or given attribute.
Definition: IH5.cpp:162
void loadFromH5(H5obj &h5obj, const std::string &datasetPath, T &v)
Load scalar dataset from HDF5 file.
Definition: Serialization.h:52
void read(T &v, const std::string &att="")
Reading scalar (non string) dataset or attributes.
Definition: IH5.icc:130
Data structure to store date time to nano-sec precision.
Definition: DateTime.h:18
void createAttribute(const std::string &name, const T &data)
Creating and writing a scalar as an attribute.
Definition: IH5.icc:824
std::vector< int > getImageDims(H5obj &h5obj, const std::string &datasetPath)
Get dimensions of complex imagery from HDF5 file.
Definition: Serialization.h:236
cell_t * data()
Access to data buffer.
Definition: Matrix.icc:219
size_t length() const
Get matrix length.
Definition: Matrix.icc:357
Our derived dataset structure that includes utility functions.
Definition: IH5.h:44
void saveToH5(H5obj &h5obj, const std::string &datasetPath, const T &val, const std::string &units="")
Write scalar dataset to HDF5 file.
Definition: Serialization.h:108
bool exists(H5obj &h5obj, const std::string &name, const std::string start=".", const std::string type="BOTH")
Check existence of a dataset or group.
Definition: Serialization.h:34
isce3::core::DateTime getRefEpoch(H5obj &h5obj, const std::string &datasetPath)
Parse time units in a dataset attribute to get a reference epoch.
Definition: Serialization.h:249
void setRefEpoch(H5obj &h5obj, const std::string &datasetPath, const isce3::core::DateTime &refEpoch)
Save reference epoch DateTime as an attribute.
Definition: Serialization.h:276
size_t width() const
Get matrix width.
Definition: Matrix.icc:349

Generated for ISCE3.0 by doxygen 1.8.5.