diff --git a/.circleci/config.yml b/.circleci/config.yml index d59f8d47f..4ba466776 100644 --- a/.circleci/config.yml +++ b/.circleci/config.yml @@ -61,13 +61,19 @@ jobs: ${MINTPY_HOME}/tests/smallbaselineApp.py --dir ${HOME}/data --dset SanFranSenDT42 - run: - name: Integration Test 3/4 - WellsEnvD2T399 (Gamma) + name: Integration Test 3/5 - RidgecrestSenDT71 (HyP3) + command: | + mkdir -p ${HOME}/data + ${MINTPY_HOME}/tests/smallbaselineApp.py --dir ${HOME}/data --dset RidgecrestSenDT71 + + - run: + name: Integration Test 4/5 - WellsEnvD2T399 (Gamma) command: | mkdir -p ${HOME}/data ${MINTPY_HOME}/tests/smallbaselineApp.py --dir ${HOME}/data --dset WellsEnvD2T399 - run: - name: Integration Test 4/4 - WCapeSenAT29 (SNAP) + name: Integration Test 5/5 - WCapeSenAT29 (SNAP) command: | mkdir -p ${HOME}/data ${MINTPY_HOME}/tests/smallbaselineApp.py --dir ${HOME}/data --dset WCapeSenAT29 diff --git a/docs/api/attributes.md b/docs/api/attributes.md index 40cc16a73..18508f31e 100644 --- a/docs/api/attributes.md +++ b/docs/api/attributes.md @@ -9,7 +9,7 @@ If using ROI_PAC as the InSAR processor, both **baseline parameter RSC** file (i + X/Y_FIRST = (for geocoded product) Longitude/easting/X and latitude/northing/Y coordinate in degrees/meters of the upper left corner of the first pixel. + X/Y_STEP = (for geocoded product) Ground resolution in degrees/meters in X/Y direction. + X/Y_UNIT = (for geocoded product) Coordinate unit in X/Y direction: degrees or meters. -+ LAT/LON_REF1/2/3/4 = Latitude/longitude at corner 1/2/3/4 (degree), used in save_unavco, PyAPS (DEM file in radar coord), not accurate; number named in order of first line near/far range, last line near/far range. ++ LAT/LON_REF1/2/3/4 = Latitude/northing and longitude/easting at corner 1/2/3/4 (in degrees or meters), used in save_unavco, PyAPS (DEM file in radar coord), not accurate; number named in order of first line near/far range, last line near/far range. + WAVELENGTH = Radar wavelength in meters. + RANGE_PIXEL_SIZE = Slant range pixel size (search for pixel_ratio to convert to ground size, in m), used in dem_error, incidence_angle, multilook, transect. + EARTH_RADIUS = Best fitting spheroid radius in meters, used in dem_error, incidence_angle, convert2mat. @@ -30,7 +30,7 @@ The following attributes vary for each interferogram: + ANTENNA_SIDE = -1 for right looking radar, used in save_unavco + AZIMUTH_PIXEL_SIZE = Azimuth pixel size at orbital altitude (multiply by Re/(Re+h) for ground size (m), where Re is the local earth radius), used in baseline_error/trop and multilook. -+ HEADING = Spacecraft heading at peg point (degree), measured from the north with clock-wise as positive, used in asc_desc ++ HEADING = Spacecraft heading at peg point (degrees), measured from the north with clock-wise as positive, used in asc_desc + PRF = Pulse repetition frequency (Hz), used in save_unavco ### Self-generated attributes ### @@ -47,7 +47,8 @@ The following attributes vary for each interferogram: + NO_DATA_VALUE = No data value, value that should be ignored. + UNIT = data unit, i.e. m, m/yr, radian, and 1 for file without unit, such as coherence [[source]](https://github.com/insarlab/MintPy/blob/main/src/mintpy/objects/stack.py#L75) + REF_DATE = reference date -+ REF_X/Y/LAT/LON = column/row/latitude/longitude of reference point ++ REF_X/Y = column/row of the reference point ++ REF_LAT/LON = latitude/northing and longitude/easting of the reference point (in degrees or meters) + SUBSET_XMIN/XMAX/YMIN/YMAX = start/end column/row number of subset in the original coverage + MODIFICATION_TIME = dataset modification time, exists in ifgramStack.h5 file for 3D dataset, used for "--update" option of unwrap error corrections. + NCORRLOOKS = number of independent looks, as explained in [SNAPHU](https://web.stanford.edu/group/radar/softwareandlinks/sw/snaphu/snaphu.conf.full) diff --git a/src/mintpy/prep_hyp3.py b/src/mintpy/prep_hyp3.py index 2817d1b1b..dbd1fc874 100644 --- a/src/mintpy/prep_hyp3.py +++ b/src/mintpy/prep_hyp3.py @@ -67,10 +67,6 @@ def add_hyp3_metadata(fname, meta, is_ifg=True): S = N + float(meta['Y_STEP']) * int(meta['LENGTH']) E = W + float(meta['X_STEP']) * int(meta['WIDTH']) - # convert UTM to lat/lon - N, W = ut.utm2latlon(meta, W, N) - S, E = ut.utm2latlon(meta, E, S) - if meta['ORBIT_DIRECTION'] == 'ASCENDING': meta['LAT_REF1'] = str(S) meta['LAT_REF2'] = str(S) @@ -110,6 +106,11 @@ def add_hyp3_metadata(fname, meta, is_ifg=True): meta['P_BASELINE_TOP_HDR'] = hyp3_meta['Baseline'] meta['P_BASELINE_BOTTOM_HDR'] = hyp3_meta['Baseline'] + # HDF-EOS5 metadata + if hyp3_meta['ReferenceGranule'].startswith('S1'): + meta['beam_mode'] = 'IW' + meta['unwrap_method'] = hyp3_meta['Unwrappingtype'] + return(meta) diff --git a/src/mintpy/save_hdfeos5.py b/src/mintpy/save_hdfeos5.py index 369cf0cdb..f60ac7c44 100644 --- a/src/mintpy/save_hdfeos5.py +++ b/src/mintpy/save_hdfeos5.py @@ -115,20 +115,24 @@ def metadata_mintpy2unavco(meta_in, dateList, geom_file): unavco_meta['last_date'] = dt.datetime.strptime(dateList[-1], '%Y%m%d').isoformat()[0:10] # footprint - lons = [meta['LON_REF1'], - meta['LON_REF3'], - meta['LON_REF4'], - meta['LON_REF2'], - meta['LON_REF1']] - - lats = [meta['LAT_REF1'], - meta['LAT_REF3'], - meta['LAT_REF4'], - meta['LAT_REF2'], - meta['LAT_REF1']] + lons = [float(meta['LON_REF1']), + float(meta['LON_REF3']), + float(meta['LON_REF4']), + float(meta['LON_REF2']), + float(meta['LON_REF1'])] + + lats = [float(meta['LAT_REF1']), + float(meta['LAT_REF3']), + float(meta['LAT_REF4']), + float(meta['LAT_REF2']), + float(meta['LAT_REF1'])] + + # convert UTM to lat/lon + if not meta_in.get('Y_UNIT', 'degrees').lower().startswith('deg'): + lats, lons = ut.utm2latlon(meta_in, easting=lons, northing=lats) unavco_meta['scene_footprint'] = "POLYGON((" + ",".join( - [lon+' '+lat for lon, lat in zip(lons, lats)]) + "))" + [f'{lon} {lat}' for lon, lat in zip(lons, lats)]) + "))" unavco_meta['history'] = dt.datetime.utcnow().isoformat()[0:10] @@ -239,6 +243,10 @@ def get_output_filename(metadata, suffix=None, update_mode=False, subset_mode=Fa lat0 = lat1 + float(metadata['Y_STEP']) * int(metadata['LENGTH']) lon1 = lon0 + float(metadata['X_STEP']) * int(metadata['WIDTH']) + # convert UTM to lat/lon + if not metadata.get('Y_UNIT', 'degrees').lower().startswith('deg'): + [lat0, lat1], [lon0, lon1] = ut.utm2latlon(metadata, easting=[lon0, lon1], northing=[lat0, lat1]) + lat0Str = f'N{round(lat0*1e3):05d}' lat1Str = f'N{round(lat1*1e3):05d}' lon0Str = f'E{round(lon0*1e3):06d}' diff --git a/src/mintpy/tropo_pyaps3.py b/src/mintpy/tropo_pyaps3.py index 49867cd0d..13c0b7475 100644 --- a/src/mintpy/tropo_pyaps3.py +++ b/src/mintpy/tropo_pyaps3.py @@ -275,10 +275,8 @@ def get_bounding_box(meta, geom_file=None): lat1 = lat0 + lat_step * (length - 1) lon1 = lon0 + lon_step * (width - 1) - # 'Y_FIRST' not in 'degree' - # e.g. meters for UTM projection from ASF HyP3 - y_unit = meta.get('Y_UNIT', 'degrees').lower() - if not y_unit.startswith('deg'): + # for UTM projection, e.g. ASF HyP3 + if not meta.get('Y_UNIT', 'degrees').lower().startswith('deg'): lat0, lon0 = ut.utm2latlon(meta, easting=lon0, northing=lat0) lat1, lon1 = ut.utm2latlon(meta, easting=lon1, northing=lat1) @@ -300,8 +298,12 @@ def get_bounding_box(meta, geom_file=None): lon1 = np.nanmax(lons) else: + # use the rough (not accurate) lat/lon info of the four corners lats = [float(meta[f'LAT_REF{i}']) for i in [1,2,3,4]] lons = [float(meta[f'LON_REF{i}']) for i in [1,2,3,4]] + # for UTM projection, e.g. ASF HyP3 + if not meta.get('Y_UNIT', 'degrees').lower().startswith('deg'): + lats, lons = ut.utm2latlon(meta, easting=lons, northing=lats) lat0 = np.mean(lats[0:2]) lat1 = np.mean(lats[2:4]) lon0 = np.mean(lons[0:3:2]) diff --git a/src/mintpy/utils/utils0.py b/src/mintpy/utils/utils0.py index f390e72a6..1e96c38c1 100644 --- a/src/mintpy/utils/utils0.py +++ b/src/mintpy/utils/utils0.py @@ -345,10 +345,10 @@ def utm2latlon(meta, easting, northing): Parameters: meta - dict, mintpy attributes that includes: UTM_ZONE - easting - scalar or 1/2D np.ndarray, UTM coordinates in x direction - northing - scalar or 1/2D np.ndarray, UTM coordinates in y direction - Returns: lat - scalar or 1/2D np.ndarray, WGS 84 coordinates in y direction - lon - scalar or 1/2D np.ndarray, WGS 84 coordinates in x direction + easting - scalar/list/tuple/1-2D np.ndarray, UTM coordinates in x direction + northing - scalar/list/tuple/1-2D np.ndarray, UTM coordinates in y direction + Returns: lat - scalar/list/tuple/1-2D np.ndarray, WGS 84 coordinates in y direction + lon - scalar/list/tuple/1-2D np.ndarray, WGS 84 coordinates in x direction """ import utm zone_num = int(meta['UTM_ZONE'][:-1]) @@ -357,20 +357,33 @@ def utm2latlon(meta, easting, northing): # which can be common for large area analysis, e.g. the Norwegian mapping authority # publishes a height data in UTM zone 33 coordinates for the whole country, even though # most of it is technically outside zone 33. - lat, lon = utm.to_latlon(easting, northing, zone_num, northern=northern, strict=False) + lat, lon = utm.to_latlon(np.array(easting), np.array(northing), zone_num, + northern=northern, strict=False) + + # output format + if any(isinstance(x, (list, tuple)) for x in [easting, northing]): + lat = lat.tolist() + lon = lon.tolist() + return lat, lon def latlon2utm(lat, lon): """Convert latitude/longitude in degrees to UTM easting/northing in meters. - Parameters: lat - scalar or 1/2D np.ndarray, WGS 84 coordinates in y direction - lon - scalar or 1/2D np.ndarray, WGS 84 coordinates in x direction - Returns: easting - scalar or 1/2D np.ndarray, UTM coordinates in x direction - northing - scalar or 1/2D np.ndarray, UTM coordinates in y direction + Parameters: lat - scalar/list/tuple/1-2D np.ndarray, WGS 84 coordinates in y direction + lon - scalar/list/tuple/1-2D np.ndarray, WGS 84 coordinates in x direction + Returns: easting - scalar/list/tuple/1-2D np.ndarray, UTM coordinates in x direction + northing - scalar/list/tuple/1-2D np.ndarray, UTM coordinates in y direction """ import utm - easting, northing = utm.from_latlon(lat, lon)[:2] + easting, northing = utm.from_latlon(np.array(lat), np.array(lon))[:2] + + # output format + if any(isinstance(x, (list, tuple)) for x in [lat, lon]): + easting = easting.tolist() + northing = northing.tolist() + return northing, easting diff --git a/tests/configs/RidgecrestSenDT71.txt b/tests/configs/RidgecrestSenDT71.txt index 0dc5fdb56..1d9dc34f5 100644 --- a/tests/configs/RidgecrestSenDT71.txt +++ b/tests/configs/RidgecrestSenDT71.txt @@ -19,4 +19,9 @@ mintpy.troposphericDelay.method = pyaps mintpy.topographicResidual = yes mintpy.topographicResidual.stepFuncDate = 20190706T0320 mintpy.topographicResidual.pixelwiseGeometry = no +mintpy.save.hdfEos5 = yes +mintpy.save.hdfEos5.subset = yes mintpy.plot.maxMemory = 2 + +##————————————————————————————— HDF-EOS5 Attributes -—————————————————————————## +relative_orbit = 71