Mercurial > octave-antonio
annotate scripts/geometry/griddata3.m @ 20158:7503499a252b stable
doc: Update docstrings to have one sentence summary as first line.
Update scripts in audio, elfun, general, geometry, and image directories.
* scripts/audio/@audioplayer/__get_properties__.m,
scripts/audio/@audioplayer/audioplayer.m, scripts/audio/@audioplayer/get.m,
scripts/audio/@audioplayer/isplaying.m, scripts/audio/@audioplayer/play.m,
scripts/audio/@audioplayer/playblocking.m, scripts/audio/@audioplayer/set.m,
scripts/audio/@audioplayer/subsasgn.m, scripts/audio/@audioplayer/subsref.m,
scripts/audio/@audiorecorder/audiorecorder.m,
scripts/audio/@audiorecorder/get.m,
scripts/audio/@audiorecorder/getaudiodata.m,
scripts/audio/@audiorecorder/getplayer.m,
scripts/audio/@audiorecorder/isrecording.m,
scripts/audio/@audiorecorder/play.m, scripts/audio/@audiorecorder/record.m,
scripts/audio/@audiorecorder/recordblocking.m,
scripts/audio/@audiorecorder/set.m, scripts/audio/@audiorecorder/stop.m,
scripts/audio/@audiorecorder/subsasgn.m,
scripts/audio/@audiorecorder/subsref.m, scripts/audio/lin2mu.m,
scripts/audio/mu2lin.m, scripts/audio/record.m, scripts/audio/sound.m,
scripts/audio/soundsc.m, scripts/audio/wavread.m, scripts/audio/wavwrite.m,
scripts/elfun/cosd.m, scripts/elfun/sind.m, scripts/elfun/tand.m,
scripts/general/accumarray.m, scripts/general/accumdim.m,
scripts/general/bitcmp.m, scripts/general/bitget.m, scripts/general/bitset.m,
scripts/general/blkdiag.m, scripts/general/cart2pol.m,
scripts/general/cart2sph.m, scripts/general/cell2mat.m,
scripts/general/celldisp.m, scripts/general/chop.m,
scripts/general/circshift.m, scripts/general/common_size.m,
scripts/general/cplxpair.m, scripts/general/cumtrapz.m,
scripts/general/dblquad.m, scripts/general/deal.m, scripts/general/del2.m,
scripts/general/display.m, scripts/general/divergence.m,
scripts/general/fieldnames.m, scripts/general/flip.m,
scripts/general/flipdim.m, scripts/general/fliplr.m, scripts/general/flipud.m,
scripts/general/gradient.m, scripts/general/interp3.m,
scripts/general/interpft.m, scripts/general/interpn.m,
scripts/general/loadobj.m, scripts/general/logspace.m,
scripts/general/methods.m, scripts/general/nargchk.m,
scripts/general/narginchk.m, scripts/general/nargoutchk.m,
scripts/general/nextpow2.m, scripts/general/nthargout.m,
scripts/general/num2str.m, scripts/general/pol2cart.m,
scripts/general/polyarea.m, scripts/general/postpad.m,
scripts/general/prepad.m, scripts/general/profile.m, scripts/general/quadgk.m,
scripts/general/quadl.m, scripts/general/quadv.m, scripts/general/randi.m,
scripts/general/rat.m, scripts/general/repmat.m, scripts/general/rot90.m,
scripts/general/rotdim.m, scripts/general/saveobj.m, scripts/general/shift.m,
scripts/general/shiftdim.m, scripts/general/sortrows.m,
scripts/general/sph2cart.m, scripts/general/structfun.m,
scripts/general/subsindex.m, scripts/general/trapz.m,
scripts/general/triplequad.m, scripts/geometry/delaunayn.m,
scripts/geometry/dsearch.m, scripts/geometry/dsearchn.m,
scripts/geometry/griddata.m, scripts/geometry/griddata3.m,
scripts/geometry/griddatan.m, scripts/geometry/inpolygon.m,
scripts/geometry/rectint.m, scripts/geometry/tsearchn.m,
scripts/geometry/voronoi.m, scripts/geometry/voronoin.m,
scripts/help/__unimplemented__.m, scripts/help/doc.m,
scripts/help/doc_cache_create.m, scripts/help/get_first_help_sentence.m,
scripts/help/help.m, scripts/help/lookfor.m, scripts/help/print_usage.m,
scripts/help/type.m, scripts/help/which.m, scripts/image/autumn.m,
scripts/image/bone.m, scripts/image/brighten.m, scripts/image/cmpermute.m,
scripts/image/colorcube.m, scripts/image/contrast.m, scripts/image/cool.m,
scripts/image/copper.m, scripts/image/cubehelix.m, scripts/image/flag.m,
scripts/image/gmap40.m, scripts/image/gray.m, scripts/image/gray2ind.m,
scripts/image/hot.m, scripts/image/hsv.m, scripts/image/image.m,
scripts/image/imagesc.m, scripts/image/imfinfo.m, scripts/image/imformats.m,
scripts/image/imread.m, scripts/image/imshow.m, scripts/image/imwrite.m,
scripts/image/iscolormap.m, scripts/image/jet.m, scripts/image/lines.m,
scripts/image/ntsc2rgb.m, scripts/image/ocean.m, scripts/image/pink.m,
scripts/image/prism.m, scripts/image/rainbow.m, scripts/image/rgb2ntsc.m,
scripts/image/spinmap.m, scripts/image/spring.m, scripts/image/summer.m,
scripts/image/white.m, scripts/image/winter.m:
Update docstrings to have one sentence summary as first line.
Re-structure to have line lengths <= 80 chars.
author | Rik <rik@octave.org> |
---|---|
date | Sun, 03 May 2015 09:36:20 -0700 |
parents | 4197fc428c7d |
children |
rev | line source |
---|---|
19697
4197fc428c7d
maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents:
17744
diff
changeset
|
1 ## Copyright (C) 2007-2015 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. |
20158
7503499a252b
doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
25 ## |
14362
cb4f1915db92
fix docstring in griddata3
Carlo de Falco <kingcrimson@tiscali.it>
parents:
14138
diff
changeset
|
26 ## 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
|
27 ## The interpolation points are specified by @var{xi}, @var{yi}, @var{zi}. |
6823 | 28 ## |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
14384
diff
changeset
|
29 ## 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:
14384
diff
changeset
|
30 ## If method is omitted it defaults to @qcode{"linear"}. |
14371
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
31 ## |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
32 ## 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
|
33 ## 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
|
34 ## @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
|
35 ## values. |
14362
cb4f1915db92
fix docstring in griddata3
Carlo de Falco <kingcrimson@tiscali.it>
parents:
14138
diff
changeset
|
36 ## @seealso{griddata, griddatan, delaunayn} |
6823 | 37 ## @end deftypefn |
38 | |
7017 | 39 ## Author: David Bateman <dbateman@free.fr> |
40 | |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
41 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
|
42 |
6826 | 43 if (nargin < 7) |
44 print_usage (); | |
6823 | 45 endif |
46 | |
14371
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
47 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
|
48 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
|
49 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
|
50 endif |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
51 elseif (! size_equal (x, y, z, v)) |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
52 error ("griddata: X, Y, Z, and V must have equal sizes"); |
6823 | 53 endif |
54 | |
14371
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
55 ## 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
|
56 ## they are vectors of the same length. |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
57 if (isvector (xi) && isvector (yi) && isvector (zi)) |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
58 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
|
59 [xi, yi, zi] = meshgrid (xi, yi, zi); |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
60 else |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
61 ## Otherwise, convert to column vectors |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
62 xi = xi(:); |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
63 yi = yi(:); |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
64 zi = zi(:); |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
65 endif |
6823 | 66 endif |
67 | |
14371
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
68 if (! size_equal (xi, yi, zi)) |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
69 error ("griddata3: XI, YI, and ZI must be vectors or matrices of the same size"); |
6823 | 70 endif |
71 | |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
72 vi = griddatan ([x(:), y(:), z(:)], v(:), [xi(:), yi(:), zi(:)], varargin{:}); |
6826 | 73 vi = reshape (vi, size (xi)); |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
74 |
6823 | 75 endfunction |
76 | |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
77 |
8153
ec0a13863eb7
Only run tests that depend on HDF5 and QHull if Octave was actually
Soren Hauberg <hauberg@gmail.com>
parents:
7585
diff
changeset
|
78 %!testif HAVE_QHULL |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
79 %! old_state = rand ("state"); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
80 %! 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
|
81 %! rand ("state", 0); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
82 %! x = 2 * rand (1000, 1) - 1; |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
83 %! y = 2 * rand (1000, 1) - 1; |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
84 %! z = 2 * rand (1000, 1) - 1; |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
85 %! v = x.^2 + y.^2 + z.^2; |
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
86 %! [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
|
87 %! 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
|
88 %! 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
|
89 %! assert (max (abs (vv(:))), 0, 0.1); |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
90 |
8153
ec0a13863eb7
Only run tests that depend on HDF5 and QHull if Octave was actually
Soren Hauberg <hauberg@gmail.com>
parents:
7585
diff
changeset
|
91 %!testif HAVE_QHULL |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
92 %! old_state = rand ("state"); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
93 %! 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
|
94 %! rand ("state", 0); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
95 %! x = 2 * rand (1000, 1) - 1; |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
96 %! y = 2 * rand (1000, 1) - 1; |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
97 %! z = 2 * rand (1000, 1) - 1; |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
98 %! v = x.^2 + y.^2 + z.^2; |
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
99 %! [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
|
100 %! 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
|
101 %! 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
|
102 %! 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
|
103 |