interp module

Interpolation for Ocean and Atmospheric Data

Routines and methods for interpolation of ocean and atmospheric model field, either on regular grids or rotated, multi-pole grids. For the ocean, is set for the CMCC Ocean model at the nominal resolution of 0.25 respectively.

The staggering requires different interpolators operators for scalar points at T-points and vector quantities carried at (u,v) points. For the moment a simple interpolation is carried out but a more accurate vector interpolation is under development.

The interpolation is obtained by triangulation of the starting grid and seaprate interpolation to the new grid. The weights are preserved and they can be used for repeated application of the same set of grids.

Classes

Atmosphere_Interpolator
Ocean_Interpolator
class interp.Atmosphere_Interpolator(grid, option='linear')[source]

Bases: object

This class creates weights for interpolation of atmospheric fields. No mask is used.

This class create an interpolator operator to a target grid tgt_grd. The interpolator then can be used to perform the actual interpolation.

The target grid must be a xarray DataArray with variables lat and lon

Parameters:
  • grid (str) -- Choice of output grids

    • 1x1 -- Regular 1 degree

    • 025x025 -- Coupled model grid, nominally 0.25,

  • option (str) -- 'linear', interpolation method

Variables:
  • name (str) -- Name of the interpolator

  • grid (str) -- Option for the grid

  • option (str) -- Interpolation method

Notes

It is a thin wrapper around xarray interp_like method.

Examples

Create the weights for interpolation

>>> w= zint.Atmosphere_Interpolator('1x1','linear')

Interpolate temperature

>>> target_xarray=w.interp_f(src_xarray)
choice

Interpolation method selected linear or nearest

Type:

str

name

Name of the Interpolator

Type:

str

tgt

Target grid

Type:

xarray

interp_scalar(xdata)[source]

Perform interpolation to the target grid. This methods can be used for scalar quantities.

Parameters:

xdata (xarray) -- 2D array to be interpolated, it must be on the src_grid

Returns:

out -- Interpolated xarray on the target grid

Return type:

xarray

class interp.Ocean_Interpolator(src_grid_name, tgt_grid_name, level=1, verbose=False, window=3, period=1)[source]

Bases: object

This class creates weights for interpolation of ocean fields.

This class create an interpolator operator from a source grid src_grid to a target grid tgt_grd. The interpolator then can be used to perform the actual interpolation. The model uses an Arakawa C-grid, that is shown in the following picture. The f-points correspond to the points where the Coriolis terms are carried.

_images/NEMOgrid.png

The Arakawa C-grid used in the ocean model show also the ordering of the points, indicating which points correspond to the (i,j) index.

The source grid must be a xarray DataSet containing coordinates latitude and longitude.

The target grid must be a xarray DataArray with variables lat and lon

Border land points at all levels are covered by a convolution value using a window that can be changed in sea_over_land.

Works only on single DataArray

Parameters:
  • src_grid (xarray) -- Source grid

  • tgt_grid (xarray) -- Target grid

  • level (float) -- Depth to generate the interpolator

  • window (int) -- Window for sea over land

  • period (int) -- Minimum number of points in the sea-over-land process

  • verbose (bool) -- Lots of output

Variables:
  • tgt_grid -- Target grid

  • mask -- Mask of the target grid

  • vt -- Weights

  • wt -- Weights

  • mdir -- Directory for masks files

  • ingrid -- Input Grid, src_grid_name

  • outgrid -- Output Grid tgt_grid_name

  • window -- Window for convolution Sea-Over-Land (default 3)

  • period -- Minimum number of points into the window (default 1)

  • T_lon -- Longitudes of input T-mask

  • T_lat -- Latitudes of input T-mask

  • U_lon -- Longitudes of input U-mask

  • U_lat -- Latitudes of input U-mask

  • V_lon -- Longitudes of input V-mask

  • V_lat -- Latitudes of input V-mask

  • tangle -- Angles of the T points of the input grid

  • mask_reg -- Mask of the Target grid

  • cent_long -- Central Longitude of the Target grid

  • name -- Name of the Interpolator Object

  • level -- Level of the Interpolator Object

Interp_T :

Interpolate Scalar quantities at T points

Interp_UV :

Interpolate Vector Velocities at (U,V)

mask_sea_over_land :

Mask border point for Sea over land

UV_sea_over_land :

Fill U,V values over land

to_file :

Writes interpolator object to file (pickled format)

Examples

Create the weights for interpolation

>>> w= zint.Ocean_Interpolator(src_grid,tgt_grid)

Interpolate temperature

>>> target_xarray=w.interp_T(src_xarray,method='linear')

Interpolate U,V

>>> target_xarray=w.interp_UV(U_xarray,V_xarray,method='linear')
mdir
ingrid
outgrid
window
period
T_lon
T_lat
U_lon
U_lat
V_lon
V_lat
name
level
tangle
mask
masktb
maskub
masku
maskvb
maskv
latlon
sea_index
masT_vec
tri_sea_T
sea_index_U
tri_sea_U
sea_index_V
tri_sea_V
mask_reg
cent_long
latlon_reg
sea_index_reg
regmask_vec
interp_T(xdata, method='linear')[source]

Perform interpolation for T Grid point to the target grid. This methods can be used for scalar quantities.

Parameters:
  • xdata (xarray) -- 2D array to be interpolated, it must be on the src_grid

  • method (str) --

    Method for interpolation
    • 'linear' , Use linear interpolation

    • 'nearest' , use nearest interpolation

Returns:

out -- Interpolated xarray on the target grid

Return type:

xarray

interp_UV(udata, vdata, method='linear')[source]

Perform interpolation for U,V Grid point to the target grid. This methods can be used for vector quantities.

The present method interpolates the U,V points to the T points, rotates them and then interpolates to the target grid.

Parameters:

udata,vdata (xarray) -- 2D array to be interpolated, it must be on the src_grid

Returns:

out -- Interpolated xarray on the target grid

Return type:

xarray

to_file(filename)[source]

This method writes to file the interpolating object in pickled format

mask_sea_over_land(mask)[source]

Mask border point for Sea over land.

Sea over land is obtained by forward and backward filling NaN land value with an arbitrary value, then athey are masked to reveal only the coastal points.

Parameters:

mask -- original mask

Returns:

border mask

Return type:

border

fill_sea_over_land(U, u_mask)[source]

Put values Sea over land.

Using the mask of the border points, the border points are filled with the convolution in 2D, using a window of width window, here The min_periods value is controlling the minimum number of points within the window that is necessary to yield a result.

They can be fixed as attributes of the interpolator

Parameters:
  • U -- Field to be treated

  • u_mask -- Mask for the field

Returns:

Filled array

Return type:

U

interp.get_sea(maskT)[source]

Obtain indexed coordinates for sea points

interp.from_file(file)[source]

Read interpolator object from file