Radial Basis Function

Module focused on the implementation of the Radial Basis Functions interpolation technique. This technique is still based on the use of a set of parameters, the so-called control points, as for FFD, but RBF is interpolatory. Another important key point of RBF strategy relies in the way we can locate the control points: in fact, instead of FFD where control points need to be placed inside a regular lattice, with RBF we have no more limitations. So we have the possibility to perform localized control points refinements. The module is analogous to the freeform one.

Theoretical Insight

As reference please consult M.D. Buhmann, Radial Basis Functions, volume 12 of Cambridge monographs on applied and computational mathematics. Cambridge University Press, UK, 2003. This implementation follows D. Forti and G. Rozza, Efficient geometrical parametrization techniques of interfaces for reduced order modelling: application to fluid-structure interaction coupling problems, International Journal of Computational Fluid Dynamics.

RBF shape parametrization technique is based on the definition of a map, \mathcal{M}(\boldsymbol{x}) : \mathbb{R}^n \rightarrow \mathbb{R}^n, that allows the possibility of transferring data across non-matching grids and facing the dynamic mesh handling. The map introduced is defines as follows

\mathcal{M}(\boldsymbol{x}) = p(\boldsymbol{x}) + \sum_{i=1}^{\mathcal{N}_C} \gamma_i \varphi(\| \boldsymbol{x} - \boldsymbol{x_{C_i}} \|)

where p(\boldsymbol{x}) is a low_degree polynomial term, \gamma_i is the weight, corresponding to the a-priori selected \mathcal{N}_C control points, associated to the i-th basis function, and \varphi(\| \boldsymbol{x} - \boldsymbol{x_{C_i}} \|) a radial function based on the Euclidean distance between the control points position \boldsymbol{x_{C_i}} and \boldsymbol{x}. A radial basis function, generally, is a real-valued function whose value depends only on the distance from the origin, so that \varphi(\boldsymbol{x}) = \tilde{\varphi}(\| \boldsymbol{x} \|).

The matrix version of the formula above is:

\mathcal{M}(\boldsymbol{x}) = \boldsymbol{c} + \boldsymbol{Q}\boldsymbol{x} + \boldsymbol{W^T}\boldsymbol{d}(\boldsymbol{x})

The idea is that after the computation of the weights and the polynomial terms from the coordinates of the control points before and after the deformation, we can deform all the points of the mesh accordingly. Among the most common used radial basis functions for modelling 2D and 3D shapes, we consider Gaussian splines, Multi-quadratic biharmonic splines, Inverted multi-quadratic biharmonic splines, Thin-plate splines, Beckert and Wendland C^2 basis and Polyharmonic splines all defined and implemented below.

class RBF(original_control_points=None, deformed_control_points=None, func='gaussian_spline', radius=0.5, extra_parameter=None)[source]

Bases: pygem.deformation.Deformation

Class that handles the Radial Basis Functions interpolation on the mesh points.

Parameters
  • 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.

  • func – the basis function to use in the transformation. Several basis function are already implemented and they are available through the RBF by passing the name of the right function (see class documentation for the updated list of basis function). A callable object can be passed as basis function.

  • radius (float) – the scaling parameter r that affects the shape of the basis functions. For details see the class RBF. The default value is 0.5.

  • extra_parameter (dict) – the additional parameters that may be passed to the kernel function. Default is None.

Variables
  • weights (numpy.ndarray) – the matrix formed by the weights corresponding to the a-priori selected N control points, associated to the basis functions and c and Q terms that describe the polynomial of order one p(x) = c + Qx. The shape is (n_control_points+1+3, 3). It is computed internally.

  • 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.

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

  • basis (callable) – the basis functions to use in the transformation.

  • radius (float) – the scaling parameter that affects the shape of the basis functions.

  • extra (dict) – the additional parameters that may be passed to the kernel function.

Example
>>> from pygem import RBF
>>> import numpy as np
>>> rbf = RBF(func='gaussian_spline')
>>> xv = np.linspace(0, 1, 20)
>>> yv = np.linspace(0, 1, 20)
>>> zv = np.linspace(0, 1, 20)
>>> z, y, x = np.meshgrid(zv, yv, xv)
>>> mesh = np.array([x.ravel(), y.ravel(), z.ravel()])
>>> deformed_mesh = rbf(mesh)
property n_control_points

Total number of control points.

Return type

int

property basis

The kernel to use in the deformation.

Getter

Returns the callable kernel

Setter

Sets the kernel. It is possible to pass the name of the function (check the list of all implemented functions in the pygem.rbf_factory.RBFFactory class) or directly the callable function.

Type

callable

_get_weights(X, Y)[source]

This private method, given the original control points and the deformed ones, returns the matrix with the weights and the polynomial terms, that is W, c^T and Q^T. The shape is (n_control_points+1+3, 3).

Parameters
  • X (numpy.ndarray) – it is an n_control_points-by-3 array with the coordinates of the original interpolation control points before the deformation.

  • Y (numpy.ndarray) – it is an n_control_points-by-3 array with the coordinates of the interpolation control points after the deformation.

Returns

weights: the 2D array with the weights and the polynomial terms.

Return type

numpy.ndarray

read_parameters(filename='parameters_rbf.prm')[source]

Reads in the parameters file and fill the self structure.

Parameters

filename (string) – parameters file to be read in. Default value is parameters_rbf.prm.

write_parameters(filename='parameters_rbf.prm')[source]

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

Parameters

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

plot_points(filename=None)[source]

Method to plot the control points. It is possible to save the resulting figure.

Parameters

filename (str) – if None the figure is shown, otherwise it is saved on the specified filename. Default is None.

compute_weights()[source]

This method compute the weights according to the original_control_points and deformed_control_points arrays.

__call__(src_pts)[source]

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