Mercurial > hg > octave-nkf
annotate scripts/sparse/treelayout.m @ 14868:5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
* lin2mu.m, loadaudio.m, wavread.m, accumarray.m, bicubic.m, celldisp.m,
colon.m, cplxpair.m, dblquad.m, divergence.m, genvarname.m, gradient.m,
int2str.m, interp1.m, interp1q.m, interp2.m, interpn.m, loadobj.m, nthargout.m,
__isequal__.m, __splinen__.m, quadgk.m, quadl.m, quadv.m, rat.m, rot90.m,
rotdim.m, saveobj.m, subsindex.m, triplequad.m, delaunay3.m, griddata.m,
inpolygon.m, tsearchn.m, voronoi.m, get_first_help_sentence.m, which.m,
gray2ind.m, pink.m, dlmwrite.m, strread.m, textread.m, textscan.m, housh.m,
ishermitian.m, issymmetric.m, krylov.m, logm.m, null.m, rref.m,
compare_versions.m, copyfile.m, dump_prefs.m, edit.m, fileparts.m,
getappdata.m, isappdata.m, movefile.m, orderfields.m, parseparams.m,
__xzip__.m, rmappdata.m, setappdata.m, swapbytes.m, unpack.m, ver.m, fminbnd.m,
fminunc.m, fsolve.m, glpk.m, lsqnonneg.m, qp.m, sqp.m, configure_make.m,
copy_files.m, describe.m, get_description.m, get_forge_pkg.m, install.m,
installed_packages.m, is_architecture_dependent.m, load_package_dirs.m,
print_package_description.m, rebuild.m, repackage.m, save_order.m, shell.m,
allchild.m, ancestor.m, area.m, axes.m, axis.m, clabel.m, close.m, colorbar.m,
comet.m, comet3.m, contour.m, cylinder.m, ezmesh.m, ezsurf.m, findobj.m,
fplot.m, hist.m, isocolors.m, isonormals.m, isosurface.m, isprop.m, legend.m,
mesh.m, meshz.m, pareto.m, pcolor.m, peaks.m, plot3.m, plotmatrix.m, plotyy.m,
polar.m, print.m, __add_datasource__.m, __add_default_menu__.m,
__axes_limits__.m, __bar__.m, __clabel__.m, __contour__.m, __errcomm__.m,
__errplot__.m, __ezplot__.m, __file_filter__.m, __fltk_print__.m,
__ghostscript__.m, __gnuplot_print__.m, __go_draw_axes__.m,
__go_draw_figure__.m, __interp_cube__.m, __marching_cube__.m, __patch__.m,
__pie__.m, __plt__.m, __print_parse_opts__.m, __quiver__.m, __scatter__.m,
__stem__.m, __tight_eps_bbox__.m, __uigetdir_fltk__.m, __uigetfile_fltk__.m,
__uiputfile_fltk__.m, quiver.m, quiver3.m, rectangle.m, refreshdata.m,
ribbon.m, scatter.m, semilogy.m, shading.m, slice.m, subplot.m, surface.m,
surfl.m, surfnorm.m, text.m, uigetfile.m, uiputfile.m, whitebg.m, deconv.m,
mkpp.m, pchip.m, polyaffine.m, polyder.m, polygcd.m, polyout.m, polyval.m,
ppint.m, ppjumps.m, ppval.m, residue.m, roots.m, spline.m, splinefit.m,
addpref.m, getpref.m, setpref.m, ismember.m, setxor.m, arch_fit.m, arch_rnd.m,
arch_test.m, autoreg_matrix.m, diffpara.m, fftconv.m, filter2.m, hanning.m,
hurst.m, periodogram.m, triangle_sw.m, sinc.m, spectral_xdf.m, spencer.m,
stft.m, synthesis.m, unwrap.m, yulewalker.m, bicgstab.m, gmres.m, pcg.m, pcr.m,
__sprand_impl__.m, speye.m, spfun.m, sprandn.m, spstats.m, svds.m,
treelayout.m, treeplot.m, bessel.m, factor.m, legendre.m, perms.m, primes.m,
magic.m, toeplitz.m, corr.m, cov.m, mean.m, median.m, mode.m, qqplot.m,
quantile.m, ranks.m, zscore.m, logistic_regression_likelihood.m,
bartlett_test.m, chisquare_test_homogeneity.m, chisquare_test_independence.m,
kolmogorov_smirnov_test.m, run_test.m, u_test.m, wilcoxon_test.m, z_test.m,
z_test_2.m, bin2dec.m, dec2base.m, mat2str.m, strcat.m, strchr.m, strjust.m,
strtok.m, substr.m, untabify.m, assert.m, demo.m, example.m, fail.m, speed.m,
test.m, now.m: Use Octave coding conventions for cuddling parentheses in
scripts directory.
author | Rik <octave@nomad.inbox5.com> |
---|---|
date | Tue, 17 Jul 2012 07:08:39 -0700 |
parents | f3d52523cde1 |
children | b81b9d079515 |
rev | line source |
---|---|
14138
72c96de7a403
maint: update copyright notices for 2012
John W. Eaton <jwe@octave.org>
parents:
13060
diff
changeset
|
1 ## Copyright (C) 2008-2012 Ivana Varekova & Radek Salac |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
2 ## |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
3 ## This file is part of Octave. |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
4 ## |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
5 ## Octave is free software; you can redistribute it and/or modify it |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
6 ## under the terms of the GNU General Public License as published by |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
7 ## the Free Software Foundation; either version 3 of the License, or (at |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
8 ## your option) any later version. |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
9 ## |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
10 ## Octave is distributed in the hope that it will be useful, but |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
11 ## WITHOUT ANY WARRANTY; without even the implied warranty of |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
12 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
13 ## General Public License for more details. |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
14 ## |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
15 ## You should have received a copy of the GNU General Public License |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
16 ## along with Octave; see the file COPYING. If not, see |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
17 ## <http://www.gnu.org/licenses/>. |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
18 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
19 ## -*- texinfo -*- |
11469
c776f063fefe
Overhaul m-script files to use common variable name between code and documentation.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
20 ## @deftypefn {Function File} {} treelayout (@var{tree}) |
c776f063fefe
Overhaul m-script files to use common variable name between code and documentation.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
21 ## @deftypefnx {Function File} {} treelayout (@var{tree}, @var{permutation}) |
c776f063fefe
Overhaul m-script files to use common variable name between code and documentation.
Rik <octave@nomad.inbox5.com>
parents:
10793
diff
changeset
|
22 ## treelayout lays out a tree or a forest. The first argument @var{tree} is a |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
23 ## vector of |
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
24 ## predecessors, optional parameter @var{permutation} is an optional postorder |
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
25 ## permutation. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
26 ## The complexity of the algorithm is O(n) in |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
27 ## terms of time and memory requirements. |
12575
d0b799dafede
Grammarcheck files for 3.4.1 release.
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
28 ## @seealso{etreeplot, gplot, treeplot} |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
29 ## @end deftypefn |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
30 |
8507 | 31 function [x_coordinate, y_coordinate, height, s] = treelayout (tree, permutation) |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
32 if (nargin < 1 || nargin > 2 || nargout > 4) |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
33 print_usage (); |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
34 elseif (! isvector (tree) || rows (tree) != 1 || ! isnumeric (tree) |
8507 | 35 || any (tree > length (tree)) || any (tree < 0)) |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
36 error ("treelayout: the first input argument must be a vector of predecessors"); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
37 else |
8506 | 38 ## Make it a row vector. |
8507 | 39 tree = tree(:)'; |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
40 |
8506 | 41 ## The count of nodes of the graph. |
8507 | 42 num_nodes = length (tree); |
8506 | 43 ## The number of children. |
8507 | 44 num_children = zeros (1, num_nodes + 1); |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
45 |
8506 | 46 ## Checking vector of predecessors. |
8507 | 47 for i = 1 : num_nodes |
48 if (tree(i) < i) | |
10549 | 49 ## This part of graph was checked before. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
50 continue; |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
51 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
52 |
8506 | 53 ## Try to find cicle in this part of graph using modified Floyd's |
54 ## cycle-finding algorithm. | |
8507 | 55 tortoise = tree(i); |
56 hare = tree(tortoise); | |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
57 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
58 while (tortoise != hare) |
10549 | 59 ## End after finding a cicle or reaching a checked part of graph. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
60 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
61 if (hare < i) |
8506 | 62 ## This part of graph was checked before. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
63 break |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
64 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
65 |
8507 | 66 tortoise = tree(tortoise); |
10549 | 67 ## Hare will move faster than tortoise so in cicle hare must |
68 ## reach tortoise. | |
8507 | 69 hare = tree(tree(hare)); |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
70 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
71 endwhile |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
72 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
73 if (tortoise == hare) |
10549 | 74 ## If hare reach tortoise we found circle. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
75 error ("treelayout: vector of predecessors has bad format"); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
76 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
77 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
78 endfor |
8506 | 79 ## Vector of predecessors has right format. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
80 |
8507 | 81 for i = 1:num_nodes |
82 ## vec_of_child is helping vector which is used to speed up the | |
8506 | 83 ## choice of descendant nodes. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
84 |
8507 | 85 num_children(tree(i)+1) = num_children(tree(i)+1) + 1; |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
86 endfor |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
87 |
8507 | 88 pos = 1; |
89 start = zeros (1, num_nodes+1); | |
90 xhelp = zeros (1, num_nodes+1); | |
91 stop = zeros (1, num_nodes+1); | |
92 for i = 1 : num_nodes + 1 | |
93 start(i) = pos; | |
94 xhelp(i) = pos; | |
95 pos += num_children(i); | |
96 stop(i) = pos; | |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
97 endfor |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
98 |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
99 if (nargin == 1) |
8507 | 100 for i = 1:num_nodes |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
101 vec_of_child(xhelp(tree(i)+1)) = i; |
8507 | 102 xhelp(tree(i)+1) = xhelp(tree(i)+1) + 1; |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
103 endfor |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
104 else |
8507 | 105 vec_of_child = permutation; |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
106 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
107 |
8506 | 108 ## The number of "parent" (actual) node (it's descendants will be |
109 ## browse in the next iteration). | |
8507 | 110 par_number = 0; |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
111 |
8506 | 112 ## The x-coordinate of the left most descendant of "parent node" |
113 ## this value is increased in each leaf. | |
8507 | 114 left_most = 0; |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
115 |
8507 | 116 ## The level of "parent" node (root level is num_nodes). |
117 level = num_nodes; | |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
118 |
8507 | 119 ## num_nodes - max_ht is the height of this graph. |
120 max_ht = num_nodes; | |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
121 |
8506 | 122 ## Main stack - each item consists of two numbers - the number of |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
123 ## node and the number it's of parent node on the top of stack |
8506 | 124 ## there is "parent node". |
8507 | 125 stk = [-1, 0]; |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
126 |
8506 | 127 ## Number of vertices s in the top-level separator. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
128 s = 0; |
8506 | 129 ## Flag which says if we are in top level separator. |
8507 | 130 top_level = 1; |
8506 | 131 ## The top of the stack. |
8507 | 132 while (par_number != -1) |
133 if (start(par_number+1) < stop(par_number+1)) | |
134 idx = vec_of_child(start(par_number+1) : stop(par_number+1) - 1); | |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
135 else |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
136 idx = zeros (1, 0); |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
137 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
138 |
8506 | 139 ## Add to idx the vector of parent descendants. |
8507 | 140 stk = [stk; [idx', ones(fliplr(size(idx))) * par_number]]; |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
141 |
8506 | 142 ## We are in top level separator when we have one child and the |
143 ## flag is 1 | |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
144 if (columns (idx) == 1 && top_level == 1) |
8507 | 145 s++; |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
146 else |
8506 | 147 # We aren't in top level separator now. |
8507 | 148 top_level = 0; |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
149 endif |
8506 | 150 ## If there is not any descendant of "parent node": |
8507 | 151 if (stk(end,2) != par_number) |
152 left_most++; | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
153 x_coordinate_r(par_number) = left_most; |
8507 | 154 max_ht = min (max_ht, level); |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
155 if (length (stk) > 1 && find ((shift (stk,1) - stk) == 0) > 1 |
10549 | 156 && stk(end,2) != stk(end-1,2)) |
157 ## Return to the nearest branching the position to return | |
158 ## position is the position on the stack, where should be | |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
159 ## started further search (there are two nodes which has the |
8506 | 160 ## same parent node). |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
161 |
8507 | 162 position = (find ((shift (stk(:,2), 1) - stk(:,2)) == 0))(end) + 1; |
163 par_number_vec = stk(position:end,2); | |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
164 |
8506 | 165 ## The vector of removed nodes (the content of stack form |
166 ## position to end). | |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
167 |
8507 | 168 level += length (par_number_vec); |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
169 |
10549 | 170 ## The level have to be decreased. |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
171 |
8507 | 172 x_coordinate_r(par_number_vec) = left_most; |
173 stk(position:end,:) = []; | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
174 endif |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
175 |
8506 | 176 ## Remove the next node from "searched branch". |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
177 |
8507 | 178 stk(end,:) = []; |
10549 | 179 ## Choose new "parent node". |
8507 | 180 par_number = stk(end,1); |
10549 | 181 ## If there is another branch start to search it. |
182 if (par_number != -1) | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
183 y_coordinate(par_number) = level; |
8507 | 184 x_coordinate_l(par_number) = left_most + 1; |
10549 | 185 endif |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
186 else |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
187 |
8506 | 188 ## There were descendants of "parent nod" choose the last of |
189 ## them and go on through it. | |
8507 | 190 level--; |
191 par_number = stk(end,1); | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
192 y_coordinate(par_number) = level; |
8507 | 193 x_coordinate_l(par_number) = left_most + 1; |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
194 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
195 endwhile |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
196 |
8506 | 197 ## Calculate the x coordinates (the known values are the position |
198 ## of most left and most right descendants). | |
8507 | 199 x_coordinate = (x_coordinate_l + x_coordinate_r) / 2; |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
200 |
8507 | 201 height = num_nodes - max_ht - 1; |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
202 endif |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
203 endfunction |
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
204 |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
205 |
13060
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
206 %!test |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
207 %! % Compute a simple tree layout |
13060
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
208 %! [x, y, h, s] = treelayout ([0, 1, 2, 2]); |
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
209 %! assert (x, [1.5, 1.5, 2, 1]); |
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
210 %! assert (y, [3, 2, 1, 1]); |
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
211 %! assert (h, 2); |
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
212 %! assert (s, 2); |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
213 |
13060
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
214 %!test |
8338
a35bf360b919
Add the cgs and treelayout functions
Radek Salac <salac.r@gmail.com>
parents:
diff
changeset
|
215 %! % Compute a simple tree layout with defined postorder permutation |
13060
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
216 %! [x, y, h, s] = treelayout ([0, 1, 2, 2], [1, 2, 4, 3]); |
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
217 %! assert (x, [1.5, 1.5, 1, 2]); |
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
218 %! assert (y, [3, 2, 1, 1]); |
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
219 %! assert (h, 2); |
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
220 %! assert (s, 2); |
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
221 |
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
222 %!test |
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
223 %! % Compute a simple tree layout with defined postorder permutation |
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
224 %! [x, y, h, s] = treelayout ([0, 1, 2, 2], [4, 2, 3, 1]); |
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
225 %! assert (x, [0, 0, 0, 1]); |
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
226 %! assert (y, [0, 0, 0, 3]); |
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
227 %! assert (h, 0); |
85dd509673e7
codesprint: tests for treelayout
John W. Eaton <jwe@octave.org>
parents:
12575
diff
changeset
|
228 %! assert (s, 1); |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
229 |