Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1

2

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

84

85

86

87

88

89

90

91

92

93

94

95

96

97

98

99

100

101

102

103

104

105

106

107

108

109

110

111

112

113

114

115

116

117

118

119

120

121

122

123

124

125

126

127

128

129

130

131

132

133

134

135

136

137

138

139

140

141

142

143

144

145

146

147

148

149

150

151

152

153

154

155

156

157

158

159

160

161

162

163

164

165

166

167

168

169

170

171

172

173

174

175

176

177

178

179

180

181

182

183

184

185

186

187

188

189

190

191

192

193

194

195

196

197

198

199

200

201

202

203

204

205

206

207

208

209

210

211

212

213

214

215

216

217

218

219

220

221

222

223

224

225

226

227

228

229

230

231

232

233

234

235

236

237

238

239

240

241

242

243

244

245

246

247

248

249

250

251

252

253

254

255

256

257

258

259

260

261

"""rbf - Radial basis functions for interpolation/smoothing scattered Nd data. 

 

Written by John Travers <jtravs@gmail.com>, February 2007 

Based closely on Matlab code by Alex Chirokov 

Additional, large, improvements by Robert Hetland 

Some additional alterations by Travis Oliphant 

 

Permission to use, modify, and distribute this software is given under the 

terms of the SciPy (BSD style) license. See LICENSE.txt that came with 

this distribution for specifics. 

 

NO WARRANTY IS EXPRESSED OR IMPLIED. USE AT YOUR OWN RISK. 

 

Copyright (c) 2006-2007, Robert Hetland <hetland@tamu.edu> 

Copyright (c) 2007, John Travers <jtravs@gmail.com> 

 

Redistribution and use in source and binary forms, with or without 

modification, are permitted provided that the following conditions are 

met: 

 

* Redistributions of source code must retain the above copyright 

notice, this list of conditions and the following disclaimer. 

 

* Redistributions in binary form must reproduce the above 

copyright notice, this list of conditions and the following 

disclaimer in the documentation and/or other materials provided 

with the distribution. 

 

* Neither the name of Robert Hetland nor the names of any 

contributors may be used to endorse or promote products derived 

from this software without specific prior written permission. 

 

THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS 

"AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT 

LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR 

A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT 

OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, 

SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT 

LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, 

DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY 

THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT 

(INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE 

OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 

""" 

from __future__ import division, print_function, absolute_import 

 

import sys 

import numpy as np 

 

from scipy import linalg 

from scipy._lib.six import callable, get_method_function, get_function_code 

from scipy.special import xlogy 

from scipy.spatial.distance import cdist, pdist, squareform 

 

__all__ = ['Rbf'] 

 

 

class Rbf(object): 

""" 

Rbf(*args) 

 

A class for radial basis function approximation/interpolation of 

n-dimensional scattered data. 

 

Parameters 

---------- 

*args : arrays 

x, y, z, ..., d, where x, y, z, ... are the coordinates of the nodes 

and d is the array of values at the nodes 

function : str or callable, optional 

The radial basis function, based on the radius, r, given by the norm 

(default is Euclidean distance); the default is 'multiquadric':: 

 

'multiquadric': sqrt((r/self.epsilon)**2 + 1) 

'inverse': 1.0/sqrt((r/self.epsilon)**2 + 1) 

'gaussian': exp(-(r/self.epsilon)**2) 

'linear': r 

'cubic': r**3 

'quintic': r**5 

'thin_plate': r**2 * log(r) 

 

If callable, then it must take 2 arguments (self, r). The epsilon 

parameter will be available as self.epsilon. Other keyword 

arguments passed in will be available as well. 

 

epsilon : float, optional 

Adjustable constant for gaussian or multiquadrics functions 

- defaults to approximate average distance between nodes (which is 

a good start). 

smooth : float, optional 

Values greater than zero increase the smoothness of the 

approximation. 0 is for interpolation (default), the function will 

always go through the nodal points in this case. 

norm : str, callable, optional 

A function that returns the 'distance' between two points, with 

inputs as arrays of positions (x, y, z, ...), and an output as an 

array of distance. E.g., the default: 'euclidean', such that the result 

is a matrix of the distances from each point in ``x1`` to each point in 

``x2``. For more options, see documentation of 

`scipy.spatial.distances.cdist`. 

 

Attributes 

---------- 

N : int 

The number of data points (as determined by the input arrays). 

di : ndarray 

The 1-D array of data values at each of the data coordinates `xi`. 

xi : ndarray 

The 2-D array of data coordinates. 

function : str or callable 

The radial basis function. See description under Parameters. 

epsilon : float 

Parameter used by gaussian or multiquadrics functions. See Parameters. 

smooth : float 

Smoothing parameter. See description under Parameters. 

norm : str or callable 

The distance function. See description under Parameters. 

nodes : ndarray 

A 1-D array of node values for the interpolation. 

A : internal property, do not use 

 

Examples 

-------- 

>>> from scipy.interpolate import Rbf 

>>> x, y, z, d = np.random.rand(4, 50) 

>>> rbfi = Rbf(x, y, z, d) # radial basis function interpolator instance 

>>> xi = yi = zi = np.linspace(0, 1, 20) 

>>> di = rbfi(xi, yi, zi) # interpolated values 

>>> di.shape 

(20,) 

 

""" 

# Available radial basis functions that can be selected as strings; 

# they all start with _h_ (self._init_function relies on that) 

def _h_multiquadric(self, r): 

return np.sqrt((1.0/self.epsilon*r)**2 + 1) 

 

def _h_inverse_multiquadric(self, r): 

return 1.0/np.sqrt((1.0/self.epsilon*r)**2 + 1) 

 

def _h_gaussian(self, r): 

return np.exp(-(1.0/self.epsilon*r)**2) 

 

def _h_linear(self, r): 

return r 

 

def _h_cubic(self, r): 

return r**3 

 

def _h_quintic(self, r): 

return r**5 

 

def _h_thin_plate(self, r): 

return xlogy(r**2, r) 

 

# Setup self._function and do smoke test on initial r 

def _init_function(self, r): 

if isinstance(self.function, str): 

self.function = self.function.lower() 

_mapped = {'inverse': 'inverse_multiquadric', 

'inverse multiquadric': 'inverse_multiquadric', 

'thin-plate': 'thin_plate'} 

if self.function in _mapped: 

self.function = _mapped[self.function] 

 

func_name = "_h_" + self.function 

if hasattr(self, func_name): 

self._function = getattr(self, func_name) 

else: 

functionlist = [x[3:] for x in dir(self) 

if x.startswith('_h_')] 

raise ValueError("function must be a callable or one of " + 

", ".join(functionlist)) 

self._function = getattr(self, "_h_"+self.function) 

elif callable(self.function): 

allow_one = False 

if hasattr(self.function, 'func_code') or \ 

hasattr(self.function, '__code__'): 

val = self.function 

allow_one = True 

elif hasattr(self.function, "im_func"): 

val = get_method_function(self.function) 

elif hasattr(self.function, "__call__"): 

val = get_method_function(self.function.__call__) 

else: 

raise ValueError("Cannot determine number of arguments to " 

"function") 

 

argcount = get_function_code(val).co_argcount 

if allow_one and argcount == 1: 

self._function = self.function 

elif argcount == 2: 

if sys.version_info[0] >= 3: 

self._function = self.function.__get__(self, Rbf) 

else: 

import new 

self._function = new.instancemethod(self.function, self, 

Rbf) 

else: 

raise ValueError("Function argument must take 1 or 2 " 

"arguments.") 

 

a0 = self._function(r) 

if a0.shape != r.shape: 

raise ValueError("Callable must take array and return array of " 

"the same shape") 

return a0 

 

def __init__(self, *args, **kwargs): 

# `args` can be a variable number of arrays; we flatten them and store 

# them as a single 2-D array `xi` of shape (n_args-1, array_size), 

# plus a 1-D array `di` for the values. 

# All arrays must have the same number of elements 

self.xi = np.asarray([np.asarray(a, dtype=np.float_).flatten() 

for a in args[:-1]]) 

self.N = self.xi.shape[-1] 

self.di = np.asarray(args[-1]).flatten() 

 

if not all([x.size == self.di.size for x in self.xi]): 

raise ValueError("All arrays must be equal length.") 

 

self.norm = kwargs.pop('norm', 'euclidean') 

self.epsilon = kwargs.pop('epsilon', None) 

if self.epsilon is None: 

# default epsilon is the "the average distance between nodes" based 

# on a bounding hypercube 

ximax = np.amax(self.xi, axis=1) 

ximin = np.amin(self.xi, axis=1) 

edges = ximax - ximin 

edges = edges[np.nonzero(edges)] 

self.epsilon = np.power(np.prod(edges)/self.N, 1.0/edges.size) 

 

self.smooth = kwargs.pop('smooth', 0.0) 

self.function = kwargs.pop('function', 'multiquadric') 

 

# attach anything left in kwargs to self for use by any user-callable 

# function or to save on the object returned. 

for item, value in kwargs.items(): 

setattr(self, item, value) 

 

self.nodes = linalg.solve(self.A, self.di) 

 

@property 

def A(self): 

# this only exists for backwards compatibility: self.A was available 

# and, at least technically, public. 

r = squareform(pdist(self.xi.T, self.norm)) # Pairwise norm 

return self._init_function(r) - np.eye(self.N)*self.smooth 

 

def _call_norm(self, x1, x2): 

return cdist(x1.T, x2.T, self.norm) 

 

def __call__(self, *args): 

args = [np.asarray(x) for x in args] 

if not all([x.shape == y.shape for x in args for y in args]): 

raise ValueError("Array lengths must be equal") 

 

shp = args[0].shape 

xa = np.asarray([a.flatten() for a in args], dtype=np.float_) 

r = self._call_norm(xa, self.xi) 

return np.dot(self._function(r), self.nodes).reshape(shp)