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-2019
6
10
11#pragma once
12
13#include <string>
14
15// isce3::core
16#include <isce3/core/Constants.h>
17#include <isce3/core/Ellipsoid.h>
19
20// isce3::io
21#include <isce3/io/IH5.h>
23
24// isce3::product
25#include <isce3/product/Metadata.h>
26#include <isce3/product/Swath.h>
27#include <isce3/product/Grid.h>
28
29#include <string>
30
32namespace isce3 {
34 namespace product {
35
40 inline void loadFromH5(isce3::io::IGroup & group, ProcessingInformation & proc) {
41
42 // Load effective velocity LUT
44 if (isce3::io::exists(group, "effectiveVelocity")) {
45 isce3::core::loadCalGrid(group, "effectiveVelocity", lut);
46 proc.effectiveVelocity(lut);
47 }
48
49 std::vector<std::string> frequencies{ "A", "B" };
50
51 for (auto frequency : frequencies) {
52
53 std::string frequency_str = "frequency" + frequency;
54 // Check for existence of given frequency group
55 if (isce3::io::exists(group, frequency_str)) {
56
57 // Get processing information subgroup
58
59 // First, try to read the coordinate vectors "zeroDopplerTime"
60 // and "slantRange" from the same group as the LUT.
61 // Check for the existence of the "zeroDopplerTime" H5 dataset in
62 // "frequency{frequency}"
63 if (isce3::io::exists(group, frequency_str+"/zeroDopplerTime")) {
64 isce3::io::IGroup freqGroup = group.openGroup(frequency_str);
65 if (isce3::io::exists(freqGroup, "azimuthFMRate")) {
66 // Load azimuth FM rate LUT (if available)
67 isce3::core::loadCalGrid(freqGroup, "azimuthFMRate", lut);
68 proc.azimuthFMRate(lut, frequency[0]);
69 }
70 if (isce3::io::exists(freqGroup, "dopplerCentroid")) {
71 // Load azimuth Doppler centroid LUT (if available)
72 isce3::core::loadCalGrid(freqGroup, "dopplerCentroid", lut);
73 proc.dopplerCentroid(lut, frequency[0]);
74 }
75 }
76
77 // If not found, read the common coordinate vectors from the parent group.
78 // This supports legacy RSLC products (created prior to ~2020) and products
79 // created by nonofficial tools which may not support the current spec.
80 else {
81 if (isce3::io::exists(group, frequency_str+"/azimuthFMRate")) {
82 // Load azimuth FM rate LUT (if available)
83 isce3::core::loadCalGrid(group, frequency_str+"/azimuthFMRate", lut);
84 proc.azimuthFMRate(lut, frequency[0]);
85 }
86 if (isce3::io::exists(group, frequency_str+"/dopplerCentroid")) {
87 // Load azimuth Doppler centroid LUT (if available)
88 isce3::core::loadCalGrid(group, frequency_str+"/dopplerCentroid", lut);
89 proc.dopplerCentroid(lut, frequency[0]);
90 }
91 }
92 }
93 }
94
95 }
96
102 inline void loadFromH5(isce3::io::IGroup & group, Swath & swath, char freq) {
103
104 // Open appropriate frequency group
105 std::string freqString("frequency");
106 freqString.push_back(freq);
107 isce3::io::IGroup fgroup = group.openGroup(freqString);
108
109 // Load slant range
110 std::valarray<double> s_array;
111 isce3::io::loadFromH5(fgroup, "slantRange", s_array);
112 swath.slantRange(s_array);
113
114 // Load zero Doppler time
115 std::valarray<double> t_array;
116 isce3::io::loadFromH5(group, "zeroDopplerTime", t_array);
117 swath.zeroDopplerTime(t_array);
118
119 // Get reference epoch
120 isce3::core::DateTime refEpoch = isce3::io::getRefEpoch(group, "zeroDopplerTime");
121 swath.refEpoch(refEpoch);
122
123 // Sub-swaths
124 int num_of_sub_swaths = 1;
125 if (isce3::io::exists(fgroup, "numberOfSubSwaths")) {
126 isce3::io::loadFromH5(fgroup, "numberOfSubSwaths", num_of_sub_swaths);
127 }
128
129 swath.subSwaths().numSubSwaths(num_of_sub_swaths);
130
131 for (int i=1; i<=num_of_sub_swaths; ++i) {
132 if (!isce3::io::exists(fgroup, "validSamplesSubSwath" +
133 std::to_string(i))) {
134 break;
135 }
136
137 std::vector<int> image_dims = getImageDims(fgroup,
138 "validSamplesSubSwath" + std::to_string(i));
139
140 if (image_dims[0] != t_array.size()) {
141 std::string error_msg = "ERROR the valid-samples";
142 error_msg += " arrays for sub-swath " + std::to_string(i);
143 error_msg += " has length ";
144 error_msg += std::to_string(image_dims[0]);
145 error_msg += " whereas dataset zeroDopplerTime has";
146 error_msg += " length ";
147 error_msg += std::to_string(t_array.size());
148 throw isce3::except::RuntimeError(ISCE_SRCINFO(), error_msg);
149 }
150
151 isce3::core::Matrix<int> valid_samples_sub_swath_array(
152 image_dims[0], image_dims[1]);
153
155 fgroup, "validSamplesSubSwath" + std::to_string(i),
156 valid_samples_sub_swath_array);
157
158 swath.subSwaths().setValidSamplesArray(i, valid_samples_sub_swath_array);
159 }
160
161 // Load other parameters
162 double value;
163 isce3::io::loadFromH5(fgroup, "acquiredCenterFrequency", value);
164 swath.acquiredCenterFrequency(value);
165
166 isce3::io::loadFromH5(fgroup, "processedCenterFrequency", value);
167 swath.processedCenterFrequency(value);
168
169 isce3::io::loadFromH5(fgroup, "acquiredRangeBandwidth", value);
170 swath.acquiredRangeBandwidth(value);
171
172 isce3::io::loadFromH5(fgroup, "processedRangeBandwidth", value);
173 swath.processedRangeBandwidth(value);
174
175 isce3::io::loadFromH5(fgroup, "nominalAcquisitionPRF", value);
176 swath.nominalAcquisitionPRF(value);
177
178 isce3::io::loadFromH5(group, "zeroDopplerTimeSpacing", value);
179 swath.zeroDopplerTimeSpacing(value);
180
181 isce3::io::loadFromH5(fgroup, "sceneCenterGroundRangeSpacing", value);
183
184 isce3::io::loadFromH5(fgroup, "sceneCenterAlongTrackSpacing", value);
185 swath.sceneCenterAlongTrackSpacing(value);
186
187 isce3::io::loadFromH5(fgroup, "processedAzimuthBandwidth", value);
188 swath.processedAzimuthBandwidth(value);
189 }
190
195 inline void loadFromH5(isce3::io::IGroup & group, std::map<char, Swath> & swaths) {
196 if (isce3::io::exists(group, "frequencyA")) {
197 loadFromH5(group, swaths['A'], 'A');
198 }
199 if (isce3::io::exists(group, "frequencyB")) {
200 loadFromH5(group, swaths['B'], 'B');
201 }
202 }
203
209 inline void populateGridSwathParameterFromH5(isce3::io::IGroup & group, Grid & grid, const char freq) {
210
211 // Open appropriate frequency group
212 std::string freqString("frequency");
213 freqString.push_back(freq);
214
215 isce3::io::IGroup fgroup = group.openGroup(freqString);
216
217 // Read other parameters
218 double value;
219
220 if (isce3::io::exists(fgroup, "rangeBandwidth")) {
221 isce3::io::loadFromH5(fgroup, "rangeBandwidth", value);
222 grid.rangeBandwidth(value);
223 }
224
225 if (isce3::io::exists(fgroup, "azimuthBandwidth")) {
226 isce3::io::loadFromH5(fgroup, "azimuthBandwidth", value);
227 grid.azimuthBandwidth(value);
228 }
229
230 if (isce3::io::exists(fgroup, "centerFrequency")) {
231 isce3::io::loadFromH5(fgroup, "centerFrequency", value);
232 grid.centerFrequency(value);
233 }
234
235 if (isce3::io::exists(fgroup, "slantRangeSpacing")) {
236 isce3::io::loadFromH5(fgroup, "slantRangeSpacing", value);
237 grid.slantRangeSpacing(value);
238 }
239
240 // Search for a dataset called "zeroDopplerTimeSpacing" in the "frequency{A,B}"
241 // group.
242 auto zero_dop_freq_vect = fgroup.find("zeroDopplerTimeSpacing",
243 ".", "DATASET");
244
245 /*
246 GSLC products are expected to have 'zeroDopplerTimeSpacing' in the frequency group
247 */
248 if (zero_dop_freq_vect.size() > 0) {
249 isce3::io::loadFromH5(fgroup, "zeroDopplerTimeSpacing", value);
250 grid.zeroDopplerTimeSpacing(value);
251 } else {
252
253 /*
254 Look for zeroDopplerTimeSpacing within parent group
255 (GCOV products)
256 */
257 auto zero_dop_vect = group.find("zeroDopplerTimeSpacing",
258 ".", "DATASET");
259 if (zero_dop_vect.size() > 0) {
260 isce3::io::loadFromH5(group, "zeroDopplerTimeSpacing", value);
261 grid.zeroDopplerTimeSpacing(value);
262 }
263 // If we don't find zeroDopplerTimeSpacing, we intentionally don't
264 // set a value to not overwrite any existing value.
265 }
266 }
267
273 inline void loadFromH5(isce3::io::IGroup & group, Grid & grid, const char freq) {
274
275 // Open appropriate frequency group
276 std::string freqString("frequency");
277 freqString.push_back(freq);
278 isce3::io::IGroup fgroup = group.openGroup(freqString);
279
280 // Get X and Y coordinates spacing
281 double dx, dy;
282 isce3::io::loadFromH5(fgroup, "xCoordinateSpacing", dx);
283 grid.spacingX(dx);
284
285 isce3::io::loadFromH5(fgroup, "yCoordinateSpacing", dy);
286 grid.spacingY(dy);
287
288 // Load X-coordinates
289 std::valarray<double> x_array;
290 isce3::io::loadFromH5(fgroup, "xCoordinates", x_array);
291 grid.startX(x_array[0] - 0.5 * dx);
292 grid.width(x_array.size());
293
294 // Load Y-coordinates
295 std::valarray<double> y_array;
296 isce3::io::loadFromH5(fgroup, "yCoordinates", y_array);
297 grid.startY(y_array[0] - 0.5 * dy);
298 grid.length(y_array.size());
299
300 auto epsg_freq_vect = fgroup.find("epsg", ".", "DATASET");
301
302 int epsg = -1;
303 // Look for epsg in frequency group
304 if (epsg_freq_vect.size() > 0) {
305 isce3::io::loadFromH5(fgroup, "epsg", epsg);
306 grid.epsg(epsg);
307 } else {
308
309 // Look for epsg in dataset projection
310 auto projection_vect = fgroup.find("projection",
311 ".", "DATASET");
312 if (projection_vect.size() > 0) {
313 for (auto projection_str: projection_vect) {
314 auto projection_obj = fgroup.openDataSet(projection_str);
315 if (projection_obj.attrExists("epsg_code")) {
316 auto attr = projection_obj.openAttribute("epsg_code");
317 attr.read(getH5Type<int>(), &epsg);
318 grid.epsg(epsg);
319 }
320 }
321 }
322 }
323
324 if (epsg == -1) {
325 throw isce3::except::RuntimeError(ISCE_SRCINFO(),
326 "ERROR could not infer EPSG code from input HDF5 file");
327 }
328
329 populateGridSwathParameterFromH5(group, grid, freq);
330
331 }
332
337 inline void loadFromH5(isce3::io::IGroup & group, std::map<char, Grid> & grids) {
338 if (isce3::io::exists(group, "frequencyA")) {
339 loadFromH5(group, grids['A'], 'A');
340 }
341 if (isce3::io::exists(group, "frequencyB")) {
342 loadFromH5(group, grids['B'], 'B');
343 }
344 }
345
350 inline void loadFromH5(isce3::io::IGroup & group, Metadata & meta,
351 const std::string& product_level) {
352
353 // Get orbit subgroup
354 isce3::io::IGroup orbGroup = group.openGroup("orbit");
355 // Configure a temporary orbit
356 isce3::core::Orbit orbit;
357 isce3::core::loadFromH5(orbGroup, orbit);
358 // Save to metadata
359 meta.orbit(orbit);
360
361 // Get attitude subgroup
362 isce3::io::IGroup attGroup = group.openGroup("attitude");
363 // Configure a temporary Euler angles object
364 isce3::core::Attitude attitude;
365 isce3::core::loadFromH5(attGroup, attitude);
366 // Save to metadata
367 meta.attitude(attitude);
368
369 // Get processing information subgroup
370 if (product_level == "L1") {
371 isce3::io::IGroup procGroup = group.openGroup("processingInformation/parameters");
372 // Configure ProcessingInformation
373 loadFromH5(procGroup, meta.procInfo());
374 }
375
376 // Get processing information subgroup for GSLC and GCOV producs
377 if (product_level == "L2" &&
378 isce3::io::exists(group, "sourceData/processingInformation/parameters")) {
379 isce3::io::IGroup procGroup = group.openGroup("sourceData/processingInformation/parameters");
380 // Configure ProcessingInformation
381 loadFromH5(procGroup, meta.procInfo());
382 }
383
384 }
385
386 }
387}
Store and interpolate attitude measurements.
Definition Attitude.h:17
Data structure to store date time to nano-sec precision.
Definition DateTime.h:18
Data structure to store 2D Lookup table.
Definition LUT2d.h:20
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
Definition IH5.h:304
IGroup openGroup(const H5std_string &name)
Open a given group.
Definition IH5.cpp:720
IDataSet openDataSet(const H5std_string &name)
Open a given dataset.
Definition IH5.cpp:712
std::vector< std::string > find(const std::string name, const std::string start=".", const std::string type="BOTH", const std::string path="FULL")
Search function for given name in the group.
Definition IH5.cpp:662
A class for representing Grid metadata originally based on NISAR L2 products.
Definition Grid.h:29
double azimuthBandwidth() const
Get acquired azimuth bandwidth in Hz.
Definition Grid.h:41
size_t epsg() const
Get epsg code for geocoded grid.
Definition Grid.h:112
double centerFrequency() const
Get processed center frequency in Hz.
Definition Grid.h:46
double rangeBandwidth() const
Get acquired range bandwidth in Hz.
Definition Grid.h:36
double spacingY() const
Get the y-coordinate spacing.
Definition Grid.h:97
size_t width() const
Get number of pixels in east-west/x direction for geocoded grid.
Definition Grid.h:102
double slantRangeSpacing() const
Get scene center ground range spacing in meters.
Definition Grid.h:56
double startX() const
Get the X-coordinate start.
Definition Grid.h:82
double startY() const
Get the y-coordinate start.
Definition Grid.h:87
double zeroDopplerTimeSpacing() const
Get time spacing of raster grid in seconds.
Definition Grid.h:75
size_t length() const
Get number of pixels in north-south/y direction for geocoded grid.
Definition Grid.h:107
double spacingX() const
Get the X-coordinate spacing.
Definition Grid.h:92
Definition Metadata.h:24
const ProcessingInformation & procInfo() const
Get read-only reference to ProcessingInformation.
Definition Metadata.h:52
const isce3::core::Orbit & orbit() const
Get read-only reference to orbit.
Definition Metadata.h:43
const isce3::core::Attitude & attitude() const
Get read-only reference to attitude.
Definition Metadata.h:34
Definition ProcessingInformation.h:24
const isce3::core::LUT2d< double > & effectiveVelocity() const
Get read-only look-up-table for effective velocity.
Definition ProcessingInformation.h:37
const isce3::core::LUT2d< double > & azimuthFMRate(char freq) const
Get read-only look-up-table for azimuth FM rate.
Definition ProcessingInformation.h:51
const isce3::core::LUT2d< double > & dopplerCentroid(char freq) const
Get read-only look-up-table for Doppler centroid.
Definition ProcessingInformation.h:65
void setValidSamplesArray(const int n, const isce3::core::Matrix< int > &v)
Set valid samples for a sub-swath's array indexed from 1 (1st sub-swath)
Definition SubSwaths.cpp:30
int numSubSwaths() const
Get number of sub-swaths.
Definition SubSwaths.h:61
Definition Swath.h:23
const isce3::core::DateTime & refEpoch() const
Get reference epoch.
Definition Swath.h:113
double sceneCenterAlongTrackSpacing() const
Get scene center along track spacing.
Definition Swath.h:90
double processedAzimuthBandwidth() const
Get processed azimuth bandwidth.
Definition Swath.h:108
double acquiredCenterFrequency() const
Get acquired center frequency.
Definition Swath.h:55
double processedCenterFrequency() const
Get processed center frequency.
Definition Swath.h:60
const std::valarray< double > & zeroDopplerTime() const
Get zero Doppler time array.
Definition Swath.h:41
isce3::product::SubSwaths & subSwaths()
Get SubSwaths object.
Definition Swath.h:118
double sceneCenterGroundRangeSpacing() const
Get scene center ground range spacing.
Definition Swath.h:99
double processedRangeBandwidth() const
Get processed range bandwidth.
Definition Swath.h:75
double acquiredRangeBandwidth() const
Get acquired range bandwidth.
Definition Swath.h:70
double zeroDopplerTimeSpacing() const
Get time spacing of raster grid.
Definition Swath.h:85
double nominalAcquisitionPRF() const
Get nominal acquisition PRF.
Definition Swath.h:80
const std::valarray< double > & slantRange() const
Get slant range array.
Definition Swath.h:30
Serialization functions for isce3::core objects.
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 loadFromH5(isce3::io::IGroup &group, Ellipsoid &ellps)
Load Ellipsoid parameters from HDF5.
Definition Serialization.h:41
Serialization utilities using HDF5 API.
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
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
The isce3::product namespace.
Definition forward.h:3
void populateGridSwathParameterFromH5(isce3::io::IGroup &group, Grid &grid, const char freq)
Populate swath-related parameters of the Grid object from HDF5.
Definition Serialization.h:209
void loadFromH5(isce3::io::IGroup &group, ProcessingInformation &proc)
Load ProcessingInformation from HDF5.
Definition Serialization.h:40
base interpolator is an abstract base class
Definition BinarySearchFunc.cpp:5

Generated for ISCE3.0 by doxygen 1.13.2.