estimateZinbwaveParams {digitalDLSorteR} | R Documentation |
Estimate the parameters of the ZINB-WaVE model to simulate new single-cell RNA-Seq expression profiles
Description
Estimate the parameters of the ZINB-WaVE model using a real single-cell
RNA-Seq data set as reference to simulate new single-cell profiles and
increase the signal of underrepresented cell types. This step is optional,
only is needed if the size of you dataset is too small or there are
underrepresented cell types in order to train the Deep Neural Network model
in a more balanced way. After this step, the simSCProfiles
function will use the estimated parameters to simulate new single-cell
profiles. See ?simSCProfiles
for more information.
Usage
estimateZinbwaveParams(
object,
cell.type.column,
cell.ID.column,
gene.ID.column,
cell.cov.columns,
gene.cov.columns,
subset.cells = NULL,
proportional = TRUE,
set.type = "All",
threads = 1,
verbose = TRUE
)
Arguments
object |
|
cell.type.column |
Name or column number corresponding to the cell type of each cell in cells metadata. |
cell.ID.column |
Name or column number corresponding to the cell names of expression matrix in cells metadata. |
gene.ID.column |
Name or column number corresponding to the notation used for features/genes in genes metadata. |
cell.cov.columns |
Name or column number(s) in cells metadata to be used
as covariates during model fitting (if no covariates are used, set to empty
or |
gene.cov.columns |
Name or column number(s) in genes metadata that will
be used as covariates during model fitting (if no covariates are used, set
to empty or |
subset.cells |
Number of cells to fit the ZINB-WaVE model. Useful when
the original data set is too large to fit the model. Set a value according
to the original data set and the resources available on your computer. If
|
proportional |
If |
set.type |
Cell type(s) to evaluate ( |
threads |
Number of threads used for estimation (1 by default). To set up the parallel environment, the BiocParallel package must be installed. |
verbose |
Show informative messages during the execution ( |
Details
ZINB-WaVE is a flexible model for zero-inflated count data. This function
carries out the model fit to real single-cell data modeling Y_{ij}
(the
count of feature j
for sample i
) as a random variable following a
zero-inflated negative binomial (ZINB) distribution. The estimated parameters
will be used for the simulation of new single-cell expression profiles by
sampling a negative binomial distribution and inserting dropouts from a
binomial distribution. To do so, digitalDLSorteR uses the
zinbFit
function from the zinbwave package
(Risso et al., 2018). For more details about the model, see Risso et al.,
2018.
Value
A DigitalDLSorter
object with zinb.params
slot containing a ZinbParametersModel
object. This
object contains a slot with the estimated ZINB-WaVE parameters from the
real single-cell RNA-Se'q data.
References
Risso, D., Perraudeau, F., Gribkova, S. et al. (2018). A general and flexible method for signal extraction from single-cell RNA-seq data. Nat Commun 9, 284. doi: doi:10.1038/s41467-017-02554-5.
Torroja, C. and Sánchez-Cabo, F. (2019). digitalDLSorter: A Deep Learning algorithm to quantify immune cell populations based on scRNA-Seq data. Frontiers in Genetics 10, 978. doi: doi:10.3389/fgene.2019.00978.
See Also
Examples
set.seed(123) # reproducibility
sce <- SingleCellExperiment::SingleCellExperiment(
assays = list(
counts = matrix(
rpois(30, lambda = 5), nrow = 15, ncol = 10,
dimnames = list(paste0("Gene", seq(15)), paste0("RHC", seq(10)))
)
),
colData = data.frame(
Cell_ID = paste0("RHC", seq(10)),
Cell_Type = sample(x = paste0("CellType", seq(2)), size = 10,
replace = TRUE)
),
rowData = data.frame(
Gene_ID = paste0("Gene", seq(15))
)
)
DDLS <- createDDLSobject(
sc.data = sce,
sc.cell.ID.column = "Cell_ID",
sc.gene.ID.column = "Gene_ID",
sc.filt.genes.cluster = FALSE,
sc.log.FC = FALSE
)
DDLS <- estimateZinbwaveParams(
object = DDLS,
cell.type.column = "Cell_Type",
cell.ID.column = "Cell_ID",
gene.ID.column = "Gene_ID",
subset.cells = 2,
verbose = TRUE
)