Mercurial > hg > octave-nkf
annotate scripts/plot/draw/shrinkfaces.m @ 19898:4197fc428c7d
maint: Update copyright notices for 2015.
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Wed, 11 Feb 2015 14:19:08 -0500 |
parents | 0e1f5a750d00 |
children | 9fc020886ae9 |
rev | line source |
---|---|
19898
4197fc428c7d
maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents:
19793
diff
changeset
|
1 ## Copyright (C) 2012-2015 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 | |
19790
446c46af4b42
strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents:
17744
diff
changeset
|
74 |
14574 | 75 sf = 0.3; |
18402 | 76 colors = []; |
14574 | 77 p = varargin{1}; |
78 | |
19256
d3a223128efc
shrinkfaces.m: improve input checking
Pantxo Diribarne <pantxo.diribarne@gmail.com>
parents:
18402
diff
changeset
|
79 if (isscalar (p) && ishandle (p) && nargin < 3 && |
d3a223128efc
shrinkfaces.m: improve input checking
Pantxo Diribarne <pantxo.diribarne@gmail.com>
parents:
18402
diff
changeset
|
80 strcmp (get (p, "type"), "patch")) |
14574 | 81 faces = get (p, "Faces"); |
82 vertices = get (p, "Vertices"); | |
83 colors = get (p, "FaceVertexCData"); | |
84 if (nargin == 2) | |
85 sf = varargin{2}; | |
86 endif | |
87 elseif (isstruct (p) && nargin < 3) | |
88 faces = p.faces; | |
89 vertices = p.vertices; | |
90 if (isfield (p, "facevertexcdata")) | |
91 colors = p.facevertexcdata; | |
92 endif | |
93 if (nargin == 2) | |
94 sf = varargin{2}; | |
95 endif | |
96 elseif (ismatrix (p) && nargin >= 2 && ismatrix (varargin{2})) | |
97 faces = p; | |
98 vertices = varargin{2}; | |
99 if (nargin == 3) | |
100 sf = varargin{3}; | |
101 endif | |
102 else | |
103 print_usage (); | |
104 endif | |
19790
446c46af4b42
strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents:
17744
diff
changeset
|
105 |
14574 | 106 if (! isscalar (sf) || sf <= 0) |
15202
f3b5cadfd6d5
fix missing semicolons in various .m files
John W. Eaton <jwe@octave.org>
parents:
14872
diff
changeset
|
107 error ("shrinkfaces: scale factor must be a positive scalar"); |
14574 | 108 endif |
109 | |
18402 | 110 nc = columns (vertices); |
111 if (nc < 2 || nc > 3) | |
17122
eaab03308c0b
doc: Rewrite docstrings for most plot functions.
Rik <rik@octave.org>
parents:
16828
diff
changeset
|
112 error ("shrinkfaces: only 2-D and 3-D patches are supported"); |
14574 | 113 endif |
114 | |
14872
c2dbdeaa25df
maint: use rows() and columns() to clarify m-files.
Rik <octave@nomad.inbox5.com>
parents:
14621
diff
changeset
|
115 m = columns (faces); |
14574 | 116 if (m < 3) |
15202
f3b5cadfd6d5
fix missing semicolons in various .m files
John W. Eaton <jwe@octave.org>
parents:
14872
diff
changeset
|
117 error ("shrinkfaces: faces must consist of at least 3 vertices"); |
14574 | 118 endif |
119 | |
120 v = vertices(faces'(:), :); | |
14872
c2dbdeaa25df
maint: use rows() and columns() to clarify m-files.
Rik <octave@nomad.inbox5.com>
parents:
14621
diff
changeset
|
121 if (isempty (colors) || rows (colors) == rows (faces)) |
14574 | 122 c = colors; |
14872
c2dbdeaa25df
maint: use rows() and columns() to clarify m-files.
Rik <octave@nomad.inbox5.com>
parents:
14621
diff
changeset
|
123 elseif (rows (colors) == rows (vertices)) |
14574 | 124 c = colors(faces'(:), :); |
125 else | |
18402 | 126 c = []; # Discard inconsistent color data. |
14574 | 127 endif |
14872
c2dbdeaa25df
maint: use rows() and columns() to clarify m-files.
Rik <octave@nomad.inbox5.com>
parents:
14621
diff
changeset
|
128 sv = rows (v); |
18402 | 129 ## We have to deal with a possibly very large number of vertices, so use |
130 ## sparse as midpoint (1/m, ..., 1/m) in generalized barycentric coordinates. | |
131 midpoints = full (kron (speye (sv / m), ones (m, m) / m) * sparse (v)); | |
14574 | 132 v = sqrt (sf) * (v - midpoints) + midpoints; |
133 f = reshape (1:sv, m, sv / m)'; | |
19790
446c46af4b42
strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents:
17744
diff
changeset
|
134 |
14574 | 135 switch (nargout) |
136 case 0 | |
137 if (ishandle (p)) | |
15202
f3b5cadfd6d5
fix missing semicolons in various .m files
John W. Eaton <jwe@octave.org>
parents:
14872
diff
changeset
|
138 ## avoid exceptions |
f3b5cadfd6d5
fix missing semicolons in various .m files
John W. Eaton <jwe@octave.org>
parents:
14872
diff
changeset
|
139 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
|
140 set (p, "Vertices", v, "Faces", f, "FaceVertexCData", c); |
14574 | 141 else |
142 nf = struct ("faces", f, "vertices", v, "facevertexcdata", c); | |
143 endif | |
144 case 1 | |
145 nf = struct ("faces", f, "vertices", v, "facevertexcdata", c); | |
146 case 2 | |
147 nf = f; | |
148 nv = v; | |
149 endswitch | |
150 | |
151 endfunction | |
152 | |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
153 |
14574 | 154 %!demo |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
155 %! clf; |
14574 | 156 %! faces = [1 2 3; 1 3 4]; |
157 %! 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
|
158 %! patch ('Faces', faces, 'Vertices', vertices, 'FaceColor', 'none'); |
14574 | 159 %! 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
|
160 %! patch (fv); |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
161 %! axis equal; |
14574 | 162 |
163 %!demo | |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
164 %! clf; |
14574 | 165 %! faces = [1 2 3 4; 5 6 7 8]; |
166 %! 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
|
167 %! patch ('Faces', faces, 'Vertices', vertices, 'FaceColor', 'none'); |
14574 | 168 %! 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
|
169 %! patch (fv); |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
170 %! axis equal; |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
171 %! grid on; |
14574 | 172 |
173 %!demo | |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
174 %! clf; |
14574 | 175 %! faces = [1 2 3 4]; |
176 %! 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
|
177 %! patch ('Faces', faces, 'Vertices', vertices, 'FaceColor', 'none'); |
14574 | 178 %! 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
|
179 %! patch (fv); |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
180 %! axis equal; |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
181 %! grid on; |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
182 %! title 'faces which are not convex are clearly not allowed' |
14574 | 183 |
184 %!demo | |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
185 %! clf; |
14574 | 186 %! [phi r] = meshgrid (linspace (0, 1.5*pi, 16), linspace (1, 2, 4)); |
187 %! tri = delaunay (phi(:), r(:)); | |
188 %! 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
|
189 %! p = patch ('Faces', tri, 'Vertices', v, 'FaceColor', 'none'); |
14574 | 190 %! 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
|
191 %! patch (fv); |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
192 %! axis equal; |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
193 %! grid on; |
14574 | 194 |
195 %!demo | |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
196 %! clf; |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
197 %! N = 10; % N intervals per axis |
14574 | 198 %! [x, y, z] = meshgrid (linspace (-4,4,N+1)); |
199 %! val = x.^3 + y.^3 + z.^3; | |
200 %! fv = isosurface (x, y, z, val, 3, z); | |
201 %! | |
16828
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
202 %! 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
|
203 %! 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
|
204 %! axis equal; |
ddac88d32d6a
Make demos in plot m-files compatible with Matlab for running comparison script.
Rik <rik@octave.org>
parents:
15332
diff
changeset
|
205 %! 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
|
206 %! drawnow; |
14574 | 207 %! shrinkfaces (p, 0.6); |
208 | |
209 %!shared faces, vertices, nfv, nfv2 | |
210 %! faces = [1 2 3]; | |
211 %! vertices = [0 0 0; 1 0 0; 1 1 0]; | |
212 %! nfv = shrinkfaces (faces, vertices, 0.7); | |
213 %! 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
|
214 %!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
|
215 %!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
|
216 %!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
|
217 %!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
|
218 %!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
|
219 |
18402 | 220 %% Test input validation |
221 %!error shrinkfaces () | |
222 %!error shrinkfaces (1,2,3,4) | |
223 %!error [a,b,c] = shrinkfaces (1) | |
224 %!error <scale factor must be a positive scalar> shrinkfaces (nfv, ones (2)) | |
225 %!error <scale factor must be a positive scalar> shrinkfaces (nfv, 0) | |
226 %!error <only 2-D and 3-D patches are supported> shrinkfaces (faces, ones (3,1)) | |
227 %!error <only 2-D and 3-D patches are supported> shrinkfaces (faces, ones (3,4)) | |
228 %!error <faces must consist of at least 3 vertices> shrinkfaces (faces(1:2), vertices) | |
229 |