Mercurial > octave
annotate scripts/geometry/tsearchn.m @ 31191:bb9d776eafac stable
Fix wrong color in PDF printout of some latex strings (bug #62884)
* octave-svgconvert (draw): For "rect" elements only set brush color if
necessary and eventually restore to previous color.
author | Pantxo Diribarne <pantxo.diribarne@gmail.com> |
---|---|
date | Sun, 14 Aug 2022 18:24:07 +0200 |
parents | 796f54d4ddbf |
children | 72ef3d097059 597f3ee61a48 |
rev | line source |
---|---|
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
1 ######################################################################## |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
2 ## |
30564
796f54d4ddbf
update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents:
30379
diff
changeset
|
3 ## Copyright (C) 2007-2022 The Octave Project Developers |
27918
b442ec6dda5c
use centralized file for copyright info for individual contributors
John W. Eaton <jwe@octave.org>
parents:
26376
diff
changeset
|
4 ## |
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
5 ## See the file COPYRIGHT.md in the top-level directory of this |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
6 ## distribution or <https://octave.org/copyright/>. |
6823 | 7 ## |
8 ## This file is part of Octave. | |
9 ## | |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
10 ## Octave is free software: you can redistribute it and/or modify it |
6823 | 11 ## under the terms of the GNU General Public License as published by |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
12 ## the Free Software Foundation, either version 3 of the License, or |
22755
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
13 ## (at your option) any later version. |
6823 | 14 ## |
15 ## Octave is distributed in the hope that it will be useful, but | |
16 ## WITHOUT ANY WARRANTY; without even the implied warranty of | |
22755
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
17 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
3a2b891d0b33
maint: Standardize Copyright formatting.
Rik <rik@octave.org>
parents:
22323
diff
changeset
|
18 ## GNU General Public License for more details. |
6823 | 19 ## |
20 ## You should have received a copy of the GNU General Public License | |
7016 | 21 ## along with Octave; see the file COPYING. If not, see |
24534
194eb4bd202b
maint: Update punctuation for GPL v3 license text.
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
22 ## <https://www.gnu.org/licenses/>. |
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
23 ## |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
24 ######################################################################## |
6823 | 25 |
26 ## -*- texinfo -*- | |
20852
516bb87ea72e
2015 Code Sprint: remove class of function from docstring for all m-files.
Rik <rik@octave.org>
parents:
20158
diff
changeset
|
27 ## @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
|
28 ## @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
|
29 ## 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
|
30 ## |
7503499a252b
doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
31 ## 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
|
32 ## 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
|
33 ## @var{idx} is NaN. |
7503499a252b
doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
34 ## |
21546
f7f97d7e9294
doc: Wrap m-file docstrings to 79 characters + newline (80 total).
Rik <rik@octave.org>
parents:
20852
diff
changeset
|
35 ## 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
|
36 ## @var{p} of the enclosing triangles. |
6823 | 37 ## @seealso{delaunay, delaunayn} |
38 ## @end deftypefn | |
39 | |
40 function [idx, p] = tsearchn (x, t, xi) | |
21758
ffad2baa90f7
maint: Use newlines to make code more readable.
Rik <rik@octave.org>
parents:
21751
diff
changeset
|
41 |
6823 | 42 if (nargin != 3) |
43 print_usage (); | |
44 endif | |
45 | |
14872
c2dbdeaa25df
maint: use rows() and columns() to clarify m-files.
Rik <octave@nomad.inbox5.com>
parents:
14868
diff
changeset
|
46 nt = rows (t); |
6823 | 47 [m, n] = size (x); |
14872
c2dbdeaa25df
maint: use rows() and columns() to clarify m-files.
Rik <octave@nomad.inbox5.com>
parents:
14868
diff
changeset
|
48 mi = rows (xi); |
10548
479536c5bb10
Replace lowercase nan with NaN for visual cue in scripts
Rik <code@nomad.inbox5.com>
parents:
9245
diff
changeset
|
49 idx = NaN (mi, 1); |
479536c5bb10
Replace lowercase nan with NaN for visual cue in scripts
Rik <code@nomad.inbox5.com>
parents:
9245
diff
changeset
|
50 p = NaN (mi, n + 1); |
6823 | 51 |
52 ni = [1:mi].'; | |
53 for i = 1 : nt | |
54 ## Only calculate the Barycentric coordinates for points that have not | |
55 ## already been found in a triangle. | |
56 b = cart2bary (x (t (i, :), :), xi(ni,:)); | |
57 | |
58 ## Our points xi are in the current triangle if | |
21751
b571fc85953f
maint: Use two spaces after period to indicate sentence break.
Rik <rik@octave.org>
parents:
21546
diff
changeset
|
59 ## (all (b >= 0) && all (b <= 1)). However as we impose that |
b571fc85953f
maint: Use two spaces after period to indicate sentence break.
Rik <rik@octave.org>
parents:
21546
diff
changeset
|
60 ## sum (b,2) == 1 we only need to test all(b>=0). Note need to add |
6823 | 61 ## a small margin for rounding errors |
62 intri = all (b >= -1e-12, 2); | |
63 idx(ni(intri)) = i; | |
64 p(ni(intri),:) = b(intri, :); | |
65 ni(intri) = []; | |
66 endfor | |
21758
ffad2baa90f7
maint: Use newlines to make code more readable.
Rik <rik@octave.org>
parents:
21751
diff
changeset
|
67 |
6823 | 68 endfunction |
69 | |
70 function Beta = cart2bary (T, P) | |
30379
363fb10055df
maint: Style check m-files ahead of 7.1 release.
Rik <rik@octave.org>
parents:
29358
diff
changeset
|
71 |
6823 | 72 ## Conversion of Cartesian to Barycentric coordinates. |
24560
1c1adf6ab75d
doc: Fix issues in geometry, polynomial, and interpolation chapters (bug #52835).
Rik <rik@octave.org>
parents:
23219
diff
changeset
|
73 ## Given a reference simplex in N dimensions represented by an |
1c1adf6ab75d
doc: Fix issues in geometry, polynomial, and interpolation chapters (bug #52835).
Rik <rik@octave.org>
parents:
23219
diff
changeset
|
74 ## N+1-by-N matrix, an arbitrary point P in Cartesian coordinates, |
1c1adf6ab75d
doc: Fix issues in geometry, polynomial, and interpolation chapters (bug #52835).
Rik <rik@octave.org>
parents:
23219
diff
changeset
|
75 ## represented by an N-by-1 column vector can be written as |
6823 | 76 ## |
77 ## P = Beta * T | |
78 ## | |
24560
1c1adf6ab75d
doc: Fix issues in geometry, polynomial, and interpolation chapters (bug #52835).
Rik <rik@octave.org>
parents:
23219
diff
changeset
|
79 ## Where Beta is an N+1 vector of the barycentric coordinates. A criteria |
6823 | 80 ## on Beta is that |
81 ## | |
82 ## sum (Beta) == 1 | |
83 ## | |
84 ## and therefore we can write the above as | |
85 ## | |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
86 ## P - T(end, :) = Beta(1:end-1) * (T(1:end-1,:) - ones (N,1) * T(end,:)) |
6823 | 87 ## |
88 ## and then we can solve for Beta as | |
89 ## | |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
90 ## 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
|
91 ## Beta(end) = sum (Beta) |
6823 | 92 ## |
24560
1c1adf6ab75d
doc: Fix issues in geometry, polynomial, and interpolation chapters (bug #52835).
Rik <rik@octave.org>
parents:
23219
diff
changeset
|
93 ## Note code below is generalized for multiple values of P, one per row. |
6823 | 94 [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
|
95 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
|
96 Beta (:,end+1) = 1 - sum (Beta, 2); |
21758
ffad2baa90f7
maint: Use newlines to make code more readable.
Rik <rik@octave.org>
parents:
21751
diff
changeset
|
97 |
6823 | 98 endfunction |
99 | |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
100 |
6823 | 101 %!shared x, tri |
102 %! x = [-1,-1;-1,1;1,-1]; | |
103 %! tri = [1, 2, 3]; | |
104 %!test | |
105 %! [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
|
106 %! assert (idx, 1); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
107 %! assert (p, [1,0,0], 1e-12); |
6823 | 108 %!test |
109 %! [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
|
110 %! assert (idx, 1); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
111 %! assert (p, [0,1,0], 1e-12); |
6823 | 112 %!test |
113 %! [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
|
114 %! assert (idx, 1); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
115 %! assert (p, [0,0,1], 1e-12); |
6823 | 116 %!test |
117 %! [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
|
118 %! assert (idx, 1); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
119 %! assert (p, [1/3,1/3,1/3], 1e-12); |
6823 | 120 %!test |
121 %! [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
|
122 %! assert (idx, NaN); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
123 %! assert (p, [NaN, NaN, NaN]); |