Mercurial > octave
view scripts/geometry/delaunayn.m @ 21546:f7f97d7e9294
doc: Wrap m-file docstrings to 79 characters + newline (80 total).
* isrecording.m, soundsc.m, delaunay3.m, cell2mat.m, cumtrapz.m, del2.m,
inputParser.m, interp1.m, interp3.m, narginchk.m, profile.m,
validateattributes.m, delaunayn.m, tsearchn.m, voronoin.m, brighten.m,
cmunique.m, colorcube.m, imfinfo.m, imshow.m, edit.m, orderfields.m, run.m,
warning_ids.m, ode23.m, ode45.m, odeget.m, integrate_adaptive.m, kahan.m,
ode_struct_value_check.m, runge_kutta_23.m, fminunc.m, fsolve.m, fzero.m,
pkg.m, build.m, specular.m, view.m, bar.m, barh.m, contour3.m, isosurface.m,
line.m, pie.m, pie3.m, quiver3.m, scatter.m, scatter3.m, stem3.m, stemleaf.m,
surfl.m, tetramesh.m, isfigure.m, mkpp.m, pchip.m, residue.m, splinefit.m,
rmpref.m, unique.m, eigs.m, ilu.m, factor.m, factorial.m, gallery.m, hankel.m,
histc.m, ols.m, finv.m, fpdf.m, kruskal_wallis_test.m, weekday.m:
Wrap m-file docstrings to 79 characters + newline (80 total).
author | Rik <rik@octave.org> |
---|---|
date | Sun, 27 Mar 2016 15:50:01 -0700 |
parents | 516bb87ea72e |
children | ecce63c99c3f |
line wrap: on
line source
## Copyright (C) 2007-2015 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 {} {@var{T} =} delaunayn (@var{pts}) ## @deftypefnx {} {@var{T} =} delaunayn (@var{pts}, @var{options}) ## Compute the Delaunay triangulation for an N-dimensional set of points. ## ## The Delaunay triangulation is a tessellation of the convex hull of a set of ## points such that no N-sphere defined by the N-triangles contains any other ## points from the set. ## ## The input matrix @var{pts} of size [n, dim] contains n points in a space of ## dimension dim. The return matrix @var{T} has size [m, dim+1]. Each row of ## @var{T} contains a set of indices back into the original set of points ## @var{pts} which describes a simplex of dimension dim. For example, a 2-D ## simplex is a triangle and 3-D simplex is a tetrahedron. ## ## An optional second 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}. ## The default options depend on the dimension of the input: ## ## @itemize ## @item 2-D and 3-D: @var{options} = @code{@{"Qt", "Qbb", "Qc", "Qz"@}} ## ## @item 4-D and higher: @var{options} = @code{@{"Qt", "Qbb", "Qc", "Qx"@}} ## @end itemize ## ## If @var{options} is not present or @code{[]} then the default arguments are ## used. Otherwise, @var{options} replaces the default argument list. ## To append user options to the defaults it is necessary to repeat the ## default arguments in @var{options}. Use a null string to pass no arguments. ## ## @seealso{delaunay, convhulln, voronoin, trimesh, tetramesh} ## @end deftypefn function T = delaunayn (pts, varargin) if (nargin < 1) print_usage (); endif T = __delaunayn__ (pts, varargin{:}); if (isa (pts, "single")) tol = 1e3 * eps ("single"); else tol = 1e3 * eps; endif ## Try to remove the zero volume simplices. The volume of the i-th simplex is ## given by abs(det(pts(T(i,1:end-1),:)-pts(T(i,2:end),:)))/prod(1:n) ## (reference http://en.wikipedia.org/wiki/Simplex). Any simplex with a ## relative volume less than some arbitrary criteria is rejected. The ## criteria we use is the volume of the simplex corresponding to an ## orthogonal simplex is equal edge length all equal to the edge length of ## the original simplex. If the relative volume is 1e3*eps then the simplex ## is rejected. Note division of the two volumes means that the factor ## prod(1:n) is dropped. idx = []; [nt, n] = size (T); ## FIXME: Vectorize this for loop or convert delaunayn to .oct function for i = 1:nt X = pts(T(i,1:end-1),:) - pts(T(i,2:end),:); if (abs (det (X)) / sqrt (sumsq (X, 2)) < tol) idx(end+1) = i; endif endfor T(idx,:) = []; endfunction %!testif HAVE_QHULL %! x = [-1, 0; 0, 1; 1, 0; 0, -1; 0, 0]; %! assert (sortrows (sort (delaunayn (x), 2)), [1,2,5;1,4,5;2,3,5;3,4,5]); ## Test 3-D input %!testif HAVE_QHULL %! x = [-1, -1, 1, 0, -1]; y = [-1, 1, 1, 0, -1]; z = [0, 0, 0, 1, 1]; %! assert (sortrows (sort (delaunayn ([x(:) y(:) z(:)]), 2)), [1,2,3,4;1,2,4,5]) ## FIXME: Need tests for delaunayn ## Input validation tests %!error delaunayn ()