File: GAlignmentsList-class.Rd

package info (click to toggle)
r-bioc-genomicalignments 1.0.6-1
  • links: PTS, VCS
  • area: main
  • in suites: jessie, jessie-kfreebsd
  • size: 2,980 kB
  • ctags: 54
  • sloc: ansic: 1,493; makefile: 4; sh: 3
file content (402 lines) | stat: -rw-r--r-- 13,471 bytes parent folder | download
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
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
\name{GAlignmentsList-class}
\docType{class}

% Class
\alias{class:GAlignmentsList}
\alias{GAlignmentsList-class}
\alias{GAlignmentsList}

% Constructors:
\alias{GAlignmentsList}
\alias{updateObject,GAlignmentsList-method}

% Accessors:
\alias{names,GAlignmentsList-method}
\alias{names<-,GAlignmentsList-method}
\alias{seqnames,GAlignmentsList-method}
\alias{seqnames<-,GAlignmentsList-method}
\alias{rname,GAlignmentsList-method}
\alias{rname<-,GAlignmentsList-method}
\alias{strand,GAlignmentsList-method}
\alias{strand<-,GAlignmentsList-method}
\alias{cigar,GAlignmentsList-method}
\alias{qwidth,GAlignmentsList-method}
\alias{njunc,GAlignmentsList-method}
\alias{elementMetadata,GAlignmentsList-method}
\alias{elementMetadata<-,GAlignmentsList-method}
\alias{seqinfo,GAlignmentsList-method}
\alias{seqinfo<-,GAlignmentsList-method}
\alias{start,GAlignmentsList-method}
\alias{end,GAlignmentsList-method}
\alias{width,GAlignmentsList-method}

% Coercion:
\alias{as.data.frame,GAlignmentsList-method}
\alias{coerce,GAlignmentsList,GRangesList-method}
\alias{coerce,GAlignmentsList,GRanges-method}
\alias{coerce,GAlignmentsList,RangesList-method}
\alias{coerce,GAlignmentsList,Ranges-method}
\alias{coerce,GAlignmentPairs,GAlignmentsList-method}
\alias{grglist,GAlignmentsList-method}
\alias{granges,GAlignmentsList-method}
\alias{rglist,GAlignmentsList-method}
\alias{ranges,GAlignmentsList-method}

% Combining:
\alias{c,GAlignmentsList-method}

% extractList() and family:
\alias{relistToClass,GAlignments-method}

% show:
\alias{show,GAlignmentsList-method}

% Old stuff (deprecated & defunct):
\alias{makeGAlignmentsListFromFeatureFragments}
\alias{ngap,GAlignmentsList-method}

\title{GAlignmentsList objects}

\description{
  The GAlignmentsList class is a container for storing a collection of
  \link{GAlignments} objects.
}

\details{
  A GAlignmentsList object contains a list of \link{GAlignments} objects.
  The majority of operations on this page are described in more detail
  on the GAlignments man page, see ?\code{GAlignments}.
} 

\section{Constructor}{
  \describe{
    \item{}{
      \code{GAlignmentsList(...)}:
      Creates a GAlignmentsList from a list of \link{GAlignments} objects.
    }
  }
}

\section{Accessors}{
  In the code snippets below, \code{x} is a GAlignmentsList object.

  \describe{
    \item{}{
      \code{length(x)}:
      Return the number of elements in \code{x}.
    }
    \item{}{
      \code{names(x)}, \code{names(x) <- value}:
      Get or set the names of the elements of \code{x}.
    }
    \item{}{
      \code{seqnames(x)}, \code{seqnames(x) <- value}:
      Get or set the name of the reference sequences of the
      alignments in each element of \code{x}.
    }
    \item{}{
      \code{rname(x)}, \code{rname(x) <- value}:
      Same as \code{seqnames(x)} and \code{seqnames(x) <- value}.
    }
    \item{}{
      \code{strand(x)}, \code{strand(x) <- value}:
      Get or set the strand of the alignments in each element 
      of \code{x}.
    }
    \item{}{
      \code{cigar(x)}:
      Returns a character list of length \code{length(x)}
      containing the CIGAR string for the alignments in
      each element of \code{x}.
    }
    \item{}{
      \code{qwidth(x)}:
      Returns an integer list of length \code{length(x)}
      containing the length of the alignments in each element of
      \code{x} *after* hard clipping (i.e. the length of the 
      query sequence that is stored in the corresponding SAM/BAM record).
    }
    \item{}{
      \code{start(x)}, \code{end(x)}:
      Returns an integer list of length \code{length(x)}
      containing the "start" and "end" (respectively) of the 
      alignments in each element of \code{x}. 
    }
    \item{}{
      \code{width(x)}:
      Returns an integer list of length \code{length(x)} containing
      the "width" of the alignments in each element of \code{x}.
    }
    \item{}{
      \code{njunc(x)}:
      Returns an integer list of length \code{x} containing the number
      of junctions (i.e. N operations in the CIGAR) for the alignments
      in each element of \code{x}.
    }
    \item{}{
      \code{seqinfo(x)}, \code{seqinfo(x) <- value}:
      Get or set the information about the underlying sequences in each
      element of \code{x}. \code{value} must be a list of \link{Seqinfo} 
      objects.
    }
    \item{}{
      \code{seqlevels(x)}, \code{seqlevels(x) <- value}:
      Get or set the sequence levels of the alignments in each element
      of \code{x}.
    }
    \item{}{
      \code{seqlengths(x)}, \code{seqlengths(x) <- value}:
      Get or set the sequence lengths for each element of \code{x}.
      \code{seqlengths(x)} is equivalent to \code{seqlengths(seqinfo(x))}.
      \code{value} can be a named non-negative integer or numeric vector
      eventually with NAs.
    }
    \item{}{
      \code{isCircular(x)}, \code{isCircular(x) <- value}:
      Get or set the circularity flags for the alignments in each
      element in \code{x}. \code{value} must be a named logical list 
      eventually with NAs.
    }
    \item{}{
      \code{genome(x)}, \code{genome(x) <- value}:
      Get or set the genome identifier or assembly name for the alignments 
      in each element of \code{x}. \code{value} must be a named character 
      list eventually with NAs.
    }
    \item{}{
      \code{seqnameStyle(x)}:
      Get or set the seqname style for alignments in each element of \code{x}.
    }
  }
}

\section{Coercion}{
  In the code snippets below, \code{x} is a GAlignmentsList object.

  \describe{
    \item{}{
      \code{granges(x, use.mcols=FALSE, ignore.strand=FALSE)}, 
      \code{ranges(x)}:
      Return either a \link{GRanges} or a \link[IRanges]{IRanges}
      object of length \code{length(x)}. Note this coercion IGNORES 
      the cigar information. The resulting ranges span the entire
      range, including any gaps or spaces between paired-end reads.

      If \code{use.mcols} is TRUE and \code{x} has metadata columns on it
      (accessible with \code{mcols(x)}), they're propagated to the returned
      object.

      \code{granges} coercion supports \code{ignore.strand} to allow 
      ranges of opposite strand to be combined (see examples). All
      ranges in the resulting GRanges will have strand \sQuote{*}.
    }
    \item{}{
      \code{grglist(x, use.mcols=FALSE, ignore.strand=FALSE)}, 
      \code{rglist(x, use.mcols=FALSE)}:
      Return either a \link{GRangesList} or a \link[IRanges]{IRangesList}
      object of length \code{length(x)}. This coercion RESPECTS the cigar 
      information. The resulting ranges are fragments of the original ranges 
      that do not include gaps or spaces between paired-end reads.
 
      If \code{use.mcols} is TRUE and \code{x} has metadata columns on it
      (accessible with \code{mcols(x)}), they're propagated to the returned
      object.

      \code{grglist} coercion supports \code{ignore.strand} to allow 
      ranges of opposite strand to be combined (see examples). All 
      ranges in the resulting GRangesList will have strand \sQuote{*}.
    }
    \item{}{
      \code{as(x, "GRangesList")}, \code{as(x, "GRanges")},
      \code{as(x, "RangesList")}, \code{as(x, "Ranges")}:
      Alternate ways of doing \code{grglist(x, use.mcols=TRUE)},
      \code{granges(x, use.mcols=TRUE)}, \code{rglist(x, use.mcols=TRUE)},
      and \code{ranges(x)}, respectively.
    }
    \item{}{
      \code{as(x, "GALignmentsList")}: Here \code{x} is a
      \link{GAlignmentPairs} object. Return a GAlignmentsList object of length
      \code{length(x)} where the i-th list element represents the ranges of
      the i-th alignment pair in \code{x}.
    }
  }
}

\section{Subsetting and related operations}{
  In the code snippets below, \code{x} is a GAlignmentsList object.

  \describe{
    \item{}{
      \code{x[i]}, \code{x[i] <- value}:
      Get or set list elements \code{i}. \code{i} can be a numeric 
      or logical vector. \code{value} must be a GAlignments.
    }
    \item{}{
      \code{x[[i]]}, \code{x[[i]] <- value}:
      Same as \code{x[i]}, \code{x[i] <- value}. 
    }
    \item{}{
      \code{x[i, j]}, \code{x[i, j] <- value}:
      Get or set list elements \code{i} with optional metadata columns
      \code{j}. \code{i} can be a numeric, logical or missing. 
      \code{value} must be a GAlignments.
    }
  }
}

\section{Combining}{
  \describe{
    \item{}{
      \code{c(...)}:
      Concatenates the GAlignmentsList objects in \code{...}.
    }
  }
}

\references{
  \url{http://samtools.sourceforge.net/}
}

\author{Valerie Obenchain <vobencha@fhcrc.org}

\seealso{
  \itemize{
    \item \code{\link{readGAlignmentsList}} for reading genomic alignments
          from a file (typically a BAM file) into a GAlignmentsList object.

    \item \link{GAlignments} and \link{GAlignmentPairs} objects for handling
          aligned single- and paired-end reads, respectively.

    \item \link{junctions-methods} for extracting and summarizing junctions
          from a GAlignmentsList object.

    \item \link[GenomicAlignments]{findOverlaps-methods} for finding range
          overlaps between a GAlignmentsList object and another range-based
          object.

    \item \code{\link[GenomicRanges]{seqinfo}} in the \pkg{GenomicRanges}
          package for getting/setting/modifying the sequence information
          stored in an object.

    \item The \link[GenomicRanges]{GRanges} and
          \link[GenomicRanges]{GRangesList} classes defined and documented
          in the \pkg{GenomicRanges} package.
  }
}

\examples{
gal1 <- GAlignments(
    seqnames=Rle(factor(c("chr1", "chr2", "chr1", "chr3")),
        c(1, 3, 2, 4)),
    pos=1:10, cigar=paste0(10:1, "M"),
    strand=Rle(strand(c("-", "+", "*", "+", "-")), c(1, 2, 2, 3, 2)),
    names=head(letters, 10), score=1:10)

gal2 <- GAlignments(
    seqnames=Rle(factor(c("chr2", "chr4")), c(3, 4)), pos=1:7,
    cigar=c("5M", "3M2N3M2N3M", "5M", "10M", "5M1N4M", "8M2N1M", "5M"),
    strand=Rle(strand(c("-", "+")), c(4, 3)),
    names=tail(letters, 7), score=1:7)

galist <- GAlignmentsList(noGaps=gal1, Gaps=gal2)


## ---------------------------------------------------------------------
## A. BASIC MANIPULATION
## ---------------------------------------------------------------------

length(galist)
names(galist)
seqnames(galist)
strand(galist)
head(cigar(galist))
head(qwidth(galist))
head(start(galist))
head(end(galist))
head(width(galist))
head(njunc(galist))
seqlevels(galist)

## Rename the reference sequences:
seqlevels(galist) <- sub("chr", "seq", seqlevels(galist))
seqlevels(galist)

grglist(galist)  # a GRangesList object
rglist(galist)   # an IRangesList object

## ---------------------------------------------------------------------
## B. SUBSETTING
## ---------------------------------------------------------------------

galist[strand(galist) == "-"]
has_junctions <- sapply(galist,
                        function(x) any(grepl("N", cigar(x), fixed=TRUE)))
galist[has_junctions]

## Different ways to subset:
galist[2]             # a GAlignments object of length 1
galist[[2]]           # a GAlignments object of length 1
grglist(galist[2])  # a GRangesList object of length 1
rglist(galist[2])   # a NormalIRangesList object of length 1

## ---------------------------------------------------------------------
## C. mcols()/elementMetadata()
## ---------------------------------------------------------------------

## Metadata can be defined on the individual GAlignment elements
## and the overall GAlignmentsList object. By default, 'level=between' 
## extracts the GALignmentsList metadata. Using 'level=within' 
## will extract the metadata on the individual GAlignments objects.

mcols(galist) ## no metadata on the GAlignmentsList object
mcols(galist, level="within")


## ---------------------------------------------------------------------
## D. readGAlignmentsList()
## ---------------------------------------------------------------------

library(pasillaBamSubset)

## 'file' as character.
fl <- untreated3_chr4()
galist1 <- readGAlignmentsList(fl)

galist1[1:3]
length(galist1)
table(elementLengths(galist1))

## When 'file' is a BamFile, 'asMates' must be TRUE. If FALSE,
## the data are treated as single-end and each list element of the
## GAlignmentsList will be of length 1. For single-end data 
## use readGAlignments() instead of readGAlignmentsList().
bf <- BamFile(fl, yieldSize=3, asMates=TRUE)
readGAlignmentsList(bf)

## Use a 'param' to fine tune the results.
param <- ScanBamParam(flag=scanBamFlag(isProperPair=TRUE))
galist2 <- readGAlignmentsList(fl, param=param)
length(galist2)


## ---------------------------------------------------------------------
## E. COERCION 
## ---------------------------------------------------------------------

## The granges() and grlist() coercions support 'ignore.strand' to 
## allow ranges from different strand to be combined. In this example 
## paired-end reads aligned to opposite strands were read into a 
## GAlignmentsList. If the desired operation is to combine these ranges, 
## regardless of gaps or the space between pairs, 'ignore.strand' must be TRUE.
granges(galist[1])
granges(galist[1], ignore.strand=TRUE)

## grglist() splits ranges by gap and the space between list elements.
galist <- GAlignmentsList(noGaps=gal1, Gaps=gal2)
grglist(galist)
grglist(galist, ignore.strand=TRUE)
}

\keyword{methods}
\keyword{classes}