import_sq {tidysq} | R Documentation |
Import sq objects from other objects
Description
Creates sq
object from object of class
from another package. Currently supported packages are ape,
bioseq, Bioconductor and seqinr. For exact list of
supported classes and resulting types, see details.
Usage
import_sq(object, ...)
Arguments
object |
[ |
... |
further arguments to be passed from or to other methods. |
Details
Currently supported classes are as follows:
-
ape
:-
AAbin
- imported as ami_bsc -
DNAbin
- imported as dna_bsc -
alignment
- exact type is guessed withinsq
function
-
-
bioseq
:-
bioseq_aa
- imported as ami_ext -
bioseq_dna
- imported as dna_ext -
bioseq_rna
- imported as rna_ext
-
-
Biostrings
:-
AAString
- imported as ami_ext with exactly one sequence -
AAStringSet
- imported as ami_ext -
DNAString
- imported as dna_ext with exactly one sequence -
DNAStringSet
- imported as dna_ext -
RNAString
- imported as rna_ext with exactly one sequence -
RNAStringSet
- imported as rna_ext -
BString
- imported as unt with exactly one sequence -
BStringSet
- imported as unt -
XStringSetList
- each element of a list can be imported as a separatetibble
, resulting in a list of tibbles; if passed argumentseparate = FALSE
, these tibbles are bound into one bigger tibble
-
-
seqinr
:-
SeqFastaAA
- imported as ami_bsc -
SeqFastadna
- imported as dna_bsc
-
Providing object of class other than specified will result in an error.
Value
A tibble
with sq
column of
sq
type representing the same sequences as given
object; the object has a type corresponding to the input type; if given
sequences have names, output tibble
will also have
another column name
with those names
See Also
Functions from input module:
random_sq()
,
read_fasta()
,
sq()
Examples
# ape example
library(ape)
ape_dna <- as.DNAbin(list(one = c("C", "T", "C", "A"), two = c("T", "G", "A", "G", "G")))
import_sq(ape_dna)
# bioseq example
library(bioseq)
bioseq_rna <- new_rna(c(one = "ANBRY", two = "YUTUGGN"))
import_sq(bioseq_rna)
# Biostrings example
library(Biostrings)
Biostrings_ami <- AAStringSet(c(one = "FEAPQLIWY", two = "EGITENAK"))
import_sq(Biostrings_ami)
# seqinr example
library(seqinr)
seqinr_dna <- as.SeqFastadna(c("C", "T", "C", "A"), name = "one")
import_sq(seqinr_dna)