Compute the Automorphisms of Mutational Events Between two Codon Sequences Represented in Z5^3.
Source:R/aut3D.R
aut3D.Rd
Given two codon sequences represented in the Z5^3 Abelian group, this function computes the automorphisms describing codon mutational events.
Usage
aut3D(
seq = NULL,
filepath = NULL,
cube = c("ACGT", "TGCA"),
cube_alt = c("CATG", "GTAC"),
field = "GF5",
start = NA,
end = NA,
chr = 1L,
strand = "+",
genetic_code = getGeneticCode("1"),
num.cores = multicoreWorkers(),
tasks = 0L,
verbose = TRUE
)
Arguments
- seq
An object from a
DNAStringSet
orDNAMultipleAlignment
class carrying the DNA pairwise alignment of two sequences. The pairwise alignment provided in argument seq or the 'fasta' file filepath must correspond to codon sequences.- filepath
A character vector containing the path to a file in fasta format to be read. This argument must be given if codon & base arguments are not provided.
- cube, cube_alt
A character string denoting pairs of the 24 Genetic-code cubes, as given in references (2-3). That is, the base pairs from the given cubes must be complementary each other. Such a cube pair are call dual cubes and, as shown in reference (3), each pair integrates group.
- field
A character string denoting the Galois field where the 3D automorphisms are estimated. This can be 'GF(4)' or 'GF(5)', but only 'GF(5)' is implemented so far.
- start, end, chr, strand
Optional parameters required to build a
GRanges-class
. If not provided the default values given for the function definition will be used.- genetic_code
The named character vector returned by
getGeneticCode
or similar. The translation of codon into aminoacids is a valuable information useful for downstream statistical analysis. The standard genetic code is the default argument value applied in the translation of codons into aminoacids (seeGENETIC_CODE_TABLE
.- num.cores, tasks
Parameters for parallel computation using package
BiocParallel-package
: the number of cores to use, i.e. at most how many child processes will be run simultaneously (seebplapply
and the number of tasks per job (only for Linux OS).- verbose
If TRUE, prints the progress bar.
Value
An object Automorphism-class
with four columns on its
metacolumn named: seq1, seq2, autm, and cube.
Details
Automorphisms in Z5^3' are described as functions \(f(x) = A x mod Z5\), where A is diagonal matrix, as noticed in reference (4).
References
Sanchez R, Morgado E, Grau R. Gene algebra from a genetic code algebraic structure. J Math Biol. 2005 Oct;51(4):431-57. doi: 10.1007/s00285-005-0332-8. Epub 2005 Jul 13. PMID: 16012800. ( PDF).
Robersy Sanchez, Jesus Barreto (2021) Genomic Abelian Finite Groups. https://doi.org/10.1101/2021.06.01.446543.
M. V Jose, E.R. Morgado, R. Sanchez, T. Govezensky, The 24 possible algebraic representations of the standard genetic code in six or in three dimensions, Adv. Stud. Biol. 4 (2012) 119-152.PDF.
R. Sanchez. Symmetric Group of the Genetic-Code Cubes. Effect of the Genetic-Code Architecture on the Evolutionary Process MATCH Commun. Math. Comput. Chem. 79 (2018) 527-560. PDF.
Author
Robersy Sanchez (https://genomaths.com).
Examples
## Load a pairwise alignment
data("aln", package = "GenomAutomorphism")
aln
#> DNAStringSet object of length 2:
#> width seq
#> [1] 51 ACCTATGTTGGTATT---GCGCTCCAACTCCTTGGCTCTAGCTCACTACAT
#> [2] 51 ATCTATGTTGGTATTACGACGCTCCAATTCCTTGGGTCC------CTCCTT
## Automorphism on Z5^3
autms <- aut3D(seq = aln)
#>
|
| | 0%
|
|======= | 10%
|
|============== | 20%
|
|===================== | 30%
|
|============================ | 40%
|
|=================================== | 50%
|
|========================================== | 60%
|
|================================================= | 70%
|
|======================================================== | 80%
|
|=============================================================== | 90%
|
|======================================================================| 100%
#>
autms
#> Automorphism object with 17 ranges and 8 metadata columns:
#> seqnames ranges strand | seq1 seq2 aa1
#> <Rle> <IRanges> <Rle> | <character> <character> <character>
#> [1] 1 1 + | ACC ATC T
#> [2] 1 2 + | TAT TAT Y
#> [3] 1 3 + | GTT GTT V
#> [4] 1 4 + | GGT GGT G
#> [5] 1 5 + | ATT ATT I
#> ... ... ... ... . ... ... ...
#> [13] 1 13 + | TCT TCC S
#> [14] 1 14 + | AGC --- S
#> [15] 1 15 + | TCA --- S
#> [16] 1 16 + | CTA CTC L
#> [17] 1 17 + | CAT CTT H
#> aa2 coord1 coord2 autm cube
#> <character> <matrix> <matrix> <character> <character>
#> [1] I 1:2:2 1:4:2 1,2,1 ACGT
#> [2] Y 4:1:4 4:1:4 1,1,1 ACGT
#> [3] V 3:4:4 3:4:4 1,1,1 ACGT
#> [4] G 3:3:4 3:3:4 1,1,1 ACGT
#> [5] I 1:4:4 1:4:4 1,1,1 ACGT
#> ... ... ... ... ... ...
#> [13] S 4:2:4 4:2:2 1,1,3 ACGT
#> [14] - 1:3:2 0:0:0 0 Trnl
#> [15] - 4:2:1 0:0:0 0 Trnl
#> [16] L 2:4:1 2:4:2 1,1,2 ACGT
#> [17] L 2:1:4 2:4:4 1,4,1 ACGT
#> -------
#> seqinfo: 1 sequence from an unspecified genome; no seqlengths