Mercurial > octave
annotate scripts/geometry/tsearchn.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 | b571fc85953f |
rev | line source |
---|---|
19697
4197fc428c7d
maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents:
17744
diff
changeset
|
1 ## Copyright (C) 2007-2015 David Bateman |
6823 | 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 | |
7016 | 7 ## the Free Software Foundation; either version 3 of the License, or (at |
8 ## your option) any later version. | |
6823 | 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 | |
7016 | 16 ## along with Octave; see the file COPYING. If not, see |
17 ## <http://www.gnu.org/licenses/>. | |
6823 | 18 |
19 ## -*- texinfo -*- | |
20852
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20158
diff
changeset
|
20 ## @deftypefn {} {@var{idx} =} tsearchn (@var{x}, @var{t}, @var{xi}) |
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20158
diff
changeset
|
21 ## @deftypefnx {} {[@var{idx}, @var{p}] =} tsearchn (@var{x}, @var{t}, @var{xi}) |
20158
7503499a252b
doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
22 ## Search for the enclosing Delaunay convex hull. |
7503499a252b
doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
23 ## |
7503499a252b
doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
24 ## For @code{@var{t} = delaunayn (@var{x})}, finds the index in @var{t} |
7503499a252b
doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
25 ## containing the points @var{xi}. For points outside the convex hull, |
7503499a252b
doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
26 ## @var{idx} is NaN. |
7503499a252b
doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
27 ## |
21546
f7f97d7e9294
doc: Wrap m-file docstrings to 79 characters + newline (80 total).
Rik <rik@octave.org>
parents:
20852
diff
changeset
|
28 ## If requested @code{tsearchn} also returns the Barycentric coordinates |
f7f97d7e9294
doc: Wrap m-file docstrings to 79 characters + newline (80 total).
Rik <rik@octave.org>
parents:
20852
diff
changeset
|
29 ## @var{p} of the enclosing triangles. |
6823 | 30 ## @seealso{delaunay, delaunayn} |
31 ## @end deftypefn | |
32 | |
33 function [idx, p] = tsearchn (x, t, xi) | |
34 if (nargin != 3) | |
35 print_usage (); | |
36 endif | |
37 | |
14872
c2dbdeaa25df
maint: use rows() and columns() to clarify m-files.
Rik <octave@nomad.inbox5.com>
parents:
14868
diff
changeset
|
38 nt = rows (t); |
6823 | 39 [m, n] = size (x); |
14872
c2dbdeaa25df
maint: use rows() and columns() to clarify m-files.
Rik <octave@nomad.inbox5.com>
parents:
14868
diff
changeset
|
40 mi = rows (xi); |
10548
479536c5bb10
Replace lowercase nan with NaN for visual cue in scripts
Rik <code@nomad.inbox5.com>
parents:
9245
diff
changeset
|
41 idx = NaN (mi, 1); |
479536c5bb10
Replace lowercase nan with NaN for visual cue in scripts
Rik <code@nomad.inbox5.com>
parents:
9245
diff
changeset
|
42 p = NaN (mi, n + 1); |
6823 | 43 |
44 ni = [1:mi].'; | |
45 for i = 1 : nt | |
46 ## Only calculate the Barycentric coordinates for points that have not | |
47 ## already been found in a triangle. | |
48 b = cart2bary (x (t (i, :), :), xi(ni,:)); | |
49 | |
50 ## 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
|
51 ## (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
|
52 ## sum (b,2) == 1 we only need to test all(b>=0). Note need to add |
6823 | 53 ## a small margin for rounding errors |
54 intri = all (b >= -1e-12, 2); | |
55 idx(ni(intri)) = i; | |
56 p(ni(intri),:) = b(intri, :); | |
57 ni(intri) = []; | |
58 endfor | |
59 endfunction | |
60 | |
61 function Beta = cart2bary (T, P) | |
62 ## Conversion of Cartesian to Barycentric coordinates. | |
63 ## Given a reference simplex in N dimensions represented by a | |
64 ## (N+1)-by-(N) matrix, and arbitrary point P in cartesion coordinates, | |
65 ## represented by a N-by-1 row vector can be written as | |
66 ## | |
67 ## P = Beta * T | |
68 ## | |
69 ## Where Beta is a N+1 vector of the barycentric coordinates. A criteria | |
70 ## on Beta is that | |
71 ## | |
72 ## sum (Beta) == 1 | |
73 ## | |
74 ## and therefore we can write the above as | |
75 ## | |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
76 ## P - T(end, :) = Beta(1:end-1) * (T(1:end-1,:) - ones (N,1) * T(end,:)) |
6823 | 77 ## |
78 ## and then we can solve for Beta as | |
79 ## | |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
80 ## 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
|
81 ## Beta(end) = sum (Beta) |
6823 | 82 ## |
83 ## Note below is generalize for multiple values of P, one per row. | |
84 [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
|
85 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
|
86 Beta (:,end+1) = 1 - sum (Beta, 2); |
6823 | 87 endfunction |
88 | |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
89 |
6823 | 90 %!shared x, tri |
91 %! x = [-1,-1;-1,1;1,-1]; | |
92 %! tri = [1, 2, 3]; | |
93 %!test | |
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, [1,0,0], 1e-12); |
6823 | 97 %!test |
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,1,0], 1e-12); |
6823 | 101 %!test |
102 %! [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
|
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, [0,0,1], 1e-12); |
6823 | 105 %!test |
106 %! [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
|
107 %! assert (idx, 1); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
108 %! assert (p, [1/3,1/3,1/3], 1e-12); |
6823 | 109 %!test |
110 %! [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
|
111 %! assert (idx, NaN); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
112 %! 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
|
113 |