CallCNVs,ExomeDepth-method | R Documentation |
Call CNV data from an ExomeDepth object.
## S4 method for signature 'ExomeDepth' CallCNVs( x, chromosome, start, end, name, transition.probability = 1e-04, expected.CNV.length = 50000 )
x |
An |
chromosome |
Chromosome information for each exon (factor). |
start |
Start (physical position) of each exon (numeric, must have the same length as the chromosome argument). |
end |
End (physical position) of each exome (numeric, must have the same length as the chromosome argument). |
name |
Name of each exon (character or factor). |
transition.probability |
Transition probability of the hidden Markov Chain from the normal copy number state to either a deletion or a duplication. The default (0.0001) expect approximately 20 CNVs genome-wide. |
expected.CNV.length |
The expectation for the length of a CNV. This value factors into the Viterbi algorithm that is used to compte the transition from one state to the next, which depends on the distance between exons. |
Must be called on an ExomeDepth object and fits a hidden Markov model to the read depth data with three hidden states (normal, deletion, duplication). Likelihood data must have been pre-computed which should have been done by default when the ExomeDepth object was created.
The same ExomeDepth object provided as input but with the slot CNVcalls containing a data frame with the output of the calling.
data(ExomeCount) ExomeCount <- ExomeCount[1:500,] ## small for the purpose of this test ref_counts <- ExomeCount$Exome2 + ExomeCount$Exome3 + ExomeCount$Exome4 ## creates a simple ExomeDepth object test_object <- new('ExomeDepth', test = ExomeCount$Exome1, reference = ref_counts) ## Call CNVs called_object <- CallCNVs(x = test_object, transition.probability = 10^-4, chromosome = GenomicRanges::seqnames(ExomeCount), start = GenomicRanges::start(ExomeCount), end = GenomicRanges::end(ExomeCount), name = ExomeCount$names) print(called_object@CNV.calls)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.