annotate scripts/geometry/griddatan.m @ 28216:f5644ccd1df5

griddata.m, griddatan.m: Small tweak in input validation. * griddata.m: Only call tolower on METHOD input when it is required. * griddatan.m: Only call tolower on METHOD input when it is required. Add note to documentation that query values outside the convex hull will return NaN.
author Rik <rik@octave.org>
date Mon, 13 Apr 2020 18:19:41 -0700
parents bb50407f9c60
children 90fea9cc9caa
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
27923
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
1 ########################################################################
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
2 ##
27919
1891570abac8 update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents: 27918
diff changeset
3 ## Copyright (C) 2007-2020 The Octave Project Developers
27918
b442ec6dda5c use centralized file for copyright info for individual contributors
John W. Eaton <jwe@octave.org>
parents: 27845
diff changeset
4 ##
27923
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
5 ## See the file COPYRIGHT.md in the top-level directory of this
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
6 ## distribution or <https://octave.org/copyright/>.
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
7 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
8 ## This file is part of Octave.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
9 ##
24534
194eb4bd202b maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents: 23220
diff changeset
10 ## Octave is free software: you can redistribute it and/or modify it
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
11 ## under the terms of the GNU General Public License as published by
24534
194eb4bd202b maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents: 23220
diff changeset
12 ## the Free Software Foundation, either version 3 of the License, or
22755
3a2b891d0b33 maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents: 22323
diff changeset
13 ## (at your option) any later version.
6823
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 ## 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
16 ## WITHOUT ANY WARRANTY; without even the implied warranty of
22755
3a2b891d0b33 maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents: 22323
diff changeset
17 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
3a2b891d0b33 maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents: 22323
diff changeset
18 ## GNU General Public License for more details.
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
19 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
20 ## 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
21 ## along with Octave; see the file COPYING. If not, see
24534
194eb4bd202b maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents: 23220
diff changeset
22 ## <https://www.gnu.org/licenses/>.
27923
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
23 ##
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
24 ########################################################################
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
25
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
26 ## -*- texinfo -*-
20852
516bb87ea72e 2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents: 20158
diff changeset
27 ## @deftypefn {} {@var{yi} =} griddatan (@var{x}, @var{y}, @var{xi})
516bb87ea72e 2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents: 20158
diff changeset
28 ## @deftypefnx {} {@var{yi} =} griddatan (@var{x}, @var{y}, @var{xi}, @var{method})
516bb87ea72e 2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents: 20158
diff changeset
29 ## @deftypefnx {} {@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
30 ##
28211
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
31 ## Interpolate irregular source data @var{x}, @var{y} at points specified by
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
32 ## @var{xi}.
20158
7503499a252b doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
33 ##
28211
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
34 ## The input @var{x} is an MxN matrix representing M points in an N-dimensional
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
35 ## space. The input @var{y} is a single-valued column vector (Mx1)
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
36 ## representing a function evaluated at the points @var{x}, i.e.,
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
37 ## @code{@var{y} = fcn (@var{x})}. The input @var{xi} is a list of points
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
38 ## for which the function output @var{yi} should be approximated through
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
39 ## interpolation. @var{xi} must have the same number of columns (@var{N})
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
40 ## as @var{x} so that the dimensionality matches.
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
41 ##
28211
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
42 ## The optional input interpolation @var{method} can be @qcode{"nearest"} or
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
43 ## @qcode{"linear"}. When the method is @qcode{"nearest"}, the output @var{yi}
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
44 ## will be the closest point in the original data @var{x} to the query point
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
45 ## @var{xi}. When the method is @qcode{"linear"}, the output @var{yi} will
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
46 ## be a linear interpolation between the two closest points in the original
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
47 ## source data. If @var{method} is omitted or empty, it defaults to
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
48 ## @qcode{"linear"}.
19593
446c46af4b42 strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents: 17744
diff changeset
49 ##
14370
fbdee844550c griddatan.m: Restore ability to pass options to underlying Qhull
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
50 ## 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
51 ## 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
52 ## @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
53 ## values.
28211
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
54 ##
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
55 ## Example
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
56 ##
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
57 ## @example
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
58 ## @group
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
59 ## ## Evaluate sombrero() function at irregular data points
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
60 ## x = 16*gallery ("uniformdata", [200,1], 1) - 8;
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
61 ## y = 16*gallery ("uniformdata", [200,1], 11) - 8;
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
62 ## z = sin (sqrt (x.^2 + y.^2)) ./ sqrt (x.^2 + y.^2);
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
63 ## ## Create a regular grid and interpolate data
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
64 ## [xi, yi] = ndgrid (linspace (-8, 8, 50));
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
65 ## zi = griddatan ([x, y], z, [xi(:), yi(:)]);
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
66 ## zi = reshape (zi, size (xi));
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
67 ## ## Plot results
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
68 ## clf ();
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
69 ## plot3 (x, y, z, "or");
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
70 ## hold on
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
71 ## surf (xi, yi, zi);
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
72 ## legend ("Original Data", "Interpolated Data");
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
73 ## @end group
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
74 ## @end example
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
75 ##
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
76 ## Programming Notes: If the input is complex the real and imaginary parts
28216
f5644ccd1df5 griddata.m, griddatan.m: Small tweak in input validation.
Rik <rik@octave.org>
parents: 28213
diff changeset
77 ## are interpolated separately. Interpolation is based on a Delaunay
f5644ccd1df5 griddata.m, griddatan.m: Small tweak in input validation.
Rik <rik@octave.org>
parents: 28213
diff changeset
78 ## triangulation and any query values outside the convex hull of the input
f5644ccd1df5 griddata.m, griddatan.m: Small tweak in input validation.
Rik <rik@octave.org>
parents: 28213
diff changeset
79 ## points will return @code{NaN}. For 2-D and 3-D data additional
f5644ccd1df5 griddata.m, griddatan.m: Small tweak in input validation.
Rik <rik@octave.org>
parents: 28213
diff changeset
80 ## interpolation methods are available by using the @code{griddata} function.
14370
fbdee844550c griddatan.m: Restore ability to pass options to underlying Qhull
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
81 ## @seealso{griddata, griddata3, delaunayn}
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
82 ## @end deftypefn
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
83
14370
fbdee844550c griddatan.m: Restore ability to pass options to underlying Qhull
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
84 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
85
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
86 if (nargin < 3)
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
87 print_usage ();
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
88 endif
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
89
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
90 [m, n] = size (x);
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
91 [mi, ni] = size (xi);
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
92
28211
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
93 if (m < n + 1)
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
94 error ("griddatan: number of points in X (rows of X) must be greater than dimensionality of data + 1 (columns of X + 1)");
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
95 endif
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
96 if (! iscolumn (y) || rows (y) != m)
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
97 error ("griddatan: Y must be a column vector with the same number of points (rows) as X");
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
98 endif
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
99 if (n != ni)
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
100 error ("griddatan: dimension of query data XI (columns) must match X");
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
101 endif
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
102
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
103 if (nargin > 3)
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
104 if (isempty (method))
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
105 method = "linear";
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
106 elseif (! ischar (method))
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
107 error ("griddatan: METHOD must be a string");
28216
f5644ccd1df5 griddata.m, griddatan.m: Small tweak in input validation.
Rik <rik@octave.org>
parents: 28213
diff changeset
108 else
f5644ccd1df5 griddata.m, griddatan.m: Small tweak in input validation.
Rik <rik@octave.org>
parents: 28213
diff changeset
109 method = tolower (method);
28211
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
110 endif
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
111
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
112 if (strcmp (method, "linear") || strcmp (method, "nearest"))
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
113 ## Do nothing, these are implemented methods
28213
bb50407f9c60 griddatan.m: Issue informational error if "cubic" METHOD is used.
Rik <rik@octave.org>
parents: 28211
diff changeset
114 elseif (strcmp (method, "v4"))
28211
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
115 error ('griddatan: "%s" METHOD is available for 2-D inputs by using "griddata"', method);
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
116
28213
bb50407f9c60 griddatan.m: Issue informational error if "cubic" METHOD is used.
Rik <rik@octave.org>
parents: 28211
diff changeset
117 elseif (any (strcmp (method, {"cubic", "natural"})))
bb50407f9c60 griddatan.m: Issue informational error if "cubic" METHOD is used.
Rik <rik@octave.org>
parents: 28211
diff changeset
118 ## FIXME: Remove when griddata.m supports these methods.
bb50407f9c60 griddatan.m: Issue informational error if "cubic" METHOD is used.
Rik <rik@octave.org>
parents: 28211
diff changeset
119 error ('griddatan: "%s" interpolation METHOD not yet implemented', method);
28211
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
120
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
121 else
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
122 error ('griddatan: unknown interpolation METHOD: "%s"', method);
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
123 endif
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
124
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
125 endif
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
126
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
127 ## triangulate data
14370
fbdee844550c griddatan.m: Restore ability to pass options to underlying Qhull
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
128 tri = delaunayn (x, varargin{:});
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
129
10548
479536c5bb10 Replace lowercase nan with NaN for visual cue in scripts
Rik <code@nomad.inbox5.com>
parents: 8920
diff changeset
130 yi = NaN (mi, 1);
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11523
diff changeset
131
28211
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
132 if (strcmp (method, "linear"))
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
133 ## search for every point the enclosing triangle
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
134 [tri_list, bary_list] = tsearchn (x, tri, xi);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
135
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
136 ## only keep the points within triangles.
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
137 valid = ! isnan (tri_list);
27244
3f354ef16400 griddatan.m: Allow a single query point as input (bug #45542, bug #56515).
Rik <rik@octave.org>
parents: 27239
diff changeset
138 tri_list = tri_list(valid);
3f354ef16400 griddatan.m: Allow a single query point as input (bug #45542, bug #56515).
Rik <rik@octave.org>
parents: 27239
diff changeset
139 bary_list = bary_list(valid, :);
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
140 nr_t = rows (tri_list);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
141
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
142 ## Use barycentric coordinate of point to calculate yi
27244
3f354ef16400 griddatan.m: Allow a single query point as input (bug #45542, bug #56515).
Rik <rik@octave.org>
parents: 27239
diff changeset
143 if (isscalar (tri_list))
3f354ef16400 griddatan.m: Allow a single query point as input (bug #45542, bug #56515).
Rik <rik@octave.org>
parents: 27239
diff changeset
144 ## Special case required by orientation rules for vector/vector index.
3f354ef16400 griddatan.m: Allow a single query point as input (bug #45542, bug #56515).
Rik <rik@octave.org>
parents: 27239
diff changeset
145 yi(valid) = sum (y(tri(tri_list,:)).' .* bary_list, 2);
3f354ef16400 griddatan.m: Allow a single query point as input (bug #45542, bug #56515).
Rik <rik@octave.org>
parents: 27239
diff changeset
146 else
3f354ef16400 griddatan.m: Allow a single query point as input (bug #45542, bug #56515).
Rik <rik@octave.org>
parents: 27239
diff changeset
147 yi(valid) = sum (y(tri(tri_list,:)) .* bary_list, 2);
27239
bee80e27dcc5 griddatan.m: Fix 4-D interpolation (bug #45542, bug #56515).
Juho Iipponen <juho.iipponen@helsinki.fi>
parents: 26376
diff changeset
148 endif
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
149
28211
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
150 else
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
151 ## search index of nearest point
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
152 idx = dsearchn (x, tri, xi);
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
153 valid = ! isnan (idx);
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
154 yi(valid) = y(idx(valid));
28210
bb929d5a34cb griddata.m: Added support for "v4" biharmonic spline interpolation method (bug #33539).
Nicholas R. Jankowski <jankowskin@asme.org>
parents: 28027
diff changeset
155
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
156 endif
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
157
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
158 endfunction
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
159
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
160
8153
ec0a13863eb7 Only run tests that depend on HDF5 and QHull if Octave was actually
Soren Hauberg <hauberg@gmail.com>
parents: 7017
diff changeset
161 %!testif HAVE_QHULL
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
162 %! [xx,yy] = meshgrid (linspace (-1,1,32));
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
163 %! xi = [xx(:), yy(:)];
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
164 %! x = 2*rand (100,2) - 1;
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
165 %! 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
166 %! y = sin (2 * sum (x.^2,2));
28211
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
167 %! zz = griddatan (x,y,xi, "linear");
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
168 %! zz2 = griddata (x(:,1),x(:,2),y,xi(:,1),xi(:,2), "linear");
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
169 %! 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
170
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
171 %!testif HAVE_QHULL
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
172 %! [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
173 %! xi = [xx(:), yy(:)];
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
174 %! 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
175 %! 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
176 %! y = sin (2*sum (x.^2,2));
28211
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
177 %! zz = griddatan (x,y,xi, "nearest");
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
178 %! zz2 = griddata (x(:,1),x(:,2),y,xi(:,1),xi(:,2), "nearest");
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
179 %! assert (zz, zz2, 1e-10);
27244
3f354ef16400 griddatan.m: Allow a single query point as input (bug #45542, bug #56515).
Rik <rik@octave.org>
parents: 27239
diff changeset
180
28011
b31a118729ed update bug status for tests
John W. Eaton <jwe@octave.org>
parents: 27978
diff changeset
181 %!testif HAVE_QHULL <*56515>
27244
3f354ef16400 griddatan.m: Allow a single query point as input (bug #45542, bug #56515).
Rik <rik@octave.org>
parents: 27239
diff changeset
182 %! x = [ 0, 0; 1, 1; 0, 1; 1, 0 ];
3f354ef16400 griddatan.m: Allow a single query point as input (bug #45542, bug #56515).
Rik <rik@octave.org>
parents: 27239
diff changeset
183 %! y = [ 1; 2; 3; 4 ];
3f354ef16400 griddatan.m: Allow a single query point as input (bug #45542, bug #56515).
Rik <rik@octave.org>
parents: 27239
diff changeset
184 %! xi = [ .5, .5 ];
3f354ef16400 griddatan.m: Allow a single query point as input (bug #45542, bug #56515).
Rik <rik@octave.org>
parents: 27239
diff changeset
185 %! yi = griddatan (x, y, xi);
28210
bb929d5a34cb griddata.m: Added support for "v4" biharmonic spline interpolation method (bug #33539).
Nicholas R. Jankowski <jankowskin@asme.org>
parents: 28027
diff changeset
186
bb929d5a34cb griddata.m: Added support for "v4" biharmonic spline interpolation method (bug #33539).
Nicholas R. Jankowski <jankowskin@asme.org>
parents: 28027
diff changeset
187 ## Test input validation
bb929d5a34cb griddata.m: Added support for "v4" biharmonic spline interpolation method (bug #33539).
Nicholas R. Jankowski <jankowskin@asme.org>
parents: 28027
diff changeset
188 %!error griddatan ()
bb929d5a34cb griddata.m: Added support for "v4" biharmonic spline interpolation method (bug #33539).
Nicholas R. Jankowski <jankowskin@asme.org>
parents: 28027
diff changeset
189 %!error griddatan (1)
bb929d5a34cb griddata.m: Added support for "v4" biharmonic spline interpolation method (bug #33539).
Nicholas R. Jankowski <jankowskin@asme.org>
parents: 28027
diff changeset
190 %!error griddatan (1,2)
28211
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
191 %!error <number of points in X> griddatan (1,2,3)
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
192 %!error <Y must be a column vector> griddatan ([1;2],[3,4], 1)
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
193 %!error <Y must .* same number of points .* as X> griddatan ([1;2],[3;4;5], 1)
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
194 %!error <dimension of .* XI .* must match X> griddatan ([1;2],[3;4], [1, 2])
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
195 %!error <METHOD must be a string> griddatan ([1;2],[3;4], 1, 5)
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
196 %!error <"v4" METHOD is available for 2-D> griddatan ([1;2],[3;4], 1, "v4")
28213
bb50407f9c60 griddatan.m: Issue informational error if "cubic" METHOD is used.
Rik <rik@octave.org>
parents: 28211
diff changeset
197 %!error <"cubic" .* not yet implemented> griddatan ([1;2],[3;4], 1, "cubic")
28211
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
198 %!error <"natural" .* not yet implemented> griddatan ([1;2],[3;4], 1, "natural")
289882040316 griddatan.m: Rewrite docstring and redo input validation.
Rik <rik@octave.org>
parents: 28210
diff changeset
199 %!error <unknown .* METHOD: "foobar"> griddatan ([1;2],[3;4], 1, "foobar")