Title: | Modeling networks for individual samples using LIONESS |
---|---|
Description: | LIONESS, or Linear Interpolation to Obtain Network Estimates for Single Samples, can be used to reconstruct single-sample networks (https://arxiv.org/abs/1505.06440). This code implements the LIONESS equation in the lioness function in R to reconstruct single-sample networks. The default network reconstruction method we use is based on Pearson correlation. However, lionessR can run on any network reconstruction algorithms that returns a complete, weighted adjacency matrix. lionessR works for both unipartite and bipartite networks. |
Authors: | Marieke Lydia Kuijjer [aut] , Ping-Han Hsieh [cre] |
Maintainer: | Ping-Han Hsieh <[email protected]> |
License: | MIT + file LICENSE |
Version: | 1.21.0 |
Built: | 2024-11-29 06:40:43 UTC |
Source: | https://github.com/bioc/lionessR |
This function uses the LIONESS equation to estimate single-sample networks. The input supports numeric matrix or a SummerizedExperiment object.
lioness(x, f = netFun)
lioness(x, f = netFun)
x |
Numeric matrix with samples in columns or SummerizedExperiment object |
f |
Network reconstruction function. Defaults to Pearson correlation. |
SummerizedExperiment object for single-sample association network. The rowData contains the information about the regulators and their targets, while the colData contains the information about the samples. The edge weights of sample specific networks can be accessed through the lioness assay of the object.
exp <- matrix(sample(1000,1000)/1000, 100, 10) genes <- paste("gene", c(1:nrow(exp)), sep="_") samples <- paste("sample", c(1:ncol(exp)), sep="_") rowData <- S4Vectors::DataFrame(row.names = genes, gene = genes) colData <- S4Vectors::DataFrame(col.names = samples, sample = samples) se <- SummarizedExperiment::SummarizedExperiment(assays = list(counts = exp), colData = colData, rowData = rowData) lionessResults <- lioness(se, netFun)
exp <- matrix(sample(1000,1000)/1000, 100, 10) genes <- paste("gene", c(1:nrow(exp)), sep="_") samples <- paste("sample", c(1:ncol(exp)), sep="_") rowData <- S4Vectors::DataFrame(row.names = genes, gene = genes) colData <- S4Vectors::DataFrame(col.names = samples, sample = samples) se <- SummarizedExperiment::SummarizedExperiment(assays = list(counts = exp), colData = colData, rowData = rowData) lionessResults <- lioness(se, netFun)
This is the network reconstruction function that will be used to build aggregate networks.
netFun(x)
netFun(x)
x |
Numeric matrix with samples in columns. |
Numeric matrix of association network.
exp <- matrix(sample(1000, 1000)/1000, 100, 10) netFun(exp)
exp <- matrix(sample(1000, 1000)/1000, 100, 10) netFun(exp)
Pre-processed gene expression data from high-grade osteosarcoma biopsies and sample characteristics were downloaded from GEO (GSE42352). We converted nuIDs to gene symbols using the annotation platform GPL10295. For genes with duplicate gene symbols, we selected the gene with the highest variance. Finally, we subsetted the data to the 53 patients for which 5 year metastasis information was available. Finally, we extracted genes with top 10,000 variance to reduce the size of the dataset.
data(OSdata)
data(OSdata)
Variable "exp": Data frame containing expression data for 10000 genes and 53 samples. Variable "targets": Data frame containing information on whether patients developed metastases within 5 years or not, 53 samples and 2 columns.