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 <memory>
15 #include <cereal/types/memory.hpp>
16 #include <cereal/types/vector.hpp>
17 #include <cereal/archives/xml.hpp>
18 #include <sstream>
19 #include <stdexcept>
20 #include <vector>
21 
22 // pyre
23 #include <pyre/journal.h>
24 
25 // isce3::core
26 #include <isce3/core/DateTime.h>
27 #include <isce3/core/EulerAngles.h>
28 #include <isce3/core/Ellipsoid.h>
29 #include <isce3/core/Metadata.h>
30 #include <isce3/core/Orbit.h>
31 #include <isce3/core/Poly2d.h>
32 #include <isce3/core/Quaternion.h>
33 #include <isce3/core/LUT1d.h>
34 #include <isce3/core/LUT2d.h>
35 #include <isce3/core/StateVector.h>
36 #include <isce3/core/TimeDelta.h>
37 
38 // isce3::io
39 #include <isce3/io/IH5.h>
40 #include <isce3/io/Serialization.h>
41 
43 namespace isce3 {
45 namespace core {
46 
47 // Archiving any isce3::core object by pointer
48 template <typename T>
49 inline void load_archive(std::string metadata, char * objectTag, T * object)
50 {
51  std::stringstream metastream;
52  metastream << metadata;
53  cereal::XMLInputArchive archive(metastream);
54  archive(cereal::make_nvp(objectTag, (*object)));
55 }
56 
57 // Archiving any isce3::core object by reference
58 template <typename T>
59 inline void load_archive_reference(std::string metadata, char * objectTag, T & object)
60 {
61  std::stringstream metastream;
62  metastream << metadata;
63  cereal::XMLInputArchive archive(metastream);
64  archive(cereal::make_nvp(objectTag, object));
65 }
66 
67 // ------------------------------------------------------------------------
68 // Serialization for Ellipsoid
69 // ------------------------------------------------------------------------
70 
71 template<class Archive>
72 inline void save(Archive & archive, const Ellipsoid & ellps) {
73  archive(cereal::make_nvp("a", ellps.a()), cereal::make_nvp("e2", ellps.e2()));
74 }
75 
76 template<class Archive>
77 inline void load(Archive & archive, Ellipsoid & ellps)
78 {
79  double a, e2;
80  archive(cereal::make_nvp("a", a), cereal::make_nvp("e2", e2));
81  ellps.a(a);
82  ellps.e2(e2);
83 }
84 
91 inline void loadFromH5(isce3::io::IGroup & group, Ellipsoid & ellps)
92 {
93  // Read data
94  std::vector<double> ellpsData;
95  isce3::io::loadFromH5(group, "ellipsoid", ellpsData);
96  // Set ellipsoid properties
97  ellps.a(ellpsData[0]);
98  ellps.e2(ellpsData[1]);
99 }
100 
101 // ------------------------------------------------------------------------
102 // Serialization for Orbit
103 // ------------------------------------------------------------------------
104 
105 template <class Archive>
106 inline void save(Archive & archive, const Orbit & orbit) {
107  archive(cereal::make_nvp("StateVectors", orbit.getStateVectors()),
108  cereal::make_nvp("InterpMethod", orbit.interpMethod()));
109 }
110 
111 template <class Archive>
112 inline void load(Archive & archive, Orbit & orbit)
113 {
114  std::vector<StateVector> statevecs;
115  OrbitInterpMethod interp_method;
116  archive(cereal::make_nvp("StateVectors", statevecs),
117  cereal::make_nvp("InterpMethod", interp_method));
118  orbit.referenceEpoch(statevecs.at(0).datetime);
119  orbit.setStateVectors(statevecs);
120  orbit.interpMethod(interp_method);
121 }
122 
129 inline void loadFromH5(isce3::io::IGroup & group, Orbit & orbit)
130 {
131  // open datasets
132  isce3::io::IDataSet time_ds = group.openDataSet("time");
133  isce3::io::IDataSet pos_ds = group.openDataSet("position");
134  isce3::io::IDataSet vel_ds = group.openDataSet("velocity");
135 
136  // get dataset dimensions
137  std::vector<int> time_dims = time_ds.getDimensions();
138  std::vector<int> pos_dims = pos_ds.getDimensions();
139  std::vector<int> vel_dims = vel_ds.getDimensions();
140 
141  // check validity
142  if (time_dims.size() != 1 || pos_dims.size() != 2 || vel_dims.size() != 2) {
143  throw std::runtime_error("unexpected orbit state vector dims");
144  }
145  if (pos_dims[1] != 3 || vel_dims[1] != 3) {
146  throw std::runtime_error("unexpected orbit position/velocity vector size");
147  }
148  std::size_t size = time_dims[0];
149  if (pos_dims[0] != size || vel_dims[0] != size) {
150  throw std::runtime_error("mismatched orbit state vector component sizes");
151  }
152 
153  // read orbit data
154  std::vector<double> time(size);
155  std::vector<double> pos(3 * size);
156  std::vector<double> vel(3 * size);
157  time_ds.read(time);
158  pos_ds.read(pos);
159  vel_ds.read(vel);
160 
161  // get reference epoch
162  DateTime reference_epoch = isce3::io::getRefEpoch(group, "time");
163 
164  // convert to state vectors
165  std::vector<StateVector> statevecs(size);
166  for (std::size_t i = 0; i < size; ++i) {
167  statevecs[i].datetime = reference_epoch + TimeDelta(time[i]);
168  statevecs[i].position = { pos[3*i+0], pos[3*i+1], pos[3*i+2] };
169  statevecs[i].velocity = { vel[3*i+0], vel[3*i+1], vel[3*i+2] };
170  }
171 
172  // build orbit
173  orbit.referenceEpoch(reference_epoch);
174  orbit.setStateVectors(statevecs);
175 
176  // get interp method
177  std::string interp_method = "Hermite";
178  if (isce3::io::exists(group, "interpMethod")) {
179  isce3::io::loadFromH5(group, "interpMethod", interp_method);
180  }
181 
182  if (interp_method == "Hermite") {
184  }
185  else if (interp_method == "Legendre") {
187  }
188  else {
189  throw std::runtime_error("unexpected orbit interpolation method '" + interp_method + "'");
190  }
191 }
192 
199 inline void saveToH5(isce3::io::IGroup & group, const Orbit & orbit)
200 {
201  // convert times to vector, get flattened position, velocity
202  std::vector<double> time(orbit.size());
203  std::vector<double> pos(3 * orbit.size());
204  std::vector<double> vel(3 * orbit.size());
205 
206  for (int i = 0; i < orbit.size(); ++i) {
207  time[i] = orbit.time(i);
208 
209  pos[3*i+0] = orbit.position(i)[0];
210  pos[3*i+1] = orbit.position(i)[1];
211  pos[3*i+2] = orbit.position(i)[2];
212 
213  vel[3*i+0] = orbit.velocity(i)[0];
214  vel[3*i+1] = orbit.velocity(i)[1];
215  vel[3*i+2] = orbit.velocity(i)[2];
216  }
217 
218  // position/velocity data dims
219  std::size_t size = orbit.size();
220  std::array<std::size_t, 2> dims = {size, 3};
221 
222  // interp method
223  std::string interp_method;
224  if (orbit.interpMethod() == OrbitInterpMethod::Hermite) {
225  interp_method = "Hermite";
226  }
227  else if (orbit.interpMethod() == OrbitInterpMethod::Legendre) {
228  interp_method = "Legendre";
229  }
230  else {
231  throw std::runtime_error("unexpected orbit interpolation method");
232  }
233 
234  // serialize
235  isce3::io::saveToH5(group, "time", time);
236  isce3::io::setRefEpoch(group, "time", orbit.referenceEpoch());
237  isce3::io::saveToH5(group, "position", pos, dims, "meters");
238  isce3::io::saveToH5(group, "velocity", vel, dims, "meters per second");
239  isce3::io::saveToH5(group, "interpMethod", interp_method);
240 }
241 
242 // ------------------------------------------------------------------------
243 // Serialization for EulerAngles
244 // ------------------------------------------------------------------------
245 
252 inline void loadFromH5(isce3::io::IGroup & group, EulerAngles & euler)
253 {
254  // Create temporary data
255  std::vector<double> time, angles, yaw, pitch, roll;
256  isce3::core::DateTime refEpoch;
257 
258  // Load angles
259  isce3::io::loadFromH5(group, "eulerAngles", angles);
260 
261  // Load time
262  isce3::io::loadFromH5(group, "time", time);
263 
264  // Get the reference epoch
265  refEpoch = isce3::io::getRefEpoch(group, "time");
266 
267  // Unpack the angles
268  const double rad = M_PI / 180.0;
269  yaw.resize(time.size());
270  pitch.resize(time.size());
271  roll.resize(time.size());
272  for (size_t i = 0; i < time.size(); ++i) {
273  yaw[i] = rad * angles[i*3 + 0];
274  pitch[i] = rad * angles[i*3 + 1];
275  roll[i] = rad * angles[i*3 + 2];
276  }
277 
278  // Save to EulerAngles object
279  euler.data(time, yaw, pitch, roll);
280  euler.refEpoch(refEpoch);
281 }
282 
289 inline void saveToH5(isce3::io::IGroup & group, const EulerAngles & euler)
290 {
291  // Create vector to store all data (convert angles to degrees)
292  const double deg = 180.0 / M_PI;
293  std::vector<double> angles(euler.nVectors() * 3);
294  for (size_t i = 0; i < euler.nVectors(); ++i) {
295  angles[i*3 + 0] = deg * euler.yaw()[i];
296  angles[i*3 + 1] = deg * euler.pitch()[i];
297  angles[i*3 + 2] = deg * euler.roll()[i];
298  }
299 
300  // Save angles
301  std::array<size_t, 2> dims = {euler.nVectors(), 3};
302  isce3::io::saveToH5(group, "eulerAngles", angles, dims, "degrees");
303 
304  // Save time and reference epoch attribute
305  isce3::io::saveToH5(group, "time", euler.time());
306  isce3::io::setRefEpoch(group, "time", euler.refEpoch());
307 }
308 
309 // ------------------------------------------------------------------------
310 // Serialization for Quaternion
311 // ------------------------------------------------------------------------
312 
313 inline void loadFromH5(isce3::io::IGroup & group, Quaternion & qs)
314 {
315  std::vector<double> time, qvec;
316  isce3::io::loadFromH5(group, "time", time);
317  isce3::io::loadFromH5(group, "quaternions", qvec);
318  qs.data(time, qvec);
319 }
320 
321 inline void saveToH5(isce3::io::IGroup & group, const Quaternion & qs)
322 {
323  std::array<size_t, 2> dims = {qs.nVectors(), 4};
324  isce3::io::saveToH5(group, "time", qs.time());
325  isce3::io::saveToH5(group, "quaternions", qs.qvec(), dims);
326 }
327 
328 // ------------------------------------------------------------------------
329 // Serialization for Metadata
330 // ------------------------------------------------------------------------
331 
332 template <class Archive>
333 inline void save(Archive & archive, const Metadata & meta)
334 {
335  archive(cereal::make_nvp("width", meta.width),
336  cereal::make_nvp("length", meta.length),
337  cereal::make_nvp("numberRangeLooks", meta.numberRangeLooks),
338  cereal::make_nvp("numberAzimuthLooks", meta.numberAzimuthLooks),
339  cereal::make_nvp("slantRangePixelSpacing", meta.slantRangePixelSpacing),
340  cereal::make_nvp("rangeFirstSample", meta.rangeFirstSample),
341  cereal::make_nvp("lookSide", meta.lookSide),
342  cereal::make_nvp("prf", meta.prf),
343  cereal::make_nvp("radarWavelength", meta.radarWavelength),
344  cereal::make_nvp("pegHeading", meta.pegHeading),
345  cereal::make_nvp("pegLatitude", meta.pegLatitude),
346  cereal::make_nvp("pegLongitude", meta.pegLongitude),
347  cereal::make_nvp("chirpSlope", meta.chirpSlope),
348  cereal::make_nvp("antennaLength", meta.antennaLength),
349  cereal::make_nvp("sensingStart", meta.sensingStart.isoformat()));
350 }
351 
352 template <class Archive>
353 inline void load(Archive & archive, Metadata & meta)
354 {
355  std::string sensingStart;
356  archive(cereal::make_nvp("width", meta.width),
357  cereal::make_nvp("length", meta.length),
358  cereal::make_nvp("numberRangeLooks", meta.numberRangeLooks),
359  cereal::make_nvp("numberAzimuthLooks", meta.numberAzimuthLooks),
360  cereal::make_nvp("slantRangePixelSpacing", meta.slantRangePixelSpacing),
361  cereal::make_nvp("rangeFirstSample", meta.rangeFirstSample),
362  cereal::make_nvp("lookSide", meta.lookSide),
363  cereal::make_nvp("prf", meta.prf),
364  cereal::make_nvp("radarWavelength", meta.radarWavelength),
365  cereal::make_nvp("pegHeading", meta.pegHeading),
366  cereal::make_nvp("pegLatitude", meta.pegLatitude),
367  cereal::make_nvp("pegLongitude", meta.pegLongitude),
368  cereal::make_nvp("chirpSlope", meta.chirpSlope),
369  cereal::make_nvp("antennaLength", meta.antennaLength),
370  cereal::make_nvp("sensingStart", sensingStart));
371  meta.sensingStart = sensingStart;
372 }
373 
374 // ------------------------------------------------------------------------
375 // Serialization for Poly2d
376 // ------------------------------------------------------------------------
377 
378 // Definition for Poly2d
379 template <class Archive>
380 inline void serialize(Archive & archive, Poly2d & poly)
381 {
382  archive(cereal::make_nvp("rangeOrder", poly.rangeOrder),
383  cereal::make_nvp("azimuthOrder", poly.azimuthOrder),
384  cereal::make_nvp("rangeMean", poly.rangeMean),
385  cereal::make_nvp("azimuthMean", poly.azimuthMean),
386  cereal::make_nvp("rangeNorm", poly.rangeNorm),
387  cereal::make_nvp("azimuthNorm", poly.azimuthNorm),
388  cereal::make_nvp("coeffs", poly.coeffs));
389 }
390 
398 inline void loadFromH5(isce3::io::IGroup & group, Poly2d & poly, std::string name)
399 {
400  // Configure the polynomial coefficients
401  isce3::io::loadFromH5(group, name, poly.coeffs);
402 
403  // Set other polynomial properties
404  poly.rangeOrder = poly.coeffs.size() - 1;
405  poly.azimuthOrder = 0;
406  poly.rangeMean = 0.0;
407  poly.azimuthMean = 0.0;
408  poly.rangeNorm = 1.0;
409  poly.azimuthNorm = 1.0;
410 }
411 
412 // ------------------------------------------------------------------------
413 // Serialization for LUT2d (specifically for calibration grids)
414 // ------------------------------------------------------------------------
415 
423 template <typename T>
424 inline void loadCalGrid(isce3::io::IGroup & group, const std::string & dsetName,
425  isce3::core::LUT2d<T> & lut)
426 {
427  // Load coordinates
428  std::valarray<double> slantRange, zeroDopplerTime;
429  isce3::io::loadFromH5(group, "slantRange", slantRange);
430  isce3::io::loadFromH5(group, "zeroDopplerTime", zeroDopplerTime);
431 
432  // Load LUT2d data in matrix
433  isce3::core::Matrix<T> matrix(zeroDopplerTime.size(), slantRange.size());
434  isce3::io::loadFromH5(group, dsetName, matrix);
435 
436  // Set in lut
437  lut.setFromData(slantRange, zeroDopplerTime, matrix);
438 }
439 
448 template <typename T>
449 inline void saveCalGrid(isce3::io::IGroup & group,
450  const std::string & dsetName,
451  const isce3::core::LUT2d<T> & lut,
452  const isce3::core::DateTime & refEpoch,
453  const std::string & units = "")
454 {
455  // Generate uniformly spaced X (slant range) coordinates
456  const double x0 = lut.xStart();
457  const double x1 = x0 + lut.xSpacing() * (lut.width() - 1.0);
458  const std::vector<double> x = isce3::core::linspace(x0, x1, lut.width());
459 
460  // Generate uniformly spaced Y (zero Doppler time) coordinates
461  const double y0 = lut.yStart();
462  const double y1 = y0 + lut.ySpacing() * (lut.length() - 1.0);
463  const std::vector<double> y = isce3::core::linspace(y0, y1, lut.length());
464 
465  // Save coordinates
466  isce3::io::saveToH5(group, "slantRange", x, "meters");
467  isce3::io::saveToH5(group, "zeroDopplerTime", y);
468  isce3::io::setRefEpoch(group, "zeroDopplerTime", refEpoch);
469 
470  // Save LUT2d data
471  isce3::io::saveToH5(group, dsetName, lut.data(), units);
472 }
473 
474 // ------------------------------------------------------------------------
475 // Serialization for LUT1d
476 // ------------------------------------------------------------------------
477 
478 // Serialization save method
479 template <class Archive, typename T>
480 inline void save(Archive & archive, LUT1d<T> const & lut)
481 {
482  // Copy LUT data from valarrays to vectors
483  std::vector<double> coords(lut.size());
484  std::vector<T> values(lut.size());
485  auto v_coords = lut.coords();
486  auto v_values = lut.values();
487  coords.assign(std::begin(v_coords), std::end(v_coords));
488  values.assign(std::begin(v_values), std::end(v_values));
489  // Archive
490  archive(cereal::make_nvp("Coords", coords),
491  cereal::make_nvp("Values", values));
492 }
493 
494 // Serialization load method
495 template<class Archive, typename T>
496 inline void load(Archive & archive, LUT1d<T> & lut)
497 {
498  // Create vector for loading results
499  std::vector<double> coords;
500  std::vector<T> values;
501  // Load the archive
502  archive(cereal::make_nvp("Coords", coords),
503  cereal::make_nvp("Values", values));
504  // Copy vector to LUT valarrays
505  std::valarray<double> v_coords(coords.data(), coords.size());
506  std::valarray<T> v_values(values.data(), values.size());
507  lut.coords(v_coords);
508  lut.values(v_values);
509 }
510 
518 template <typename T>
519 inline void loadFromH5(isce3::io::IGroup & group, LUT1d<T> & lut,
520  std::string name_coords, std::string name_values)
521 {
522  // Valarrays for storing results
523  std::valarray<double> x, y;
524  // Load the LUT values
525  isce3::io::loadFromH5(group, name_values, y);
526  // Load the LUT coordinates
527  isce3::io::loadFromH5(group, name_coords, x);
528  // Set LUT data
529  lut.coords(x);
530  lut.values(y);
531 }
532 
533 // ------------------------------------------------------------------------
534 // Serialization for StateVector
535 // ------------------------------------------------------------------------
536 
537 // Serialization save method
538 template<class Archive>
539 inline void save(Archive & archive, const StateVector & sv)
540 {
541  // Serialize position vector to string as whitespace-delimited values
542  std::stringstream pos_stream;
543  pos_stream << sv.position[0] << " " << sv.position[1] << " " << sv.position[2];
544  std::string position_string = pos_stream.str();
545  // Serialize velocity vector to string as whitespace-delimited values
546  std::stringstream vel_stream;
547  vel_stream << sv.velocity[0] << " " << sv.velocity[1] << " " << sv.velocity[2];
548  std::string velocity_string = vel_stream.str();
549  // Archive
550  archive(cereal::make_nvp("Time", sv.datetime.isoformat()),
551  cereal::make_nvp("Position", position_string),
552  cereal::make_nvp("Velocity", velocity_string));
553 }
554 
555 // Serialization load method
556 template<class Archive>
557 inline void load(Archive & archive, StateVector & sv)
558 {
559  // Make strings for position, velocity, and datetime
560  std::string position_string, velocity_string, datetime_string;
561  // Load the archive
562  archive(cereal::make_nvp("Time", datetime_string),
563  cereal::make_nvp("Position", position_string),
564  cereal::make_nvp("Velocity", velocity_string));
565  // Send datetime string to datetime object parser
566  sv.datetime = datetime_string;
567  // De-serialize position vector from stringstream
568  std::stringstream pos_stream (position_string);
569  pos_stream >> sv.position[0] >> sv.position[1] >> sv.position[2];
570  // De-serialize velocity vector from stringstream
571  std::stringstream vel_stream (velocity_string);
572  vel_stream >> sv.velocity[0] >> sv.velocity[1] >> sv.velocity[2];
573 }
574 
575 }
576 }
Eighth-order Legendre polynomial interpolation.
void loadCalGrid(isce3::io::IGroup &group, const std::string &dsetName, isce3::core::LUT2d< T > &lut)
Load LUT2d data from HDF5 product.
Definition: Serialization.h:424
Data structure for representing 1D polynomials.
Definition: Poly2d.h:25
std::valarray< T > & values()
Get a reference to the coordinates.
Definition: LUT1d.h:89
std::vector< int > getDimensions(const std::string &v="")
Get the size of each dimension of the dataset or given attribute.
Definition: IH5.cpp:162
CUDA_HOSTDEV double e2() const
Return eccentricity^2.
Definition: Ellipsoid.h:52
Data structure to store Ellipsoid information.
Definition: Ellipsoid.h:20
IDataSet openDataSet(const H5std_string &name)
Open a given dataset.
Definition: IH5.cpp:854
const std::vector< double > & qvec() const
Get the quaternion elements (packed in size N*4 vector)
Definition: Quaternion.h:54
void loadFromH5(H5obj &h5obj, const std::string &datasetPath, T &v)
Load scalar dataset from HDF5 file.
Definition: Serialization.h:52
const DateTime & referenceEpoch() const
Reference epoch (UTC)
Definition: Orbit.h:77
int rangeOrder
Order of polynomial in range or x.
Definition: Poly2d.h:29
void read(T &v, const std::string &att="")
Reading scalar (non string) dataset or attributes.
Definition: IH5.icc:130
const std::vector< double > & roll() const
Return data vector of roll.
Definition: EulerAngles.h:55
int azimuthOrder
Order of polynomial in azimuth or y.
Definition: Poly2d.h:31
Definition: StateVector.h:10
Data structure to store date time to nano-sec precision.
Definition: DateTime.h:18
OrbitInterpMethod
Orbit interpolation method.
Definition: Orbit.h:20
std::vector< T > linspace(T low, T high, std::size_t number)
Function to return a Matlab/Python style linspace vector.
Definition: Utilities.h:81
CUDA_HOSTDEV double a() const
Return semi-major axis.
Definition: Ellipsoid.h:44
const std::vector< Vec3 > & position() const
Get state vector positions in ECEF coordinates (m)
Definition: Orbit.h:116
double azimuthNorm
Norm in azimuth or y direction.
Definition: Poly2d.h:39
std::vector< double > coeffs
Linearized vector of coefficients in row-major format.
Definition: Poly2d.h:41
void saveToH5(isce3::io::IGroup &group, const Orbit &orbit)
Save orbit data to HDF5 product.
Definition: Serialization.h:199
size_t size() const
Get size of LUT.
Definition: LUT1d.h:115
double rangeNorm
Norm in range or x direction.
Definition: Poly2d.h:37
Data structure to hold a 1D Lookup table.
Definition: forward.h:29
OrbitInterpMethod interpMethod() const
Interpolation method.
Definition: Orbit.h:83
const std::vector< double > & time() const
Return data vector of time.
Definition: Quaternion.h:51
Serialization utilities using HDF5 API.
double azimuthMean
Mean in azimuth or y direction.
Definition: Poly2d.h:35
Data structure to store TimeDelta to double precision seconds.
Definition: TimeDelta.h:16
Data structure to store 2D Lookup table.
Definition: forward.h:30
Data structure for storing basic radar geometry image metadata.
Definition: Metadata.h:17
void setFromData(const std::valarray< double > &xcoord, const std::valarray< double > &ycoord, const isce3::core::Matrix< T > &data)
Definition: LUT2d.cpp:52
void data(const std::vector< double > &time, const std::vector< double > &quaternions)
Set all quaternion elements from a vector.
Definition: Quaternion.cpp:156
std::string isoformat() const
Return date formatted as ISO-8601 string.
Definition: DateTime.cpp:452
const Linspace< double > & time() const
Get state vector times relative to reference epoch (s)
Definition: Orbit.h:113
Our derived dataset structure that includes utility functions.
Definition: IH5.h:44
Quaternion representation of attitude information.
Definition: Quaternion.h:17
int size() const
Number of state vectors in orbit.
Definition: Orbit.h:110
Definition: IH5.h:319
void data(const std::vector< double > &time, const std::vector< double > &yaw, const std::vector< double > &pitch, const std::vector< double > &roll)
Set data after construction.
Definition: EulerAngles.cpp:105
size_t nVectors() const
Return number of epochs.
Definition: Quaternion.h:57
const isce3::core::DateTime & refEpoch() const
Get reference epoch.
Definition: EulerAngles.h:84
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
size_t nVectors() const
Return number of epochs.
Definition: EulerAngles.h:89
double rangeMean
Mean in range or x direction.
Definition: Poly2d.h:33
Third-order Hermite polynomial interpolation.
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
void setStateVectors(const std::vector< StateVector > &)
Set orbit state vectors.
Definition: Orbit.cpp:42
const std::vector< double > & yaw() const
Return data vector of yaw.
Definition: EulerAngles.h:49
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 loadFromH5(isce3::io::IGroup &group, Ellipsoid &ellps)
Load Ellipsoid parameters from HDF5.
Definition: Serialization.h:91
void setRefEpoch(H5obj &h5obj, const std::string &datasetPath, const isce3::core::DateTime &refEpoch)
Save reference epoch DateTime as an attribute.
Definition: Serialization.h:276
std::valarray< double > & coords()
Get a reference to the coordinates.
Definition: LUT1d.h:77
const std::vector< Vec3 > & velocity() const
Get state vector velocities in ECEF coordinates (m/s)
Definition: Orbit.h:119
Data structure for Euler Angle representation of attitude information.
Definition: EulerAngles.h:17
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="")
Save LUT2d data to HDF5 product.
Definition: Serialization.h:449
const std::vector< double > & time() const
Return data vector of time.
Definition: EulerAngles.h:46
Sequence of platform ephemeris samples (state vectors) with uniform temporal spacing, supporting efficient lookup and interpolation.
Definition: Orbit.h:43
std::vector< StateVector > getStateVectors() const
Export list of state vectors.
Definition: Orbit.cpp:31
const std::vector< double > & pitch() const
Return data vector of pitch.
Definition: EulerAngles.h:52

Generated for ISCE3.0 by doxygen 1.8.5.