Title: | A model for network-based prioritisation of genes |
---|---|
Description: | A model for semi-supervised prioritisation of genes integrating network data, phenotypes and additional prior knowledge about TP and TN gene labels from the literature or experts. |
Authors: | Fabian Schmich |
Maintainer: | Fabian Schmich <[email protected]> |
License: | GPL-3 |
Version: | 1.33.0 |
Built: | 2024-10-30 09:02:24 UTC |
Source: | https://github.com/bioc/netprioR |
This package provides a model for semi-supervised prioritisation of genes integrating network data, phenotypes and additional prior knowledge about TP and TN gene labels.
Fabian Schmich | Computational Biology Group, ETH Zurich | [email protected]
Fabian Schmich et. al (2016).
Compute the bandwidth of a matrix
bandwidth(x)
bandwidth(x)
x |
Inpute matrix |
Bandwidth
Fabian Schmich
Class Mass Normalization (CMN) from Zhu et al., 2003
cmn(yhat, l, u)
cmn(yhat, l, u)
yhat |
Response for labeled (l) and unlabeld (u) genes |
l |
Indices of labeled genes |
u |
Indices of unlabeled genes |
Class normalized yhat
Fabian Schmich
Solves linear equation systems iteratively
conjugate_gradient(A, b, x0 = rep(0, ncol(A)), threshold = 1e-15, verbose = FALSE)
conjugate_gradient(A, b, x0 = rep(0, ncol(A)), threshold = 1e-15, verbose = FALSE)
A |
Matrix |
b |
Coefficients |
x0 |
Starting solution |
threshold |
Termination threshold |
verbose |
Show iterative progress |
Solution for equation system
Fabian Schmich
Transform sparse matrix into a band matrix
cuthill_mckee(x)
cuthill_mckee(x)
x |
Input matrix |
Band matrix
Fabian Schmich
netprioR
modelFit netprioR
model
fit(object, ...) ## S4 method for signature 'netprioR' fit(object, refit = FALSE, ...)
fit(object, ...) ## S4 method for signature 'netprioR' fit(object, refit = FALSE, ...)
object |
A |
... |
Additional arguments |
refit |
Flag whether to overwrite existing fit |
A netprioR
object with fitted model
Fabian Schmich
data(simulation) np <- netprioR(networks = simulation$networks, phenotypes = simulation$phenotypes, labels = simulation$labels.obs, model.fit = FALSE) summary(np) np <- fit(np, nrestarts = 1, verbose = FALSE) summary(np)
data(simulation) np <- netprioR(networks = simulation$networks, phenotypes = simulation$phenotypes, labels = simulation$labels.obs, model.fit = FALSE) summary(np) np <- fit(np, nrestarts = 1, verbose = FALSE) summary(np)
Compute the Laplacian matrix of a graph given its adjacency matrix
laplacian(x, norm = c("none", "sym", "asym"))
laplacian(x, norm = c("none", "sym", "asym"))
x |
Adjacency matrix |
norm |
Type of normalisation |
Laplacian matrix
Fabian Schmich
Infer parameters and hidden data using the EM algorithm of netprioR
learn(Yobs, X, G, l, u, a = 0.1, b = 0.1, sigma2 = 1, tau2 = 10, eps = 1e-11, max.iter = 500, thresh = 0.001, use.cg = TRUE, thresh.cg = 1e-05, nrestarts = 5, max.cores = detectCores(), verbose = FALSE)
learn(Yobs, X, G, l, u, a = 0.1, b = 0.1, sigma2 = 1, tau2 = 10, eps = 1e-11, max.iter = 500, thresh = 0.001, use.cg = TRUE, thresh.cg = 1e-05, nrestarts = 5, max.cores = detectCores(), verbose = FALSE)
Yobs |
Observed labels (NA, if not observed) |
X |
Phenotypes |
G |
Graph Laplacians |
l |
Indices of labelled instances |
u |
Indices of unlabelled instances |
a |
Shape parameter of Gamma prior for W |
b |
Scale parameter of Gamma prior for W |
sigma2 |
Cariance for Gaussian labels |
tau2 |
Variance for Gaussian prior for beta |
eps |
Small value added to diagonal of Q in order to make it non-singular |
max.iter |
Maximum number of iterations for EM |
thresh |
Threshold for termination of EM with respect to change in parameters |
use.cg |
Flag whether to use conjugate gradient instead of exact computation of expectations |
thresh.cg |
Threshold for the termination of the conjugate gradient solver |
nrestarts |
Number of restarts for EM |
max.cores |
Maximum number of cores to use for parallel computation |
verbose |
Print verbose output |
List containing: Predicted labels Yhat and inferred parameters W and beta
Fabian Schmich
Class that represents a netprioR model.
netprioR(networks, phenotypes, labels, ...) ## S4 method for signature 'list,matrix,factor' netprioR(networks, phenotypes, labels, fit.model = FALSE, a = 0.1, b = 0.1, sigma2 = 0.1, tau2 = 100, eps = 1e-10, max.iter = 500, thresh = 1e-06, use.cg = FALSE, thresh.cg = 1e-06, nrestarts = 5, max.cores = detectCores(), verbose = TRUE, ...)
netprioR(networks, phenotypes, labels, ...) ## S4 method for signature 'list,matrix,factor' netprioR(networks, phenotypes, labels, fit.model = FALSE, a = 0.1, b = 0.1, sigma2 = 0.1, tau2 = 100, eps = 1e-10, max.iter = 500, thresh = 1e-06, use.cg = FALSE, thresh.cg = 1e-06, nrestarts = 5, max.cores = detectCores(), verbose = TRUE, ...)
networks |
List of NxN adjacency matrices of gene-gene similarities |
phenotypes |
Matrix of dimension NxP containing covariates |
labels |
Vector of Nx1 labels for all genes (NA if no label available) |
... |
Additional arguments |
fit.model |
Indicator whether to fit the model |
a |
Shape parameter of Gamma prior for W |
b |
Scale parameter of Gamma prior for W |
sigma2 |
Cariance for Gaussian labels |
tau2 |
Variance for Gaussian prior for beta |
eps |
Small value added to diagonal of Q in order to make it non-singular |
max.iter |
Maximum number of iterations for EM |
thresh |
Threshold for termination of EM with respect to change in parameters |
use.cg |
Flag whether to use conjugate gradient instead of exact computation of expectations |
thresh.cg |
Threshold for the termination of the conjugate gradient solver |
nrestarts |
Number of restarts for EM |
max.cores |
Maximum number of cores to use for parallel computation |
verbose |
Print verbose output |
A netprioR
object
networks
List of NxN adjacency matrices of gene-gene similarities
phenotypes
Matrix of dimension NxP containing covariates
labels
Vector of Nx1 labels for all genes. NA if no label available.
is.fitted
Flag indicating if model is fitted
model
List containing estimated parameters and imputed missing data
Fabian Schmich
# runs long-ish data(simulation) np <- netprioR(networks = simulation$networks, phenotypes = simulation$phenotypes, labels = simulation$labels.obs, fit.model = TRUE) summary(np)
# runs long-ish data(simulation) np <- netprioR(networks = simulation$networks, phenotypes = simulation$phenotypes, labels = simulation$labels.obs, fit.model = TRUE) summary(np)
adopted from GeneMania, Mostafavi et al, 2009
norm_kern(x)
norm_kern(x)
x |
kernel |
Normalised kernel
Fabian Schmich
netprioR
objectsPlot method for netprioR
objects
## S3 method for class 'netprioR' plot(x, which = c("all", "weights", "lik", "scores"), ...)
## S3 method for class 'netprioR' plot(x, which = c("all", "weights", "lik", "scores"), ...)
x |
A |
which |
Flag for which plot should be shown, options: weights, lik, scores, all |
... |
Additional paramters for plot |
Plot of the weights, likelihood, ranks, or all three
Fabian Schmich
data(simulation) plot(simulation$model)
data(simulation) plot(simulation$model)
Retrieve ranked prioritisation list
ranks(object) ## S4 method for signature 'netprioR' ranks(object)
ranks(object) ## S4 method for signature 'netprioR' ranks(object)
object |
A |
Ranked list of prioritised genes
Fabian Schmich
data(simulation) ranks(simulation$model)
data(simulation) ranks(simulation$model)
Compute ROC curve from netprioR model and true labels
ROC(object, ...) ## S4 method for signature 'netprioR' ROC(object, true.labels, plot = FALSE, ...)
ROC(object, ...) ## S4 method for signature 'netprioR' ROC(object, true.labels, plot = FALSE, ...)
object |
A |
... |
Additional arguments |
true.labels |
True full set of underlying labels |
plot |
Flag whether to plot the AUC curve |
ROC curve with AUC
Fabian Schmich
data(simulation) ROC(simulation$model, true.labels = simulation$labels.true)
data(simulation) ROC(simulation$model, true.labels = simulation$labels.true)
Simulate labels
simulate_labels(values, sizes, nobs)
simulate_labels(values, sizes, nobs)
values |
Vector of labels for groups |
sizes |
Vector of group sizes |
nobs |
Vector of number of observed labels per group |
List of Y, Yobs and indices for labeled instances
Fabian Schmich
labels <- simulate_labels(values = c("Positive", "Negative"), sizes = c(10, 10), nobs = c(5, 5))
labels <- simulate_labels(values = c("Positive", "Negative"), sizes = c(10, 10), nobs = c(5, 5))
Simulate random networks with predefined number of members for each of the two groups and the number of neighbours for each node
simulate_network_random(nmemb, nnei = 1)
simulate_network_random(nmemb, nnei = 1)
nmemb |
Vector of number of members for each group |
nnei |
Number of neighbours for each node |
Adjacency matrix of graph
Fabian Schmich
network <- simulate_network_random(nmemb = c(10, 10), nnei = 1)
network <- simulate_network_random(nmemb = c(10, 10), nnei = 1)
Simulate scale free networks for predefined number of members for each of two groups and a parameter pclus that determines how strictly distinct the groups are
simulate_network_scalefree(nmemb, pclus = 1)
simulate_network_scalefree(nmemb, pclus = 1)
nmemb |
Vector of numbers of members per group |
pclus |
Scalar in [0, 1] determining how strictly distinct groups are |
Adjacency matrix
Fabian Schmich
network <- simulate_network_scalefree(nmemb = c(10, 10), pclus = 0.8)
network <- simulate_network_scalefree(nmemb = c(10, 10), pclus = 0.8)
Simulate phenotypes correlated to labels pivoted into two groups
simulate_phenotype(labels.true, meandiff, sd)
simulate_phenotype(labels.true, meandiff, sd)
labels.true |
Vector of labels |
meandiff |
difference of means between positive and negative groups |
sd |
Standard deviation of the phenotype |
Simulated phenotype
Fabian Schmich
data(simulation) phenotypes <- simulate_phenotype(labels.true = simulation$labels.true, meandiff = 0.5, sd = 1)
data(simulation) phenotypes <- simulate_phenotype(labels.true = simulation$labels.true, meandiff = 0.5, sd = 1)
The data set contains simulated data for N = 1000 genes and P = 1 (univariate) phenotypes. The list of networks contains 2 low noise networks and two high noise networks. The class labels are "Positive" and "Negative".
data(simulation)
data(simulation)
The code used to simluate the data can be found in system.file("example", "data_simulation.R", package = "netprioR")
List of simulated networks, phenotypes and labels for 1000 genes
Retrieve network weights
weights(object, ...) ## S4 method for signature 'netprioR' weights(object)
weights(object, ...) ## S4 method for signature 'netprioR' weights(object)
object |
A |
... |
Additional arguments |
Estimated network weights
Fabian Schmich
data(simulation) weights(simulation$model)
data(simulation) weights(simulation$model)