view scripts/geometry/voronoi.m @ 14868:5d3a684236b0

maint: Use Octave coding conventions for cuddling parentheses in scripts directory * lin2mu.m, loadaudio.m, wavread.m, accumarray.m, bicubic.m, celldisp.m, colon.m, cplxpair.m, dblquad.m, divergence.m, genvarname.m, gradient.m, int2str.m, interp1.m, interp1q.m, interp2.m, interpn.m, loadobj.m, nthargout.m, __isequal__.m, __splinen__.m, quadgk.m, quadl.m, quadv.m, rat.m, rot90.m, rotdim.m, saveobj.m, subsindex.m, triplequad.m, delaunay3.m, griddata.m, inpolygon.m, tsearchn.m, voronoi.m, get_first_help_sentence.m, which.m, gray2ind.m, pink.m, dlmwrite.m, strread.m, textread.m, textscan.m, housh.m, ishermitian.m, issymmetric.m, krylov.m, logm.m, null.m, rref.m, compare_versions.m, copyfile.m, dump_prefs.m, edit.m, fileparts.m, getappdata.m, isappdata.m, movefile.m, orderfields.m, parseparams.m, __xzip__.m, rmappdata.m, setappdata.m, swapbytes.m, unpack.m, ver.m, fminbnd.m, fminunc.m, fsolve.m, glpk.m, lsqnonneg.m, qp.m, sqp.m, configure_make.m, copy_files.m, describe.m, get_description.m, get_forge_pkg.m, install.m, installed_packages.m, is_architecture_dependent.m, load_package_dirs.m, print_package_description.m, rebuild.m, repackage.m, save_order.m, shell.m, allchild.m, ancestor.m, area.m, axes.m, axis.m, clabel.m, close.m, colorbar.m, comet.m, comet3.m, contour.m, cylinder.m, ezmesh.m, ezsurf.m, findobj.m, fplot.m, hist.m, isocolors.m, isonormals.m, isosurface.m, isprop.m, legend.m, mesh.m, meshz.m, pareto.m, pcolor.m, peaks.m, plot3.m, plotmatrix.m, plotyy.m, polar.m, print.m, __add_datasource__.m, __add_default_menu__.m, __axes_limits__.m, __bar__.m, __clabel__.m, __contour__.m, __errcomm__.m, __errplot__.m, __ezplot__.m, __file_filter__.m, __fltk_print__.m, __ghostscript__.m, __gnuplot_print__.m, __go_draw_axes__.m, __go_draw_figure__.m, __interp_cube__.m, __marching_cube__.m, __patch__.m, __pie__.m, __plt__.m, __print_parse_opts__.m, __quiver__.m, __scatter__.m, __stem__.m, __tight_eps_bbox__.m, __uigetdir_fltk__.m, __uigetfile_fltk__.m, __uiputfile_fltk__.m, quiver.m, quiver3.m, rectangle.m, refreshdata.m, ribbon.m, scatter.m, semilogy.m, shading.m, slice.m, subplot.m, surface.m, surfl.m, surfnorm.m, text.m, uigetfile.m, uiputfile.m, whitebg.m, deconv.m, mkpp.m, pchip.m, polyaffine.m, polyder.m, polygcd.m, polyout.m, polyval.m, ppint.m, ppjumps.m, ppval.m, residue.m, roots.m, spline.m, splinefit.m, addpref.m, getpref.m, setpref.m, ismember.m, setxor.m, arch_fit.m, arch_rnd.m, arch_test.m, autoreg_matrix.m, diffpara.m, fftconv.m, filter2.m, hanning.m, hurst.m, periodogram.m, triangle_sw.m, sinc.m, spectral_xdf.m, spencer.m, stft.m, synthesis.m, unwrap.m, yulewalker.m, bicgstab.m, gmres.m, pcg.m, pcr.m, __sprand_impl__.m, speye.m, spfun.m, sprandn.m, spstats.m, svds.m, treelayout.m, treeplot.m, bessel.m, factor.m, legendre.m, perms.m, primes.m, magic.m, toeplitz.m, corr.m, cov.m, mean.m, median.m, mode.m, qqplot.m, quantile.m, ranks.m, zscore.m, logistic_regression_likelihood.m, bartlett_test.m, chisquare_test_homogeneity.m, chisquare_test_independence.m, kolmogorov_smirnov_test.m, run_test.m, u_test.m, wilcoxon_test.m, z_test.m, z_test_2.m, bin2dec.m, dec2base.m, mat2str.m, strcat.m, strchr.m, strjust.m, strtok.m, substr.m, untabify.m, assert.m, demo.m, example.m, fail.m, speed.m, test.m, now.m: Use Octave coding conventions for cuddling parentheses in scripts directory.
author Rik <octave@nomad.inbox5.com>
date Tue, 17 Jul 2012 07:08:39 -0700
parents ce2b59a6d0e5
children 02952657182e
line wrap: on
line source

## Copyright (C) 2000-2012 Kai Habel
##
## 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} {} voronoi (@var{x}, @var{y})
## @deftypefnx {Function File} {} voronoi (@var{x}, @var{y}, @var{options})
## @deftypefnx {Function File} {} voronoi (@dots{}, "linespec")
## @deftypefnx {Function File} {} voronoi (@var{hax}, @dots{})
## @deftypefnx {Function File} {@var{h} =} voronoi (@dots{})
## @deftypefnx {Function File} {[@var{vx}, @var{vy}] =} voronoi (@dots{})
## Plot the Voronoi diagram of points @code{(@var{x}, @var{y})}.
## The Voronoi facets with points at infinity are not drawn.
## 
## If "linespec" is given it is used to set the color and line style of the
## plot.  If an axis graphics handle @var{hax} is supplied then the Voronoi
## diagram is drawn on the specified axis rather than in a new figure.
##
## The @var{options} argument, which must be a string or cell array of strings,
## contains options passed to the underlying qhull command.
## See the documentation for the Qhull library for details
## @url{http://www.qhull.org/html/qh-quick.htm#options}.
##
## If a single output argument is requested then the Voronoi diagram will be
## plotted and a graphics handle @var{h} to the plot is returned.
## [@var{vx}, @var{vy}] = voronoi (@dots{}) returns the Voronoi vertices
## instead of plotting the diagram.
##
## @example
## @group
## x = rand (10, 1);
## y = rand (size (x));
## h = convhull (x, y);
## [vx, vy] = voronoi (x, y);
## plot (vx, vy, "-b", x, y, "o", x(h), y(h), "-g");
## legend ("", "points", "hull");
## @end group
## @end example
##
## @seealso{voronoin, delaunay, convhull}
## @end deftypefn

## Author: Kai Habel <kai.habel@gmx.de>
## First Release: 20/08/2000

## 2002-01-04 Paul Kienzle <pkienzle@users.sf.net>
## * limit the default graph to the input points rather than the whole diagram
## * provide example
## * use unique(x,"rows") rather than __unique_rows__

## 2003-12-14 Rafael Laboissiere <rafael@laboissiere.net>
## Added optional fourth argument to pass options to the underlying
## qhull command

function [vx, vy] = voronoi (varargin)

  if (nargin < 1)
    print_usage ();
  endif

  narg = 1;
  if (isscalar (varargin{1}) && ishandle (varargin{1}))
    handl = varargin{1};
    if (! strcmp (get (handl, "type"), "axes"))
      error ("voronoi: expecting first argument to be an axes object");
    endif
    narg++;
  elseif (nargout < 2)
    handl = gca ();
  endif

  if (nargin < 1 + narg || nargin > 3 + narg)
    print_usage ();
  endif

  x = varargin{narg++};
  y = varargin{narg++};

  opts = {};
  if (narg <= nargin)
    if (iscell (varargin{narg}))
      opts = varargin(narg++);
    elseif (isnumeric (varargin{narg}))
      ## Accept, but ignore, the triangulation
      narg++;
    endif
  endif

  linespec = {"b"};
  if (narg <= nargin && ischar (varargin{narg}))
    linespec = varargin(narg);
  endif

  lx = length (x);
  ly = length (y);

  if (lx != ly)
    error ("voronoi: X and Y must be vectors of the same length");
  endif

  ## Add box to approximate rays to infinity. For Voronoi diagrams the
  ## box can (and should) be close to the points themselves. To make the
  ## job of finding the exterior edges it should be at least two times the
  ## delta below however
  xmax = max (x(:));
  xmin = min (x(:));
  ymax = max (y(:));
  ymin = min (y(:));
  xdelta = xmax - xmin;
  ydelta = ymax - ymin;
  scale = 2;

  xbox = [xmin - scale * xdelta; xmin - scale * xdelta; ...
          xmax + scale * xdelta; xmax + scale * xdelta];
  ybox = [xmin - scale * xdelta; xmax + scale * xdelta; ...
          xmax + scale * xdelta; xmin - scale * xdelta];

  [p, c, infi] = __voronoi__ ("voronoi",
                              [[x(:) ; xbox(:)], [y(:); ybox(:)]],
                              opts{:});

  idx = find (! infi);
  ll = length (idx);
  c = c(idx).';
  k = sum (cellfun ("length", c));
  edges = cell2mat (cellfun (@(x) [x ; [x(end), x(1:end-1)]], c,
                             "uniformoutput", false));

  ## Identify the unique edges of the Voronoi diagram
  edges = sortrows (sort (edges).').';
  edges = edges (:, [(edges(1, 1: end - 1) != edges(1, 2 : end) | ...
                      edges(2, 1 :end - 1) != edges(2, 2 : end)), true]);

  ## Eliminate the edges of the diagram representing the box
  poutside = (1 : rows (p)) ...
      (p (:, 1) < xmin - xdelta | p (:, 1) > xmax + xdelta | ...
       p (:, 2) < ymin - ydelta | p (:, 2) > ymax + ydelta);
  edgeoutside = ismember (edges (1, :), poutside) & ...
      ismember (edges (2, :), poutside);
  edges (:, edgeoutside) = [];

  ## Get points of the diagram
  Vvx = reshape (p(edges, 1), size (edges));
  Vvy = reshape (p(edges, 2), size (edges));

  if (nargout < 2)
    lim = [xmin, xmax, ymin, ymax];
    h = plot (handl, Vvx, Vvy, linespec{:}, x, y, '+');
    axis (lim + 0.1 * [[-1, 1] * (lim (2) - lim (1)), ...
                       [-1, 1] * (lim (4) - lim (3))]);
    if (nargout == 1)
      vx = h;
    endif
  else
    vx = Vvx;
    vy = Vvy;
  endif

endfunction


%!demo
%! voronoi (rand (10,1), rand (10,1));

%!testif HAVE_QHULL
%! phi = linspace (-pi, 3/4*pi, 8);
%! [x,y] = pol2cart (phi, 1);
%! [vx,vy] = voronoi (x,y);
%! assert (vx(2,:), zeros (1, columns (vx)), eps);
%! assert (vy(2,:), zeros (1, columns (vy)), eps);

%% FIXME: Need input validation tests