View source: R/karyotypeMeasures.R
karyotypeMeasures | R Documentation |
Computes measures for karyotype heterogeneity. See the Details section for how these measures are defined.
karyotypeMeasures(hmms, normalChromosomeNumbers = NULL, regions = NULL,
exclude.regions = NULL)
hmms |
A list with |
normalChromosomeNumbers |
A named integer vector or matrix with physiological copy numbers, where each element (vector) or column (matrix) corresponds to a chromosome. This is useful to specify male or female samples, e.g. |
regions |
A |
exclude.regions |
A |
We define x
as the vector of copy number states for each position. The number of HMMs is S
. The measures are computed for each bin as follows:
D = mean( abs(x-P) )
, where P is the physiological number of chromosomes at that position.
H = sum( table(x) * 0:(length(table(x))-1) ) / S
A list
with two data.frame
s, containing the karyotype measures $genomewide and $per.chromosome. If region
was specified, a third list entry $regions will contain the regions with karyotype measures.
Aaron Taudt
### Example 1 ###
## Get results from a small-cell-lung-cancer
lung.folder <- system.file("extdata", "primary-lung", "hmms", package="AneuFinderData")
lung.files <- list.files(lung.folder, full.names=TRUE)
## Get results from the liver metastasis of the same patient
liver.folder <- system.file("extdata", "metastasis-liver", "hmms", package="AneuFinderData")
liver.files <- list.files(liver.folder, full.names=TRUE)
## Compare karyotype measures between the two cancers
normal.chrom.numbers <- rep(2, 23)
names(normal.chrom.numbers) <- c(1:22,'X')
lung <- karyotypeMeasures(lung.files, normalChromosomeNumbers=normal.chrom.numbers)
liver <- karyotypeMeasures(liver.files, normalChromosomeNumbers=normal.chrom.numbers)
print(lung$genomewide)
print(liver$genomewide)
### Example 2 ###
## Construct a matrix with physiological copy numbers for a mix of 5 male and 5 female samples
normal.chrom.numbers <- matrix(2, nrow=10, ncol=24,
dimnames=list(sample=c(paste('male', 1:5), paste('female', 6:10)),
chromosome=c(1:22,'X','Y')))
normal.chrom.numbers[1:5,c('X','Y')] <- 1
normal.chrom.numbers[6:10,c('Y')] <- 0
print(normal.chrom.numbers)
### Example 3 ###
## Exclude artifact regions with high variance
consensus <- consensusSegments(c(lung.files, liver.files))
variance <- apply(consensus$copy.number, 1, var)
exclude.regions <- consensus[variance > quantile(variance, 0.999)]
## Compare karyotype measures between the two cancers
normal.chrom.numbers <- rep(2, 23)
names(normal.chrom.numbers) <- c(1:22,'X')
lung <- karyotypeMeasures(lung.files, normalChromosomeNumbers=normal.chrom.numbers,
exclude.regions = exclude.regions)
liver <- karyotypeMeasures(liver.files, normalChromosomeNumbers=normal.chrom.numbers,
exclude.regions = exclude.regions)
print(lung$genomewide)
print(liver$genomewide)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.