ConScape.jl

A Julia package for landscape connectivity.

Example

In the following, a small example of the functionality is given. First, we canstruct an artificial landscape with a wall and two corridors. The landscape is stored in a Grid struct.

using ConScape
g = ConScape.perm_wall_sim(30, 60, corridorwidths=(3,2), costs=ConScape.MinusLog()) # Generate an artificial landscape
ConScape.plot_outdegrees(g)

From a Grid, we can now create a GridRSP which we can use to compute the randomized shortest path based quality weighted betweenness with the temperature parameter θ=0.2.

h = ConScape.GridRSP(g, θ=0.2)
bet_q = ConScape.betweenness_qweighted(h)
ConScape.heatmap(bet_q, yflip=true)

Details

The section provides more details about the functions of this package including exact call signatures.

Grid

ConScape.GridType
Grid(nrows::Integer,
          ncols::Integer;
          affinities=nothing,
          qualities::Matrix=ones(nrows, ncols),
          source_qualities::Matrix=qualities,
          target_qualities::AbstractMatrix=qualities,
          costs::Union{Transformation,SparseMatrixCSC{Float64,Int}}=MinusLog(),
          prune=true)::Grid

Construct a Grid from an affinities matrix of type SparseMatrixCSC. It is possible to also supply matrices of source_qualities and target_qualities as well as a costs function that maps the affinities matrix to a costs matrix. Alternatively, it is possible to supply a matrix to costs directly. If prune=true (the default), the affinity and cost matrices will be pruned to exclude unreachable nodes.

ConScape.free_energy_distanceType
free_energy_distance(
    g::Grid;
    target::Union{Tuple{Int,Int},Nothing}=nothing,
    θ::Union{Real,Nothing}=nothing,
    approx::Bool=false
)

Compute the directed free energy distance from all source nodes to all target nodes in the graph defined by g using the inverse temperature parameter θ. The computation can either continue until convergence when setting approx=false (the default) or return an approximate result based on just a single iteration of the Bellman-Ford algorithm when approx=true.

LightGraphs.is_strongly_connectedFunction
is_strongly_connected(g::Grid)::Bool

Test if graph defined by Grid is fully connected.

Examples

julia> affinities = [1/4 0 1/4 1/4
                     1/4 0 1/4 1/4
                     1/4 0 1/4 1/4
                     1/4 0 1/4 1/4];

julia> grid = ConScape.Grid(size(affinities)..., affinities=ConScape.graph_matrix_from_raster(affinities), prune=false)
ConScape.Grid of size 4x4

julia> ConScape.is_strongly_connected(grid)
false
ConScape.largest_subgraphFunction
largest_subgraph(g::Grid)::Grid

Extract the largest fully connected subgraph of the Grid. The returned Grid will have the same size as the input Grid but only nodes associated with the largest subgraph of the affinities will be active.

ConScape.least_cost_distanceType
least_cost_distance(g::Grid)::Matrix{Float64}

Compute the least cost distance from all the cells in the grid to all target cells.

Examples

julia> affinities = [1/4 0 1/2 1/4
                     1/4 0 1/2 1/4
                     1/4 0 1/2 1/4
                     1/4 0 1/2 1/4];

julia> grid = ConScape.Grid(size(affinities)..., affinities=ConScape.graph_matrix_from_raster(affinities))
[ Info: cost graph contains 6 strongly connected subgraphs
[ Info: removing 8 nodes from affinity and cost graphs
ConScape.Grid of size 4x4

julia> ConScape.least_cost_distance(grid)
8×8 Matrix{Float64}:
 0.0       0.693147  1.38629   2.07944   0.693147  1.03972   1.73287   2.42602
 0.693147  0.0       0.693147  1.38629   1.03972   0.693147  1.03972   1.73287
 1.38629   0.693147  0.0       0.693147  1.73287   1.03972   0.693147  1.03972
 2.07944   1.38629   0.693147  0.0       2.42602   1.73287   1.03972   0.693147
 1.38629   1.73287   2.42602   3.11916   0.0       1.38629   2.77259   3.46574
 1.73287   1.38629   1.73287   2.42602   1.38629   0.0       1.38629   2.77259
 2.42602   1.73287   1.38629   1.73287   2.77259   1.38629   0.0       1.38629
 3.11916   2.42602   1.73287   1.38629   3.46574   2.77259   1.38629   0.0
ConScape.sum_neighborhoodFunction
sum_neighborhood(g::Grid, rc::Tuple{Int,Int}, npix::Integer)::Float64

A helper-function, used by coarse_graining, that computes the sum of pixels within a npix neighborhood around the target rc.

ConScape.coarse_grainingFunction
coarse_graining(g::Grid, npix::Integer)::Array

Creates a sparse matrix of target qualities for the landmarks based on merging npix pixels into the center pixel.

GridRSP (Randomized Shortest Path)

ConScape.GridRSPType
GridRSP(g::Grid; θ=nothing)::GridRSP

Construct a GridRSP from a g::Grid based on the inverse temperature parameter θ::Real.

ConScape.expected_costType
free_energy_distance(
    g::Grid;
    θ::Union{Real,Nothing}=nothing,
    approx::Bool=false
)

Compute the randomized shorted path based expected costs from all source nodes to all target nodes in the graph defined by g using the inverse temperature parameter θ. The computation can either continue until convergence when setting approx=false (the default) or return an approximate result based on just a single iteration of the Bellman-Ford algorithm when approx=true.

expected_cost(grsp::GridRSP)::Matrix{Float64}

Compute RSP expected costs from all nodes.

ConScape.betweenness_qweightedFunction
betweenness_qweighted(grsp::GridRSP)::Matrix{Float64}

Compute RSP betweenness of all nodes weighted by source and target qualities.

ConScape.betweenness_kweightedFunction
betweenness_kweighted(grsp::GridRSP;
    connectivity_function=expected_cost,
    distance_transformation=inv(grsp.g.costfunction),
    diagvalue=nothing])::SparseMatrixCSC{Float64,Int}

Compute RSP betweenness of all nodes weighted with proximities computed with respect to the distance/proximity measure defined by connectivity_function. Optionally, an inverse cost function can be passed. The function will be applied elementwise to the matrix of distances to convert it to a matrix of proximities. If no inverse cost function is passed the the inverse of the cost function is used for the conversion of distances.

The optional diagvalue element specifies which value to use for the diagonal of the matrix of proximities, i.e. after applying the inverse cost function to the matrix of distances. When nothing is specified, the diagonal elements won't be adjusted.

ConScape.edge_betweenness_qweightedFunction
edge_betweenness_qweighted(grsp::GridRSP)::Matrix{Float64}

Compute RSP betweenness of all edges weighted by source and target qualities. Returns a sparse matrix where element (i,j) is the betweenness of edge (i,j).

ConScape.edge_betweenness_kweightedFunction
edge_betweenness_kweighted(grsp::GridRSP; [distance_transformation=inv(grsp.g.costfunction), diagvalue=nothing])::SparseMatrixCSC{Float64,Int}

Compute RSP betweenness of all edges weighted by qualities of source s and target t and the proximity between s and t. Returns a
sparse matrix where element (i,j) is the betweenness of edge (i,j).

The optional `diagvalue` element specifies which value to use for the diagonal of the matrix
of proximities, i.e. after applying the inverse cost function to the matrix of expected costs.
When nothing is specified, the diagonal elements won't be adjusted.
ConScape.mean_kl_divergenceFunction
mean_kl_divergence(grsp::GridRSP)::Float64

Compute the mean Kullback–Leibler divergence between the free energy distances and the RSP expected costs for grsp::GridRSP.

ConScape.mean_lc_kl_divergenceFunction
mean_lc_kl_divergence(grsp::GridRSP)::Float64

Compute the mean Kullback–Leibler divergence between the least-cost path and the random path distribution for grsp::GridRSP, weighted by the qualities of the source and target node.

ConScape.least_cost_kl_divergenceFunction
least_cost_kl_divergence(grsp::GridRSP, target::Tuple{Int,Int})

Compute the least cost Kullback-Leibler divergence from each cell in the g in h to the target cell.

ConScape.connected_habitatFunction
connected_habitat(grsp::Union{Grid,GridRSP};
    connectivity_function=expected_cost,
    distance_transformation=nothing,
    diagvalue=nothing,
    θ::Union{Nothing,Real}=nothing,
    approx::Bool=false)::Matrix{Float64}

Compute RSP connected_habitat of all nodes. An inverse cost function must be passed for a Grid argument but is optional for GridRSP. The function will be applied elementwise to the matrix of distances to convert it to a matrix of proximities. If no inverse cost function is passed the the inverse of the cost function is used for the conversion of the proximities.

The optional diagvalue element specifies which value to use for the diagonal of the matrix of proximities, i.e. after applying the inverse cost function to the matrix of distances. When nothing is specified, the diagonal elements won't be adjusted.

connectivity_function determines which function is used for computing the matrix of proximities. If connectivity_function is a DistanceFunction, then it is used for computing distances, which is converted to proximities using distance_transformation. If connectivity_function is a ProximityFunction, then proximities are computed directly using it. The default is expected_cost.

For Grid objects, the inverse temperature parameter θ must be passed when the connectivity_function requires it such as expected_cost. Also for Grid objects, the approx Boolean argument can be set to true to switch to a cheaper approximate solution of the connectivity_function. The default value is false.

LinearAlgebra.eigmaxFunction
eigmax(grsp::GridRSP;
    connectivity_function=expected_cost,
    distance_transformation=nothing,
    diagvalue=nothing,
    tol=1e-14)

Compute the largest eigenvalue triple (left vector, value, and right vector) of the quality scaled proximities with respect to the distance/proximity measure defined by connectivity_function. If connectivity_function is a distance measure then the distances are transformed to proximities by distance_transformation which defaults to the inverse of the costfunction in the underlying Grid (if defined). Optionally, the diagonal values of the proximity matrix may be set to diagvalue. The tol argument specifies the convergence tolerance in the Arnoldi based eigensolver.

ConScape.criticalityFunction
criticality(grsp::GridRSP[;
            distance_transformation=inv(grsp.g.costfunction),
            diagvalue=nothing,
            avalue=floatmin(),
            qˢvalue=0.0,
            qᵗvalue=0.0])

Compute the landscape criticality for each target cell by setting setting affinities for the cell to avalue as well as the source and target qualities associated with the cell to qˢvalue and qᵗvalue respectively. It is required that avalue is positive to avoid that the graph becomes disconnected.

Utility functions

ConScape.graph_matrix_from_rasterFunction
graph_matrix_from_raster(R::Matrix[, type=AffinityMatrix, neighbors::Tuple=N8, weight=TargetWeight])::SparseMatrixCSC

Compute a graph matrix, i.e. an affinity or cost matrix of the raster image R of cell affinities or cell costs. The values are computed as either the value of the target cell (TargetWeight) or as harmonic (arithmetic) means of the cell affinities (costs) weighted by the grid distance (AverageWeight). The values can be computed with respect to eight neighbors (N8) or four neighbors (N4).

ConScape.mapnzFunction
mapnz(f, A::SparseMatrixCSC)::SparseMatrixCSC

Map the non-zero values of a sparse matrix A with the function f.

ConScape.readascFunction
readasc(::Union{String,IO}; nodata_value=0.0)::Tuple{Matrix{Float64}, Dict{String,Int}}

Read asc file of raster data and return tuple of a raster matrix and a dictionary containing the metadata information.

ConScape.writeascFunction
writeasc(fn::String, m::Matrix{<:Real}; kwargs...)

Write the raster matrix m to a file named fn. It's possible to pass metadata arguments xllcorner, yllcorner, cellsize, nodata_value as keywords or as a single dictionary.

This package is derived from the Python package reindeers.