Estimate parameters to be used to simulate repeated measures scRNA-Seq data
Usage
estRescueSimParams(
sce,
paramObj = NULL,
sampleVariable = NULL,
subjectVariable = NULL,
groupVariable = NULL,
timepointVariable = NULL,
nonDEGs = NULL,
cellParamsByCondition = FALSE
)Arguments
- sce
SingleCellExperimentobject containing empirical data with a counts matrix in thecountsslot and cell level data (ex. sample/subject/timepoint labels) in thecolDataslot.- paramObj
RescueSimParams-classobject with empty slots for parameters which need to be estimated and parameters. IfNULL, all parameter values will be estimated.- sampleVariable
String denoting name of sample identifier variable in the
colDataofSingleCellExperimentobject. IfNULL, data is assumed to contain only one sample parameterssampleFacVarMeanandsampleFacVarSDparameters cannot be estimated.- subjectVariable
String denoting name of subject identifier variable in the
colDataofSingleCellExperimentobject. IfNULL, parameterssubjectFacVarMeanandsubjectFacVarSDparameters cannot be estimated.- groupVariable
String denoting name of group identifier variable in the
colDataofSingleCellExperimentobject. IfNULL, a single group design is assumed.- timepointVariable
String denoting name of timepoint identifier variable in the
colDataofSingleCellExperimentobject. IfNULL, the number of timepoints for the simulated data is set to 1.- nonDEGs
Vector containing names or row indices of genes to be used to estimate sample and subject factor parameters. Using genes that are differentially expressed across timepoints may lead to inaccurate estimates.
- cellParamsByCondition
Logical value indicating whether
maxCellsPerSampandminCellsPerSampshould be estimated for each condition (group/timepoint) or overall.
Details
All parameters in RescueSimParams-class object can be
estimated/extracted from empirical data except propDE and
deLog2FC which are both set to 0 (no differential expression) if not
set manually. If a paramObj is provided, only parameters with empty slots in
the object will be estimated.
Examples
## Load the data
data("RecAM_sce")
## Subset data to 40 genes
RecAM_sce_small <- RecAM_sce[1:40, ]
## Estimate all parameters from the data
myParams_estimated <- estRescueSimParams(
sce = RecAM_sce_small, paramObj = NULL,
sampleVariable = "sampleID",
subjectVariable = "subjectID",
timepointVariable = "time",
groupVariable = NULL, nonDEGs = NULL,
cellParamsByCondition = FALSE
)
## Example using Single timepoint data
## Subset the data to a single timepoint
RecAM_sce_single_time <- RecAM_sce_small[, RecAM_sce_small$time == "time0"]
## Create a params object with necessary parameters filled
myParams <- RescueSimParams(subjectFacVarMean = -4, subjectFacVarSD = 1, nSubjsPerGroup = 5)
## Estimate remaining parameters from the data
myParams_estimated <- estRescueSimParams(
sce = RecAM_sce_single_time,
paramObj = myParams,
sampleVariable = "sampleID"
)
#> Warning: Since subjectVariable is NULL, subjectFacVarMean and subjectFacVarSD cannot be estimated.