2.3.8.6. Content¶
Regridding utilities
See also
Tutorials: Regrillage 1D et 2D
-
class
CDATRegridder
(fromgrid, togrid, missing=None, order=None, mask=None, **keywords)[source]¶ Bases:
object
Regridding using CDAT regridders
Note
This code is adapted from the
regrid()
method of MV2 arrays.Params: - fromgrid: Input grid, or variable with a grid.
- togrid: Output grid.
- tool, optional: One of
"esmf"
,"libcf"
and"regrid2"
. - method, optional: One of
"linear"
,"path"
,"conservative"
and"cellave"
.
-
regrid
(vari, weidstfracs=None, check_mask=None, csvhack=False, **keywords)[source]¶ Regrid the variable
Params: vari: Input variable.
weidstfracs, optional: Specify how to apply weights computed using destination fractions. Divide if <0, else multiply. It is not used with linear and match methods.
check_mask, optional: Mask point using masked data (the algo interpolate the mask and check level 0.999). MUST BE IMPROVED! If None, mask is checked if tool!=’regrid2’
csvhack, optional: Hack to prevent a bug with conservative-like method of the ESMF regridder. Use it if you have strange result in the output most right longitude.
Warning
The algo mask the must right longitude of input data before processing.
-
class
CurvedInterpolator
(fromgrid, topts, g2g=False)[source]¶ Bases:
object
Interpolator from a curved grid to random points or another grid
Params: - fromgrid: Input grid, or variable with a grid.
- topts: Output coordinates or grid. It a tuple of coordinates, it is assumed that it refers to random points, and NOT a grid.
- g2g, optional: Force the interpretation of
topts
as a grid or axes of a grid.
Examples: >>> interpolator = CurvedInterpolator(ssti.getGrid(), (lono, lato)) # lono/lato 1D >>> interpolator = CurvedInterpolator(ssti.getGrid(), grido) # grid >>> ssto = interpolator(ssti)
-
regrid
(vari, method='bilinear')¶
-
valid_methods
= ['bilinear', 'nearest', 'dstwgt']¶
-
class
GridData
(xi, yi, ggo, nl=False, ext=False, geo=None, method='nat', sub=30, margin=0.5, compress=False, **kwargs)[source]¶ Bases:
object
2D interpolator from a randomly spaced sample on a regular grid
Possible algorithms:
Natural Neighbor using nat.Natgrid:
2D splines using css.css.Cssgrid
Parameters: xi: Input 1D X positions.
yi: Input 1D Y positions.
ggo: Output grid. Can either (xo,yo), a cdms grid or a cdms variable with a grid.
method: Interpolator type, either ‘nat’ (Natural Neighbors) or ‘css’ (=’splines’ using splines) [default: ‘nat’]
nl: Nonlinear interpolator (usually gives better results) [default: False]
ext: Extrapolate value outsite convex hull [‘nat’ only, default: False]
mask: Mask to apply to output data [default: None]
compress: If
True
, interolate only unmasked data, and thus does not try guess the best mask (that’s more efficient but very bad if data are masked!).sub: Size of blocks for subblocking [“nat” only]
margin: Margin around ouput grid (or block) to select input data. Value is relative to X and Y extent.
Other keywords are set as attribute to the interpolator instance ; to get the list of parameters:
>>> import nat ; nat.printParameterTable() >>> import css ; css.printParameterTable()
Example: >>> r = GridData(gridi, grido, method='nat', ext=False, margin=.7) >>> varo1 = r(vari1) >>> varo2 = r(vari2)
-
regrid
(zi, missing_value=None, **kwargs)¶ Interpolate zi on output grid
- zi: At least a 1D array.
-
rgrd
(zi, missing_value=None, **kwargs)¶ Interpolate zi on output grid
- zi: At least a 1D array.
-
class
GriddedMerger
(grid, id=None, long_name=None, units=None, **kwargs)[source]¶ Bases:
object
Merge several gridded variables onto a grid
Params: - grid: Output grid
- id, optional: Output id
- long_name, optional: Output long name
- units, optional: Output units
- Other keywords are set a output attributes
Example: >>> gm = GriddedMerger(mygrid) >>> gm += var1 >>> gm.append(var2) >>> gm += var3 >>> gm -= var3 >>> gm.insert(0, var3) >>> print len(gm) 3 >>> print gm .... >>> myvar = gm.merge(res_ratio=.4, pad=3)
-
merge
(res_ratio=0.5, pad=5, **kwargs)[source]¶ Merge all the variables on to a grid
- grid: Out put grid or axes.
- res_ratio: Resolution ratio for choosing between cell
averaging and bilinear interpolation (see:
regrid_method()
). - regrid_<kwparam>: <kwparam> is passed to
regrid2d()
for interpolation.
-
cargen
(xi, yi, zi, ggo, mask=None, compress=False, missing_value=None, **kwargs)[source]¶ Interpolator from IFREMER
Params: - xi: Input 1D X positions.
- yi: Input 1D Y positions.
- ggo, optional: Output grid. Can either (xo,yo), a cdms grid or a cdms variable with a grid.
- mask, optional: Mask to apply to output data [default: None]
-
cellave1d
(vari, axo, conserv=False, **kwargs)[source]¶ Cell averaging or conservative regridding along an axis
Params: - vari: Input cdms array
- axo: Output cdms axis
- axis, optional: Axis on wich to operate
- conservative, optional: If True, regridding is conservative
- Other keywords are passed to
regrid1d()
Note
This is an wrapper to
regrid1d()
usingcellave
orconservative
as a default method. See its help for more information.
-
cellave2d
(vari, ggo, **kwargs)[source]¶ Shortcut to
regrid2d()
call withmethod="cellave"
-
cellerr1d
(vari, axo, erri, errl=None, **kwargs)[source]¶ Cell averaing with weights based on errors
Params: - vari: Input cdms array
- axo: Output cdms axis
- erri: Input measurement errors
- errl, optional: Input lag error relative to lag
- axis, optional: Axis on wich to operate
- Other keywords are passed to
regrid1d()
Note
This is an wrapper to
regrid1d()
usingcellerr
method. See its help for more information.
-
cubic1d
(vari, axo, **kwargs)[source]¶ Cubic interpolation along an axes
Params: - vari: Input cdms array
- axo: Output cdms axis
- axis, optional: Axis on wich to operate
- Other keywords are passed to
regrid1d()
Note
This is an wrapper to
regrid1d()
usingcubic
as a default method. See its help for more information.
-
extend1d
(var, ext=0, mode=None, axis=-1, copy=False, num=False)[source]¶ Extrapolate data on an axis
Params: var: Array.
ext, optional: Size of extension. If a tuple, it gives left and right extensions, else the same for both.
mode, optional: Interpolation mode for boundary point outside initial positions.
None
:"extrap"
if axis else"same"
"extrap"
or"linear"
: Linear extrapolation."same"
: Constant extrapolation."cyclic"
: Cyclic extrapolation."masked"
: Masked.
axis, optional: Axis on which to operate.
copy, optional: Always input copy data.
-
extend2d
(vari, iext=0, jext=0, mode=None, copy=False)[source]¶ Interpolate data on an grid shifted by an half cell size in X and Y
X and Y are supposed to be the -1 and -2 axes of var.
Params: var: Array.
i/jext, optional: Size of extension along i/j. If a tuple,it gives left and right extensions, else the same for both.
mode, optional: Interpolation mode for boundary point outside initial positions.
None
:"extrap"
if axis else"same"
"extrap"
or"linear"
: Linear extrapolation."same"
: Constant extrapolation."masked"
: Masked.
copy, optional: Always input copy data.
-
extendgrid
(gg, iext=0, jext=0, mode='extrap')[source]¶ Extrapolate a grid
Params: gg: cdms2 grid.
i/jext: Size of extrapolation along i/j.
mode, optional: Interpolation mode for boundary point outside initial positions.
"extrap"
: Linear extrapolation."same"
: Constant extrapolation."masked"
: Masked.
-
fill1d
(vari, axis=0, method='linear', maxgap=0)[source]¶ Fill missing values of a 1D array using interpolation.
Params: - vari: Input
cdms2
variable - axis: Axis number on which filling is performed
- method: Interpolation method (see
interp1d()
) - maxgap: Maximal size of filled gaps (in steps)
Example: >>> fill1d(vari, axis=2, method='cubic', maxgap=5)
Return: Filled
cdms2
variable similar to input one- vari: Input
-
fill1d
(vari, axis=0, method='linear', maxgap=0)[source] Fill missing values of a 1D array using interpolation.
Params: - vari: Input
cdms2
variable - axis: Axis number on which filling is performed
- method: Interpolation method (see
interp1d()
) - maxgap: Maximal size of filled gaps (in steps)
Example: >>> fill1d(vari, axis=2, method='cubic', maxgap=5)
Return: Filled
cdms2
variable similar to input one- vari: Input
-
fill1d2
(vi, axis=0, k=1, padding=None, clip=False, min_padding=None, method='linear')[source]¶ Fill missing values of a 1D array using spline interpolation.
- vi: Input cdms variable
- k: Order of splines [default: 1 = linear]
- padding: Padding around an gap defining where on which part of the sample we must fit splines [default: max([min_padding,len(gap)*5])]
- min_padding: See padding [default: k]
- method: See
interp1d()
[default: linear]
Return: Filled cdms2
variable
-
fill2d
(var, xx=None, yy=None, mask=None, copy=True, **kwargs)[source]¶ Fill missing value of 2D variable using inter/extrapolation
- var: A cdms 2D variable.
- xx/yy: Substitutes for axis coordinates [default: None]
- Other keywords are passe to
griddata()
-
grid2xy
(vari, xo, yo, zo=None, to=None, zi=None, method='linear', outaxis=None, distmode='haversine')[source]¶ Interpolate gridded data to ramdom positions
Params: vari: Input cdms variable on a grid
xo: Output longitudes
yo: Output latitudes
method, optional: Interpolation method
nearest
: Nearest neighborlinear
: Linear interpolation
zo, optional: Output depths (negative in the ocean).
to, optional: Output times.
zi, optional: Input depths when variable in space.
outaxis, optional: Output spatial axis
- A cdms2 axis.
None
or'auto'
: Longitudes or latitudes depending on the range if coordinates are monotonic, else'dist'
.'lon'
or'x'
: Longitudes.'lat'
or'y'
: Latitudes.'dist'
or'd'
: Distance in km.
distmode, optional: Distance computation mode. See
get_distances()
.
-
griddata
(xi, yi, zi, ggo, method='linear', cgrid=False, cache=None, proj=True, **kwargs)[source]¶ Interpolation in one single shot using GridData
Params: - xi: 1D input x coordinates.
- yi: 1D input y coordinates (same length as xi).
- zi: 1D input values (same length as xi).
- method, optional: Method of interpolation, within (‘nat’, ‘css’, ‘carg’, ‘krig’) [default: ‘carg’]
- cgrid, optional: Output on a C-grid at U- and V-points deduced from ggo [default: False]. Not available for ‘carg’ and ‘krig’ methods.
See also:
-
interp1d
(vari, axo, method='linear', **kwargs)[source]¶ Linear or cubic interpolation along an axes
Params: - vari: Input cdms array
- axo: Output cdms axis
- axis, optional: Axis on wich to operate
- Other keywords are passed to
regrid1d()
Note
This is an wrapper to
regrid1d()
usinglinear
as a default method. See its help for more information.
-
interp2d
(vari, ggo, method='interp', **kwargs)[source]¶ Shortcut to
regrid2d()
call withmethod="interp"
by default
-
krig
(xi, yi, zi, ggo, mask=True, proj=True, missing_value=None, **kwargs)[source]¶ Kriging interpolator to a grid
Params: - xi: Input 1D X positions.
- yi: Input 1D Y positions.
- zi: Input N-D with last dim as space.
- ggo, optional: Output grid. Can either (xo,yo), a cdms grid or a cdms variable with a grid.
- mask, optional: Mask to apply to output data [default: None]
-
nearest1d
(vari, axo, **kwargs)[source]¶ Interpolation along an axes
Params: - vari: Input cdms array
- axo: Output cdms axis
- axis, optional: Axis on wich to operate
- Other keywords are passed to
regrid1d()
Note
This is an wrapper to
regrid1d()
usingnearest
as a default method. See its help for more information.
-
refine
(vari, factor, geo=True, smoothcoast=False, noaxes=False)[source]¶ Refine a variable on a grid by a factor
Params: - vari: 1D or 2D variable.
- factor: Refinement factor > 1
-
regrid1d
(vari, axo, method='auto', axis=None, axi=None, iaxo=None, iaxi=None, xmap=None, xmapper=None, mask_thres=0.5, extrap=0, erri=None, errl=None, geterr=False)[source]¶ Interpolation along one axis
Params: vari: Input cdms array.
axo: Output cdms axis or array. It can be of any dimensions.
method:
"nearest"|0
: Nearest neighbor"linear"|2
: Linear interpolation"cubic"|2
: Cubic interpolation (not used, switched to3
)"hermit"|3
: Cubic hermit interpolation"cellerr"|4
: Cell averaging based on errors"cellave"|-1
: Cell averaging"conserv"|-1
: Conservative cel averaging (likecellave
but with integral preserved)
axis, optional: Dimension (int) on which the interpolation is performed. If not specified, it is guessed from the input and output axis types, or set to
0
.axi, optional: Input axis. It defaults to the axis-th axis of
vari
. Likeaxo
, it can be of any dimensions.iaxo, optional: Dimension of
axo
on which the interpolation is performed whenaxo
has more than one dimension.iaxi, optional: Same as
iaxo
but foraxi
.mask_thres, optional: Time steps when interpolated mask is greater than this value are masked.
extrap, optional: Extrapolate outside input grid when the “nearest” method is used:
0
orFalse
: No extrapolation.-1
or"min"
, or"bottom"
, or"lower"
, or"first"
: Extrapolate toward first values of the axis.1
or"max"
, or"top"
, or"upper"
, or"last"
: Extrapolate toward last values of the axis.2
or"both"
: Extrapolate toward both first and last values.
erri, optional: Input “measurement” errors with the same shape as input variable.
errl, optional: Derivative of lag error with respect to lag. Note that the lag must expressed in days for time axes. If positive, it based on quadratic errors, else on error itself. Estimate for instance it using the slope of a linear regression. It is usually varying in space and constant in time.
geterr, optional: When method is “cellerr”, also return the error along with the variable.
Examples: >>> varo = regrid1d(vari, taxis, method='linear') # interpolation in time >>> varo = regrid1d(vari, zo, axis=1) # Z interpolation on second axis >>> varo = regrid1d(vari, zzo, iaxo=1, axi=zzi, iaxi=1) # sigma to sigma
Note
Cubic method, use “linear” interpolation when less than 4 valid points are available. Linear interpolation uses “nearest” interpolation when less than 2 points are available.
-
regrid1d_method_name
(method, raiseerr=True)[source]¶ Check the method name and return its generic name
-
regrid1dold
(vari, axo, method='auto', axis=None, xmap=None, xmapper=None, mask_thres=0.5, extrap=0)[source]¶ Interpolation along one axis
Params: vari: Input cdms array.
axo: Output cdms axis.
method:
"nearest"
: Nearest neighbor"linear"
: Linear interpolation"cubic"
: Cubic interpolation"cellave"
: Cell averaging"conserv"
: Conservative cel averaging (likecellave
but with integral preserved)
axis, optional: Axis (int) on which to operate. If not specified, it is guessed from the input and output axis types, or set to
0
.xmap, optional: Integer or tuple that specify on which axes input axis is varying.
xmapper, optional: Array that specify values of input axis along axes specified by
xmap
. It is an array of size(...,len(var.getAxis(xmap[-2])), len(var.getAxis(xmap[-1])), len(var.getAxis(axis))]
.mask_thres, optional: Time steps when interpolated mask is greater than this value are masked.
extrap, optional: Extrapolate outside input grid when the “nearest” method is used:
0
orFalse
: No extrapolation.-1
or"min"
, or"bottom"
, or"lower"
, or"first"
: Extrapolate toward first values of the axis.1
or"max"
, or"top"
, or"upper"
, or"last"
: Extrapolate toward last values of the axis.2
or"both"
: Extrapolate toward both first and last values.
Note
Cubic method, use “linear” interpolation when less than 4 valid points are available. Linear interpolation uses “nearest” interpolation when less than 2 points are available.
-
regrid2d
(vari, ggo, method='auto', tool=None, rgdr=None, getrgdr=False, erode=False, **kwargs)[source]¶ Regrid a variable from a regular grid to another
If the input or output grid is curvilinear and
method
is set to"linear"
,"cellave"
or"conserv"
, theCDATRegridder
is used.Params: vari: Variable cdms on regular grid
ggo: Tuple of (x,y) or a cdms grid or a cdms variable with a grid
method, optional: One of:
"auto"
: method guessed betweenlinear
andcellave
according to resolution of input and output grid (seeregrid_method()
)"nearest"
: nearest neighbour"linear"
or"bilinear"
: bilinear interpolation (low res. to high res.)"dstwgt"
: distance weighting between the four nearest grid points (low res. to high res.)"patch"
: patch recovery interpolation (low res. to high res.)"cellave"
: weighted regridding based on areas of cells (high res. to low res.)"conserv"
: same ascell
but conservative (high res. to low res.)
tool, optional: Regridder. One of:
"auto"
: tool guessed depending on the method, the first available tool and the grids (rectangular or curvilinear)."vacumm"
: Internal routines."esmf"
and"libcf"
: Regridders provided by UVCDAT."regrid2"
: Old regridder provided by CDAT (rectangular only).
rgdr, optional: An already set up regridder instance to speed up regridding:
CDATRegridder
instance forregrid2
,esmf
andlibcf
tools, else aCurvedInterpolator
instance forvacumm
tool with interpolation on curvilinear grids.getrgdr, optional: Also return the regridder instance if it applies, or None.
Other keywords are passed to special interpolation functions depending on method and choices :
cargen()
when “nat” or “carg” method is used
Tools/methods: Overview table of method availability for each tool.
RECT
means that it only works with rectangular grids.Met/Tool Vacumm regrid2d ESMF Libcf nearest OK bilinear OK OK OK dstwgt OK patch OK cellave RECT OK conserv RECT OK Examples: >>> regrid2d(var, (lon, lat), method='linear') >>> regrid2d(var, grid, method='cellave')
-
regrid2d_method_name
(method, raiseerr=True)[source]¶ Check the method name and return its generic name
-
regrid2dnew
(vari, ggo, method='auto', tool=None, rgdr=None, getrgdr=False, erode=False, **kwargs)¶ Regrid a variable from a regular grid to another
If the input or output grid is curvilinear and
method
is set to"linear"
,"cellave"
or"conserv"
, theCDATRegridder
is used.Params: vari: Variable cdms on regular grid
ggo: Tuple of (x,y) or a cdms grid or a cdms variable with a grid
method, optional: One of:
"auto"
: method guessed betweenlinear
andcellave
according to resolution of input and output grid (seeregrid_method()
)"nearest"
: nearest neighbour"linear"
or"bilinear"
: bilinear interpolation (low res. to high res.)"dstwgt"
: distance weighting between the four nearest grid points (low res. to high res.)"patch"
: patch recovery interpolation (low res. to high res.)"cellave"
: weighted regridding based on areas of cells (high res. to low res.)"conserv"
: same ascell
but conservative (high res. to low res.)
tool, optional: Regridder. One of:
"auto"
: tool guessed depending on the method, the first available tool and the grids (rectangular or curvilinear)."vacumm"
: Internal routines."esmf"
and"libcf"
: Regridders provided by UVCDAT."regrid2"
: Old regridder provided by CDAT (rectangular only).
rgdr, optional: An already set up regridder instance to speed up regridding:
CDATRegridder
instance forregrid2
,esmf
andlibcf
tools, else aCurvedInterpolator
instance forvacumm
tool with interpolation on curvilinear grids.getrgdr, optional: Also return the regridder instance if it applies, or None.
Other keywords are passed to special interpolation functions depending on method and choices :
cargen()
when “nat” or “carg” method is used
Tools/methods: Overview table of method availability for each tool.
RECT
means that it only works with rectangular grids.Met/Tool Vacumm regrid2d ESMF Libcf nearest OK bilinear OK OK OK dstwgt OK patch OK cellave RECT OK conserv RECT OK Examples: >>> regrid2d(var, (lon, lat), method='linear') >>> regrid2d(var, grid, method='cellave')
-
regrid2dold
(vari, ggo, method='auto', mask_thres=0.5, ext=False, bilinear_masking='dstwgt', ext_masking='poly', cdr=None, getcdr=False, usecdr=None, useoldcdr=True, check_mask=True, clipminmax=False, geo=None, **kwargs)[source]¶ Regrid a variable from a regular grid to another
If the input or output grid is curvilinear and
method
is set to"linear"
,"cellave"
or"conserv"
, theCDATRegridder
is used.Params: vari: Variable cdms on regular grid
ggo: Tuple of (x,y) or a cdms grid or a cdms variable with a grid
method, optional: One of:
"auto"
: method guessed according to resolution of input and output grid (seeregrid_method()
)"nearest"
: nearest neighbour"linear"
or"bilinear"
: bilinear interpolation (low res. to high res.)"dstwgt"
: distance weighting (low res. to high res.)"cellave"
: weighted regridding based on areas of cells (high res. to low res.)"conserv"
: same ascell
but conservative (high res. to low res.)"bining"
: simple averaging using bining (very high res. to low res.)"nat"
: Natgrid interpolation (low res. to high res.) (seeGridData
for more info)"carg"
: Interpolation with minicargen(low res. to high res.) (seecargen()
for more info)
cdr, optional:
CDATRegridder
instance.getcdr, optional: Also return the computed
CDATRegridder
instance.usecdr, optional: Force the use or not of a
CDATRegridder
instance, even for rectangular grids.useoldcdr, optional: Force the use the old CDAT regridder (rectangular grids only).
ext, optional: Perform extrapolation when possible
bilinear_masking: the way to handle interpolation near masked values
"nearest"
: brut masking using nearest neighbor"dstwgt"
: distance weight data are used where interpolated mask is lowermask_thres
mask_thres, optional: Threshold for masking points for some methods (~ land fraction) for rectangular grids only:
method="bilinear"
andbilinear_masking="dstwght"
method="cellave"
ormethod="bining"
- ext_masking, optional: Manual masking method when
ext=False
(when needed) with methods [“carg”,] (see
grid_envelop_mask()
) if input grid is not rectangular"poly"
: use the polygon defined by the input grid envelopp and check if output points are inside"nearest"
: use hack with nearest 2d interpolation
- ext_masking, optional: Manual masking method when
Other keywords are passed to special interpolation functions depending on method and choices :
cargen()
when “nat” or “carg” method is used- mask_thres, optional: Time steps when interpolated mask is greater than this value are masked.
Examples: >>> regrid2d(var, (lon, lat), method='bilinear', bilinear_masking='nearest') >>> regrid2d(var, grid, method='cellave', mask_thres=.8) >>> regrid2d(var, grid, method='nat', hor=.2)
-
regrid_method
(gridi, grido, ratio=0.55, iaxi=-1, iaxo=-1)[source]¶ Guess the best regridding method for passing from
gridi
togrido
If
resolution(gridi) <= ratio*resolution(grido)
,method="cellave"
elsemethod="interp"
Theresol()
function is used to estimate the resolution. For grids, a resolution common to X and Y axes is estimated using the following sequence:xres, yres = resol(grid) res = (xres**2+yres**2)**.5
Params: - gridi: Grid, tuple of axes or single axis or array.
- grido: Grid, tuple of axes or single axis or array.
- ratio, optional: Limit ratio of output grid resolution to input grid resolution.
- iaxi/iaxo, optional: Dimension on which to compute the resolution with
resol()
whengridi
andgrido
are single axes with several dimensions.
..note:
The resolution of the grids is checked in their attributes "_xres" and "_yres" before before trying to compute them.
Returns: 'cellave'
OR'linear'
-
regular
(vi, dx=None, verbose=True, auto_bounds=False)[source]¶ Fill a variable with missing values when step of first axis is increasing
- vi: Input array on almost regular axis
- dx: Force grid step to this. Else, auto evaluated.
-
regular_fill1d
(var, k=1, dx=None)[source]¶ Combination: fill1d(regular_fill)) (with their parameter)
-
shift1d
(var, shift=0, bmode=None, axis=-1, copy=True, shiftaxis=True, **kwargs)[source]¶ Interpolate data on an axis shifted by an half cell size
Params: var: array.
shift, optional: Shift to operate.
0
: No shift.<0
: Shilt toward bottom of the axis.>0
: Shilt toward top of the axis.
bmode, optional: Boundary mode.
None
:"extrap"
if axis else"same"
"linear"
: Linear extrapolation."same"
: Constant extrapolation."masked"
: Mask outside data."cyclic"
: Cyclic extrapolation.
axis, optional: Axis on which to operate.
copy, optional: Always input copy data.
-
shift2d
(vari, ishift=0, jshift=0, bmode=None, copy=True, **kwargs)[source]¶ Interpolate data on an grid shifted by an half cell size in X and Y
X and Y are supposed to be the -1 and -2 axes of var.
Params: var: array.
[i/j]shift, optional: Shift to operate along X/Y.
0
: No shift.<0
: Shilt toward bottom of the axis.>0
: Shilt toward top of the axis.
bmode, optional: Boundary mode.
"linear"
: Linear extrapolation."same"
: Constant extrapolation."masked"
: Mask outside data."cyclic"
: Cyclic extrapolation.
copy, optional: Always input copy data.
-
shiftgrid
(gg, ishift=0, jshift=0, bmode='linear', **kwargs)[source]¶ Shift a grid of an half cell
Params: - gg: cdms2 grid.
- i/jshift: Fraction cell to shift.
-
spline_interp1d
(old_var, new_axis, check_missing=True, k=3, **kwargs)[source]¶ Backward compatibility function
See
regrid1d()
-
transect
(var, lons, lats, depths=None, times=None, method='linear', subsamp=3, getcoords=False, outaxis=None, depth=None, **kwargs)[source]¶ Make a transect in a -[T][Z]YX variable
It calls
transect_specs()
to compute transect coordinates when not explictly specified, andgrid2xy()
to perform 4D interpolations.Example: >>> tsst = transect(sst, (1,1.6), (42.,43.), subsamp=4) >>> tmld = transect(mld, [1.,1.,2.], [42.,43.,43.], outaxis='dist') >>> tprof = transect(temp, 1., 42.)
Params: var: MV2 array of at least rank 2 (YX).
lons/lats: Specification of transect, either
- Coordinates of first and last point in degrees as
tuples in the form
(lon0,lon1)
and(lat0,lat1)
. The array of coordinates is generated usingtransect_specs()
. - Or explicit array of coordinates (as scalars, lists or arrays).
- Coordinates of first and last point in degrees as
tuples in the form
depths, optional: Output depths. If not a tuple, it must have the same size as lons and lats.
times, optional: Tuple, or time sequence or axis of the same length as resulting coordinates. If provided, the interpolation is first done in space, them onto this lagrangian time, and the final space-time trajectory is returned. If
outaxis
is None,taxis
becomes the output axis.subsamp, optional: Subsampling with respect to grid cell (only when coordinates are not explicit).
method, optional: Interpolation method (see
grid2xy()
).getcoords, optional: Also get computed coordiantes.
outaxis, optional: Output spatial axis (see
grid2xy()
).
Return: tvar
ortvar,tons,tlats
-
xy2grid
(*args, **kwargs)[source]¶ Alias for
griddata()