annotate scripts/geometry/tsearchn.m @ 14872:c2dbdeaa25df

maint: use rows() and columns() to clarify m-files. * gradient.m, interp1q.m, rat.m, tsearchn.m, image.m, imwrite.m, area.m, contourc.m, hist.m, isocolors.m, isonormals.m, meshz.m, print.m, __bar__.m, __go_draw_axes__.m, __interp_cube__.m, __marching_cube__.m, __patch__.m, __print_parse_opts__.m, __quiver__.m, rose.m, shrinkfaces.m, stairs.m, surfnorm.m, tetramesh.m, text.m, deconv.m, spline.m, intersect.m, setdiff.m, setxor.m, union.m, periodogram.m, pcg.m, perms.m: Replace size (x,1) with rows (x) and size(x,2) with columns(x).
author Rik <octave@nomad.inbox5.com>
date Tue, 17 Jul 2012 13:34:19 -0700
parents 5d3a684236b0
children d63878346099
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: 12575
diff changeset
1 ## Copyright (C) 2007-2012 David Bateman
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
2 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
3 ## This file is part of Octave.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
4 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
5 ## Octave is free software; you can redistribute it and/or modify it
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
6 ## 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
7 ## 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
8 ## your option) any later version.
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
9 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
10 ## Octave is distributed in the hope that it will be useful, but
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
11 ## WITHOUT ANY WARRANTY; without even the implied warranty of
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
12 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
13 ## General Public License for more details.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
14 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
15 ## 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
16 ## 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
17 ## <http://www.gnu.org/licenses/>.
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
18
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
19 ## -*- texinfo -*-
6846
7911a62a300d [project @ 2007-08-30 07:25:57 by dbateman]
dbateman
parents: 6823
diff changeset
20 ## @deftypefn {Function File} {[@var{idx}, @var{p}] =} tsearchn (@var{x}, @var{t}, @var{xi})
12575
d0b799dafede Grammarcheck files for 3.4.1 release.
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
21 ## Search for the enclosing Delaunay convex hull. For @code{@var{t} =
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
22 ## delaunayn (@var{x})}, finds the index in @var{t} containing the
9051
1bf0ce0930be Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents: 7017
diff changeset
23 ## points @var{xi}. For points outside the convex hull, @var{idx} is NaN.
9070
e9dc2ed2ec0f Cleanup documentation for poly.texi, interp.texi, geometry.texi
Rik <rdrider0-list@yahoo.com>
parents: 9051
diff changeset
24 ## If requested @code{tsearchn} also returns the Barycentric coordinates @var{p}
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
25 ## of the enclosing triangles.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
26 ## @seealso{delaunay, delaunayn}
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
27 ## @end deftypefn
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
28
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
29 function [idx, p] = tsearchn (x, t, xi)
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
30 if (nargin != 3)
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
31 print_usage ();
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
32 endif
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
33
14872
c2dbdeaa25df maint: use rows() and columns() to clarify m-files.
Rik <octave@nomad.inbox5.com>
parents: 14868
diff changeset
34 nt = rows (t);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
35 [m, n] = size (x);
14872
c2dbdeaa25df maint: use rows() and columns() to clarify m-files.
Rik <octave@nomad.inbox5.com>
parents: 14868
diff changeset
36 mi = rows (xi);
10548
479536c5bb10 Replace lowercase nan with NaN for visual cue in scripts
Rik <code@nomad.inbox5.com>
parents: 9245
diff changeset
37 idx = NaN (mi, 1);
479536c5bb10 Replace lowercase nan with NaN for visual cue in scripts
Rik <code@nomad.inbox5.com>
parents: 9245
diff changeset
38 p = NaN (mi, n + 1);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
39
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
40 ni = [1:mi].';
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
41 for i = 1 : nt
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
42 ## Only calculate the Barycentric coordinates for points that have not
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
43 ## already been found in a triangle.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
44 b = cart2bary (x (t (i, :), :), xi(ni,:));
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
45
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
46 ## Our points xi are in the current triangle if
14868
5d3a684236b0 maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
47 ## (all (b >= 0) && all (b <= 1)). However as we impose that
5d3a684236b0 maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
48 ## sum (b,2) == 1 we only need to test all(b>=0). Note need to add
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
49 ## a small margin for rounding errors
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
50 intri = all (b >= -1e-12, 2);
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
51 idx(ni(intri)) = i;
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
52 p(ni(intri),:) = b(intri, :);
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
53 ni(intri) = [];
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
54 endfor
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
55 endfunction
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
56
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
57 function Beta = cart2bary (T, P)
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
58 ## Conversion of Cartesian to Barycentric coordinates.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
59 ## Given a reference simplex in N dimensions represented by a
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
60 ## (N+1)-by-(N) matrix, and arbitrary point P in cartesion coordinates,
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
61 ## represented by a N-by-1 row vector can be written as
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
62 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
63 ## P = Beta * T
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
64 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
65 ## Where Beta is a N+1 vector of the barycentric coordinates. A criteria
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
66 ## on Beta is that
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
67 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
68 ## sum (Beta) == 1
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
69 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
70 ## and therefore we can write the above as
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
71 ##
14868
5d3a684236b0 maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
72 ## P - T(end, :) = Beta(1:end-1) * (T(1:end-1,:) - ones (N,1) * T(end,:))
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
73 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
74 ## and then we can solve for Beta as
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
75 ##
14868
5d3a684236b0 maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
76 ## Beta(1:end-1) = (P - T(end,:)) / (T(1:end-1,:) - ones (N,1) * T(end,:))
5d3a684236b0 maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
77 ## Beta(end) = sum (Beta)
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
78 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
79 ## Note below is generalize for multiple values of P, one per row.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
80 [M, N] = size (P);
14868
5d3a684236b0 maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
81 Beta = (P - ones (M,1) * T(end,:)) / (T(1:end-1,:) - ones (N,1) * T(end,:));
5d3a684236b0 maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
82 Beta (:,end+1) = 1 - sum (Beta, 2);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
83 endfunction
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
84
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
85
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
86 %!shared x, tri
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
87 %! x = [-1,-1;-1,1;1,-1];
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
88 %! tri = [1, 2, 3];
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
89 %!test
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
90 %! [idx, p] = tsearchn (x,tri,[-1,-1]);
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
91 %! assert (idx, 1);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
92 %! assert (p, [1,0,0], 1e-12);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
93 %!test
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
94 %! [idx, p] = tsearchn (x,tri,[-1,1]);
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
95 %! assert (idx, 1);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
96 %! assert (p, [0,1,0], 1e-12);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
97 %!test
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
98 %! [idx, p] = tsearchn (x,tri,[1,-1]);
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
99 %! assert (idx, 1);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
100 %! assert (p, [0,0,1], 1e-12);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
101 %!test
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
102 %! [idx, p] = tsearchn (x,tri,[-1/3,-1/3]);
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
103 %! assert (idx, 1);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
104 %! assert (p, [1/3,1/3,1/3], 1e-12);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
105 %!test
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
106 %! [idx, p] = tsearchn (x,tri,[1,1]);
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
107 %! assert (idx, NaN);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
108 %! assert (p, [NaN, NaN, NaN]);
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
109