Mercurial > hg > octave-nkf
annotate scripts/geometry/griddata3.m @ 14363:f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
* wavread.m, acosd.m, acot.m, acotd.m, acoth.m, acsc.m, acscd.m, acsch.m,
asec.m, asecd.m, asech.m, asind.m, atand.m, cosd.m, cot.m, cotd.m, coth.m,
csc.m, cscd.m, csch.m, sec.m, secd.m, sech.m, sind.m, tand.m, accumarray.m,
accumdim.m, bitcmp.m, bitget.m, bitset.m, blkdiag.m, cart2pol.m, cart2sph.m,
celldisp.m, chop.m, circshift.m, colon.m, common_size.m, cplxpair.m,
cumtrapz.m, curl.m, dblquad.m, deal.m, divergence.m, flipdim.m, fliplr.m,
flipud.m, genvarname.m, gradient.m, idivide.m, int2str.m, interp1.m,
interp1q.m, interp2.m, interp3.m, interpft.m, interpn.m, isa.m, isdir.m,
isequal.m, isequalwithequalnans.m, issquare.m, logspace.m, nargchk.m,
narginchk.m, nargoutchk.m, nextpow2.m, nthargout.m, num2str.m, pol2cart.m,
polyarea.m, postpad.m, prepad.m, profile.m, profshow.m, quadgk.m, quadv.m,
randi.m, rat.m, repmat.m, rot90.m, rotdim.m, shift.m, shiftdim.m, sph2cart.m,
structfun.m, trapz.m, triplequad.m, convhull.m, dsearch.m, dsearchn.m,
griddata3.m, griddatan.m, rectint.m, tsearchn.m, __makeinfo__.m, doc.m,
get_first_help_sentence.m, help.m, type.m, unimplemented.m, which.m, imread.m,
imwrite.m, dlmwrite.m, fileread.m, is_valid_file_id.m, strread.m, textread.m,
textscan.m, commutation_matrix.m, cond.m, condest.m, cross.m,
duplication_matrix.m, expm.m, housh.m, isdefinite.m, ishermitian.m,
issymmetric.m, logm.m, normest.m, null.m, onenormest.m, orth.m, planerot.m,
qzhess.m, rank.m, rref.m, trace.m, vech.m, ans.m, bincoeff.m, bug_report.m,
bzip2.m, comma.m, compare_versions.m, computer.m, edit.m, fileparts.m,
fullfile.m, getfield.m, gzip.m, info.m, inputname.m, isappdata.m, isdeployed.m,
ismac.m, ispc.m, isunix.m, list_primes.m, ls.m, mexext.m, namelengthmax.m,
news.m, orderfields.m, paren.m, recycle.m, rmappdata.m, semicolon.m,
setappdata.m, setfield.m, substruct.m, symvar.m, ver.m, version.m,
warning_ids.m, xor.m, fminbnd.m, fsolve.m, fzero.m, lsqnonneg.m, optimset.m,
pqpnonneg.m, sqp.m, matlabroot.m, __gnuplot_drawnow__.m,
__plt_get_axis_arg__.m, ancestor.m, cla.m, clf.m, close.m, colorbar.m,
colstyle.m, comet3.m, contourc.m, figure.m, gca.m, gcbf.m, gcbo.m, gcf.m,
ginput.m, graphics_toolkit.m, gtext.m, hggroup.m, hist.m, hold.m, isfigure.m,
ishghandle.m, ishold.m, isocolors.m, isonormals.m, isosurface.m, isprop.m,
legend.m, line.m, loglog.m, loglogerr.m, meshgrid.m, ndgrid.m, newplot.m,
orient.m, patch.m, plot3.m, plotyy.m, __print_parse_opts__.m, quiver3.m,
refreshdata.m, ribbon.m, semilogx.m, semilogxerr.m, semilogy.m, stem.m,
stem3.m, subplot.m, title.m, uigetfile.m, view.m, whitebg.m, compan.m, conv.m,
deconv.m, mkpp.m, mpoles.m, pchip.m, poly.m, polyaffine.m, polyder.m,
polyfit.m, polygcd.m, polyint.m, polyout.m, polyval.m, polyvalm.m, ppder.m,
ppint.m, ppjumps.m, ppval.m, residue.m, roots.m, spline.m, intersect.m,
ismember.m, powerset.m, setdiff.m, setxor.m, union.m, unique.m,
autoreg_matrix.m, bartlett.m, blackman.m, detrend.m, fftconv.m, fftfilt.m,
fftshift.m, freqz.m, hamming.m, hanning.m, ifftshift.m, sinc.m, sinetone.m,
sinewave.m, unwrap.m, bicg.m, bicgstab.m, gmres.m, gplot.m, nonzeros.m, pcg.m,
pcr.m, spaugment.m, spconvert.m, spdiags.m, speye.m, spfun.m, spones.m,
sprand.m, sprandsym.m, spstats.m, spy.m, svds.m, treelayout.m, bessel.m,
beta.m, betaln.m, factor.m, factorial.m, isprime.m, lcm.m, legendre.m,
nchoosek.m, nthroot.m, perms.m, pow2.m, primes.m, reallog.m, realpow.m,
realsqrt.m, hadamard.m, hankel.m, hilb.m, invhilb.m, magic.m, rosser.m,
vander.m, __finish__.m, center.m, cloglog.m, corr.m, cov.m, gls.m, histc.m,
iqr.m, kendall.m, kurtosis.m, logit.m, mahalanobis.m, mean.m, meansq.m,
median.m, mode.m, moment.m, ols.m, ppplot.m, prctile.m, probit.m, quantile.m,
range.m, ranks.m, run_count.m, runlength.m, skewness.m, spearman.m,
statistics.m, std.m, table.m, var.m, zscore.m, betacdf.m, betainv.m, betapdf.m,
betarnd.m, binocdf.m, binoinv.m, binopdf.m, binornd.m, cauchy_cdf.m,
cauchy_inv.m, cauchy_pdf.m, cauchy_rnd.m, chi2cdf.m, chi2inv.m, chi2pdf.m,
chi2rnd.m, discrete_cdf.m, discrete_inv.m, discrete_pdf.m, discrete_rnd.m,
empirical_cdf.m, empirical_inv.m, empirical_pdf.m, empirical_rnd.m, expcdf.m,
expinv.m, exppdf.m, exprnd.m, fcdf.m, finv.m, fpdf.m, frnd.m, gamcdf.m,
gaminv.m, gampdf.m, gamrnd.m, geocdf.m, geoinv.m, geopdf.m, geornd.m,
hygecdf.m, hygeinv.m, hygepdf.m, hygernd.m, kolmogorov_smirnov_cdf.m,
laplace_cdf.m, laplace_inv.m, laplace_pdf.m, laplace_rnd.m, logistic_cdf.m,
logistic_inv.m, logistic_pdf.m, logistic_rnd.m, logncdf.m, logninv.m,
lognpdf.m, lognrnd.m, nbincdf.m, nbininv.m, nbinpdf.m, nbinrnd.m, normcdf.m,
norminv.m, normpdf.m, normrnd.m, poisscdf.m, poissinv.m, poisspdf.m,
poissrnd.m, stdnormal_cdf.m, stdnormal_inv.m, stdnormal_pdf.m, stdnormal_rnd.m,
tcdf.m, tinv.m, tpdf.m, trnd.m, unidcdf.m, unidinv.m, unidpdf.m, unidrnd.m,
unifcdf.m, unifinv.m, unifpdf.m, unifrnd.m, wblcdf.m, wblinv.m, wblpdf.m,
wblrnd.m, kolmogorov_smirnov_test.m, kruskal_wallis_test.m, base2dec.m,
bin2dec.m, blanks.m, cstrcat.m, deblank.m, dec2base.m, dec2bin.m, dec2hex.m,
findstr.m, hex2dec.m, index.m, isletter.m, mat2str.m, rindex.m, str2num.m,
strcat.m, strjust.m, strmatch.m, strsplit.m, strtok.m, strtrim.m, strtrunc.m,
substr.m, validatestring.m, demo.m, example.m, fail.m, speed.m, addtodate.m,
asctime.m, clock.m, ctime.m, date.m, datenum.m, datetick.m, datevec.m,
eomday.m, etime.m, is_leap_year.m, now.m:
Use Octave coding conventions in all m-file %!test blocks
author | Rik <octave@nomad.inbox5.com> |
---|---|
date | Mon, 13 Feb 2012 07:29:44 -0800 |
parents | 72c96de7a403 |
children | 079e6f3a0977 |
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 -*- | |
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. |
23 ## The function is defined by @code{@var{y} = f (@var{x},@var{y},@var{z})}. | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
24 ## The interpolation points are all @var{xi}. |
6823 | 25 ## |
6826 | 26 ## The interpolation method can be @code{"nearest"} or @code{"linear"}. |
27 ## If method is omitted it defaults to @code{"linear"}. | |
6823 | 28 ## @seealso{griddata, delaunayn} |
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); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
69 %! 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
|
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); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
82 %! 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
|
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) |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
85 |