Call copy number for each cell-chromosome using Gaussian mixture models
Source:R/GaussianMixtureCalls.R
calcGMMCopyNumber.Rd
Uses control cells to simulate expected smoothed copy number distributions for all chromosomes across each of model.components
(copy number level).
Then uses the distributions to calculate posterior probabilities for each cell-chromosome belonging to each of copy number level.
Each cell-chromosome is assigned the copy number value for which its posterior probability is highest.
This is done for both whole chromosomes and chromosome arms.
Usage
calcGMMCopyNumber(
TapestriExperiment,
cell.barcodes,
control.copy.number,
model.components = 1:5,
model.priors = NULL,
...
)
Arguments
- TapestriExperiment
TapestriExperiment
object.- cell.barcodes
character, vector of cell barcodes to fit GMM. Usually corresponds to diploid control.
- control.copy.number
data.frame
with columnsarm
andcopy.number
, indicating of known copy number of cells incell.barcodes
.- model.components
numeric, vector of copy number GMM components to calculate, default
1:5
(for copy number = 1, 2, 3, 4, 5).- model.priors
numeric, relative prior probabilities for each GMM component. If
NULL
(default), assumes equal priors.- ...
Additional parameters to be passed to internal functions.
Value
TapestriExperiment
object with copy number calls based on the calculated GMMs, saved to gmmCopyNumber
slot of smoothedCopyNumberByChr
and smoothedCopyNumberByArm
altExps.
GMM parameters for each feature.id
are saved to the metadata
slot.
Examples
# \donttest{
tap.object <- newTapestriExperimentExample() # example TapestriExperiment object
#> ℹ Moving gRNA probe to `altExp` slot "grnaCounts".
#> ℹ Moving barcode probe to `altExp` slot "barcodeCounts".
#> ℹ Moving chrY probe(s) probe_231, probe_232, probe_233, probe_234, probe_235, probe_236, probe_237, probe_238, probe_239, and probe_240 to `altExp` slot "chrYCounts".
tap.object <- calcNormCounts(tap.object)
control.copy.number <- generateControlCopyNumberTemplate(tap.object,
copy.number = 2,
sample.feature.label = "cellline1"
)
tap.object <- calcCopyNumber(tap.object,
control.copy.number,
sample.feature = "test.cluster"
)
tap.object <- calcSmoothCopyNumber(tap.object)
#> ℹ Smoothing copy number by median...
#> ✔ Smoothing copy number by median... [1s]
#>
tap.object <- calcGMMCopyNumber(tap.object,
cell.barcodes = colnames(tap.object),
control.copy.number = control.copy.number,
model.components = 1:5
)
#> ℹ Calculating GMMs using 300 input cells.
#> ℹ Generating probe values for 500 simulated cells...
#> ✔ Generating probe values for 500 simulated cells... [1000ms]
#>
#> ℹ Fitting Gaussian distributions to simulated cells...
#> ✔ Fitting Gaussian distributions to simulated cells... [11.1s]
#>
#> ℹ Calculating posterior probabilities...
#> ✔ Calculating posterior probabilities... [301ms]
#>
#> ℹ Calling copy number from posterior probabilities...
#> ✔ Calling copy number from posterior probabilities... [53ms]
#>
#> ✔ Saving whole chromosome copy number calls to altExp: smoothedCopyNumberByChr,
#> assay: gmmCopyNumber...
#> ✔ Saving chromosome arm copy number calls to altExp: smoothedCopyNumberByArm,
#> assay: gmmCopyNumber...
#> ✔ Saving GMM models and metadata to `gmmParams` slot...
# }