annotate scripts/general/interp1.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 7277fe922e99
children 7ce925166af6
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
14138
72c96de7a403 maint: update copyright notices for 2012
John W. Eaton <jwe@octave.org>
parents: 13141
diff changeset
1 ## Copyright (C) 2000-2012 Paul Kienzle
9929
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
2 ## Copyright (C) 2009 VZLU Prague
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
3 ##
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
4 ## This file is part of Octave.
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
5 ##
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
6 ## Octave is free software; you can redistribute it and/or modify it
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
7 ## under the terms of the GNU General Public License as published by
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 7001
diff changeset
8 ## the Free Software Foundation; either version 3 of the License, or (at
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 7001
diff changeset
9 ## your option) any later version.
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
10 ##
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
11 ## Octave is distributed in the hope that it will be useful, but
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
12 ## WITHOUT ANY WARRANTY; without even the implied warranty of
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
13 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
14 ## General Public License for more details.
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
15 ##
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
16 ## You should have received a copy of the GNU General Public License
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 7001
diff changeset
17 ## along with Octave; see the file COPYING. If not, see
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 7001
diff changeset
18 ## <http://www.gnu.org/licenses/>.
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
19
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
20 ## -*- texinfo -*-
10846
a4f482e66b65 Grammarcheck more of the documentation.
Rik <octave@nomad.inbox5.com>
parents: 10821
diff changeset
21 ## @deftypefn {Function File} {@var{yi} =} interp1 (@var{x}, @var{y}, @var{xi})
10820
c44c786f87ba interp1.m: When absent set X equal to the inices of Y.
Ben Abbott <bpabbott@mac.com>
parents: 10793
diff changeset
22 ## @deftypefnx {Function File} {@var{yi} =} interp1 (@var{y}, @var{xi})
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
23 ## @deftypefnx {Function File} {@var{yi} =} interp1 (@dots{}, @var{method})
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
24 ## @deftypefnx {Function File} {@var{yi} =} interp1 (@dots{}, @var{extrap})
14359
7277fe922e99 doc: Use Octave preference for double quote in docstrings in scripts/
Rik <octave@nomad.inbox5.com>
parents: 14335
diff changeset
25 ## @deftypefnx {Function File} {@var{pp} =} interp1 (@dots{}, "pp")
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
26 ##
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
27 ## One-dimensional interpolation. Interpolate @var{y}, defined at the
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11536
diff changeset
28 ## points @var{x}, at the points @var{xi}. The sample points @var{x}
10846
a4f482e66b65 Grammarcheck more of the documentation.
Rik <octave@nomad.inbox5.com>
parents: 10821
diff changeset
29 ## must be monotonic. If not specified, @var{x} is taken to be the
a4f482e66b65 Grammarcheck more of the documentation.
Rik <octave@nomad.inbox5.com>
parents: 10821
diff changeset
30 ## indices of @var{y}. If @var{y} is an array, treat the columns
7001
8b0cfeb06365 [project @ 2007-10-10 18:02:59 by jwe]
jwe
parents: 6742
diff changeset
31 ## of @var{y} separately.
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
32 ##
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
33 ## Method is one of:
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
34 ##
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
35 ## @table @asis
14359
7277fe922e99 doc: Use Octave preference for double quote in docstrings in scripts/
Rik <octave@nomad.inbox5.com>
parents: 14335
diff changeset
36 ## @item "nearest"
9070
e9dc2ed2ec0f Cleanup documentation for poly.texi, interp.texi, geometry.texi
Rik <rdrider0-list@yahoo.com>
parents: 9051
diff changeset
37 ## Return the nearest neighbor.
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10820
diff changeset
38 ##
14359
7277fe922e99 doc: Use Octave preference for double quote in docstrings in scripts/
Rik <octave@nomad.inbox5.com>
parents: 14335
diff changeset
39 ## @item "linear"
9070
e9dc2ed2ec0f Cleanup documentation for poly.texi, interp.texi, geometry.texi
Rik <rdrider0-list@yahoo.com>
parents: 9051
diff changeset
40 ## Linear interpolation from nearest neighbors
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10820
diff changeset
41 ##
14359
7277fe922e99 doc: Use Octave preference for double quote in docstrings in scripts/
Rik <octave@nomad.inbox5.com>
parents: 14335
diff changeset
42 ## @item "pchip"
11536
702dbd0c53f5 Add undocumented ppder, ppint, ppjumps functions to documentation.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
43 ## Piecewise cubic Hermite interpolating polynomial
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10820
diff changeset
44 ##
14359
7277fe922e99 doc: Use Octave preference for double quote in docstrings in scripts/
Rik <octave@nomad.inbox5.com>
parents: 14335
diff changeset
45 ## @item "cubic"
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
46 ## Cubic interpolation (same as @code{pchip})
10821
693e22af08ae Grammarcheck documentation of m-files
Rik <octave@nomad.inbox5.com>
parents: 10820
diff changeset
47 ##
14359
7277fe922e99 doc: Use Octave preference for double quote in docstrings in scripts/
Rik <octave@nomad.inbox5.com>
parents: 14335
diff changeset
48 ## @item "spline"
12175
2090995ca588 Correct en-dash,em-dash instances in docstrings.
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
49 ## Cubic spline interpolation---smooth first and second derivatives
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
50 ## throughout the curve
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
51 ## @end table
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
52 ##
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
53 ## Appending '*' to the start of the above method forces @code{interp1}
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
54 ## to assume that @var{x} is uniformly spaced, and only @code{@var{x}
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 8920
diff changeset
55 ## (1)} and @code{@var{x} (2)} are referenced. This is usually faster,
14359
7277fe922e99 doc: Use Octave preference for double quote in docstrings in scripts/
Rik <octave@nomad.inbox5.com>
parents: 14335
diff changeset
56 ## and is never slower. The default method is "linear".
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
57 ##
14359
7277fe922e99 doc: Use Octave preference for double quote in docstrings in scripts/
Rik <octave@nomad.inbox5.com>
parents: 14335
diff changeset
58 ## If @var{extrap} is the string "extrap", then extrapolate values beyond
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
59 ## the endpoints. If @var{extrap} is a number, replace values beyond the
6742
ebf96cc00ee9 [project @ 2007-06-18 16:27:52 by jwe]
jwe
parents: 6721
diff changeset
60 ## endpoints with that number. If @var{extrap} is missing, assume NA.
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
61 ##
14359
7277fe922e99 doc: Use Octave preference for double quote in docstrings in scripts/
Rik <octave@nomad.inbox5.com>
parents: 14335
diff changeset
62 ## If the string argument "pp" is specified, then @var{xi} should not be
11536
702dbd0c53f5 Add undocumented ppder, ppint, ppjumps functions to documentation.
Rik <octave@nomad.inbox5.com>
parents: 11523
diff changeset
63 ## supplied and @code{interp1} returns the piecewise polynomial that
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
64 ## can later be used with @code{ppval} to evaluate the interpolation.
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
65 ## There is an equivalence, such that @code{ppval (interp1 (@var{x},
14359
7277fe922e99 doc: Use Octave preference for double quote in docstrings in scripts/
Rik <octave@nomad.inbox5.com>
parents: 14335
diff changeset
66 ## @var{y}, @var{method}, "pp"), @var{xi}) == interp1 (@var{x}, @var{y},
7277fe922e99 doc: Use Octave preference for double quote in docstrings in scripts/
Rik <octave@nomad.inbox5.com>
parents: 14335
diff changeset
67 ## @var{xi}, @var{method}, "extrap")}.
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
68 ##
10711
fbd7843974fa Periodic grammar check of documentation files to ensure common format.
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
69 ## Duplicate points in @var{x} specify a discontinuous interpolant. There
9929
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
70 ## should be at most 2 consecutive points with the same value.
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
71 ## The discontinuous interpolant is right-continuous if @var{x} is increasing,
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
72 ## left-continuous if it is decreasing.
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
73 ## Discontinuities are (currently) only allowed for "nearest" and "linear"
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
74 ## methods; in all other cases, @var{x} must be strictly monotonic.
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
75 ##
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
76 ## An example of the use of @code{interp1} is
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
77 ##
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
78 ## @example
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
79 ## @group
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7671
diff changeset
80 ## xf = [0:0.05:10];
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7671
diff changeset
81 ## yf = sin (2*pi*xf/5);
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7671
diff changeset
82 ## xp = [0:10];
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7671
diff changeset
83 ## yp = sin (2*pi*xp/5);
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7671
diff changeset
84 ## lin = interp1 (xp, yp, xf);
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7671
diff changeset
85 ## spl = interp1 (xp, yp, xf, "spline");
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7671
diff changeset
86 ## cub = interp1 (xp, yp, xf, "cubic");
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7671
diff changeset
87 ## near = interp1 (xp, yp, xf, "nearest");
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7671
diff changeset
88 ## plot (xf, yf, "r", xf, lin, "g", xf, spl, "b",
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7671
diff changeset
89 ## xf, cub, "c", xf, near, "m", xp, yp, "r*");
14327
4d917a6a858b doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
90 ## legend ("original", "linear", "spline", "cubic", "nearest");
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
91 ## @end group
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
92 ## @end example
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
93 ##
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
94 ## @seealso{interpft}
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
95 ## @end deftypefn
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
96
5838
376e02b2ce70 [project @ 2006-06-01 20:23:53 by jwe]
jwe
parents: 5837
diff changeset
97 ## Author: Paul Kienzle
376e02b2ce70 [project @ 2006-06-01 20:23:53 by jwe]
jwe
parents: 5837
diff changeset
98 ## Date: 2000-03-25
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
99 ## added 'nearest' as suggested by Kai Habel
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
100 ## 2000-07-17 Paul Kienzle
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
101 ## added '*' methods and matrix y
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
102 ## check for proper table lengths
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
103 ## 2002-01-23 Paul Kienzle
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
104 ## fixed extrapolation
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
105
5838
376e02b2ce70 [project @ 2006-06-01 20:23:53 by jwe]
jwe
parents: 5837
diff changeset
106 function yi = interp1 (x, y, varargin)
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
107
10820
c44c786f87ba interp1.m: When absent set X equal to the inices of Y.
Ben Abbott <bpabbott@mac.com>
parents: 10793
diff changeset
108 if (nargin < 2 || nargin > 6)
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
109 print_usage ();
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
110 endif
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
111
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
112 method = "linear";
6742
ebf96cc00ee9 [project @ 2007-06-18 16:27:52 by jwe]
jwe
parents: 6721
diff changeset
113 extrap = NA;
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
114 xi = [];
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
115 ispp = false;
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
116 firstnumeric = true;
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
117
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
118 if (nargin > 2)
5838
376e02b2ce70 [project @ 2006-06-01 20:23:53 by jwe]
jwe
parents: 5837
diff changeset
119 for i = 1:length (varargin)
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
120 arg = varargin{i};
5838
376e02b2ce70 [project @ 2006-06-01 20:23:53 by jwe]
jwe
parents: 5837
diff changeset
121 if (ischar (arg))
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 9929
diff changeset
122 arg = tolower (arg);
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 9929
diff changeset
123 if (strcmp ("extrap", arg))
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 9929
diff changeset
124 extrap = "extrap";
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 9929
diff changeset
125 elseif (strcmp ("pp", arg))
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
126 ispp = true;
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 9929
diff changeset
127 else
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 9929
diff changeset
128 method = arg;
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 9929
diff changeset
129 endif
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
130 else
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 9929
diff changeset
131 if (firstnumeric)
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 9929
diff changeset
132 xi = arg;
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 9929
diff changeset
133 firstnumeric = false;
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 9929
diff changeset
134 else
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 9929
diff changeset
135 extrap = arg;
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 9929
diff changeset
136 endif
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
137 endif
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
138 endfor
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
139 endif
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
140
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
141 if (isempty (xi) && firstnumeric && ! ispp)
10820
c44c786f87ba interp1.m: When absent set X equal to the inices of Y.
Ben Abbott <bpabbott@mac.com>
parents: 10793
diff changeset
142 xi = y;
c44c786f87ba interp1.m: When absent set X equal to the inices of Y.
Ben Abbott <bpabbott@mac.com>
parents: 10793
diff changeset
143 y = x;
c44c786f87ba interp1.m: When absent set X equal to the inices of Y.
Ben Abbott <bpabbott@mac.com>
parents: 10793
diff changeset
144 x = 1:numel(y);
c44c786f87ba interp1.m: When absent set X equal to the inices of Y.
Ben Abbott <bpabbott@mac.com>
parents: 10793
diff changeset
145 endif
c44c786f87ba interp1.m: When absent set X equal to the inices of Y.
Ben Abbott <bpabbott@mac.com>
parents: 10793
diff changeset
146
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
147 ## reshape matrices for convenience
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
148 x = x(:);
9754
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
149 nx = rows (x);
9769
9a1c4fe44af8 small interp1 simplification
Jaroslav Hajek <highegg@gmail.com>
parents: 9754
diff changeset
150 szx = size (xi);
9754
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
151 if (isvector (y))
5838
376e02b2ce70 [project @ 2006-06-01 20:23:53 by jwe]
jwe
parents: 5837
diff changeset
152 y = y(:);
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
153 endif
13141
e81ddf9cacd5 maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 12608
diff changeset
154
5838
376e02b2ce70 [project @ 2006-06-01 20:23:53 by jwe]
jwe
parents: 5837
diff changeset
155 szy = size (y);
9754
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
156 y = y(:,:);
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
157 [ny, nc] = size (y);
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
158 xi = xi(:);
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
159
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
160 ## determine sizes
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
161 if (nx < 2 || ny < 2)
5838
376e02b2ce70 [project @ 2006-06-01 20:23:53 by jwe]
jwe
parents: 5837
diff changeset
162 error ("interp1: table too short");
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
163 endif
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
164
9754
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
165 ## check whether x is sorted; sort if not.
9929
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
166 if (! issorted (x, "either"))
9754
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
167 [x, p] = sort (x);
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
168 y = y(p,:);
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
169 endif
5838
376e02b2ce70 [project @ 2006-06-01 20:23:53 by jwe]
jwe
parents: 5837
diff changeset
170
9754
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
171 starmethod = method(1) == "*";
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
172
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
173 if (starmethod)
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
174 dx = x(2) - x(1);
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
175 else
9929
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
176 jumps = x(1:nx-1) == x(2:nx);
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
177 have_jumps = any (jumps);
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
178 if (have_jumps)
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
179 if (any (strcmp (method, {"nearest", "linear"})))
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
180 if (any (jumps(1:nx-2) & jumps(2:nx-1)))
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
181 warning ("interp1: extra points in discontinuities");
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
182 endif
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
183 else
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
184 error ("interp1: discontinuities not supported for method %s", method);
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
185 endif
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
186 endif
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
187 endif
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
188
9754
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
189 ## Proceed with interpolating by all methods.
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
190
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
191 switch (method)
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
192 case "nearest"
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
193 pp = mkpp ([x(1); (x(1:nx-1)+x(2:nx))/2; x(nx)], shiftdim (y, 1), szy(2:end));
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
194 pp.orient = "first";
13141
e81ddf9cacd5 maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 12608
diff changeset
195
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
196 if (ispp)
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
197 yi = pp;
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
198 else
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
199 yi = ppval (pp, reshape (xi, szx));
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
200 endif
9754
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
201 case "*nearest"
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
202 pp = mkpp ([x(1), x(1)+[0.5:(nx-1)]*dx, x(nx)], shiftdim (y, 1), szy(2:end));
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
203 pp.orient = "first";
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
204 if (ispp)
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
205 yi = pp;
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
206 else
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
207 yi = ppval(pp, reshape (xi, szx));
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
208 endif
9754
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
209 case "linear"
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
210 dy = diff (y);
13141
e81ddf9cacd5 maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 12608
diff changeset
211 dx = diff (x);
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
212 dx = repmat (dx, [1 size(dy)(2:end)]);
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
213 coefs = [(dy./dx).'(:), y(1:nx-1, :).'(:)];
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
214 xx = x;
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
215
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
216 if (have_jumps)
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
217 ## Omit zero-size intervals.
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
218 coefs(jumps, :) = [];
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
219 xx(jumps) = [];
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
220 endif
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
221
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
222 pp = mkpp (xx, coefs, szy(2:end));
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
223 pp.orient = "first";
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
224
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
225 if (ispp)
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
226 yi = pp;
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
227 else
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
228 yi = ppval(pp, reshape (xi, szx));
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
229 endif
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
230
9754
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
231 case "*linear"
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
232 dy = diff (y);
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
233 coefs = [(dy/dx).'(:), y(1:nx-1, :).'(:)];
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
234 pp = mkpp (x, coefs, szy(2:end));
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
235 pp.orient = "first";
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
236
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
237 if (ispp)
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
238 yi = pp;
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
239 else
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
240 yi = ppval(pp, reshape (xi, szx));
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
241 endif
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
242
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
243 case {"pchip", "*pchip", "cubic", "*cubic"}
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11536
diff changeset
244 if (nx == 2 || starmethod)
6374
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
245 x = linspace (x(1), x(nx), ny);
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
246 endif
13141
e81ddf9cacd5 maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 12608
diff changeset
247
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
248 if (ispp)
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
249 y = shiftdim (reshape (y, szy), 1);
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
250 yi = pchip (x, y);
9754
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
251 else
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
252 y = shiftdim (y, 1);
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
253 yi = pchip (x, y, reshape (xi, szx));
9754
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
254 endif
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
255 case {"spline", "*spline"}
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11536
diff changeset
256 if (nx == 2 || starmethod)
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11536
diff changeset
257 x = linspace(x(1), x(nx), ny);
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
258 endif
13141
e81ddf9cacd5 maint: untabify and remove trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 12608
diff changeset
259
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
260 if (ispp)
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
261 y = shiftdim (reshape (y, szy), 1);
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
262 yi = spline (x, y);
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
263 else
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
264 y = shiftdim (y, 1);
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
265 yi = spline (x, y, reshape (xi, szx));
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
266 endif
9754
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
267 otherwise
5838
376e02b2ce70 [project @ 2006-06-01 20:23:53 by jwe]
jwe
parents: 5837
diff changeset
268 error ("interp1: invalid method '%s'", method);
9754
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
269 endswitch
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
270
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
271 if (! ispp)
9754
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
272 if (! ischar (extrap))
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
273 ## determine which values are out of range and set them to extrap,
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
274 ## unless extrap == "extrap".
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
275 minx = min (x(1), x(nx));
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
276 maxx = max (x(1), x(nx));
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
277
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
278 outliers = xi < minx | ! (xi <= maxx); # this catches even NaNs
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
279 if (size_equal (outliers, yi))
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
280 yi(outliers) = extrap;
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
281 yi = reshape (yi, szx);
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
282 elseif (!isvector (yi))
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
283 if (strcmp (method, "pchip") || strcmp (method, "*pchip")
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
284 ||strcmp (method, "cubic") || strcmp (method, "*cubic")
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
285 ||strcmp (method, "spline") || strcmp (method, "*spline"))
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
286 yi(:, outliers) = extrap;
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
287 yi = shiftdim(yi, 1);
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
288 else
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
289 yi(outliers, :) = extrap;
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
290 endif
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
291 else
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
292 yi(outliers.') = extrap;
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
293 endif
9754
4219e5cf773d improve interp1 and pchip
Jaroslav Hajek <highegg@gmail.com>
parents: 9070
diff changeset
294 endif
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
295 else
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
296 yi.orient = "first";
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
297 endif
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
298
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
299 endfunction
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
300
14237
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
301
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
302 %!demo
14237
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
303 %! clf;
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
304 %! xf = 0:0.05:10; yf = sin (2*pi*xf/5);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
305 %! xp = 0:10; yp = sin (2*pi*xp/5);
14237
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
306 %! lin = interp1 (xp,yp,xf, "linear");
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
307 %! spl = interp1 (xp,yp,xf, "spline");
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
308 %! cub = interp1 (xp,yp,xf, "pchip");
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
309 %! near= interp1 (xp,yp,xf, "nearest");
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
310 %! plot (xf,yf,"r",xf,near,"g",xf,lin,"b",xf,cub,"c",xf,spl,"m",xp,yp,"r*");
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
311 %! legend ("original", "nearest", "linear", "pchip", "spline");
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
312 %! %--------------------------------------------------------
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
313 %! % confirm that interpolated function matches the original
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
314
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
315 %!demo
14237
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
316 %! clf;
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
317 %! xf = 0:0.05:10; yf = sin (2*pi*xf/5);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
318 %! xp = 0:10; yp = sin (2*pi*xp/5);
14237
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
319 %! lin = interp1 (xp,yp,xf, "*linear");
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
320 %! spl = interp1 (xp,yp,xf, "*spline");
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
321 %! cub = interp1 (xp,yp,xf, "*cubic");
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
322 %! near= interp1 (xp,yp,xf, "*nearest");
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
323 %! plot (xf,yf,"r",xf,near,"g",xf,lin,"b",xf,cub,"c",xf,spl,"m",xp,yp,"r*");
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
324 %! legend ("*original", "*nearest", "*linear", "*cubic", "*spline");
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
325 %! %--------------------------------------------------------
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
326 %! % confirm that interpolated function matches the original
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
327
6721
01036667884a [project @ 2007-06-14 06:56:41 by dbateman]
dbateman
parents: 6702
diff changeset
328 %!demo
14237
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
329 %! clf;
6721
01036667884a [project @ 2007-06-14 06:56:41 by dbateman]
dbateman
parents: 6702
diff changeset
330 %! t = 0 : 0.3 : pi; dt = t(2)-t(1);
01036667884a [project @ 2007-06-14 06:56:41 by dbateman]
dbateman
parents: 6702
diff changeset
331 %! n = length (t); k = 100; dti = dt*n/k;
01036667884a [project @ 2007-06-14 06:56:41 by dbateman]
dbateman
parents: 6702
diff changeset
332 %! ti = t(1) + [0 : k-1]*dti;
01036667884a [project @ 2007-06-14 06:56:41 by dbateman]
dbateman
parents: 6702
diff changeset
333 %! y = sin (4*t + 0.3) .* cos (3*t - 0.1);
14237
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
334 %! ddyc = diff (diff (interp1 (t,y,ti, "cubic")) ./dti)./dti;
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
335 %! ddys = diff (diff (interp1 (t,y,ti, "spline"))./dti)./dti;
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
336 %! ddyp = diff (diff (interp1 (t,y,ti, "pchip")) ./dti)./dti;
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
337 %! plot (ti(2:end-1),ddyc,'g+', ti(2:end-1),ddys,'b*', ti(2:end-1),ddyp,'c^');
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
338 %! legend ("cubic", "spline", "pchip");
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
339 %! title ("Second derivative of interpolated 'sin (4*t + 0.3) .* cos (3*t - 0.1)'");
6721
01036667884a [project @ 2007-06-14 06:56:41 by dbateman]
dbateman
parents: 6702
diff changeset
340
9929
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
341 %!demo
14237
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
342 %! clf;
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
343 %! xf = 0:0.05:10; yf = sin (2*pi*xf/5) - (xf >= 5);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
344 %! xp = [0:.5:4.5,4.99,5:.5:10]; yp = sin (2*pi*xp/5) - (xp >= 5);
14237
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
345 %! lin = interp1 (xp,yp,xf, "linear");
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
346 %! near= interp1 (xp,yp,xf, "nearest");
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
347 %! plot (xf,yf,"r", xf,near,"g", xf,lin,"b", xp,yp,"r*");
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
348 %! legend ("original", "nearest", "linear");
9929
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
349 %! %--------------------------------------------------------
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
350 %! % confirm that interpolated function matches the original
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
351
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
352 ##FIXME: add test for n-d arguments here
9929
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
353
6374
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
354 ## For each type of interpolated test, confirm that the interpolated
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
355 ## value at the knots match the values at the knots. Points away
14359
7277fe922e99 doc: Use Octave preference for double quote in docstrings in scripts/
Rik <octave@nomad.inbox5.com>
parents: 14335
diff changeset
356 ## from the knots are requested, but only "nearest" and "linear"
6374
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
357 ## confirm they are the correct values.
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
358
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
359 %!shared xp, yp, xi, style
14237
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
360 %! xp = 0:2:10;
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
361 %! yp = sin (2*pi*xp/5);
6374
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
362 %! xi = [-1, 0, 2.2, 4, 6.6, 10, 11];
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
363
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
364 ## The following BLOCK/ENDBLOCK section is repeated for each style
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
365 ## nearest, linear, cubic, spline, pchip
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
366 ## The test for ppval of cubic has looser tolerance, but otherwise
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
367 ## the tests are identical.
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
368 ## Note that the block checks style and *style; if you add more tests
14237
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
369 ## be sure to add them to both sections of each block. One test,
6374
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
370 ## style vs. *style, occurs only in the first section.
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
371 ## There is an ENDBLOCKTEST after the final block
14237
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
372
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
373 %!test style = "nearest";
6374
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
374 ## BLOCK
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
375 %!assert (interp1 (xp, yp, [min(xp)-1, max(xp)+1],style), [NA, NA])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
376 %!assert (interp1 (xp,yp,xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
377 %!assert (interp1 (xp,yp,xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
378 %!assert (interp1 (xp',yp',xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
379 %!assert (interp1 (xp',yp',xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
380 %!assert (isempty (interp1 (xp',yp',[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
381 %!assert (isempty (interp1 (xp,yp,[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
382 %!assert (interp1 (xp,[yp',yp'],xi(:),style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
383 %! [interp1(xp,yp,xi(:),style),interp1(xp,yp,xi(:),style)])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
384 %!assert (interp1 (xp,yp,xi,style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
385 %! interp1 (fliplr (xp),fliplr (yp),xi,style),100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
386 %!assert (ppval (interp1 (xp,yp,style,"pp"),xi),
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
387 %! interp1 (xp,yp,xi,style,"extrap"),10*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
388 %!error interp1 (1,1,1, style)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
389 %!assert (interp1 (xp,[yp',yp'],xi,style),
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
390 %! interp1 (xp,[yp',yp'],xi,["*",style]),100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
391 %!test style = ["*",style];
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
392 %!assert (interp1 (xp, yp, [min(xp)-1, max(xp)+1],style), [NA, NA])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
393 %!assert (interp1 (xp,yp,xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
394 %!assert (interp1 (xp,yp,xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
395 %!assert (interp1 (xp',yp',xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
396 %!assert (interp1 (xp',yp',xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
397 %!assert (isempty (interp1 (xp',yp',[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
398 %!assert (isempty (interp1 (xp,yp,[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
399 %!assert (interp1 (xp,[yp',yp'],xi(:),style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
400 %! [interp1(xp,yp,xi(:),style),interp1(xp,yp,xi(:),style)])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
401 %!assert (interp1 (xp,yp,xi,style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
402 %! interp1 (fliplr (xp),fliplr (yp),xi,style),100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
403 %!assert (ppval (interp1 (xp,yp,style,"pp"),xi),
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
404 %! interp1 (xp,yp,xi,style,"extrap"),10*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
405 %!error interp1 (1,1,1, style)
6374
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
406 ## ENDBLOCK
14237
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
407
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
408 %!test style = "linear";
6374
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
409 ## BLOCK
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
410 %!assert (interp1 (xp, yp, [min(xp)-1, max(xp)+1],style), [NA, NA])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
411 %!assert (interp1 (xp,yp,xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
412 %!assert (interp1 (xp,yp,xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
413 %!assert (interp1 (xp',yp',xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
414 %!assert (interp1 (xp',yp',xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
415 %!assert (isempty (interp1 (xp',yp',[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
416 %!assert (isempty (interp1 (xp,yp,[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
417 %!assert (interp1 (xp,[yp',yp'],xi(:),style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
418 %! [interp1(xp,yp,xi(:),style),interp1(xp,yp,xi(:),style)])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
419 %!assert (interp1 (xp,yp,xi,style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
420 %! interp1 (fliplr (xp),fliplr (yp),xi,style),100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
421 %!assert (ppval (interp1 (xp,yp,style,"pp"),xi),
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
422 %! interp1 (xp,yp,xi,style,"extrap"),10*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
423 %!error interp1 (1,1,1, style)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
424 %!assert (interp1 (xp,[yp',yp'],xi,style),
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
425 %! interp1 (xp,[yp',yp'],xi,["*",style]),100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
426 %!test style = ['*',style];
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
427 %!assert (interp1 (xp, yp, [min(xp)-1, max(xp)+1],style), [NA, NA])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
428 %!assert (interp1 (xp,yp,xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
429 %!assert (interp1 (xp,yp,xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
430 %!assert (interp1 (xp',yp',xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
431 %!assert (interp1 (xp',yp',xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
432 %!assert (isempty (interp1 (xp',yp',[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
433 %!assert (isempty (interp1 (xp,yp,[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
434 %!assert (interp1 (xp,[yp',yp'],xi(:),style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
435 %! [interp1(xp,yp,xi(:),style),interp1(xp,yp,xi(:),style)])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
436 %!assert (interp1 (xp,yp,xi,style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
437 %! interp1 (fliplr (xp),fliplr (yp),xi,style),100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
438 %!assert (ppval (interp1 (xp,yp,style,"pp"),xi),
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
439 %! interp1 (xp,yp,xi,style,"extrap"),10*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
440 %!error interp1 (1,1,1, style)
6374
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
441 ## ENDBLOCK
14237
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
442
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
443 %!test style = "cubic";
6374
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
444 ## BLOCK
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
445 %!assert (interp1 (xp, yp, [min(xp)-1, max(xp)+1],style), [NA, NA])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
446 %!assert (interp1 (xp,yp,xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
447 %!assert (interp1 (xp,yp,xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
448 %!assert (interp1 (xp',yp',xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
449 %!assert (interp1 (xp',yp',xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
450 %!assert (isempty (interp1 (xp',yp',[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
451 %!assert (isempty (interp1 (xp,yp,[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
452 %!assert (interp1 (xp,[yp',yp'],xi(:),style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
453 %! [interp1(xp,yp,xi(:),style),interp1(xp,yp,xi(:),style)])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
454 %!assert (interp1 (xp,yp,xi,style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
455 %! interp1 (fliplr (xp),fliplr (yp),xi,style),100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
456 %!assert (ppval (interp1 (xp,yp,style,"pp"),xi),
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
457 %! interp1 (xp,yp,xi,style,"extrap"),100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
458 %!error interp1 (1,1,1, style)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
459 %!assert (interp1 (xp,[yp',yp'],xi,style),
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
460 %! interp1 (xp,[yp',yp'],xi,["*",style]),100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
461 %!test style = ["*",style];
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
462 %!assert (interp1 (xp, yp, [min(xp)-1, max(xp)+1],style), [NA, NA])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
463 %!assert (interp1 (xp,yp,xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
464 %!assert (interp1 (xp,yp,xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
465 %!assert (interp1 (xp',yp',xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
466 %!assert (interp1 (xp',yp',xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
467 %!assert (isempty (interp1 (xp',yp',[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
468 %!assert (isempty (interp1 (xp,yp,[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
469 %!assert (interp1 (xp,[yp',yp'],xi(:),style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
470 %! [interp1(xp,yp,xi(:),style),interp1(xp,yp,xi(:),style)])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
471 %!assert (interp1 (xp,yp,xi,style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
472 %! interp1 (fliplr (xp),fliplr (yp),xi,style),100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
473 %!assert (ppval (interp1 (xp,yp,style,"pp"),xi),
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
474 %! interp1 (xp,yp,xi,style,"extrap"),100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
475 %!error interp1 (1,1,1, style)
6374
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
476 ## ENDBLOCK
14237
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
477
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
478 %!test style = "pchip";
6374
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
479 ## BLOCK
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
480 %!assert (interp1 (xp, yp, [min(xp)-1, max(xp)+1],style), [NA, NA])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
481 %!assert (interp1 (xp,yp,xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
482 %!assert (interp1 (xp,yp,xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
483 %!assert (interp1 (xp',yp',xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
484 %!assert (interp1 (xp',yp',xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
485 %!assert (isempty (interp1 (xp',yp',[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
486 %!assert (isempty (interp1 (xp,yp,[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
487 %!assert (interp1 (xp,[yp',yp'],xi(:),style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
488 %! [interp1(xp,yp,xi(:),style),interp1(xp,yp,xi(:),style)])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
489 %!assert (interp1 (xp,yp,xi,style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
490 %! interp1 (fliplr (xp),fliplr (yp),xi,style),100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
491 %!assert (ppval (interp1 (xp,yp,style,"pp"),xi),
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
492 %! interp1 (xp,yp,xi,style,"extrap"),10*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
493 %!error interp1 (1,1,1, style)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
494 %!assert (interp1 (xp,[yp',yp'],xi,style),
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
495 %! interp1 (xp,[yp',yp'],xi,["*",style]),100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
496 %!test style = ["*",style];
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
497 %!assert (interp1 (xp, yp, [min(xp)-1, max(xp)+1],style), [NA, NA])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
498 %!assert (interp1 (xp,yp,xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
499 %!assert (interp1 (xp,yp,xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
500 %!assert (interp1 (xp',yp',xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
501 %!assert (interp1 (xp',yp',xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
502 %!assert (isempty (interp1 (xp',yp',[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
503 %!assert (isempty (interp1 (xp,yp,[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
504 %!assert (interp1 (xp,[yp',yp'],xi(:),style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
505 %! [interp1(xp,yp,xi(:),style),interp1(xp,yp,xi(:),style)])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
506 %!assert (interp1 (xp,yp,xi,style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
507 %! interp1 (fliplr (xp),fliplr (yp),xi,style),100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
508 %!assert (ppval (interp1 (xp,yp,style,"pp"),xi),
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
509 %! interp1 (xp,yp,xi,style,"extrap"),10*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
510 %!error interp1 (1,1,1, style)
6374
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
511 ## ENDBLOCK
14237
11949c9795a0 Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
512
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
513 %!test style = "spline";
6374
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
514 ## BLOCK
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
515 %!assert (interp1 (xp, yp, [min(xp)-1, max(xp)+1],style), [NA, NA])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
516 %!assert (interp1 (xp,yp,xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
517 %!assert (interp1 (xp,yp,xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
518 %!assert (interp1 (xp',yp',xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
519 %!assert (interp1 (xp',yp',xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
520 %!assert (isempty (interp1 (xp',yp',[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
521 %!assert (isempty (interp1 (xp,yp,[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
522 %!assert (interp1 (xp,[yp',yp'],xi(:),style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
523 %! [interp1(xp,yp,xi(:),style),interp1(xp,yp,xi(:),style)])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
524 %!assert (interp1 (xp,yp,xi,style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
525 %! interp1 (fliplr (xp),fliplr (yp),xi,style),100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
526 %!assert (ppval (interp1 (xp,yp,style,"pp"),xi),
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
527 %! interp1 (xp,yp,xi,style,"extrap"),10*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
528 %!error interp1 (1,1,1, style)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
529 %!assert (interp1 (xp,[yp',yp'],xi,style),
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
530 %! interp1 (xp,[yp',yp'],xi,["*",style]),100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
531 %!test style = ["*",style];
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
532 %!assert (interp1 (xp, yp, [min(xp)-1, max(xp)+1],style), [NA, NA])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
533 %!assert (interp1 (xp,yp,xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
534 %!assert (interp1 (xp,yp,xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
535 %!assert (interp1 (xp',yp',xp',style), yp', 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
536 %!assert (interp1 (xp',yp',xp,style), yp, 100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
537 %!assert (isempty (interp1 (xp',yp',[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
538 %!assert (isempty (interp1 (xp,yp,[],style)))
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
539 %!assert (interp1 (xp,[yp',yp'],xi(:),style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
540 %! [interp1(xp,yp,xi(:),style),interp1(xp,yp,xi(:),style)])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
541 %!assert (interp1 (xp,yp,xi,style),...
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
542 %! interp1 (fliplr (xp),fliplr (yp),xi,style),100*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
543 %!assert (ppval (interp1 (xp,yp,style,"pp"),xi),
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
544 %! interp1 (xp,yp,xi,style,"extrap"),10*eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
545 %!error interp1 (1,1,1, style)
6374
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
546 ## ENDBLOCK
419017274c1e [project @ 2007-03-01 15:57:50 by jwe]
jwe
parents: 6366
diff changeset
547 ## ENDBLOCKTEST
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
548
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
549 %!# test linear extrapolation
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
550 %!assert (interp1 ([1:5],[3:2:11],[0,6],"linear","extrap"), [1, 13], eps)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
551 %!assert (interp1 (xp, yp, [-1, max(xp)+1],"linear",5), [5, 5])
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
552
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
553 %!assert (interp1 (1:2,1:2,1.4,"nearest"), 1)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
554 %!assert (interp1 (1:2,1:2,1.4,"linear"), 1.4)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
555 %!assert (interp1 (1:4,1:4,1.4,"cubic"), 1.4)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
556 %!assert (interp1 (1:2,1:2,1.1, "spline"), 1.1)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
557 %!assert (interp1 (1:3,1:3,1.4,"spline"), 1.4)
5837
55404f3b0da1 [project @ 2006-06-01 19:05:31 by jwe]
jwe
parents:
diff changeset
558
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
559 %!assert (interp1 (1:2:4,1:2:4,1.4,"*nearest"), 1)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
560 %!assert (interp1 (1:2:4,1:2:4,[0,1,1.4,3,4],"*linear"), [NA,1,1.4,3,NA])
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
561 %!assert (interp1 (1:2:8,1:2:8,1.4,"*cubic"), 1.4)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
562 %!assert (interp1 (1:2,1:2,1.3, "*spline"), 1.3)
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
563 %!assert (interp1 (1:2:6,1:2:6,1.4,"*spline"), 1.4)
7671
4fbaba9abec1 implement compiled binary lookup
Jaroslav Hajek <highegg@gmail.com>
parents: 7017
diff changeset
564
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
565 %!assert (interp1 ([3,2,1],[3,2,2],2.5), 2.5)
9929
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
566
45c08d7c2c79 allow discontinuous interpolant in interp1
Jaroslav Hajek <highegg@gmail.com>
parents: 9769
diff changeset
567 %!assert (interp1 ([1,2,2,3,4],[0,1,4,2,1],[-1,1.5,2,2.5,3.5], "linear", "extrap"), [-2,0.5,4,3,1.5])
12608
59e2460acae1 make piecewise polynomial (pp) functions more compatible
Kai Habel <kai.habel@gmx.de>
parents: 12459
diff changeset
568 %!assert (interp1 ([4,4,3,2,0],[0,1,4,2,1],[1.5,4,4.5], "linear"), [1.75,1,NA])
10820
c44c786f87ba interp1.m: When absent set X equal to the inices of Y.
Ben Abbott <bpabbott@mac.com>
parents: 10793
diff changeset
569 %!assert (interp1 (0:4, 2.5), 1.5)
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
570
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
571 %!error interp1 ()
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
572 %!error interp1 (1,1,1, "linear")
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
573 %!error interp1 (1,1,1, "*nearest")
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
574 %!error interp1 (1,1,1, "*linear")
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
575 %!error interp1 (1:2,1:2,1, "bogus")
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14359
diff changeset
576