cv.glmSparseNet: Calculate cross validating GLM model with network-based...

View source: R/glm_sparse_net.R

cv.glmDegreeR Documentation

Calculate cross validating GLM model with network-based regularization

Description

network parameter accepts:

Usage

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

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

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

cv.glmSparseNet(
  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::cv.glmnet() accepts.

Details

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

  • matrix representing the network

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

Value

an object just as cv.glmnet

Functions

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

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

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

See Also

Model with the same penalizations glmSparseNet().

Examples

# Degree penalization

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

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

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

# Gaussian model
xdata <- matrix(rnorm(500), ncol = 5)
cv.glmSparseNet(
    xdata, rnorm(nrow(xdata)), "correlation",
    family = "gaussian"
)
cv.glmSparseNet(
    xdata, rnorm(nrow(xdata)), "covariance",
    family = "gaussian"
)


#
#
# Using MultiAssayExperiment with survival model
library(MultiAssayExperiment)
data("miniACC", package = "MultiAssayExperiment")

xdata <- miniACC

#
# 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")

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


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