Mercurial > hg > octave-nkf
annotate scripts/plot/draw/shrinkfaces.m @ 18402:6165053c56b3
shrinkfaces.m: Overhaul function.
* shrinkfaces.m: Improve docstring. Add %!tests for input validation.
author | Rik <rik@octave.org> |
---|---|
date | Wed, 22 Jan 2014 17:20:40 -0800 |
parents | d63878346099 |
children | d3a223128efc |
rev | line source |
---|---|
17744
d63878346099
maint: Update copyright notices for release.
John W. Eaton <jwe@octave.org>
parents:
17572
diff
changeset
|
1 ## Copyright (C) 2012-2013 Martin Helm |
14574 | 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 -*- | |
14620
cd375519eab0
doc: Periodic grammar check of documentation
Rik <octave@nomad.inbox5.com>
parents:
14574
diff
changeset
|
20 ## @deftypefn {Function File} {} shrinkfaces (@var{p}, @var{sf}) |
14574 | 21 ## @deftypefnx {Function File} {@var{nfv} =} shrinkfaces (@var{p}, @var{sf}) |
22 ## @deftypefnx {Function File} {@var{nfv} =} shrinkfaces (@var{fv}, @var{sf}) | |
23 ## @deftypefnx {Function File} {@var{nfv} =} shrinkfaces (@var{f}, @var{v}, @var{sf}) | |
24 ## @deftypefnx {Function File} {[@var{nf}, @var{nv}] =} shrinkfaces (@dots{}) | |
25 ## | |
18402 | 26 ## Reduce the size of faces in a patch by the shrink factor @var{sf}. |
27 ## | |
28 ## The patch object can be specified by a graphics handle (@var{p}), a patch | |
29 ## structure (@var{fv}) with the fields @qcode{"faces"} and @qcode{"vertices"}, | |
30 ## or as two separate matrices (@var{f}, @var{v}) of faces and vertices. | |
14574 | 31 ## |
18402 | 32 ## The shrink factor @var{sf} is a positive number specifying the percentage |
33 ## of the original area the new face will occupy. If no factor is given the | |
34 ## default is 0.3 (a reduction to 30% of the original size). A factor greater | |
35 ## than 1.0 will result in the expansion of faces. | |
36 ## | |
37 ## Given a patch handle as the first input argument and no output parameters, | |
38 ## perform the shrinking of the patch faces in place and redraw the patch. | |
14574 | 39 ## |
40 ## If called with one output argument, return a structure with fields | |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17122
diff
changeset
|
41 ## @qcode{"faces"}, @qcode{"vertices"}, and @qcode{"facevertexcdata"} |
18402 | 42 ## containing the data after shrinking. This structure can be used directly |
43 ## as an input argument to the @code{patch} function. | |
14574 | 44 ## |
18402 | 45 ## @strong{Caution:}: Performing the shrink operation on faces which are not |
46 ## convex can lead to undesirable results. | |
14574 | 47 ## |
18402 | 48 ## Example: a triangulated 3/4 circle and the corresponding shrunken version. |
14574 | 49 ## |
50 ## @example | |
51 ## @group | |
52 ## [phi r] = meshgrid (linspace (0, 1.5*pi, 16), linspace (1, 2, 4)); | |
53 ## tri = delaunay (phi(:), r(:)); | |
54 ## v = [r(:).*sin(phi(:)) r(:).*cos(phi(:))]; | |
55 ## clf () | |
56 ## p = patch ("Faces", tri, "Vertices", v, "FaceColor", "none"); | |
57 ## fv = shrinkfaces (p); | |
58 ## patch (fv) | |
59 ## axis equal | |
60 ## grid on | |
61 ## @end group | |
62 ## @end example | |
63 ## | |
64 ## @seealso{patch} | |
65 ## @end deftypefn | |
66 | |
67 ## Author: Martin Helm <martin@mhelm.de> | |
68 | |
69 function [nf, nv] = shrinkfaces (varargin) | |
70 | |
71 if (nargin < 1 || nargin > 3 || nargout > 2) | |
72 print_usage (); | |
73 endif | |
74 | |
75 sf = 0.3; | |
18402 | 76 colors = []; |
14574 | 77 p = varargin{1}; |
78 | |
79 if (ishandle (p) && nargin < 3) | |
80 faces = get (p, "Faces"); | |
81 vertices = get (p, "Vertices"); | |
82 colors = get (p, "FaceVertexCData"); | |
83 if (nargin == 2) | |
84 sf = varargin{2}; | |
85 endif | |
86 elseif (isstruct (p) && nargin < 3) | |
87 faces = p.faces; | |
88 vertices = p.vertices; | |
89 if (isfield (p, "facevertexcdata")) | |
90 colors = p.facevertexcdata; | |
91 endif | |
92 if (nargin == 2) | |
93 sf = varargin{2}; | |
94 endif | |
95 elseif (ismatrix (p) && nargin >= 2 && ismatrix (varargin{2})) | |
96 faces = p; | |
97 vertices = varargin{2}; | |
98 if (nargin == 3) | |
99 sf = varargin{3}; | |
100 endif | |
101 else | |
102 print_usage (); | |
103 endif | |
104 | |
105 if (! isscalar (sf) || sf <= 0) | |
15202
f3b5cadfd6d5
fix missing semicolons in various .m files
John W. Eaton <jwe@octave.org>
parents:
14872
diff
changeset
|
106 error ("shrinkfaces: scale factor must be a positive scalar"); |
14574 | 107 endif |
108 | |
18402 | 109 nc = columns (vertices); |
110 if (nc < 2 || nc > 3) | |
17122
eaab03308c0b
doc: Rewrite docstrings for most plot functions.
Rik <rik@octave.org>
parents:
16828
diff
changeset
|
111 error ("shrinkfaces: only 2-D and 3-D patches are supported"); |
14574 | 112 endif |
113 | |
14872
c2dbdeaa25df
maint: use rows() and columns() to clarify m-files.
Rik <octave@nomad.inbox5.com>
parents:
14621
diff
changeset
|
114 m = columns (faces); |
14574 | 115 if (m < 3) |
15202
f3b5cadfd6d5
fix missing semicolons in various .m files
John W. Eaton <jwe@octave.org>
parents:
14872
diff
changeset
|
116 error ("shrinkfaces: faces must consist of at least 3 vertices"); |
14574 | 117 endif |
118 | |
119 v = vertices(faces'(:), :); | |
14872
c2dbdeaa25df
maint: use rows() and columns() to clarify m-files.
Rik <octave@nomad.inbox5.com>
parents:
14621
diff
changeset
|
120 if (isempty (colors) || rows (colors) == rows (faces)) |
14574 | 121 c = colors; |
14872
c2dbdeaa25df
maint: use rows() and columns() to clarify m-files.
Rik <octave@nomad.inbox5.com>
parents:
14621
diff
changeset
|
122 elseif (rows (colors) == rows (vertices)) |
14574 | 123 c = colors(faces'(:), :); |
124 else | |
18402 | 125 c = []; # Discard inconsistent color data. |
14574 | 126 endif |
14872
c2dbdeaa25df
maint: use rows() and columns() to clarify m-files.
Rik <octave@nomad.inbox5.com>
parents:
14621
diff
changeset
|
127 sv = rows (v); |
18402 | 128 ## We have to deal with a possibly very large number of vertices, so use |
129 ## sparse as midpoint (1/m, ..., 1/m) in generalized barycentric coordinates. | |
130 midpoints = full (kron (speye (sv / m), ones (m, m) / m) * sparse (v)); | |
14574 | 131 v = sqrt (sf) * (v - midpoints) + midpoints; |
132 f = reshape (1:sv, m, sv / m)'; | |
133 | |
134 switch (nargout) | |
135 case 0 | |
136 if (ishandle (p)) | |
15202
f3b5cadfd6d5
fix missing semicolons in various .m files
John W. Eaton <jwe@octave.org>
parents:
14872
diff
changeset
|
137 ## avoid exceptions |
f3b5cadfd6d5
fix missing semicolons in various .m files
John W. Eaton <jwe@octave.org>
parents:
14872
diff
changeset
|
138 set (p, "FaceVertexCData", [], "CData", []); |
15332
b613757ff5be
Fix typo in shrinkfaces.m (thanks Marco Atzeri)
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
15202
diff
changeset
|
139 set (p, "Vertices", v, "Faces", f, "FaceVertexCData", c); |
14574 | 140 else |
141 nf = struct ("faces", f, "vertices", v, "facevertexcdata", c); | |
142 endif | |
143 case 1 | |
144 nf = struct ("faces", f, "vertices", v, "facevertexcdata", c); | |
145 case 2 | |
146 nf = f; | |
147 nv = v; | |
148 endswitch | |
149 | |
150 endfunction | |
151 | |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
152 |
14574 | 153 %!demo |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
154 %! clf; |
14574 | 155 %! faces = [1 2 3; 1 3 4]; |
156 %! vertices = [0 0; 1 0; 1 1; 0 1]; | |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
157 %! patch ('Faces', faces, 'Vertices', vertices, 'FaceColor', 'none'); |
14574 | 158 %! fv = shrinkfaces (faces, vertices, 0.25); |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
159 %! patch (fv); |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
160 %! axis equal; |
14574 | 161 |
162 %!demo | |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
163 %! clf; |
14574 | 164 %! faces = [1 2 3 4; 5 6 7 8]; |
165 %! vertices = [0 0; 1 0; 2 1; 1 1; 2 0; 3 0; 4 1; 3.5 1]; | |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
166 %! patch ('Faces', faces, 'Vertices', vertices, 'FaceColor', 'none'); |
14574 | 167 %! fv = shrinkfaces (faces, vertices, 0.25); |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
168 %! patch (fv); |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
169 %! axis equal; |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
170 %! grid on; |
14574 | 171 |
172 %!demo | |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
173 %! clf; |
14574 | 174 %! faces = [1 2 3 4]; |
175 %! vertices = [-1 2; 0 0; 1 2; 0 1]; | |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
176 %! patch ('Faces', faces, 'Vertices', vertices, 'FaceColor', 'none'); |
14574 | 177 %! fv = shrinkfaces (faces, vertices, 0.25); |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
178 %! patch (fv); |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
179 %! axis equal; |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
180 %! grid on; |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
181 %! title 'faces which are not convex are clearly not allowed' |
14574 | 182 |
183 %!demo | |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
184 %! clf; |
14574 | 185 %! [phi r] = meshgrid (linspace (0, 1.5*pi, 16), linspace (1, 2, 4)); |
186 %! tri = delaunay (phi(:), r(:)); | |
187 %! v = [r(:).*sin(phi(:)) r(:).*cos(phi(:))]; | |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
188 %! p = patch ('Faces', tri, 'Vertices', v, 'FaceColor', 'none'); |
14574 | 189 %! fv = shrinkfaces (p); |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
190 %! patch (fv); |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
191 %! axis equal; |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
192 %! grid on; |
14574 | 193 |
194 %!demo | |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
195 %! clf; |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
196 %! N = 10; % N intervals per axis |
14574 | 197 %! [x, y, z] = meshgrid (linspace (-4,4,N+1)); |
198 %! val = x.^3 + y.^3 + z.^3; | |
199 %! fv = isosurface (x, y, z, val, 3, z); | |
200 %! | |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
201 %! p = patch ('Faces', fv.faces, 'Vertices', fv.vertices, 'FaceVertexCData', ... |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
202 %! fv.facevertexcdata, 'FaceColor', 'interp', 'EdgeColor', 'black'); |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
203 %! axis equal; |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
204 %! view (115, 30); |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
205 %! drawnow; |
14574 | 206 %! shrinkfaces (p, 0.6); |
207 | |
208 %!shared faces, vertices, nfv, nfv2 | |
209 %! faces = [1 2 3]; | |
210 %! vertices = [0 0 0; 1 0 0; 1 1 0]; | |
211 %! nfv = shrinkfaces (faces, vertices, 0.7); | |
212 %! nfv2 = shrinkfaces (nfv, 1/0.7); | |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
213 %!assert (isfield (nfv, "faces")) |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
214 %!assert (isfield (nfv, "vertices")) |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
215 %!assert (size (nfv.faces), [1 3]) |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
216 %!assert (size (nfv.vertices), [3 3]) |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
217 %!assert (norm (nfv2.vertices - vertices), 0, 2*eps) |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
218 |
18402 | 219 %% Test input validation |
220 %!error shrinkfaces () | |
221 %!error shrinkfaces (1,2,3,4) | |
222 %!error [a,b,c] = shrinkfaces (1) | |
223 %!error <scale factor must be a positive scalar> shrinkfaces (nfv, ones (2)) | |
224 %!error <scale factor must be a positive scalar> shrinkfaces (nfv, 0) | |
225 %!error <only 2-D and 3-D patches are supported> shrinkfaces (faces, ones (3,1)) | |
226 %!error <only 2-D and 3-D patches are supported> shrinkfaces (faces, ones (3,4)) | |
227 %!error <faces must consist of at least 3 vertices> shrinkfaces (faces(1:2), vertices) | |
228 |