Download this notebook from github.


GMFs and LUTs

This notebook will show how to use xsarsea.windspeed.models.Model

Models are functions (GMF) or lookup table (LUT) that returns a simulated sigma0 from incidence angle, wind speed, and wind direction relative to the antenna.

Models are used for wind invertion with xsarsea.windspeed.invert_from_model, but they also be used independently.

[1]:
import xsar
from xsarsea import windspeed
import xsarsea
import xarray as xr
import numpy as np
import holoviews as hv
hv.extension('bokeh')
[2]:
# optional debug messages
import logging
logging.basicConfig()
logging.getLogger('xsarsea.windspeed').setLevel(logging.DEBUG) # or .setLevel(logging.INFO)

Available models

Available models can be retrieved with xsarsea.windspeed.available_models. By default, analytical models are already in the dataframe.

[3]:
windspeed.available_models()
[3]:
alias pol model
gmf_cmod5 cmod5 VV <GmfModel('gmf_cmod5') pol=VV>
gmf_cmod5n cmod5n VV <GmfModel('gmf_cmod5n') pol=VV>
gmf_cmodifr2 cmodifr2 VV <GmfModel('gmf_cmodifr2') pol=VV>
gmf_rs2_v2 rs2_v2 VH <GmfModel('gmf_rs2_v2') pol=VH>
gmf_s1_v2 s1_v2 VH <GmfModel('gmf_s1_v2') pol=VH>
gmf_rcm_noaa rcm_noaa VH <GmfModel('gmf_rcm_noaa') pol=VH>

Add all models

replace paths by your own path containing all luts

[4]:
nc_luts_path = xsarsea.get_test_file('nc_luts_reduce')
path_cmod7 = xsarsea.get_test_file("cmod7_and_python_script")
windspeed.register_luts(nc_luts_path, path_cmod7)
windspeed.available_models()
DEBUG:xsarsea.windspeed:register model nc_lut_gmf_cmod5n_Rlow_hh_mouche1 pol=HH units=linear
 inc_range=[17.0, 50.0] wspd_range=[0.2, 50.0] phi_range=[0.0, 180.0]
 inc_step=0.1 wspd_step=0.1 phi_step=1.0
 inc_step_lr=1.0 wspd_step_lr=0.2 phi_step_lr=2.5
DEBUG:xsarsea.windspeed:register model nc_lut_cmodms1ahw pol=VH units=dB
 inc_range=[17.0, 50.0] wspd_range=[3.0, 80.0] phi_range=None
 inc_step=0.1 wspd_step=0.1 phi_step=1.0
 inc_step_lr=1.0 wspd_step_lr=0.2 phi_step_lr=2.5
DEBUG:xsarsea.windspeed:register model nc_lut_gmf_cmod5n_Rlow_hh_zhangB pol=HH units=linear
 inc_range=[17.0, 50.0] wspd_range=[0.2, 50.0] phi_range=[0.0, 180.0]
 inc_step=0.1 wspd_step=0.1 phi_step=1.0
 inc_step_lr=1.0 wspd_step_lr=0.2 phi_step_lr=2.5
DEBUG:xsarsea.windspeed:register model nc_lut_gmf_cmod5n_Rlow_hh_zhangA pol=HH units=linear
 inc_range=[17.0, 50.0] wspd_range=[0.2, 50.0] phi_range=[0.0, 180.0]
 inc_step=0.1 wspd_step=0.1 phi_step=1.0
 inc_step_lr=1.0 wspd_step_lr=0.2 phi_step_lr=2.5
DEBUG:xsarsea.windspeed:register model nc_lut_gmf_cmod7_Rlow_hh_zhangB pol=HH units=linear
 inc_range=[16.0, 66.0] wspd_range=[0.2, 50.0] phi_range=[0.0, 180.0]
 inc_step=0.1 wspd_step=0.1 phi_step=1.0
 inc_step_lr=1 wspd_step_lr=0.2 phi_step_lr=2.5
DEBUG:xsarsea.windspeed:register model nc_lut_gmf_cmod7_Rlow_hh_zhangA pol=HH units=linear
 inc_range=[16.0, 66.0] wspd_range=[0.2, 50.0] phi_range=[0.0, 180.0]
 inc_step=0.1 wspd_step=0.1 phi_step=1.0
 inc_step_lr=1 wspd_step_lr=0.2 phi_step_lr=2.5
DEBUG:xsarsea.windspeed:register model nc_lut_gmf_cmod7_Rlow_hh_mouche1 pol=HH units=linear
 inc_range=[16.0, 66.0] wspd_range=[0.2, 50.0] phi_range=[0.0, 180.0]
 inc_step=0.1 wspd_step=0.1 phi_step=1.0
 inc_step_lr=1 wspd_step_lr=0.2 phi_step_lr=2.5
DEBUG:xsarsea.windspeed:register model gmf_cmod7 pol=VV units=None
 inc_range=[16.0, 66.0] wspd_range=None phi_range=None
 inc_step=0.1 wspd_step=0.1 phi_step=1.0
 inc_step_lr=1.0 wspd_step_lr=0.2 phi_step_lr=2.5
[4]:
alias pol model
gmf_cmod7 cmod7 VV <Cmod7Model('gmf_cmod7') pol=VV>
gmf_cmod5 cmod5 VV <GmfModel('gmf_cmod5') pol=VV>
gmf_cmod5n cmod5n VV <GmfModel('gmf_cmod5n') pol=VV>
gmf_cmodifr2 cmodifr2 VV <GmfModel('gmf_cmodifr2') pol=VV>
gmf_rs2_v2 rs2_v2 VH <GmfModel('gmf_rs2_v2') pol=VH>
gmf_s1_v2 s1_v2 VH <GmfModel('gmf_s1_v2') pol=VH>
gmf_rcm_noaa rcm_noaa VH <GmfModel('gmf_rcm_noaa') pol=VH>
nc_lut_gmf_cmod5n_Rlow_hh_mouche1 cmod5n_Rlow_hh_mouche1 HH <NcLutModel('nc_lut_gmf_cmod5n_Rlow_hh_mouche1...
nc_lut_cmodms1ahw cmodms1ahw VH <NcLutModel('nc_lut_cmodms1ahw') pol=VH>
nc_lut_gmf_cmod5n_Rlow_hh_zhangB cmod5n_Rlow_hh_zhangB HH <NcLutModel('nc_lut_gmf_cmod5n_Rlow_hh_zhangB'...
nc_lut_gmf_cmod5n_Rlow_hh_zhangA cmod5n_Rlow_hh_zhangA HH <NcLutModel('nc_lut_gmf_cmod5n_Rlow_hh_zhangA'...
nc_lut_gmf_cmod7_Rlow_hh_zhangB cmod7_Rlow_hh_zhangB HH <NcLutModel('nc_lut_gmf_cmod7_Rlow_hh_zhangB')...
nc_lut_gmf_cmod7_Rlow_hh_zhangA cmod7_Rlow_hh_zhangA HH <NcLutModel('nc_lut_gmf_cmod7_Rlow_hh_zhangA')...
nc_lut_gmf_cmod7_Rlow_hh_mouche1 cmod7_Rlow_hh_mouche1 HH <NcLutModel('nc_lut_gmf_cmod7_Rlow_hh_mouche1'...

Add models by type

Adding analytical models

they are already available in the dataframe by default

Adding netcdf models (LUT)

Netcdf models are not available by default, because they needs to be loaded from external file with xsarsea.windspeed.register_nc_luts

[5]:
nc_luts_path = xsarsea.get_test_file('nc_luts_reduce')
windspeed.register_nc_luts(nc_luts_path)
windspeed.available_models()
DEBUG:xsarsea.windspeed:register model nc_lut_gmf_cmod5n_Rlow_hh_mouche1 pol=HH units=linear
 inc_range=[17.0, 50.0] wspd_range=[0.2, 50.0] phi_range=[0.0, 180.0]
 inc_step=0.1 wspd_step=0.1 phi_step=1.0
 inc_step_lr=1.0 wspd_step_lr=0.2 phi_step_lr=2.5
DEBUG:xsarsea.windspeed:register model nc_lut_cmodms1ahw pol=VH units=dB
 inc_range=[17.0, 50.0] wspd_range=[3.0, 80.0] phi_range=None
 inc_step=0.1 wspd_step=0.1 phi_step=1.0
 inc_step_lr=1.0 wspd_step_lr=0.2 phi_step_lr=2.5
DEBUG:xsarsea.windspeed:register model nc_lut_gmf_cmod5n_Rlow_hh_zhangB pol=HH units=linear
 inc_range=[17.0, 50.0] wspd_range=[0.2, 50.0] phi_range=[0.0, 180.0]
 inc_step=0.1 wspd_step=0.1 phi_step=1.0
 inc_step_lr=1.0 wspd_step_lr=0.2 phi_step_lr=2.5
DEBUG:xsarsea.windspeed:register model nc_lut_gmf_cmod5n_Rlow_hh_zhangA pol=HH units=linear
 inc_range=[17.0, 50.0] wspd_range=[0.2, 50.0] phi_range=[0.0, 180.0]
 inc_step=0.1 wspd_step=0.1 phi_step=1.0
 inc_step_lr=1.0 wspd_step_lr=0.2 phi_step_lr=2.5
DEBUG:xsarsea.windspeed:register model nc_lut_gmf_cmod7_Rlow_hh_zhangB pol=HH units=linear
 inc_range=[16.0, 66.0] wspd_range=[0.2, 50.0] phi_range=[0.0, 180.0]
 inc_step=0.1 wspd_step=0.1 phi_step=1.0
 inc_step_lr=1 wspd_step_lr=0.2 phi_step_lr=2.5
DEBUG:xsarsea.windspeed:register model nc_lut_gmf_cmod7_Rlow_hh_zhangA pol=HH units=linear
 inc_range=[16.0, 66.0] wspd_range=[0.2, 50.0] phi_range=[0.0, 180.0]
 inc_step=0.1 wspd_step=0.1 phi_step=1.0
 inc_step_lr=1 wspd_step_lr=0.2 phi_step_lr=2.5
DEBUG:xsarsea.windspeed:register model nc_lut_gmf_cmod7_Rlow_hh_mouche1 pol=HH units=linear
 inc_range=[16.0, 66.0] wspd_range=[0.2, 50.0] phi_range=[0.0, 180.0]
 inc_step=0.1 wspd_step=0.1 phi_step=1.0
 inc_step_lr=1 wspd_step_lr=0.2 phi_step_lr=2.5
[5]:
alias pol model
gmf_cmod7 cmod7 VV <Cmod7Model('gmf_cmod7') pol=VV>
gmf_cmod5 cmod5 VV <GmfModel('gmf_cmod5') pol=VV>
gmf_cmod5n cmod5n VV <GmfModel('gmf_cmod5n') pol=VV>
gmf_cmodifr2 cmodifr2 VV <GmfModel('gmf_cmodifr2') pol=VV>
gmf_rs2_v2 rs2_v2 VH <GmfModel('gmf_rs2_v2') pol=VH>
gmf_s1_v2 s1_v2 VH <GmfModel('gmf_s1_v2') pol=VH>
gmf_rcm_noaa rcm_noaa VH <GmfModel('gmf_rcm_noaa') pol=VH>
nc_lut_gmf_cmod5n_Rlow_hh_mouche1 cmod5n_Rlow_hh_mouche1 HH <NcLutModel('nc_lut_gmf_cmod5n_Rlow_hh_mouche1...
nc_lut_cmodms1ahw cmodms1ahw VH <NcLutModel('nc_lut_cmodms1ahw') pol=VH>
nc_lut_gmf_cmod5n_Rlow_hh_zhangB cmod5n_Rlow_hh_zhangB HH <NcLutModel('nc_lut_gmf_cmod5n_Rlow_hh_zhangB'...
nc_lut_gmf_cmod5n_Rlow_hh_zhangA cmod5n_Rlow_hh_zhangA HH <NcLutModel('nc_lut_gmf_cmod5n_Rlow_hh_zhangA'...
nc_lut_gmf_cmod7_Rlow_hh_zhangB cmod7_Rlow_hh_zhangB HH <NcLutModel('nc_lut_gmf_cmod7_Rlow_hh_zhangB')...
nc_lut_gmf_cmod7_Rlow_hh_zhangA cmod7_Rlow_hh_zhangA HH <NcLutModel('nc_lut_gmf_cmod7_Rlow_hh_zhangA')...
nc_lut_gmf_cmod7_Rlow_hh_mouche1 cmod7_Rlow_hh_mouche1 HH <NcLutModel('nc_lut_gmf_cmod7_Rlow_hh_mouche1'...

Adding cmod7

[6]:
try :
    path_cmod7 = xsar.get_test_file("cmod7_and_python_script")
    windspeed.register_cmod7(path_cmod7)
except Exception as e:
    print(e)
DEBUG:xsarsea.windspeed:register model gmf_cmod7 pol=VV units=None
 inc_range=[16.0, 66.0] wspd_range=None phi_range=None
 inc_step=0.1 wspd_step=0.1 phi_step=1.0
 inc_step_lr=1.0 wspd_step_lr=0.2 phi_step_lr=2.5

Using models

Model instance can be retrieved with xsarsea.windspeed.get_model

[7]:
cmod5 = windspeed.get_model('cmod5')
cmod5
[7]:
<GmfModel('gmf_cmod5') pol=VV>

Models can be used as a regular function because they have a __call__ method.

Argument to __call__ are (incidence, wspd, phi) (phi is optionnal for crosspol lut).

If arguments are 1d arrays, output shape will be (incidence,.size wspd.size, phi.size)

2d arrays input are only implemented for GmfModel. The output will have the same shape as the input.

[8]:
incidence = np.array([25,35,45])
wspd = np.array([5, 40])
phi = np.array([0, 45, 90])
cmod5(incidence, wspd, phi)
DEBUG:xsarsea.windspeed:timing _gmf_function : 0.40s. mem: +17.0Mb
DEBUG:xsarsea.windspeed:timing _gmf_function : 0.76s. mem: +20.6Mb
DEBUG:xsarsea.windspeed:timing _get_function_for_args : 0.76s. mem: +20.6Mb
[8]:
<xarray.DataArray (incidence: 3, wspd: 2, phi: 3)>
array([[[0.14506564, 0.12403658, 0.10445991],
        [0.75708097, 0.66227162, 0.57208342]],

       [[0.02992838, 0.02286184, 0.01591826],
        [0.29075383, 0.26333267, 0.23685913]],

       [[0.01090204, 0.00765933, 0.00460774],
        [0.1568922 , 0.14702434, 0.13737013]]])
Coordinates:
  * incidence  (incidence) int64 25 35 45
  * wspd       (wspd) int64 5 40
  * phi        (phi) int64 0 45 90
Attributes:
    units:    linear

The full lut can also be retrieved with Model.to_lut

[9]:
cmod5.to_lut()
DEBUG:xsarsea.windspeed:_raw_lut gmf at res low, inc_step = 1.0 , wspd_step = 0.2, phi_step = 2.5
DEBUG:xsarsea.windspeed:timing _gmf_function : 0.00s. mem: +0.0Mb
DEBUG:xsarsea.windspeed:timing _get_function_for_args : 0.00s. mem: +0.0Mb
DEBUG:xsarsea.windspeed:timing _raw_lut : 0.13s. mem: +7.1Mb
DEBUG:xsarsea.windspeed:lut_resolution gmf_cmod5 from _raw_lut : low
DEBUG:xsarsea.windspeed:desired gmf_cmod5 resolution : high
DEBUG:xsarsea.windspeed:interp lut gmf_cmod5 to high res
[9]:
<xarray.DataArray 'sigma0_model' (incidence: 501, wspd: 499, phi: 181)>
array([[[2.76625219e-01, 2.76612056e-01, 2.76598894e-01, ...,
         2.75679466e-01, 2.75692262e-01, 2.75705057e-01],
        [3.33053413e-01, 3.33037643e-01, 3.33021874e-01, ...,
         3.32363669e-01, 3.32379177e-01, 3.32394684e-01],
        [3.89481608e-01, 3.89463230e-01, 3.89444853e-01, ...,
         3.89047872e-01, 3.89066092e-01, 3.89084311e-01],
        ...,
        [2.39831755e+00, 2.39798890e+00, 2.39766026e+00, ...,
         2.39759365e+00, 2.39792226e+00, 2.39825087e+00],
        [2.39541012e+00, 2.39508353e+00, 2.39475693e+00, ...,
         2.39469156e+00, 2.39501813e+00, 2.39534469e+00],
        [2.39250269e+00, 2.39217815e+00, 2.39185360e+00, ...,
         2.39178948e+00, 2.39211400e+00, 2.39243852e+00]],

       [[2.65686830e-01, 2.65673854e-01, 2.65660877e-01, ...,
         2.64687769e-01, 2.64700357e-01, 2.64712944e-01],
        [3.20443379e-01, 3.20427792e-01, 3.20412205e-01, ...,
         3.19669526e-01, 3.19684816e-01, 3.19700106e-01],
        [3.75199927e-01, 3.75181730e-01, 3.75163532e-01, ...,
         3.74651283e-01, 3.74669276e-01, 3.74687268e-01],
...
        [7.34475443e-02, 7.34465877e-02, 7.34456310e-02, ...,
         7.34473770e-02, 7.34483344e-02, 7.34492917e-02],
        [7.35197851e-02, 7.35188405e-02, 7.35178959e-02, ...,
         7.35195927e-02, 7.35205380e-02, 7.35214832e-02],
        [7.35920260e-02, 7.35910934e-02, 7.35901608e-02, ...,
         7.35918083e-02, 7.35927416e-02, 7.35936748e-02]],

       [[7.49902187e-04, 7.49576182e-04, 7.49250177e-04, ...,
         5.88705553e-04, 5.88944985e-04, 5.89184416e-04],
        [7.82970021e-04, 7.82629638e-04, 7.82289255e-04, ...,
         6.15317447e-04, 6.15567787e-04, 6.15818126e-04],
        [8.16037855e-04, 8.15683094e-04, 8.15328333e-04, ...,
         6.41929340e-04, 6.42190588e-04, 6.42451836e-04],
        ...,
        [7.31446797e-02, 7.31437288e-02, 7.31427779e-02, ...,
         7.31445397e-02, 7.31454913e-02, 7.31464429e-02],
        [7.32168718e-02, 7.32159328e-02, 7.32149938e-02, ...,
         7.32167060e-02, 7.32176456e-02, 7.32185852e-02],
        [7.32890638e-02, 7.32881368e-02, 7.32872098e-02, ...,
         7.32888722e-02, 7.32897999e-02, 7.32907275e-02]]])
Coordinates:
  * incidence  (incidence) float64 16.0 16.1 16.2 16.3 ... 65.7 65.8 65.9 66.0
  * wspd       (wspd) float64 0.2 0.3 0.4 0.5 0.6 ... 49.6 49.7 49.8 49.9 50.0
  * phi        (phi) float64 0.0 1.0 2.0 3.0 4.0 ... 177.0 178.0 179.0 180.0
Attributes:
    units:       linear
    resolution:  high
    model:       gmf_cmod5
    pol:         VV

Man can play with **kwargs of GmfModel(Model).

[10]:
# this will directly generate a LUT at high resolution
cmod5.to_lut(**{'wspd_step' : 0.1, 'phi_step' : 1.5, 'inc_step' : 0.5, 'resolution' : 'high'})

DEBUG:xsarsea.windspeed:_raw_lut gmf at res high, inc_step = 0.5 , wspd_step = 0.1, phi_step = 1.5
DEBUG:xsarsea.windspeed:timing _gmf_function : 0.00s. mem: +0.0Mb
DEBUG:xsarsea.windspeed:timing _get_function_for_args : 0.00s. mem: +0.0Mb
DEBUG:xsarsea.windspeed:timing _raw_lut : 0.88s. mem: +46.5Mb
DEBUG:xsarsea.windspeed:lut_resolution gmf_cmod5 from _raw_lut : high
DEBUG:xsarsea.windspeed:desired gmf_cmod5 resolution : high
DEBUG:xsarsea.windspeed:lut gmf_cmod5 already at desired resolution high with exact same params : no interpolation needed
[10]:
<xarray.DataArray 'sigma0_model' (incidence: 101, wspd: 499, phi: 121)>
array([[[2.76625219e-01, 2.76613368e-01, 2.76577849e-01, ...,
         2.75659007e-01, 2.75693536e-01, 2.75705057e-01],
        [3.38003758e-01, 3.37989339e-01, 3.37946121e-01, ...,
         3.37212494e-01, 3.37254922e-01, 3.37269079e-01],
        [3.89481608e-01, 3.89465061e-01, 3.89415467e-01, ...,
         3.89018740e-01, 3.89067907e-01, 3.89084311e-01],
        ...,
        [2.39831755e+00, 2.39802164e+00, 2.39713480e+00, ...,
         2.39706823e+00, 2.39795499e+00, 2.39825087e+00],
        [2.39540775e+00, 2.39511369e+00, 2.39423240e+00, ...,
         2.39416708e+00, 2.39504830e+00, 2.39534233e+00],
        [2.39250269e+00, 2.39221047e+00, 2.39133469e+00, ...,
         2.39127061e+00, 2.39214632e+00, 2.39243852e+00]],

       [[2.14509772e-01, 2.14498628e-01, 2.14465229e-01, ...,
         2.13172772e-01, 2.13204767e-01, 2.13215442e-01],
        [2.65734356e-01, 2.65720599e-01, 2.65679369e-01, ...,
         2.64392328e-01, 2.64432160e-01, 2.64445450e-01],
        [3.09207356e-01, 3.09191405e-01, 3.09143595e-01, ...,
         3.08011650e-01, 3.08058229e-01, 3.08073771e-01],
...
        [7.46565063e-02, 7.46556250e-02, 7.46529837e-02, ...,
         7.46546677e-02, 7.46573107e-02, 7.46581926e-02],
        [7.47289450e-02, 7.47280749e-02, 7.47254670e-02, ...,
         7.47271029e-02, 7.47297125e-02, 7.47305832e-02],
        [7.48014304e-02, 7.48005713e-02, 7.47979963e-02, ...,
         7.47995855e-02, 7.48021621e-02, 7.48030218e-02]],

       [[7.49902187e-04, 7.49608634e-04, 7.48728979e-04, ...,
         5.88322798e-04, 5.88968806e-04, 5.89184416e-04],
        [7.82366939e-04, 7.82060688e-04, 7.81142980e-04, ...,
         6.14415565e-04, 6.15090410e-04, 6.15315645e-04],
        [8.16037855e-04, 8.15718408e-04, 8.14761160e-04, ...,
         6.41511708e-04, 6.42216580e-04, 6.42451836e-04],
        ...,
        [7.31446797e-02, 7.31438235e-02, 7.31412574e-02, ...,
         7.31430181e-02, 7.31455861e-02, 7.31464429e-02],
        [7.32168468e-02, 7.32160015e-02, 7.32134678e-02, ...,
         7.32151782e-02, 7.32177136e-02, 7.32185596e-02],
        [7.32890638e-02, 7.32882291e-02, 7.32857275e-02, ...,
         7.32873889e-02, 7.32898923e-02, 7.32907275e-02]]])
Coordinates:
  * incidence  (incidence) float64 16.0 16.5 17.0 17.5 ... 64.5 65.0 65.5 66.0
  * wspd       (wspd) float64 0.2 0.3 0.4 0.5 0.6 ... 49.6 49.7 49.8 49.9 50.0
  * phi        (phi) float64 0.0 1.5 3.0 4.5 6.0 ... 175.5 177.0 178.5 180.0
Attributes:
    units:       linear
    resolution:  high
    model:       gmf_cmod5
    pol:         VV
[11]:
# this will generate a lut at low resolution and then interpolate at high spedified resolution ('resolution'=None)
cmod5.to_lut(**{'wspd_step' : 0.1, 'phi_step' : 1.5, 'inc_step' : 0.5, 'resolution' : None})
DEBUG:xsarsea.windspeed:_raw_lut gmf at res low, inc_step = 1.0 , wspd_step = 0.2, phi_step = 2.5
DEBUG:xsarsea.windspeed:timing _gmf_function : 0.00s. mem: +0.0Mb
DEBUG:xsarsea.windspeed:timing _get_function_for_args : 0.00s. mem: +0.0Mb
DEBUG:xsarsea.windspeed:timing _raw_lut : 0.14s. mem: +7.0Mb
DEBUG:xsarsea.windspeed:lut_resolution gmf_cmod5 from _raw_lut : low
DEBUG:xsarsea.windspeed:desired gmf_cmod5 resolution : high
DEBUG:xsarsea.windspeed:interp lut gmf_cmod5 to high res
[11]:
<xarray.DataArray 'sigma0_model' (incidence: 101, wspd: 499, phi: 121)>
array([[[2.76625219e-01, 2.76605475e-01, 2.76572622e-01, ...,
         2.75653926e-01, 2.75685864e-01, 2.75705057e-01],
        [3.33053413e-01, 3.33029759e-01, 3.32990396e-01, ...,
         3.32332716e-01, 3.32371423e-01, 3.32394684e-01],
        [3.89481608e-01, 3.89454042e-01, 3.89408169e-01, ...,
         3.89011505e-01, 3.89056982e-01, 3.89084311e-01],
        ...,
        [2.39831755e+00, 2.39782458e+00, 2.39700436e+00, ...,
         2.39693780e+00, 2.39775795e+00, 2.39825087e+00],
        [2.39541012e+00, 2.39492023e+00, 2.39410511e+00, ...,
         2.39403980e+00, 2.39485485e+00, 2.39534469e+00],
        [2.39250269e+00, 2.39201588e+00, 2.39120587e+00, ...,
         2.39114181e+00, 2.39195174e+00, 2.39243852e+00]],

       [[2.21933278e-01, 2.21914926e-01, 2.21884387e-01, ...,
         2.20697518e-01, 2.20726860e-01, 2.20744494e-01],
        [2.70003241e-01, 2.69980957e-01, 2.69943875e-01, ...,
         2.68864172e-01, 2.68900164e-01, 2.68921794e-01],
        [3.18073204e-01, 3.18046988e-01, 3.18003363e-01, ...,
         3.17030826e-01, 3.17073468e-01, 3.17099094e-01],
...
        [7.46590026e-02, 7.46575333e-02, 7.46550882e-02, ...,
         7.46567697e-02, 7.46592164e-02, 7.46606867e-02],
        [7.47314386e-02, 7.47299877e-02, 7.47275734e-02, ...,
         7.47292075e-02, 7.47316235e-02, 7.47330753e-02],
        [7.48038746e-02, 7.48024422e-02, 7.48000585e-02, ...,
         7.48016453e-02, 7.48040306e-02, 7.48054639e-02]],

       [[7.49902187e-04, 7.49413180e-04, 7.48599714e-04, ...,
         5.88227958e-04, 5.88825269e-04, 5.89184416e-04],
        [7.82970021e-04, 7.82459447e-04, 7.81610104e-04, ...,
         6.14818091e-04, 6.15442617e-04, 6.15818126e-04],
        [8.16037855e-04, 8.15505713e-04, 8.14620494e-04, ...,
         6.41408225e-04, 6.42059964e-04, 6.42451836e-04],
        ...,
        [7.31446797e-02, 7.31432533e-02, 7.31408797e-02, ...,
         7.31426402e-02, 7.31450155e-02, 7.31464429e-02],
        [7.32168718e-02, 7.32154633e-02, 7.32131195e-02, ...,
         7.32148303e-02, 7.32171758e-02, 7.32185852e-02],
        [7.32890638e-02, 7.32876733e-02, 7.32853593e-02, ...,
         7.32870205e-02, 7.32893361e-02, 7.32907275e-02]]])
Coordinates:
  * incidence  (incidence) float64 16.0 16.5 17.0 17.5 ... 64.5 65.0 65.5 66.0
  * wspd       (wspd) float64 0.2 0.3 0.4 0.5 0.6 ... 49.6 49.7 49.8 49.9 50.0
  * phi        (phi) float64 0.0 1.5 3.0 4.5 6.0 ... 175.5 177.0 178.5 180.0
Attributes:
    units:       linear
    resolution:  high
    model:       gmf_cmod5
    pol:         VV

It won’t have the same impact on NcLutModel or cmod7Model.

Indeed, these are saved at at desired format with a certain resolution.

By specifying **kwargs, it forced to interpolate the gmf at the desired resolution.

[12]:
# here we specify the exact same params than the saved LUT have
windspeed.get_model('nc_lut_cmodms1ahw').to_lut(**{'wspd_step' : 0.1, 'phi_step' : 1, 'inc_step' : 0.1, 'resolution' : 'high'})
DEBUG:xsarsea.windspeed:lut_resolution nc_lut_cmodms1ahw from _raw_lut : high
DEBUG:xsarsea.windspeed:desired nc_lut_cmodms1ahw resolution : high
DEBUG:xsarsea.windspeed:lut nc_lut_cmodms1ahw already at desired resolution high with exact same params : no interpolation needed
[12]:
<xarray.DataArray 'sigma0_model' (incidence: 331, wspd: 771)>
array([[2.01394183e-04, 2.12451029e-04, 2.23734904e-04, ...,
        2.34376562e-02, 2.34646757e-02, 2.34916925e-02],
       [2.01394183e-04, 2.12451029e-04, 2.23734904e-04, ...,
        2.34376562e-02, 2.34646757e-02, 2.34916925e-02],
       [2.01394183e-04, 2.12451029e-04, 2.23734904e-04, ...,
        2.34376562e-02, 2.34646757e-02, 2.34916925e-02],
       ...,
       [5.04205334e-05, 5.37849910e-05, 5.72563977e-05, ...,
        2.48140448e-02, 2.48560327e-02, 2.48980389e-02],
       [5.04205334e-05, 5.37849910e-05, 5.72563977e-05, ...,
        2.48140448e-02, 2.48560327e-02, 2.48980389e-02],
       [5.04205334e-05, 5.37849910e-05, 5.72563977e-05, ...,
        2.48140448e-02, 2.48560327e-02, 2.48980389e-02]])
Coordinates:
  * incidence  (incidence) float64 17.0 17.1 17.2 17.3 ... 49.7 49.8 49.9 50.0
  * wspd       (wspd) float64 3.0 3.1 3.2 3.3 3.4 ... 79.6 79.7 79.8 79.9 80.0
Attributes:
    units:    linear
    model:    nc_lut_cmodms1ahw
    pol:      VH
[13]:
# here we specify different params than the saved LUT have : interpolation is made
windspeed.get_model('nc_lut_cmodms1ahw').to_lut(**{'wspd_step' : 0.1, 'phi_step' : 1, 'inc_step' : 0.2, 'resolution' : 'high'})
DEBUG:xsarsea.windspeed:lut_resolution nc_lut_cmodms1ahw from _raw_lut : high
DEBUG:xsarsea.windspeed:desired nc_lut_cmodms1ahw resolution : high
DEBUG:xsarsea.windspeed:Even if lut_resolution is already set to high, lut nc_lut_cmodms1ahw needs interpolation to match your desired resolution
DEBUG:xsarsea.windspeed:interp lut nc_lut_cmodms1ahw to high res
[13]:
<xarray.DataArray 'sigma0_model' (incidence: 166, wspd: 771)>
array([[2.01394183e-04, 2.12451029e-04, 2.23734904e-04, ...,
        2.34376562e-02, 2.34646757e-02, 2.34916925e-02],
       [2.01394183e-04, 2.12451029e-04, 2.23734904e-04, ...,
        2.34376562e-02, 2.34646757e-02, 2.34916925e-02],
       [2.01394183e-04, 2.12451029e-04, 2.23734904e-04, ...,
        2.34376562e-02, 2.34646757e-02, 2.34916925e-02],
       ...,
       [5.04205334e-05, 5.37849910e-05, 5.72563977e-05, ...,
        2.48140448e-02, 2.48560327e-02, 2.48980389e-02],
       [5.04205334e-05, 5.37849910e-05, 5.72563977e-05, ...,
        2.48140448e-02, 2.48560327e-02, 2.48980389e-02],
       [5.04205334e-05, 5.37849910e-05, 5.72563977e-05, ...,
        2.48140448e-02, 2.48560327e-02, 2.48980389e-02]])
Coordinates:
  * incidence  (incidence) float64 17.0 17.2 17.4 17.6 ... 49.4 49.6 49.8 50.0
  * wspd       (wspd) float64 3.0 3.1 3.2 3.3 3.4 ... 79.6 79.7 79.8 79.9 80.0
Attributes:
    units:    linear
    model:    nc_lut_cmodms1ahw
    pol:      VH

In practice for wind retrieval

xsarsea.windspeed.invert_from_model can be called with **kwargs.

We can use kwargs to force the use of high resolution Luts.

If possible, use analytical_luts (gmf_…).

Else, better LUTs that has not be interpolated.

We use kwargs = {“wspd_step”: 0.1, “inc_step”: 0.1, “phi_step”: 0.1, “resolution”: “high”}

Adding models (writing your own GMF)

A Geophysical Modeling Function (GMF) is a function that return a simulated sigma0 from wind condition and instrument incidence angle.

A new gmf function is registered with the decorator @xsarsea.windspeed.gmfs.GmfModel.register

To register a new GMF with xsarsea.windspeed, you have to follow the following rules:

  • parameters are (incidence, windspeed, phi)

    • incidence is in degrees

    • windspeed is in m/s

    • phi is wind direction, in degrees, relative to antenna look (0 is downwind in the antenna direction)

      note that phi is mandatory. If the gmf doesn’t need phi, you have to explicitely set phi=None kwarg.

  • all parameters must be float. numpy array are not allowed. xsarsea.windspeed will vectorize the function with numba to allow numpy arrays.

  • allowed units are linear or dB

  • function name must start with gmf_

[14]:
@windspeed.gmfs.GmfModel.register(pol='VH', units='linear', defer=False)
def gmf_dummy(inc, wspd, phi=None):
    a0 = 0.00013106836021008122
    a1 = -4.530598283705591e-06
    a2 = 4.429277425062766e-08
    b0 = 1.3925444179360706
    b1 = 0.004157838450541205
    b2 = 3.4735809771069953e-05


    a = a0 + a1 * inc + a2 * inc ** 2
    b = b0 + b1 * inc + b2 * inc ** 2

    sig = a * wspd ** b

    return sig
DEBUG:xsarsea.windspeed:gmf_dummy doesn't needs phi
DEBUG:xsarsea.windspeed:register model gmf_dummy pol=VH units=linear
 inc_range=[16.0, 66.0] wspd_range=[3.0, 80.0] phi_range=None
 inc_step=0.1 wspd_step=0.1 phi_step=1.0
 inc_step_lr=1.0 wspd_step_lr=0.2 phi_step_lr=2.5

Note : if defer = True, you will have to use again xsarsea.windspeed.gmfs.GmfModel.activate_gmfs_impl() to activate the new model

[15]:
gmf_dummy
[15]:
<function __main__.gmf_dummy(inc, wspd, phi=None)>
[16]:
windspeed.get_model('gmf_dummy').to_lut()
DEBUG:xsarsea.windspeed:_raw_lut gmf at res low, inc_step = 1.0 , wspd_step = 0.2, phi_step = 2.5
DEBUG:xsarsea.windspeed:timing _gmf_function : 0.13s. mem: +0.0Mb
DEBUG:xsarsea.windspeed:timing _gmf_function : 0.36s. mem: +0.5Mb
DEBUG:xsarsea.windspeed:timing _get_function_for_args : 0.36s. mem: +0.5Mb
DEBUG:xsarsea.windspeed:timing _raw_lut : 0.37s. mem: +0.5Mb
DEBUG:xsarsea.windspeed:lut_resolution gmf_dummy from _raw_lut : low
DEBUG:xsarsea.windspeed:desired gmf_dummy resolution : high
DEBUG:xsarsea.windspeed:interp lut gmf_dummy to high res
[16]:
<xarray.DataArray 'sigma0_model' (incidence: 501, wspd: 771)>
array([[0.00035074, 0.00036817, 0.00038559, ..., 0.04331674, 0.04339647,
        0.0434762 ],
       [0.00034939, 0.00036676, 0.00038413, ..., 0.04322394, 0.04330352,
        0.04338311],
       [0.00034805, 0.00036536, 0.00038266, ..., 0.04313113, 0.04321058,
        0.04329002],
       ...,
       [0.00018198, 0.0001933 , 0.00020462, ..., 0.07055381, 0.07071459,
        0.07087537],
       [0.00018309, 0.00019448, 0.00020587, ..., 0.07117419, 0.07133645,
        0.07149872],
       [0.00018419, 0.00019566, 0.00020713, ..., 0.07179456, 0.07195831,
        0.07212207]])
Coordinates:
  * incidence  (incidence) float64 16.0 16.1 16.2 16.3 ... 65.7 65.8 65.9 66.0
  * wspd       (wspd) float64 3.0 3.1 3.2 3.3 3.4 ... 79.6 79.7 79.8 79.9 80.0
Attributes:
    units:       linear
    resolution:  high
    model:       gmf_dummy
    pol:         VH

HH Luts

We also created HH LUTS using CMODs and Polarization Ratio (PR) Models :

  • PR “mouche1” from Mouche, A., Hauser, D., Kudryavtsev, V., and Daloze, J.-F. (2005). Multi-polarisation ocean radar cross-section from envisat asar observations, airborne polarimetric radar measurements and empirical or semiempirical models &

  • PR “zhang” from Zhang, B., Perrie, W., and He, Y. (2011). Wind speed retrieval from radarsat-2 quad-polarization images using a new polarization ratio model. Journal of Geophysical Research: Oceans.

We simply used this equation and created the NcLutModels

\[nrcs_{HH} = \frac{nrcs_{VV}}{PR}\]

For CMOD5n, we created high resolution LUTS (0.1m/s, 0.1°, 1°).

For CMOD7, we created low resolution LUTS (base cmod7 resolution from files) and high resolution (from file + interpolation).

Then we can directly use these luts

[17]:
# loading low resolution & interpolating
windspeed.get_model('nc_lut_gmf_cmod7_Rlow_hh_mouche1').to_lut(**{'wspd_step' : 0.1, 'phi_step' : 1, 'inc_step' : 0.1, 'resolution' : 'high'})
# or
# loading high resolution [not computed in the doc cause nc_lut_gmf_cmod7_Rhigh_hh_mouche1 is too big]
#windspeed.get_model('nc_lut_gmf_cmod7_Rhigh_hh_mouche1').to_lut(**{'wspd_step' : 0.1, 'phi_step' : 1, 'inc_step' : 0.1, 'resolution' : 'high'})
DEBUG:xsarsea.windspeed:lut_resolution nc_lut_gmf_cmod7_Rlow_hh_mouche1 from _raw_lut : low
DEBUG:xsarsea.windspeed:desired nc_lut_gmf_cmod7_Rlow_hh_mouche1 resolution : high
DEBUG:xsarsea.windspeed:interp lut nc_lut_gmf_cmod7_Rlow_hh_mouche1 to high res
[17]:
<xarray.DataArray 'sigma0_model' (incidence: 501, wspd: 499, phi: 181)>
array([[[1.23932622e-01, 1.23926588e-01, 1.23920554e-01, ...,
         1.23092110e-01, 1.23097779e-01, 1.23103447e-01],
        [1.46205354e-01, 1.46198246e-01, 1.46191139e-01, ...,
         1.45058429e-01, 1.45065097e-01, 1.45071765e-01],
        [1.68478086e-01, 1.68469905e-01, 1.68461724e-01, ...,
         1.67024748e-01, 1.67032415e-01, 1.67040082e-01],
        ...,
        [2.43469425e+00, 2.43433781e+00, 2.43398138e+00, ...,
         2.43389735e+00, 2.43425406e+00, 2.43461077e+00],
        [2.41908574e+00, 2.41874062e+00, 2.41839550e+00, ...,
         2.41831848e+00, 2.41866374e+00, 2.41900900e+00],
        [2.40347723e+00, 2.40314343e+00, 2.40280962e+00, ...,
         2.40273960e+00, 2.40307341e+00, 2.40340722e+00]],

       [[1.18065426e-01, 1.18059547e-01, 1.18053668e-01, ...,
         1.17228098e-01, 1.17233613e-01, 1.17239128e-01],
        [1.39484671e-01, 1.39477733e-01, 1.39470795e-01, ...,
         1.38345545e-01, 1.38352045e-01, 1.38358544e-01],
        [1.60903916e-01, 1.60895919e-01, 1.60887922e-01, ...,
         1.59462993e-01, 1.59470477e-01, 1.59477960e-01],
...
        [2.46566096e-03, 2.46563232e-03, 2.46560369e-03, ...,
         2.46563538e-03, 2.46566417e-03, 2.46569296e-03],
        [2.47024784e-03, 2.47021951e-03, 2.47019118e-03, ...,
         2.47022218e-03, 2.47025072e-03, 2.47027925e-03],
        [2.47483472e-03, 2.47480669e-03, 2.47477867e-03, ...,
         2.47480899e-03, 2.47483726e-03, 2.47486554e-03]],

       [[3.10420750e-06, 3.10294843e-06, 3.10168936e-06, ...,
         3.13727881e-06, 3.13855226e-06, 3.13982572e-06],
        [5.30036492e-06, 5.29834257e-06, 5.29632021e-06, ...,
         5.29234646e-06, 5.29436955e-06, 5.29639265e-06],
        [7.49652234e-06, 7.49373670e-06, 7.49095107e-06, ...,
         7.44741411e-06, 7.45018684e-06, 7.45295958e-06],
        ...,
        [2.41565417e-03, 2.41562627e-03, 2.41559837e-03, ...,
         2.41562971e-03, 2.41565785e-03, 2.41568600e-03],
        [2.42016627e-03, 2.42013862e-03, 2.42011097e-03, ...,
         2.42014167e-03, 2.42016959e-03, 2.42019752e-03],
        [2.42467837e-03, 2.42465097e-03, 2.42462358e-03, ...,
         2.42465363e-03, 2.42468133e-03, 2.42470904e-03]]])
Coordinates:
  * incidence  (incidence) float64 16.0 16.1 16.2 16.3 ... 65.7 65.8 65.9 66.0
  * wspd       (wspd) float64 0.2 0.3 0.4 0.5 0.6 ... 49.6 49.7 49.8 49.9 50.0
  * phi        (phi) float64 0.0 1.0 2.0 3.0 4.0 ... 177.0 178.0 179.0 180.0
Attributes:
    units:       linear
    model:       nc_lut_gmf_cmod7_Rlow_hh_mouche1
    resolution:  high
    pol:         HH

Model comparison

This example function can be used to compare models.

[18]:

def model_compare(compare_models): luts = [ windspeed.get_model(name).to_lut(units='dB') for name in compare_models] if 'phi' not in luts[0].dims: kdims=['incidence'] dim_range=dict(incidence=(17,50)) else: kdims=['incidence', 'phi'] dim_range=dict(incidence=(17,50), phi=(0,360)) def model_curve(incidence, phi=None): if 'phi' not in luts[0].dims: sel = dict(incidence=incidence) else: sel = dict(phi=phi, incidence=incidence) return hv.Overlay( [ hv.Curve(lut.sel(**sel, method='nearest'),'wspd','sigma0', label=lut.attrs['model']) for lut in luts ] ) dmap = hv.DynamicMap(model_curve, kdims=kdims).opts(height=600, width=600) return dmap.redim.range(**dim_range)
[19]:
windspeed.available_models(pol='VH')
[19]:
alias pol model
gmf_rs2_v2 rs2_v2 VH <GmfModel('gmf_rs2_v2') pol=VH>
gmf_s1_v2 s1_v2 VH <GmfModel('gmf_s1_v2') pol=VH>
gmf_rcm_noaa rcm_noaa VH <GmfModel('gmf_rcm_noaa') pol=VH>
gmf_dummy dummy VH <GmfModel('gmf_dummy') pol=VH>
nc_lut_cmodms1ahw cmodms1ahw VH <NcLutModel('nc_lut_cmodms1ahw') pol=VH>
[20]:
model_compare([ 'gmf_dummy', 'nc_lut_cmodms1ahw'] )
DEBUG:xsarsea.windspeed:_raw_lut gmf at res low, inc_step = 1.0 , wspd_step = 0.2, phi_step = 2.5
DEBUG:xsarsea.windspeed:timing _gmf_function : 0.00s. mem: +0.0Mb
DEBUG:xsarsea.windspeed:timing _get_function_for_args : 0.00s. mem: +0.0Mb
DEBUG:xsarsea.windspeed:timing _raw_lut : 0.00s. mem: +0.0Mb
DEBUG:xsarsea.windspeed:lut_resolution gmf_dummy from _raw_lut : low
DEBUG:xsarsea.windspeed:desired gmf_dummy resolution : high
DEBUG:xsarsea.windspeed:interp lut gmf_dummy to high res
DEBUG:xsarsea.windspeed:lut_resolution nc_lut_cmodms1ahw from _raw_lut : high
DEBUG:xsarsea.windspeed:desired nc_lut_cmodms1ahw resolution : high
DEBUG:xsarsea.windspeed:lut nc_lut_cmodms1ahw already at desired resolution high with exact same params : no interpolation needed
[20]:
[21]:
windspeed.available_models(pol='VV')
[21]:
alias pol model
gmf_cmod7 cmod7 VV <Cmod7Model('gmf_cmod7') pol=VV>
gmf_cmod5 cmod5 VV <GmfModel('gmf_cmod5') pol=VV>
gmf_cmod5n cmod5n VV <GmfModel('gmf_cmod5n') pol=VV>
gmf_cmodifr2 cmodifr2 VV <GmfModel('gmf_cmodifr2') pol=VV>
[22]:
model_compare([ 'gmf_cmod7', 'gmf_cmod5n', 'gmf_cmod5'])
INFO:xsarsea.windspeed:load gmf lut from /tmp/cmod7_and_python_script
DEBUG:xsarsea.windspeed:lut_resolution gmf_cmod7 from _raw_lut : low
DEBUG:xsarsea.windspeed:desired gmf_cmod7 resolution : high
DEBUG:xsarsea.windspeed:interp lut gmf_cmod7 to high res
DEBUG:xsarsea.windspeed:_raw_lut gmf at res low, inc_step = 1.0 , wspd_step = 0.2, phi_step = 2.5
DEBUG:xsarsea.windspeed:timing _gmf_function : 0.28s. mem: +2.9Mb
DEBUG:xsarsea.windspeed:timing _gmf_function : 0.52s. mem: +5.1Mb
DEBUG:xsarsea.windspeed:timing _get_function_for_args : 0.52s. mem: +5.1Mb
DEBUG:xsarsea.windspeed:timing _raw_lut : 0.64s. mem: +5.1Mb
DEBUG:xsarsea.windspeed:lut_resolution gmf_cmod5n from _raw_lut : low
DEBUG:xsarsea.windspeed:desired gmf_cmod5n resolution : high
DEBUG:xsarsea.windspeed:interp lut gmf_cmod5n to high res
DEBUG:xsarsea.windspeed:_raw_lut gmf at res low, inc_step = 1.0 , wspd_step = 0.2, phi_step = 2.5
DEBUG:xsarsea.windspeed:timing _gmf_function : 0.00s. mem: +0.0Mb
DEBUG:xsarsea.windspeed:timing _get_function_for_args : 0.00s. mem: +0.0Mb
DEBUG:xsarsea.windspeed:timing _raw_lut : 0.13s. mem: +0.0Mb
DEBUG:xsarsea.windspeed:lut_resolution gmf_cmod5 from _raw_lut : low
DEBUG:xsarsea.windspeed:desired gmf_cmod5 resolution : high
DEBUG:xsarsea.windspeed:interp lut gmf_cmod5 to high res
[22]: