Mercurial > octave-antonio
annotate scripts/geometry/griddata3.m @ 14362:cb4f1915db92 stable
fix docstring in griddata3
* griddata3.m: some fixes in the docstring
author | Carlo de Falco <kingcrimson@tiscali.it> |
---|---|
date | Mon, 13 Feb 2012 21:34:47 +0100 |
parents | 72c96de7a403 |
children | 4e8f1d1b0d75 |
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 -*- | |
14362
cb4f1915db92
fix docstring in griddata3
Carlo de Falco <kingcrimson@tiscali.it>
parents:
14138
diff
changeset
|
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}) |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
21 ## |
6823 | 22 ## 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
|
23 ## 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
|
24 ## The interpolation points are specified by @var{xi}, @var{yi}, @var{zi}. |
6823 | 25 ## |
6826 | 26 ## The interpolation method can be @code{"nearest"} or @code{"linear"}. |
27 ## If method is omitted it defaults to @code{"linear"}. | |
14362
cb4f1915db92
fix docstring in griddata3
Carlo de Falco <kingcrimson@tiscali.it>
parents:
14138
diff
changeset
|
28 ## @seealso{griddata, griddatan, delaunayn} |
6823 | 29 ## @end deftypefn |
30 | |
7017 | 31 ## Author: David Bateman <dbateman@free.fr> |
32 | |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
33 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
|
34 |
6826 | 35 if (nargin < 7) |
36 print_usage (); | |
6823 | 37 endif |
38 | |
6826 | 39 if (!all (size (x) == size (y) & size (x) == size(z) & size(x) == size (v))) |
11472
1740012184f9
Use uppercase for variable names in error() strings to match Info documentation. Only m-files done.
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
40 error ("griddata3: X, Y, Z, and V must be vectors of same length"); |
6823 | 41 endif |
42 | |
43 ## meshgrid xi, yi and zi if they are vectors unless they | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
44 ## are vectors of the same length |
6826 | 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); | |
6823 | 48 endif |
49 | |
6826 | 50 if (any (size(xi) != size(yi)) || any (size(xi) != size(zi))) |
11472
1740012184f9
Use uppercase for variable names in error() strings to match Info documentation. Only m-files done.
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
51 error ("griddata3: XI, YI and ZI must be vectors or matrices of same size"); |
6823 | 52 endif |
53 | |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
54 vi = griddatan ([x(:), y(:), z(:)], v(:), [xi(:), yi(:), zi(:)], varargin{:}); |
6826 | 55 vi = reshape (vi, size (xi)); |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
56 |
6823 | 57 endfunction |
58 | |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
59 |
8153
ec0a13863eb7
Only run tests that depend on HDF5 and QHull if Octave was actually
Soren Hauberg <hauberg@gmail.com>
parents:
7585
diff
changeset
|
60 %!testif HAVE_QHULL |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
61 %! old_state = rand ("state"); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
62 %! 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
|
63 %! rand ("state", 0); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
64 %! x = 2 * rand (1000, 1) - 1; |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
65 %! y = 2 * rand (1000, 1) - 1; |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
66 %! z = 2 * rand (1000, 1) - 1; |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
67 %! v = x.^2 + y.^2 + z.^2; |
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
68 %! [xi, yi, zi] = meshgrid (-0.8:0.2:0.8); |
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
69 %! vi = griddata3 (x, y, z, v, xi, yi, zi, 'linear'); |
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
70 %! 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
|
71 %! assert (max (abs (vv(:))), 0, 0.1); |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
72 |
8153
ec0a13863eb7
Only run tests that depend on HDF5 and QHull if Octave was actually
Soren Hauberg <hauberg@gmail.com>
parents:
7585
diff
changeset
|
73 %!testif HAVE_QHULL |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
74 %! old_state = rand ("state"); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
75 %! 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
|
76 %! rand ("state", 0); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
77 %! x = 2 * rand (1000, 1) - 1; |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
78 %! y = 2 * rand (1000, 1) - 1; |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
79 %! z = 2 * rand (1000, 1) - 1; |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
80 %! v = x.^2 + y.^2 + z.^2; |
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
81 %! [xi, yi, zi] = meshgrid (-0.8:0.2:0.8); |
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
82 %! vi = griddata3 (x, y, z, v, xi, yi, zi, 'nearest'); |
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
83 %! 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
|
84 %! assert (max (abs (vv(:))), 0, 0.1) |