igraph Reference Manual

For using the igraph C library

Search the manual:

Chapter 20. Generating layouts for graph drawing

1. 2D layout generators

Layout generator functions (or at least most of them) try to place the vertices and edges of a graph on a 2D plane or in 3D space in a way which visually pleases the human eye.

They take a graph object and a number of parameters as arguments and return an igraph_matrix_t, in which each row gives the coordinates of a vertex.

1.1. igraph_layout_random — Places the vertices uniform randomly on a plane.

igraph_error_t igraph_layout_random(const igraph_t *graph, igraph_matrix_t *res);

Arguments: 

graph:

Pointer to an initialized graph object.

res:

Pointer to an initialized matrix object. This will contain the result and will be resized as needed.

Returns: 

Error code. The current implementation always returns with success.

Time complexity: O(|V|), the number of vertices.

1.2. igraph_layout_circle — Places the vertices uniformly on a circle in arbitrary order.

igraph_error_t igraph_layout_circle(const igraph_t *graph, igraph_matrix_t *res,
                         igraph_vs_t order);

Arguments: 

graph:

Pointer to an initialized graph object.

res:

Pointer to an initialized matrix object. This will contain the result and will be resized as needed.

order:

The order of the vertices on the circle. The vertices not included here, will be placed at (0,0). Supply igraph_vss_all() here to place vertices in the order of their vertex IDs.

Returns: 

Error code.

Time complexity: O(|V|), the number of vertices.

1.3. igraph_layout_star — Generates a star-like layout.

igraph_error_t igraph_layout_star(const igraph_t *graph, igraph_matrix_t *res,
                       igraph_integer_t center, const igraph_vector_int_t *order);

Arguments: 

graph:

The input graph. Its edges are ignored by this function.

res:

Pointer to an initialized matrix object. This will contain the result and will be resized as needed.

center:

The id of the vertex to put in the center. You can set it to any arbitrary value for the special case when the input graph has no vertices; otherwise it must be between 0 and the number of vertices minus 1.

order:

A numeric vector giving the order of the vertices (including the center vertex!). If a null pointer, then the vertices are placed in increasing vertex ID order.

Returns: 

Error code.

Time complexity: O(|V|), linear in the number of vertices.

See also: 

igraph_layout_circle() and other layout generators.

1.4. igraph_layout_grid — Places the vertices on a regular grid on the plane.

igraph_error_t igraph_layout_grid(const igraph_t *graph, igraph_matrix_t *res, igraph_integer_t width);

Arguments: 

graph:

Pointer to an initialized graph object.

res:

Pointer to an initialized matrix object. This will contain the result and will be resized as needed.

width:

The number of vertices in a single row of the grid. When zero or negative, the width of the grid will be the square root of the number of vertices, rounded up if needed.

Returns: 

Error code. The current implementation always returns with success.

Time complexity: O(|V|), the number of vertices.

1.5. igraph_layout_graphopt — Optimizes vertex layout via the graphopt algorithm.

igraph_error_t igraph_layout_graphopt(const igraph_t *graph, igraph_matrix_t *res,
                           igraph_integer_t niter,
                           igraph_real_t node_charge, igraph_real_t node_mass,
                           igraph_real_t spring_length,
                           igraph_real_t spring_constant,
                           igraph_real_t max_sa_movement,
                           igraph_bool_t use_seed);

This is a port of the graphopt layout algorithm by Michael Schmuhl. graphopt version 0.4.1 was rewritten in C, the support for layers was removed and the code was reorganized to avoid some unnecessary steps when the node charge (see below) is zero.

Graphopt uses physical analogies for defining attracting and repelling forces among the vertices and then the physical system is simulated until it reaches an equilibrium. (There is no simulated annealing or anything like that, so a stable fixed point is not guaranteed.)

See also http://www.schmuhl.org/graphopt/ for the original graphopt.

Arguments: 

graph:

The input graph.

res:

Pointer to an initialized matrix, the result will be stored here and its initial contents are used as the starting point of the simulation if the use_seed argument is true. Note that in this case the matrix should have the proper size, otherwise a warning is issued and the supplied values are ignored. If no starting positions are given (or they are invalid) then a random starting position is used. The matrix will be resized if needed.

niter:

Integer constant, the number of iterations to perform. Should be a couple of hundred in general. If you have a large graph then you might want to only do a few iterations and then check the result. If it is not good enough you can feed it in again in the res argument. The original graphopt default is 500.

node_charge:

The charge of the vertices, used to calculate electric repulsion. The original graphopt default is 0.001.

node_mass:

The mass of the vertices, used for the spring forces. The original graphopt defaults to 30.

spring_length:

The length of the springs. The original graphopt defaults to zero.

spring_constant:

The spring constant, the original graphopt defaults to one.

max_sa_movement:

Real constant, it gives the maximum amount of movement allowed in a single step along a single axis. The original graphopt default is 5.

use_seed:

Logical scalar, whether to use the positions in res as a starting configuration. See also res above.

Returns: 

Error code.

Time complexity: O(n (|V|^2+|E|) ), n is the number of iterations, |V| is the number of vertices, |E| the number of edges. If node_charge is zero then it is only O(n|E|).

1.6. igraph_layout_bipartite — Simple layout for bipartite graphs.

igraph_error_t igraph_layout_bipartite(const igraph_t *graph,
                            const igraph_vector_bool_t *types,
                            igraph_matrix_t *res, igraph_real_t hgap,
                            igraph_real_t vgap, igraph_integer_t maxiter);

The layout is created by first placing the vertices in two rows, according to their types. Then the positions within the rows are optimized to minimize edge crossings, by calling igraph_layout_sugiyama().

Arguments: 

graph:

The input graph.

types:

A boolean vector containing ones and zeros, the vertex types. Its length must match the number of vertices in the graph.

res:

Pointer to an initialized matrix, the result, the x and y coordinates are stored here.

hgap:

The preferred minimum horizontal gap between vertices in the same layer (i.e. vertices of the same type).

vgap:

The distance between layers.

maxiter:

Maximum number of iterations in the crossing minimization stage. 100 is a reasonable default; if you feel that you have too many edge crossings, increase this.

Returns: 

Error code.

See also: 

1.7. The DrL layout generator

DrL is a sophisticated layout generator developed and implemented by Shawn Martin et al. As of October 2012 the original DrL homepage is unfortunately not available. You can read more about this algorithm in the following technical report: Martin, S., Brown, W.M., Klavans, R., Boyack, K.W., DrL: Distributed Recursive (Graph) Layout. SAND Reports, 2008. 2936: p. 1-10.

Only a subset of the complete DrL functionality is included in igraph, parallel runs and recursive, multi-level layouting is not supported.

The parameters of the layout are stored in an igraph_layout_drl_options_t structure, this can be initialized by calling the function igraph_layout_drl_options_init(). The fields of this structure can then be adjusted by hand if needed. The layout is calculated by an igraph_layout_drl() call.

1.7.1. igraph_layout_drl_options_t — Parameters for the DrL layout generator

typedef struct igraph_layout_drl_options_t {
    igraph_real_t    edge_cut;
    igraph_integer_t init_iterations;
    igraph_real_t    init_temperature;
    igraph_real_t    init_attraction;
    igraph_real_t    init_damping_mult;
    igraph_integer_t liquid_iterations;
    igraph_real_t    liquid_temperature;
    igraph_real_t    liquid_attraction;
    igraph_real_t    liquid_damping_mult;
    igraph_integer_t expansion_iterations;
    igraph_real_t    expansion_temperature;
    igraph_real_t    expansion_attraction;
    igraph_real_t    expansion_damping_mult;
    igraph_integer_t cooldown_iterations;
    igraph_real_t    cooldown_temperature;
    igraph_real_t    cooldown_attraction;
    igraph_real_t    cooldown_damping_mult;
    igraph_integer_t crunch_iterations;
    igraph_real_t    crunch_temperature;
    igraph_real_t    crunch_attraction;
    igraph_real_t    crunch_damping_mult;
    igraph_integer_t simmer_iterations;
    igraph_real_t    simmer_temperature;
    igraph_real_t    simmer_attraction;
    igraph_real_t    simmer_damping_mult;
} igraph_layout_drl_options_t;

Values: 

edge_cut:

The edge cutting parameter. Edge cutting is done in the late stages of the algorithm in order to achieve less dense layouts. Edges are cut if there is a lot of stress on them (a large value in the objective function sum). The edge cutting parameter is a value between 0 and 1 with 0 representing no edge cutting and 1 representing maximal edge cutting. The default value is 32/40.

init_iterations:

Number of iterations, initial phase.

init_temperature:

Start temperature, initial phase.

init_attraction:

Attraction, initial phase.

init_damping_mult:

Damping factor, initial phase.

liquid_iterations:

Number of iterations in the liquid phase.

liquid_temperature:

Start temperature in the liquid phase.

liquid_attraction:

Attraction in the liquid phase.

liquid_damping_mult:

Multiplicatie damping factor, liquid phase.

expansion_iterations:

Number of iterations in the expansion phase.

expansion_temperature:

Start temperature in the expansion phase.

expansion_attraction:

Attraction, expansion phase.

expansion_damping_mult:

Damping factor, expansion phase.

cooldown_iterations:

Number of iterations in the cooldown phase.

cooldown_temperature:

Start temperature in the cooldown phase.

cooldown_attraction:

Attraction in the cooldown phase.

cooldown_damping_mult:

Damping fact int the cooldown phase.

crunch_iterations:

Number of iterations in the crunch phase.

crunch_temperature:

Start temperature in the crunch phase.

crunch_attraction:

Attraction in the crunch phase.

crunch_damping_mult:

Damping factor in the crunch phase.

simmer_iterations:

Number of iterations in the simmer phase.

simmer_temperature:

Start temperature in te simmer phase.

simmer_attraction:

Attraction in the simmer phase.

simmer_damping_mult:

Multiplicative damping factor in the simmer phase.

1.7.2. igraph_layout_drl_default_t — Predefined parameter templates for the DrL layout generator

typedef enum { IGRAPH_LAYOUT_DRL_DEFAULT = 0,
               IGRAPH_LAYOUT_DRL_COARSEN,
               IGRAPH_LAYOUT_DRL_COARSEST,
               IGRAPH_LAYOUT_DRL_REFINE,
               IGRAPH_LAYOUT_DRL_FINAL
             } igraph_layout_drl_default_t;

These constants can be used to initialize a set of DrL parameters. These can then be modified according to the user's needs.

Values: 

IGRAPH_LAYOUT_DRL_DEFAULT:

The deafult parameters.

IGRAPH_LAYOUT_DRL_COARSEN:

Slightly modified parameters to get a coarser layout.

IGRAPH_LAYOUT_DRL_COARSEST:

An even coarser layout.

IGRAPH_LAYOUT_DRL_REFINE:

Refine an already calculated layout.

IGRAPH_LAYOUT_DRL_FINAL:

Finalize an already refined layout.

1.7.3. igraph_layout_drl_options_init — Initialize parameters for the DrL layout generator

igraph_error_t igraph_layout_drl_options_init(igraph_layout_drl_options_t *options,
                                   igraph_layout_drl_default_t templ);

This function can be used to initialize the struct holding the parameters for the DrL layout generator. There are a number of predefined templates available, it is a good idea to start from one of these by modifying some parameters.

Arguments: 

options:

The struct to initialize.

templ:

The template to use. Currently the following templates are supplied: IGRAPH_LAYOUT_DRL_DEFAULT, IGRAPH_LAYOUT_DRL_COARSEN, IGRAPH_LAYOUT_DRL_COARSEST, IGRAPH_LAYOUT_DRL_REFINE and IGRAPH_LAYOUT_DRL_FINAL.

Returns: 

Error code.

Time complexity: O(1).

1.7.4. igraph_layout_drl — The DrL layout generator

igraph_error_t igraph_layout_drl(const igraph_t *graph, igraph_matrix_t *res,
                      igraph_bool_t use_seed,
                      const igraph_layout_drl_options_t *options,
                      const igraph_vector_t *weights);

This function implements the force-directed DrL layout generator. Please see more in the following technical report: Martin, S., Brown, W.M., Klavans, R., Boyack, K.W., DrL: Distributed Recursive (Graph) Layout. SAND Reports, 2008. 2936: p. 1-10.

Arguments: 

graph:

The input graph.

use_seed:

Logical scalar, if true, then the coordinates supplied in the res argument are used as starting points.

res:

Pointer to a matrix, the result layout is stored here. It will be resized as needed.

options:

The parameters to pass to the layout generator.

weights:

Edge weights, pointer to a vector. If this is a null pointer then every edge will have the same weight.

Returns: 

Error code.

Time complexity: ???.

1.7.5. igraph_layout_drl_3d — The DrL layout generator, 3d version.

igraph_error_t igraph_layout_drl_3d(const igraph_t *graph, igraph_matrix_t *res,
                         igraph_bool_t use_seed,
                         const igraph_layout_drl_options_t *options,
                         const igraph_vector_t *weights);

This function implements the force-directed DrL layout generator. Please see more in the technical report: Martin, S., Brown, W.M., Klavans, R., Boyack, K.W., DrL: Distributed Recursive (Graph) Layout. SAND Reports, 2008. 2936: p. 1-10.

This function uses a modified DrL generator that does the layout in three dimensions.

Arguments: 

graph:

The input graph.

use_seed:

Logical scalar, if true, then the coordinates supplied in the res argument are used as starting points.

res:

Pointer to a matrix, the result layout is stored here. It will be resized as needed.

options:

The parameters to pass to the layout generator.

weights:

Edge weights, pointer to a vector. If this is a null pointer then every edge will have the same weight.

Returns: 

Error code.

Time complexity: ???.

See also: 

igraph_layout_drl() for the standard 2d version.

1.8. igraph_layout_fruchterman_reingold — Places the vertices on a plane according to the Fruchterman-Reingold algorithm.

igraph_error_t igraph_layout_fruchterman_reingold(const igraph_t *graph,
                                       igraph_matrix_t *res,
                                       igraph_bool_t use_seed,
                                       igraph_integer_t niter,
                                       igraph_real_t start_temp,
                                       igraph_layout_grid_t grid,
                                       const igraph_vector_t *weights,
                                       const igraph_vector_t *minx,
                                       const igraph_vector_t *maxx,
                                       const igraph_vector_t *miny,
                                       const igraph_vector_t *maxy);

This is a force-directed layout that simulates an attractive force f_a between connected vertex pairs and a repulsive force f_r between all vertex pairs. The forces are computed as a function of the distance d between the two vertices as

f_a(d) = -w * d^2 and f_r(d) = 1/d,

where w represents the edge weight. The equilibrium distance of two connected vertices is thus 1/w^3, assuming no other forces acting on them.

In disconnected graphs, igraph effectively inserts a weak connection of weight n^(-3/2) between all pairs of vertices, where n is the vertex count. This ensures that components are kept near each other.

Reference:

Fruchterman, T.M.J. and Reingold, E.M.: Graph Drawing by Force-directed Placement. Software -- Practice and Experience, 21/11, 1129--1164, 1991. https://doi.org/10.1002/spe.4380211102

Arguments: 

graph:

Pointer to an initialized graph object.

res:

Pointer to an initialized matrix object. This will contain the result and will be resized as needed.

use_seed:

Logical, if true the supplied values in the res argument are used as an initial layout, if false a random initial layout is used.

niter:

The number of iterations to do. A reasonable default value is 500.

start_temp:

Start temperature. This is the maximum amount of movement allowed along one axis, within one step, for a vertex. Currently it is decreased linearly to zero during the iteration.

grid:

Whether to use the (fast but less accurate) grid based version of the algorithm. Possible values: IGRAPH_LAYOUT_GRID, IGRAPH_LAYOUT_NOGRID, IGRAPH_LAYOUT_AUTOGRID. The last one uses the grid based version only for large graphs, currently the ones with more than 1000 vertices.

weights:

Pointer to a vector containing edge weights, the attraction along the edges will be multiplied by these. Weights must be positive. It will be ignored if it is a null-pointer.

minx:

Pointer to a vector, or a NULL pointer. If not a NULL pointer then the vector gives the minimum x coordinate for every vertex.

maxx:

Same as minx, but the maximum x coordinates.

miny:

Pointer to a vector, or a NULL pointer. If not a NULL pointer then the vector gives the minimum y coordinate for every vertex.

maxy:

Same as miny, but the maximum y coordinates.

Returns: 

Error code.

Time complexity: O(|V|^2) in each iteration, |V| is the number of vertices in the graph.

1.9. igraph_layout_kamada_kawai — Places the vertices on a plane according to the Kamada-Kawai algorithm.

igraph_error_t igraph_layout_kamada_kawai(const igraph_t *graph, igraph_matrix_t *res,
                               igraph_bool_t use_seed, igraph_integer_t maxiter,
                               igraph_real_t epsilon, igraph_real_t kkconst,
                               const igraph_vector_t *weights,
                               const igraph_vector_t *minx, const igraph_vector_t *maxx,
                               const igraph_vector_t *miny, const igraph_vector_t *maxy);

This is a force-directed layout. A spring is inserted between all pairs of vertices, both those which are directly connected and those that are not. The unstretched length of springs is chosen based on the undirected graph distance between the corresponding pair of vertices. Thus, in a weighted graph, increasing the weight between two vertices pushes them apart. The Young modulus of springs is inversely proportional to the graph distance, ensuring that springs between far-apart veritces will have a smaller effect on the layout.

Disconnected graphs are handled by assuming that the graph distance between vertices in different components is the same as the graph diameter.

This layout works particularly well for locally connected spatial networks such as lattices.

This layout algorithm is not suitable for large graphs. The memory requirements are of the order O(|V|^2).

Reference:

Kamada, T. and Kawai, S.: An Algorithm for Drawing General Undirected Graphs. Information Processing Letters, 31/1, 7--15, 1989. https://doi.org/10.1016/0020-0190(89)90102-6

Arguments: 

graph:

A graph object.

res:

Pointer to an initialized matrix object. This will contain the result (x-positions in column zero and y-positions in column one) and will be resized if needed.

use_seed:

Boolean, whether to use the values supplied in the res argument as the initial configuration. If zero and there are any limits on the X or Y coordinates, then a random initial configuration is used. Otherwise the vertices are placed on a circle of radius 1 as the initial configuration.

maxiter:

The maximum number of iterations to perform. A reasonable default value is at least ten (or more) times the number of vertices.

epsilon:

Stop the iteration, if the maximum delta value of the algorithm is smaller than still. It is safe to leave it at zero, and then maxiter iterations are performed.

kkconst:

The Kamada-Kawai vertex attraction constant. Typical value: number of vertices.

weights:

Edge weights, larger values will result longer edges. Weights must be positive. Pass NULL to assume unit weights for all edges.

minx:

Pointer to a vector, or a NULL pointer. If not a NULL pointer then the vector gives the minimum x coordinate for every vertex.

maxx:

Same as minx, but the maximum x coordinates.

miny:

Pointer to a vector, or a NULL pointer. If not a NULL pointer then the vector gives the minimum y coordinate for every vertex.

maxy:

Same as miny, but the maximum y coordinates.

Returns: 

Error code.

Time complexity: O(|V|) for each iteration, after an O(|V|^2 log|V|) initialization step. |V| is the number of vertices in the graph.

1.10. igraph_layout_gem — Layout graph according to GEM algorithm.

igraph_error_t igraph_layout_gem(const igraph_t *graph, igraph_matrix_t *res,
                      igraph_bool_t use_seed, igraph_integer_t maxiter,
                      igraph_real_t temp_max, igraph_real_t temp_min,
                      igraph_real_t temp_init);

The GEM layout algorithm, as described in Arne Frick, Andreas Ludwig, Heiko Mehldau: A Fast Adaptive Layout Algorithm for Undirected Graphs, Proc. Graph Drawing 1994, LNCS 894, pp. 388-403, 1995.

Arguments: 

graph:

The input graph. Edge directions are ignored in directed graphs.

res:

The result is stored here. If the use_seed argument is true (non-zero), then this matrix is also used as the starting point of the algorithm.

use_seed:

Boolean, whether to use the supplied coordinates in res as the starting point. If false (zero), then a uniform random starting point is used.

maxiter:

The maximum number of iterations to perform. Updating a single vertex counts as an iteration. A reasonable default is 40 * n * n, where n is the number of vertices. The original paper suggests 4 * n * n, but this usually only works if the other parameters are set up carefully.

temp_max:

The maximum allowed local temperature. A reasonable default is the number of vertices.

temp_min:

The global temperature at which the algorithm terminates (even before reaching maxiter iterations). A reasonable default is 1/10.

temp_init:

Initial local temperature of all vertices. A reasonable default is the square root of the number of vertices.

Returns: 

Error code.

Time complexity: O(t * n * (n+e)), where n is the number of vertices, e is the number of edges and t is the number of time steps performed.

1.11. igraph_layout_davidson_harel — Davidson-Harel layout algorithm.

igraph_error_t igraph_layout_davidson_harel(const igraph_t *graph, igraph_matrix_t *res,
                                 igraph_bool_t use_seed, igraph_integer_t maxiter,
                                 igraph_integer_t fineiter, igraph_real_t cool_fact,
                                 igraph_real_t weight_node_dist, igraph_real_t weight_border,
                                 igraph_real_t weight_edge_lengths,
                                 igraph_real_t weight_edge_crossings,
                                 igraph_real_t weight_node_edge_dist);

This function implements the algorithm by Davidson and Harel, see Ron Davidson, David Harel: Drawing Graphs Nicely Using Simulated Annealing. ACM Transactions on Graphics 15(4), pp. 301-331, 1996. https://doi.org/10.1145/234535.234538

The algorithm uses simulated annealing and a sophisticated energy function, which is unfortunately hard to parameterize for different graphs. The original publication did not disclose any parameter values, and the ones below were determined by experimentation.

The algorithm consists of two phases, an annealing phase, and a fine-tuning phase. There is no simulated annealing in the second phase.

Our implementation tries to follow the original publication, as much as possible. The only major difference is that coordinates are explicitly kept within the bounds of the rectangle of the layout.

Arguments: 

graph:

The input graph, edge directions are ignored.

res:

A matrix, the result is stored here. It can be used to supply start coordinates, see use_seed.

use_seed:

Boolean, whether to use the supplied res as start coordinates.

maxiter:

The maximum number of annealing iterations. A reasonable value for smaller graphs is 10.

fineiter:

The number of fine tuning iterations. A reasonable value is max(10, log2(n)) where n is the number of vertices.

cool_fact:

Cooling factor. A reasonable value is 0.75.

weight_node_dist:

Weight for the node-node distances component of the energy function. Reasonable value: 1.0.

weight_border:

Weight for the distance from the border component of the energy function. It can be set to zero, if vertices are allowed to sit on the border.

weight_edge_lengths:

Weight for the edge length component of the energy function, a reasonable value is the density of the graph divided by 10.

weight_edge_crossings:

Weight for the edge crossing component of the energy function, a reasonable default is 1 minus the square root of the density of the graph.

weight_node_edge_dist:

Weight for the node-edge distance component of the energy function. A reasonable value is 1 minus the density, divided by 5.

Returns: 

Error code.

Time complexity: one first phase iteration has time complexity O(n^2+m^2), one fine tuning iteration has time complexity O(mn). Time complexity might be smaller if some of the weights of the components of the energy function are set to zero.

1.12. igraph_layout_mds — Place the vertices on a plane using multidimensional scaling.

igraph_error_t igraph_layout_mds(const igraph_t* graph, igraph_matrix_t *res,
                      const igraph_matrix_t *dist, igraph_integer_t dim);

This layout requires a distance matrix, where the intersection of row i and column j specifies the desired distance between vertex i and vertex j. The algorithm will try to place the vertices in a space having a given number of dimensions in a way that approximates the distance relations prescribed in the distance matrix. igraph uses the classical multidimensional scaling by Torgerson; for more details, see Cox & Cox: Multidimensional Scaling (1994), Chapman and Hall, London.

If the input graph is disconnected, igraph will decompose it first into its subgraphs, lay out the subgraphs one by one using the appropriate submatrices of the distance matrix, and then merge the layouts using igraph_layout_merge_dla. Since igraph_layout_merge_dla works for 2D layouts only, you cannot run the MDS layout on disconnected graphs for more than two dimensions.

Warning: if the graph is symmetric to the exchange of two vertices (as is the case with leaves of a tree connecting to the same parent), classical multidimensional scaling may assign the same coordinates to these vertices.

Arguments: 

graph:

A graph object.

res:

Pointer to an initialized matrix object. This will contain the result and will be resized if needed.

dist:

The distance matrix. It must be symmetric and this function does not check whether the matrix is indeed symmetric. Results are unspecified if you pass a non-symmetric matrix here. You can set this parameter to null; in this case, the undirected shortest path lengths between vertices will be used as distances.

dim:

The number of dimensions in the embedding space. For 2D layouts, supply 2 here.

Returns: 

Error code.

Added in version 0.6.

Time complexity: usually around O(|V|^2 dim).

1.13. igraph_layout_lgl — Force based layout algorithm for large graphs.

igraph_error_t igraph_layout_lgl(const igraph_t *graph, igraph_matrix_t *res,
                      igraph_integer_t maxit, igraph_real_t maxdelta,
                      igraph_real_t area, igraph_real_t coolexp,
                      igraph_real_t repulserad, igraph_real_t cellsize,
                      igraph_integer_t proot);

This is a layout generator similar to the Large Graph Layout algorithm and program (http://lgl.sourceforge.net/). But unlike LGL, this version uses a Fruchterman-Reingold style simulated annealing algorithm for placing the vertices. The speedup is achieved by placing the vertices on a grid and calculating the repulsion only for vertices which are closer to each other than a limit.

Arguments: 

graph:

The (initialized) graph object to place. It must be connnected; disconnected graphs are not handled by the algorithm.

res:

Pointer to an initialized matrix object to hold the result. It will be resized if needed.

maxit:

The maximum number of cooling iterations to perform for each layout step. A reasonable default is 150.

maxdelta:

The maximum length of the move allowed for a vertex in a single iteration. A reasonable default is the number of vertices.

area:

This parameter gives the area of the square on which the vertices will be placed. A reasonable default value is the number of vertices squared.

coolexp:

The cooling exponent. A reasonable default value is 1.5.

repulserad:

Determines the radius at which vertex-vertex repulsion cancels out attraction of adjacent vertices. A reasonable default value is area times the number of vertices.

cellsize:

The size of the grid cells, one side of the square. A reasonable default value is the fourth root of area (or the square root of the number of vertices if area is also left at its default value).

proot:

The root vertex, this is placed first, its neighbors in the first iteration, second neighbors in the second, etc. If negative then a random vertex is chosen.

Returns: 

Error code.

Added in version 0.2.

Time complexity: ideally O(dia*maxit*(|V|+|E|)), |V| is the number of vertices, dia is the diameter of the graph, worst case complexity is still O(dia*maxit*(|V|^2+|E|)), this is the case when all vertices happen to be in the same grid cell.

2. Layouts for trees and acyclic graphs

2.1. igraph_layout_reingold_tilford — Reingold-Tilford layout for tree graphs.

igraph_error_t igraph_layout_reingold_tilford(const igraph_t *graph,
                                   igraph_matrix_t *res,
                                   igraph_neimode_t mode,
                                   const igraph_vector_int_t *roots,
                                   const igraph_vector_int_t *rootlevel);

Arranges the nodes in a tree where the given node is used as the root. The tree is directed downwards and the parents are centered above its children. For the exact algorithm, see:

Reingold, E and Tilford, J: Tidier drawing of trees. IEEE Trans. Softw. Eng., SE-7(2):223--228, 1981. https://doi.org/10.1109/TSE.1981.234519

If the given graph is not a tree, a breadth-first search is executed first to obtain a possible spanning tree.

Arguments: 

graph:

The graph object.

res:

The result, the coordinates in a matrix. The parameter should point to an initialized matrix object and will be resized.

mode:

Specifies which edges to consider when building the tree. If it is IGRAPH_OUT then only the outgoing, if it is IGRAPH_IN then only the incoming edges of a parent are considered. If it is IGRAPH_ALL then all edges are used (this was the behavior in igraph 0.5 and before). This parameter also influences how the root vertices are calculated, if they are not given. See the roots parameter.

roots:

The index of the root vertex or root vertices. The set of roots should be specified so that all vertices of the graph are reachable from them. Simply put, in the undirected case, one root should be given from each connected component. If roots is NULL or a pointer to an empty vector, then the roots will be selected automatically. Currently, automatic root selection prefers low eccentricity vertices in graphs with fewer than 500 vertices, and high degree vertices (according to mode) in larger graphs. The root selection heuristic may change without notice. To ensure a consistent output, please specify the roots manually. The igraph_roots_for_tree_layout() function gives more control over automatic root selection.

rootlevel:

This argument can be useful when drawing forests which are not trees (i.e. they are unconnected and have tree components). It specifies the level of the root vertices for every tree in the forest. It is only considered if not a null pointer and the roots argument is also given (and it is not a null pointer of an empty vector).

Returns: 

Error code.

Added in version 0.2.

See also: 

Example 20.1.  File examples/simple/igraph_layout_reingold_tilford.c

#include <igraph.h>
#include <math.h>

int main(void) {

    igraph_t g;
    FILE *f;
    igraph_matrix_t coords;
    /* igraph_integer_t i, n; */

    f = fopen("igraph_layout_reingold_tilford.in", "r");
    igraph_read_graph_edgelist(&g, f, 0, IGRAPH_DIRECTED);
    igraph_matrix_init(&coords, 0, 0);
    igraph_layout_reingold_tilford(&g, &coords, IGRAPH_IN, 0, 0);

    /*n=igraph_vcount(&g);
    for (i=0; i<n; i++) {
      printf("%6.3f %6.3f\n", MATRIX(coords, i, 0), MATRIX(coords, i, 1));
    }*/

    igraph_matrix_destroy(&coords);
    igraph_destroy(&g);
    return 0;
}


2.2. igraph_layout_reingold_tilford_circular — Circular Reingold-Tilford layout for trees.

igraph_error_t igraph_layout_reingold_tilford_circular(const igraph_t *graph,
        igraph_matrix_t *res,
        igraph_neimode_t mode,
        const igraph_vector_int_t *roots,
        const igraph_vector_int_t *rootlevel);

This layout is almost the same as igraph_layout_reingold_tilford(), but the tree is drawn in a circular way, with the root vertex in the center.

Arguments: 

graph:

The graph object.

res:

The result, the coordinates in a matrix. The parameter should point to an initialized matrix object and will be resized.

mode:

Specifies which edges to consider when building the tree. If it is IGRAPH_OUT then only the outgoing, if it is IGRAPH_IN then only the incoming edges of a parent are considered. If it is IGRAPH_ALL then all edges are used (this was the behavior in igraph 0.5 and before). This parameter also influences how the root vertices are calculated, if they are not given. See the roots parameter.

roots:

The index of the root vertex or root vertices. The set of roots should be specified so that all vertices of the graph are reachable from them. Simply put, in the undirected case, one root should be given from each connected component. If roots is NULL or a pointer to an empty vector, then the roots will be selected automatically. Currently, automatic root selection prefers low eccentricity vertices in graphs with fewer than 500 vertices, and high degree vertices (according to mode) in larger graphs. The root selection heuristic may change without notice. To ensure a consistent output, please specify the roots manually.

rootlevel:

This argument can be useful when drawing forests which are not trees (i.e. they are unconnected and have tree components). It specifies the level of the root vertices for every tree in the forest. It is only considered if not a null pointer and the roots argument is also given (and it is not a null pointer or an empty vector).

Returns: 

Error code.

See also: 

2.3. igraph_roots_for_tree_layout — Roots suitable for a nice tree layout.

igraph_error_t igraph_roots_for_tree_layout(
        const igraph_t *graph,
        igraph_neimode_t mode,
        igraph_vector_int_t *roots,
        igraph_root_choice_t heuristic);

This function chooses a root, or a set of roots suitable for visualizing a tree, or a tree-like graph. It is typically used with igraph_layout_reingold_tilford(). The principle is to select a minimal set of roots so that all other vertices will be reachable from them.

In the undirected case, one root is chosen from each connected component. In the directed case, one root is chosen from each strongly connected component that has no incoming (or outgoing) edges (depending on 'mode'). When more than one root choice is possible, vertices are prioritized based on the given heuristic.

Arguments: 

graph:

The graph, typically a tree, but any graph is accepted.

mode:

Whether to interpret the input as undirected, a directed out-tree or in-tree.

roots:

An initialized integer vector, the roots will be returned here.

heuristic:

The heuristic to use for breaking ties when multiple root choices are possible.

IGRAPH_ROOT_CHOICE_DEGREE

Choose the vertices with the highest degree (out- or in-degree in directed mode). This simple heuristic is fast even in large graphs.

IGRAPH_ROOT_CHOICE_ECCENTRICITY

Choose the vertices with the lowest eccentricity. This usually results in a "wide and shallow" tree layout. While this heuristic produces high-quality results, it is slow for large graphs: computing the eccentricities has quadratic complexity in the number of vertices.

Returns: 

Error code.

Time complexity: depends on the heuristic.

2.4. igraph_layout_sugiyama — Sugiyama layout algorithm for layered directed acyclic graphs.

igraph_error_t igraph_layout_sugiyama(const igraph_t *graph, igraph_matrix_t *res,
                           igraph_t *extd_graph, igraph_vector_int_t *extd_to_orig_eids,
                           const igraph_vector_int_t* layers, igraph_real_t hgap, igraph_real_t vgap,
                           igraph_integer_t maxiter, const igraph_vector_t *weights);

This layout algorithm is designed for directed acyclic graphs where each vertex is assigned to a layer. Layers are indexed from zero, and vertices of the same layer will be placed on the same horizontal line. The X coordinates of vertices within each layer are decided by the heuristic proposed by Sugiyama et al to minimize edge crossings.

You can also try to lay out undirected graphs, graphs containing cycles, or graphs without an a priori layered assignment with this algorithm. igraph will try to eliminate cycles and assign vertices to layers, but there is no guarantee on the quality of the layout in such cases.

The Sugiyama layout may introduce "bends" on the edges in order to obtain a visually more pleasing layout. This is achieved by adding dummy nodes to edges spanning more than one layer. The resulting layout assigns coordinates not only to the nodes of the original graph but also to the dummy nodes. The layout algorithm will also return the extended graph with the dummy nodes. An edge in the original graph may either be mapped to a single edge in the extended graph or a path that starts and ends in the original source and target vertex and passes through multiple dummy vertices. In such cases, the user may also request the mapping of the edges of the extended graph back to the edges of the original graph.

For more details, see K. Sugiyama, S. Tagawa and M. Toda, "Methods for Visual Understanding of Hierarchical Systems". IEEE Transactions on Systems, Man and Cybernetics 11(2):109-125, 1981.

Arguments: 

graph:

Pointer to an initialized graph object.

res:

Pointer to an initialized matrix object. This will contain the result and will be resized as needed. The first |V| rows of the layout will contain the coordinates of the original graph, the remaining rows contain the positions of the dummy nodes. Therefore, you can use the result both with graph or with extended_graph.

extended_graph:

Pointer to an uninitialized graph object or NULL. The extended graph with the added dummy nodes will be returned here. In this graph, each edge points downwards to lower layers, spans exactly one layer and the first |V| vertices coincide with the vertices of the original graph.

extd_to_orig_eids:

Pointer to a vector or NULL. If not NULL, the mapping from the edge IDs of the extended graph back to the edge IDs of the original graph will be stored here.

layers:

The layer index for each vertex or NULL if the layers should be determined automatically by igraph.

hgap:

The preferred minimum horizontal gap between vertices in the same layer.

vgap:

The distance between layers.

maxiter:

Maximum number of iterations in the crossing minimization stage. 100 is a reasonable default; if you feel that you have too many edge crossings, increase this.

weights:

Weights of the edges. These are used only if the graph contains cycles; igraph will tend to reverse edges with smaller weights when breaking the cycles.

2.5. igraph_layout_umap — Layout using Uniform Manifold Approximation and Projection (UMAP).

igraph_error_t igraph_layout_umap(const igraph_t *graph,
                                  igraph_matrix_t *res,
                                  igraph_bool_t use_seed,
                                  const igraph_vector_t *distances,
                                  igraph_real_t min_dist,
                                  igraph_integer_t epochs,
                                  igraph_bool_t distances_are_weights);

Warning

This function is experimental and its signature is not considered final yet. We reserve the right to change the function signature without changing the major version of igraph. Use it at your own risk.

UMAP is mostly used to embed high-dimensional vectors in a low-dimensional space (most commonly 2D). The algorithm is probabilistic and introduces nonlinearities, unlike e.g. PCA and similar to T-distributed Stochastic Neighbor Embedding (t-SNE). Nonlinearity helps "cluster" very similar vectors together without imposing a global geometry on the embedded space (e.g. a rigid rotation + compression in PCA). UMAP uses graphs as intermediate data structures, hence it can be used as a graph layout algorithm as well.

The general UMAP workflow is to start from vectors, compute a sparse distance graph that only contains edges between simiar points (e.g. a k-nearest neighbors graph), and then convert these distances into exponentially decaying weights between 0 and 1 that are larger for points that are closest neighbors in the distance graph. If a graph without any distances associated to the edges is used, all weights will be set to 1.

If you are trying to use this function to embed high-dimensional vectors, you should first compute a k-nearest neighbors graph between your vectors and compute the associated distances, and then call this function on that graph. If you already have a distance graph, or you have a graph with no distances, you can call this function directly. If you already have a graph with meaningful weights associated to each edge, you can also call this function, but set the argument distances_are_weights to true. To compute weights from distances without computing the layout, see igraph_layout_umap_compute_weights().

References:

Leland McInnes, John Healy, and James Melville. https://arxiv.org/abs/1802.03426

Arguments: 

graph:

Pointer to the graph to find a layout for (i.e. to embed). This is typically a sparse graph with only edges for the shortest distances stored, e.g. a k-nearest neighbors graph.

res:

Pointer to the n by 2 matrix where the layout coordinates will be stored.

use_seed:

Logical, if true the supplied values in the res argument are used as an initial layout, if false a random initial layout is used.

distances:

Pointer to a vector of distances associated with the graph edges. If this argument is NULL, all weights will be set to 1.

min_dist:

A fudge parameter that decides how close two unconnected vertices can be in the embedding before feeling a repulsive force. It must not be negative. Typical values are between 0 and 1.

epochs:

Number of iterations of the main stochastic gradient descent loop on the cross-entropy. Typical values are between 30 and 500.

distances_are_weights:

Whether to use precomputed weights. If true, the "distances" vector contains precomputed weights. If false (the typical use case), this function will compute weights from distances and then use them to compute the layout.

Returns: 

Error code.

2.6. igraph_layout_umap_compute_weights — Compute weights for a UMAP layout starting from distances.

igraph_error_t igraph_layout_umap_compute_weights(
        const igraph_t *graph,
        const igraph_vector_t *distances,
        igraph_vector_t *weights);

Warning

This function is experimental and its signature is not considered final yet. We reserve the right to change the function signature without changing the major version of igraph. Use it at your own risk.

UMAP is used to embed high-dimensional vectors in a low-dimensional space (most commonly 2D). It uses a distance graph as an intermediate data structure, making it also a useful graph layout algorithm. See igraph_layout_umap() for more information.

An early step in UMAP is to compute exponentially decaying "weights" from the distance graph. Connectivities can also be viewed as edge weights that quantify similarity between two vertices. This function computes weights from the distance graph. To compute the layout from precomputed weights, call igraph_layout_umap() with the distances_are_weights argument set to true.

While the distance graph can be directed (e.g. in a k-nearest neighbors, it is clear *who* you are a neighbor of), the weights are usually undirected. Whenever two vertices are doubly connected in the distance graph, the resulting weight W is set as: W = W1 + W2 - W1 * W2 Because UMAP weights are interpreted as probabilities, this is just the probability that either edge is present, without double counting. It is called "fuzzy union" in the original UMAP implementation and is the default. One could also require that both edges are there, i.e. W = W1 * W2: this would represent the fuzzy intersection and is not implemented in igraph. As a consequence of this symmetrization, information is lost, i.e. one needs fewer weights than one had distances. To keep things efficient, here we set the weight for one of the two edges as above and the weight for its opposite edge as 0, so that it will be skipped in the UMAP gradient descent later on.

Technical note: For each vertex, this function computes its scale factor (sigma), its connectivity correction (rho), and finally the weights themselves.

References:

Leland McInnes, John Healy, and James Melville. https://arxiv.org/abs/1802.03426

Arguments: 

graph:

Pointer to the distance graph. This can be directed (e.g. connecting each vertex to its neighbors in a k-nearest neighbor) or undirected, but must have no loops nor parallel edges. The only exception is: if the graph is directed, having pairs of edges with opposite direction is accepted.

distances:

Pointer to the vector with the vertex-to-vertex distance associated with each edge. This argument can be NULL, in which case all edges are assumed to have the same distance.

weights:

Pointer to an initialized vector where the result will be stored. If the input graph is directed, the weights represent a symmetrized version which contains less information. Therefore, whenever two edges between the same vertices and opposite direction are present in the input graph, only one of the weights is set and the other is fixed to zero. That format is accepted by igraph_layout_umap(), which skips all zero-weight edges from the layout optimization.

Returns: 

Error code.

3. 3D layout generators

3.1. igraph_layout_random_3d — Places the vertices uniform randomly in a cube.

igraph_error_t igraph_layout_random_3d(const igraph_t *graph, igraph_matrix_t *res);

Vertex coordinates range from -1 to 1, and are placed in 3 columns of a matrix, with a row for each vertex.

Arguments: 

graph:

The graph to place.

res:

Pointer to an initialized matrix object. It will be resized to hold the result.

Returns: 

Error code. The current implementation always returns with success.

Added in version 0.2.

Time complexity: O(|V|), the number of vertices.

3.2. igraph_layout_sphere — Places vertices (more or less) uniformly on a sphere.

igraph_error_t igraph_layout_sphere(const igraph_t *graph, igraph_matrix_t *res);

The vertices are placed with approximately equal spacing on a spiral wrapped around a sphere, in the order of their vertex IDs. Vertices with consecutive vertex IDs are placed near each other.

The algorithm was described in the following paper:

Distributing many points on a sphere by E.B. Saff and A.B.J. Kuijlaars, Mathematical Intelligencer 19.1 (1997) 5--11. https://doi.org/10.1007/BF03024331

Arguments: 

graph:

Pointer to an initialized graph object.

res:

Pointer to an initialized matrix object. This will contain the result and will be resized as needed.

Returns: 

Error code. The current implementation always returns with success.

Added in version 0.2.

Time complexity: O(|V|), the number of vertices in the graph.

3.3. igraph_layout_grid_3d — Places the vertices on a regular grid in the 3D space.

igraph_error_t igraph_layout_grid_3d(const igraph_t *graph, igraph_matrix_t *res,
                          igraph_integer_t width, igraph_integer_t height);

Arguments: 

graph:

Pointer to an initialized graph object.

res:

Pointer to an initialized matrix object. This will contain the result and will be resized as needed.

width:

The number of vertices in a single row of the grid. When zero or negative, the width is determined automatically.

height:

The number of vertices in a single column of the grid. When zero or negative, the height is determined automatically.

Returns: 

Error code. The current implementation always returns with success.

Time complexity: O(|V|), the number of vertices.

3.4. igraph_layout_fruchterman_reingold_3d — 3D Fruchterman-Reingold algorithm.

igraph_error_t igraph_layout_fruchterman_reingold_3d(const igraph_t *graph,
        igraph_matrix_t *res,
        igraph_bool_t use_seed,
        igraph_integer_t niter,
        igraph_real_t start_temp,
        const igraph_vector_t *weights,
        const igraph_vector_t *minx,
        const igraph_vector_t *maxx,
        const igraph_vector_t *miny,
        const igraph_vector_t *maxy,
        const igraph_vector_t *minz,
        const igraph_vector_t *maxz);

This is the 3D version of the force based Fruchterman-Reingold layout. See igraph_layout_fruchterman_reingold() for the 2D version.

Arguments: 

graph:

Pointer to an initialized graph object.

res:

Pointer to an initialized matrix object. This will contain the result and will be resized as needed.

use_seed:

Logical, if true the supplied values in the res argument are used as an initial layout, if false a random initial layout is used.

niter:

The number of iterations to do. A reasonable default value is 500.

start_temp:

Start temperature. This is the maximum amount of movement alloved along one axis, within one step, for a vertex. Currently it is decreased linearly to zero during the iteration.

weights:

Pointer to a vector containing edge weights, the attraction along the edges will be multiplied by these. Weights must be positive. It will be ignored if it is a null-pointer.

minx:

Pointer to a vector, or a NULL pointer. If not a NULL pointer then the vector gives the minimum x coordinate for every vertex.

maxx:

Same as minx, but the maximum x coordinates.

miny:

Pointer to a vector, or a NULL pointer. If not a NULL pointer then the vector gives the minimum y coordinate for every vertex.

maxy:

Same as miny, but the maximum y coordinates.

minz:

Pointer to a vector, or a NULL pointer. If not a NULL pointer then the vector gives the minimum z coordinate for every vertex.

maxz:

Same as minz, but the maximum z coordinates.

Returns: 

Error code.

Added in version 0.2.

Time complexity: O(|V|^2) in each iteration, |V| is the number of vertices in the graph.

3.5. igraph_layout_kamada_kawai_3d — 3D version of the Kamada-Kawai layout generator.

igraph_error_t igraph_layout_kamada_kawai_3d(const igraph_t *graph, igraph_matrix_t *res,
                                  igraph_bool_t use_seed, igraph_integer_t maxiter,
                                  igraph_real_t epsilon, igraph_real_t kkconst,
                                  const igraph_vector_t *weights,
                                  const igraph_vector_t *minx, const igraph_vector_t *maxx,
                                  const igraph_vector_t *miny, const igraph_vector_t *maxy,
                                  const igraph_vector_t *minz, const igraph_vector_t *maxz);

This is the 3D version of igraph_layout_kamada_kawai(). See the documentation of that function for more information.

This layout algorithm is not suitable for large graphs. The memory requirements are of the order O(|V|^2).

Arguments: 

graph:

A graph object.

res:

Pointer to an initialized matrix object. This will contain the result (x-, y- and z-positions in columns one through three) and will be resized if needed.

use_seed:

Boolean, whether to use the values supplied in the res argument as the initial configuration. If zero and there are any limits on the z, y or z coordinates, then a random initial configuration is used. Otherwise the vertices are placed uniformly on a sphere of radius 1 as the initial configuration.

maxiter:

The maximum number of iterations to perform. A reasonable default value is at least ten (or more) times the number of vertices.

epsilon:

Stop the iteration, if the maximum delta value of the algorithm is smaller than still. It is safe to leave it at zero, and then maxiter iterations are performed.

kkconst:

The Kamada-Kawai vertex attraction constant. Typical value: number of vertices.

weights:

Edge weights, larger values will result longer edges. Weights must be positive. Pass NULL to assume unit weights for all edges.

minx:

Pointer to a vector, or a NULL pointer. If not a NULL pointer then the vector gives the minimum x coordinate for every vertex.

maxx:

Same as minx, but the maximum x coordinates.

miny:

Pointer to a vector, or a NULL pointer. If not a NULL pointer then the vector gives the minimum y coordinate for every vertex.

maxy:

Same as miny, but the maximum y coordinates.

minz:

Pointer to a vector, or a NULL pointer. If not a NULL pointer then the vector gives the minimum z coordinate for every vertex.

maxz:

Same as minz, but the maximum z coordinates.

Returns: 

Error code.

Time complexity: O(|V|) for each iteration, after an O(|V|^2 log|V|) initialization step. |V| is the number of vertices in the graph.

3.6. igraph_layout_umap_3d — 3D layout using UMAP.

igraph_error_t igraph_layout_umap_3d(const igraph_t *graph,
                                     igraph_matrix_t *res,
                                     igraph_bool_t use_seed,
                                     const igraph_vector_t *distances,
                                     igraph_real_t min_dist,
                                     igraph_integer_t epochs,
                                     igraph_bool_t distances_are_weights);

Warning

This function is experimental and its signature is not considered final yet. We reserve the right to change the function signature without changing the major version of igraph. Use it at your own risk.

This is the 3D version of the UMAP algorithm (see igraph_layout_umap() for the 2D version).

Arguments: 

graph:

Pointer to the graph to find a layout for (i.e. to embed). This is typically a directed, sparse graph with only edges for the shortest distances stored, e.g. a k-nearest neighbors graph with the edges going from each focal vertex to its neighbors. However, it can also be an undirected graph. If the distances_are_weights is true, this is treated as an undirected graph.

res:

Pointer to the n by 3 matrix where the layout coordinates will be stored.

use_seed:

Logical, if true the supplied values in the res argument are used as an initial layout, if false a random initial layout is used.

distances:

Pointer to a vector of distances associated with the graph edges. If this argument is NULL, all edges are assumed to have the same distance.

min_dist:

A fudge parameter that decides how close two unconnected vertices can be in the embedding before feeling a repulsive force. It must not be negative. Typical values are between 0 and 1.

epochs:

Number of iterations of the main stochastic gradient descent loop on the cross-entropy. Typical values are between 30 and 500.

distances_are_weights:

Whether to use precomputed weights. If false (the typical use case), this function will compute weights from distances and then use them to compute the layout. If true, the "distances" vector contains precomputed weights, including possibly some weights equal to zero that are inconsequential for the layout optimization.

Returns: 

Error code.

4. Merging layouts

4.1. igraph_layout_merge_dla — Merges multiple layouts by using a DLA algorithm.

igraph_error_t igraph_layout_merge_dla(
    const igraph_vector_ptr_t *thegraphs, const igraph_matrix_list_t *coords,
    igraph_matrix_t *res
);

Warning

This function is experimental and its signature is not considered final yet. We reserve the right to change the function signature without changing the major version of igraph. Use it at your own risk.

First each layout is covered by a circle. Then the layout of the largest graph is placed at the origin. Then the other layouts are placed by the DLA algorithm, larger ones first and smaller ones last.

Arguments: 

thegraphs:

Pointer vector containing the graph objects of which the layouts will be merged.

coords:

List of matrices with the 2D layouts of the graphs in thegraphs.

res:

Pointer to an initialized matrix object, the result will be stored here. It will be resized if needed.

Returns: 

Error code.

Added in version 0.2.

Time complexity: TODO.