mfsetup.interpolate module

class mfsetup.interpolate.Interpolator(xyz, uvw, d=2, source_values_mask=None)[source]

Bases: object

Speed up barycentric interpolation similar to scipy.interpolate.griddata (method=’linear’), by computing the weights once and then re-using them for successive interpolation with the same source and destination points.

Parameters:
xyzndarray or tuple

x, y, z, … locations of source data. (shape n source points x ndims)

uvwndarray or tuple

x, y, z, … locations of where source data will be interpolated (shape n destination points x ndims)

dint

Number of dimensions (2 for 2D, 3 for 3D, etc.)

source_values_maskboolean array

Boolean array of same structure as the source_values array input to the interpolate() method, with the same number of active values as the size of xyz.

Notes

The methods employed are based on this Stack Overflow post: https://stackoverflow.com/questions/20915502/speedup-scipy-griddata-for-multiple-interpolations-between-two-irregular-grids

property interp_weights

Calculate the interpolation weights.

interpolate(source_values, method='linear')[source]

Interpolate values in source_values to the destination points in the uvw attribute. using modelgrid instances attached to the source and destination models.

Parameters:
source_valuesndarray

Values to be interpolated to destination points. Array must be the same size as the number of source points, or the number of active points within source points, as defined by the source_values_mask array input to the Interpolator.

methodstr (‘linear’, ‘nearest’)

Interpolation method. With ‘linear’ a triangular mesh is discretized around the source points, and barycentric weights representing the influence of the d +1 source points on each destination point (where d is the number of dimensions), are computed. With ‘nearest’, the input is simply passed to scipy.interpolate.griddata().

Returns:
interpolated1D numpy array

Array of interpolated values at the destination locations.

property source_values_mask
mfsetup.interpolate.get_source_dest_model_xys(source_model, dest_model, source_mask=None)[source]

Get the xyz and uvw inputs to the interp_weights function.

Parameters:
source_modelflopy.modeflow.Modflow, flopy.mf6.MFModel, or MFsetupGrid instance
dest_modelmfsetup.MFnwtModel, mfsetup.MF6model instance
mfsetup.interpolate.interp_weights(xyz, uvw, d=2, mask=None)[source]

Speed up interpolation vs scipy.interpolate.griddata (method=’linear’), by only computing the weights once: https://stackoverflow.com/questions/20915502/speedup-scipy-griddata-for-multiple-interpolations-between-two-irregular-grids

Parameters:
xyzndarray or tuple

x, y, z, … locations of source data. (shape n source points x ndims)

uvwndarray or tuple

x, y, z, … locations of where source data will be interpolated (shape n destination points x ndims)

dint

Number of dimensions (2 for 2D, 3 for 3D, etc.)

Returns:
indicesndarray of shape n destination points x 3

Index positions in flattened (1D) xyz array

weightsndarray of shape n destination points x 3

Fractional weights for each row position in indices. Weights in each row sum to 1 across the 3 columns.

mfsetup.interpolate.interpolate(values, vtx, wts, fill_value='mean')[source]

Apply the interpolation weights to a set of values.

Parameters:
values1D array of length n source points (same as xyz in interp_weights)
vtxindices returned by interp_weights
wtsweights returned by interp_weights
fill_valuefloat

Value used to fill in for requested points outside of the convex hull of the input points (i.e., those with at least one negative weight). If not provided, then the default is nan.

Returns
——-
interpolated values
mfsetup.interpolate.regrid(arr, grid, grid2, mask1=None, mask2=None, method='linear')[source]

Interpolate array values from one model grid to another, using scipy.interpolate.griddata.

Parameters:
arr2D numpy array

Source data

gridflopy.discretization.StructuredGrid instance

Source grid

grid2flopy.discretization.StructuredGrid instance

Destination grid (to interpolate onto)

mask1boolean array

mask for source grid. Areas that are masked will be converted to nans, and not included in the interpolation.

mask2boolean array

mask denoting active area for destination grid. The mean value will be applied to inactive areas if linear interpolation is used (not for integer/categorical arrays).

methodstr

interpolation method (‘nearest’, ‘linear’, or ‘cubic’)

mfsetup.interpolate.regrid3d(arr, grid, grid2, mask1=None, mask2=None, method='linear')[source]

Interpolate array values from one model grid to another, using scipy.interpolate.griddata.

Parameters:
arr3D numpy array

Source data

gridflopy.discretization.StructuredGrid instance

Source grid

grid2flopy.discretization.StructuredGrid instance

Destination grid (to interpolate onto)

mask1boolean array

mask for source grid. Areas that are masked will be converted to nans, and not included in the interpolation.

mask2boolean array

mask denoting active area for destination grid. The mean value will be applied to inactive areas if linear interpolation is used (not for integer/categorical arrays).

methodstr

interpolation method (‘nearest’, ‘linear’, or ‘cubic’)

Returns:
arr3D numpy array

Interpolated values at the x, y, z locations in grid2.