View source: R/export_exportProteoDiscography.R
exportProteoDiscography | R Documentation |
Exports (mutant) peptide sequences to FASTA using several fields as identifier. In addition, users can specify to only export (mutant) sequences containing at a min. number of proteotypic fragments.
exportProteoDiscography( ProteoDiscography, outFile = NULL, minProteotypicFragments = NULL, aggregateSamples = TRUE )
ProteoDiscography |
(ProteoDiscography): ProteoDiscography object which stores the annotation and genomic sequences. |
outFile |
(character): Filepath to output FASTA, will append the samplenames if aggregateSamples is FALSE. If left NULL, it will return an AAStringSet with the given records. |
minProteotypicFragments |
(integer): Only output mutant protein-isoforms from incorporated genomic variants with at least this many proteotypic fragments (if 'checkProteotypicFragments()' was performed). |
aggregateSamples |
(logical): Should samples be aggregated into the same output database (TRUE) or should a seperate FASTA-file be generated per sample (FALSE). |
Writes a FASTA file containing the mutant protein isoforms if outFile is given. Otherwise, will return an AAStringSet.
Job van Riet j.vanriet@erasmusmc.nl
Wesley van de Geer w.vandegeer@erasmusmc.nl
# Import example ProteoDiscography (hg19) data('ProteoDiscographyExample.hg19', package = 'ProteoDisco') ProteoDiscographyExample.hg19 <- setTxDb(ProteoDiscographyExample.hg19, TxDb.Hsapiens.UCSC.hg19.knownGene::TxDb.Hsapiens.UCSC.hg19.knownGene) ProteoDiscographyExample.hg19 <- setGenomicSequences(ProteoDiscographyExample.hg19, BSgenome.Hsapiens.UCSC.hg19::BSgenome.Hsapiens.UCSC.hg19) # Export peptide sequences to FASTA file. Optionally, only export those with at least a min. # number of proteotypic fragments. exportProteoDiscography(ProteoDiscographyExample.hg19, outFile = 'out.fasta')
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.