Milo.Rd
The Milo class extends the SingleCellExperiment class and is designed to work with neighbourhoods of cells. Therefore, it inherits from the SingleCellExperiment class and follows the same usage conventions. There is additional support for cell-to-cell distances via distance, and the KNN-graph used to define the neighbourhoods.
Milo( ..., graph = list(), nhoodDistances = Matrix(0L, sparse = TRUE), nhoods = list(), nhoodCounts = Matrix(0L, sparse = TRUE), nhoodIndex = list(), nhoodExpression = Matrix(0L, sparse = TRUE) )
... | Arguments passed to the Milo constructor to fill the slots of the
base class. This should be either a |
---|---|
graph | An igraph object or list of adjacent vertices that represents the KNN-graph |
nhoodDistances | A list containing sparse matrices of cell-to-cell distances for cells in the same neighbourhoods, one list entry per neighbourhood. |
nhoods | A list of graph vertices, each containing the indices of the constiuent graph vertices in the respective neighbourhood |
nhoodCounts | A matrix of neighbourhood X sample counts of the number of cells in each neighbourhood derived from the respective samples |
nhoodIndex | A list of cells that are the neighborhood index cells. |
nhoodExpression | A matrix of gene X neighbourhood expression. |
a Milo object
In this class the underlying structure is the gene/feature X cell expression data. The additional slots provide a link between these single cells and the neighbourhood representation. This can be further extended by the use of an abstracted graph for visualisation that preserves the structure of the single-cell KNN-graph
A Milo object can also be constructed by inputting a feature X cell gene expression matrix. In this case it simply constructs a SingleCellExperiment and fills the relevant slots, such as reducedDims.
library(SingleCellExperiment)#>#>#>#>#>#> #>#>#> #> #> #>#>#> #>#>#> #>#>#> #> #> #> #> #> #>#>#> #>#>#> #>#>#>#>#>#> #> #> #>#>#>#> #>#>#> #>#>#> #>#>#> #>#>#> #>#> #>#>#> #>ux <- matrix(rpois(12000, 5), ncol=200) vx <- log2(ux + 1) pca <- prcomp(t(vx)) sce <- SingleCellExperiment(assays=list(counts=ux, logcounts=vx), reducedDims=SimpleList(PCA=pca$x)) milo <- Milo(sce) milo#> class: Milo #> dim: 60 200 #> metadata(0): #> assays(2): counts logcounts #> rownames: NULL #> rowData names(0): #> colnames: NULL #> colData names(0): #> reducedDimNames(1): PCA #> spikeNames(0): #> altExpNames(0): #> nhoods dimensions(1): 0 #> nhoodCounts dimensions(2): 1 1 #> nhoodDistances dimension(1): 0 #> graph names(0): #> nhoodIndex names(1): 0 #> nhoodExpression dimension(2): 1 1 #> nhoodReducedDim names(0): #> nhoodGraph names(0): #> nhoodAdjacency dimension(0):