Mercurial > hg > octave-lyh
annotate scripts/geometry/griddatan.m @ 17531:f9abc8e5fc2e
-Werror at JIT related files
author | LYH <lyh.kernel@gmail.com> |
---|---|
date | Fri, 27 Sep 2013 17:21:33 +0800 |
parents | bc924baa2c4e |
children |
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 | 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 | |
7016 | 7 ## the Free Software Foundation; either version 3 of the License, or (at |
8 ## your option) any later version. | |
6823 | 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 | |
7016 | 16 ## along with Octave; see the file COPYING. If not, see |
17 ## <http://www.gnu.org/licenses/>. | |
6823 | 18 |
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 | 24 ## Generate a regular mesh from irregular data using interpolation. |
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 | 27 ## |
17289
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
14370
diff
changeset
|
28 ## The interpolation method can be @qcode{"nearest"} or @qcode{"linear"}. |
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
14370
diff
changeset
|
29 ## If method is omitted it defaults to @qcode{"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 | 36 ## @end deftypefn |
37 | |
7017 | 38 ## Author: David Bateman <dbateman@free.fr> |
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 | 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 | 43 print_usage (); |
6823 | 44 endif |
45 | |
6826 | 46 if (ischar (method)) |
47 method = tolower (method); | |
6823 | 48 endif |
49 | |
6826 | 50 [m, n] = size (x); |
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 | 54 error ("griddatan: dimensional mismatch"); |
55 endif | |
56 | |
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 | 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 | 62 if (strcmp (method, "nearest")) |
6823 | 63 ## search index of nearest point |
6826 | 64 idx = dsearchn (x, tri, xi); |
65 valid = !isnan (idx); | |
6823 | 66 yi(valid) = y(idx(valid)); |
67 | |
6826 | 68 elseif (strcmp (method, "linear")) |
6823 | 69 ## search for every point the enclosing triangle |
6826 | 70 [tri_list, bary_list] = tsearchn (x, tri, xi); |
6823 | 71 |
72 ## only keep the points within triangles. | |
6826 | 73 valid = !isnan (tri_list); |
74 tri_list = tri_list(!isnan (tri_list)); | |
75 bary_list = bary_list(!isnan (tri_list), :); | |
76 nr_t = rows (tri_list); | |
6823 | 77 |
78 ## assign x,y for each point of simplex | |
6826 | 79 xt = reshape (x(tri(tri_list,:),:), [nr_t, n+1, n]); |
6823 | 80 yt = y(tri(tri_list,:)); |
81 | |
82 ## Use barycentric coordinate of point to calculate yi | |
83 yi(valid) = sum (y(tri(tri_list,:)) .* bary_list, 2); | |
84 | |
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 | 87 endif |
88 | |
89 endfunction | |
90 | |
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 | 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 | 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 |