Inverse Distance Weighting

Module focused on the Inverse Distance Weighting interpolation technique. The IDW algorithm is an average moving interpolation that is usually applied to highly variable data. The main idea of this interpolation strategy lies in fact that it is not desirable to honour local high/low values but rather to look at a moving average of nearby data points and estimate the local trends. The node value is calculated by averaging the weighted sum of all the points. Data points that lie progressively farther from the node inuence much less the computed value than those lying closer to the node.

Theoretical Insight

This implementation is based on the simplest form of inverse distance weighting interpolation, proposed by D. Shepard, A two-dimensional interpolation function for irregularly-spaced data, Proceedings of the 23 rd ACM National Conference.

The interpolation value u of a given point \mathrm{x} from a set of samples u_k = u(\mathrm{x}_k), with k = 1,2,\dotsc,\mathcal{N}, is given by:

u(\mathrm{x}) = \displaystyle\sum_{k=1}^\mathcal{N} \frac{w(\mathrm{x},\mathrm{x}_k)} {\displaystyle\sum_{j=1}^\mathcal{N} w(\mathrm{x},\mathrm{x}_j)} u_k

where, in general, w(\mathrm{x}, \mathrm{x}_i) represents the weighting function:

w(\mathrm{x}, \mathrm{x}_i) = \| \mathrm{x} - \mathrm{x}_i \|^{-p}

being \| \mathrm{x} - \mathrm{x}_i \|^{-p} \ge 0 is the Euclidean distance between \mathrm{x} and data point \mathrm{x}_i and p is a power parameter, typically equal to 2.

class IDW(original_control_points=None, deformed_control_points=None, power=2)[source]

Bases: pygem.deformation.Deformation

Class that perform the Inverse Distance Weighting (IDW).

Parameters
  • power (int) – the power parameter. The default value is 2.

  • original_control_points (numpy.ndarray) – it is an (n_control_points, 3) array with the coordinates of the original interpolation control points before the deformation. The default is the vertices of the unit cube.

  • deformed_control_points (numpy.ndarray) – it is an (n_control_points, 3) array with the coordinates of the interpolation control points after the deformation. The default is the vertices of the unit cube.

Variables
  • power (int) – the power parameter. The default value is 2.

  • original_control_points (numpy.ndarray) – it is an (n_control_points, 3) array with the coordinates of the original interpolation control points before the deformation. The default is the vertices of the unit cube.

  • deformed_control_points (numpy.ndarray) – it is an (n_control_points, 3) array with the coordinates of the interpolation control points after the deformation. The default is the vertices of the unit cube.

Example

>>> from pygem import IDW
>>> import numpy as np
>>> nx, ny, nz = (20, 20, 20)
>>> mesh = np.zeros((nx * ny * nz, 3))
>>> xv = np.linspace(0, 1, nx)
>>> yv = np.linspace(0, 1, ny)
>>> zv = np.linspace(0, 1, nz)
>>> z, y, x = np.meshgrid(zv, yv, xv)
>>> mesh_points = np.array([x.ravel(), y.ravel(), z.ravel()])
>>> idw = IDW()
>>> idw.read_parameters('tests/test_datasets/parameters_idw_cube.prm')
>>> new_mesh_points = idw(mesh_points.T)
__call__(src_pts)[source]

This method performs the deformation of the mesh points. After the execution it sets self.modified_mesh_points.

read_parameters(filename)[source]

Reads in the parameters file and fill the self structure.

Parameters

filename (string) – parameters file to be read in.

write_parameters(filename)[source]

This method writes a parameters file (.prm) called filename and fills it with all the parameters class members.

Parameters

filename (string) – parameters file to be written out.