ScpModel | R Documentation |
An ScpModel
object must be always stored in the metadata()
of
an object that inherits from the SingleCellExperiment
class. The
ScpModel
object should never be accessed directly by the
user. Instead, we provide several setter function to retrieve
information that may be useful to the user.The ScpModel
class
contains several slots:
scpModelFormula
: a formula
object controlling which
variables are to be modelled.
scpModelInputIndex
: a numeric(1)
, selecting the assay to use
in the SingleCellExperiment
object as input matrix. Note that
this slot serves as a pointer, meaning that the quantitative
data is not duplicated. Any change to the assay in the
SingleCellExperiment
will impact the estimation of the
ScpModel
object.
scpModelFilterThreshold
: A numeric(1)
indicating the minimal
n/p ratio required for a feature to be included in further model
exploration. n/p is the number of measured values for a features
divided by the number of coefficients to estimate. n/p cannot be
smaller than 1 because this would lead to over-specified models.
scpModelFitList
: A List
that contains the model results for
each feature. Each element is a ScpModelFit
object (see
ScpModelFit
)
scpModelFormula(object, name)
scpModelInput(object, name, filtered = TRUE)
scpModelFilterThreshold(object, name)
scpModelFilterNPRatio(object, name, filtered = TRUE)
scpModelResiduals(object, name, join = TRUE, filtered = TRUE)
scpModelEffects(object, name, join = TRUE, filtered = TRUE)
scpModelNames(object)
scpModelFilterThreshold(object, name) <- value
object |
An object that inherits from the
|
name |
A |
filtered |
A |
join |
A |
value |
An |
Each slot has a getter function associated:
scpModelNames()
: returns a vector of names of ScpModel
objects stored in the SingleCellExperiment
object.
scpModelFormula()
: returns the formula
slot of the ScpModel
within an object that inherits from the SummarizedExperiment
class.
scpModelFilterThreshold()
: returns the n/p ration threshold
used for feature filtering.
scpModelInput()
: returns a matrix
with the quantitative
values used as input of the model. Hence, the matrix contains
the data before modelling. If filtered = TRUE
, the feature of
the matrix are restricted to the features that satisfy the n/p
ratio threshold.
scpModelFilterNPRatio()
: returns the computed n/p ratio for
each feature. If filtered = TRUE
, the function returns only
the n/p of the features that satisfy the n/p ratio threshold.
scpModelResiduals()
: when join = FALSE
, the function returns
a list where each element corresponds to a feature and contains
the estimated residuals. When join = TRUE
(default), the function
combines the list into a matrix with features in rows and cells
in columns, and filling the gaps with NA
. If filtered = TRUE
,
the feature of the matrix are restricted to the features that
satisfy the n/p ratio threshold.
scpModelEffects()
: when join = FALSE
, the function return a
list where each element of the list corresponds to a feature.
Each element contains another list with as many elements as
variable in the model and each element contains the data effect vector
for that vector. When join = TRUE
(default), each element of the list is
a matrix with features in rows and cells in columns where gaps
are filled with NA
. If filtered = TRUE
, the feature of the
matrix are restricted to the features that satisfy the n/p
ratio threshold.
Setter:
scpModelFilterThreshold<-()
: the function changes the n/p
ratio threshold used for filtering features.
Christophe Vanderaa, Laurent Gatto
ScpModelFit for a description of the class that store modelling results
ScpModel-Workflow that uses the class to store the estimated model.
data("leduc_minimal")
####---- Getters ----####
scpModelNames(leduc_minimal)
scpModelFormula(leduc_minimal)
dim(leduc_minimal)
dim(scpModelInput(leduc_minimal))
dim(scpModelInput(leduc_minimal, filtered = FALSE))
head(scpModelFilterNPRatio(leduc_minimal))
dim(scpModelResiduals(leduc_minimal))
dim(scpModelResiduals(leduc_minimal, filtered = FALSE))
scpModelResiduals(leduc_minimal, join = FALSE)
scpModelEffects(leduc_minimal)
dim(scpModelEffects(leduc_minimal)$Set)
dim(scpModelEffects(leduc_minimal, filtered = FALSE)$Set)
scpModelEffects(leduc_minimal, join = FALSE)[[1]]
scpModelFilterThreshold(leduc_minimal)
####---- Setter ----####
scpModelFilterThreshold(leduc_minimal) <- 2
scpModelFilterThreshold(leduc_minimal)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.