codon_coord | R Documentation |
Given a string denoting a codon or base from the DNA (or RNA) alphabet and a genetic-code Abelian group as given in reference (1).
codon_coord(codon = NULL, ...)
## S4 method for signature 'BaseGroup'
codon_coord(codon, group = NULL)
## S4 method for signature 'DNAStringSet_OR_NULL'
codon_coord(
codon = NULL,
filepath = NULL,
cube = c("ACGT", "AGCT", "TCGA", "TGCA", "CATG", "GTAC", "CTAG", "GATC", "ACTG",
"ATCG", "GTCA", "GCTA", "CAGT", "TAGC", "TGAC", "CGAT", "AGTC", "ATGC", "CGTA",
"CTGA", "GACT", "GCAT", "TACG", "TCAG"),
group = c("Z4", "Z5", "Z64", "Z125", "Z4^3", "Z5^3"),
start = NA,
end = NA,
chr = 1L,
strand = "+"
)
## S4 method for signature 'matrix_OR_data_frame'
codon_coord(
codon,
cube = c("ACGT", "AGCT", "TCGA", "TGCA", "CATG", "GTAC", "CTAG", "GATC", "ACTG",
"ATCG", "GTCA", "GCTA", "CAGT", "TAGC", "TGAC", "CGAT", "AGTC", "ATGC", "CGTA",
"CTGA", "GACT", "GCAT", "TACG", "TCAG"),
group = c("Z64", "Z125", "Z4^3", "Z5^3")
)
codon |
An object from |
... |
Not in use. |
group |
A character string denoting the group representation for the given base or codon as shown in reference (2-3). |
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 |
A character string denoting one of the 24 Genetic-code cubes, as given in references (2-3). |
start , end , chr , strand |
Optional parameters required to build a
|
Symbols "-" and "N" usually found in DNA sequence alignments to denote gaps and missing/unknown bases are represented by the number: '-1' on Z4 and '0' on Z5. In Z64 the symbol 'NA' will be returned for codons including symbols "-" and "N".
This function returns a GRanges-class
object
carrying the codon sequence(s) and their respective coordinates in the
requested Abelian group or simply, when group = 'Z5^3'
3D-coordinates, which are derive from Z5 as indicated in reference (3).
Notice that the coordinates can be 3D or just one-dimension ("Z64" or
"Z125"). Hence, the pairwise alignment provided in argument
codon must correspond to codon sequences.
A CodonGroup-class
object.
Robersy Sanchez https://genomaths.com
Robersy Sanchez, Jesus Barreto (2021) Genomic Abelian Finite Groups. doi: 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.
Symmetric Group of the Genetic-Code Cubes.
codon_matrix, base_coord and base2int.
## Load a pairwise alignment
data("aln", package = "GenomAutomorphism")
aln
## DNA base representation in the Abelian group Z5
bs_cor <- codon_coord(
codon = aln,
cube = "ACGT",
group = "Z5"
)
bs_cor ## 3-D coordinates
## DNA base representation in the Abelian group Z64
bs_cor <- codon_coord(
codon = aln,
cube = "ACGT",
group = "Z64"
)
bs_cor
## Giving a matrix of codons
codon_coord(base2codon(x = aln))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.