This file is indexed.

/usr/lib/R/site-library/GenomicRanges/doc/GenomicRangesHOWTOs.Rnw is in r-bioc-genomicranges 1.16.4-2.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

   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
 403
 404
 405
 406
 407
 408
 409
 410
 411
 412
 413
 414
 415
 416
 417
 418
 419
 420
 421
 422
 423
 424
 425
 426
 427
 428
 429
 430
 431
 432
 433
 434
 435
 436
 437
 438
 439
 440
 441
 442
 443
 444
 445
 446
 447
 448
 449
 450
 451
 452
 453
 454
 455
 456
 457
 458
 459
 460
 461
 462
 463
 464
 465
 466
 467
 468
 469
 470
 471
 472
 473
 474
 475
 476
 477
 478
 479
 480
 481
 482
 483
 484
 485
 486
 487
 488
 489
 490
 491
 492
 493
 494
 495
 496
 497
 498
 499
 500
 501
 502
 503
 504
 505
 506
 507
 508
 509
 510
 511
 512
 513
 514
 515
 516
 517
 518
 519
 520
 521
 522
 523
 524
 525
 526
 527
 528
 529
 530
 531
 532
 533
 534
 535
 536
 537
 538
 539
 540
 541
 542
 543
 544
 545
 546
 547
 548
 549
 550
 551
 552
 553
 554
 555
 556
 557
 558
 559
 560
 561
 562
 563
 564
 565
 566
 567
 568
 569
 570
 571
 572
 573
 574
 575
 576
 577
 578
 579
 580
 581
 582
 583
 584
 585
 586
 587
 588
 589
 590
 591
 592
 593
 594
 595
 596
 597
 598
 599
 600
 601
 602
 603
 604
 605
 606
 607
 608
 609
 610
 611
 612
 613
 614
 615
 616
 617
 618
 619
 620
 621
 622
 623
 624
 625
 626
 627
 628
 629
 630
 631
 632
 633
 634
 635
 636
 637
 638
 639
 640
 641
 642
 643
 644
 645
 646
 647
 648
 649
 650
 651
 652
 653
 654
 655
 656
 657
 658
 659
 660
 661
 662
 663
 664
 665
 666
 667
 668
 669
 670
 671
 672
 673
 674
 675
 676
 677
 678
 679
 680
 681
 682
 683
 684
 685
 686
 687
 688
 689
 690
 691
 692
 693
 694
 695
 696
 697
 698
 699
 700
 701
 702
 703
 704
 705
 706
 707
 708
 709
 710
 711
 712
 713
 714
 715
 716
 717
 718
 719
 720
 721
 722
 723
 724
 725
 726
 727
 728
 729
 730
 731
 732
 733
 734
 735
 736
 737
 738
 739
 740
 741
 742
 743
 744
 745
 746
 747
 748
 749
 750
 751
 752
 753
 754
 755
 756
 757
 758
 759
 760
 761
 762
 763
 764
 765
 766
 767
 768
 769
 770
 771
 772
 773
 774
 775
 776
 777
 778
 779
 780
 781
 782
 783
 784
 785
 786
 787
 788
 789
 790
 791
 792
 793
 794
 795
 796
 797
 798
 799
 800
 801
 802
 803
 804
 805
 806
 807
 808
 809
 810
 811
 812
 813
 814
 815
 816
 817
 818
 819
 820
 821
 822
 823
 824
 825
 826
 827
 828
 829
 830
 831
 832
 833
 834
 835
 836
 837
 838
 839
 840
 841
 842
 843
 844
 845
 846
 847
 848
 849
 850
 851
 852
 853
 854
 855
 856
 857
 858
 859
 860
 861
 862
 863
 864
 865
 866
 867
 868
 869
 870
 871
 872
 873
 874
 875
 876
 877
 878
 879
 880
 881
 882
 883
 884
 885
 886
 887
 888
 889
 890
 891
 892
 893
 894
 895
 896
 897
 898
 899
 900
 901
 902
 903
 904
 905
 906
 907
 908
 909
 910
 911
 912
 913
 914
 915
 916
 917
 918
 919
 920
 921
 922
 923
 924
 925
 926
 927
 928
 929
 930
 931
 932
 933
 934
 935
 936
 937
 938
 939
 940
 941
 942
 943
 944
 945
 946
 947
 948
 949
 950
 951
 952
 953
 954
 955
 956
 957
 958
 959
 960
 961
 962
 963
 964
 965
 966
 967
 968
 969
 970
 971
 972
 973
 974
 975
 976
 977
 978
 979
 980
 981
 982
 983
 984
 985
 986
 987
 988
 989
 990
 991
 992
 993
 994
 995
 996
 997
 998
 999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
%\VignetteIndexEntry{GenomicRanges HOWTOs}
%\VignetteDepends{GenomicRanges, Rsamtools, GenomicAlignments, pasillaBamSubset, TxDb.Dmelanogaster.UCSC.dm3.ensGene, TxDb.Athaliana.BioMart.plantsmart21, AnnotationHub, DESeq, edgeR, TxDb.Hsapiens.UCSC.hg19.knownGene, GenomicFeatures, Biostrings, BSgenome.Hsapiens.UCSC.hg19, KEGG.db, KEGGgraph, BSgenome.Scerevisiae.UCSC.sacCer2}
%\VignetteKeywords{sequence, sequencing, alignments}
%\VignettePackage{GenomicRanges}

\documentclass{article}

\usepackage[authoryear,round]{natbib}

<<style, echo=FALSE, results=tex>>=
BiocStyle::latex(use.unsrturl=FALSE)
@

\title{\Biocpkg{GenomicRanges} HOWTOs}
\author{Bioconductor Team}
\date{Edited: July 2014; Compiled: \today}

\begin{document}

\maketitle

\tableofcontents

<<options, echo=FALSE>>=
options(width=72)
options(showHeadLines=3)
options(showTailLines=3)
.precomputed_results_path <- "precomputed_results"
@



%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\section{Introduction}

\subsection{Purpose of this document}

This document is a collection of {\it HOWTOs}. Each {\it HOWTO} is
a short section that demonstrates how to use the containers and
operations implemented in the \Biocpkg{GenomicRanges} and related
packages (\Biocpkg{IRanges}, \Biocpkg{Biostrings}, \Biocpkg{Rsamtools},
\Biocpkg{GenomicAlignments}, \Biocpkg{BSgenome}, and
\Biocpkg{GenomicFeatures}) to perform a task typically found in the context
of a high throughput sequence analysis.

Unless stated otherwise, the {\it HOWTOs} are self contained, independent
of each other, and can be studied and reproduced in any order.

\subsection{Prerequisites and additional recommended reading}

We assume the reader has some previous experience with \R{} and
with basic manipulation of \Rcode{GRanges}, \Rcode{GRangesList}, \Rcode{Rle},
\Rcode{RleList}, and \Rcode{DataFrame} objects. See the ``An Introduction
to Genomic Ranges Classes'' vignette located in the \Biocpkg{GenomicRanges}
package (in the same folder as this document) for an introduction to these
containers.

Additional recommended readings after this document are the ``Software for
Computing and Annotating Genomic Ranges'' paper[\citet{Lawrence2013ranges}]
and the ``Counting reads with \Rfunction{summarizeOverlaps}'' vignette
located in the \Biocpkg{GenomicAlignments} package.

To display the list of vignettes available in the \Biocpkg{GenomicRanges}
package, use \Rcode{browseVignettes("GenomicRanges")}.

\subsection{Input data and terminology used across the HOWTOs}

In order to avoid repetition, input data, concepts and terms used in more
than one {\it HOWTO} are described here:

\begin{itemize}
  \item {\bf The \Biocpkg{pasillaBamSubset} data package}: contains both a BAM
        file with single-end reads (untreated1\_chr4) and a BAM file with
        paired-end reads (untreated3\_chr4). Each file is a subset of chr4
        from the "Pasilla" experiment.

<<pasillaBamSubset>>=
library(pasillaBamSubset)
untreated1_chr4()
untreated3_chr4()
@

        See \Rcode{?pasillaBamSubset} for more information.

<<pasillaBamSubset_help, eval=FALSE>>=
?pasillaBamSubset
@

  \item {\bf Gene models and \Rclass{TranscriptDb} objects}: A \textit{gene model}
        is essentially a set of annotations that describes the genomic
        locations of the known genes, transcripts, exons, and CDS, for a
        given organism. In \Bioconductor{} it is typically represented as
        a \Rclass{TranscriptDb} object but also sometimes as a \Rclass{GRanges}
        or \Rclass{GRangesList} object.
        The \Biocpkg{GenomicFeatures} package contains tools for making and
        manipulating \Rclass{TranscriptDb} objects.
\end{itemize}



%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\section{HOWTOs}


%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\subsection{How to read single-end reads from a BAM file}

As sample data we use the \Biocpkg{pasillaBamSubset} data package
described in the introduction.

<<readGAlignments_1, results=hide>>=
library(pasillaBamSubset)
un1 <- untreated1_chr4()  # single-end reads
@

Several functions are available for reading BAM files into \R{}:

\begin{verbatim}
  readGAlignments()
  readGAlignmentPairs()
  readGAlignmentsList()
  scanBam()
\end{verbatim}

\Rfunction{scanBam} is a low-level function that returns a list of lists
and is not discussed further here. See \Rcode{?scanBam} in the
\Biocpkg{Rsamtools} package for more information.

Single-end reads can be loaded with the \Rfunction{readGAlignments} function
from the \Biocpkg{GenomicAlignments} package.

<<readGAlignments_2>>=
library(GenomicAlignments)
gal <- readGAlignments(un1)
@

Data subsets can be specified by genomic position, field names, or flag
criteria in the \Rcode{ScanBamParam}. Here we input records that overlap 
position 1 to 5000 on the negative strand with \Rcode{flag} and 
\Rcode{cigar} as metadata columns.
 
<<readGAlignments_3>>=
what <- c("flag", "cigar") 
which <- GRanges("chr4", IRanges(1, 5000)) 
flag <- scanBamFlag(isMinusStrand = TRUE)
param <- ScanBamParam(which=which, what=what, flag=flag)
neg <- readGAlignments(un1, param=param)
neg
@

Another approach to subsetting the data is to use \Rfunction{filterBam}.
This function creates a new BAM file of records passing user-defined 
criteria. See \Rcode{?filterBam} in the \Biocpkg{Rsamtools} package for
more information.


%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\subsection{How to read paired-end reads from a BAM file}

As sample data we use the \Biocpkg{pasillaBamSubset} data package
described in the introduction.

<<readGAlignmentPairs_1>>=
library(pasillaBamSubset)
un3 <- untreated3_chr4()  # paired-end reads
@

Paired-end reads can be loaded with the \Rfunction{readGAlignmentPairs}
or \Rfunction{readGAlignmentsList} function from the
\Biocpkg{GenomicAlignments} package. These functions use the same
mate paring algorithm but output different objects.

Let's start with \Rfunction{readGAlignmentPairs}:

<<readGAlignmentPairs_2>>=
un3 <- untreated3_chr4()
gapairs <- readGAlignmentPairs(un3)
@

The \Robject{GAlignmentPairs} class holds only pairs; reads with no
mate or with ambiguous pairing are discarded.
Each list element holds exactly 2 records (a mated pair). Records
can be accessed as the \Rcode{first} and\Rcode{last} segments in
a template or as \Rcode{left} and \Rcode{right} alignments.
See \Rcode{?GAlignmentPairs} in the \Biocpkg{GenomicAlignments} package
for more information.

<<readGAlignmentPairs_3>>=
gapairs
@ 

For \Rcode{readGAlignmentsList}, mate pairing is performed when \Rcode{asMates}
is set to \Rcode{TRUE} on the \Rcode{BamFile} object, otherwise records are
treated as single-end. 

<<readGAlignmentsList_1>>=
galist <- readGAlignmentsList(BamFile(un3, asMates=TRUE))
@

\Robject{GAlignmentsList} is a more general `list-like' structure
that holds mate pairs as well as non-mates (i.e., singletons, records 
with unmapped mates etc.) A \Rcode{mates} metadata column (accessed
with \Rfunction{mcols}) indicates which records were paired and is set 
on both the individual \Robject{GAlignments} and the outer list elements.

<<readGAlignmentsList_2>>=
galist
@

Non-mated reads are returned as groups by QNAME and contain any number 
of records. Here the non-mate groups range in size from 1 to 9.

<<readGAlignmentsList_3>>=
non_mates <- galist[unlist(mcols(galist)$mates) == FALSE]
table(elementLengths(non_mates))
@


%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\subsection{How to read and process a big BAM file by chunks in order to
            reduce memory usage}

A large BAM file can be iterated through in chunks by setting a
\Rcode{yieldSize} on the \Rclass{BamFile} object.
As sample data we use the \Biocpkg{pasillaBamSubset} data package
described in the introduction.

<<yieldSize>>=
library(pasillaBamSubset)
un1 <- untreated1_chr4()
bf <- BamFile(un1, yieldSize=100000)
@

Iteration through a BAM file requires that the file be opened, repeatedly
queried inside a loop, then closed. Repeated calls to 
\Rfunction{readGAlignments} without opening the file first result
in the same 100000 records returned each time.

<<readGAlignments_by_chunk>>=
open(bf)
cvg <- NULL
repeat {
    chunk <- readGAlignments(bf)
    if (length(chunk) == 0L)
        break
    chunk_cvg <- coverage(chunk)
    if (is.null(cvg)) {
        cvg <- chunk_cvg
    } else {
        cvg <- cvg + chunk_cvg
    }
}
close(bf)
cvg
@


%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\subsection{How to compute read coverage}

The ``read coverage'' is the number of reads that cover a given genomic
position. Computing the read coverage generally consists in computing
the coverage at each position in the genome. This can be done with the
\Rcode{coverage()} function.

As sample data we use the \Biocpkg{pasillaBamSubset} data package
described in the introduction.

<<coverage_1>>=
library(pasillaBamSubset)
un1 <- untreated1_chr4()  # single-end reads
library(GenomicAlignments)
reads1 <- readGAlignments(un1)
cvg1 <- coverage(reads1)
cvg1
@

Coverage on chr4:

<<coverage_2>>=
cvg1$chr4
@

Average and max coverage:

<<coverage_3>>=
mean(cvg1$chr4)
max(cvg1$chr4)
@

Note that \Rcode{coverage()} is a generic function with methods for
different types of objects. See \Rcode{?coverage} for more information.


%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\subsection{How to find peaks in read coverage}

ChIP-Seq analysis usually involves finding peaks in read coverage.
This process is sometimes called ``peak calling'' or ``peak detection''.
Here we're only showing a naive way to find peaks in the object returned
by the \Rcode{coverage()} function. \Bioconductor{} packages
\Biocpkg{BayesPeak}, \Biocpkg{bumphunter}, \Biocpkg{Starr}, \Biocpkg{CexoR},
\Biocpkg{exomePeak}, \Biocpkg{RIPSeeker}, and others, provide sophisticated
peak calling tools for ChIP-Seq, RIP-Seq, and other kind of high throughput
sequencing data.

Let's assume \Rcode{cvg1} is the object returned by \Rcode{coverage()}
(see previous {\it HOWTO} for how to compute it). We can use the
\Rcode{slice()} function to find the genomic regions where the coverage
is greater or equal to a given threshold.

<<peaks_1>>=
chr4_peaks <- slice(cvg1$chr4, lower=500)
chr4_peaks
length(chr4_peaks)  # nb of peaks
@

The weight of a given peak can be defined as the number of aligned
nucleotides that belong to the peak (a.k.a. the area under the peak in
mathematics). It can be obtained with \Rcode{sum()}:

<<peaks_2>>=
sum(chr4_peaks)
@


%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\subsection{How to retrieve a gene model from the UCSC genome browser}

See introduction for a quick description of what \textit{gene models} and
\Rclass{TranscriptDb} objects are.
We can use the \Rcode{make\-Transcript\-Db\-From\-UCSC()} function from
the \Biocpkg{GenomicFeatures} package to import a UCSC genome browser
track as a \Rclass{TranscriptDb} object.

<<makeTranscriptDbFromUCSC_1, eval=FALSE>>=
library(GenomicFeatures)
### Internet connection required! Can take several minutes...
txdb <- makeTranscriptDbFromUCSC(genome="sacCer2", tablename="ensGene")
@

See \Rcode{?makeTranscriptDbFromUCSC} in the \Biocpkg{GenomicFeatures}
package for more information.

Note that some of the most frequently used gene models are available
as TxDb packages. A TxDb package consists of a pre-made \Rclass{TranscriptDb}
object wrapped into an annotation data package. Go to
\url{http://bioconductor.org/packages/release/BiocViews.html#\_\_\_TranscriptDb}
to browse the list of available TxDb packages.

<<TxDb.Hsapiens.UCSC.hg19.knownGene_1>>=
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
txdb <- TxDb.Hsapiens.UCSC.hg19.knownGene
txdb
@

Extract the transcript coordinates from this gene model:

<<TxDb.Hsapiens.UCSC.hg19.knownGene_2>>=
transcripts(txdb)
@


%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\subsection{How to retrieve a gene model from Ensembl}

See introduction for a quick description of what \textit{gene models} and
\Rclass{TranscriptDb} objects are.
We can use the \Rcode{make\-Transcript\-Db\-From\-Biomart()} function from
the \Biocpkg{GenomicFeatures} package to retrieve a gene model from the
Ensembl Mart.

<<makeTranscriptDbFromBiomart_1, eval=FALSE>>=
library(GenomicFeatures)
### Internet connection required! Can take several minutes...
txdb <- makeTranscriptDbFromBiomart(biomart="ensembl",
                                    dataset="hsapiens_gene_ensembl")
@

See \Rcode{?makeTranscriptDbFromBiomart} in the \Biocpkg{GenomicFeatures}
package for more information.

Note that some of the most frequently used gene models are available
as TxDb packages. A TxDb package consists of a pre-made \Rclass{TranscriptDb}
object wrapped into an annotation data package. Go to
\url{http://bioconductor.org/packages/release/BiocViews.html#\_\_\_TranscriptDb}
to browse the list of available TxDb packages.

<<TxDb.Athaliana.BioMart.plantsmart21_2>>=
library(TxDb.Athaliana.BioMart.plantsmart21)
txdb <- TxDb.Athaliana.BioMart.plantsmart21
txdb
@

Extract the exon coordinates from this gene model:

<<TxDb.Athaliana.BioMart.plantsmart21_2>>=
exons(txdb)
@


%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\subsection{How to load a gene model from a GFF or GTF file}

See introduction for a quick description of what \textit{gene models} and
\Rclass{TranscriptDb} objects are.
We can use the \Rcode{make\-Transcript\-Db\-From\-GFF()} function from
the \Biocpkg{GenomicFeatures} package to import a GFF or GTF file as a
\Rclass{TranscriptDb} object.

<<makeTranscriptDbFromGFF_1>>=
library(GenomicFeatures)
gff_file <- system.file("extdata", "a.gff3", package="GenomicFeatures")
txdb <- makeTranscriptDbFromGFF(gff_file, format="gff3")
txdb
@

See \Rcode{?makeTranscriptDbFromGFF} in the \Biocpkg{GenomicFeatures}
package for more information.

Extract the exon coordinates grouped by gene from this gene model:

<<makeTranscriptDbFromGFF_2>>=
exonsBy(txdb, by="gene")
@


%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\subsection{How to retrieve a gene model from \Biocpkg{AnnotationHub}}

When a gene model is not available as a \Rclass{GRanges} or
\Rclass{GRangesList} object or as a \Bioconductor{} data package, it may
be available on \Biocpkg{AnnotationHub}.
In this {\it HOWTO}, will look for a gene model for Drosophila melanogaster
on \Biocpkg{AnnotationHub}.
Create a `hub' and filter on Drosophila melanogaster:

<<hub_1>>=
library(AnnotationHub)
### Internet connection required!
hub <- AnnotationHub()
filters(hub) <- list(Species="Drosophila melanogaster")
@

There are 87 files that match Drosophila melanogaster.

<<hub_2>>=
length(hub)
head(names(hub))
@

Retrieve a dm3 file as a \Rcode{GRanges}. 

<<hub_3>>=
gr <- hub$goldenpath.dm3.database.ensGene_0.0.1.RData
summary(gr)
@

The metadata fields contain the details of file origin and content.

<<hub_4>>=
names(metadata(gr)[[2]])
metadata(gr)[[2]]$Tags
@

Split the \Rclass{GRanges} object by gene name to get a \Rclass{GRangesList}
object of transcript ranges grouped by gene.

<<hub_5>>= 
txbygn <- split(gr, gr$name)
@

You can now use \Rcode{txbygn} with the \Rcode{summarizeOverlaps} function
to prepare a table of read counts for RNA-Seq differential gene expression.

Note that before passing \Rcode{txbygn} to \Rcode{summarizeOverlaps},
you should confirm that the seqlevels (chromosome names) in it match those
in the BAM file. See \Rcode{?renameSeqlevels}, \Rcode{?keepSeqlevels}
and \Rcode{?seqlevels} for examples of renaming seqlevels.

 
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\subsection{How to annotate peaks in read coverage}

[coming soon...]


%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\subsection{How to prepare a table of read counts for RNA-Seq differential
            gene expression}

Methods for RNA-Seq gene expression analysis generally require a table of
counts that summarize the number of reads that overlap or `hit' a 
particular gene. In this {\it HOWTO} we count with the
\Rcode{summarizeOverlaps} function from the \Biocpkg{GenomicAlignments}
package and create a count table from the results. 

Other packages that provide read counting are \Biocpkg{Rsubread} and 
\Biocpkg{easyRNASeq}. The \Biocpkg{parathyroidSE} package vignette 
contains a workflow on counting and other common operations required for 
differential expression analysis. 

As sample data we use the \Biocpkg{pasillaBamSubset} data package
described in the introduction.

<<count_1>>=
library(pasillaBamSubset)
un1 <- untreated1_chr4()  # single-end reads
@

\Rcode{summarizeOverlaps} requires the name of a BAM file(s) and a
{\textit gene model} to count against. See introduction for a quick
description of what a \textit{gene models} is.
The gene model must match the genome build the reads in the BAM file were
aligned to. For the pasilla data this is dm3 Dmelanogaster which is
available as a \Bioconductor{} package. Load the package and extract
the exon ranges grouped by gene:

<<count_2>>=
library(TxDb.Dmelanogaster.UCSC.dm3.ensGene)
exbygene <- exonsBy(TxDb.Dmelanogaster.UCSC.dm3.ensGene, "gene")
@

\Rcode{exbygene} is a \Rclass{GRangesList} object with one list
element per gene in the gene model.

\Rcode{summarizeOverlaps} automatically sets a \Rcode{yieldSize}
on large BAM files and iterates over them in chunks. When reading
paired-end data set the \Rcode{singleEnd} argument to FALSE.
See ?\Rfunction{summarizeOverlaps} for details reguarding the
count \Rcode{modes} and additional arguments. 

<<count_3>>=
library(GenomicAlignments)
se <- summarizeOverlaps(exbygene, un1, mode="IntersectionNotEmpty")
@

The return object is a \Rcode{SummarizedExperiment} with counts in 
the \Rcode{assays} slot.

<<count_4>>=
class(se)
head(table(assays(se)$counts))
@

The count vector is the same length as \Rcode{exbygene}:

<<count_5>>=
identical(length(exbygene), length(assays(se)$counts))
@

A copy of \Rcode{exbygene} is stored in the \Rcode{rowData} slot:

<<count_6>>=
rowData(se)
@

Two popular packages for RNA-Seq differential gene expression are
\Biocpkg{DESeq} and \Biocpkg{edgeR}. Tables of counts per gene are required
for both and can be easily created with a vector of counts. Here we use the
counts from our \Rclass{SummarizedExperiment} object:

<<count_table>>=
library(DESeq)
deseq <- newCountDataSet(assays(se)$counts, rownames(colData(se)))
library(edgeR)
edger <- DGEList(assays(se)$counts, group=rownames(colData(se)))
@


%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\subsection{How to summarize junctions from a BAM file containing RNA-Seq
            reads}

As sample data we use the \Biocpkg{pasillaBamSubset} data package
described in the introduction.

<<summarizeJunctions_1>>=
library(pasillaBamSubset)
un1 <- untreated1_chr4()  # single-end reads
library(GenomicAlignments)
reads1 <- readGAlignments(un1)
reads1
@

For each alignment, the aligner generated a CIGAR string that describes
its "geometry", that is, the locations of insertions, deletions and
junctions in the alignment. See the SAM Spec available on the SAMtools
website for the details (\url{http://samtools.sourceforge.net/}).

The \Rcode{summarizeJunctions()} function from the
\Biocpkg{GenomicAlignments} package can be used to summarize the junctions
in \Rcode{reads1}.

<<summarizeJunctions_2>>=
junc_summary <- summarizeJunctions(reads1)
junc_summary
@

See \Rcode{?summarizeJunctions} in the \Biocpkg{GenomicAlignments}
package for more information.


%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\subsection{How to get the exon and intron sequences of a given gene}

The exon and intron sequences of a gene are essentially the DNA
sequences of the introns and exons of all known transcripts of the gene.
The first task is to identify all transcripts associated with the gene of
interest. Our sample gene is the human TRAK2 which is involved in
regulation of endosome-to-lysosome trafficking of membrane cargo.
The Entrez gene id is `66008'. 

<<trak_1>>=
trak2 <- "66008"
@

The \Biocpkg{TxDb.Hsapiens.UCSC.hg19.knownGene} data package contains the
gene model corresponding to the UCSC `Known Genes' track.

<<trak_2>>=
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
txdb <- TxDb.Hsapiens.UCSC.hg19.knownGene
@

The transcript ranges for all the genes in the gene model can be extracted
with the \Rfunction{transcriptsBy} function from the \Biocpkg{GenomicFeatures}
package. They will be returned in a named \Rclass{GRangesList} object
containing all the transcripts grouped by gene. In order to keep only the
transcripts of the TRAK2 gene we will subset the \Rclass{GRangesList}
object using the \Rcode{[[} operator.

<<trak_3>>=
library(GenomicFeatures)
trak2_txs <- transcriptsBy(txdb, by="gene")[[trak2]]
trak2_txs
@

\Rcode{trak2\_txs} is a \Rclass{GRanges} object with one range per
transcript in the TRAK2 gene. The transcript names are stored in the 
\Rcode{tx\_name} metadata column. We will need them to subset the
extracted intron and exon regions:

<<trak_4>>=
trak2_tx_names <- mcols(trak2_txs)$tx_name
trak2_tx_names
@

The exon and intron genomic ranges for all the transcripts in the gene model
can be extracted with the \Rfunction{exonsBy} and
\Rfunction{intronsByTranscript} functions, respectively. Both functions
return a \Rclass{GRangesList} object. Then we keep only the exon and intron
for the transcripts of the TRAK2 gene by subsetting each \Rclass{GRangesList}
object by the TRAK2 transcript names.

Extract the exon regions:

<<trak_5>>=
trak2_exbytx <- exonsBy(txdb, "tx", use.names=TRUE)[trak2_tx_names]
elementLengths(trak2_exbytx)
@

... and the intron regions:

<<trak_7>>=
trak2_inbytx <- intronsByTranscript(txdb, use.names=TRUE)[trak2_tx_names]
elementLengths(trak2_inbytx)
@

Next we want the DNA sequences for these exons and introns.
The \Rfunction{getSeq} function from the \Biocpkg{Biostrings} package can
be used to query a \Biocpkg{BSgenome} object with a set of genomic ranges
and retrieve the corresponding DNA sequences. 

<<trak_8>>=
library(BSgenome.Hsapiens.UCSC.hg19)
@

Extract the exon sequences:

<<trak_9>>=
trak2_ex_seqs <- getSeq(Hsapiens, trak2_exbytx)
trak2_ex_seqs
trak2_ex_seqs[["uc002uyb.4"]]
trak2_ex_seqs[["uc002uyc.2"]]
@

... and the intron sequences:

<<trak_10>>=
trak2_in_seqs <- getSeq(Hsapiens, trak2_inbytx)
trak2_in_seqs
trak2_in_seqs[["uc002uyb.4"]]
trak2_in_seqs[["uc002uyc.2"]]
@


%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\subsection{How to get the CDS and UTR sequences of genes associated 
            with colorectal cancer}

In this {\it HOWTO} we extract the CDS and UTR sequences of genes involved
in colorectal cancer. The workflow extends the ideas presented in the
previous {\it HOWTO} and suggests an approach for identifying
disease-related genes.

\subsubsection{Build a gene list}

We start with a list of gene or transcript ids. If you do not have 
pre-defined list one can be created with the \Biocpkg{KEGG.db} and 
\Biocpkg{KEGGgraph} packages. Updates to the data in the \Biocpkg{KEGG.db} 
package are no longer available, however, the resource is still useful for 
identifying pathway names and ids. 

Create a table of KEGG pathways and ids and search on the term `cancer'.

<<cancer_1>>=
library(KEGG.db)
pathways <- toTable(KEGGPATHNAME2ID)
pathways[grepl("cancer", pathways$path_name, fixed=TRUE),] 
@

Use the "05210" id to query the KEGG web resource (accesses the currently
maintained data).

<<cancer_2>>=
library(KEGGgraph)
dest <- tempfile()
retrieveKGML("05200", "hsa", dest, "internal")
@

The suffix of the KEGG id is the Entrez gene id. The 
\Rfunction{translateKEGGID2GeneID} simply removes the prefix leaving 
just the Entrez gene ids.

<<cancer_3>>=
crids <- as.character(parseKGML2DataFrame(dest)[,1])
crgenes <- unique(translateKEGGID2GeneID(crids))
head(crgenes)
@

\subsubsection{Identify genomic coordinates}

The list of gene ids is used to extract genomic positions of the regions 
of interest. The Known Gene table from UCSC will be the annotation and
is available as a \Bioconductor{} package.

<<cancer_4>>=
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
txdb <- TxDb.Hsapiens.UCSC.hg19.knownGene
@

If an annotation is not available as a \Bioconductor{} annotation package
it may be available in \Biocpkg{AnnotationHub}. Additionally, there are
functions in \Biocpkg{GenomicFeatures} which can retrieve data from UCSC and
Ensembl to create a \Robject{TranscriptDb}. See 
\Rcode{?makeTranscriptDbFromUCSC} for more information.

As in the previous {\it HOWTO} we need to identify the transcripts 
corresponding to each gene. The transcript id (or name) is used
to isolate the UTR and coding regions of interest. This grouping of
transcript by gene is also used to re-group the final sequence results.

The \Rcode{transcriptsBy} function outputs both the gene and transcript
identifiers which we use to create a map between the two. The 
\Rcode{map} is a \Robject{CharacterList} with gene ids as names and 
transcript ids as the list elements.

<<cancer_5>>=
txbygene <- transcriptsBy(txdb, "gene")[crgenes] ## subset on colorectal genes
map <- relist(unlist(txbygene, use.names=FALSE)$tx_id, txbygene)
map
@

Extract the UTR and coding regions.

<<cancer_6>>=
cds <- cdsBy(txdb, "tx")
threeUTR <- threeUTRsByTranscript(txdb)
fiveUTR <- fiveUTRsByTranscript(txdb)
@

Coding and UTR regions may not be present for all transcripts specified 
in \Rcode{map}. Consequently, the subset results will not be the same 
length. This length discrepancy must be taken into account when re-listing 
the final results by gene.

<<cancer_7>>=
txid <- unlist(map, use.names=FALSE)
cds <- cds[names(cds) %in% txid]
threeUTR <- threeUTR[names(threeUTR) %in% txid]
fiveUTR <- fiveUTR[names(fiveUTR) %in% txid]
@

Note the different lengths of the subset regions.

<<cancer_8>>=
length(txid) ## all possible transcripts
length(cds)
length(threeUTR)
length(fiveUTR)
@

These objects are \Robject{GRangesList}s with the transcript id as the 
outer list element. 

<<cancer_9>>=
cds
@

\subsubsection{Extract sequences from BSgenome}

The \Rcode{BSgenome} packages contain complete genome sequences
for a given organism.

Load the \Rcode{BSgenome} package for homo sapiens.

<<cancer_10>>=
library(BSgenome.Hsapiens.UCSC.hg19)
genome <- BSgenome.Hsapiens.UCSC.hg19
@

Use \Rfunction{extractTranscriptSeqs} to extract the UTR and coding 
regions from the \Rcode{BSgenome}. This function retrieves the sequences 
for an any \Robject{GRanges} or \Robject{GRangesList} (i.e., not just 
transcripts like the name implies).

<<cancer_11>>=
threeUTR_seqs <- extractTranscriptSeqs(genome, threeUTR) 
fiveUTR_seqs <- extractTranscriptSeqs(genome, fiveUTR) 
cds_seqs <- extractTranscriptSeqs(genome, cds) 
@

The return values are \Robject{DNAStringSet} objects.

<<cancer_12>>=
cds_seqs
@

Our final step is to collect the coding and UTR regions (currently 
organzied by transcript) into groups by gene id. The \Rfunction{relist} 
function groups the sequences of a \Robject{DNAStringSet} object into
a \Robject{DNAStringSetList} object, based on the specified \Rcode{skeleton}
argument. The \Rcode{skeleton} must be a list-like object and only its shape
(i.e. its element lengths) matters (its exact content is ignored). A simple
form of \Rcode{skeleton} is to use a partitioning object that we make by
specifying the size of each partition. The partitioning objects are different
for each type of region because not all transcripts had a coding or 3' or 5'
UTR region defined. 

<<cancer_13>>=
lst3 <- relist(threeUTR_seqs, PartitioningByWidth(sum(map %in% names(threeUTR))))
lst5 <- relist(fiveUTR_seqs, PartitioningByWidth(sum(map %in% names(fiveUTR))))
lstc <- relist(cds_seqs, PartitioningByWidth(sum(map %in% names(cds))))
@

There are 239 genes in \Rcode{map} each of which have 1 or more transcripts. 
The table of element lengths shows how many genes have each number of
transcripts. For example, 47 genes have 1 transcript, 48 genes have 2 etc.

<<cancer_14>>=
length(map)
table(elementLengths(map))
@

The lists of DNA sequences all have the same length as \Rcode{map} but one or
more of the element lengths may be zero. This would indicate that data were
not available for that gene. The tables below show that there was at least
1 coding region available for all genes (i.e., none of the element lengths
are 0). However, both the 3' and 5' UTR results have element lengths of 0 
which indicates no UTR data were available for that gene.

<<cancer_15>>=
table(elementLengths(lstc))
table(elementLengths(lst3))
names(lst3)[elementLengths(lst3) == 0L] ## genes with no 3' UTR data
table(elementLengths(lst5))
names(lst5)[elementLengths(lst5) == 0L] ## genes with no 5' UTR data
@


%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\subsection{How to create DNA consensus sequences for read group `families'}

The motivation for this {\it HOWTO} comes from a study which explored the
dynamics of point mutations. The mutations of interest exist with a range 
of frequencies in the control group (e.g., 0.1\% - 50\%). PCR and sequencing 
error rates make it difficult to identify low frequency events 
(e.g., < 20\%).

When a library is prepared with Nextera, random fragments are generated 
followed by a few rounds of PCR. When the genome is large enough, reads 
aligning to the same start position are likely descendant from the same 
template fragment and should have identical sequences. 

The goal is to elimininate noise by grouping the reads by common start 
position and discarding those that do not exceed a certain threshold within 
each family. A new consensus sequence will be created for each read group
family.

\subsubsection{Sort reads into groups by start position}

Load the BAM file into a GAlignments object.

<<cseq_1>>=
library(Rsamtools)
bamfile <- system.file("extdata", "ex1.bam", package="Rsamtools")
param <- ScanBamParam(what=c("seq", "qual"))
library(GenomicAlignments)
gal <- readGAlignments(bamfile, use.names=TRUE, param=param)
@

Use the \Rfunction{sequenceLayer} function to {\it lay} the query sequences
and quality strings on the reference.

<<cseq_2>>=
qseq <- setNames(mcols(gal)$seq, names(gal))
qual <- setNames(mcols(gal)$qual, names(gal))
qseq_on_ref <- sequenceLayer(qseq, cigar(gal),
                             from="query", to="reference")
qual_on_ref <- sequenceLayer(qual, cigar(gal),
                             from="query", to="reference")
@

Split by chromosome.

<<cseq_3>>=
qseq_on_ref_by_chrom <- splitAsList(qseq_on_ref, seqnames(gal))
qual_on_ref_by_chrom <- splitAsList(qual_on_ref, seqnames(gal))
pos_by_chrom <- splitAsList(start(gal), seqnames(gal))
@

For each chromosome generate one GRanges object that contains
unique alignment start positions and attach 3 metadata columns
to it: the number of reads, the query sequences, and the quality
strings.

<<cseq_4>>=
gr_by_chrom <- lapply(seqlevels(gal),
  function(seqname)
  {
    qseq_on_ref2 <- qseq_on_ref_by_chrom[[seqname]]
    qual_on_ref2 <- qual_on_ref_by_chrom[[seqname]]
    pos2 <- pos_by_chrom[[seqname]]
    qseq_on_ref_per_pos <- split(qseq_on_ref2, pos2)
    qual_on_ref_per_pos <- split(qual_on_ref2, pos2)
    nread <- elementLengths(qseq_on_ref_per_pos)
    gr_mcols <- DataFrame(nread=unname(nread),
                          qseq_on_ref=unname(qseq_on_ref_per_pos),
                          qual_on_ref=unname(qual_on_ref_per_pos))
    gr <- GRanges(Rle(seqname, nrow(gr_mcols)),
                  IRanges(as.integer(names(nread)), width=1))
    mcols(gr) <- gr_mcols
    seqlevels(gr) <- seqlevels(gal)
    gr
  })
@

Combine all the GRanges objects obtained in (4) in 1 big GRanges
object:

<<cseq_5>>=
gr <- do.call(c, gr_by_chrom)
seqinfo(gr) <- seqinfo(gal)
@

`gr' is a GRanges object that contains unique alignment start positions:

<<cseq_6>>=
gr[1:6]
@

Look at qseq\_on\_ref and qual\_on\_ref.

<<cseq_7>>= 
qseq_on_ref
qual_on_ref
@

2 reads align to start position 13. Let's have a close look at their 
sequences:

<<cseq_8>>=
mcols(gr)$qseq_on_ref[[6]]
@

and their qualities:

<<cseq_9>>=
mcols(gr)$qual_on_ref[[6]]
@

Note that the sequence and quality strings are those projected to the 
reference so the first letter in those strings are on top of start 
position 13, the 2nd letter on top of position 14, etc...

\subsubsection{Remove low frequency reads}

For each start position, remove reads with and under-represented sequence 
(e.g. threshold = 20\% for the data used here which is low coverage).
A unique number is assigned to each unique sequence. This will make
future calculations easier and a little bit faster.

<<cseq_10>>=
qseq_on_ref <- mcols(gr)$qseq_on_ref
tmp <- unlist(qseq_on_ref, use.names=FALSE)
qseq_on_ref_id <- relist(match(tmp, tmp), qseq_on_ref)
@

Quick look at `qseq\_on\_ref\_id':
It's an IntegerList object with the same length and "shape"
as `qseq\_on\_ref'.

<<cseq_11>>=
qseq_on_ref_id
@

Remove the under represented ids from each list element of `qseq\_on\_ref\_id':

<<cseq_12>>=
qseq_on_ref_id2 <- endoapply(qseq_on_ref_id,
    function(ids) ids[countMatches(ids, ids) >= 0.2 * length(ids)])
@

Remove corresponding sequences from `qseq\_on\_ref':

<<cseq_13>>=
tmp <- unlist(qseq_on_ref_id2, use.names=FALSE)
qseq_on_ref2 <- relist(unlist(qseq_on_ref, use.names=FALSE)[tmp],
                       qseq_on_ref_id2)
@

\subsubsection{Create a consensus sequence for each read group family}

Compute 1 consensus matrix per chromosome:

<<cseq_14>>=
split_factor <- rep.int(seqnames(gr), elementLengths(qseq_on_ref2))
qseq_on_ref2 <- unlist(qseq_on_ref2, use.names=FALSE)
qseq_on_ref2_by_chrom <- splitAsList(qseq_on_ref2, split_factor)
qseq_pos_by_chrom <- splitAsList(start(gr), split_factor)

cm_by_chrom <- lapply(names(qseq_pos_by_chrom),
    function(seqname)
        consensusMatrix(qseq_on_ref2_by_chrom[[seqname]],
                        as.prob=TRUE,
                        shift=qseq_pos_by_chrom[[seqname]]-1,
                        width=seqlengths(gr)[[seqname]]))
names(cm_by_chrom) <- names(qseq_pos_by_chrom)
@

'cm\_by\_chrom' is a list of consensus matrices. Each matrix has 17 rows 
(1 per letter in the DNA alphabet) and 1 column per chromosome position.

<<cseq_15>>=
lapply(cm_by_chrom, dim)
@

Compute the consensus string from each consensus matrix. We'll put "+" 
in the strings wherever there is no coverage for that position, and "N" 
where there is coverage but no consensus.

<<cseq_16>>=
cs_by_chrom <- lapply(cm_by_chrom,
    function(cm) {
        ## need to "fix" 'cm' because consensusString()
        ## doesn't like consensus matrices with columns
        ## that contain only zeroes (e.g., chromosome
        ## positions with no coverage)
        idx <- colSums(cm) == 0L
        cm["+", idx] <- 1
        DNAString(consensusString(cm, ambiguityMap="N"))
    })
@

The new consensus strings.

<<cseq_17>>=
cs_by_chrom
@


%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\subsection{How to compute binned averages along a genome}

In some applications, there is the need to compute the average of a 
variable along a genome for a set of predefined fixed-width regions 
(sometimes called "bins"). One such example is coverage.
Coverage is an \Robject{RleList} with one list element per
chromosome. Here we simulate a coverage list.

<<bim_1>>=
library(BSgenome.Scerevisiae.UCSC.sacCer2)
set.seed(22)
cov <- RleList(
    lapply(seqlengths(Scerevisiae),
           function(len) Rle(sample(-10:10, len, replace=TRUE))),
    compress=FALSE)
head(cov, 3)
@

Use the \Rfunction{tileGenome} function to create a set of bins along
the genome.

<<bin_2>>=
bins1 <- tileGenome(seqinfo(Scerevisiae), tilewidth=100,
                    cut.last.tile.in.chrom=TRUE)
@

We define the following function to compute the binned average of a
numerical variable defined along a genome.

\begin{verbatim}
Arguments:
  'bins': a GRanges object representing the genomic bins.
       Typically obtained by calling tileGenome() with
       'cut.last.tile.in.chrom=TRUE'.
  'numvar': a named RleList object representing a numerical
       variable defined along the genome covered by 'bins', which
       is the genome described by 'seqinfo(bins)'.
  'mcolname': the name to give to the metadata column that will
       contain the binned average in the returned object.
\end{verbatim}

The function returns `bins' with an additional metadata column named 
`mcolname' containing the binned average.

<<bin_3>>=
binnedAverage <- function(bins, numvar, mcolname)
{
    stopifnot(is(bins, "GRanges"))
    stopifnot(is(numvar, "RleList"))
    stopifnot(identical(seqlevels(bins), names(numvar)))
    bins_per_chrom <- split(ranges(bins), seqnames(bins))
    means_list <- lapply(names(numvar),
        function(seqname) {
            views <- Views(numvar[[seqname]],
                           bins_per_chrom[[seqname]])
            viewMeans(views)
        })
    new_mcol <- unsplit(means_list, as.factor(seqnames(bins)))
    mcols(bins)[[mcolname]] <- new_mcol
    bins
}
@

Compute the binned average for `cov':

<<bin_4>>=
bins1 <- binnedAverage(bins1, cov, "binned_cov")
bins1
@

The bin size can be modified with the \Rcode{tilewidth} argument
to \Rfunction{tileGenome}. For additional examples see
?\Rfunction{tileGenome}.



%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\section{Session Information}

<<SessionInfo, echo=FALSE>>=
sessionInfo()
@



%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

\bibliography{GenomicRanges}
\bibliographystyle{plainnat}

\end{document}