Mercurial > hg > octave-nkf
annotate scripts/geometry/griddata3.m @ 15107:03381a36f70d
generate convenience libraries for new parse-tree and interpfcn subdirectories
* src/Makefile.am (liboctinterp_la_SOURCES): Include octave.cc in the
list, not $(DIST_SRC).
(liboctinterp_la_LIBADD): Include octave-value/liboctave-value.la,
parse-tree/libparse-tree.la, interp-core/libinterp-core.la,
interpfcn/libinterpfcn.la, and corefcn/libcorefcn.la in the list.
* src/interp-core/module.mk (noinst_LTLIBRARIES): Add
interp-core/libinterp-core.la to the list.
(interp_core_libinterp_core_la_SOURCES): New variable.
* src/interpfcn/module.mk (noinst_LTLIBRARIES): Add
interpfcn/libinterpfcn.la to the list.
(interpfcn_libinterpfcn_la_SOURCES): New variable.
* src/parse-tree/module.mk (noinst_LTLIBRARIES): Add
parse-tree/libparse-tree.la to the list.
(parse_tree_libparse_tree_la_SOURCES): New variable.
* src/octave-value/module.mk (noinst_LTLIBRARIES): Add
octave-value/liboctave-value.la to the list.
(octave_value_liboctave_value_la_SOURCES): New variable.
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Sun, 05 Aug 2012 09:04:30 -0400 |
parents | 4e8f1d1b0d75 |
children | bc924baa2c4e |
rev | line source |
---|---|
14138
72c96de7a403
maint: update copyright notices for 2012
John W. Eaton <jwe@octave.org>
parents:
13747
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 -*- | |
14384
4e8f1d1b0d75
maint: periodic merge of stable to default.
Rik <octave@nomad.inbox5.com>
diff
changeset
|
20 ## @deftypefn {Function File} {@var{vi} =} griddata3 (@var{x}, @var{y}, @var{z}, @var{v}, @var{xi}, @var{yi}, @var{zi}) |
4e8f1d1b0d75
maint: periodic merge of stable to default.
Rik <octave@nomad.inbox5.com>
diff
changeset
|
21 ## @deftypefnx {Function File} {@var{vi} =} griddata3 (@var{x}, @var{y}, @var{z}, @var{v}, @var{xi}, @var{yi}, @var{zi}, @var{method}) |
4e8f1d1b0d75
maint: periodic merge of stable to default.
Rik <octave@nomad.inbox5.com>
diff
changeset
|
22 ## @deftypefnx {Function File} {@var{vi} =} griddata3 (@var{x}, @var{y}, @var{z}, @var{v}, @var{xi}, @var{yi}, @var{zi}, @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. |
14362
cb4f1915db92
fix docstring in griddata3
Carlo de Falco <kingcrimson@tiscali.it>
parents:
14138
diff
changeset
|
25 ## The function is defined by @code{@var{v} = f (@var{x}, @var{y}, @var{z})}. |
cb4f1915db92
fix docstring in griddata3
Carlo de Falco <kingcrimson@tiscali.it>
parents:
14138
diff
changeset
|
26 ## The interpolation points are specified by @var{xi}, @var{yi}, @var{zi}. |
6823 | 27 ## |
6826 | 28 ## The interpolation method can be @code{"nearest"} or @code{"linear"}. |
29 ## If method is omitted it defaults to @code{"linear"}. | |
14371
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
30 ## |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
31 ## The optional argument @var{options} is passed directly to Qhull when |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
32 ## computing the Delaunay triangulation used for interpolation. See |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
33 ## @code{delaunayn} for information on the defaults and how to pass different |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
34 ## values. |
14362
cb4f1915db92
fix docstring in griddata3
Carlo de Falco <kingcrimson@tiscali.it>
parents:
14138
diff
changeset
|
35 ## @seealso{griddata, griddatan, delaunayn} |
6823 | 36 ## @end deftypefn |
37 | |
7017 | 38 ## Author: David Bateman <dbateman@free.fr> |
39 | |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
40 function vi = griddata3 (x, y, z, v, xi, yi, zi, method, varargin) |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
41 |
6826 | 42 if (nargin < 7) |
43 print_usage (); | |
6823 | 44 endif |
45 | |
14371
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
46 if (isvector (x) && isvector (y) && isvector (z) && isvector (v)) |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
47 if (! isequal (length (x), length (y), length (z), length (v))) |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
48 error ("griddata: X, Y, Z, and V must be vectors of the same length"); |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
49 endif |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
50 elseif (! size_equal (x, y, z, v)) |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
51 error ("griddata: X, Y, Z, and V must have equal sizes"); |
6823 | 52 endif |
53 | |
14371
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
54 ## meshgrid xi, yi and zi if they are vectors unless |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
55 ## they are vectors of the same length. |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
56 if (isvector (xi) && isvector (yi) && isvector (zi)) |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
57 if (! isequal (length (xi), length (yi), length (zi))) |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
58 [xi, yi, zi] = meshgrid (xi, yi, zi); |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
59 else |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
60 ## Otherwise, convert to column vectors |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
61 xi = xi(:); |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
62 yi = yi(:); |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
63 zi = zi(:); |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
64 endif |
6823 | 65 endif |
66 | |
14371
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
67 if (! size_equal (xi, yi, zi)) |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
68 error ("griddata3: XI, YI, and ZI must be vectors or matrices of the same size"); |
6823 | 69 endif |
70 | |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
71 vi = griddatan ([x(:), y(:), z(:)], v(:), [xi(:), yi(:), zi(:)], varargin{:}); |
6826 | 72 vi = reshape (vi, size (xi)); |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
73 |
6823 | 74 endfunction |
75 | |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
76 |
8153
ec0a13863eb7
Only run tests that depend on HDF5 and QHull if Octave was actually
Soren Hauberg <hauberg@gmail.com>
parents:
7585
diff
changeset
|
77 %!testif HAVE_QHULL |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
78 %! old_state = rand ("state"); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
79 %! restore_state = onCleanup (@() rand ("state", old_state)); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
80 %! rand ("state", 0); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
81 %! x = 2 * rand (1000, 1) - 1; |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
82 %! y = 2 * rand (1000, 1) - 1; |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
83 %! z = 2 * rand (1000, 1) - 1; |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
84 %! v = x.^2 + y.^2 + z.^2; |
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
85 %! [xi, yi, zi] = meshgrid (-0.8:0.2:0.8); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
86 %! vi = griddata3 (x, y, z, v, xi, yi, zi, "linear"); |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
87 %! vv = vi - xi.^2 - yi.^2 - zi.^2; |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
88 %! assert (max (abs (vv(:))), 0, 0.1); |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
89 |
8153
ec0a13863eb7
Only run tests that depend on HDF5 and QHull if Octave was actually
Soren Hauberg <hauberg@gmail.com>
parents:
7585
diff
changeset
|
90 %!testif HAVE_QHULL |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
91 %! old_state = rand ("state"); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
92 %! restore_state = onCleanup (@() rand ("state", old_state)); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
93 %! rand ("state", 0); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
94 %! x = 2 * rand (1000, 1) - 1; |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
95 %! y = 2 * rand (1000, 1) - 1; |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
96 %! z = 2 * rand (1000, 1) - 1; |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
97 %! v = x.^2 + y.^2 + z.^2; |
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
98 %! [xi, yi, zi] = meshgrid (-0.8:0.2:0.8); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
99 %! vi = griddata3 (x, y, z, v, xi, yi, zi, "nearest"); |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
100 %! vv = vi - xi.^2 - yi.^2 - zi.^2; |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
101 %! assert (max (abs (vv(:))), 0, 0.1) |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
102 |