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 <memory>
15#include <sstream>
16#include <stdexcept>
17#include <vector>
18
19#include <isce3/core/Attitude.h>
20#include <isce3/core/DateTime.h>
21#include <isce3/core/Ellipsoid.h>
22#include <isce3/core/Metadata.h>
23#include <isce3/core/Orbit.h>
24#include <isce3/core/Poly2d.h>
25#include <isce3/core/LUT1d.h>
26#include <isce3/core/LUT2d.h>
27#include <isce3/core/StateVector.h>
28#include <isce3/core/TimeDelta.h>
29#include <isce3/io/IH5.h>
31
32namespace isce3 {
33namespace core {
34
41inline void loadFromH5(isce3::io::IGroup & group, Ellipsoid & ellps)
42{
43 // Read data
44 std::vector<double> ellpsData;
45 isce3::io::loadFromH5(group, "ellipsoid", ellpsData);
46 // Set ellipsoid properties
47 ellps.a(ellpsData[0]);
48 ellps.e2(ellpsData[1]);
49}
50
57inline void loadFromH5(isce3::io::IGroup & group, Orbit & orbit)
58{
59 // open datasets
60 isce3::io::IDataSet time_ds = group.openDataSet("time");
61 isce3::io::IDataSet pos_ds = group.openDataSet("position");
62 isce3::io::IDataSet vel_ds = group.openDataSet("velocity");
63
64 // get dataset dimensions
65 std::vector<int> time_dims = time_ds.getDimensions();
66 std::vector<int> pos_dims = pos_ds.getDimensions();
67 std::vector<int> vel_dims = vel_ds.getDimensions();
68
69 // check validity
70 if (time_dims.size() != 1 || pos_dims.size() != 2 || vel_dims.size() != 2) {
71 throw std::runtime_error("unexpected orbit state vector dims");
72 }
73 if (pos_dims[1] != 3 || vel_dims[1] != 3) {
74 throw std::runtime_error("unexpected orbit position/velocity vector size");
75 }
76 std::size_t size = time_dims[0];
77 if (pos_dims[0] != size || vel_dims[0] != size) {
78 throw std::runtime_error("mismatched orbit state vector component sizes");
79 }
80
81 // read orbit data
82 std::vector<double> time(size);
83 std::vector<double> pos(3 * size);
84 std::vector<double> vel(3 * size);
85 time_ds.read(time);
86 pos_ds.read(pos);
87 vel_ds.read(vel);
88
89 // get reference epoch
90 DateTime reference_epoch = isce3::io::getRefEpoch(group, "time");
91
92 // convert to state vectors
93 std::vector<StateVector> statevecs(size);
94 for (std::size_t i = 0; i < size; ++i) {
95 statevecs[i].datetime = reference_epoch + TimeDelta(time[i]);
96 statevecs[i].position = { pos[3*i+0], pos[3*i+1], pos[3*i+2] };
97 statevecs[i].velocity = { vel[3*i+0], vel[3*i+1], vel[3*i+2] };
98 }
99
100 // build orbit
101 orbit.referenceEpoch(reference_epoch);
102 orbit.setStateVectors(statevecs);
103
104 // get interp method
105 std::string interp_method = "Hermite";
106 if (isce3::io::exists(group, "interpMethod")) {
107 isce3::io::loadFromH5(group, "interpMethod", interp_method);
108 }
109
110 if (interp_method == "Hermite") {
111 orbit.interpMethod(OrbitInterpMethod::Hermite);
112 }
113 else if (interp_method == "Legendre") {
114 orbit.interpMethod(OrbitInterpMethod::Legendre);
115 }
116 else {
117 throw std::runtime_error("unexpected orbit interpolation method '" + interp_method + "'");
118 }
119
120 // get orbit ephemeris precision type
121 std::string orbit_type = "Undefined";
122 if (isce3::io::exists(group, "orbitType")) {
123 isce3::io::loadFromH5(group, "orbitType", orbit_type);
124 }
125 orbit.type(orbit_type);
126
127}
128
138inline void saveToH5(isce3::io::IGroup & group, const Orbit & orbit,
139 const bool ensureEpochIntegerSeconds = true)
140{
141 // convert times to vector, get flattened position, velocity
142 std::vector<double> time(orbit.size());
143 std::vector<double> pos(3 * orbit.size());
144 std::vector<double> vel(3 * orbit.size());
145
146 for (int i = 0; i < orbit.size(); ++i) {
147 time[i] = orbit.time(i);
148
149 pos[3*i+0] = orbit.position(i)[0];
150 pos[3*i+1] = orbit.position(i)[1];
151 pos[3*i+2] = orbit.position(i)[2];
152
153 vel[3*i+0] = orbit.velocity(i)[0];
154 vel[3*i+1] = orbit.velocity(i)[1];
155 vel[3*i+2] = orbit.velocity(i)[2];
156 }
157
158 // position/velocity data dims
159 std::size_t size = orbit.size();
160 std::array<std::size_t, 2> dims = {size, 3};
161
162 // interp method
163 std::string interp_method;
164 if (orbit.interpMethod() == OrbitInterpMethod::Hermite) {
165 interp_method = "Hermite";
166 }
167 else if (orbit.interpMethod() == OrbitInterpMethod::Legendre) {
168 interp_method = "Legendre";
169 }
170 else {
171 throw std::runtime_error("unexpected orbit interpolation method");
172 }
173
174 std::string orbit_type = orbit.type();
175 if (orbit_type == "") {
176 orbit_type = "Undefined";
177 }
178
179 // serialize
180 isce3::io::saveToH5(group, "time", time);
181 isce3::io::setRefEpoch(group, "time", orbit.referenceEpoch(),
182 ensureEpochIntegerSeconds);
183 isce3::io::saveToH5(group, "position", pos, dims, "meters");
184 isce3::io::saveToH5(group, "velocity", vel, dims, "meters / second");
185 isce3::io::saveToH5(group, "orbitType", orbit_type);
186 isce3::io::saveToH5(group, "interpMethod", interp_method);
187}
188
214
215// ------------------------------------------------------------------------
216// Serialization for Attitude
217// ------------------------------------------------------------------------
218
219inline void loadFromH5(isce3::io::IGroup& group, Attitude& att)
220{
221 // load data from file
222 DateTime epoch = isce3::io::getRefEpoch(group, "time");
223 std::vector<double> time, packed_quat;
224 isce3::io::loadFromH5(group, "time", time);
225 isce3::io::loadFromH5(group, "quaternions", packed_quat);
226 // convert quaternion representation
227 int n = packed_quat.size() / 4;
228 std::vector<isce3::core::Quaternion> quat(n);
229 for (int i = 0; i < n; ++i) {
230 const double* q = &packed_quat[i * 4];
231 // XXX Careful to use 4-arg ctor (not array) to avoid mixing up order of
232 // elements. Want real part followed by bivector part.
233 quat[i] = isce3::core::Quaternion(q[0], q[1], q[2], q[3]);
234 }
235 // use ctor for remaining checks
236 att = isce3::core::Attitude(time, quat, epoch);
237}
238
248inline void saveToH5(isce3::io::IGroup & group, const Attitude& att,
249 const bool ensureEpochIntegerSeconds = true)
250{
251 // Flatten quaternion vector.
252 int n = att.size();
253 std::vector<double> qflat(n * 4);
254 auto qvec = att.quaternions();
255 for (int i = 0; i < n; ++i) {
256 double *q = &qflat[i * 4];
257 // XXX Don't use internal Quaternion array since it uses a
258 // different storage order!
259 q[0] = qvec[i].w();
260 q[1] = qvec[i].x();
261 q[2] = qvec[i].y();
262 q[3] = qvec[i].z();
263 }
264 // Save to disk.
265 isce3::io::saveToH5(group, "time", att.time());
266 isce3::io::setRefEpoch(group, "time", att.referenceEpoch(),
267 ensureEpochIntegerSeconds);
268 std::array<size_t, 2> dims = {static_cast<size_t>(n), 4};
269 isce3::io::saveToH5(group, "quaternions", qflat, dims);
270 // TODO convert and save EulerAngles
271}
272
280inline void loadFromH5(isce3::io::IGroup & group, Poly2d & poly, std::string name)
281{
282 // Configure the polynomial coefficients
283 isce3::io::loadFromH5(group, name, poly.coeffs);
284
285 // Set other polynomial properties
286 poly.xOrder = poly.coeffs.size() - 1;
287 poly.yOrder = 0;
288 poly.xMean = 0.0;
289 poly.yMean = 0.0;
290 poly.xNorm = 1.0;
291 poly.yNorm = 1.0;
292}
293
294// ------------------------------------------------------------------------
295// Serialization for LUT2d (specifically for calibration grids)
296// ------------------------------------------------------------------------
297
305template <typename T>
306inline void loadCalGrid(isce3::io::IGroup & group, const std::string & dsetName,
308{
309 // Load coordinates
310 std::valarray<double> slantRange, zeroDopplerTime;
311 isce3::io::loadFromH5(group, "slantRange", slantRange);
312 isce3::io::loadFromH5(group, "zeroDopplerTime", zeroDopplerTime);
313
314 // Load LUT2d data in matrix
315 isce3::core::Matrix<T> matrix(zeroDopplerTime.size(), slantRange.size());
316 isce3::io::loadFromH5(group, dsetName, matrix);
317
318 // Set in lut
319 lut.setFromData(slantRange, zeroDopplerTime, matrix);
320}
321
333template <typename T>
334inline void saveCalGrid(isce3::io::IGroup & group,
335 const std::string & dsetName,
336 const isce3::core::LUT2d<T> & lut,
337 const isce3::core::DateTime & refEpoch,
338 const std::string & units = "",
339 const bool ensureEpochIntegerSeconds = true)
340{
341 // Generate uniformly spaced X (slant range) coordinates
342 const double x0 = lut.xStart();
343 const double x1 = lut.xEnd();
344 const std::vector<double> x = isce3::core::linspace(x0, x1, lut.width());
345
346 // Generate uniformly spaced Y (zero Doppler time) coordinates
347 const double y0 = lut.yStart();
348 const double y1 = lut.yEnd();
349 const std::vector<double> y = isce3::core::linspace(y0, y1, lut.length());
350
351 // Save coordinates
352 isce3::io::saveToH5(group, "slantRange", x, "meters");
353 isce3::io::saveToH5(group, "zeroDopplerTime", y);
354 isce3::io::setRefEpoch(group, "zeroDopplerTime", refEpoch,
355 ensureEpochIntegerSeconds);
356
357 // Save LUT2d data
358 isce3::io::saveToH5(group, dsetName, lut.data(), units);
359}
360
368template <typename T>
369inline void loadFromH5(isce3::io::IGroup & group, LUT1d<T> & lut,
370 std::string name_coords, std::string name_values)
371{
372 // Valarrays for storing results
373 std::valarray<double> x, y;
374 // Load the LUT values
375 isce3::io::loadFromH5(group, name_values, y);
376 // Load the LUT coordinates
377 isce3::io::loadFromH5(group, name_coords, x);
378 // Set LUT data
379 lut.coords(x);
380 lut.values(y);
381}
382
383}} // namespace isce3::core
Store and interpolate attitude measurements.
Definition Attitude.h:17
const DateTime & referenceEpoch() const
Get reference epoch (UTC) for time tags.
Definition Attitude.h:46
int size() const
Return number of epochs.
Definition Attitude.h:43
const std::vector< double > & time() const
Return data vector of time.
Definition Attitude.h:37
const std::vector< Quaternion > & quaternions() const
Return data vector of quaternions.
Definition Attitude.h:40
Data structure to store date time to nano-sec precision.
Definition DateTime.h:18
Data structure to store Ellipsoid information.
Definition Ellipsoid.h:20
CUDA_HOSTDEV double a() const
Return semi-major axis.
Definition Ellipsoid.h:38
CUDA_HOSTDEV double e2() const
Return eccentricity^2.
Definition Ellipsoid.h:46
Data structure to hold a 1D Lookup table.
Definition LUT1d.h:15
std::valarray< T > & values()
Get a reference to the coordinates.
Definition LUT1d.h:82
std::valarray< double > & coords()
Get a reference to the coordinates.
Definition LUT1d.h:70
Data structure to store 2D Lookup table.
Definition LUT2d.h:20
void setFromData(const std::valarray< double > &xcoord, const std::valarray< double > &ycoord, const isce3::core::Matrix< T > &data)
Definition LUT2d.cpp:53
Data structure for a 2D row-major matrix.
Definition Matrix.h:23
Sequence of platform ephemeris samples (state vectors) with uniform temporal spacing,...
Definition Orbit.h:44
const Linspace< double > & time() const
Get state vector times relative to reference epoch (s)
Definition Orbit.h:157
const std::vector< Vec3 > & velocity() const
Get state vector velocities in ECEF coordinates (m/s)
Definition Orbit.h:163
const std::vector< Vec3 > & position() const
Get state vector positions in ECEF coordinates (m)
Definition Orbit.h:160
const DateTime & referenceEpoch() const
Reference epoch (UTC)
Definition Orbit.h:105
void setStateVectors(const std::vector< StateVector > &)
Set orbit state vectors.
Definition Orbit.cpp:52
OrbitInterpMethod interpMethod() const
Interpolation method.
Definition Orbit.h:116
const std::string & type() const
Orbit ephemeris precision type.
Definition Orbit.h:122
int size() const
Number of state vectors in orbit.
Definition Orbit.h:154
Data structure for representing 2D polynomials.
Definition Poly2d.h:20
int yOrder
Order of polynomial in azimuth or y.
Definition Poly2d.h:26
double yNorm
Norm in azimuth or y direction.
Definition Poly2d.h:34
double yMean
Mean in azimuth or y direction.
Definition Poly2d.h:30
double xNorm
Norm in range or x direction.
Definition Poly2d.h:32
int xOrder
Order of polynomial in range or x.
Definition Poly2d.h:24
std::vector< double > coeffs
Linearized vector of coefficients in row-major format.
Definition Poly2d.h:36
double xMean
Mean in range or x direction.
Definition Poly2d.h:28
Quaternion representation of rotations, based on double precision Eigen::Quaterniond.
Definition Quaternion.h:25
Data structure to store TimeDelta to double precision seconds.
Definition TimeDelta.h:16
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
Definition IH5.h:304
IDataSet openDataSet(const H5std_string &name)
Open a given dataset.
Definition IH5.cpp:712
void loadCalGrid(isce3::io::IGroup &group, const std::string &dsetName, isce3::core::LUT2d< T > &lut)
Load LUT2d data from HDF5 product.
Definition Serialization.h:306
void saveToH5(isce3::io::IGroup &group, const Orbit &orbit, const bool ensureEpochIntegerSeconds=true)
Save orbit data to HDF5 product.
Definition Serialization.h:138
void loadFromH5(isce3::io::IGroup &group, Ellipsoid &ellps)
Load Ellipsoid parameters from HDF5.
Definition Serialization.h:41
void saveCalGrid(isce3::io::IGroup &group, const std::string &dsetName, const isce3::core::LUT2d< T > &lut, const isce3::core::DateTime &refEpoch, const std::string &units="", const bool ensureEpochIntegerSeconds=true)
Save LUT2d data to HDF5 product.
Definition Serialization.h:334
Serialization utilities using HDF5 API.
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
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.