Mercurial > octave
annotate scripts/geometry/griddata3.m @ 33625:d213a148b3f1 default tip @
ensure exp. terminal widget has focus at startup
* main-window.cc (main_window): call focus_command_window only if
event loop is idle by using a single shot timer
* main-window.h: make focus_command_window a public slot
author | Torsten Lilge <ttl-octave@mailbox.org> |
---|---|
date | Sun, 26 May 2024 02:29:44 +0200 |
parents | 2e484f9f1f18 |
children |
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 ## |
32632
2e484f9f1f18
maint: update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents:
31706
diff
changeset
|
3 ## Copyright (C) 2007-2024 The Octave Project Developers |
27918
b442ec6dda5c
use centralized file for copyright info for individual contributors
John W. Eaton <jwe@octave.org>
parents:
26376
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 | 7 ## |
8 ## This file is part of Octave. | |
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 | 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 | 14 ## |
15 ## Octave is distributed in the hope that it will be useful, but | |
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 | 19 ## |
20 ## You should have received a copy of the GNU General Public License | |
7016 | 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 | 25 |
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{vi} =} griddata3 (@var{x}, @var{y}, @var{z}, @var{v}, @var{xi}, @var{yi}, @var{zi}) |
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{vi} =} griddata3 (@var{x}, @var{y}, @var{z}, @var{v}, @var{xi}, @var{yi}, @var{zi}, @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{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
|
30 ## |
28215 | 31 ## Interpolate irregular 3-D source data at specified points. |
20158
7503499a252b
doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
32 ## |
28215 | 33 ## The inputs @var{x}, @var{y}, and @var{z} define the points where the |
34 ## function @code{@var{v} = f (@var{x}, @var{y}, @var{z})} is evaluated. The | |
35 ## inputs @var{x}, @var{y}, @var{z} are either vectors of the same length, or | |
36 ## if they are of unequal length, then they are expanded to a 3-D grid with | |
37 ## @code{meshgrid}. The size of the input @var{v} must match the size of the | |
38 ## original data, either as a vector or a matrix. | |
39 ## | |
14362
cb4f1915db92
fix docstring in griddata3
Carlo de Falco <kingcrimson@tiscali.it>
parents:
14138
diff
changeset
|
40 ## The interpolation points are specified by @var{xi}, @var{yi}, @var{zi}. |
6823 | 41 ## |
28215 | 42 ## The optional input interpolation @var{method} can be @qcode{"nearest"} or |
43 ## @qcode{"linear"}. When the method is @qcode{"nearest"}, the output @var{vi} | |
44 ## will be the closest point in the original data (@var{x}, @var{y}, @var{z}) | |
45 ## to the query point (@var{xi}, @var{yi}, @var{zi}). When the method is | |
46 ## @qcode{"linear"}, the output @var{vi} will be a linear interpolation between | |
47 ## the two closest points in the original source data in each dimension. | |
48 ## If @var{method} is omitted or empty, it defaults to @qcode{"linear"}. | |
14371
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
49 ## |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
50 ## 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
|
51 ## 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
|
52 ## @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
|
53 ## values. |
28215 | 54 ## |
55 ## Programming Notes: If the input is complex the real and imaginary parts | |
56 ## are interpolated separately. Interpolation is based on a Delaunay | |
57 ## triangulation and any query values outside the convex hull of the input | |
58 ## points will return @code{NaN}. | |
14362
cb4f1915db92
fix docstring in griddata3
Carlo de Falco <kingcrimson@tiscali.it>
parents:
14138
diff
changeset
|
59 ## @seealso{griddata, griddatan, delaunayn} |
6823 | 60 ## @end deftypefn |
61 | |
28085
739045a86cfd
griddata3.m: Fix ignored input "method" (bug #57835).
Nicholas R. Jankowski <jankowskin@asme.org>
parents:
27978
diff
changeset
|
62 function vi = griddata3 (x, y, z, v, xi, yi, zi, method = "linear", varargin) |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
63 |
6826 | 64 if (nargin < 7) |
65 print_usage (); | |
6823 | 66 endif |
67 | |
14371
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
68 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
|
69 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
|
70 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
|
71 endif |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
72 elseif (! size_equal (x, y, z, v)) |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
73 error ("griddata: X, Y, Z, and V must have equal sizes"); |
6823 | 74 endif |
75 | |
14371
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
76 ## 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
|
77 ## they are vectors of the same length. |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
78 if (isvector (xi) && isvector (yi) && isvector (zi)) |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
79 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
|
80 [xi, yi, zi] = meshgrid (xi, yi, zi); |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
81 else |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
82 ## Otherwise, convert to column vectors |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
83 xi = xi(:); |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
84 yi = yi(:); |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
85 zi = zi(:); |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
86 endif |
6823 | 87 endif |
88 | |
14371
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
89 if (! size_equal (xi, yi, zi)) |
079e6f3a0977
griddata3.m: Accept vectors of any orientation.
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
90 error ("griddata3: XI, YI, and ZI must be vectors or matrices of the same size"); |
6823 | 91 endif |
92 | |
30330
01de0045b2e3
maint: Shorten some long lines to <= 80 characters (bug #57599)
Rik <rik@octave.org>
parents:
29359
diff
changeset
|
93 vi = griddatan ([x(:), y(:), z(:)], v(:), [xi(:), yi(:), zi(:)], method, ... |
01de0045b2e3
maint: Shorten some long lines to <= 80 characters (bug #57599)
Rik <rik@octave.org>
parents:
29359
diff
changeset
|
94 varargin{:}); |
6826 | 95 vi = reshape (vi, size (xi)); |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
96 |
6823 | 97 endfunction |
98 | |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
99 |
8153
ec0a13863eb7
Only run tests that depend on HDF5 and QHull if Octave was actually
Soren Hauberg <hauberg@gmail.com>
parents:
7585
diff
changeset
|
100 %!testif HAVE_QHULL |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
101 %! old_state = rand ("state"); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
102 %! 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
|
103 %! rand ("state", 0); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
104 %! x = 2 * rand (1000, 1) - 1; |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
105 %! y = 2 * rand (1000, 1) - 1; |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
106 %! z = 2 * rand (1000, 1) - 1; |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
107 %! v = x.^2 + y.^2 + z.^2; |
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
108 %! [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
|
109 %! 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
|
110 %! 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
|
111 %! assert (max (abs (vv(:))), 0, 0.1); |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
112 |
8153
ec0a13863eb7
Only run tests that depend on HDF5 and QHull if Octave was actually
Soren Hauberg <hauberg@gmail.com>
parents:
7585
diff
changeset
|
113 %!testif HAVE_QHULL |
13747
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
114 %! old_state = rand ("state"); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
115 %! 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
|
116 %! rand ("state", 0); |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
117 %! x = 2 * rand (1000, 1) - 1; |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
118 %! y = 2 * rand (1000, 1) - 1; |
e8564e8b0043
Restore random number state after %!demos or %!tests
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
119 %! z = 2 * rand (1000, 1) - 1; |
7585
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
120 %! v = x.^2 + y.^2 + z.^2; |
522433b05f45
Fix griddata3 and add test code
David Bateman <dbateman@free.fr>
parents:
7418
diff
changeset
|
121 %! [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
|
122 %! 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
|
123 %! vv = vi - xi.^2 - yi.^2 - zi.^2; |
28085
739045a86cfd
griddata3.m: Fix ignored input "method" (bug #57835).
Nicholas R. Jankowski <jankowskin@asme.org>
parents:
27978
diff
changeset
|
124 %! assert (max (abs (vv(:))), 0.385, 0.1); |
28215 | 125 |
126 ## FIXME: Ideally, there should be BIST tests for input validation. | |
127 ## However, this function is deprecated in Matlab and it probably isn't worth | |
128 ## the coding time to work on a function that will be removed soon. |