Mercurial > octave-nkf
view scripts/general/interpn.m @ 17281:bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Macro handles options ("on") or properties ("position") more elegantly
than @code{"text"}.
* doc/interpreter/macros.texi: Add new @qcode macro.
* doc/interpreter/tips.txi: Add documentation about @qcode macro.
* doc/interpreter/basics.txi, doc/interpreter/container.txi,
doc/interpreter/emacs.txi, doc/interpreter/errors.txi,
doc/interpreter/eval.txi, doc/interpreter/expr.txi,
doc/interpreter/external.txi, doc/interpreter/func.txi,
doc/interpreter/grammar.txi, doc/interpreter/image.txi,
doc/interpreter/install.txi, doc/interpreter/interp.txi,
doc/interpreter/io.txi, doc/interpreter/matrix.txi,
doc/interpreter/numbers.txi, doc/interpreter/oop.txi,
doc/interpreter/package.txi, doc/interpreter/plot.txi,
doc/interpreter/quad.txi, doc/interpreter/sparse.txi,
doc/interpreter/strings.txi, doc/interpreter/system.txi,
doc/interpreter/vectorize.txi, libinterp/corefcn/balance.cc,
libinterp/corefcn/bitfcns.cc, libinterp/corefcn/cellfun.cc,
libinterp/corefcn/conv2.cc, libinterp/corefcn/data.cc,
libinterp/corefcn/debug.cc, libinterp/corefcn/defaults.cc,
libinterp/corefcn/dirfns.cc, libinterp/corefcn/dlmread.cc,
libinterp/corefcn/error.cc, libinterp/corefcn/file-io.cc,
libinterp/corefcn/find.cc, libinterp/corefcn/gammainc.cc,
libinterp/corefcn/graphics.cc, libinterp/corefcn/help.cc,
libinterp/corefcn/hex2num.cc, libinterp/corefcn/input.cc,
libinterp/corefcn/load-path.cc, libinterp/corefcn/load-save.cc,
libinterp/corefcn/ls-oct-ascii.cc, libinterp/corefcn/lu.cc,
libinterp/corefcn/luinc.cc, libinterp/corefcn/matrix_type.cc,
libinterp/corefcn/oct-hist.cc, libinterp/corefcn/pager.cc,
libinterp/corefcn/pr-output.cc, libinterp/corefcn/pt-jit.cc,
libinterp/corefcn/qz.cc, libinterp/corefcn/rand.cc,
libinterp/corefcn/regexp.cc, libinterp/corefcn/schur.cc,
libinterp/corefcn/sighandlers.cc, libinterp/corefcn/sparse.cc,
libinterp/corefcn/spparms.cc, libinterp/corefcn/str2double.cc,
libinterp/corefcn/svd.cc, libinterp/corefcn/symtab.cc,
libinterp/corefcn/syscalls.cc, libinterp/corefcn/toplev.cc,
libinterp/corefcn/tril.cc, libinterp/corefcn/typecast.cc,
libinterp/corefcn/utils.cc, libinterp/corefcn/variables.cc,
libinterp/dldfcn/__init_fltk__.cc, libinterp/dldfcn/chol.cc,
libinterp/dldfcn/colamd.cc, libinterp/dldfcn/fftw.cc, libinterp/dldfcn/qr.cc,
libinterp/dldfcn/symbfact.cc, libinterp/octave-value/ov-base.cc,
libinterp/octave-value/ov-fcn-handle.cc,
libinterp/octave-value/ov-fcn-inline.cc, libinterp/octave-value/ov-java.cc,
libinterp/octave-value/ov-range.cc, libinterp/octave-value/ov-struct.cc,
libinterp/octave-value/ov-usr-fcn.cc, libinterp/parse-tree/oct-parse.in.yy,
libinterp/parse-tree/pt-binop.cc, libinterp/parse-tree/pt-eval.cc,
libinterp/parse-tree/pt-mat.cc, scripts/@ftp/ftp.m,
scripts/deprecated/java_convert_matrix.m, scripts/deprecated/java_debug.m,
scripts/deprecated/java_unsigned_conversion.m, scripts/deprecated/shell_cmd.m,
scripts/general/dblquad.m, scripts/general/display.m,
scripts/general/genvarname.m, scripts/general/idivide.m,
scripts/general/interp1.m, scripts/general/interp2.m,
scripts/general/interp3.m, scripts/general/interpn.m, scripts/general/isa.m,
scripts/general/profexplore.m, scripts/general/profile.m,
scripts/general/quadgk.m, scripts/general/randi.m, scripts/general/structfun.m,
scripts/general/subsindex.m, scripts/general/triplequad.m,
scripts/geometry/griddata.m, scripts/geometry/griddata3.m,
scripts/geometry/griddatan.m, scripts/geometry/voronoi.m, scripts/help/help.m,
scripts/help/lookfor.m, scripts/image/cmpermute.m, scripts/image/colormap.m,
scripts/image/image.m, scripts/image/imagesc.m, scripts/image/imfinfo.m,
scripts/image/imformats.m, scripts/image/imread.m, scripts/image/imshow.m,
scripts/image/imwrite.m, scripts/image/ind2gray.m, scripts/image/lines.m,
scripts/image/rgb2ind.m, scripts/image/spinmap.m, scripts/io/dlmwrite.m,
scripts/io/strread.m, scripts/io/textread.m, scripts/io/textscan.m,
scripts/java/javaclasspath.m, scripts/java/usejava.m,
scripts/miscellaneous/bzip2.m, scripts/miscellaneous/computer.m,
scripts/miscellaneous/copyfile.m, scripts/miscellaneous/debug.m,
scripts/miscellaneous/dos.m, scripts/miscellaneous/edit.m,
scripts/miscellaneous/gzip.m, scripts/miscellaneous/license.m,
scripts/miscellaneous/mkoctfile.m, scripts/miscellaneous/movefile.m,
scripts/miscellaneous/parseparams.m, scripts/miscellaneous/unix.m,
scripts/optimization/fminbnd.m, scripts/optimization/fminsearch.m,
scripts/optimization/fminunc.m, scripts/optimization/fsolve.m,
scripts/optimization/fzero.m, scripts/optimization/glpk.m,
scripts/optimization/lsqnonneg.m, scripts/optimization/optimset.m,
scripts/optimization/pqpnonneg.m, scripts/pkg/pkg.m, scripts/plot/allchild.m,
scripts/plot/ancestor.m, scripts/plot/area.m, scripts/plot/axis.m,
scripts/plot/bar.m, scripts/plot/barh.m, scripts/plot/box.m,
scripts/plot/caxis.m, scripts/plot/cla.m, scripts/plot/clabel.m,
scripts/plot/clf.m, scripts/plot/close.m, scripts/plot/colorbar.m,
scripts/plot/daspect.m, scripts/plot/ezmesh.m, scripts/plot/ezmeshc.m,
scripts/plot/ezsurf.m, scripts/plot/ezsurfc.m, scripts/plot/findall.m,
scripts/plot/findobj.m, scripts/plot/gcbo.m, scripts/plot/gcf.m,
scripts/plot/gco.m, scripts/plot/grid.m, scripts/plot/guihandles.m,
scripts/plot/hdl2struct.m, scripts/plot/hidden.m, scripts/plot/hold.m,
scripts/plot/isonormals.m, scripts/plot/isosurface.m, scripts/plot/legend.m,
scripts/plot/mesh.m, scripts/plot/meshc.m, scripts/plot/meshz.m,
scripts/plot/newplot.m, scripts/plot/orient.m, scripts/plot/pareto.m,
scripts/plot/patch.m, scripts/plot/pbaspect.m, scripts/plot/pcolor.m,
scripts/plot/plot.m, scripts/plot/print.m,
scripts/plot/private/__add_default_menu__.m, scripts/plot/quiver.m,
scripts/plot/quiver3.m, scripts/plot/refreshdata.m, scripts/plot/saveas.m,
scripts/plot/scatter.m, scripts/plot/scatter3.m, scripts/plot/shading.m,
scripts/plot/shrinkfaces.m, scripts/plot/slice.m, scripts/plot/stem.m,
scripts/plot/stem3.m, scripts/plot/struct2hdl.m, scripts/plot/subplot.m,
scripts/plot/surf.m, scripts/plot/surfc.m, scripts/plot/surfl.m,
scripts/plot/tetramesh.m, scripts/plot/uigetfile.m, scripts/plot/uimenu.m,
scripts/plot/uiputfile.m, scripts/plot/waterfall.m, scripts/plot/whitebg.m,
scripts/plot/xlim.m, scripts/plot/ylim.m, scripts/plot/zlim.m,
scripts/polynomial/conv.m, scripts/polynomial/polyout.m,
scripts/polynomial/splinefit.m, scripts/set/ismember.m, scripts/set/powerset.m,
scripts/set/setdiff.m, scripts/set/union.m, scripts/set/unique.m,
scripts/signal/detrend.m, scripts/signal/filter2.m, scripts/signal/freqz.m,
scripts/signal/periodogram.m, scripts/signal/spectral_adf.m,
scripts/signal/spectral_xdf.m, scripts/sparse/eigs.m, scripts/sparse/svds.m,
scripts/specfun/legendre.m, scripts/special-matrix/gallery.m,
scripts/statistics/base/mean.m, scripts/statistics/base/moment.m,
scripts/statistics/tests/cor_test.m,
scripts/statistics/tests/kolmogorov_smirnov_test.m,
scripts/statistics/tests/kolmogorov_smirnov_test_2.m,
scripts/statistics/tests/kruskal_wallis_test.m,
scripts/statistics/tests/prop_test_2.m, scripts/statistics/tests/sign_test.m,
scripts/statistics/tests/t_test.m, scripts/statistics/tests/t_test_2.m,
scripts/statistics/tests/t_test_regression.m,
scripts/statistics/tests/u_test.m, scripts/statistics/tests/var_test.m,
scripts/statistics/tests/welch_test.m,
scripts/statistics/tests/wilcoxon_test.m, scripts/statistics/tests/z_test.m,
scripts/statistics/tests/z_test_2.m, scripts/strings/base2dec.m,
scripts/strings/index.m, scripts/strings/isstrprop.m,
scripts/strings/mat2str.m, scripts/strings/regexptranslate.m,
scripts/strings/rindex.m, scripts/strings/str2num.m, scripts/strings/strcat.m,
scripts/strings/strjust.m, scripts/strings/strmatch.m,
scripts/strings/validatestring.m, scripts/testfun/demo.m,
scripts/testfun/example.m, scripts/testfun/test.m, scripts/time/addtodate.m,
scripts/time/asctime.m, scripts/time/datestr.m, scripts/time/datetick.m,
scripts/time/weekday.m, scripts/ui/errordlg.m, scripts/ui/helpdlg.m,
scripts/ui/inputdlg.m, scripts/ui/listdlg.m, scripts/ui/msgbox.m,
scripts/ui/questdlg.m, scripts/ui/warndlg.m: Use new @qcode macro.
author | Rik <rik@octave.org> |
---|---|
date | Mon, 19 Aug 2013 20:46:38 -0700 |
parents | 486c3e2731ff |
children | d63878346099 |
line wrap: on
line source
## Copyright (C) 2007-2012 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);