forked from BIMSBbioinfo/genomation
-
Notifications
You must be signed in to change notification settings - Fork 0
/
scoreMatrixList.R
367 lines (327 loc) · 13.2 KB
/
scoreMatrixList.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
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
#######################
# S3 #
#######################
.plotXYlab = function(lab, pos){
if(!pos %in% c('x','y'))
stop('pos can only be x or y')
if(pos == 'x'){
.plotLab(lab, srt=0)
}
if(pos == 'y'){
.plotLab(lab, srt=90)
}
}
.plotLab = function(lab, srt, xlim, ylim, x, y){
plot.new()
plot.window(xlim=c(1,5), ylim=c(1,5))
text(3, 3, labels=lab, cex=2, srt=srt)
}
# ---------------------------------------------------------------------------- #
#' Make ScoreMatrixList from multiple targets
#'
#' The function constructs a list of \code{ScoreMatrix} objects in the form
#' of \code{ScoreMatrixList} object. This object can be visualized using
#' \code{multiHeatMatrix}, \code{heatMeta} or \code{plotMeta}
#'
#' @param targets can be a list of \code{scoreMatrix} objects, that are coerced
#' to the \code{ScoreMatrixList}, a list of \code{RleList} objects, or a
#' character vector specifying the locations of mulitple bam files that
#' are used to construct the \code{scoreMatrixList}. If l is either a
#' RleList object or a character vector of files, it is obligatory to
#' give a granges argument.
#' @param windows \code{GenomicRanges} containing viewpoints for the scoreMatrix
#' or ScoreMatrixList functions
#' @param bin.num an integer telling the number of bins to bin the score matrix
#' @param bin.op an name of the function that will be used for smoothing windows of ranges
#' @param strand.aware a boolean telling the function whether to reverse the
#' coverage of ranges that come from - strand (e.g. when plotting
#' enrichment around transcription start sites)
#' @param weight.col if the object is \code{GRanges} object a numeric column
#' in meta data part can be used as weights. This is particularly
#' useful when genomic regions have scores other than their
#' coverage values, such as percent methylation, conservation
#' scores, GC content, etc.
#' @param is.noCovNA (Default:FALSE)
#' if TRUE,and if 'targets' is a GRanges object with 'weight.col'
#' provided, the bases that are uncovered will be preserved as
#' NA in the returned object. This useful for situations where
#' you can not have coverage all over the genome, such as CpG
#' methylation values.
#'
#' @param type if \code{targets} is a character vector of file paths, then type
#' designates the type of the corresponding files (bam or bigWig)
#' @param rpm boolean telling whether to normalize the coverage to per milion reads.
#' FALSE by default.
#' @param unique boolean which tells the function to remove duplicated reads
#' based on chr, start, end and strand
#' @param extend numeric which tells the function to extend the features
#' ( i.e aligned reads) to total
#' length ofwidth+extend
#' @param param ScanBamParam object
#'
#' @return returns a \code{ScoreMatrixList} object
#'
#' @examples
#'
#' # visualize the distribution of cage clusters and cpg islands around promoters
#' data(cage)
#' data(cpgi)
#' data(promoters)
#'
#' cage$tpm = NULL
#' targets = GRangesList(cage=cage, cpgi=cpgi)
#' sml = ScoreMatrixList(targets, promoters, bin.num=10, strand.aware=TRUE)
#' multiHeatMatrix(sml)
#'
#' @export
#' @docType methods
#' @rdname ScoreMatrixList-methods
ScoreMatrixList = function(targets, windows=NULL, bin.num=NULL,
bin.op='mean', strand.aware=FALSE, weight.col=NULL,
is.noCovNA=FALSE, type='', rpm=FALSE, unique=FALSE,
extend=0, param=NULL){
len = length(targets)
if(len == 0L)
stop('targets argument is empty')
# this checks whether we can work with the corresponding targets object class set
list.ind = grepl('list', class(targets)) | grepl('List', class(targets))
if(len > 1L & !list.ind){
if(all(is.character(targets)) && is.null(type) && all(file.exists(targets)))
stop('targets argument is neither a list like object (e.g. GRangesList),
nor a set of files')
}
# ----------------------------------------------------------------- #
# checks whether the list argument contains only scoreMatrix objects
if(all(unlist(lapply(targets, class)) == 'ScoreMatrix'))
return(new("ScoreMatrixList",targets))
# ----------------------------------------------------------------- #
if(is.null(windows))
stop("windows of class GRanges must be defined")
# Given a list of RleList objects and a granges object, returns the scoreMatrix list Object
if(list.ind & !all(unlist(lapply(targets, class)) %in% c('SimpleRleList', 'RleList','GRanges'))){
stop('targets should be one of the following: an RleList, list of Rle,
GRangesList, a list of GRanges objects')
}
if(!list.ind && all(file.exists(targets)) && is.null(type))
stop('When providing a file, it is necessary to give the type of the file')
# gets the names for the resulting list
if(all(is.character(targets))){
names = basename(targets)
}else{
names = names(targets)
}
sml = list()
for(i in 1:length(targets)){
message(paste('working on:',names[i]))
if(is.null(bin.num) && all(width(windows) == unique(width(windows)))){
sml[[i]] = ScoreMatrix(targets[[i]], windows=windows,
strand.aware=strand.aware,
weight.col=weight.col,
is.noCovNA=is.noCovNA,
type=type,
rpm=rpm,
unique=unique,
extend=extend,
param=param)
} else{
if(is.null(bin.num))
bin.num = 10
sml[[i]] = ScoreMatrixBin(targets[[i]], windows=windows,
bin.num=bin.num,
bin.op=bin.op,
strand.aware=strand.aware,
weight.col=weight.col,
is.noCovNA = is.noCovNA,
type=type,
rpm=rpm,
unique=unique,
extend=extend,
param=param)
}
}
names(sml) = names
return(new("ScoreMatrixList",sml))
}
# ---------------------------------------------------------------------------- #
# Validator
.valid.ScoreMatrixList = function(l){
errors = character()
# checks whether all matrices are of class scoreMatrix
if(!all(unlist(lapply(l, function(x)class(x) == 'scoreMatrix'))))
errors = paste(errors,
'All elements for ScoreMatrixList
need to be of class scoreMatrix',
sep='\n')
# checks whether all matrices are numeric
if(!all(unlist(lapply(l, function(x)all(is.integer(x) | is.numeric(x))))))
errors = paste(errors, '
Not all matrices are of type integer or numeric',
sep='\n')
if(length(errors) == 0) TRUE else errors
}
# ---------------------------------------------------------------------------- #
# show Methods
#' @rdname show-methods
setMethod("show", "ScoreMatrixList",
function(object){
dims = lapply(object, dim)
len = length(object)
widths = apply(do.call(rbind, dims),2, function(x)max(nchar(x)))
message('scoreMatrixlist of length:', len, '\n')
for(i in 1:len){
s=sprintf(paste('%d%s ','%',widths[1],'d %',widths[2],'d', sep=''),
i, '. scoreMatrix with dims:', dims[[i]][1], dims[[i]][2])
message(s)
}
}
)
# ---------------------------------------------------------------------------- #
#' Scale the ScoreMatrixList
#'
#' Scales each ScoreMatrix in the ScoreMatrixList object, by rows and/or columns
#'
#' @param sml a \code{ScoreMatrixList} object
#' @param columns a \code{columns} whether to scale the matrix by columns.
#' Set by default to FALSE
#' @param rows a \code{rows} Whether to scale the matrix by rows. Set by default
#' to TRUE
#' @param scalefun a function object that takes as input a matrix and returns a
#' matrix.
#' @param ... other argments that be passed to the function
#' By default the argument is set to the R scale function with center=TRUE and
#' scale=TRUE
#'
#' @usage scaleScoreMatrixList(sml, columns, rows, scalefun, ...)
#' @return \code{ScoreMatrixList} object
#' @examples
#'
#' data(cage)
#' data(cpgi)
#' data(promoters)
#'
#' cage$tpm = NULL
#' targets = GRangesList(cage=cage, cpgi=cpgi)
#' sml = ScoreMatrixList(targets, promoters, bin.num=10, strand.aware=TRUE)
#' sml.scaled = scaleScoreMatrixList(sml, rows=TRUE)
#' multiHeatMatrix(sml)
#'
#'
#' @docType methods
#' @rdname scaleScoreMatrixList
#' @export
setGeneric("scaleScoreMatrixList",
function(sml,
columns=FALSE, rows=TRUE,
scalefun=NULL,
...)
standardGeneric("scaleScoreMatrixList") )
#' @aliases scaleScoreMatrixList,ScoreMatrixList-method
#' @rdname scaleScoreMatrixList
setMethod("scaleScoreMatrixList", signature("ScoreMatrixList"),
function(sml, columns, rows, scalefun, ...){
sml = lapply(sml, function(x)
scaleScoreMatrix(x,
columns=columns,
rows=rows,
scalefun=scalefun))
sml = as(sml,'ScoreMatrixList')
return (sml)
}
)
# ---------------------------------------------------------------------------- #
#' Get common rows from all matrices in a ScoreMatrixList object
#'
#' Returns a intersection of rows for each matrix in a ScoreMatrixList object.
#' This is done using the rownames of each element in the list.
#'
#' @param sml a \code{ScoreMatrixList} object
#' @param reorder if TRUE \code{ScoreMatrix} objects in the list are sorted
#' based on their common row ids.
#'
#' @return \code{ScoreMatrixList} object
#' @examples
#' target = GRanges(rep(c(1,2),each=7),
#' IRanges(rep(c(1,1,2,3,7,8,9), times=2), width=5),
#' weight = rep(c(1,2),each=7))
#'
#' windows1 = GRanges(rep(c(1,2),each=2),
#' IRanges(rep(c(1,2), times=2), width=5),
#' strand=c('-','+','-','+'))
#' windows2 = windows1[c(1,3)]
#' sml = as(list(ScoreMatrix(target, windows1),
#' ScoreMatrix(target, windows2)), 'ScoreMatrixList')
#' sml
#' intersectScoreMatrixList(sml)
#'
#' @docType methods
#' @rdname intersectScoreMatrixList-methods
#' @export
setGeneric("intersectScoreMatrixList",
function(sml,reorder=FALSE)
standardGeneric("intersectScoreMatrixList") )
#' @aliases intersectScoreMatrixList,ScoreMatrixList-method
#' @rdname intersectScoreMatrixList-methods
setMethod("intersectScoreMatrixList", signature("ScoreMatrixList"),
function(sml,reorder){
rnames = Reduce('intersect' ,lapply(sml, rownames))
if(reorder){
sml = as(lapply(sml, function(x){
x=x[rownames(x) %in% rnames,]
x[order(rownames(x)),]
}),
'ScoreMatrixList')
}else{
sml = as(lapply(sml, function(x)x[rownames(x) %in% rnames,]),
'ScoreMatrixList')
}
return (sml)
}
)
# ---------------------------------------------------------------------------- #
#' Reorder all elements of a ScoreMatrixList to a given ordering vector
#' @param sml \code{ScoreMatrixList} object
#' @param ord.vec an integer vector
#' @return \code{ScoreMatrixList} object
#'
#' @examples
#' data(cage)
#' data(cpgi)
#' data(promoters)
#'
#' cage$tpm = NULL
#' targets = GRangesList(cage=cage, cpgi=cpgi)
#' sml = ScoreMatrixList(targets, promoters, bin.num=10)
#' kmeans.clust = kmeans(sml$cage,3)
#'
#' sml.ordered = orderBy(sml, kmeans.clust$cluster)
#' multiHeatMatrix(sml.ordered)
#'
#' @docType methods
#' @rdname orderBy-methods
#' @export
setGeneric("orderBy",
function(sml,ord.vec)
standardGeneric("orderBy") )
#' @aliases orderBy,ScoreMatrixList-method
#' @rdname orderBy-methods
setMethod("orderBy", signature("ScoreMatrixList"),
function(sml, ord.vec){
if(is.null(ord.vec))
return(sml)
if(!is.integer(ord.vec))
stop('ord.vec needs to be of an integer vector')
sml = lapply(sml, function(x)x[ord.vec,])
return (as(sml,'ScoreMatrixList'))
}
)
# ---------------------------------------------------------------------------- #
#' @aliases binMatrix,ScoreMatrixList-method
#' @rdname binMatrix-methods
setMethod("binMatrix", signature("ScoreMatrixList"),
function(x, bin.num=NULL, fun='mean'){
if(is.null(bin.num))
return(x)
return(new("ScoreMatrix",
lapply(x, function(y)binMatrix(y, bin.num=bin.num, fun=fun))))
}
)