changelog shortlog tags changeset files revisions annotate raw

scripts/geometry/griddata3.m

changeset 10289: 4b124317dc38
parent:ec0a13863eb7
author: John W. Eaton <jwe@octave.org>
date: Tue Feb 09 20:58:55 2010 -0500 (31 minutes ago)
permissions: -rw-r--r--
description: base_properties::set_children: account for hidden children
1## Copyright (C) 2007, 2008 David Bateman
2##
3## This file is part of Octave.
4##
5## Octave is free software; you can redistribute it and/or modify it
6## under the terms of the GNU General Public License as published by
7## the Free Software Foundation; either version 3 of the License, or (at
8## your option) any later version.
9##
10## Octave is distributed in the hope that it will be useful, but
11## WITHOUT ANY WARRANTY; without even the implied warranty of
12## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13## General Public License for more details.
14##
15## You should have received a copy of the GNU General Public License
16## along with Octave; see the file COPYING. If not, see
17## <http://www.gnu.org/licenses/>.
18
19## -*- texinfo -*-
20## @deftypefn {Function File} {@var{vi} =} griddata3 (@var{x}, @var{y}, @var{z}, @var{v} @var{xi}, @var{yi}, @var{zi}, @var{method}, @var{options})
21##
22## Generate a regular mesh from irregular data using interpolation.
23## The function is defined by @code{@var{y} = f (@var{x},@var{y},@var{z})}.
24## The interpolation points are all @var{xi}.
25##
26## The interpolation method can be @code{"nearest"} or @code{"linear"}.
27## If method is omitted it defaults to @code{"linear"}.
28## @seealso{griddata, delaunayn}
29## @end deftypefn
30
31## Author: David Bateman <dbateman@free.fr>
32
33function vi = griddata3 (x, y, z, v, xi, yi, zi, method, varargin)
34
35 if (nargin < 7)
36 print_usage ();
37 endif
38
39 if (!all (size (x) == size (y) & size (x) == size(z) & size(x) == size (v)))
40 error ("griddata3: x, y, z, and v must be vectors of same length");
41 endif
42
43 ## meshgrid xi, yi and zi if they are vectors unless they
44 ## are vectors of the same length
45 if (isvector (xi) && isvector (yi) && isvector (zi)
46 && (numel (xi) != numel (yi) || numel (xi) != numel (zi)))
47 [xi, yi, zi] = meshgrid (xi, yi, zi);
48 endif
49
50 if (any (size(xi) != size(yi)) || any (size(xi) != size(zi)))
51 error ("griddata3: xi, yi and zi must be vectors or matrices of same size");
52 endif
53
54 vi = griddatan ([x(:), y(:), z(:)], v(:), [xi(:), yi(:), zi(:)], varargin{:});
55 vi = reshape (vi, size (xi));
56endfunction
57
58%!testif HAVE_QHULL
59%! rand('state', 0);
60%! x = 2 * rand(1000, 1) - 1;
61%! y = 2 * rand(1000, 1) - 1;
62%! z = 2 * rand(1000, 1) - 1;
63%! v = x.^2 + y.^2 + z.^2;
64%! [xi, yi, zi] = meshgrid (-0.8:0.2:0.8);
65%! ##vi = reshape (griddatan([x(:), y(:), z(:)], v, [xi(:), yi(:), zi(:)], 'linear'), size (xi));
66%! vi = griddata3 (x, y, z, v, xi, yi, zi, 'linear');
67%! vv = vi - xi.^2 - yi.^2 - zi.^2;
68%! assert (max(abs(vv(:))), 0, 0.1)
69
70%!testif HAVE_QHULL
71%! rand('state', 0);
72%! x = 2 * rand(1000, 1) - 1;
73%! y = 2 * rand(1000, 1) - 1;
74%! z = 2 * rand(1000, 1) - 1;
75%! v = x.^2 + y.^2 + z.^2;
76%! [xi, yi, zi] = meshgrid (-0.8:0.2:0.8);
77%! ##vi = reshape (griddatan([x(:), y(:), z(:)], v, [xi(:), yi(:), zi(:)], 'linear'), size (xi));
78%! vi = griddata3 (x, y, z, v, xi, yi, zi, 'nearest');
79%! vv = vi - xi.^2 - yi.^2 - zi.^2;
80%! assert (max(abs(vv(:))), 0, 0.1)