-
Notifications
You must be signed in to change notification settings - Fork 1
/
coerce.R
199 lines (173 loc) · 6.16 KB
/
coerce.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
187
188
189
190
191
192
193
194
195
196
197
198
#' @title Coercing functions
#' @name as
#' @rdname as
#'
#' @description
#' Coercing functions available for HiCExperiment objects.
#'
#' @aliases coerce,HiCExperiment,GInteractions-method
#' @aliases coerce,HiCExperiment,InteractionSet-method
#' @aliases coerce,HiCExperiment,ContactMatrix-method
#' @aliases coerce,HiCExperiment,matrix-method
#' @aliases coerce,HiCExperiment,data.frame-method
#' @aliases as.matrix,HiCExperiment-method
#' @aliases as.data.frame,HiCExperiment-method
#'
#' @import InteractionSet
#' @importFrom GenomicRanges mcols
#' @param gi GInteractions object
#' @param x HiCExperiment object
#' @param use.scores Which scores to use to inflate GInteractions
#' @param cm A `ContactMatrix` object
#' @param replace_NA Replace NA values
#' @param sparse Whether to return the contact matrix as a sparse matrix
#' @param seqnames1,start1,end1,seqnames2,start2,end2 Names (as strings) of
#' columns containing corresponding information in a data.frame parsed into
#' GInteractions
#' (default: FALSE)
#' @param df A `data.frame` object
#' @importFrom Matrix as.matrix
#' @export
#' @examples
#' mcoolPath <- HiContactsData::HiContactsData('yeast_wt', 'mcool')
#' contacts <- import(mcoolPath, focus = 'XVI', resolution = 16000, format = 'cool')
#' gis <- interactions(contacts)
#' cm <- gi2cm(gis, 'balanced')
#' cm
#' cm2matrix(cm)[1:10, 1:10]
#' df2gi(data.frame(
#' chr1 = 'I', start1 = 10, end1 = 100,
#' chr2 = 'I', start2 = 40, end2 = 1000,
#' score = 12,
#' weight = 0.234,
#' filtered = TRUE
#' ), seqnames1 = 'chr1', seqnames2 = 'chr2')
NULL
################################################################################
################################################################################
############### ###############
############### HICEXPERIMENT COERCING ###############
############### ###############
################################################################################
################################################################################
#' @export
#' @name as
setAs("HiCExperiment", "GInteractions", function(from) interactions(from))
#' @export
#' @name as
setAs("HiCExperiment", "InteractionSet", function(from) {
gi <- interactions(from)
is <- InteractionSet::InteractionSet(
gi,
assays = lapply(
names(scores(from)), function(n) as.matrix(scores(from, n))
)
)
return(is)
})
#' @export
#' @name as
setAs("HiCExperiment", "ContactMatrix", function(from) {
x <- interactions(from, fillout.regions = TRUE)
if ('balanced' %in% names(scores(from))) {
x$score <- scores(from, 'balanced')
gi2cm(x)
}
else {
x$score <- scores(from, 1)
gi2cm(x)
}
})
#' @export
#' @name as
setAs("HiCExperiment", "matrix", function(from) {
as(from, "ContactMatrix") |> cm2matrix(sparse = FALSE)
})
#' @export
#' @name as
setAs("HiCExperiment", "dgTMatrix", function(from) {
as(from, "ContactMatrix") |> cm2matrix(sparse = TRUE)
})
#' @export
#' @name as
setAs("HiCExperiment", "data.frame", function(from) {
x <- interactions(from)
x <- as.data.frame(x)
x <- x[, !colnames(x) %in% c("chr1", "chr2", "bin_id1.1", "bin_id2.1")]
for (n in names(scores(from))) {
x[[n]] <- scores(from, n)
}
return(x)
})
#' @export
#' @name as
setMethod("as.matrix", "HiCExperiment", function(x, use.scores = "balanced", sparse = FALSE) {
interactions(x) |> gi2cm(use.scores = use.scores) |> cm2matrix(sparse = sparse)
})
#' @export
#' @name as
setMethod("as.data.frame", "HiCExperiment", function(x) {
as(x, 'data.frame')
})
################################################################################
################################################################################
############### ###############
############### OTHER COERCING ###############
############### ###############
################################################################################
################################################################################
#' @export
#' @name as
gi2cm <- function(gi, use.scores = 'score') {
if (!use.scores %in% colnames(S4Vectors::mcols(gi)))
stop("`use.scores` argument not found in the provided interactions")
# seqnames <- unique(GenomeInfoDb::seqnames(regions(gi)))
# if (length(seqnames == 1)) {
# gr <- GenomicRanges::GRanges(seqinfo(gi))
# res <- GenomicRanges::width(regions(gi))[1]
# gr <- GenomicRanges::tileGenome(
# seqinfo(gi),
# tilewidth = res,
# cut.last.tile.in.chrom = TRUE
# )
# gr <- gr[GenomicRanges::seqnames(gr) == seqnames]
# InteractionSet::replaceRegions(gi) <- gr
# }
an <- anchors(gi)
allTrans <- all(GenomicRanges::seqnames(an[[1]]) != GenomicRanges::seqnames(an[[2]]))
if (allTrans) {
stop("Only trans contacts are found.")
}
InteractionSet::inflate(
gi,
rows = InteractionSet::regions(gi),
columns = InteractionSet::regions(gi),
fill = GenomicRanges::mcols(gi)[[use.scores]],
sparse = TRUE
)
}
#' @name as
#' @export
cm2matrix <- function(cm, replace_NA = NA, sparse = FALSE) {
m <- Matrix::as.matrix(cm)
m <- as(m, 'TsparseMatrix')
if (!is.na(replace_NA)) m[is.na(m)] <- replace_NA
if (!sparse) m <- base::as.matrix(m)
m
}
#' @name as
#' @export
df2gi <- function(df, seqnames1 = 'seqnames1', start1 = 'start1', end1 = 'end1', seqnames2 = 'seqnames2', start2 = 'start2', end2 = 'end2') {
gi <- InteractionSet::GInteractions(
anchor1 = GenomicRanges::GRanges(
df[[seqnames1]], IRanges::IRanges(df[[start1]], df[[end1]])
),
anchor2 = GenomicRanges::GRanges(
df[[seqnames2]], IRanges::IRanges(df[[start2]], df[[end2]])
)
)
S4Vectors::mcols(gi) <- dplyr::select(
df, !dplyr::any_of(c(seqnames1, start1, end1, seqnames2, start2, end2))
)
return(gi)
}