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.

import_sq(object, ...)

Arguments

object

[any(1)]
An object of one of supported classes.

...

further arguments to be passed from or to other methods.

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

Details

Currently supported classes are as follows:

  • ape:

    • AAbin - imported as ami_bsc

    • DNAbin - imported as dna_bsc

    • alignment - exact type is guessed within sq 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 separate tibble, resulting in a list of tibbles; if passed argument separate = 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.

See also

sq class

Functions from input module: random_sq(), read_fasta(), sq()

Examples

# ape example
library(ape)
#> 
#> Attaching package: ‘ape’
#> The following object is masked from ‘package:tidysq’:
#> 
#>     complement
ape_dna <- as.DNAbin(list(one = c("C", "T", "C", "A"), two = c("T", "G", "A", "G", "G")))
import_sq(ape_dna)
#> # A tibble: 2 × 2
#>   name  sq       
#>   <chr> <dna_bsc>
#> 1 one   CTCA  <4>
#> 2 two   TGAGG <5>

# bioseq example
library(bioseq)
#> 
#> Attaching package: ‘bioseq’
#> The following objects are masked from ‘package:tidysq’:
#> 
#>     read_fasta, write_fasta
bioseq_rna <- new_rna(c(one = "ANBRY", two = "YUTUGGN"))
#> Warning: Non-standard IUPAC symbols detected for RNA: 1 characters were converted to N.
import_sq(bioseq_rna)
#> # A tibble: 2 × 2
#>   name  sq         
#>   <chr> <rna_ext>  
#> 1 one   ANBRY   <5>
#> 2 two   YUNUGGN <7>

# Biostrings example
library(Biostrings)
#> Loading required package: BiocGenerics
#> Loading required package: parallel
#> 
#> Attaching package: ‘BiocGenerics’
#> The following objects are masked from ‘package:parallel’:
#> 
#>     clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
#>     clusterExport, clusterMap, parApply, parCapply, parLapply,
#>     parLapplyLB, parRapply, parSapply, parSapplyLB
#> The following objects are masked from ‘package:dplyr’:
#> 
#>     combine, intersect, setdiff, union
#> The following object is masked from ‘package:tidysq’:
#> 
#>     paste
#> The following objects are masked from ‘package:stats’:
#> 
#>     IQR, mad, sd, var, xtabs
#> The following objects are masked from ‘package:base’:
#> 
#>     Filter, Find, Map, Position, Reduce, anyDuplicated, append,
#>     as.data.frame, basename, cbind, colnames, dirname, do.call,
#>     duplicated, eval, evalq, get, grep, grepl, intersect, is.unsorted,
#>     lapply, mapply, match, mget, order, paste, pmax, pmax.int, pmin,
#>     pmin.int, rank, rbind, rownames, sapply, setdiff, sort, table,
#>     tapply, union, unique, unsplit, which.max, which.min
#> Loading required package: S4Vectors
#> Loading required package: stats4
#> 
#> Attaching package: ‘S4Vectors’
#> The following objects are masked from ‘package:dplyr’:
#> 
#>     first, rename
#> The following objects are masked from ‘package:base’:
#> 
#>     I, expand.grid, unname
#> Loading required package: IRanges
#> 
#> Attaching package: ‘IRanges’
#> The following objects are masked from ‘package:dplyr’:
#> 
#>     collapse, desc, slice
#> The following objects are masked from ‘package:tidysq’:
#> 
#>     collapse, reverse
#> Loading required package: XVector
#> Loading required package: GenomeInfoDb
#> 
#> Attaching package: ‘Biostrings’
#> The following object is masked from ‘package:ape’:
#> 
#>     complement
#> The following objects are masked from ‘package:tidysq’:
#> 
#>     alphabet, complement, translate
#> The following object is masked from ‘package:base’:
#> 
#>     strsplit
Biostrings_ami <- AAStringSet(c(one = "FEAPQLIWY", two = "EGITENAK"))
import_sq(Biostrings_ami)
#> # A tibble: 2 × 2
#>   name  sq           
#>   <chr> <ami_ext>    
#> 1 one   FEAPQLIWY <9>
#> 2 two   EGITENAK  <8>

# seqinr example
library(seqinr)
#> 
#> Attaching package: ‘seqinr’
#> The following object is masked from ‘package:Biostrings’:
#> 
#>     translate
#> The following objects are masked from ‘package:ape’:
#> 
#>     as.alignment, consensus
#> The following object is masked from ‘package:dplyr’:
#> 
#>     count
#> The following object is masked from ‘package:tidysq’:
#> 
#>     translate
seqinr_dna <- as.SeqFastadna(c("C", "T", "C", "A"), name = "one")
import_sq(seqinr_dna)
#> # A tibble: 1 × 2
#>   name  sq       
#>   <chr> <dna_bsc>
#> 1 one   CTCA <4>