annotate scripts/geometry/griddatan.m @ 15063:36cbcc37fdb8

Refactor configure.ac to make it more understandable. Use common syntax for messages in config.h Correct typos, refer to libraries in all caps, use two spaces after period. Follow Autoconf guidelines and place general tests before specific tests. * configure.ac, m4/acinclude.m4: Use common syntax for messages in config.h Correct typos, refer to libraries in all caps, use two spaces after period. Follow Autoconf guidelines and place general tests before specific tests.
author Rik <rik@octave.org>
date Tue, 31 Jul 2012 10:28:51 -0700
parents fbdee844550c
children bc924baa2c4e
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
14138
72c96de7a403 maint: update copyright notices for 2012
John W. Eaton <jwe@octave.org>
parents: 11587
diff changeset
1 ## Copyright (C) 2007-2012 David Bateman
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
2 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
3 ## This file is part of Octave.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
4 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
5 ## Octave is free software; you can redistribute it and/or modify it
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
6 ## under the terms of the GNU General Public License as published by
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6826
diff changeset
7 ## the Free Software Foundation; either version 3 of the License, or (at
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6826
diff changeset
8 ## your option) any later version.
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
9 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
10 ## Octave is distributed in the hope that it will be useful, but
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
11 ## WITHOUT ANY WARRANTY; without even the implied warranty of
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
12 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
13 ## General Public License for more details.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
14 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
15 ## You should have received a copy of the GNU General Public License
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6826
diff changeset
16 ## along with Octave; see the file COPYING. If not, see
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6826
diff changeset
17 ## <http://www.gnu.org/licenses/>.
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
18
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
19 ## -*- texinfo -*-
14370
fbdee844550c griddatan.m: Restore ability to pass options to underlying Qhull
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
20 ## @deftypefn {Function File} {@var{yi} =} griddatan (@var{x}, @var{y}, @var{xi})
fbdee844550c griddatan.m: Restore ability to pass options to underlying Qhull
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
21 ## @deftypefnx {Function File} {@var{yi} =} griddatan (@var{x}, @var{y}, @var{xi}, @var{method})
fbdee844550c griddatan.m: Restore ability to pass options to underlying Qhull
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
22 ## @deftypefnx {Function File} {@var{yi} =} griddatan (@var{x}, @var{y}, @var{xi}, @var{method}, @var{options})
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
23 ##
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
24 ## Generate a regular mesh from irregular data using interpolation.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
25 ## The function is defined by @code{@var{y} = f (@var{x})}.
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
26 ## The interpolation points are all @var{xi}.
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
27 ##
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
28 ## The interpolation method can be @code{"nearest"} or @code{"linear"}.
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
29 ## If method is omitted it defaults to @code{"linear"}.
14370
fbdee844550c griddatan.m: Restore ability to pass options to underlying Qhull
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
30 ##
fbdee844550c griddatan.m: Restore ability to pass options to underlying Qhull
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
31 ## The optional argument @var{options} is passed directly to Qhull when
fbdee844550c griddatan.m: Restore ability to pass options to underlying Qhull
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
32 ## computing the Delaunay triangulation used for interpolation. See
fbdee844550c griddatan.m: Restore ability to pass options to underlying Qhull
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
33 ## @code{delaunayn} for information on the defaults and how to pass different
fbdee844550c griddatan.m: Restore ability to pass options to underlying Qhull
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
34 ## values.
fbdee844550c griddatan.m: Restore ability to pass options to underlying Qhull
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
35 ## @seealso{griddata, griddata3, delaunayn}
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
36 ## @end deftypefn
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
37
7017
a1dbe9d80eee [project @ 2007-10-12 21:27:11 by jwe]
jwe
parents: 7016
diff changeset
38 ## Author: David Bateman <dbateman@free.fr>
a1dbe9d80eee [project @ 2007-10-12 21:27:11 by jwe]
jwe
parents: 7016
diff changeset
39
14370
fbdee844550c griddatan.m: Restore ability to pass options to underlying Qhull
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
40 function yi = griddatan (x, y, xi, method = "linear", varargin)
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
41
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
42 if (nargin < 3)
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
43 print_usage ();
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
44 endif
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
45
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
46 if (ischar (method))
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
47 method = tolower (method);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
48 endif
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
49
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
50 [m, n] = size (x);
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
51 [mi, ni] = size (xi);
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
52
14370
fbdee844550c griddatan.m: Restore ability to pass options to underlying Qhull
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
53 if (n != ni || rows (y) != m || columns (y) != 1)
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
54 error ("griddatan: dimensional mismatch");
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
55 endif
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
56
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
57 ## triangulate data
14370
fbdee844550c griddatan.m: Restore ability to pass options to underlying Qhull
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
58 tri = delaunayn (x, varargin{:});
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
59
10548
479536c5bb10 Replace lowercase nan with NaN for visual cue in scripts
Rik <code@nomad.inbox5.com>
parents: 8920
diff changeset
60 yi = NaN (mi, 1);
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
61
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
62 if (strcmp (method, "nearest"))
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
63 ## search index of nearest point
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
64 idx = dsearchn (x, tri, xi);
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
65 valid = !isnan (idx);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
66 yi(valid) = y(idx(valid));
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
67
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
68 elseif (strcmp (method, "linear"))
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
69 ## search for every point the enclosing triangle
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
70 [tri_list, bary_list] = tsearchn (x, tri, xi);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
71
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
72 ## only keep the points within triangles.
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
73 valid = !isnan (tri_list);
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
74 tri_list = tri_list(!isnan (tri_list));
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
75 bary_list = bary_list(!isnan (tri_list), :);
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
76 nr_t = rows (tri_list);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
77
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
78 ## assign x,y for each point of simplex
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
79 xt = reshape (x(tri(tri_list,:),:), [nr_t, n+1, n]);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
80 yt = y(tri(tri_list,:));
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
81
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
82 ## Use barycentric coordinate of point to calculate yi
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
83 yi(valid) = sum (y(tri(tri_list,:)) .* bary_list, 2);
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
84
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
85 else
11472
1740012184f9 Use uppercase for variable names in error() strings to match Info documentation. Only m-files done.
Rik <octave@nomad.inbox5.com>
parents: 10548
diff changeset
86 error ("griddatan: unknown interpolation METHOD");
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
87 endif
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
88
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
89 endfunction
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
90
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
91
8153
ec0a13863eb7 Only run tests that depend on HDF5 and QHull if Octave was actually
Soren Hauberg <hauberg@gmail.com>
parents: 7017
diff changeset
92 %!testif HAVE_QHULL
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
93 %! [xx,yy] = meshgrid (linspace (-1,1,32));
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
94 %! xi = [xx(:), yy(:)];
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
95 %! x = 2*rand (100,2) - 1;
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
96 %! x = [x;1,1;1,-1;-1,-1;-1,1];
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
97 %! y = sin (2 * sum (x.^2,2));
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
98 %! zz = griddatan (x,y,xi,"linear");
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
99 %! zz2 = griddata (x(:,1),x(:,2),y,xi(:,1),xi(:,2),"linear");
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
100 %! assert (zz, zz2, 1e-10);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
101
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
102 %!testif HAVE_QHULL
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
103 %! [xx,yy] = meshgrid (linspace (-1,1,32));
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
104 %! xi = [xx(:), yy(:)];
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
105 %! x = 2*rand (100,2) - 1;
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
106 %! x = [x;1,1;1,-1;-1,-1;-1,1];
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
107 %! y = sin (2*sum (x.^2,2));
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
108 %! zz = griddatan (x,y,xi,"nearest");
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
109 %! zz2 = griddata (x(:,1),x(:,2),y,xi(:,1),xi(:,2),"nearest");
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
110 %! assert (zz, zz2, 1e-10);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
111