Download this notebook from github.


Advanced explanation for RCM

Contrary to Sentinel-1, Radarsat Constellation Mission doesn’t have the notion of multi dataset

[1]:
import xsar
import geoviews as gv
import holoviews as hv
import geoviews.feature as gf
hv.extension('bokeh')
path = xsar.get_test_file('RCM1_OK1050603_PK1050605_1_SC50MB_20200214_115905_HH_HV_Z010')
[2]:
import warnings
warnings.filterwarnings('ignore')

Access metadata from a product

Raw information is stocked in different files such as tiff ones (for digital numbers). A file named product.xml is constitued of the main information (geolocation grid, orbit attitude, noise look up tables…). Calibration look up tables are located in xml files. All the useful data is grouped in a datatree thanks to dependencie xarray-safe-rcm. This datatree is than used as an attribute of RcmMeta.

[3]:
#Instanciate a RadarSat2Meta object
meta = xsar.RcmMeta(name=path)
[4]:
#Access the datatree extracted from the reader
meta.dt
[4]:
<xarray.DatasetView>
Dimensions:  ()
Data variables:
    *empty*
Attributes:
    copyright:               RCM Data and Products (c) CSA, 2020 - All Rights...
    productId:               1050605_1
    documentIdentifier:      RCM-SP-53-0419, Issue 2/5
    securityClassification:  Non classifié / Unclassified

Examples of alias to datasets (from the datatree above)

[5]:
#geolocation grid (low resolution)
meta.geoloc
[5]:
<xarray.Dataset>
Dimensions:    (line: 13, pixel: 13)
Coordinates:
  * line       (line) float64 0.0 1.594e+03 3.188e+03 ... 1.753e+04 1.913e+04
  * pixel      (pixel) float64 0.0 1.503e+03 3.007e+03 ... 1.654e+04 1.804e+04
Data variables:
    latitude   (line, pixel) float64 48.72 48.77 48.82 ... 45.78 45.82 45.86
    longitude  (line, pixel) float64 -84.63 -85.03 -85.44 ... -89.82 -90.21
    height     (line, pixel) float64 229.9 229.9 229.9 ... 229.9 229.9 229.9
Attributes:
    footprint:  POLYGON ((-84.63397773377159 48.71534850877886, -89.495130111...

In the metadata class : noise lut, calibration lut and incidence are processed because samples are expressed thanks to a firstPixelValue, a step and a nbOfValues. It is made internally to the reader’s datatree thanks to the method assign_index See Documentation : RCM Image Product Format Definition (7.5.1)

[6]:
#Calibration look up tables in range
meta.lut
[6]:
<xarray.DatasetView>
Dimensions:             (sarCalibrationType: 3, pole: 2, pixel: 786)
Coordinates:
  * sarCalibrationType  (sarCalibrationType) <U12 'Beta Nought' ... 'Sigma No...
  * pole                (pole) <U2 'HH' 'HV'
  * pixel               (pixel) int64 -14 9 32 55 78 ... 17972 17995 18018 18041
Data variables:
    lookup_tables       (pixel, sarCalibrationType, pole) float64 2.731e+09 ....
[7]:
#noise look up tables in range
meta.noise_lut
[7]:
<xarray.DatasetView>
Dimensions:               (pole: 2, sarCalibrationType: 3, pixel: 786)
Coordinates:
  * sarCalibrationType    (sarCalibrationType) <U12 'Beta Nought' ... 'Sigma ...
    pixelFirstNoiseValue  int64 18041
    stepSize              int64 -23
  * pole                  (pole) <U2 'HH' 'HV'
  * pixel                 (pixel) int64 -14 9 32 55 ... 17972 17995 18018 18041
Data variables:
    noiseLevelValues      (pole, sarCalibrationType, pixel) float64 -33.18 .....
Attributes:
    numberOfValues:  786
[8]:
#incidence angles
meta.incidence
[8]:
<xarray.DatasetView>
Dimensions:  (pixel: 786)
Coordinates:
  * pixel    (pixel) int64 -14 9 32 55 78 101 ... 17949 17972 17995 18018 18041
Data variables:
    angles   (pixel) float64 26.87 26.91 26.95 26.99 ... 50.89 50.91 50.93 50.95
Attributes:
    pixelFirstAnglesValue:  18041
    stepSize:               -23
    numberOfValues:         786

Open a dataset

[9]:
# Define the resolution to load the dataset at a lower resolution (if not specified or None, the dataset is loaded at high resolution)
resolution = '1000m'

# Instanciate a RadarSatDataset object
rcmds = xsar.RcmDataset(dataset_id=path, resolution=resolution)

Get the Dataset object

[10]:
rcmds
[10]:
<RcmDataset full coverage object>

Access the metadata object from the Dataset object

[11]:
rcmds.sar_meta
[11]:
<BlockingActorProxy: RcmMeta>

Access the dataset

In this dataset, we can find variables like latitude, longitude, look up tables (before and after denoising), incidence…

[12]:
rcmds.dataset
[12]:
<xarray.Dataset>
Dimensions:          (pol: 2, line: 382, sample: 360)
Coordinates:
  * pol              (pol) <U2 'HH' 'HV'
  * line             (line) float64 24.5 74.5 124.5 ... 1.902e+04 1.907e+04
  * sample           (sample) float64 24.5 74.5 124.5 ... 1.792e+04 1.797e+04
Data variables: (12/21)
    digital_number   (pol, line, sample) uint16 dask.array<chunksize=(1, 382, 360), meta=np.ndarray>
    lines_flipped    bool False
    samples_flipped  bool True
    sampleSpacing    float64 1e+03
    lineSpacing      float64 1e+03
    sigma0_raw       (pol, line, sample) float64 dask.array<chunksize=(1, 382, 360), meta=np.ndarray>
    ...               ...
    ground_heading   (line, sample) float32 dask.array<chunksize=(382, 360), meta=np.ndarray>
    latitude         (line, sample) float64 dask.array<chunksize=(382, 360), meta=np.ndarray>
    longitude        (line, sample) float64 dask.array<chunksize=(382, 360), meta=np.ndarray>
    altitude         (line, sample) float64 dask.array<chunksize=(382, 360), meta=np.ndarray>
    incidence        (line, sample) float32 dask.array<chunksize=(382, 360), meta=np.ndarray>
    elevation        (line, sample) float64 dask.array<chunksize=(382, 360), meta=np.ndarray>
Attributes: (12/18)
    name:                    RCM_DS:/tmp/RCM1_OK1050603_PK1050605_1_SC50MB_20...
    short_name:              RCM_DS:RCM1_OK1050603_PK1050605_1_SC50MB_2020021...
    product:                 Z010
    safe:                    RCM1_OK1050603_PK1050605_1_SC50MB_20200214_11590...
    swath:                   SC50MB
    multidataset:            False
    ...                      ...
    stop_date:               2020-02-14 12:00:02.000000
    footprint:               POLYGON ((-84.63397773377159 48.71534850877886, ...
    coverage:                381km * 359km (line * sample )
    pixel_line_m:            20.0
    pixel_sample_m:          20.0
    approx_transform:        |-0.00,-0.00,-84.71|\n|-0.00, 0.00, 48.73|\n| 0....

Variables lines_flippedand samples_flipped are added to the dataset to know if these have been flipped (in order to follow xsar convention).

See RCM Image Product Format Definition (4.2.1) for more explication about the different cases when lines or samples are flipped

Alternatives solutions to open dataset and datatree

[13]:
# Open dataset
xsar.open_dataset(path, resolution=resolution)
[13]:
<xarray.Dataset>
Dimensions:          (pol: 2, line: 382, sample: 360)
Coordinates:
  * pol              (pol) <U2 'HH' 'HV'
  * line             (line) float64 24.5 74.5 124.5 ... 1.902e+04 1.907e+04
  * sample           (sample) float64 24.5 74.5 124.5 ... 1.792e+04 1.797e+04
Data variables: (12/21)
    digital_number   (pol, line, sample) uint16 dask.array<chunksize=(1, 382, 360), meta=np.ndarray>
    lines_flipped    bool False
    samples_flipped  bool True
    sampleSpacing    float64 1e+03
    lineSpacing      float64 1e+03
    sigma0_raw       (pol, line, sample) float64 dask.array<chunksize=(1, 382, 360), meta=np.ndarray>
    ...               ...
    ground_heading   (line, sample) float32 dask.array<chunksize=(382, 360), meta=np.ndarray>
    latitude         (line, sample) float64 dask.array<chunksize=(382, 360), meta=np.ndarray>
    longitude        (line, sample) float64 dask.array<chunksize=(382, 360), meta=np.ndarray>
    altitude         (line, sample) float64 dask.array<chunksize=(382, 360), meta=np.ndarray>
    incidence        (line, sample) float32 dask.array<chunksize=(382, 360), meta=np.ndarray>
    elevation        (line, sample) float64 dask.array<chunksize=(382, 360), meta=np.ndarray>
Attributes: (12/18)
    name:                    RCM_DS:/tmp/RCM1_OK1050603_PK1050605_1_SC50MB_20...
    short_name:              RCM_DS:RCM1_OK1050603_PK1050605_1_SC50MB_2020021...
    product:                 Z010
    safe:                    RCM1_OK1050603_PK1050605_1_SC50MB_20200214_11590...
    swath:                   SC50MB
    multidataset:            False
    ...                      ...
    stop_date:               2020-02-14 12:00:02.000000
    footprint:               POLYGON ((-84.63397773377159 48.71534850877886, ...
    coverage:                381km * 359km (line * sample )
    pixel_line_m:            20.0
    pixel_sample_m:          20.0
    approx_transform:        |-0.00,-0.00,-84.71|\n|-0.00, 0.00, 48.73|\n| 0....
[14]:
# Open datatree
xsar.open_datatree(path, resolution=resolution)
[14]:
<xarray.DatasetView>
Dimensions:  ()
Data variables:
    *empty*
Attributes: (12/18)
    name:                    RCM_DS:/tmp/RCM1_OK1050603_PK1050605_1_SC50MB_20...
    short_name:              RCM_DS:RCM1_OK1050603_PK1050605_1_SC50MB_2020021...
    product:                 Z010
    safe:                    RCM1_OK1050603_PK1050605_1_SC50MB_20200214_11590...
    swath:                   SC50MB
    multidataset:            False
    ...                      ...
    stop_date:               2020-02-14 12:00:02.000000
    footprint:               POLYGON ((-84.63397773377159 48.71534850877886, ...
    coverage:                381km * 359km (line * sample )
    pixel_line_m:            20.0
    pixel_sample_m:          20.0
    approx_transform:        |-0.00,-0.00,-84.71|\n|-0.00, 0.00, 48.73|\n| 0....

How to apply calibration?

All the operation below are already performed by default for GRD products. So, what is following is a simple explanation about how is made calibration.

Load digital numbers

load_digital_number is a function that allows to load digital numbers from tiff files at chosen resolution and return it as a DataArray. Resampling is made thanks to rasterio.open(tiff_file).read. For dual pol products, there is a tiff file for each pol. So that digital numbers are computed for both pol. Posting of lines and samples is computed thanks to Affine.translation and Affine.scale.

[15]:
import rasterio

#Define resampling method (here it is the root mean square from rasterio)
resampling = rasterio.enums.Resampling.rms

#Define the chunks size for line and samples
chunks = {'line': 5000, 'sample': 5000}
[16]:
dn_low_res = rcmds.load_digital_number(resolution=resolution, resampling=resampling, chunks=chunks)
dn_low_res
[16]:
<xarray.Dataset>
Dimensions:         (pol: 2, line: 382, sample: 360)
Coordinates:
  * pol             (pol) <U2 'HH' 'HV'
  * line            (line) float64 24.5 74.5 124.5 ... 1.902e+04 1.907e+04
  * sample          (sample) float64 24.5 74.5 124.5 ... 1.792e+04 1.797e+04
Data variables:
    digital_number  (pol, line, sample) uint16 dask.array<chunksize=(1, 382, 360), meta=np.ndarray>

Get the raw normalized radar cross section

First the calibration look up tables are loaded at the good resolution in xsar.RcmDataset._luts thanks to the method lazy_load_luts. The resampling is made thanks to an interpolation with the method ̀xsar.RcmDataset._interpolate_var.

This last method uses RectBivariateSplinefor the interpolation, but it is necessary that data is expressed as a 2D vector. Here, calibration/noise look up tables and incidence are expressed as 1D vector (range sample). Consequently, we need to convert these in 2D (adding an azimuth dimension dependency) before applying the interpolation. Conversion is made thanks to numpy.tile, using the low resolution lines expressed in the geolocation grid part of the reader; reducing the calculation. A template of a DataArray that uses the posting of digital numbers (with applied resolution) is given on this interpolation function so the result is now at the right resolution.

_apply_calibration_lut is a method that applies a calibration look up table to digital numbers to compute gamma0, beta0 and sigma0 (depending on the variable name passed in argument) and return the result as a DataArray. It applies the following formula :

\[\frac{(digitalNumbers^2)+offset}{Gain}\]

Reference : RCM Image Product Format Definition (7.5.1)

Different resampling method were tried such as scipy.ndimage.gaussian_filter1d that had the convenience to accept 1d vectors. Data was computed with this function and the chosen posting was this of digital numbers. But in order to be homogenous with other operations made in xsar, we chose to keep the solution with RectBivariateSpline.

[17]:
sigma0_raw = rcmds._apply_calibration_lut('sigma0').sigma0_raw
sigma0_raw
[17]:
<xarray.DataArray 'sigma0_raw' (pol: 2, line: 382, sample: 360)>
dask.array<where, shape=(2, 382, 360), dtype=float64, chunksize=(1, 382, 360), chunktype=numpy.ndarray>
Coordinates:
  * line     (line) float64 24.5 74.5 124.5 ... 1.897e+04 1.902e+04 1.907e+04
  * sample   (sample) float64 24.5 74.5 124.5 ... 1.787e+04 1.792e+04 1.797e+04
  * pol      (pol) <U2 'HH' 'HV'
Attributes:
    pixelFirstLutValue:  18041
    stepSize:            -23
    numberOfValues:      786
    offset:              0.0
[18]:
import matplotlib.pyplot as plt
[19]:
plt.figure(figsize=(13, 6))

plt.subplot(1, 2, 1)
sigma0_raw_cross = sigma0_raw.sel(pol='HV')
plt.pcolor(sigma0_raw_cross, vmax=0.02, cmap='Greys_r')
plt.title('sigma0_raw HV')
plt.xlabel('samples')
plt.ylabel('lines')
plt.colorbar()

plt.subplot(1, 2, 2)
sigma0_raw_co = sigma0_raw.sel(pol='HH')
plt.pcolor(sigma0_raw_co, vmax=0.7, cmap='Greys_r')
plt.title('sigma0_raw HH')
plt.xlabel('samples')
plt.ylabel('lines')
plt.colorbar()
[19]:
<matplotlib.colorbar.Colorbar at 0x7fe788ed2110>
../_images/examples_rcm_38_1.png

How to apply denoising ?

All the operation below are already performed by default for GRD products. So, what is following is a simple explanation about how is made denoising.

How to get the Noise Equivalent Sigma Zero ?

The noise look up tables are loaded at the good resolution in xsar.RcmDataset._noise_luts thanks to the method lazy_load_noise_luts. The resampling is made thanks to an interpolation with the method ̀xsar.RcmDataset._interpolate_var.

Notes: - noise luts are already calibrated so we don’t have to apply a calibration on these - interpolate_var method explained above converts the NoiseLevelValues in linear beacause these are expressed in dB in the reader. The formula used is :

\[NoiseLevelValues_{linear} = 10^\frac{NoiseLevelValues_{dB}}{10}\]
[20]:
nesz_low_res = rcmds.lazy_load_noise_luts().sigma0
nesz_low_res
[20]:
<xarray.DataArray 'sigma0' (pol: 2, line: 382, sample: 360)>
dask.array<concatenate, shape=(2, 382, 360), dtype=float32, chunksize=(1, 382, 360), chunktype=numpy.ndarray>
Coordinates:
  * line                  (line) float64 24.5 74.5 124.5 ... 1.902e+04 1.907e+04
  * sample                (sample) float64 24.5 74.5 ... 1.792e+04 1.797e+04
    pixelFirstNoiseValue  int64 18041
    stepSize              int64 -23
  * pol                   (pol) <U2 'HH' 'HV'
[21]:
plt.pcolor(nesz_low_res.sel(pol='HH'), cmap='Greys_r')
plt.title('nesz (pol = HH)')
plt.xlabel('samples')
plt.ylabel('lines')
plt.colorbar()
[21]:
<matplotlib.colorbar.Colorbar at 0x7fe77f283670>
../_images/examples_rcm_45_1.png

How to get the noise substracted Sigma0

Right now we only have to substract the noise_lut to the raw normalized radar cross section. It is made with the function _add_denoised, that add the variables to the RcmDataset.dataset

[22]:
sigma0 = sigma0_raw - nesz_low_res
sigma0
[22]:
<xarray.DataArray (pol: 2, line: 382, sample: 360)>
dask.array<sub, shape=(2, 382, 360), dtype=float64, chunksize=(1, 382, 360), chunktype=numpy.ndarray>
Coordinates:
  * line                  (line) float64 24.5 74.5 124.5 ... 1.902e+04 1.907e+04
  * sample                (sample) float64 24.5 74.5 ... 1.792e+04 1.797e+04
  * pol                   (pol) <U2 'HH' 'HV'
    pixelFirstNoiseValue  int64 18041
    stepSize              int64 -23

Comparison between noised sigma0 and noised substracted sigma0

[23]:
plt.figure(figsize=(26, 12))

sigma0_cross = sigma0.sel(pol='HV')
sigma0_co = sigma0.sel(pol='HH')

plt.subplot(2,2,1)
plt.pcolor(sigma0_raw_cross, vmax=0.02, cmap='Greys_r')
plt.title('Sigma0 HV with noise')
plt.xlabel('samples')
plt.ylabel('lines')
plt.colorbar()

plt.subplot(2,2,3)
plt.pcolor(sigma0_cross, vmax=0.02, cmap='Greys_r')
plt.title('Sigma0 HV without noise')
plt.xlabel('samples')
plt.ylabel('lines')
plt.colorbar()

plt.subplot(2,2,2)
plt.pcolor(sigma0_raw_co, vmax=0.7, cmap='Greys_r')
plt.title('Sigma0 HH with noise')
plt.xlabel('samples')
plt.ylabel('lines')
plt.colorbar()

plt.subplot(2,2,4)
plt.pcolor(sigma0_co, vmax=0.7, cmap='Greys_r')
plt.title('Sigma0 HH without noise')
plt.xlabel('samples')
plt.ylabel('lines')
plt.colorbar()
[23]:
<matplotlib.colorbar.Colorbar at 0x7fe765d335b0>
../_images/examples_rcm_50_1.png

How to get the incidence ?

Such as noise/calibration look up tables, incidence depends on samples, which are expressed in the reader thanks to step, firstPixelValue and nbOfValues. Samples are computed in the meta class directly on the datatree thanks to the method assign_index.

_load_incidence_from_lut is a function that applies an interpolation with _interpolate_var and then return the incidence at the good resolution

[24]:
incidence = rcmds._load_incidence_from_lut()
incidence
[24]:
<xarray.DataArray "Delayed('RectBivariateSpline-56b4b549-5360-441d-aa-b0467323585fa1c326b419a819e4fe93" (
                                                                                                         line: 382,
                                                                                                         sample: 360)>
dask.array<Delayed('RectBivariateSpline-56b4b549-5360-441d-aa, shape=(382, 360), dtype=float32, chunksize=(382, 360), chunktype=numpy.ndarray>
Coordinates:
  * line     (line) float64 24.5 74.5 124.5 ... 1.897e+04 1.902e+04 1.907e+04
  * sample   (sample) float64 24.5 74.5 124.5 ... 1.787e+04 1.792e+04 1.797e+04

How to get the elevation ?

To get the incidence, we apply a formula :

\[\theta = \arcsin{[\sin{(Incidence)} . \frac{r}{r + h}]}\]
\[( r \text{ is the earth radius} , h \text{ is the orbit altitude} )\]

2 variables give orbit altitude so we considered theSatelliteHeight (and not thealtitude).

RcmDataset._load_elevation_from_lut permit to calculate the elevation (in degrees).

[25]:
elevation = rcmds._load_elevation_from_lut()
elevation
[25]:
<xarray.DataArray (line: 382, sample: 360)>
dask.array<degrees, shape=(382, 360), dtype=float64, chunksize=(382, 360), chunktype=numpy.ndarray>
Coordinates:
  * line     (line) float64 24.5 74.5 124.5 ... 1.897e+04 1.902e+04 1.907e+04
  * sample   (sample) float64 24.5 74.5 124.5 ... 1.787e+04 1.792e+04 1.797e+04