changelog shortlog tags changeset files revisions annotate raw

scripts/geometry/tsearchn.m

changeset 10289: 4b124317dc38
parent:e9dc2ed2ec0f
author: John W. Eaton <jwe@octave.org>
date: Tue Feb 09 20:58:55 2010 -0500 (42 minutes ago)
permissions: -rw-r--r--
description: base_properties::set_children: account for hidden children
1## Copyright (C) 2007, 2009 David Bateman
2##
3## This file is part of Octave.
4##
5## Octave is free software; you can redistribute it and/or modify it
6## under the terms of the GNU General Public License as published by
7## the Free Software Foundation; either version 3 of the License, or (at
8## your option) any later version.
9##
10## Octave is distributed in the hope that it will be useful, but
11## WITHOUT ANY WARRANTY; without even the implied warranty of
12## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13## General Public License for more details.
14##
15## You should have received a copy of the GNU General Public License
16## along with Octave; see the file COPYING. If not, see
17## <http://www.gnu.org/licenses/>.
18
19## -*- texinfo -*-
20## @deftypefn {Function File} {[@var{idx}, @var{p}] =} tsearchn (@var{x}, @var{t}, @var{xi})
21## Searches for the enclosing Delaunay convex hull. For @code{@var{t} =
22## delaunayn (@var{x})}, finds the index in @var{t} containing the
23## points @var{xi}. For points outside the convex hull, @var{idx} is NaN.
24## If requested @code{tsearchn} also returns the Barycentric coordinates @var{p}
25## of the enclosing triangles.
26## @seealso{delaunay, delaunayn}
27## @end deftypefn
28
29function [idx, p] = tsearchn (x, t, xi)
30 if (nargin != 3)
31 print_usage ();
32 endif
33
34 nt = size (t, 1);
35 [m, n] = size (x);
36 mi = size (xi, 1);
37 idx = nan (mi, 1);
38 p = nan (mi, n + 1);
39
40 ni = [1:mi].';
41 for i = 1 : nt
42 ## Only calculate the Barycentric coordinates for points that have not
43 ## already been found in a triangle.
44 b = cart2bary (x (t (i, :), :), xi(ni,:));
45
46 ## Our points xi are in the current triangle if
47 ## (all(b >= 0) && all (b <= 1)). However as we impose that
48 ## sum(b,2) == 1 we only need to test all(b>=0). Note need to add
49 ## a small margin for rounding errors
50 intri = all (b >= -1e-12, 2);
51 idx(ni(intri)) = i;
52 p(ni(intri),:) = b(intri, :);
53 ni(intri) = [];
54 endfor
55endfunction
56
57function Beta = cart2bary (T, P)
58 ## Conversion of Cartesian to Barycentric coordinates.
59 ## Given a reference simplex in N dimensions represented by a
60 ## (N+1)-by-(N) matrix, and arbitrary point P in cartesion coordinates,
61 ## represented by a N-by-1 row vector can be written as
62 ##
63 ## P = Beta * T
64 ##
65 ## Where Beta is a N+1 vector of the barycentric coordinates. A criteria
66 ## on Beta is that
67 ##
68 ## sum (Beta) == 1
69 ##
70 ## and therefore we can write the above as
71 ##
72 ## P - T(end, :) = Beta(1:end-1) * (T(1:end-1,:) - ones(N,1) * T(end,:))
73 ##
74 ## and then we can solve for Beta as
75 ##
76 ## Beta(1:end-1) = (P - T(end,:)) / (T(1:end-1,:) - ones(N,1) * T(end,:))
77 ## Beta(end) = sum(Beta)
78 ##
79 ## Note below is generalize for multiple values of P, one per row.
80 [M, N] = size (P);
81 Beta = (P - ones (M,1) * T(end,:)) / (T(1:end-1,:) - ones(N,1) * T(end,:));
82 Beta (:,end+1) = 1 - sum(Beta, 2);
83endfunction
84
85%!shared x, tri
86%! x = [-1,-1;-1,1;1,-1];
87%! tri = [1, 2, 3];
88%!test
89%! [idx, p] = tsearchn (x,tri,[-1,-1]);
90%! assert (idx, 1)
91%! assert (p, [1,0,0], 1e-12)
92%!test
93%! [idx, p] = tsearchn (x,tri,[-1,1]);
94%! assert (idx, 1)
95%! assert (p, [0,1,0], 1e-12)
96%!test
97%! [idx, p] = tsearchn (x,tri,[1,-1]);
98%! assert (idx, 1)
99%! assert (p, [0,0,1], 1e-12)
100%!test
101%! [idx, p] = tsearchn (x,tri,[-1/3,-1/3]);
102%! assert (idx, 1)
103%! assert (p, [1/3,1/3,1/3], 1e-12)
104%!test
105%! [idx, p] = tsearchn (x,tri,[1,1]);
106%! assert (idx, NaN)
107%! assert (p, [NaN, NaN, NaN])