/
conserved_regions.R
186 lines (173 loc) · 5.68 KB
/
conserved_regions.R
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
## Copyright (C) 2021 Robersy Sanchez <https://genomaths.com/>
## Author: Robersy Sanchez This file is part of the R package
## 'GenomAutomorphism'. 'GenomAutomorphism' is a free
## software: you can redistribute it and/or modify it under the
## terms of the GNU General Public License as published by the Free
## Software Foundation, either version 3 of the License, or (at
## your option) any later version. This program is distributed in
## the hope that it will be useful, but WITHOUT ANY WARRANTY;
## without even the implied warranty of MERCHANTABILITY or FITNESS
## FOR A PARTICULAR PURPOSE. See the GNU General Public License for
## more details. You should have received a copy of the GNU
## General Public License along with this program; if not, see
## <http://www.gnu.org/licenses/>.
#' @rdname conserved_regions
#' @aliases conserved_regions
#' @title Conserved and Non-conserved Regions from a MSA
#' @description Returns the Conserved or the Non-conserved Regions from a MSA.
#' @param x A \code{\link{Automorphism-class}}, a
#' \code{\link{AutomorphismList-class}},
#' a \code{\link{AutomorphismByCoef}} or a
#' \code{\link{AutomorphismByCoefList}} class object.
#' @param conserved Logical, Whether to return the \emph{conserved} or the
#' \emph{non-conserved regions}.
#' @param output A character string. Type of output.
#' @param ... Not in use.
#' @return A \code{\link{AutomorphismByCoef}} class object containing the
#' requested regions.
#' @import S4Vectors
#' @export
#' @examples
#' ## Load dataset
#' data("autm", package = "GenomAutomorphism")
#' conserved_regions(autm[1:3])
setGeneric(
"conserved_regions",
function(x,
...) {
standardGeneric("conserved_regions")
}
)
#' @rdname conserved_regions
#' @aliases conserved_regions
#' @import GenomicRanges
#' @export
setMethod("conserved_regions",
signature = "Automorphism",
function(x,
conserved = TRUE,
output = c("all_pairs", "unique_pairs", "unique")) {
output <- match.arg(output)
x <- automorphism_bycoef(x)
x <- conserved_regions(
x,
conserved = conserved,
output = output
)
return(x)
}
)
#' @rdname conserved_regions
#' @aliases conserved_regions
#' @import GenomicRanges
#' @param num.cores,tasks Integers. Argument \emph{num.cores} denotes the
#' number of cores to use, i.e. at most how many child processes will be run
#' simultaneously (see \code{\link[BiocParallel]{bplapply}} function from
#' BiocParallel package). Argument \emph{tasks} denotes the number of tasks per
#' job. value must be a scalar integer >= 0L. In this documentation a job is
#' defined as a single call to a function, such as
#' \code{\link[BiocParallel]{bplapply}}. A task is the division of the \eqn{X}
#' argument into chunks. When tasks == 0 (default), \eqn{X} is divided as
#' evenly as possible over the number of workers (see
#' \code{\link[BiocParallel]{MulticoreParam}} from BiocParallel package).
#' @param verbose logic(1). If TRUE, enable progress bar.
#' @importFrom BiocParallel multicoreWorkers
#' @export
#' @examples
#' ## Load automorphism found COVID datatset
#' data("covid_autm", package = "GenomAutomorphism")
#'
#' ## Conserved regions in the first 100 codons
#' conserv <- conserved_regions(covid_autm[1:100], output = "unique")
#' conserv
setMethod("conserved_regions",
signature = "AutomorphismList",
function(x,
conserved = TRUE,
output = c("all_pairs", "unique_pairs", "unique"),
num.cores = multicoreWorkers(),
tasks = 0L,
verbose = FALSE) {
output <- match.arg(output)
x <- automorphism_bycoef(x,
num.cores = num.cores,
tasks = tasks,
verbose = verbose
)
x <- unlist(x)
x <- conserved_regions(
x,
conserved = conserved,
output = output
)
return(x)
}
)
#' @rdname conserved_regions
#' @aliases conserved_regions
#' @import GenomicRanges
#' @export
setMethod("conserved_regions",
signature = "AutomorphismByCoef",
function(x,
conserved = TRUE,
output = c("all_pairs", "unique_pairs", "unique")) {
autm <- end <- NULL
output <- match.arg(output)
if (inherits(x$autm, "numeric")) {
if (conserved) {
x <- x[x$autm == 1]
} else {
x <- x[x$autm != 1]
}
}
else {
if (conserved) {
x <- x[x$autm == "1,1,1"]
} else {
x <- x[x$autm != "1,1,1"]
}
}
x <- sortByChromAndEnd(x)
x <- switch(output,
all_pairs = x,
unique_pairs = unique(x),
unique = {
x <- unique(x)
x <- data.table(data.frame(x))
x <- x[, list(
seqnames = unique(seqnames),
end = min(end),
strand = unique(strand),
autm = unique(autm)
),
by = c("start", "cube")
]
x <- makeGRangesFromDataFrame(x,
keep.extra.columns = TRUE
)
x <- x[, c("autm", "cube")]
}
)
return(x)
}
)
#' @rdname conserved_regions
#' @aliases conserved_regions
#' @import GenomicRanges
#' @export
setMethod("conserved_regions",
signature = "AutomorphismByCoefList",
function(x,
conserved = TRUE,
output = c("all_pairs", "unique_pairs", "unique")) {
output <- match.arg(output)
x <- unlist(x)
x <- conserved_regions(
x,
conserved = conserved,
output = output
)
return(x)
}
)