View source: R/import_manualSequences.R
importTranscriptSequences | R Documentation |
The given DataFrame (transcripts) should include the following columns:
sample: Names of the samples. (character)
identifier: The identifier which will be used in downstream analysis. (character)
Tx.SequenceMut: The mutant mRNA sequence. (DNAStringSet)
gene: Name of the gene. (character, optional)
The supplied mutant sequences do not need to have a respective record within the used TxDb as these are handled as independent sequences.
importTranscriptSequences( ProteoDiscography, transcripts, removeExisting = FALSE, overwriteDuplicateSamples = FALSE )
ProteoDiscography |
(ProteoDiscography): ProteoDiscography object which stores the annotation and genomic sequences. |
transcripts |
(DataFrame): DataFrame containing metadata and the transcript sequences. |
removeExisting |
(logical): Should previous mutations within the ProteoDiscography be removed? |
overwriteDuplicateSamples |
(logical): Replace duplicate samples (TRUE) or throw an error if duplicate samples are found. |
ProteoDiscography with added mutant transcript sequences from manual input.
Job van Riet j.vanriet@erasmusmc.nl
Wesley van de Geer w.vandegeer@erasmusmc.nl
ProteoDiscography.hg19 <- ProteoDisco::generateProteoDiscography( TxDb = TxDb.Hsapiens.UCSC.hg19.knownGene::TxDb.Hsapiens.UCSC.hg19.knownGene, genomeSeqs = BSgenome.Hsapiens.UCSC.hg19::BSgenome.Hsapiens.UCSC.hg19 ) manualSeq <- S4Vectors::DataFrame( sample = rep('Example', 2), identifier = c('Example 1', 'Example 2'), gene = c('GeneA', 'GeneB'), Tx.SequenceMut = Biostrings::DNAStringSet(list(Biostrings::DNAString('ATCGGGCCCGACGTT'), Biostrings::DNAString('GCTAGCGATCAGGGA'))) ) # Add to ProteoDiscography. ProteoDiscography.hg19 <- ProteoDisco::importTranscriptSequences( ProteoDiscography.hg19, transcripts = manualSeq )
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.