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-2019
6 
11 #pragma once
12 
13 // isce3::core
14 #include <isce3/core/Constants.h>
15 #include <isce3/core/Ellipsoid.h>
17 
18 // isce3::io
19 #include <isce3/io/IH5.h>
20 #include <isce3/io/Serialization.h>
21 
22 // isce3::product
23 #include <isce3/product/Metadata.h>
24 #include <isce3/product/Swath.h>
25 
27 namespace isce3 {
29  namespace product {
30 
35  inline void loadFromH5(isce3::io::IGroup & group, ProcessingInformation & proc) {
36 
37  // Load slant range
38  std::valarray<double> values;
39  isce3::io::loadFromH5(group, "slantRange", values);
40  proc.slantRange(values);
41 
42  // Load zero Doppler time
43  isce3::io::loadFromH5(group, "zeroDopplerTime", values);
44  proc.zeroDopplerTime(values);
45 
46  // Load effective velocity LUT
48  isce3::core::loadCalGrid(group, "effectiveVelocity", lut);
49  proc.effectiveVelocity(lut);
50 
51  // Load azimuth FM rate and Doppler centroid for primary frequency (A)
52  isce3::core::loadCalGrid(group, "frequencyA/azimuthFMRate", lut);
53  proc.azimuthFMRate(lut, 'A');
54  isce3::core::loadCalGrid(group, "frequencyA/dopplerCentroid", lut);
55  proc.dopplerCentroid(lut, 'A');
56 
57  // Check for existence of secondary frequencies
58  if (isce3::io::exists(group, "frequencyB")) {
59  isce3::core::loadCalGrid(group, "frequencyB/azimuthFMRate", lut);
60  proc.azimuthFMRate(lut, 'B');
61  isce3::core::loadCalGrid(group, "frequencyB/dopplerCentroid", lut);
62  proc.dopplerCentroid(lut, 'B');
63  }
64  }
65 
71  inline void loadFromH5(isce3::io::IGroup & group, Swath & swath, char freq) {
72 
73  // Open appropriate frequency group
74  std::string freqString("frequency");
75  freqString.push_back(freq);
76  isce3::io::IGroup fgroup = group.openGroup(freqString);
77 
78  // Load slant range
79  std::valarray<double> s_array;
80  isce3::io::loadFromH5(fgroup, "slantRange", s_array);
81  swath.slantRange(s_array);
82 
83  // Load zero Doppler time
84  std::valarray<double> t_array;
85  isce3::io::loadFromH5(group, "zeroDopplerTime", t_array);
86  swath.zeroDopplerTime(t_array);
87 
88  // Get reference epoch
89  isce3::core::DateTime refEpoch = isce3::io::getRefEpoch(group, "zeroDopplerTime");
90  swath.refEpoch(refEpoch);
91 
92  // Load other parameters
93  double value;
94  isce3::io::loadFromH5(fgroup, "acquiredCenterFrequency", value);
95  swath.acquiredCenterFrequency(value);
96 
97  isce3::io::loadFromH5(fgroup, "processedCenterFrequency", value);
98  swath.processedCenterFrequency(value);
99 
100  isce3::io::loadFromH5(fgroup, "acquiredRangeBandwidth", value);
101  swath.acquiredRangeBandwidth(value);
102 
103  isce3::io::loadFromH5(fgroup, "processedRangeBandwidth", value);
104  swath.processedRangeBandwidth(value);
105 
106  isce3::io::loadFromH5(fgroup, "nominalAcquisitionPRF", value);
107  swath.nominalAcquisitionPRF(value);
108 
109  isce3::io::loadFromH5(fgroup, "sceneCenterGroundRangeSpacing", value);
110  swath.sceneCenterGroundRangeSpacing(value);
111 
112  isce3::io::loadFromH5(fgroup, "processedAzimuthBandwidth", value);
113  swath.processedAzimuthBandwidth(value);
114  }
115 
120  inline void loadFromH5(isce3::io::IGroup & group, std::map<char, Swath> & swaths) {
121  loadFromH5(group, swaths['A'], 'A');
122  if (isce3::io::exists(group, "frequencyB")) {
123  loadFromH5(group, swaths['B'], 'B');
124  }
125  }
126 
131  inline void loadFromH5(isce3::io::IGroup & group, Metadata & meta) {
132 
133  // Get orbit subgroup
134  isce3::io::IGroup orbGroup = group.openGroup("orbit");
135  // Configure a temporary orbit
136  isce3::core::Orbit orbit;
137  isce3::core::loadFromH5(orbGroup, orbit);
138  // Save to metadata
139  meta.orbit(orbit);
140 
141  // Get attitude subgroup
142  isce3::io::IGroup attGroup = group.openGroup("attitude");
143  // Configure a temporary Euler angles object
145  isce3::core::loadFromH5(attGroup, euler);
146  // Save to metadata
147  meta.attitude(euler);
148 
149  // Get processing information subgroup
150  isce3::io::IGroup procGroup = group.openGroup("processingInformation/parameters");
151  // Configure ProcessingInformation
152  loadFromH5(procGroup, meta.procInfo());
153  }
154 
155  }
156 }
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
double processedCenterFrequency() const
Get processed center frequency.
Definition: Swath.h:56
Serialization functions for isce3::core objects.
double sceneCenterGroundRangeSpacing() const
Get scene center ground range spacing.
Definition: Swath.h:90
Definition: Swath.h:25
double processedRangeBandwidth() const
Get processed range bandwidth.
Definition: Swath.h:71
const isce3::core::LUT2d< double > & dopplerCentroid(char freq) const
Get read-only look-up-table for Doppler centroid.
Definition: ProcessingInformation.h:75
const std::valarray< double > & slantRange() const
Get slant range array.
Definition: Swath.h:32
Definition: Metadata.h:24
void loadFromH5(H5obj &h5obj, const std::string &datasetPath, T &v)
Load scalar dataset from HDF5 file.
Definition: Serialization.h:52
Data structure to store date time to nano-sec precision.
Definition: DateTime.h:18
double processedAzimuthBandwidth() const
Get processed azimuth bandwidth.
Definition: Swath.h:99
Serialization utilities using HDF5 API.
const std::valarray< double > & slantRange() const
Get read-only slant range coordinates.
Definition: ProcessingInformation.h:37
double acquiredCenterFrequency() const
Get acquired center frequency.
Definition: Swath.h:51
const isce3::core::Orbit & orbit() const
Get read-only reference to orbit.
Definition: Metadata.h:43
Definition: IH5.h:319
const isce3::core::DateTime & refEpoch() const
Get reference epoch.
Definition: Swath.h:104
const ProcessingInformation & procInfo() const
Get read-only reference to ProcessingInformation.
Definition: Metadata.h:52
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 loadFromH5(isce3::io::IGroup &group, Ellipsoid &ellps)
Load Ellipsoid parameters from HDF5.
Definition: Serialization.h:91
const isce3::core::LUT2d< double > & effectiveVelocity() const
Get read-only look-up-table for effective velocity.
Definition: ProcessingInformation.h:47
Data structure for Euler Angle representation of attitude information.
Definition: EulerAngles.h:17
const isce3::core::EulerAngles & attitude() const
Get read-only reference to attitude.
Definition: Metadata.h:34
Sequence of platform ephemeris samples (state vectors) with uniform temporal spacing, supporting efficient lookup and interpolation.
Definition: Orbit.h:43
IGroup openGroup(const H5std_string &name)
Open a given group.
Definition: IH5.cpp:862
const isce3::core::LUT2d< double > & azimuthFMRate(char freq) const
Get read-only look-up-table for azimuth FM rate.
Definition: ProcessingInformation.h:61
double acquiredRangeBandwidth() const
Get acquired range bandwidth.
Definition: Swath.h:66
double nominalAcquisitionPRF() const
Get nominal acquisition PRF.
Definition: Swath.h:76
const std::valarray< double > & zeroDopplerTime() const
Get zero Doppler time array.
Definition: Swath.h:40
void loadFromH5(isce3::io::IGroup &group, ProcessingInformation &proc)
Load ProcessingInformation from HDF5.
Definition: Serialization.h:35
Definition: ProcessingInformation.h:24
const std::valarray< double > & zeroDopplerTime() const
Get read-only zero Doppler time coordinates.
Definition: ProcessingInformation.h:42

Generated for ISCE3.0 by doxygen 1.8.5.