scran
C++ library for basic single-cell RNA-seq analyses
Loading...
Searching...
No Matches
Classes | Public Types | Public Member Functions | List of all members
scran::BuildSnnGraph Class Reference

Build a shared nearest-neighbor graph with cells as nodes. More...

#include <BuildSnnGraph.hpp>

Classes

struct  Defaults
 Default parameter settings. More...
 
struct  Results
 Results of SNN graph construction. More...
 

Public Types

enum  Scheme { RANKED , NUMBER , JACCARD }
 

Public Member Functions

BuildSnnGraphset_neighbors (int k=Defaults::neighbors)
 
BuildSnnGraphset_approximate (bool a=Defaults::approximate)
 
BuildSnnGraphset_weighting_scheme (Scheme w=Defaults::weighting_scheme)
 
BuildSnnGraphset_num_threads (int n=Defaults::num_threads)
 
Results run (size_t ndims, size_t ncells, const double *mat) const
 
template<class Algorithm >
Results run (const Algorithm *search) const
 
template<typename Index_ , typename Distance_ >
Results run (const knncolle::NeighborList< Index_, Distance_ > &neighbors) const
 
template<class Indices_ >
Results run (const std::vector< Indices_ > &indices) const
 

Detailed Description

Build a shared nearest-neighbor graph with cells as nodes.

In a shared nearest neighbor graph, pairs of cells are connected to each other by an edge with weight determined from their shared nearest neighbors. If two cells are close together but have distinct sets of neighbors, the corresponding edge is downweighted as the two cells are unlikely to be part of the same neighborhood. In this manner, highly weighted edges will form within highly interconnected neighborhoods where many cells share the same neighbors. This provides a more sophisticated definition of similarity between cells compared to a simpler (unweighted) nearest neighbor graph that just focuses on immediate proximity.

A key parameter in the construction of the graph is the number of nearest neighbors $k$ to consider. Larger values increase the connectivity of the graph and reduce the granularity of any subsequent community detection steps (see scran::ClusterSnnGraph) at the cost of speed. The nearest neighbor search can be performed using either vantage point trees (exact) or with the Annoy algorithm (approximate) - see the knncolle library for details.

For the edges, a variety of weighting schemes are possible:

See the ClusterSNNGraph class to perform community detection on the graph returned by run().

See also
Xu C and Su Z (2015). Identification of cell types from single-cell transcriptomes using a novel clustering method. Bioinformatics 31, 1974-80

Member Enumeration Documentation

◆ Scheme

Choices for the edge weighting scheme during graph construction.

Member Function Documentation

◆ set_neighbors()

BuildSnnGraph & scran::BuildSnnGraph::set_neighbors ( int  k = Defaults::neighbors)
inline
Parameters
kNumber of neighbors to use in the nearest neighbor search.
Returns
A reference to this BuildSnnGraph object.

◆ set_approximate()

BuildSnnGraph & scran::BuildSnnGraph::set_approximate ( bool  a = Defaults::approximate)
inline
Parameters
aWhether to perform an approximate nearest neighbor search.
Returns
A reference to this BuildSnnGraph object.

◆ set_weighting_scheme()

BuildSnnGraph & scran::BuildSnnGraph::set_weighting_scheme ( Scheme  w = Defaults::weighting_scheme)
inline
Parameters
wThe edge weighting scheme to use.
Returns
A reference to this BuildSnnGraph object.

◆ set_num_threads()

BuildSnnGraph & scran::BuildSnnGraph::set_num_threads ( int  n = Defaults::num_threads)
inline
Parameters
nNumber of threads to use.
Returns
A reference to this BuildSnnGraph object.

◆ run() [1/4]

Results scran::BuildSnnGraph::run ( size_t  ndims,
size_t  ncells,
const double *  mat 
) const
inline
Parameters
ndimsNumber of dimensions.
ncellsNumber of cells.
matPointer to an array of expression values or a low-dimensional representation thereof. Rows should be dimensions while columns should be cells. Data should be stored in column-major format.
Returns
The edges and weights of the constructed SNN graph.

◆ run() [2/4]

template<class Algorithm >
Results scran::BuildSnnGraph::run ( const Algorithm *  search) const
inline
Template Parameters
AlgorithmAny instance of a knncolle::Base subclass.
Parameters
searchPointer to a knncolle::Base instance to use for the nearest-neighbor search.
Returns
The edges and weights of the constructed SNN graph.

◆ run() [3/4]

template<typename Index_ , typename Distance_ >
Results scran::BuildSnnGraph::run ( const knncolle::NeighborList< Index_, Distance_ > &  neighbors) const
inline
Template Parameters
Index_Integer type for the indices.
Distance_Floating-point type for the distances.
Parameters
neighborsVector of indices and distances of the neighbors for each cell, sorted by increasing distance.
Returns
The edges and weights of the constructed SNN graph.

Distances are ignored here; this overload is only provided to enable convenient usage with pre-computed neighbors from knncolle.

◆ run() [4/4]

template<class Indices_ >
Results scran::BuildSnnGraph::run ( const std::vector< Indices_ > &  indices) const
inline
Template Parameters
Indices_Vector-like class containing integer indices. This should provide the [ and size() methods.
Parameters
indicesVector of indices of the neighbors for each cell, sorted by increasing distance.
Returns
The edges and weights of the constructed SNN graph.

The documentation for this class was generated from the following file: