Converts object of class sq
to a class
from another package. Currently supported packages are ape,
bioseq, Bioconductor and seqinr. For exact list of
supported classes and resulting types, see details.
export_sq(x, export_format, name = NULL, ...)
[sq
]
An object this function is applied to.
[character(1)
]
A string indicating desired class (with specified package for clarity).
[character
]
Vector of sequence names. Must be of the same length as sq
object.
Can be NULL
.
further arguments to be passed from or to other methods.
An object with the format specified in the parameter. To find information about the detailed structure of this object, see documentation of these objects.
Currently supported formats are as follows (grouped by sq
types):
ami:
"ape::AAbin"
"bioseq::bioseq_aa"
"Biostrings::AAString"
"Biostrings::AAStringSet"
"seqinr::SeqFastaAA"
dna:
"ape::DNAbin"
"bioseq::bioseq_dna"
"Biostrings::DNAString"
"Biostrings::DNAStringSet"
"seqinr::SeqFastadna"
rna:
"bioseq::bioseq_rna"
"Biostrings::RNAString"
"Biostrings::RNAStringSet"
Functions from output module:
as.character.sq()
,
as.matrix.sq()
,
as.sq()
,
write_fasta()
# DNA and amino acid sequences can be exported to most packages
sq_ami <- sq(c("MVVGL", "LAVPP"), alphabet = "ami_bsc")
export_sq(sq_ami, "ape::AAbin")
#> 2 amino acid sequences in a list
#>
#> All sequences of the same length: 5
#>
export_sq(sq_ami, "bioseq::bioseq_aa")
#> AA vector of 2 sequences
#> > MVVGL
#> > LAVPP
export_sq(sq_ami, "Biostrings::AAStringSet", c("one", "two"))
#> AAStringSet object of length 2:
#> width seq names
#> [1] 5 MVVGL one
#> [2] 5 LAVPP two
export_sq(sq_ami, "seqinr::SeqFastaAA")
#> [[1]]
#> [1] "M" "V" "V" "G" "L"
#> attr(,"class")
#> [1] "SeqFastaAA"
#>
#> [[2]]
#> [1] "L" "A" "V" "P" "P"
#> attr(,"class")
#> [1] "SeqFastaAA"
#>
sq_dna <- sq(c("TGATGAAGCGCA", "TTGATGGGAA"), alphabet = "dna_bsc")
export_sq(sq_dna, "ape::DNAbin", name = c("one", "two"))
#> 2 DNA sequences in binary format stored in a list.
#>
#> Mean sequence length: 11
#> Shortest sequence: 10
#> Longest sequence: 12
#>
#> Labels:
#> one
#> two
#>
#> Base composition:
#> a c g t
#> 0.318 0.091 0.364 0.227
#> (Total: 22 bases)
export_sq(sq_dna, "bioseq::bioseq_dna")
#> DNA vector of 2 sequences
#> > TGATGAAGCGCA
#> > TTGATGGGAA
export_sq(sq_dna, "Biostrings::DNAStringSet")
#> DNAStringSet object of length 2:
#> width seq
#> [1] 12 TGATGAAGCGCA
#> [2] 10 TTGATGGGAA
export_sq(sq_dna, "seqinr::SeqFastadna")
#> [[1]]
#> [1] "T" "G" "A" "T" "G" "A" "A" "G" "C" "G" "C" "A"
#> attr(,"class")
#> [1] "SeqFastadna"
#>
#> [[2]]
#> [1] "T" "T" "G" "A" "T" "G" "G" "G" "A" "A"
#> attr(,"class")
#> [1] "SeqFastadna"
#>
# RNA sequences are limited to Biostrings and bioseq
sq_rna <- sq(c("NUARYGCB", "", "DRKCNYBAU"), alphabet = "rna_ext")
export_sq(sq_rna, "bioseq::bioseq_rna")
#> RNA vector of 3 sequences
#> > NUARYGCB
#> >
#> > DRKCNYBAU
export_sq(sq_rna, "Biostrings::RNAStringSet")
#> RNAStringSet object of length 3:
#> width seq
#> [1] 8 NUARYGCB
#> [2] 0
#> [3] 9 DRKCNYBAU
# Biostrings can export single sequences to simple strings as well
export_sq(sq_dna[1], "Biostrings::DNAString")
#> 12-letter DNAString object
#> seq: TGATGAAGCGCA