14#include <pyre/journal.h>
18#include <isce3/core/DateTime.h>
21#include <isce3/io/IH5.h>
35template<
typename H5obj>
36inline bool exists(H5obj& h5obj,
const std::string& name,
37 const std::string start =
".",
38 const std::string type =
"BOTH")
41 std::vector<std::string> objs = h5obj.find(name, start, type,
"FULL");
43 if (objs.size() > 0) {
55template<
typename H5obj,
typename T>
56inline void loadFromH5(H5obj& h5obj,
const std::string& datasetPath, T& v)
69template<
typename H5obj,
typename T>
70inline void loadFromH5(H5obj& h5obj,
const std::string& datasetPath,
84template<
typename H5obj,
typename T>
85inline void loadFromH5(H5obj& h5obj,
const std::string& datasetPath,
99template<
typename H5obj,
typename T>
100inline void loadFromH5(H5obj& h5obj,
const std::string& datasetPath,
106 dataset.
read(m.data());
115template<
typename H5obj,
typename T>
116inline void saveToH5(H5obj& h5obj,
const std::string& datasetPath,
const T& val,
117 const std::string& units =
"")
120 if (
exists(h5obj, datasetPath)) {
126 if (units.length() > 0) {
137template<
typename H5obj,
typename T>
138inline void saveToH5(H5obj& h5obj,
const std::string& datasetPath,
139 const std::vector<T>& v,
const std::string& units =
"")
142 if (
exists(h5obj, datasetPath)) {
148 if (units.length() > 0) {
159template<
typename H5obj,
typename T>
160inline void saveToH5(H5obj& h5obj,
const std::string& datasetPath,
161 const std::valarray<T>& v,
const std::string& units =
"")
164 if (
exists(h5obj, datasetPath)) {
170 if (units.length() > 0) {
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 =
"")
187 if (
exists(h5obj, datasetPath)) {
193 if (units.length() > 0) {
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 =
"")
210 if (
exists(h5obj, datasetPath)) {
216 if (units.length() > 0) {
227template<
typename H5obj,
typename T>
228inline void saveToH5(H5obj& h5obj,
const std::string& datasetPath,
230 const std::string& units =
"")
233 if (
exists(h5obj, datasetPath)) {
237 std::array<size_t, 2> dims {mat.
length(), mat.
width()};
239 h5obj.createDataSet(datasetPath, mat.data(), dims);
241 if (units.length() > 0) {
250template<
typename H5obj>
252 const std::string& datasetPath)
268template<
typename H5obj>
270 const std::string& datasetPath)
277 std::string unitAttr;
278 std::string attribute(
"units");
279 dset.
read(unitAttr, attribute);
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);
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);
306template<
typename H5obj>
307inline void setRefEpoch(H5obj& h5obj,
const std::string& datasetPath,
309 const bool ensureEpochIntegerSeconds =
true)
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;
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);
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);
335 std::string unitsAttr {buffer};
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