PageRenderTime 52ms CodeModel.GetById 19ms RepoModel.GetById 0ms app.codeStats 0ms

/octave-3.6.2/scripts/general/interp2.m

#
Objective C | 610 lines | 547 code | 63 blank | 0 comment | 81 complexity | 7b378071460a9ccb2d26057d4b279319 MD5 | raw file
Possible License(s): GPL-3.0
  1. ## Copyright (C) 2000-2012 Kai Habel
  2. ## Copyright (C) 2009 Jaroslav Hajek
  3. ##
  4. ## This file is part of Octave.
  5. ##
  6. ## Octave is free software; you can redistribute it and/or modify it
  7. ## under the terms of the GNU General Public License as published by
  8. ## the Free Software Foundation; either version 3 of the License, or (at
  9. ## your option) any later version.
  10. ##
  11. ## Octave is distributed in the hope that it will be useful, but
  12. ## WITHOUT ANY WARRANTY; without even the implied warranty of
  13. ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
  14. ## General Public License for more details.
  15. ##
  16. ## You should have received a copy of the GNU General Public License
  17. ## along with Octave; see the file COPYING. If not, see
  18. ## <http://www.gnu.org/licenses/>.
  19. ## -*- texinfo -*-
  20. ## @deftypefn {Function File} {@var{zi} =} interp2 (@var{x}, @var{y}, @var{z}, @var{xi}, @var{yi})
  21. ## @deftypefnx {Function File} {@var{zi} =} interp2 (@var{Z}, @var{xi}, @var{yi})
  22. ## @deftypefnx {Function File} {@var{zi} =} interp2 (@var{Z}, @var{n})
  23. ## @deftypefnx {Function File} {@var{zi} =} interp2 (@dots{}, @var{method})
  24. ## @deftypefnx {Function File} {@var{zi} =} interp2 (@dots{}, @var{method}, @var{extrapval})
  25. ##
  26. ## Two-dimensional interpolation. @var{x}, @var{y} and @var{z} describe a
  27. ## surface function. If @var{x} and @var{y} are vectors their length
  28. ## must correspondent to the size of @var{z}. @var{x} and @var{y} must be
  29. ## monotonic. If they are matrices they must have the @code{meshgrid}
  30. ## format.
  31. ##
  32. ## @table @code
  33. ## @item interp2 (@var{x}, @var{y}, @var{Z}, @var{xi}, @var{yi}, @dots{})
  34. ## Returns a matrix corresponding to the points described by the
  35. ## matrices @var{xi}, @var{yi}.
  36. ##
  37. ## If the last argument is a string, the interpolation method can
  38. ## be specified. The method can be 'linear', 'nearest' or 'cubic'.
  39. ## If it is omitted 'linear' interpolation is assumed.
  40. ##
  41. ## @item interp2 (@var{z}, @var{xi}, @var{yi})
  42. ## Assumes @code{@var{x} = 1:rows (@var{z})} and @code{@var{y} =
  43. ## 1:columns (@var{z})}
  44. ##
  45. ## @item interp2 (@var{z}, @var{n})
  46. ## Interleaves the matrix @var{z} n-times. If @var{n} is omitted a value
  47. ## of @code{@var{n} = 1} is assumed.
  48. ## @end table
  49. ##
  50. ## The variable @var{method} defines the method to use for the
  51. ## interpolation. It can take one of the following values
  52. ##
  53. ## @table @asis
  54. ## @item 'nearest'
  55. ## Return the nearest neighbor.
  56. ##
  57. ## @item 'linear'
  58. ## Linear interpolation from nearest neighbors.
  59. ##
  60. ## @item 'pchip'
  61. ## Piecewise cubic Hermite interpolating polynomial.
  62. ##
  63. ## @item 'cubic'
  64. ## Cubic interpolation from four nearest neighbors.
  65. ##
  66. ## @item 'spline'
  67. ## Cubic spline interpolation---smooth first and second derivatives
  68. ## throughout the curve.
  69. ## @end table
  70. ##
  71. ## If a scalar value @var{extrapval} is defined as the final value, then
  72. ## values outside the mesh as set to this value. Note that in this case
  73. ## @var{method} must be defined as well. If @var{extrapval} is not
  74. ## defined then NA is assumed.
  75. ##
  76. ## @seealso{interp1}
  77. ## @end deftypefn
  78. ## Author: Kai Habel <kai.habel@gmx.de>
  79. ## 2005-03-02 Thomas Weber <weber@num.uni-sb.de>
  80. ## * Add test cases
  81. ## 2005-03-02 Paul Kienzle <pkienzle@users.sf.net>
  82. ## * Simplify
  83. ## 2005-04-23 Dmitri A. Sergatskov <dasergatskov@gmail.com>
  84. ## * Modified demo and test for new gnuplot interface
  85. ## 2005-09-07 Hoxide <hoxide_dirac@yahoo.com.cn>
  86. ## * Add bicubic interpolation method
  87. ## * Fix the eat line bug when the last element of XI or YI is
  88. ## negative or zero.
  89. ## 2005-11-26 Pierre Baldensperger <balden@libertysurf.fr>
  90. ## * Rather big modification (XI,YI no longer need to be
  91. ## "meshgridded") to be consistent with the help message
  92. ## above and for compatibility.
  93. function ZI = interp2 (varargin)
  94. Z = X = Y = XI = YI = n = [];
  95. method = "linear";
  96. extrapval = NA;
  97. switch (nargin)
  98. case 1
  99. Z = varargin{1};
  100. n = 1;
  101. case 2
  102. if (ischar (varargin{2}))
  103. [Z, method] = deal (varargin{:});
  104. n = 1;
  105. else
  106. [Z, n] = deal (varargin{:});
  107. endif
  108. case 3
  109. if (ischar (varargin{3}))
  110. [Z, n, method] = deal (varargin{:});
  111. else
  112. [Z, XI, YI] = deal (varargin{:});
  113. endif
  114. case 4
  115. if (ischar (varargin{4}))
  116. [Z, XI, YI, method] = deal (varargin{:});
  117. else
  118. [Z, n, method, extrapval] = deal (varargin{:});
  119. endif
  120. case 5
  121. if (ischar (varargin{4}))
  122. [Z, XI, YI, method, extrapval] = deal (varargin{:});
  123. else
  124. [X, Y, Z, XI, YI] = deal (varargin{:});
  125. endif
  126. case 6
  127. [X, Y, Z, XI, YI, method] = deal (varargin{:});
  128. case 7
  129. [X, Y, Z, XI, YI, method, extrapval] = deal (varargin{:});
  130. otherwise
  131. print_usage ();
  132. endswitch
  133. ## Type checking.
  134. if (!ismatrix (Z))
  135. error ("interp2: Z must be a matrix");
  136. endif
  137. if (!isempty (n) && !isscalar (n))
  138. error ("interp2: N must be a scalar");
  139. endif
  140. if (!ischar (method))
  141. error ("interp2: METHOD must be a string");
  142. endif
  143. if (ischar (extrapval) || strcmp (extrapval, "extrap"))
  144. extrapval = [];
  145. elseif (!isscalar (extrapval))
  146. error ("interp2: EXTRAPVAL must be a scalar");
  147. endif
  148. ## Define X, Y, XI, YI if needed
  149. [zr, zc] = size (Z);
  150. if (isempty (X))
  151. X = 1:zc;
  152. Y = 1:zr;
  153. endif
  154. if (! isnumeric (X) || ! isnumeric (Y))
  155. error ("interp2: X, Y must be numeric matrices");
  156. endif
  157. if (! isempty (n))
  158. ## Calculate the interleaved input vectors.
  159. p = 2^n;
  160. XI = (p:p*zc)/p;
  161. YI = (p:p*zr)'/p;
  162. endif
  163. if (! isnumeric (XI) || ! isnumeric (YI))
  164. error ("interp2: XI, YI must be numeric");
  165. endif
  166. if (strcmp (method, "linear") || strcmp (method, "nearest") ...
  167. || strcmp (method, "pchip"))
  168. ## If X and Y vectors produce a grid from them
  169. if (isvector (X) && isvector (Y))
  170. X = X(:); Y = Y(:);
  171. elseif (size_equal (X, Y))
  172. X = X(1,:)'; Y = Y(:,1);
  173. else
  174. error ("interp2: X and Y must be matrices of same size");
  175. endif
  176. if (columns (Z) != length (X) || rows (Z) != length (Y))
  177. error ("interp2: X and Y size must match the dimensions of Z");
  178. endif
  179. ## If Xi and Yi are vectors of different orientation build a grid
  180. if ((rows (XI) == 1 && columns (YI) == 1)
  181. || (columns (XI) == 1 && rows (YI) == 1))
  182. [XI, YI] = meshgrid (XI, YI);
  183. elseif (! size_equal (XI, YI))
  184. error ("interp2: XI and YI must be matrices of equal size");
  185. endif
  186. ## if XI, YI are vectors, X and Y should share their orientation.
  187. if (rows (XI) == 1)
  188. if (rows (X) != 1)
  189. X = X.';
  190. endif
  191. if (rows (Y) != 1)
  192. Y = Y.';
  193. endif
  194. elseif (columns (XI) == 1)
  195. if (columns (X) != 1)
  196. X = X.';
  197. endif
  198. if (columns (Y) != 1)
  199. Y = Y.';
  200. endif
  201. endif
  202. xidx = lookup (X, XI, "lr");
  203. yidx = lookup (Y, YI, "lr");
  204. if (strcmp (method, "linear"))
  205. ## each quad satisfies the equation z(x,y)=a+b*x+c*y+d*xy
  206. ##
  207. ## a-b
  208. ## | |
  209. ## c-d
  210. a = Z(1:(zr - 1), 1:(zc - 1));
  211. b = Z(1:(zr - 1), 2:zc) - a;
  212. c = Z(2:zr, 1:(zc - 1)) - a;
  213. d = Z(2:zr, 2:zc) - a - b - c;
  214. ## scale XI, YI values to a 1-spaced grid
  215. Xsc = (XI - X(xidx)) ./ (diff (X)(xidx));
  216. Ysc = (YI - Y(yidx)) ./ (diff (Y)(yidx));
  217. ## Get 2D index.
  218. idx = sub2ind (size (a), yidx, xidx);
  219. ## We can dispose of the 1D indices at this point to save memory.
  220. clear xidx yidx;
  221. ## apply plane equation
  222. ZI = a(idx) + b(idx).*Xsc + c(idx).*Ysc + d(idx).*Xsc.*Ysc;
  223. elseif (strcmp (method, "nearest"))
  224. ii = (XI - X(xidx) >= X(xidx + 1) - XI);
  225. jj = (YI - Y(yidx) >= Y(yidx + 1) - YI);
  226. idx = sub2ind (size (Z), yidx+jj, xidx+ii);
  227. ZI = Z(idx);
  228. elseif (strcmp (method, "pchip"))
  229. if (length (X) < 2 || length (Y) < 2)
  230. error ("interp2: pchip2 requires at least 2 points in each dimension");
  231. endif
  232. ## first order derivatives
  233. DX = __pchip_deriv__ (X, Z, 2);
  234. DY = __pchip_deriv__ (Y, Z, 1);
  235. ## Compute mixed derivatives row-wise and column-wise, use the average.
  236. DXY = (__pchip_deriv__ (X, DY, 2) + __pchip_deriv__ (Y, DX, 1))/2;
  237. ## do the bicubic interpolation
  238. hx = diff (X); hx = hx(xidx);
  239. hy = diff (Y); hy = hy(yidx);
  240. tx = (XI - X(xidx)) ./ hx;
  241. ty = (YI - Y(yidx)) ./ hy;
  242. ## construct the cubic hermite base functions in x, y
  243. ## formulas:
  244. ## b{1,1} = ( 2*t.^3 - 3*t.^2 + 1);
  245. ## b{2,1} = h.*( t.^3 - 2*t.^2 + t );
  246. ## b{1,2} = (-2*t.^3 + 3*t.^2 );
  247. ## b{2,2} = h.*( t.^3 - t.^2 );
  248. ## optimized equivalents of the above:
  249. t1 = tx.^2;
  250. t2 = tx.*t1 - t1;
  251. xb{2,2} = hx.*t2;
  252. t1 = t2 - t1;
  253. xb{2,1} = hx.*(t1 + tx);
  254. t2 += t1;
  255. xb{1,2} = -t2;
  256. xb{1,1} = t2 + 1;
  257. t1 = ty.^2;
  258. t2 = ty.*t1 - t1;
  259. yb{2,2} = hy.*t2;
  260. t1 = t2 - t1;
  261. yb{2,1} = hy.*(t1 + ty);
  262. t2 += t1;
  263. yb{1,2} = -t2;
  264. yb{1,1} = t2 + 1;
  265. ZI = zeros (size (XI));
  266. for i = 1:2
  267. for j = 1:2
  268. zidx = sub2ind (size (Z), yidx+(j-1), xidx+(i-1));
  269. ZI += xb{1,i} .* yb{1,j} .* Z(zidx);
  270. ZI += xb{2,i} .* yb{1,j} .* DX(zidx);
  271. ZI += xb{1,i} .* yb{2,j} .* DY(zidx);
  272. ZI += xb{2,i} .* yb{2,j} .* DXY(zidx);
  273. endfor
  274. endfor
  275. endif
  276. if (! isempty (extrapval))
  277. ## set points outside the table to 'extrapval'
  278. if (X (1) < X (end))
  279. if (Y (1) < Y (end))
  280. ZI (XI < X(1,1) | XI > X(end) | YI < Y(1,1) | YI > Y(end)) = ...
  281. extrapval;
  282. else
  283. ZI (XI < X(1) | XI > X(end) | YI < Y(end) | YI > Y(1)) = ...
  284. extrapval;
  285. endif
  286. else
  287. if (Y (1) < Y (end))
  288. ZI (XI < X(end) | XI > X(1) | YI < Y(1) | YI > Y(end)) = ...
  289. extrapval;
  290. else
  291. ZI (XI < X(1,end) | XI > X(1) | YI < Y(end) | YI > Y(1)) = ...
  292. extrapval;
  293. endif
  294. endif
  295. endif
  296. else
  297. ## Check dimensions of X and Y
  298. if (isvector (X) && isvector (Y))
  299. X = X(:).';
  300. Y = Y(:);
  301. if (!isequal ([length(Y), length(X)], size(Z)))
  302. error ("interp2: X and Y size must match the dimensions of Z");
  303. endif
  304. elseif (!size_equal (X, Y))
  305. error ("interp2: X and Y must be matrices of equal size");
  306. if (! size_equal (X, Z))
  307. error ("interp2: X and Y size must match the dimensions of Z");
  308. endif
  309. endif
  310. ## Check dimensions of XI and YI
  311. if (isvector (XI) && isvector (YI) && ! size_equal (XI, YI))
  312. XI = XI(:).';
  313. YI = YI(:);
  314. [XI, YI] = meshgrid (XI, YI);
  315. elseif (! size_equal (XI, YI))
  316. error ("interp2: XI and YI must be matrices of equal size");
  317. endif
  318. if (strcmp (method, "cubic"))
  319. if (isgriddata (XI) && isgriddata (YI'))
  320. ZI = bicubic (X, Y, Z, XI (1, :), YI (:, 1), extrapval);
  321. elseif (isgriddata (X) && isgriddata (Y'))
  322. ## Allocate output
  323. ZI = zeros (size (X));
  324. ## Find inliers
  325. inside = !(XI < X (1) | XI > X (end) | YI < Y (1) | YI > Y (end));
  326. ## Scale XI and YI to match indices of Z
  327. XI = (columns (Z) - 1) * (XI - X (1)) / (X (end) - X (1)) + 1;
  328. YI = (rows (Z) - 1) * (YI - Y (1)) / (Y (end) - Y (1)) + 1;
  329. ## Start the real work
  330. K = floor (XI);
  331. L = floor (YI);
  332. ## Coefficients
  333. AY1 = bc ((YI - L + 1));
  334. AX1 = bc ((XI - K + 1));
  335. AY0 = bc ((YI - L + 0));
  336. AX0 = bc ((XI - K + 0));
  337. AY_1 = bc ((YI - L - 1));
  338. AX_1 = bc ((XI - K - 1));
  339. AY_2 = bc ((YI - L - 2));
  340. AX_2 = bc ((XI - K - 2));
  341. ## Perform interpolation
  342. sz = size(Z);
  343. ZI = AY_2 .* AX_2 .* Z (sym_sub2ind (sz, L+2, K+2)) ...
  344. + AY_2 .* AX_1 .* Z (sym_sub2ind (sz, L+2, K+1)) ...
  345. + AY_2 .* AX0 .* Z (sym_sub2ind (sz, L+2, K)) ...
  346. + AY_2 .* AX1 .* Z (sym_sub2ind (sz, L+2, K-1)) ...
  347. + AY_1 .* AX_2 .* Z (sym_sub2ind (sz, L+1, K+2)) ...
  348. + AY_1 .* AX_1 .* Z (sym_sub2ind (sz, L+1, K+1)) ...
  349. + AY_1 .* AX0 .* Z (sym_sub2ind (sz, L+1, K)) ...
  350. + AY_1 .* AX1 .* Z (sym_sub2ind (sz, L+1, K-1)) ...
  351. + AY0 .* AX_2 .* Z (sym_sub2ind (sz, L, K+2)) ...
  352. + AY0 .* AX_1 .* Z (sym_sub2ind (sz, L, K+1)) ...
  353. + AY0 .* AX0 .* Z (sym_sub2ind (sz, L, K)) ...
  354. + AY0 .* AX1 .* Z (sym_sub2ind (sz, L, K-1)) ...
  355. + AY1 .* AX_2 .* Z (sym_sub2ind (sz, L-1, K+2)) ...
  356. + AY1 .* AX_1 .* Z (sym_sub2ind (sz, L-1, K+1)) ...
  357. + AY1 .* AX0 .* Z (sym_sub2ind (sz, L-1, K)) ...
  358. + AY1 .* AX1 .* Z (sym_sub2ind (sz, L-1, K-1));
  359. ZI (!inside) = extrapval;
  360. else
  361. error ("interp2: input data must have `meshgrid' format");
  362. endif
  363. elseif (strcmp (method, "spline"))
  364. if (isgriddata (XI) && isgriddata (YI'))
  365. ZI = __splinen__ ({Y(:,1).', X(1,:)}, Z, {YI(:,1), XI(1,:)}, extrapval,
  366. "spline");
  367. else
  368. error ("interp2: input data must have `meshgrid' format");
  369. endif
  370. else
  371. error ("interp2: interpolation METHOD not recognized");
  372. endif
  373. endif
  374. endfunction
  375. function b = isgriddata (X)
  376. d1 = diff (X, 1, 1);
  377. b = all (d1 (:) == 0);
  378. endfunction
  379. ## Compute the bicubic interpolation coefficients
  380. function o = bc(x)
  381. x = abs(x);
  382. o = zeros(size(x));
  383. idx1 = (x < 1);
  384. idx2 = !idx1 & (x < 2);
  385. o(idx1) = 1 - 2.*x(idx1).^2 + x(idx1).^3;
  386. o(idx2) = 4 - 8.*x(idx2) + 5.*x(idx2).^2 - x(idx2).^3;
  387. endfunction
  388. ## This version of sub2ind behaves as if the data was symmetrically padded
  389. function ind = sym_sub2ind(sz, Y, X)
  390. Y (Y < 1) = 1 - Y (Y < 1);
  391. while (any (Y (:) > 2 * sz (1)))
  392. Y (Y > 2 * sz (1)) = round (Y (Y > 2 * sz (1)) / 2);
  393. endwhile
  394. Y (Y > sz (1)) = 1 + 2 * sz (1) - Y (Y > sz (1));
  395. X (X < 1) = 1 - X (X < 1);
  396. while (any (X (:) > 2 * sz (2)))
  397. X (X > 2 * sz (2)) = round (X (X > 2 * sz (2)) / 2);
  398. endwhile
  399. X (X > sz (2)) = 1 + 2 * sz (2) - X (X > sz (2));
  400. ind = sub2ind(sz, Y, X);
  401. endfunction
  402. %!demo
  403. %! A=[13,-1,12;5,4,3;1,6,2];
  404. %! x=[0,1,4]; y=[10,11,12];
  405. %! xi=linspace(min(x),max(x),17);
  406. %! yi=linspace(min(y),max(y),26)';
  407. %! mesh(xi,yi,interp2(x,y,A,xi,yi,'linear'));
  408. %! [x,y] = meshgrid(x,y);
  409. %! hold on; plot3(x(:),y(:),A(:),"b*"); hold off;
  410. %!demo
  411. %! [x,y,A] = peaks(10);
  412. %! x = x(1,:)'; y = y(:,1);
  413. %! xi=linspace(min(x),max(x),41);
  414. %! yi=linspace(min(y),max(y),41)';
  415. %! mesh(xi,yi,interp2(x,y,A,xi,yi,'linear'));
  416. %! [x,y] = meshgrid(x,y);
  417. %! hold on; plot3(x(:),y(:),A(:),"b*"); hold off;
  418. %!demo
  419. %! A=[13,-1,12;5,4,3;1,6,2];
  420. %! x=[0,1,4]; y=[10,11,12];
  421. %! xi=linspace(min(x),max(x),17);
  422. %! yi=linspace(min(y),max(y),26)';
  423. %! mesh(xi,yi,interp2(x,y,A,xi,yi,'nearest'));
  424. %! [x,y] = meshgrid(x,y);
  425. %! hold on; plot3(x(:),y(:),A(:),"b*"); hold off;
  426. %!demo
  427. %! [x,y,A] = peaks(10);
  428. %! x = x(1,:)'; y = y(:,1);
  429. %! xi=linspace(min(x),max(x),41);
  430. %! yi=linspace(min(y),max(y),41)';
  431. %! mesh(xi,yi,interp2(x,y,A,xi,yi,'nearest'));
  432. %! [x,y] = meshgrid(x,y);
  433. %! hold on; plot3(x(:),y(:),A(:),"b*"); hold off;
  434. %!demo
  435. %! A=[13,-1,12;5,4,3;1,6,2];
  436. %! x=[0,1,2]; y=[10,11,12];
  437. %! xi=linspace(min(x),max(x),17);
  438. %! yi=linspace(min(y),max(y),26)';
  439. %! mesh(xi,yi,interp2(x,y,A,xi,yi,'pchip'));
  440. %! [x,y] = meshgrid(x,y);
  441. %! hold on; plot3(x(:),y(:),A(:),"b*"); hold off;
  442. %!demo
  443. %! [x,y,A] = peaks(10);
  444. %! x = x(1,:)'; y = y(:,1);
  445. %! xi=linspace(min(x),max(x),41);
  446. %! yi=linspace(min(y),max(y),41)';
  447. %! mesh(xi,yi,interp2(x,y,A,xi,yi,'pchip'));
  448. %! [x,y] = meshgrid(x,y);
  449. %! hold on; plot3(x(:),y(:),A(:),"b*"); hold off;
  450. %!demo
  451. %! A=[13,-1,12;5,4,3;1,6,2];
  452. %! x=[0,1,2]; y=[10,11,12];
  453. %! xi=linspace(min(x),max(x),17);
  454. %! yi=linspace(min(y),max(y),26)';
  455. %! mesh(xi,yi,interp2(x,y,A,xi,yi,'cubic'));
  456. %! [x,y] = meshgrid(x,y);
  457. %! hold on; plot3(x(:),y(:),A(:),"b*"); hold off;
  458. %!demo
  459. %! [x,y,A] = peaks(10);
  460. %! x = x(1,:)'; y = y(:,1);
  461. %! xi=linspace(min(x),max(x),41);
  462. %! yi=linspace(min(y),max(y),41)';
  463. %! mesh(xi,yi,interp2(x,y,A,xi,yi,'cubic'));
  464. %! [x,y] = meshgrid(x,y);
  465. %! hold on; plot3(x(:),y(:),A(:),"b*"); hold off;
  466. %!demo
  467. %! A=[13,-1,12;5,4,3;1,6,2];
  468. %! x=[0,1,2]; y=[10,11,12];
  469. %! xi=linspace(min(x),max(x),17);
  470. %! yi=linspace(min(y),max(y),26)';
  471. %! mesh(xi,yi,interp2(x,y,A,xi,yi,'spline'));
  472. %! [x,y] = meshgrid(x,y);
  473. %! hold on; plot3(x(:),y(:),A(:),"b*"); hold off;
  474. %!demo
  475. %! [x,y,A] = peaks(10);
  476. %! x = x(1,:)'; y = y(:,1);
  477. %! xi=linspace(min(x),max(x),41);
  478. %! yi=linspace(min(y),max(y),41)';
  479. %! mesh(xi,yi,interp2(x,y,A,xi,yi,'spline'));
  480. %! [x,y] = meshgrid(x,y);
  481. %! hold on; plot3(x(:),y(:),A(:),"b*"); hold off;
  482. %!test % simple test
  483. %! x = [1,2,3];
  484. %! y = [4,5,6,7];
  485. %! [X, Y] = meshgrid(x,y);
  486. %! Orig = X.^2 + Y.^3;
  487. %! xi = [1.2,2, 1.5];
  488. %! yi = [6.2, 4.0, 5.0]';
  489. %!
  490. %! Expected = ...
  491. %! [243, 245.4, 243.9;
  492. %! 65.6, 68, 66.5;
  493. %! 126.6, 129, 127.5];
  494. %! Result = interp2(x,y,Orig, xi, yi);
  495. %!
  496. %! assert(Result, Expected, 1000*eps);
  497. %!test % 2^n form
  498. %! x = [1,2,3];
  499. %! y = [4,5,6,7];
  500. %! [X, Y] = meshgrid(x,y);
  501. %! Orig = X.^2 + Y.^3;
  502. %! xi = [1:0.25:3]; yi = [4:0.25:7]';
  503. %! Expected = interp2(x,y,Orig, xi, yi);
  504. %! Result = interp2(Orig,2);
  505. %!
  506. %! assert(Result, Expected, 10*eps);
  507. %!test % matrix slice
  508. %! A = eye(4);
  509. %! assert(interp2(A,[1:4],[1:4]),[1,1,1,1]);
  510. %!test % non-gridded XI,YI
  511. %! A = eye(4);
  512. %! assert(interp2(A,[1,2;3,4],[1,3;2,4]),[1,0;0,1]);
  513. %!test % for values outside of boundaries
  514. %! x = [1,2,3];
  515. %! y = [4,5,6,7];
  516. %! [X, Y] = meshgrid(x,y);
  517. %! Orig = X.^2 + Y.^3;
  518. %! xi = [0,4];
  519. %! yi = [3,8]';
  520. %! assert(interp2(x,y,Orig, xi, yi),[NA,NA;NA,NA]);
  521. %! assert(interp2(x,y,Orig, xi, yi,'linear', 0),[0,0;0,0]);
  522. %!test % for values at boundaries
  523. %! A=[1,2;3,4];
  524. %! x=[0,1];
  525. %! y=[2,3]';
  526. %! assert(interp2(x,y,A,x,y,'linear'), A);
  527. %! assert(interp2(x,y,A,x,y,'nearest'), A);
  528. %!test % for Matlab-compatible rounding for 'nearest'
  529. %! X = meshgrid (1:4);
  530. %! assert (interp2 (X, 2.5, 2.5, 'nearest'), 3);
  531. %!shared z, zout, tol
  532. %! z = [1 3 5; 3 5 7; 5 7 9];
  533. %! zout = [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];
  534. %! tol = 2 * eps;
  535. %!assert (interp2 (z), zout, tol);
  536. %!assert (interp2 (z, "linear"), zout, tol);
  537. %!assert (interp2 (z, "pchip"), zout, tol);
  538. %!assert (interp2 (z, "cubic"), zout, 10 * tol);
  539. %!assert (interp2 (z, "spline"), zout, tol);
  540. %!assert (interp2 (z, [2 3 1], [2 2 2]', "linear"), repmat ([5, 7, 3], [3, 1]), tol)
  541. %!assert (interp2 (z, [2 3 1], [2 2 2]', "pchip"), repmat ([5, 7, 3], [3, 1]), tol)
  542. %!assert (interp2 (z, [2 3 1], [2 2 2]', "cubic"), repmat ([5, 7, 3], [3, 1]), 10 * tol)
  543. %!assert (interp2 (z, [2 3 1], [2 2 2]', "spline"), repmat ([5, 7, 3], [3, 1]), tol)
  544. %!assert (interp2 (z, [2 3 1], [2 2 2], "linear"), [5 7 3], tol);
  545. %!assert (interp2 (z, [2 3 1], [2 2 2], "pchip"), [5 7 3], tol);
  546. %!assert (interp2 (z, [2 3 1], [2 2 2], "cubic"), [5 7 3], 10 * tol);
  547. %!assert (interp2 (z, [2 3 1], [2 2 2], "spline"), [5 7 3], tol);