isce3 0.25.0
Loading...
Searching...
No Matches
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
10
11#pragma once
12
13#include <iostream>
14#include <pyre/journal.h>
15#include <sstream>
16
17// isce3::core
18#include <isce3/core/DateTime.h>
19
20// isce3::io
21#include <isce3/io/IH5.h>
22
24namespace isce3 {
26namespace io {
27
35template<typename H5obj>
36inline bool exists(H5obj& h5obj, const std::string& name,
37 const std::string start = ".",
38 const std::string type = "BOTH")
39{
40 // Search for objects
41 std::vector<std::string> objs = h5obj.find(name, start, type, "FULL");
42 // Check size of vector
43 if (objs.size() > 0) {
44 return true;
45 } else {
46 return false;
47 }
48}
49
55template<typename H5obj, typename T>
56inline void loadFromH5(H5obj& h5obj, const std::string& datasetPath, T& v)
57{
58 // Open dataset
59 isce3::io::IDataSet dataset = h5obj.openDataSet(datasetPath);
60 // Read the scalar dataset
61 dataset.read(v);
62}
63
69template<typename H5obj, typename T>
70inline void loadFromH5(H5obj& h5obj, const std::string& datasetPath,
71 std::vector<T>& v)
72{
73 // Open dataset
74 isce3::io::IDataSet dataset = h5obj.openDataSet(datasetPath);
75 // Read the vector dataset
76 dataset.read(v);
77}
78
84template<typename H5obj, typename T>
85inline void loadFromH5(H5obj& h5obj, const std::string& datasetPath,
86 std::valarray<T>& v)
87{
88 // Open dataset
89 isce3::io::IDataSet dataset = h5obj.openDataSet(datasetPath);
90 // Read the valarray dataset
91 dataset.read(v);
92}
93
99template<typename H5obj, typename T>
100inline void loadFromH5(H5obj& h5obj, const std::string& datasetPath,
102{
103 // Open dataset
104 isce3::io::IDataSet dataset = h5obj.openDataSet(datasetPath);
105 // Read the Matrix dataset using raw pointer interface
106 dataset.read(m.data());
107}
108
115template<typename H5obj, typename T>
116inline void saveToH5(H5obj& h5obj, const std::string& datasetPath, const T& val,
117 const std::string& units = "")
118{
119 // Check for existence of dataset
120 if (exists(h5obj, datasetPath)) {
121 return;
122 }
123 // Create dataset
124 isce3::io::IDataSet dset = h5obj.createDataSet(datasetPath, val);
125 // Update units attribute if long enough
126 if (units.length() > 0) {
127 dset.createAttribute("units", units);
128 }
129}
130
137template<typename H5obj, typename T>
138inline void saveToH5(H5obj& h5obj, const std::string& datasetPath,
139 const std::vector<T>& v, const std::string& units = "")
140{
141 // Check for existence of dataset
142 if (exists(h5obj, datasetPath)) {
143 return;
144 }
145 // Create dataset
146 isce3::io::IDataSet dset = h5obj.createDataSet(datasetPath, v);
147 // Update units attribute if long enough
148 if (units.length() > 0) {
149 dset.createAttribute("units", units);
150 }
151}
152
159template<typename H5obj, typename T>
160inline void saveToH5(H5obj& h5obj, const std::string& datasetPath,
161 const std::valarray<T>& v, const std::string& units = "")
162{
163 // Check for existence of dataset
164 if (exists(h5obj, datasetPath)) {
165 return;
166 }
167 // Create dataset
168 isce3::io::IDataSet dset = h5obj.createDataSet(datasetPath, v);
169 // Update units attribute if long enough
170 if (units.length() > 0) {
171 dset.createAttribute("units", units);
172 }
173}
174
181template<typename H5obj, typename T, size_t S>
182inline void saveToH5(H5obj& h5obj, const std::string& datasetPath,
183 const std::vector<T>& v, std::array<size_t, S> dims,
184 const std::string& units = "")
185{
186 // Check for existence of dataset
187 if (exists(h5obj, datasetPath)) {
188 return;
189 }
190 // Create dataset
191 isce3::io::IDataSet dset = h5obj.createDataSet(datasetPath, v, dims);
192 // Update units attribute if long enough
193 if (units.length() > 0) {
194 dset.createAttribute("units", units);
195 }
196}
197
204template<typename H5obj, typename T, size_t S>
205inline void saveToH5(H5obj& h5obj, const std::string& datasetPath,
206 const std::valarray<T>& v, std::array<size_t, S> dims,
207 const std::string& units = "")
208{
209 // Check for existence of dataset
210 if (exists(h5obj, datasetPath)) {
211 return;
212 }
213 // Create dataset
214 isce3::io::IDataSet dset = h5obj.createDataSet(datasetPath, v, dims);
215 // Update units attribute if long enough
216 if (units.length() > 0) {
217 dset.createAttribute("units", units);
218 }
219}
220
227template<typename H5obj, typename T>
228inline void saveToH5(H5obj& h5obj, const std::string& datasetPath,
229 const isce3::core::Matrix<T>& mat,
230 const std::string& units = "")
231{
232 // Check for existence of dataset
233 if (exists(h5obj, datasetPath)) {
234 return;
235 }
236 // Create dataset
237 std::array<size_t, 2> dims {mat.length(), mat.width()};
239 h5obj.createDataSet(datasetPath, mat.data(), dims);
240 // Update units attribute if long enough
241 if (units.length() > 0) {
242 dset.createAttribute("units", units);
243 }
244}
245
250template<typename H5obj>
251inline std::vector<int> getImageDims(H5obj& h5obj,
252 const std::string& datasetPath)
253{
254 // Open dataset
255 isce3::io::IDataSet dataset = h5obj.openDataSet(datasetPath);
256 // Get dimensions
257 return dataset.getDimensions();
258}
259
267
268template<typename H5obj>
270 const std::string& datasetPath)
271{
272
273 // Open the dataset
274 isce3::io::IDataSet dset = h5obj.openDataSet(datasetPath);
275
276 // Get the attribute
277 std::string unitAttr;
278 std::string attribute("units");
279 dset.read(unitAttr, attribute);
280
281 // Parse it
282 // find UTC date-time string
283 std::string utc_ref {};
284 const auto pos_iso {unitAttr.find_last_of('-')};
285 if (pos_iso != std::string::npos) {
286 utc_ref = unitAttr.substr(pos_iso - 7, 30);
287 // get rid of any trailing white space
288 auto pos_last_wspc = utc_ref.find_last_not_of(" ");
289 if (pos_last_wspc != std::string::npos)
290 utc_ref = utc_ref.substr(0, pos_last_wspc + 1);
291 }
292 isce3::core::DateTime epoch(utc_ref);
293 // Done
294 return epoch;
295}
296
306template<typename H5obj>
307inline void setRefEpoch(H5obj& h5obj, const std::string& datasetPath,
308 const isce3::core::DateTime& refEpoch,
309 const bool ensureEpochIntegerSeconds = true)
310{
311 if (refEpoch.frac != 0.0 && ensureEpochIntegerSeconds) {
312 auto log = pyre::journal::error_t("isce.io.Serialization.setRefEpoch");
313 log << pyre::journal::at(__HERE__) << "Truncating fractional part of "
314 "reference epoch " << std::string(refEpoch) << " while serializing "
315 "dataset " << datasetPath << " to HDF5" << pyre::journal::endl;
316 }
317
318 char buffer[50];
319
320 // Need to create string representation of DateTime manually
321 if (ensureEpochIntegerSeconds) {
322 snprintf(buffer, sizeof(buffer), "seconds since %04d-%02d-%02dT%02d:%02d:%02d",
323 refEpoch.year, refEpoch.months, refEpoch.days, refEpoch.hours,
324 refEpoch.minutes, refEpoch.seconds);
325
326 } else {
327 snprintf(buffer, sizeof(buffer), "seconds since %04d-%02d-%02dT%02d:%02d:%012.9f",
328 refEpoch.year, refEpoch.months, refEpoch.days, refEpoch.hours,
329 refEpoch.minutes, refEpoch.seconds + refEpoch.frac);
330 }
331
332 // Open the dataset
333 isce3::io::IDataSet dset = h5obj.openDataSet(datasetPath);
334
335 std::string unitsAttr {buffer};
336
337 // Save buffer to attribute
338 dset.createAttribute("units", unitsAttr);
339}
340
341} // namespace io
342} // namespace isce3
Data structure to store date time to nano-sec precision.
Definition DateTime.h:18
Data structure for a 2D row-major matrix.
Definition Matrix.h:23
size_t length() const
Get matrix length.
Definition Matrix.h:77
size_t width() const
Get matrix width.
Definition Matrix.h:74
Our derived dataset structure that includes utility functions.
Definition IH5.h:41
std::vector< int > getDimensions(const std::string &v="")
Get the size of each dimension of the dataset or given attribute.
Definition IH5.cpp:172
void read(T &v, const std::string &att="")
Reading scalar (non string) dataset or attributes.
Definition IH5.icc:143
void createAttribute(const std::string &name, const T &data)
Creating and writing a scalar as an attribute.
Definition IH5.icc:837
The isce3::io namespace.
Definition Constants.h:14
void saveToH5(H5obj &h5obj, const std::string &datasetPath, const T &val, const std::string &units="")
Write scalar dataset to HDF5 file.
Definition Serialization.h:116
std::vector< int > getImageDims(H5obj &h5obj, const std::string &datasetPath)
Get dimensions of complex imagery from HDF5 file.
Definition Serialization.h:251
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:36
void setRefEpoch(H5obj &h5obj, const std::string &datasetPath, const isce3::core::DateTime &refEpoch, const bool ensureEpochIntegerSeconds=true)
Save reference epoch DateTime as an attribute.
Definition Serialization.h:307
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:269
void loadFromH5(H5obj &h5obj, const std::string &datasetPath, T &v)
Load scalar dataset from HDF5 file.
Definition Serialization.h:56
base interpolator is an abstract base class
Definition BinarySearchFunc.cpp:5

Generated for ISCE3.0 by doxygen 1.13.2.