view scripts/plot/triplot.m @ 13929:9cae456085c2

Grammarcheck of documentation before 3.6.0 release. * accumarray.m, blkdiag.m, nargoutchk.m, nthargout.m, profexplore.m, profile.m, computer.m, orderfields.m, recycle.m, version.m, sqp.m, matlabroot.m, __plt_get_axis_arg__.m, isonormals.m, isosurface.m, __fltk_file_filter__.m, __is_function__.m, __uigetdir_fltk__.m, __uigetfile_fltk__.m, __uiobject_split_args__.m, __uiputfile_fltk__.m, uicontextmenu.m, uiresume.m, uiwait.m, mkpp.m, ppder.m, residue.m, addpref.m, getpref.m, ispref.m, loadprefs.m, prefsfile.m, saveprefs.m, rmpref.m, setpref.m, fftshift.m, bicg.m, bicgstab.m, cgs.m, gmres.m, __sprand_impl__.m, quantile.m, deblank.m, strsplit.m, addtodate.m, bsxfun.cc, kron.cc, regexp.cc, data.cc, file-io.cc, graphics.cc, load-save.cc, mappers.cc: Grammarcheck of documentation before 3.6.0 release.
author Rik <octave@nomad.inbox5.com>
date Wed, 23 Nov 2011 08:38:19 -0800
parents e8564e8b0043
children 5f0bb45e615c
line wrap: on
line source

## Copyright (C) 2007-2011 David Bateman
##
## This file is part of Octave.
##
## Octave is free software; you can redistribute it and/or modify it
## under the terms of the GNU General Public License as published by
## the Free Software Foundation; either version 3 of the License, or (at
## your option) any later version.
##
## Octave is distributed in the hope that it will be useful, but
## WITHOUT ANY WARRANTY; without even the implied warranty of
## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
## General Public License for more details.
##
## You should have received a copy of the GNU General Public License
## along with Octave; see the file COPYING.  If not, see
## <http://www.gnu.org/licenses/>.

## -*- texinfo -*-
## @deftypefn  {Function File} {} triplot (@var{tri}, @var{x}, @var{y})
## @deftypefnx {Function File} {} triplot (@var{tri}, @var{x}, @var{y}, @var{linespec})
## @deftypefnx {Function File} {@var{h} =} triplot (@dots{})
## Plot a triangular mesh in 2D@.  The variable @var{tri} is the triangular
## meshing of the points @code{(@var{x}, @var{y})} which is returned from
## @code{delaunay}.  If given, @var{linespec} determines the properties
## to use for the lines.  The output argument @var{h} is the graphic handle
## of the plot.
## @seealso{plot, trimesh, trisurf, delaunay}
## @end deftypefn

function h = triplot (tri, x, y, varargin)

  if (nargin < 3)
    print_usage ();
  endif

  idx = tri(:, [1, 2, 3, 1]).';
  nt = rows (tri);
  handle = plot ([x(idx); NaN(1, nt)](:),
                 [y(idx); NaN(1, nt)](:), varargin{:});

  if (nargout > 0)
    h = handle;
  endif

endfunction


%!demo
%! old_state = rand ("state");
%! restore_state = onCleanup (@() rand ("state", old_state));
%! rand ("state", 2);
%! N = 20;
%! x = rand (N, 1);
%! y = rand (N, 1);
%! tri = delaunay (x, y);
%! triplot (tri, x, y);