glmSparseNet: Calculate GLM model with network-based regularization

View source: R/glm_sparse_net.R

glmSparseNetR Documentation

Calculate GLM model with network-based regularization

Description

network parameter accepts:

  • string to calculate network based on data (correlation, covariance)

  • matrix representing the network

  • vector with already calculated penalty weights (can also be used directly with glmnet)

Usage

glmSparseNet(
  xdata,
  ydata,
  network,
  options = networkOptions(),
  experiment = NULL,
  network.options = deprecated(),
  experiment.name = deprecated(),
  ...
)

glmDegree(
  xdata,
  ydata,
  network,
  options = networkOptions(),
  experiment = NULL,
  network.options = deprecated(),
  experiment.name = deprecated(),
  ...
)

glmHub(
  xdata,
  ydata,
  network,
  options = networkOptions(),
  experiment = NULL,
  network.options = deprecated(),
  experiment.name = deprecated(),
  ...
)

glmOrphan(
  xdata,
  ydata,
  network,
  options = networkOptions(),
  experiment = NULL,
  network.options = deprecated(),
  experiment.name = deprecated(),
  ...
)

Arguments

xdata

input data, can be a matrix or MultiAssayExperiment.

ydata

response data compatible with glmnet.

network

type of network, see below.

options

options to calculate network.

experiment

name of experiment to use as input in MultiAssayExperiment object (only if xdata is an object of this class).

network.options

[Deprecated]

experiment.name

[Deprecated]

...

parameters that glmnet::glmnet() accepts.

Value

an object just as glmnet

Functions

  • glmDegree(): penalizes nodes with small degree (inversion penalization h(x) = 1 / x).

  • glmHub(): Penalizes nodes with small degree (normalized heuristic that promotes nodes with many edges).

  • glmOrphan(): Penalizes nodes with high degree (normalized heuristic that promotes nodes with few edges).

See Also

Cross-validation functions cv.glmSparseNet().

Examples

xdata <- matrix(rnorm(100), ncol = 20)
glmSparseNet(xdata, rnorm(nrow(xdata)), "correlation", family = "gaussian")
glmSparseNet(xdata, rnorm(nrow(xdata)), "covariance", family = "gaussian")


#
#
# Using MultiAssayExperiment
# load data
library(MultiAssayExperiment)
data("miniACC", package = "MultiAssayExperiment")

xdata <- miniACC
# TODO aking out x individuals missing values
# build valid data with days of last follow up or to event
event.ix <- which(!is.na(xdata$days_to_death))
cens.ix <- which(!is.na(xdata$days_to_last_followup))

xdata$surv_event_time <- array(NA, nrow(colData(xdata)))
xdata$surv_event_time[event.ix] <- xdata$days_to_death[event.ix]
xdata$surv_event_time[cens.ix] <- xdata$days_to_last_followup[cens.ix]

# Keep only valid individuals
valid.ix <- as.vector(!is.na(xdata$surv_event_time) &
    !is.na(xdata$vital_status) &
    xdata$surv_event_time > 0)
xdata.valid <- xdata[, rownames(colData(xdata))[valid.ix]]
ydata.valid <- colData(xdata.valid)[, c("surv_event_time", "vital_status")]
colnames(ydata.valid) <- c("time", "status")

glmSparseNet(
    xdata.valid,
    ydata.valid,
    family = "cox",
    network = "correlation",
    experiment = "RNASeq2GeneNorm"
)

# Degree penalization

xdata <- matrix(rnorm(100), ncol = 5)
glmDegree(
    xdata,
    rnorm(nrow(xdata)),
    "correlation",
    family = "gaussian",
    options = networkOptions(minDegree = .2)
)
xdata <- matrix(rnorm(100), ncol = 5)
glmHub(
    xdata,
    rnorm(nrow(xdata)),
    "correlation",
    family = "gaussian",
    options = networkOptions(minDegree = .2)
)
# Orphan penalization

xdata <- matrix(rnorm(100), ncol = 5)
glmOrphan(
    xdata,
    rnorm(nrow(xdata)),
    "correlation",
    family = "gaussian",
    options = networkOptions(minDegree = .2)
)

sysbiomed/glmSparseNet documentation built on Feb. 17, 2024, 1:38 p.m.