view scripts/sparse/treelayout.m @ 30875:5d3faba0342e

doc: Ensure documentation lists output argument when it exists for all m-files. For new users of Octave it is best to show explicit calling forms in the documentation and to show a return argument when it exists. * bp-table.cc, shift.m, accumarray.m, accumdim.m, bincoeff.m, bitcmp.m, bitget.m, bitset.m, blkdiag.m, celldisp.m, cplxpair.m, dblquad.m, flip.m, fliplr.m, flipud.m, idivide.m, int2str.m, interpft.m, logspace.m, num2str.m, polyarea.m, postpad.m, prepad.m, randi.m, repmat.m, rng.m, rot90.m, rotdim.m, structfun.m, triplequad.m, uibuttongroup.m, uicontrol.m, uipanel.m, uipushtool.m, uitoggletool.m, uitoolbar.m, waitforbuttonpress.m, help.m, __additional_help_message__.m, hsv.m, im2double.m, im2frame.m, javachk.m, usejava.m, argnames.m, char.m, formula.m, inline.m, __vectorize__.m, findstr.m, flipdim.m, strmatch.m, vectorize.m, commutation_matrix.m, cond.m, cross.m, duplication_matrix.m, expm.m, orth.m, rank.m, rref.m, trace.m, vech.m, cast.m, compare_versions.m, delete.m, dir.m, fileattrib.m, grabcode.m, gunzip.m, inputname.m, license.m, list_primes.m, ls.m, mexext.m, movefile.m, namelengthmax.m, nargoutchk.m, nthargout.m, substruct.m, swapbytes.m, ver.m, verLessThan.m, what.m, fminunc.m, fsolve.m, fzero.m, optimget.m, __fdjac__.m, matlabroot.m, savepath.m, campos.m, camroll.m, camtarget.m, camup.m, camva.m, camzoom.m, clabel.m, diffuse.m, legend.m, orient.m, rticks.m, specular.m, thetaticks.m, xlim.m, xtickangle.m, xticklabels.m, xticks.m, ylim.m, ytickangle.m, yticklabels.m, yticks.m, zlim.m, ztickangle.m, zticklabels.m, zticks.m, ellipsoid.m, isocolors.m, isonormals.m, stairs.m, surfnorm.m, __actual_axis_position__.m, __pltopt__.m, close.m, graphics_toolkit.m, pan.m, print.m, printd.m, __ghostscript__.m, __gnuplot_print__.m, __opengl_print__.m, rotate3d.m, subplot.m, zoom.m, compan.m, conv.m, poly.m, polyaffine.m, polyder.m, polyint.m, polyout.m, polyreduce.m, polyvalm.m, roots.m, prefdir.m, prefsfile.m, profexplore.m, profexport.m, profshow.m, powerset.m, unique.m, arch_rnd.m, arma_rnd.m, autoreg_matrix.m, bartlett.m, blackman.m, detrend.m, durbinlevinson.m, fftconv.m, fftfilt.m, fftshift.m, fractdiff.m, hamming.m, hanning.m, hurst.m, ifftshift.m, rectangle_lw.m, rectangle_sw.m, triangle_lw.m, sinc.m, sinetone.m, sinewave.m, spectral_adf.m, spectral_xdf.m, spencer.m, ilu.m, __sprand__.m, sprand.m, sprandn.m, sprandsym.m, treelayout.m, beta.m, betainc.m, betaincinv.m, betaln.m, cosint.m, expint.m, factorial.m, gammainc.m, gammaincinv.m, lcm.m, nthroot.m, perms.m, reallog.m, realpow.m, realsqrt.m, sinint.m, hadamard.m, hankel.m, hilb.m, invhilb.m, magic.m, pascal.m, rosser.m, toeplitz.m, vander.m, wilkinson.m, center.m, corr.m, cov.m, discrete_cdf.m, discrete_inv.m, discrete_pdf.m, discrete_rnd.m, empirical_cdf.m, empirical_inv.m, empirical_pdf.m, empirical_rnd.m, kendall.m, kurtosis.m, mad.m, mean.m, meansq.m, median.m, mode.m, moment.m, range.m, ranks.m, run_count.m, skewness.m, spearman.m, statistics.m, std.m, base2dec.m, bin2dec.m, blanks.m, cstrcat.m, deblank.m, dec2base.m, dec2bin.m, dec2hex.m, hex2dec.m, index.m, regexptranslate.m, rindex.m, strcat.m, strjust.m, strtrim.m, strtrunc.m, substr.m, untabify.m, __have_feature__.m, __prog_output_assert__.m, __run_test_suite__.m, example.m, fail.m, asctime.m, calendar.m, ctime.m, date.m, etime.m: Add return arguments to @deftypefn macros where they were missing. Rename variables in functions (particularly generic "retval") to match documentation. Rename some return variables for (hopefully) better clarity (e.g., 'ax' to 'hax' to indicate it is a graphics handle to an axes object).
author Rik <rik@octave.org>
date Wed, 30 Mar 2022 20:40:27 -0700
parents 397d29f7135c
children 597f3ee61a48
line wrap: on
line source

########################################################################
##
## Copyright (C) 2008-2022 The Octave Project Developers
##
## See the file COPYRIGHT.md in the top-level directory of this
## distribution or <https://octave.org/copyright/>.
##
## This file is part of Octave.
##
## Octave is free software: you can redistribute it and/or modify it
## under the terms of the GNU General Public License as published by
## the Free Software Foundation, either version 3 of the License, or
## (at your option) any later version.
##
## Octave is distributed in the hope that it will be useful, but
## WITHOUT ANY WARRANTY; without even the implied warranty of
## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
## GNU General Public License for more details.
##
## You should have received a copy of the GNU General Public License
## along with Octave; see the file COPYING.  If not, see
## <https://www.gnu.org/licenses/>.
##
########################################################################

## -*- texinfo -*-
## @deftypefn  {} {[@var{x}, @var{y}] =} treelayout (@var{tree})
## @deftypefnx {} {[@var{x}, @var{y}] =} treelayout (@var{tree}, @var{permutation})
## @deftypefnx {} {[@var{x}, @var{y}, @var{h}, @var{s}] =} treelayout (@dots{})
## treelayout lays out a tree or a forest.
##
## The first argument @var{tree} is a vector of predecessors.
##
## The optional parameter @var{permutation} is a postorder permutation.
##
## The complexity of the algorithm is O(n) in terms of time and memory
## requirements.
## @seealso{etreeplot, gplot, treeplot}
## @end deftypefn

function [x, y, h, s] = treelayout (tree, permutation)

  if (nargin < 1)
    print_usage ();
  elseif (! isvector (tree) || rows (tree) != 1 || ! isnumeric (tree)
          || any (tree > length (tree)) || any (tree < 0))
    error ("treelayout: the first input argument must be a vector of predecessors");
  endif

  ## Make it a row vector.
  tree = tree(:)';

  ## The count of nodes of the graph.
  num_nodes = length (tree);
  ## The number of children.
  num_children = zeros (1, num_nodes + 1);

  ## Checking vector of predecessors.
  for i = 1 : num_nodes
    if (tree(i) < i)
      ## This part of graph was checked before.
      continue;
    endif

    ## Try to find cicle in this part of graph using modified Floyd's
    ## cycle-finding algorithm.
    tortoise = tree(i);
    hare = tree(tortoise);

    while (tortoise != hare)
      ## End after finding a cicle or reaching a checked part of graph.

      if (hare < i)
        ## This part of graph was checked before.
        break;
      endif

      tortoise = tree(tortoise);
      ## Hare will move faster than tortoise so in cicle hare must
      ## reach tortoise.
      hare = tree(tree(hare));

    endwhile

    if (tortoise == hare)
      ## If hare reach tortoise we found circle.
      error ("treelayout: vector of predecessors has bad format");
    endif

  endfor
  ## Vector of predecessors has right format.

  for i = 1:num_nodes
    ## vec_of_child is helping vector which is used to speed up the
    ## choice of descendant nodes.

    num_children(tree(i)+1) = num_children(tree(i)+1) + 1;
  endfor

  pos = 1;
  start = zeros (1, num_nodes+1);
  xhelp = zeros (1, num_nodes+1);
  stop = zeros (1, num_nodes+1);
  for i = 1 : num_nodes + 1
    start(i) = pos;
    xhelp(i) = pos;
    pos += num_children(i);
    stop(i) = pos;
  endfor

  if (nargin == 1)
    for i = 1:num_nodes
      vec_of_child(xhelp(tree(i)+1)) = i;
      xhelp(tree(i)+1) = xhelp(tree(i)+1) + 1;
    endfor
  else
    vec_of_child = permutation;
  endif

  ## The number of "parent" (actual) node (its descendants will be
  ## browse in the next iteration).
  par_number = 0;

  ## The x-coordinate of the left most descendant of "parent node"
  ## this value is increased in each leaf.
  left_most = 0;

  ## The level of "parent" node (root level is num_nodes).
  level = num_nodes;

  ## num_nodes - max_ht is the height of this graph.
  max_ht = num_nodes;

  ## Main stack - each item consists of two numbers - the number of
  ## node and the number it's of parent node on the top of stack
  ## there is "parent node".
  stk = [-1, 0];

  ## Number of vertices s in the top-level separator.
  s = 0;
  ## Flag which says if we are in top level separator.
  top_level = 1;
  ## The top of the stack.
  while (par_number != -1)
    if (start(par_number+1) < stop(par_number+1))
      idx = vec_of_child(start(par_number+1) : stop(par_number+1) - 1);
    else
      idx = zeros (1, 0);
    endif

    ## Add to idx the vector of parent descendants.
    stk = [stk; [idx', ones(fliplr(size(idx))) * par_number]];

    ## We are in top level separator when we have one child and the
    ## flag is 1
    if (columns (idx) == 1 && top_level == 1)
      s += 1;
    else
      ## We aren't in top level separator now.
      top_level = 0;
    endif
    ## If there is not any descendant of "parent node":
    if (stk(end,2) != par_number)
     left_most += 1;
     x_r(par_number) = left_most;
     max_ht = min (max_ht, level);
     if (length (stk) > 1 && find ((circshift (stk,1) - stk) == 0) > 1
         && stk(end,2) != stk(end-1,2))
        ## Return to the nearest branching the position to return
        ## position is the position on the stack, where should be
        ## started further search (there are two nodes which has the
        ## same parent node).

        position = (find ((circshift (stk(:,2), 1) - stk(:,2)) == 0))(end) + 1;
        par_number_vec = stk(position:end,2);

        ## The vector of removed nodes (the content of stack form
        ## position to end).

        level += length (par_number_vec);

        ## The level have to be decreased.

        x_r(par_number_vec) = left_most;
        stk(position:end,:) = [];
      endif

      ## Remove the next node from "searched branch".

      stk(end,:) = [];
      ## Choose new "parent node".
      par_number = stk(end,1);
      ## If there is another branch start to search it.
      if (par_number != -1)
        y(par_number) = level;
        x_l(par_number) = left_most + 1;
      endif
    else

      ## There were descendants of "parent nod" choose the last of
      ## them and go on through it.
      level -= 1;
      par_number = stk(end,1);
      y(par_number) = level;
      x_l(par_number) = left_most + 1;
    endif
  endwhile

  ## Calculate the x coordinates (the known values are the position
  ## of most left and most right descendants).
  x = (x_l + x_r) / 2;

  h = num_nodes - max_ht - 1;

endfunction


%!test
%! % Compute a simple tree layout
%! [x, y, h, s] = treelayout ([0, 1, 2, 2]);
%! assert (x, [1.5, 1.5, 2, 1]);
%! assert (y, [3, 2, 1, 1]);
%! assert (h, 2);
%! assert (s, 2);

%!test
%! % Compute a simple tree layout with defined postorder permutation
%! [x, y, h, s] = treelayout ([0, 1, 2, 2], [1, 2, 4, 3]);
%! assert (x, [1.5, 1.5, 1, 2]);
%! assert (y, [3, 2, 1, 1]);
%! assert (h, 2);
%! assert (s, 2);

%!test
%! % Compute a simple tree layout with defined postorder permutation
%! [x, y, h, s] = treelayout ([0, 1, 2, 2], [4, 2, 3, 1]);
%! assert (x, [0, 0, 0, 1]);
%! assert (y, [0, 0, 0, 3]);
%! assert (h, 0);
%! assert (s, 1);