File: interpn.m

package info (click to toggle)
octave 3.8.2-4
  • links: PTS, VCS
  • area: main
  • in suites: jessie, jessie-kfreebsd
  • size: 84,396 kB
  • ctags: 45,547
  • sloc: cpp: 293,356; ansic: 42,041; fortran: 23,669; sh: 13,629; objc: 7,890; yacc: 7,093; lex: 3,442; java: 2,125; makefile: 1,589; perl: 1,009; awk: 974; xml: 34
file content (320 lines) | stat: -rw-r--r-- 9,881 bytes parent folder | download | duplicates (3)
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
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
## Copyright (C) 2007-2013 David Bateman
##
## This file is part of Octave.
##
## Octave is free software; you can redistribute it and/or modify it
## under the terms of the GNU General Public License as published by
## the Free Software Foundation; either version 3 of the License, or (at
## your option) any later version.
##
## Octave is distributed in the hope that it will be useful, but
## WITHOUT ANY WARRANTY; without even the implied warranty of
## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
## General Public License for more details.
##
## You should have received a copy of the GNU General Public License
## along with Octave; see the file COPYING.  If not, see
## <http://www.gnu.org/licenses/>.

## -*- texinfo -*-
## @deftypefn  {Function File} {@var{vi} =} interpn (@var{x1}, @var{x2}, @dots{}, @var{v}, @var{y1}, @var{y2}, @dots{})
## @deftypefnx {Function File} {@var{vi} =} interpn (@var{v}, @var{y1}, @var{y2}, @dots{})
## @deftypefnx {Function File} {@var{vi} =} interpn (@var{v}, @var{m})
## @deftypefnx {Function File} {@var{vi} =} interpn (@var{v})
## @deftypefnx {Function File} {@var{vi} =} interpn (@dots{}, @var{method})
## @deftypefnx {Function File} {@var{vi} =} interpn (@dots{}, @var{method}, @var{extrapval})
##
## Perform @var{n}-dimensional interpolation, where @var{n} is at least two.
## Each element of the @var{n}-dimensional array @var{v} represents a value
## at a location given by the parameters @var{x1}, @var{x2}, @dots{}, @var{xn}.
## The parameters @var{x1}, @var{x2}, @dots{}, @var{xn} are either
## @var{n}-dimensional arrays of the same size as the array @var{v} in
## the @qcode{"ndgrid"} format or vectors.  The parameters @var{y1}, etc.
## respect a similar format to @var{x1}, etc., and they represent the points
## at which the array @var{vi} is interpolated.
##
## If @var{x1}, @dots{}, @var{xn} are omitted, they are assumed to be
## @code{x1 = 1 : size (@var{v}, 1)}, etc.  If @var{m} is specified, then
## the interpolation adds a point half way between each of the interpolation
## points.  This process is performed @var{m} times.  If only @var{v} is
## specified, then @var{m} is assumed to be @code{1}.
##
## Method is one of:
##
## @table @asis
## @item @qcode{"nearest"}
## Return the nearest neighbor.
##
## @item @qcode{"linear"}
## Linear interpolation from nearest neighbors.
##
## @item @qcode{"cubic"}
## Cubic interpolation from four nearest neighbors (not implemented yet).
##
## @item @qcode{"spline"}
## Cubic spline interpolation---smooth first and second derivatives
## throughout the curve.
## @end table
##
## The default method is @qcode{"linear"}.
##
## If @var{extrapval} is the scalar value, use it to replace the values
## beyond the endpoints with that number.  If @var{extrapval} is missing,
## assume NA.
## @seealso{interp1, interp2, spline, ndgrid}
## @end deftypefn

function vi = interpn (varargin)

  method = "linear";
  extrapval = NA;
  nargs = nargin;

  if (nargin < 1 || ! isnumeric (varargin{1}))
    print_usage ();
  endif

  if (ischar (varargin{end}))
    method = varargin{end};
    nargs -= 1;
  elseif (nargs > 1 && ischar (varargin{end - 1}))
    if (! isnumeric (varargin{end}) || ! isscalar (varargin{end}))
      error ("interpn: extrapal is expected to be a numeric scalar");
    endif
    method = varargin{end - 1};
    extrapval = varargin{end};
    nargs -= 2;
  endif

  if (nargs < 3)
    v = varargin{1};
    m = 1;
    if (nargs == 2)
      if (ischar (varargin{2}))
        method = varargin{2};
      elseif (isnumeric (m) && isscalar (m) && fix (m) == m)
        m = varargin{2};
      else
        print_usage ();
      endif
    endif
    sz = size (v);
    nd = ndims (v);
    x = cell (1, nd);
    y = cell (1, nd);
    for i = 1 : nd
      x{i} = 1 : sz(i);
      y{i} = 1 : (1 / (2 ^ m)) : sz(i);
    endfor
    y{1} = y{1}.';
    [y{:}] = ndgrid (y{:});
  elseif (! isvector (varargin{1}) && nargs == (ndims (varargin{1}) + 1))
    v = varargin{1};
    sz = size (v);
    nd = ndims (v);
    x = cell (1, nd);
    y = varargin(2 : nargs);
    for i = 1 : nd
      x{i} = 1 : sz(i);
    endfor
  elseif (rem (nargs, 2) == 1
          && nargs == (2 * ndims (varargin{ceil (nargs / 2)})) + 1)
    nv = ceil (nargs / 2);
    v = varargin{nv};
    sz = size (v);
    nd = ndims (v);
    x = varargin(1 : (nv - 1));
    y = varargin((nv + 1) : nargs);
  else
    error ("interpn: wrong number or incorrectly formatted input arguments");
  endif

  if (any (! cellfun ("isvector", x)))
    for i = 2 : nd
      if (! size_equal (x{1}, x{i}) || ! size_equal (x{i}, v))
        error ("interpn: dimensional mismatch");
      endif
      idx(1 : nd) = {1};
      idx(i) = ":";
      x{i} = x{i}(idx{:})(:);
    endfor
    idx(1 : nd) = {1};
    idx(1) = ":";
    x{1} = x{1}(idx{:})(:);
  endif

  method = tolower (method);

  all_vectors = all (cellfun ("isvector", y));
  different_lengths = numel (unique (cellfun ("numel", y))) > 1;
  if (all_vectors && different_lengths)
    [foobar(1:numel(y)).y] = ndgrid (y{:});
    y = {foobar.y};
  endif

  if (strcmp (method, "linear"))
    vi = __lin_interpn__ (x{:}, v, y{:});
    vi(isna (vi)) = extrapval;
  elseif (strcmp (method, "nearest"))
    yshape = size (y{1});
    yidx = cell (1, nd);
    for i = 1 : nd
      y{i} = y{i}(:);
      yidx{i} = lookup (x{i}, y{i}, "lr");
    endfor
    idx = cell (1,nd);
    for i = 1 : nd
      idx{i} = yidx{i} + (y{i} - x{i}(yidx{i})(:) >= x{i}(yidx{i} + 1)(:) - y{i});
    endfor
    vi = v(sub2ind (sz, idx{:}));
    idx = zeros (prod (yshape), 1);
    for i = 1 : nd
      idx |= y{i} < min (x{i}(:)) | y{i} > max (x{i}(:));
    endfor
    vi(idx) = extrapval;
    vi = reshape (vi, yshape);
  elseif (strcmp (method, "spline"))
    if (any (! cellfun ("isvector", y)))
      for i = 2 : nd
        if (! size_equal (y{1}, y{i}))
          error ("interpn: dimensional mismatch");
        endif
        idx(1 : nd) = {1};
        idx(i) = ":";
        y{i} = y{i}(idx{:});
      endfor
      idx(1 : nd) = {1};
      idx(1) = ":";
      y{1} = y{1}(idx{:});
    endif

    vi = __splinen__ (x, v, y, extrapval, "interpn");

    if (size_equal (y{:}))
      ly = length (y{1});
      idx = cell (1, ly);
      q = cell (1, nd);
      for i = 1 : ly
        q(:) = i;
        idx{i} = q;
      endfor
      vi = vi(cellfun (@(x) sub2ind (size (vi), x{:}), idx));
      vi = reshape (vi, size (y{1}));
    endif
  elseif (strcmp (method, "cubic"))
    error ("interpn: cubic interpolation not yet implemented");
  else
    error ("interpn: unrecognized interpolation METHOD");
  endif

endfunction


%!demo
%! clf;
%! colormap ("default");
%! A = [13,-1,12;5,4,3;1,6,2];
%! x = [0,1,4];  y = [10,11,12];
%! xi = linspace (min (x), max (x), 17);
%! yi = linspace (min (y), max (y), 26)';
%! mesh (xi, yi, interpn (x,y,A.',xi,yi, "linear").');
%! [x,y] = meshgrid (x,y);
%! hold on; plot3 (x(:),y(:),A(:),"b*"); hold off;

%!demo
%! clf;
%! colormap ("default");
%! A = [13,-1,12;5,4,3;1,6,2];
%! x = [0,1,4];  y = [10,11,12];
%! xi = linspace (min (x), max (x), 17);
%! yi = linspace (min (y), max (y), 26)';
%! mesh (xi, yi, interpn (x,y,A.',xi,yi, "nearest").');
%! [x,y] = meshgrid (x,y);
%! hold on; plot3 (x(:),y(:),A(:),"b*"); hold off;

%!#demo  # FIXME: Uncomment when support for "cubic" has been added
%! clf;
%! colormap ("default");
%! A = [13,-1,12;5,4,3;1,6,2];
%! x = [0,1,2];  y = [10,11,12];
%! xi = linspace (min (x), max (x), 17);
%! yi = linspace (min (y), max (y), 26)';
%! mesh (xi, yi, interpn (x,y,A.',xi,yi, "cubic").');
%! [x,y] = meshgrid (x,y);
%! hold on; plot3 (x(:),y(:),A(:),"b*"); hold off;

%!demo
%! clf;
%! colormap ("default");
%! A = [13,-1,12;5,4,3;1,6,2];
%! x = [0,1,2];  y = [10,11,12];
%! xi = linspace (min (x), max (x), 17);
%! yi = linspace (min (y), max (y), 26)';
%! mesh (xi, yi, interpn (x,y,A.',xi,yi, "spline").');
%! [x,y] = meshgrid (x,y);
%! hold on; plot3 (x(:),y(:),A(:),"b*"); hold off;

%!demo
%! clf;
%! colormap ("default");
%! x = y = z = -1:1;
%! f = @(x,y,z) x.^2 - y - z.^2;
%! [xx, yy, zz] = meshgrid (x, y, z);
%! v = f (xx,yy,zz);
%! xi = yi = zi = -1:0.1:1;
%! [xxi, yyi, zzi] = ndgrid (xi, yi, zi);
%! vi = interpn (x, y, z, v, xxi, yyi, zzi, "spline");
%! mesh (yi, zi, squeeze (vi(1,:,:)));

%!test
%! [x,y,z] = ndgrid (0:2);
%! f = x + y + z;
%! assert (interpn (x,y,z,f,[.5 1.5],[.5 1.5],[.5 1.5]), [1.5, 4.5]);
%! assert (interpn (x,y,z,f,[.51 1.51],[.51 1.51],[.51 1.51],"nearest"), [3, 6]);
%! assert (interpn (x,y,z,f,[.5 1.5],[.5 1.5],[.5 1.5],"spline"), [1.5, 4.5]);
%! assert (interpn (x,y,z,f,x,y,z), f);
%! assert (interpn (x,y,z,f,x,y,z,"nearest"), f);
%! assert (interpn (x,y,z,f,x,y,z,"spline"), f);

%!test
%! [x, y, z] = ndgrid (0:2, 1:4, 2:6);
%! f = x + y + z;
%! xi = [0.5 1.0 1.5];  yi = [1.5 2.0 2.5 3.5];  zi = [2.5 3.5 4.0 5.0 5.5];
%! fi = interpn (x, y, z, f, xi, yi, zi);
%! [xi, yi, zi] = ndgrid (xi, yi, zi);
%! assert (fi, xi + yi + zi);

%!test
%! xi = 0:2;  yi = 1:4;  zi = 2:6;
%! [x, y, z] = ndgrid (xi, yi, zi);
%! f = x + y + z;
%! fi = interpn (x, y, z, f, xi, yi, zi, "nearest");
%! assert (fi, x + y + z);

%!test
%! [x,y,z] = ndgrid (0:2);
%! f = x.^2 + y.^2 + z.^2;
%! assert (interpn (x,y,-z,f,1.5,1.5,-1.5), 7.5);

%!test  # for Matlab-compatible rounding for "nearest"
%! x = meshgrid (1:4);
%! assert (interpn (x, 2.5, 2.5, "nearest"), 3);

%!test
%! z = zeros (3, 3, 3);
%! zout = zeros (5, 5, 5);
%! z(:,:,1) = [1 3 5; 3 5 7; 5 7 9];
%! z(:,:,2) = z(:,:,1) + 2;
%! z(:,:,3) = z(:,:,2) + 2;
%! for n = 1:5
%!   zout(:,:,n) = [1 2 3 4 5;
%!                  2 3 4 5 6;
%!                  3 4 5 6 7;
%!                  4 5 6 7 8;
%!                  5 6 7 8 9] + (n-1);
%! endfor
%! tol = 10*eps;
%! assert (interpn (z), zout, tol);
%! assert (interpn (z, "linear"), zout, tol);
%! assert (interpn (z, "spline"), zout, tol);