In this manual, we will show how to use the methylKit package. methylKit is an R package for analysis and annotation of DNA methylation information obtained by high-throughput bisulfite sequencing. The package is designed to deal with sequencing data from RRBS and its variants. But, it can potentially handle whole-genome bisulfite sequencing data if proper input format is provided.
DNA methylation in vertebrates typically occurs at CpG dinucleotides, however non-CpG Cs are also methylated in certain tissues such as embryonic stem cells. DNA methylation can act as an epigenetic control mechanism for gene regulation. Methylation can hinder binding of transcription factors and/or methylated bases can be bound by methyl-binding-domain proteins which can recruit chromatin remodeling factors. In both cases, the transcription of the regulated gene will be effected. In addition, aberrant DNA methylation patterns have been associated with many human malignancies and can be used in a predictive manner. In malignant tissues, DNA is either hypo-methylated or hyper-methylated compared to the normal tissue. The location of hyper- and hypo-methylated sites gives a distinct signature to many diseases. Traditionally, hypo-methylation is associated with gene transcription (if it is on a regulatory region such as promoters) and hyper-methylation is associated with gene repression.
Bisulfite sequencing is a technique that can determine DNA methylation patterns. The major difference from regular sequencing experiments is that, in bisulfite sequencing DNA is treated with bisulfite which converts cytosine residues to uracil, but leaves 5-methylcytosine residues unaffected. By sequencing and aligning those converted DNA fragments it is possible to call methylation status of a base. Usually, the methylation status of a base determined by a high-throughput bisulfite sequencing will not be a binary score, but it will be a percentage. The percentage simply determines how many of the bases that are aligning to a given cytosine location in the genome have actual C bases in the reads. Since bisulfite treatment leaves methylated Cs intact, that percentage will give us percent methylation score on that base. The reasons why we will not get a binary response are:
We start by reading in the methylation call data from bisulfite sequencing with read
function. Reading in the data this way will return a methylRawList
object which stores methylation information per sample for each covered base. The methylation call files are basically text files that contain percent methylation score per base. Such input files may be obtained from AMP pipeline developed for aligning RRBS reads or from processBismarkAln
function. . However, “cytosineReport” and “coverage” files from Bismark aligner can be read in to methylKit as well.
A typical methylation call file looks like this:
## chrBase chr base strand coverage freqC freqT
## 1 chr21.9764539 chr21 9764539 R 12 25.00 75.00
## 2 chr21.9764513 chr21 9764513 R 12 0.00 100.00
## 3 chr21.9820622 chr21 9820622 F 13 0.00 100.00
## 4 chr21.9837545 chr21 9837545 F 11 0.00 100.00
## 5 chr21.9849022 chr21 9849022 F 124 72.58 27.42
Most of the time bisulfite sequencing experiments have test and control samples. The test samples can be from a disease tissue while the control samples can be from a healthy tissue. You can read a set of methylation call files that have test/control conditions giving treatment
vector option. For sake of subsequent analysis, file.list, sample.id and treatment option should have the same order. In the following example, first two files are have the sample ids “test1” and “test2” and as determined by treatment vector they belong to the same group. The third and fourth files have sample ids “ctrl1” and “ctrl2” and they belong to the same group as indicated by the treatment vector.
library(methylKit)
file.list=list( system.file("extdata",
"test1.myCpG.txt", package = "methylKit"),
system.file("extdata",
"test2.myCpG.txt", package = "methylKit"),
system.file("extdata",
"control1.myCpG.txt", package = "methylKit"),
system.file("extdata",
"control2.myCpG.txt", package = "methylKit") )
# read the files to a methylRawList object: myobj
myobj=methRead(file.list,
sample.id=list("test1","test2","ctrl1","ctrl2"),
assembly="hg18",
treatment=c(1,1,0,0),
context="CpG"
)
In addition to the options we mentioned above, any tab separated text file with a generic format can be read in using methylKit, such as methylation ratio files from BSMAP. See here for an example.
Sometimes, when dealing with multiple samples and increased sample sizes coming from genome wide bisulfite sequencing experiments, the memory of your computer might not be sufficient enough.
Therefore methylKit offers a new group of classes, that are basically pendants to the original methylKit classes with one important difference: The methylation information, which normally is internally stored as data.frame, is stored in an external bgzipped file and is indexed by tabix (H. Li 2011), to enable fast retrieval of records or regions. This group contains methylRawListDB
, methylRawDB
, methylBaseDB
and methylDiffDB
, let us call them methylDB objects.
We can now create a methylRawListDB
object, which stores the same content as myobj from above. But the single methylRaw
objects retrieve their data from the tabix-file linked under dbpath
.
library(methylKit)
file.list=list( system.file("extdata", "test1.myCpG.txt", package = "methylKit"),
system.file("extdata", "test2.myCpG.txt", package = "methylKit"),
system.file("extdata", "control1.myCpG.txt", package = "methylKit"),
system.file("extdata", "control2.myCpG.txt", package = "methylKit") )
# read the files to a methylRawListDB object: myobjDB
# and save in databases in folder methylDB
myobjDB=methRead(file.list,
sample.id=list("test1","test2","ctrl1","ctrl2"),
assembly="hg18",
treatment=c(1,1,0,0),
context="CpG",
dbtype = "tabix",
dbdir = "methylDB"
)
print(myobjDB[[1]]@dbpath)
## [1] "/tmp/RtmpTpOekb/Rbuild5ebe4d8329c7/methylKit/vignettes/methylDB/test1.txt.bgz"
Most if not all functions in this package will work with methylDB objects the same way as it does with normal methylKit objects. Functions that return methylKit objects, will return a methylDB object if provided, but there are a few exceptions such as the select
, the [
and the selectByOverlap
functions.
Alternatively, methylation percentage calls can be calculated from sorted SAM or BAM file(s) from Bismark aligner and read-in to the memory. Bismark is a popular aligner for bisulfite sequencing reads, available here (Krueger and Andrews 2011). processBismarkAln
function is designed to read-in Bismark SAM/BAM files as methylRaw
or methylRawList
objects which store per base methylation calls. SAM files must be sorted by chromosome and read position columns, using ‘sort’ command in unix-like machines will accomplish such a sort easily. BAM files should be sorted and indexed. This could be achieved with samtools (http://www.htslib.org/doc/samtools.html).
The following command reads a sorted SAM file and creates a methylRaw
object for CpG methylation.The user has the option to save the methylation call files to a folder given by save.folder
option. The saved files can be read-in using the read
function when needed.
my.methRaw=processBismarkAln( location =
system.file("extdata",
"test.fastq_bismark.sorted.min.sam",
package = "methylKit"),
sample.id="test1", assembly="hg18",
read.context="CpG",save.folder=getwd())
It is also possible to read multiple SAM files at the same time, check processBismarkAln
documentation.
Since we read the methylation data now, we can check the basic stats about the methylation data such as coverage and percent methylation. We now have a methylRawList
object which contains methylation information per sample. The following command prints out percent methylation statistics for second sample: “test2”
getMethylationStats(myobj[[2]],plot=FALSE,both.strands=FALSE)
## methylation statistics per base
## summary:
## Min. 1st Qu. Median Mean 3rd Qu. Max.
## 0.00 20.00 82.79 63.17 94.74 100.00
## percentiles:
## 0% 10% 20% 30% 40% 50% 60%
## 0.00000 0.00000 0.00000 48.38710 70.00000 82.78556 90.00000
## 70% 80% 90% 95% 99% 99.5% 99.9%
## 93.33333 96.42857 100.00000 100.00000 100.00000 100.00000 100.00000
## 100%
## 100.00000
The following command plots the histogram for percent methylation distribution.The figure below is the histogram and numbers on bars denote what percentage of locations are contained in that bin. Typically, percent methylation histogram should have two peaks on both ends. In any given cell, any given base are either methylated or not. Therefore, looking at many cells should yield a similar pattern where we see lots of locations with high methylation and lots of locations with low methylation.
getMethylationStats(myobj[[2]],plot=TRUE,both.strands=FALSE)
We can also plot the read coverage per base information in a similar way, again numbers on bars denote what percentage of locations are contained in that bin. Experiments that are highly suffering from PCR duplication bias will have a secondary peak towards the right hand side of the histogram.
getCoverageStats(myobj[[2]],plot=TRUE,both.strands=FALSE)
It might be useful to filter samples based on coverage. Particularly, if our samples are suffering from PCR bias it would be useful to discard bases with very high read coverage. Furthermore, we would also like to discard bases that have low read coverage, a high enough read coverage will increase the power of the statistical tests. The code below filters a methylRawList
and discards bases that have coverage below 10X and also discards the bases that have more than 99.9th percentile of coverage in each sample.
filtered.myobj=filterByCoverage(myobj,lo.count=10,lo.perc=NULL,
hi.count=NULL,hi.perc=99.9)
In order to do further analysis, we will need to get the bases covered in all samples. The following function will merge all samples to one object for base-pair locations that are covered in all samples. Setting destrand=TRUE
(the default is FALSE) will merge reads on both strands of a CpG dinucleotide. This provides better coverage, but only advised when looking at CpG methylation (for CpH methylation this will cause wrong results in subsequent analyses). In addition, setting destrand=TRUE
will only work when operating on base-pair resolution, otherwise setting this option TRUE will have no effect. The unite()
function will return a methylBase
object which will be our main object for all comparative analysis. The methylBase
object contains methylation information for regions/bases that are covered in all samples.
meth=unite(myobj, destrand=FALSE)
Let us take a look at the data content of methylBase object:
head(meth)
## chr start end strand coverage1 numCs1 numTs1 coverage2 numCs2
## 1 chr21 9853296 9853296 + 17 10 7 333 268
## 2 chr21 9853326 9853326 + 17 12 5 329 249
## 3 chr21 9860126 9860126 + 39 38 1 83 78
## 4 chr21 9906604 9906604 + 68 42 26 111 97
## 5 chr21 9906616 9906616 + 68 52 16 111 104
## 6 chr21 9906619 9906619 + 68 59 9 111 109
## numTs2 coverage3 numCs3 numTs3 coverage4 numCs4 numTs4
## 1 65 18 16 2 395 341 54
## 2 79 16 14 2 379 284 95
## 3 5 83 83 0 41 40 1
## 4 14 23 18 5 37 33 4
## 5 7 23 14 9 37 27 10
## 6 2 22 18 4 37 29 8
By default, unite
function produces bases/regions covered in all samples. That requirement can be relaxed using “min.per.group” option in unite
function.
# creates a methylBase object,
# where only CpGs covered with at least 1 sample per group will be returned
# there were two groups defined by the treatment vector,
# given during the creation of myobj: treatment=c(1,1,0,0)
meth.min=unite(myobj,min.per.group=1L)
We can check the correlation between samples using getCorrelation
. This function will either plot scatter plot and correlation coefficients or just print a correlation matrix
getCorrelation(meth,plot=TRUE)
## test1 test2 ctrl1 ctrl2
## test1 1.0000000 0.9252530 0.8767865 0.8737509
## test2 0.9252530 1.0000000 0.8791864 0.8801669
## ctrl1 0.8767865 0.8791864 1.0000000 0.9465369
## ctrl2 0.8737509 0.8801669 0.9465369 1.0000000
We can cluster the samples based on the similarity of their methylation profiles. The following function will cluster the samples and draw a dendrogram.
clusterSamples(meth, dist="correlation", method="ward", plot=TRUE)
## The "ward" method has been renamed to "ward.D"; note new "ward.D2"
##
## Call:
## hclust(d = d, method = HCLUST.METHODS[hclust.method])
##
## Cluster method : ward.D
## Distance : pearson
## Number of objects: 4
Setting the plot=FALSE
will return a dendrogram object which can be manipulated by users or fed in to other user functions that can work with dendrograms.
hc = clusterSamples(meth, dist="correlation", method="ward", plot=FALSE)
We can also do a PCA analysis on our samples. The following function will plot a scree plot for importance of components.
PCASamples(meth, screeplot=TRUE)
We can also plot PC1 and PC2 axis and a scatter plot of our samples on those axis which will reveal how they cluster.
PCASamples(meth)
We have implemented some rudimentary functionality for batch effect control. You can check which one of the principal components are statistically associated with the potential batch effects such as batch processing dates, age of subjects, sex of subjects using assocComp
. The function gets principal components from the percent methylation matrix derived from the input methylBase
object, and checks for association. The tests for association are either via Kruskal-Wallis test or Wilcoxon test for categorical attributes and correlation test for numerical attributes for samples such as age. If you are convinced that some principal components are accounting for batch effects, you can remove those principal components from your data using removeComp
.
# make some batch data frame
# this is a bogus data frame
# we don't have batch information
# for the example data
sampleAnnotation=data.frame(batch_id=c("a","a","b","b"),
age=c(19,34,23,40))
as=assocComp(mBase=meth,sampleAnnotation)
as
## $pcs
## PC1 PC2 PC3 PC4
## test1 -0.4978699 -0.5220504 0.68923849 -0.06737363
## test2 -0.4990924 -0.4805506 -0.71827964 0.06365693
## ctrl1 -0.5016543 0.4938800 0.08068700 0.70563101
## ctrl2 -0.5013734 0.5026102 -0.05014261 -0.70249091
##
## $vars
## [1] 92.271885 4.525328 1.870950 1.331837
##
## $association
## PC1 PC2 PC3 PC4
## batch_id 0.3333333 0.3333333 1.0000000 1.0000000
## age 0.5864358 0.6794346 0.3140251 0.3467957
# construct a new object by removing the first pricipal component
# from percent methylation value matrix
newObj=removeComp(meth,comp=1)
In addition to the methods described above, if you have used other ways to correct for batch effects and obtained a corrected percent methylation matrix, you can use reconstruct
function to reconstruct a corrected methylBase
object. Users have to supply a corrected percent methylation matrix and methylBase
object (where the uncorrected percent methylation matrix obtained from) to the reconstruct
function. Corrected percent methylation matrix should have the same row and column order as the original percent methylation matrix. All of these functions described in this section work on a methylBase
object that does not have missing values (that means all bases in methylBase object should have coverage in all samples).
mat=percMethylation(meth)
# do some changes in the matrix
# this is just a toy example
# ideally you want to correct the matrix
# for batch effects
mat[mat==100]=80
# reconstruct the methylBase from the corrected matrix
newobj=reconstruct(mat,meth)
For some situations, it might be desirable to summarize methylation information over tiling windows rather than doing base-pair resolution analysis. methylKit
provides functionality to do such analysis. The function below tiles the genome with windows 1000bp length and 1000bp step-size and summarizes the methylation information on those tiles. In this case, it returns a methylRawList
object which can be fed into unite
and calculateDiffMeth
functions consecutively to get differentially methylated regions. The tilling function adds up C and T counts from each covered cytosine and returns a total C and T count for each tile.
tiles=tileMethylCounts(myobj,win.size=1000,step.size=1000)
head(tiles[[1]],3)
## chr start end strand coverage numCs numTs
## 1 chr21 9764001 9765000 * 24 3 21
## 2 chr21 9820001 9821000 * 13 0 13
## 3 chr21 9837001 9838000 * 11 0 11
The calculateDiffMeth()
function is the main function to calculate differential methylation. Depending on the sample size per each set it will either use Fisher’s exact or logistic regression to calculate P-values. P-values will be adjusted to Q-values using SLIM method (Wang, Tuominen, and Tsai 2011). If you have replicates, the function will automatically use logistic regression. You can force the calculateDiffMeth()
function to use Fisher’s exact test if you pool the replicates when there is only test and control sample groups. This can be achieved with pool()
function, see FAQ for more info.
In its simplest form ,where there are no covariates, the logistic regression will try to model the the log odds ratio which is based on methylation proportion of a CpG, \(\pi_i\), using the treatment vector which denotes the sample group membership for the CpGs in the model. Below, the “Treatment” variable is used to predict the log-odds ratio of methylation proportions.
\[ \text{log}\left(\dfrac{\pi_i}{1-\pi_i}\right) =\beta_0 + \beta_1 Treatment_i \]
The logistic regression model is fitted per CpG or per region and we test if treatment vector has any effect on the outcome variable or not. In other words, we are testing if \(log(\pi_i/(1-\pi_i)) = \beta_0 + \beta_1 Treatment_i\) is a “better” model than \(log(\pi_i/(1-\pi_i)) = \beta_0\).
The following code snippet tests for differential methylation. Since the example data has replicates, the logistic regression based modeling and test will be used.
myDiff=calculateDiffMeth(meth)
After q-value calculation, we can select the differentially methylated regions/bases based on q-value and percent methylation difference cutoffs. Following bit selects the bases that have q-value<0.01 and percent methylation difference larger than 25%. If you specify type="hyper"
or type="hypo"
options, you will get hyper-methylated or hypo-methylated regions/bases.
# get hyper methylated bases
myDiff25p.hyper=getMethylDiff(myDiff,difference=25,qvalue=0.01,type="hyper")
#
# get hypo methylated bases
myDiff25p.hypo=getMethylDiff(myDiff,difference=25,qvalue=0.01,type="hypo")
#
#
# get all differentially methylated bases
myDiff25p=getMethylDiff(myDiff,difference=25,qvalue=0.01)
We can also visualize the distribution of hypo/hyper-methylated bases/regions per chromosome using the following function. In this case, the example set includes only one chromosome. The list
shows percentages of hypo/hyper methylated bases over all the covered bases in a given chromosome.
diffMethPerChr(myDiff,plot=FALSE,qvalue.cutoff=0.01, meth.cutoff=25)
## $diffMeth.per.chr
## chr number.of.hypermethylated percentage.of.hypermethylated
## 1 chr21 75 7.788162
## number.of.hypomethylated percentage.of.hypomethylated
## 1 59 6.126687
##
## $diffMeth.all
## number.of.hypermethylated percentage.of.hypermethylated
## 1 75 7.788162
## number.of.hypomethylated percentage.of.hypomethylated
## 1 59 6.126687
Overdispersion occurs when there is more variability in the data than assumed by the distribution. In the logistic regression model, the response variable \(meth_i\) (number of methylated CpGs) is expected to have a binomial distribution: \[meth_i \sim Bin(n_i, \pi_i)\] Therefore, the methylated CpGs will have the variance \(n_i \pi_i(1-\pi_i)\) and mean \(\mu_i=n_i \pi_i\). \(n_i\) is the coverage for the CpG or a region and \(\pi_i\) is the underlying methylation proportion.
Overdispersion occurs when the variance of \(meth_i\) is greater than \(n_i \hat{\pi_i}(1-\hat{\pi_i})\), where \(\hat{\pi_i}\) is the estimated methylation proportion from the model. This can be corrected by calculating a scaling parameter \(\phi\) and adjusting the variance as \(\phi n_i \hat{\pi_i}(1-\hat{\pi_i})\). calculateDiffMeth
can calculate that scaling parameter and use it in statistical tests to correct for overdispersion. The parameter is calculated as proposed by (McCullagh and Nelder 1989) as follows: \(\hat{\phi}=X^2/(N-P)\), where \(X\) is Pearson goodness-of-fit statistic, \(N\) is the number of samples, and \(P\) is the number of parameters. This scaling parameter also effects the statistical tests and if there is overdispersion correction the tests will be more stringent in general.
By default,this overdispersion correction is not applied. This can be achieved by setting overdispersion="MN"
. The Chisq-test is used by default only when no overdispersion correction is applied. If overdispersion correction is applied, the function automatically switches to the F-test. The Chisq-test can be manually chosen in this case as well, but the F-test only works with overdispersion correction switched on. In both cases, the procedure tests if the full model (the model where treatment is included as an explanatory variable) explains the data better than the null model (the model with no treatment, just intercept). If there is no effect based on samples being from different groups adding a treatment vector for sample groupings will be no better than not adding the treatment vector. Below, we simulate methylation data and use overdispersion correction for the logistic regression model.
sim.methylBase1<-dataSim(replicates=6,sites=1000,
treatment=c(rep(1,3),rep(0,3)),
sample.ids=c(paste0("test",1:3),paste0("ctrl",1:3))
)
my.diffMeth<-calculateDiffMeth(sim.methylBase1[1:,],
overdispersion="MN",test="Chisq",mc.cores=1)
Covariates can be included in the analysis. The function will then try to separate the influence of the covariates from the treatment effect via the logistic regression model. In this case, we will test if full model (model with treatment and covariates) is better than the model with the covariates only. If there is no effect due to the treatment (sample groups), the full model will not explain the data better than the model with covariates only. In calculateDiffMeth
, this is achieved by supplying the covariates
argument in the format of a data.frame
. Below, we simulate methylation data and add make a data.frame
for the age. The data frame can include more columns, and those columns can also be factor
variables. The row order of the data.frame should match the order of samples in the methylBase
object.
covariates=data.frame(age=c(30,80,34,30,80,40))
sim.methylBase<-dataSim(replicates=6,sites=1000,
treatment=c(rep(1,3),rep(0,3)),
covariates=covariates,
sample.ids=c(paste0("test",1:3),paste0("ctrl",1:3))
)
my.diffMeth3<-calculateDiffMeth(sim.methylBase,
covariates=covariates,
overdispersion="MN",test="Chisq",mc.cores=1)
The differential methylation calculation speed can be increased substantially by utilizing multiple-cores in a machine if available. Both Fisher’s Exact test and logistic regression based test are able to use multiple-core option.
The following piece of code will run differential methylation calculation using 2 cores.
myDiff=calculateDiffMeth(meth,num.cores=2)
We can annotate our differentially methylated regions/bases based on gene annotation using genomation package. In this example, we read the gene annotation from a BED file and annotate our differentially methylated regions with that information using genomation functions. Note that these functions operate on GRanges
objects ,so we first coerce methylKit objects to GRanges. This annotation operation will tell us what percentage of our differentially methylated regions are on promoters/introns/exons/intergenic region. In this case we read annotation from a BED file, similar gene annotation information can be fetched using GenomicFeatures
package or other packages available from Bioconductor.org.
library(genomation)
## Loading required package: grid
##
## Attaching package: 'genomation'
## The following objects are masked from 'package:methylKit':
##
## getFeatsWithTargetsStats, getFlanks, getMembers,
## getTargetAnnotationStats, plotTargetAnnotation
# read the gene BED file
gene.obj=readTranscriptFeatures(system.file("extdata", "refseq.hg18.bed.txt",
package = "methylKit"))
## Reading the table...
## Calculating intron coordinates...
## Calculating exon coordinates...
## Calculating TSS coordinates...
## Calculating promoter coordinates...
## Outputting the final GRangesList...
#
# annotate differentially methylated CpGs with
# promoter/exon/intron using annotation data
#
annotateWithGeneParts(as(myDiff25p,"GRanges"),gene.obj)
## Summary of target set annotation with genic parts
## Rows in target set: 133
## -----------------------
## percentage of target features overlapping with annotation:
## promoter exon intron intergenic
## 27.82 15.04 34.59 57.14
##
## percentage of target features overlapping with annotation:
## (with promoter > exon > intron precedence):
## promoter exon intron intergenic
## 27.82 0.00 15.04 57.14
##
## percentage of annotation boundaries with feature overlap:
## promoter exon intron
## 0.29 0.03 0.17
##
## summary of distances to the nearest TSS:
## Min. 1st Qu. Median Mean 3rd Qu. Max.
## 5 828 45160 52030 94640 313500
##
Similarly, we can read the CpG island annotation and annotate our differentially methylated bases/regions with them.
# read the shores and flanking regions and name the flanks as shores
# and CpG islands as CpGi
cpg.obj=readFeatureFlank(system.file("extdata", "cpgi.hg18.bed.txt",
package = "methylKit"),
feature.flank.name=c("CpGi","shores"))
#
# convert methylDiff object to GRanges and annotate
diffCpGann=annotateWithFeatureFlank(as(myDiff25p,"GRanges"),
cpg.obj$CpGi,cpg.obj$shores,
feature.name="CpGi",flank.name="shores")
We can also summarize methylation information over a set of defined regions such as promoters or CpG islands. The function below summarizes the methylation information over a given set of promoter regions and outputs a methylRaw
or methylRawList
object depending on the input. We are using the output of genomation functions used above to provide the locations of promoters. For regional summary functions, we need to provide regions of interest as GRanges object.
promoters=regionCounts(myobj,gene.obj$promoters)
head(promoters[[1]])
## chr start end strand coverage numCs numTs
## 1 chr21 10011791 10013791 - 7953 6662 1290
## 2 chr21 10119796 10121796 - 1725 1171 554
## 3 chr21 10119808 10121808 - 1725 1171 554
## 4 chr21 13903368 13905368 + 10 10 0
## 5 chr21 14273636 14275636 - 282 220 62
## 6 chr21 14509336 14511336 + 1058 55 1003
After getting the annotation of differentially methylated regions, we can get the distance to TSS and nearest gene name using the getAssociationWithTSS
function from genomation package.
diffAnn=annotateWithGeneParts(as(myDiff25p,"GRanges"),gene.obj)
# target.row is the row number in myDiff25p
head(getAssociationWithTSS(diffAnn))
## target.row dist.to.feature feature.name feature.strand
## 60 1 106111 NM_199260 -
## 60.1 2 106098 NM_199260 -
## 60.2 3 106092 NM_199260 -
## 60.3 4 105919 NM_199260 -
## 60.4 5 85265 NM_199260 -
## 60.5 6 68287 NM_199260 -
It is also desirable to get percentage/number of differentially methylated regions that overlap with intron/exon/promoters
getTargetAnnotationStats(diffAnn,percentage=TRUE,precedence=TRUE)
## promoter exon intron intergenic
## 27.82 0.00 15.04 57.14
We can also plot the percentage of differentially methylated bases overlapping with exon/intron/promoters
plotTargetAnnotation(diffAnn,precedence=TRUE,
main="differential methylation annotation")
We can also plot the CpG island annotation the same way. The plot below shows what percentage of differentially methylated bases are on CpG islands, CpG island shores and other regions.
plotTargetAnnotation(diffCpGann,col=c("green","gray","white"),
main="differential methylation annotation")
It might be also useful to get percentage of intron/exon/promoters that overlap with differentially methylated bases.
getFeatsWithTargetsStats(diffAnn,percentage=TRUE)
## promoter exon intron
## 0.29 0.03 0.17
Most methylKit
objects (methylRaw,methylBase and methylDiff), including methylDB objects (methylRawDB,methylBaseDB and methylDiffDB) can be coerced to GRanges
objects from GenomicRanges
package. Coercing methylKit objects to GRanges
will give users additional flexibility when customizing their analyses.
class(meth)
## [1] "methylBase"
## attr(,"package")
## [1] "methylKit"
as(meth,"GRanges")
## GRanges object with 963 ranges and 12 metadata columns:
## seqnames ranges strand | coverage1 numCs1
## <Rle> <IRanges> <Rle> | <integer> <numeric>
## [1] chr21 [9853296, 9853296] + | 17 10
## [2] chr21 [9853326, 9853326] + | 17 12
## [3] chr21 [9860126, 9860126] + | 39 38
## [4] chr21 [9906604, 9906604] + | 68 42
## [5] chr21 [9906616, 9906616] + | 68 52
## ... ... ... ... . ... ...
## [959] chr21 [19855690, 19855690] + | 27 26
## [960] chr21 [19855706, 19855706] + | 27 27
## [961] chr21 [19855711, 19855711] + | 18 18
## [962] chr21 [19943653, 19943653] + | 12 12
## [963] chr21 [19943695, 19943695] + | 12 11
## numTs1 coverage2 numCs2 numTs2 coverage3 numCs3
## <numeric> <integer> <numeric> <numeric> <integer> <numeric>
## [1] 7 333 268 65 18 16
## [2] 5 329 249 79 16 14
## [3] 1 83 78 5 83 83
## [4] 26 111 97 14 23 18
## [5] 16 111 104 7 23 14
## ... ... ... ... ... ... ...
## [959] 1 19 17 2 34 34
## [960] 0 19 19 0 34 34
## [961] 0 18 15 3 34 34
## [962] 0 32 30 2 26 25
## [963] 1 32 32 0 26 26
## numTs3 coverage4 numCs4 numTs4
## <numeric> <integer> <numeric> <numeric>
## [1] 2 395 341 54
## [2] 2 379 284 95
## [3] 0 41 40 1
## [4] 5 37 33 4
## [5] 9 37 27 10
## ... ... ... ... ...
## [959] 0 12 12 0
## [960] 0 12 11 1
## [961] 0 12 12 0
## [962] 1 24 22 2
## [963] 0 27 24 3
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
class(myDiff)
## [1] "methylDiff"
## attr(,"package")
## [1] "methylKit"
as(myDiff,"GRanges")
## GRanges object with 963 ranges and 2 metadata columns:
## seqnames ranges strand | qvalue
## <Rle> <IRanges> <Rle> | <numeric>
## [1] chr21 [9853296, 9853296] + | 0.0215658126063664
## [2] chr21 [9853326, 9853326] + | 0.592173028310101
## [3] chr21 [9860126, 9860126] + | 0.0697808391745445
## [4] chr21 [9906604, 9906604] + | 0.259453089661393
## [5] chr21 [9906616, 9906616] + | 0.00432220069940914
## ... ... ... ... . ...
## [959] chr21 [19855690, 19855690] + | 0.0660274910679262
## [960] chr21 [19855706, 19855706] + | 0.282958728725395
## [961] chr21 [19855711, 19855711] + | 0.0446545923565475
## [962] chr21 [19943653, 19943653] + | 0.592173028310101
## [963] chr21 [19943695, 19943695] + | 0.396045532748492
## meth.diff
## <numeric>
## [1] -7.01210653753027
## [2] 0.209135938359928
## [3] -4.11158117398203
## [4] -7.34636871508379
## [5] 18.8175046554935
## ... ...
## [959] -6.52173913043478
## [960] 2.17391304347826
## [961] -8.33333333333334
## [962] 1.45454545454546
## [963] 3.38765008576329
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
methylDB objects (methylRawDB
, methylBaseDB
and methylDiffDB
) can be coerced to normal methylKit
objects. This might speed up the analysis if sufficient computing resources are available. This can be done via “as()” function.
class(myobjDB[[1]])
## [1] "methylRawDB"
## attr(,"package")
## [1] "methylKit"
as(myobjDB[[1]],"methylRaw")
## chr start end strand coverage numCs numTs
## 1 chr21 9764513 9764513 - 12 0 12
## 2 chr21 9764539 9764539 - 12 3 9
## 3 chr21 9820622 9820622 + 13 0 13
## 4 chr21 9837545 9837545 + 11 0 11
## 5 chr21 9849022 9849022 + 124 90 34
## 6 chr21 9853296 9853296 + 17 10 7
## 7 chr21 9853326 9853326 + 17 12 5
## 8 chr21 9860126 9860126 + 39 38 1
## 9 chr21 9874553 9874553 + 12 12 0
## 10 chr21 9882409 9882409 + 10 7 3
## 11 chr21 9882434 9882434 + 10 0 10
## 12 chr21 9906604 9906604 + 68 42 26
## 13 chr21 9906616 9906616 + 68 52 16
## 14 chr21 9906619 9906619 + 68 59 9
## 15 chr21 9906634 9906634 + 68 56 12
## 16 chr21 9906635 9906635 - 37 24 13
## 17 chr21 9906640 9906640 + 67 6 61
## 18 chr21 9906641 9906641 - 36 4 32
## 19 chr21 9906644 9906644 + 68 63 5
## 20 chr21 9906645 9906645 - 36 30 6
## 21 chr21 9906655 9906655 - 37 26 11
## 22 chr21 9906657 9906657 - 37 4 33
## 23 chr21 9906660 9906660 - 37 24 13
## 24 chr21 9906663 9906663 - 37 25 12
## 25 chr21 9906676 9906676 + 21 17 4
## 26 chr21 9906677 9906677 - 37 28 9
## 27 chr21 9906681 9906681 + 21 12 9
## 28 chr21 9906694 9906694 + 21 9 12
## 29 chr21 9906700 9906700 + 13 6 7
## 30 chr21 9906714 9906714 + 14 3 11
## 31 chr21 9906846 9906846 + 14 13 1
## 32 chr21 9906864 9906864 + 14 10 4
## 33 chr21 9906869 9906869 + 14 4 10
## 34 chr21 9906873 9906873 + 12 8 4
## 35 chr21 9906884 9906884 + 14 1 13
## 36 chr21 9906889 9906889 + 11 0 11
## 37 chr21 9906933 9906933 + 20 17 3
## 38 chr21 9906949 9906949 + 20 14 6
## 39 chr21 9906952 9906952 + 20 15 5
## 40 chr21 9906954 9906954 + 20 15 5
## 41 chr21 9906957 9906957 + 20 17 3
## 42 chr21 9911433 9911433 - 19 7 12
## 43 chr21 9911447 9911447 - 19 17 2
## 44 chr21 9913543 9913543 + 36 26 10
## 45 chr21 9913575 9913575 + 35 30 5
## 46 chr21 9927527 9927527 + 17 5 12
## 47 chr21 9944505 9944505 + 37 2 35
## 48 chr21 9944663 9944663 - 61 19 42
## 49 chr21 9959407 9959407 + 44 17 27
## 50 chr21 9959541 9959541 - 26 12 14
## 51 chr21 9959569 9959569 - 25 17 8
## 52 chr21 9959577 9959577 - 25 25 0
## 53 chr21 9959644 9959644 - 21 0 21
## 54 chr21 9959650 9959650 - 21 6 15
## 55 chr21 9967634 9967634 - 10 0 10
## 56 chr21 10011833 10011833 + 174 173 1
## 57 chr21 10011841 10011841 + 173 164 9
## 58 chr21 10011855 10011855 + 175 175 0
## 59 chr21 10011858 10011858 + 175 131 44
## 60 chr21 10011861 10011861 + 174 147 27
## 61 chr21 10011872 10011872 + 167 160 7
## 62 chr21 10011876 10011876 + 160 148 12
## 63 chr21 10011878 10011878 + 150 134 16
## 64 chr21 10011925 10011925 - 120 65 55
## 65 chr21 10011938 10011938 - 134 127 7
## 66 chr21 10011954 10011954 - 133 131 2
## 67 chr21 10011964 10011964 - 139 110 29
## 68 chr21 10011974 10011974 - 139 106 33
## 69 chr21 10011987 10011987 + 224 214 10
## 70 chr21 10011990 10011990 + 224 184 40
## 71 chr21 10012008 10012008 + 229 91 138
## 72 chr21 10012030 10012030 + 229 180 49
## 73 chr21 10012069 10012069 + 94 88 6
## 74 chr21 10012075 10012075 + 95 93 2
## 75 chr21 10012079 10012079 + 94 88 6
## 76 chr21 10012089 10012089 + 94 80 14
## 77 chr21 10012095 10012095 + 94 89 5
## 78 chr21 10012101 10012101 + 98 96 2
## 79 chr21 10012316 10012316 - 14 11 3
## 80 chr21 10012318 10012318 - 14 14 0
## 81 chr21 10012321 10012321 - 14 14 0
## 82 chr21 10012324 10012324 - 14 14 0
## 83 chr21 10012328 10012328 - 14 13 1
## 84 chr21 10012425 10012425 - 70 48 22
## 85 chr21 10012428 10012428 - 69 52 17
## 86 chr21 10012430 10012430 - 70 54 16
## 87 chr21 10012442 10012442 - 67 50 17
## 88 chr21 10012454 10012454 - 70 65 5
## 89 chr21 10012465 10012465 - 70 56 14
## 90 chr21 10012587 10012587 - 39 8 31
## 91 chr21 10012606 10012606 - 43 4 39
## 92 chr21 10012609 10012609 + 162 148 14
## 93 chr21 10012610 10012610 - 44 41 3
## 94 chr21 10012636 10012636 + 158 115 43
## 95 chr21 10012637 10012637 - 136 88 48
## 96 chr21 10012640 10012640 + 156 129 27
## 97 chr21 10012641 10012641 - 136 111 25
## 98 chr21 10012642 10012642 + 155 146 9
## 99 chr21 10012643 10012643 - 136 130 6
## 100 chr21 10012658 10012658 + 148 125 22
## 101 chr21 10012659 10012659 - 136 109 27
## 102 chr21 10012662 10012662 - 136 133 3
## 103 chr21 10012696 10012696 + 148 147 1
## 104 chr21 10012699 10012699 + 148 140 8
## 105 chr21 10012862 10012862 - 10 9 1
## 106 chr21 10012876 10012876 + 29 28 1
## 107 chr21 10012877 10012877 - 10 9 1
## 108 chr21 10012881 10012881 + 30 0 30
## 109 chr21 10012883 10012883 + 29 29 0
## 110 chr21 10012887 10012887 + 30 27 3
## 111 chr21 10012891 10012891 + 30 29 1
## 112 chr21 10012894 10012894 + 30 23 7
## 113 chr21 10012898 10012898 + 20 20 0
## 114 chr21 10012900 10012900 + 24 22 2
## 115 chr21 10012902 10012902 + 24 24 0
## 116 chr21 10012921 10012921 + 23 23 0
## 117 chr21 10012924 10012924 + 23 23 0
## 118 chr21 10012934 10012934 - 84 82 2
## 119 chr21 10012938 10012938 - 119 116 3
## 120 chr21 10012940 10012940 - 118 118 0
## 121 chr21 10012942 10012942 - 128 126 2
## 122 chr21 10012944 10012944 - 130 123 7
## 123 chr21 10012947 10012947 - 137 134 3
## 124 chr21 10012961 10012961 - 143 136 7
## 125 chr21 10012975 10012975 - 149 107 42
## 126 chr21 10012981 10012981 + 18 13 5
## 127 chr21 10012982 10012982 - 152 150 2
## 128 chr21 10012989 10012989 + 18 15 3
## 129 chr21 10012994 10012994 + 18 13 5
## 130 chr21 10013000 10013000 + 12 6 6
## 131 chr21 10013039 10013039 - 106 0 106
## 132 chr21 10013060 10013060 - 109 56 53
## 133 chr21 10013062 10013062 - 111 55 56
## 134 chr21 10013072 10013072 - 111 97 14
## 135 chr21 10013086 10013086 - 112 106 6
## 136 chr21 10013594 10013594 - 11 4 7
## 137 chr21 10013827 10013827 - 16 14 2
## 138 chr21 10022131 10022131 - 17 7 10
## 139 chr21 10022157 10022157 - 17 8 9
## 140 chr21 10022165 10022165 - 17 16 1
## 141 chr21 10065022 10065022 + 26 26 0
## 142 chr21 10065061 10065061 + 25 12 13
## 143 chr21 10077699 10077699 + 47 10 37
## 144 chr21 10077705 10077705 + 47 12 35
## 145 chr21 10077904 10077904 - 17 10 7
## 146 chr21 10081627 10081627 + 12 2 10
## 147 chr21 10081638 10081638 + 14 6 8
## 148 chr21 10081652 10081652 + 11 3 8
## 149 chr21 10086785 10086785 + 20 8 12
## 150 chr21 10087014 10087014 - 10 3 7
## 151 chr21 10099942 10099942 + 13 13 0
## 152 chr21 10099950 10099950 + 13 10 3
## 153 chr21 10114305 10114305 - 12 2 10
## 154 chr21 10114343 10114343 - 13 8 5
## 155 chr21 10114345 10114345 - 13 5 8
## 156 chr21 10114351 10114351 - 13 6 7
## 157 chr21 10118483 10118483 + 12 0 12
## 158 chr21 10120564 10120564 - 25 15 10
## 159 chr21 10120584 10120584 - 27 8 19
## 160 chr21 10120599 10120599 + 84 78 6
## 161 chr21 10120600 10120600 - 27 21 6
## 162 chr21 10120603 10120603 + 94 12 82
## 163 chr21 10120634 10120634 + 99 14 85
## 164 chr21 10120643 10120643 + 97 18 79
## 165 chr21 10121028 10121028 + 11 11 0
## 166 chr21 10121033 10121033 + 11 9 2
## 167 chr21 10121126 10121126 + 818 686 132
## 168 chr21 10121171 10121171 + 406 284 122
## 169 chr21 10121198 10121198 + 26 15 11
## 170 chr21 10122567 10122567 + 14 14 0
## 171 chr21 10128590 10128590 + 310 275 35
## 172 chr21 10128601 10128601 + 306 270 36
## 173 chr21 10128603 10128603 + 308 284 24
## 174 chr21 10128609 10128609 + 309 292 17
## 175 chr21 10128617 10128617 + 307 255 51
## 176 chr21 10128639 10128639 + 527 227 300
## 177 chr21 10128659 10128659 + 237 159 78
## 178 chr21 10128660 10128660 - 135 90 45
## 179 chr21 10128661 10128661 + 236 196 40
## 180 chr21 10128662 10128662 - 136 115 17
## 181 chr21 10128670 10128670 + 237 225 12
## 182 chr21 10128671 10128671 - 138 117 21
## 183 chr21 10128673 10128673 + 227 194 33
## 184 chr21 10128674 10128674 - 139 126 13
## 185 chr21 10128682 10128682 + 145 2 142
## 186 chr21 10128683 10128683 - 101 39 62
## 187 chr21 10128702 10128702 + 87 87 0
## 188 chr21 10128705 10128705 + 87 86 1
## 189 chr21 10128706 10128706 - 49 47 2
## 190 chr21 10128707 10128707 + 87 84 3
## 191 chr21 10128708 10128708 - 50 47 3
## 192 chr21 10128710 10128710 + 87 85 2
## 193 chr21 10128711 10128711 - 50 48 2
## 194 chr21 10128720 10128720 + 87 83 4
## 195 chr21 10128721 10128721 - 50 48 2
## 196 chr21 10128734 10128734 + 260 201 59
## 197 chr21 10128735 10128735 - 50 48 2
## 198 chr21 10128738 10128738 + 286 223 63
## 199 chr21 10128745 10128745 + 292 168 124
## 200 chr21 10128748 10128748 + 294 258 36
## 201 chr21 10128758 10128758 + 192 108 84
## 202 chr21 10128775 10128775 + 204 161 43
## 203 chr21 10128839 10128839 + 518 472 46
## 204 chr21 10128850 10128850 + 514 395 119
## 205 chr21 10128867 10128867 + 522 487 35
## 206 chr21 10128881 10128881 + 515 463 52
## 207 chr21 10128883 10128883 + 511 465 46
## 208 chr21 10128964 10128964 - 68 40 28
## 209 chr21 10128980 10128980 - 70 38 32
## 210 chr21 10128989 10128989 - 70 62 8
## 211 chr21 10128991 10128991 - 70 39 31
## 212 chr21 10129002 10129002 + 78 77 1
## 213 chr21 10129003 10129003 - 71 60 11
## 214 chr21 10129032 10129032 + 80 61 19
## 215 chr21 10129036 10129036 + 81 79 2
## 216 chr21 10129044 10129044 + 80 78 2
## 217 chr21 10129181 10129181 + 34 28 6
## 218 chr21 10129546 10129546 + 36 35 1
## 219 chr21 10129550 10129550 + 36 26 10
## 220 chr21 10131969 10131969 + 14 14 0
## 221 chr21 10132033 10132033 + 12 7 5
## 222 chr21 10132081 10132081 - 29 27 2
## 223 chr21 10132089 10132089 - 29 25 4
## 224 chr21 10132094 10132094 + 35 28 7
## 225 chr21 10132095 10132095 - 29 15 14
## 226 chr21 10132122 10132122 + 36 24 12
## 227 chr21 10132130 10132130 + 136 92 44
## 228 chr21 10132131 10132131 - 16 12 4
## 229 chr21 10132177 10132177 + 108 41 67
## 230 chr21 10132398 10132398 + 19 14 5
## 231 chr21 10132406 10132406 + 19 14 5
## 232 chr21 10132418 10132418 + 19 8 11
## 233 chr21 10133048 10133048 + 133 105 28
## 234 chr21 10133056 10133056 + 132 94 38
## 235 chr21 10133080 10133080 + 131 73 58
## 236 chr21 10133091 10133091 + 124 41 83
## 237 chr21 10133092 10133092 - 16 1 15
## 238 chr21 10133111 10133111 - 74 47 27
## 239 chr21 10133159 10133159 + 50 24 26
## 240 chr21 10133160 10133160 - 114 60 54
## 241 chr21 10133401 10133401 + 17 9 8
## 242 chr21 10133402 10133402 - 26 12 14
## 243 chr21 10136614 10136614 + 97 79 18
## 244 chr21 10136632 10136632 + 133 9 124
## 245 chr21 10136648 10136648 + 138 29 109
## 246 chr21 10136731 10136731 - 145 97 48
## 247 chr21 10136737 10136737 - 144 114 30
## 248 chr21 10136755 10136755 - 143 113 30
## 249 chr21 10136772 10136772 + 20 11 9
## 250 chr21 10136773 10136773 - 142 76 66
## 251 chr21 10136966 10136966 - 11 0 11
## 252 chr21 10145439 10145439 + 32 16 16
## 253 chr21 10145470 10145470 + 32 24 8
## 254 chr21 10145471 10145471 - 39 17 22
## 255 chr21 10145507 10145507 - 41 24 17
## 256 chr21 10148190 10148190 + 16 11 5
## 257 chr21 10148224 10148224 + 17 14 3
## 258 chr21 10148396 10148396 - 16 4 12
## 259 chr21 10151952 10151952 - 11 7 4
## 260 chr21 10151967 10151967 + 12 5 7
## 261 chr21 10151968 10151968 - 11 9 2
## 262 chr21 10152050 10152050 - 17 14 3
## 263 chr21 10154024 10154024 - 19 18 1
## 264 chr21 10154062 10154062 - 19 9 10
## 265 chr21 10154072 10154072 - 19 19 0
## 266 chr21 10154210 10154210 - 57 49 8
## 267 chr21 10154226 10154226 - 57 43 14
## 268 chr21 10154246 10154246 - 57 48 9
## 269 chr21 10154306 10154306 - 98 82 14
## 270 chr21 10154872 10154872 + 49 24 25
## 271 chr21 10154873 10154873 - 28 22 6
## 272 chr21 10154893 10154893 + 49 32 17
## 273 chr21 10154894 10154894 - 28 19 9
## 274 chr21 10154897 10154897 + 43 0 43
## 275 chr21 10154898 10154898 - 27 0 27
## 276 chr21 10154916 10154916 - 29 24 5
## 277 chr21 10156390 10156390 + 13 7 6
## 278 chr21 10156406 10156406 + 13 10 3
## 279 chr21 10156434 10156434 + 11 11 0
## 280 chr21 10156485 10156485 - 13 10 3
## 281 chr21 10156496 10156496 - 14 13 1
## 282 chr21 10156520 10156520 + 148 77 71
## 283 chr21 10156521 10156521 - 14 7 7
## 284 chr21 10157440 10157440 - 68 42 26
## 285 chr21 10158020 10158020 + 10 8 2
## 286 chr21 10158023 10158023 + 10 8 2
## 287 chr21 10160749 10160749 - 11 9 2
## 288 chr21 10160753 10160753 - 11 7 4
## 289 chr21 10162951 10162951 + 12 0 12
## 290 chr21 10166512 10166512 - 12 12 0
## 291 chr21 10166531 10166531 - 12 0 12
## 292 chr21 10176577 10176577 - 68 37 31
## 293 chr21 10176592 10176592 - 68 52 16
## 294 chr21 10186386 10186386 - 52 39 13
## 295 chr21 10186452 10186452 + 13 10 3
## 296 chr21 10186507 10186507 - 54 42 12
## 297 chr21 10186519 10186519 - 76 24 52
## 298 chr21 10186541 10186541 - 70 5 65
## 299 chr21 10187865 10187865 - 13 10 3
## 300 chr21 10190503 10190503 - 293 224 69
## 301 chr21 10190515 10190515 + 47 43 4
## 302 chr21 10190516 10190516 - 297 249 48
## 303 chr21 10190523 10190523 + 47 37 10
## 304 chr21 10190524 10190524 - 22 19 3
## 305 chr21 10190541 10190541 + 50 49 1
## 306 chr21 10190542 10190542 - 25 22 3
## 307 chr21 10190546 10190546 + 47 34 13
## 308 chr21 10190547 10190547 - 25 21 4
## 309 chr21 10190550 10190550 + 35 33 2
## 310 chr21 10190551 10190551 - 24 21 3
## 311 chr21 10190555 10190555 + 48 42 6
## 312 chr21 10190556 10190556 - 25 22 3
## 313 chr21 10190563 10190563 + 47 38 9
## 314 chr21 10190574 10190574 + 45 34 11
## 315 chr21 10190582 10190582 + 44 0 44
## 316 chr21 10190583 10190583 - 113 92 21
## 317 chr21 10190586 10190586 + 24 19 5
## 318 chr21 10190587 10190587 - 116 111 5
## 319 chr21 10190611 10190611 - 118 0 118
## 320 chr21 10190622 10190622 - 118 97 21
## 321 chr21 10190627 10190627 - 118 38 80
## 322 chr21 10190629 10190629 + 91 82 9
## 323 chr21 10190630 10190630 - 117 101 16
## 324 chr21 10190633 10190633 + 91 67 24
## 325 chr21 10190658 10190658 + 91 83 8
## 326 chr21 10190669 10190669 + 87 83 4
## 327 chr21 10190676 10190676 + 81 70 11
## 328 chr21 10190754 10190754 + 11 9 2
## 329 chr21 10190787 10190787 + 10 10 0
## 330 chr21 10190834 10190834 - 17 11 6
## 331 chr21 10190852 10190852 - 17 7 10
## 332 chr21 10190868 10190868 - 17 5 12
## 333 chr21 10190881 10190881 + 296 152 144
## 334 chr21 10190882 10190882 - 16 7 9
## 335 chr21 10190896 10190896 + 293 239 54
## 336 chr21 10190897 10190897 - 344 266 78
## 337 chr21 10190924 10190924 + 204 195 9
## 338 chr21 10190925 10190925 - 387 373 14
## 339 chr21 10190929 10190929 + 175 129 46
## 340 chr21 10190930 10190930 - 390 248 142
## 341 chr21 10190937 10190937 - 381 351 30
## 342 chr21 10190944 10190944 - 390 366 24
## 343 chr21 10190970 10190970 + 17 17 0
## 344 chr21 10190983 10190983 + 18 18 0
## 345 chr21 10190985 10190985 + 18 17 1
## 346 chr21 10190995 10190995 + 17 13 4
## 347 chr21 10191023 10191023 - 153 88 65
## 348 chr21 10191062 10191062 - 154 110 44
## 349 chr21 10191089 10191089 - 48 41 7
## 350 chr21 10191109 10191109 - 48 42 6
## 351 chr21 10191112 10191112 - 47 33 14
## 352 chr21 10191118 10191118 - 48 46 2
## 353 chr21 10191326 10191326 - 69 40 29
## 354 chr21 10191348 10191348 - 64 51 13
## 355 chr21 10191408 10191408 - 238 54 184
## 356 chr21 10191430 10191430 - 237 54 183
## 357 chr21 10191769 10191769 + 46 21 25
## 358 chr21 10191770 10191770 - 190 152 38
## 359 chr21 10191773 10191773 + 45 28 17
## 360 chr21 10191786 10191786 + 45 25 20
## 361 chr21 10191818 10191818 + 45 36 9
## 362 chr21 10191819 10191819 - 304 160 144
## 363 chr21 10191824 10191824 - 307 181 126
## 364 chr21 10191830 10191830 - 310 259 51
## 365 chr21 10191838 10191838 - 309 235 74
## 366 chr21 10191851 10191851 - 310 218 92
## 367 chr21 10191878 10191878 + 309 289 20
## 368 chr21 10191899 10191899 + 305 250 55
## 369 chr21 10191900 10191900 - 299 242 57
## 370 chr21 10191906 10191906 + 308 191 117
## 371 chr21 10191907 10191907 - 311 216 95
## 372 chr21 10191923 10191923 + 262 200 61
## 373 chr21 10191924 10191924 - 313 223 90
## 374 chr21 10191940 10191940 - 328 237 90
## 375 chr21 10192279 10192279 + 245 138 107
## 376 chr21 10192318 10192318 + 230 122 107
## 377 chr21 10192319 10192319 - 178 141 37
## 378 chr21 10192322 10192322 + 241 190 51
## 379 chr21 10192323 10192323 - 249 203 46
## 380 chr21 10192331 10192331 - 248 163 85
## 381 chr21 10192334 10192334 + 11 7 4
## 382 chr21 10192335 10192335 - 254 188 66
## 383 chr21 10192368 10192368 + 10 7 3
## 384 chr21 10192378 10192378 + 11 9 2
## 385 chr21 13282433 13282433 + 63 41 22
## 386 chr21 13282440 13282440 + 62 20 42
## 387 chr21 13294742 13294742 + 20 20 0
## 388 chr21 13294755 13294755 + 21 21 0
## 389 chr21 13294757 13294757 + 21 21 0
## 390 chr21 13294785 13294785 + 22 19 3
## 391 chr21 13306319 13306319 - 30 21 9
## 392 chr21 13306339 13306339 - 30 30 0
## 393 chr21 13306349 13306349 - 30 22 8
## 394 chr21 13490512 13490512 + 23 17 6
## 395 chr21 13490521 13490521 + 23 11 12
## 396 chr21 13490535 13490535 + 22 19 2
## 397 chr21 13490540 13490540 + 23 21 2
## 398 chr21 13490561 13490561 + 10 8 2
## 399 chr21 13662766 13662766 + 11 11 0
## 400 chr21 13664888 13664888 + 20 10 10
## 401 chr21 13664889 13664889 - 71 38 33
## 402 chr21 13666397 13666397 - 61 45 16
## 403 chr21 13666553 13666553 - 10 7 3
## 404 chr21 13673290 13673290 + 17 17 0
## 405 chr21 13673294 13673294 + 17 17 0
## 406 chr21 13673713 13673713 + 10 10 0
## 407 chr21 13673740 13673740 + 10 9 1
## 408 chr21 13673754 13673754 + 10 9 1
## 409 chr21 13673865 13673865 + 13 13 0
## 410 chr21 13719180 13719180 - 11 5 6
## 411 chr21 13719188 13719188 - 11 0 11
## 412 chr21 13786212 13786212 + 20 10 10
## 413 chr21 13808881 13808881 - 11 8 3
## 414 chr21 13808911 13808911 - 12 10 2
## 415 chr21 13866816 13866816 + 11 3 8
## 416 chr21 13866838 13866838 + 11 7 4
## 417 chr21 13877808 13877808 - 10 6 4
## 418 chr21 13878804 13878804 - 11 9 2
## 419 chr21 13903783 13903783 - 10 10 0
## 420 chr21 13958077 13958077 - 18 12 6
## 421 chr21 13958080 13958080 - 17 15 2
## 422 chr21 13958094 13958094 - 17 14 3
## 423 chr21 13958107 13958107 - 18 13 5
## 424 chr21 13958266 13958266 - 20 15 5
## 425 chr21 13958466 13958466 + 11 5 6
## 426 chr21 13960041 13960041 + 16 10 6
## 427 chr21 13960067 13960067 + 16 14 2
## 428 chr21 13960069 13960069 + 15 14 1
## 429 chr21 13960080 13960080 + 16 14 2
## 430 chr21 13960088 13960088 + 16 14 2
## 431 chr21 13989277 13989277 + 48 47 1
## 432 chr21 13989282 13989282 + 48 47 1
## 433 chr21 13989284 13989284 + 48 44 4
## 434 chr21 13989303 13989303 + 48 47 1
## 435 chr21 13989320 13989320 + 47 45 2
## 436 chr21 13989384 13989384 - 36 24 12
## 437 chr21 13989387 13989387 - 34 23 11
## 438 chr21 13989431 13989431 + 16 15 1
## 439 chr21 13989432 13989432 - 39 38 1
## 440 chr21 13989448 13989448 + 16 14 2
## 441 chr21 13989449 13989449 - 29 25 4
## 442 chr21 13989451 13989451 + 16 15 1
## 443 chr21 13989452 13989452 - 30 28 2
## 444 chr21 13989454 13989454 + 16 12 4
## 445 chr21 13989455 13989455 - 30 28 2
## 446 chr21 13989497 13989497 - 30 27 3
## 447 chr21 13990130 13990130 + 30 25 5
## 448 chr21 13990159 13990159 + 30 29 1
## 449 chr21 13990171 13990171 + 29 28 1
## 450 chr21 13990874 13990874 - 14 13 1
## 451 chr21 13990884 13990884 - 12 12 0
## 452 chr21 13990904 13990904 - 12 10 2
## 453 chr21 13990906 13990906 - 12 12 0
## 454 chr21 13990915 13990915 + 14 3 11
## 455 chr21 13990916 13990916 - 12 7 5
## 456 chr21 13990925 13990925 + 14 14 0
## 457 chr21 13991040 13991040 - 17 9 8
## 458 chr21 13991046 13991046 - 17 13 4
## 459 chr21 13991062 13991062 - 17 15 2
## 460 chr21 13991068 13991068 + 10 10 0
## 461 chr21 13991069 13991069 - 17 16 1
## 462 chr21 13991080 13991080 + 10 10 0
## 463 chr21 13991085 13991085 + 10 9 1
## 464 chr21 13991220 13991220 - 19 6 13
## 465 chr21 13991229 13991229 - 19 6 13
## 466 chr21 13991248 13991248 - 19 19 0
## 467 chr21 13998967 13998967 + 77 38 39
## 468 chr21 13998987 13998987 + 79 68 11
## 469 chr21 13999005 13999005 + 76 71 5
## 470 chr21 13999011 13999011 + 77 69 8
## 471 chr21 13999012 13999012 - 99 92 7
## 472 chr21 13999023 13999023 - 104 90 14
## 473 chr21 13999031 13999031 - 102 49 53
## 474 chr21 13999039 13999039 - 105 37 68
## 475 chr21 13999041 13999041 - 107 33 74
## 476 chr21 13999059 13999059 - 106 90 16
## 477 chr21 13999082 13999082 + 98 93 5
## 478 chr21 13999089 13999089 + 98 88 10
## 479 chr21 13999095 13999095 + 98 91 7
## 480 chr21 13999103 13999103 + 98 76 22
## 481 chr21 13999107 13999107 + 98 84 14
## 482 chr21 13999114 13999114 + 89 69 20
## 483 chr21 13999128 13999128 + 19 18 1
## 484 chr21 13999129 13999129 - 83 33 50
## 485 chr21 13999131 13999131 + 23 9 14
## 486 chr21 13999132 13999132 - 86 17 69
## 487 chr21 13999142 13999142 - 82 43 39
## 488 chr21 13999150 13999150 - 83 39 44
## 489 chr21 13999168 13999168 - 84 62 22
## 490 chr21 13999192 13999192 + 100 66 34
## 491 chr21 13999199 13999199 + 101 59 42
## 492 chr21 13999206 13999206 + 100 82 18
## 493 chr21 13999210 13999210 + 101 80 21
## 494 chr21 13999212 13999212 + 100 55 45
## 495 chr21 13999236 13999236 + 65 57 8
## 496 chr21 13999237 13999237 - 87 72 15
## 497 chr21 13999238 13999238 + 52 39 13
## 498 chr21 13999239 13999239 - 91 62 29
## 499 chr21 13999247 13999247 - 94 64 30
## 500 chr21 13999249 13999249 - 94 77 17
## 501 chr21 13999252 13999252 - 96 94 2
## 502 chr21 13999256 13999256 - 91 85 6
## 503 chr21 13999262 13999262 + 60 55 5
## 504 chr21 13999263 13999263 - 96 94 2
## 505 chr21 13999275 13999275 + 58 50 8
## 506 chr21 13999276 13999276 - 66 61 5
## 507 chr21 13999277 13999277 + 60 51 9
## 508 chr21 13999278 13999278 - 67 63 4
## 509 chr21 13999287 13999287 + 60 49 11
## 510 chr21 13999288 13999288 - 66 58 8
## 511 chr21 13999294 13999294 + 58 52 6
## 512 chr21 13999295 13999295 - 69 61 8
## 513 chr21 13999299 13999299 + 54 25 29
## 514 chr21 13999300 13999300 - 69 38 31
## 515 chr21 13999304 13999304 + 40 38 2
## 516 chr21 13999305 13999305 - 69 60 9
## 517 chr21 13999309 13999309 + 41 37 4
## 518 chr21 13999310 13999310 - 69 60 9
## 519 chr21 13999313 13999313 + 59 57 2
## 520 chr21 13999314 13999314 - 104 67 37
## 521 chr21 13999343 13999343 + 50 32 18
## 522 chr21 13999344 13999344 - 35 17 18
## 523 chr21 13999347 13999347 + 32 29 3
## 524 chr21 13999348 13999348 - 35 18 17
## 525 chr21 13999362 13999362 + 45 6 39
## 526 chr21 13999363 13999363 - 35 33 2
## 527 chr21 13999446 13999446 + 70 60 10
## 528 chr21 13999449 13999449 + 68 42 26
## 529 chr21 13999451 13999451 + 69 59 10
## 530 chr21 13999456 13999456 + 69 58 11
## 531 chr21 13999457 13999457 - 12 12 0
## 532 chr21 13999458 13999458 + 70 51 19
## 533 chr21 13999459 13999459 - 13 12 1
## 534 chr21 13999467 13999467 + 70 37 33
## 535 chr21 13999468 13999468 - 13 11 2
## 536 chr21 13999475 13999475 + 69 56 13
## 537 chr21 13999476 13999476 - 19 16 3
## 538 chr21 13999477 13999477 + 67 47 20
## 539 chr21 13999478 13999478 - 19 19 0
## 540 chr21 13999493 13999493 + 99 45 53
## 541 chr21 13999494 13999494 - 20 18 2
## 542 chr21 13999509 13999509 + 51 42 9
## 543 chr21 13999517 13999517 + 50 37 13
## 544 chr21 13999528 13999528 + 43 42 1
## 545 chr21 13999533 13999533 + 39 36 3
## 546 chr21 13999541 13999541 + 118 89 29
## 547 chr21 13999547 13999547 + 94 84 10
## 548 chr21 13999551 13999551 + 93 67 26
## 549 chr21 13999561 13999561 + 87 62 25
## 550 chr21 13999562 13999562 - 100 67 33
## 551 chr21 13999566 13999566 + 93 89 4
## 552 chr21 13999567 13999567 - 103 101 2
## 553 chr21 13999568 13999568 + 92 82 10
## 554 chr21 13999569 13999569 - 104 102 2
## 555 chr21 13999574 13999574 + 88 84 4
## 556 chr21 13999575 13999575 - 104 88 16
## 557 chr21 13999588 13999588 - 105 100 5
## 558 chr21 13999610 13999610 - 105 92 13
## 559 chr21 13999621 13999621 + 11 11 0
## 560 chr21 13999645 13999645 + 11 8 3
## 561 chr21 13999647 13999647 + 11 11 0
## 562 chr21 13999937 13999937 + 76 71 5
## 563 chr21 13999971 13999971 + 77 68 9
## 564 chr21 14000056 14000056 + 42 26 16
## 565 chr21 14000092 14000092 + 42 19 23
## 566 chr21 14000617 14000617 + 27 16 11
## 567 chr21 14017662 14017662 + 21 13 8
## 568 chr21 14017671 14017671 + 21 21 0
## 569 chr21 14017672 14017672 - 15 14 1
## 570 chr21 14017675 14017675 + 21 21 0
## 571 chr21 14017676 14017676 - 15 14 1
## 572 chr21 14017679 14017679 + 21 20 1
## 573 chr21 14017680 14017680 - 15 14 1
## 574 chr21 14017683 14017683 + 21 21 0
## 575 chr21 14017684 14017684 - 15 14 1
## 576 chr21 14017687 14017687 + 21 17 4
## 577 chr21 14017688 14017688 - 15 8 7
## 578 chr21 14017689 14017689 + 21 21 0
## 579 chr21 14017690 14017690 - 15 14 1
## 580 chr21 14017695 14017695 + 21 19 2
## 581 chr21 14017696 14017696 - 15 12 3
## 582 chr21 14017699 14017699 + 21 21 0
## 583 chr21 14017700 14017700 - 15 11 4
## 584 chr21 14017702 14017702 + 21 21 0
## 585 chr21 14017703 14017703 - 15 12 3
## 586 chr21 14017707 14017707 - 15 9 6
## 587 chr21 14017772 14017772 - 52 36 16
## 588 chr21 14017775 14017775 - 51 26 25
## 589 chr21 14017778 14017778 - 52 26 26
## 590 chr21 14017781 14017781 - 51 29 22
## 591 chr21 14017784 14017784 - 52 21 31
## 592 chr21 14017796 14017796 - 52 29 23
## 593 chr21 14017929 14017929 - 65 10 55
## 594 chr21 14017936 14017936 - 64 27 37
## 595 chr21 14017945 14017945 - 65 36 29
## 596 chr21 14018027 14018027 + 82 76 6
## 597 chr21 14018046 14018046 + 79 28 50
## 598 chr21 14018050 14018050 + 81 70 11
## 599 chr21 14018053 14018053 + 80 42 38
## 600 chr21 14018056 14018056 + 79 25 54
## 601 chr21 14018069 14018069 + 12 6 6
## 602 chr21 14018224 14018224 + 90 86 4
## 603 chr21 14018261 14018261 + 88 61 27
## 604 chr21 14018268 14018268 + 84 51 33
## 605 chr21 14018269 14018269 - 15 14 1
## 606 chr21 14018279 14018279 - 62 29 33
## 607 chr21 14018298 14018298 - 62 55 7
## 608 chr21 14018304 14018304 - 63 61 2
## 609 chr21 14020063 14020063 + 18 12 6
## 610 chr21 14020268 14020268 + 10 7 3
## 611 chr21 14020269 14020269 - 12 11 1
## 612 chr21 14056557 14056557 - 49 48 1
## 613 chr21 14056561 14056561 - 52 49 3
## 614 chr21 14056566 14056566 - 51 48 3
## 615 chr21 14056604 14056604 - 51 47 4
## 616 chr21 14057234 14057234 + 62 31 31
## 617 chr21 14057394 14057394 + 35 28 7
## 618 chr21 14057829 14057829 - 34 29 5
## 619 chr21 14057832 14057832 - 35 27 8
## 620 chr21 14057847 14057847 - 35 29 6
## 621 chr21 14057861 14057861 - 33 24 9
## 622 chr21 14057892 14057892 - 39 30 9
## 623 chr21 14057932 14057932 - 39 30 9
## 624 chr21 14071266 14071266 + 10 9 1
## 625 chr21 14071286 14071286 + 10 5 5
## 626 chr21 14087211 14087211 - 11 11 0
## 627 chr21 14087227 14087227 - 11 11 0
## 628 chr21 14120949 14120949 + 20 18 2
## 629 chr21 14120991 14120991 + 20 12 8
## 630 chr21 14121557 14121557 - 13 13 0
## 631 chr21 14121571 14121571 - 14 14 0
## 632 chr21 14121733 14121733 - 13 11 2
## 633 chr21 14121754 14121754 - 13 13 0
## 634 chr21 14121762 14121762 - 13 9 4
## 635 chr21 14229161 14229161 + 13 13 0
## 636 chr21 14229171 14229171 + 13 13 0
## 637 chr21 14229194 14229194 + 13 13 0
## 638 chr21 14229195 14229195 - 13 11 2
## 639 chr21 14229197 14229197 + 12 12 0
## 640 chr21 14229198 14229198 - 13 10 3
## 641 chr21 14229209 14229209 - 13 12 1
## 642 chr21 14229227 14229227 - 13 10 3
## 643 chr21 14229238 14229238 - 13 10 3
## 644 chr21 14229877 14229877 + 13 10 3
## 645 chr21 14229946 14229946 - 17 16 1
## 646 chr21 14229948 14229948 - 17 17 0
## 647 chr21 14229992 14229992 - 15 14 1
## 648 chr21 14230065 14230065 + 10 5 5
## 649 chr21 14230136 14230136 - 15 11 4
## 650 chr21 14230140 14230140 - 15 11 4
## 651 chr21 14230187 14230187 - 17 10 7
## 652 chr21 14230216 14230216 - 17 16 1
## 653 chr21 14230284 14230284 - 16 3 13
## 654 chr21 14230293 14230293 - 17 15 2
## 655 chr21 14230660 14230660 - 14 11 3
## 656 chr21 14230665 14230665 - 14 14 0
## 657 chr21 14230669 14230669 - 14 13 1
## 658 chr21 14230706 14230706 - 14 5 9
## 659 chr21 14274239 14274239 - 18 17 1
## 660 chr21 14274242 14274242 - 18 2 16
## 661 chr21 14274273 14274273 - 18 0 18
## 662 chr21 14274275 14274275 - 18 16 2
## 663 chr21 14274285 14274285 - 18 16 2
## 664 chr21 14275419 14275419 + 21 21 0
## 665 chr21 14275446 14275446 + 21 21 0
## 666 chr21 14275564 14275564 - 75 52 23
## 667 chr21 14275582 14275582 - 75 75 0
## 668 chr21 14280366 14280366 - 24 18 6
## 669 chr21 14281495 14281495 - 13 12 1
## 670 chr21 14302317 14302317 - 14 7 7
## 671 chr21 14305155 14305155 + 13 8 5
## 672 chr21 14305244 14305244 + 18 12 6
## 673 chr21 14306496 14306496 + 14 11 3
## 674 chr21 14306497 14306497 - 10 9 1
## 675 chr21 14306746 14306746 - 11 11 0
## 676 chr21 14306754 14306754 - 11 9 2
## 677 chr21 14307314 14307314 + 10 2 8
## 678 chr21 14321687 14321687 - 13 0 13
## 679 chr21 14321842 14321842 + 17 0 17
## 680 chr21 14321855 14321855 + 17 0 17
## 681 chr21 14321861 14321861 + 17 0 17
## 682 chr21 14321875 14321875 + 17 0 17
## 683 chr21 14321880 14321880 + 15 0 15
## 684 chr21 14348810 14348810 + 18 18 0
## 685 chr21 14348845 14348845 + 18 16 2
## 686 chr21 14348846 14348846 - 14 14 0
## 687 chr21 14348848 14348848 + 18 18 0
## 688 chr21 14348849 14348849 - 15 15 0
## 689 chr21 14348893 14348893 - 15 15 0
## 690 chr21 14350317 14350317 + 15 15 0
## 691 chr21 14350562 14350562 - 12 12 0
## 692 chr21 14350764 14350764 - 18 18 0
## 693 chr21 14350806 14350806 - 18 18 0
## 694 chr21 14352384 14352384 + 46 44 2
## 695 chr21 14356768 14356768 + 11 11 0
## 696 chr21 14356777 14356777 + 11 11 0
## 697 chr21 14356778 14356778 - 16 16 0
## 698 chr21 14356802 14356802 - 15 15 0
## 699 chr21 14356810 14356810 - 16 16 0
## 700 chr21 14357926 14357926 + 125 121 4
## 701 chr21 14357942 14357942 + 124 122 2
## 702 chr21 14357945 14357945 + 121 110 11
## 703 chr21 14357956 14357956 + 91 91 0
## 704 chr21 14357962 14357962 + 104 93 11
## 705 chr21 14358042 14358042 - 123 123 0
## 706 chr21 14358047 14358047 - 122 118 4
## 707 chr21 14358077 14358077 + 75 75 0
## 708 chr21 14358078 14358078 - 123 115 8
## 709 chr21 14358090 14358090 + 74 73 1
## 710 chr21 14358091 14358091 - 76 73 3
## 711 chr21 14358092 14358092 + 75 74 1
## 712 chr21 14358093 14358093 - 77 76 1
## 713 chr21 14358096 14358096 + 75 74 1
## 714 chr21 14358097 14358097 - 77 73 4
## 715 chr21 14358116 14358116 + 72 64 8
## 716 chr21 14358117 14358117 - 79 78 1
## 717 chr21 14358131 14358131 + 82 82 0
## 718 chr21 14358132 14358132 - 79 78 1
## 719 chr21 14358134 14358134 + 82 81 1
## 720 chr21 14358143 14358143 + 79 79 0
## 721 chr21 14358146 14358146 + 80 79 1
## 722 chr21 14358154 14358154 + 80 78 2
## 723 chr21 14358155 14358155 - 80 77 3
## 724 chr21 14358178 14358178 + 79 36 43
## 725 chr21 14358179 14358179 - 100 44 56
## 726 chr21 14358180 14358180 + 78 71 6
## 727 chr21 14358181 14358181 - 104 102 2
## 728 chr21 14358186 14358186 - 104 100 4
## 729 chr21 14358196 14358196 - 104 98 6
## 730 chr21 14358203 14358203 - 104 102 2
## 731 chr21 14358275 14358275 - 58 58 0
## 732 chr21 14358284 14358284 - 59 58 1
## 733 chr21 14358291 14358291 - 59 58 1
## 734 chr21 14358314 14358314 - 58 54 4
## 735 chr21 14358487 14358487 + 85 84 1
## 736 chr21 14358528 14358528 + 87 82 5
## 737 chr21 14358534 14358534 + 86 79 7
## 738 chr21 14358577 14358577 - 60 56 4
## 739 chr21 14358585 14358585 - 57 55 2
## 740 chr21 14358592 14358592 - 57 57 0
## 741 chr21 14358604 14358604 + 95 95 0
## 742 chr21 14358612 14358612 + 95 43 52
## 743 chr21 14358635 14358635 + 93 93 0
## 744 chr21 14358648 14358648 + 84 83 0
## 745 chr21 14358651 14358651 + 88 88 0
## 746 chr21 14358700 14358700 - 51 48 3
## 747 chr21 14358725 14358725 - 53 50 3
## 748 chr21 14358738 14358738 - 53 52 1
## 749 chr21 14360174 14360174 - 34 29 5
## 750 chr21 14360181 14360181 - 33 29 4
## 751 chr21 14364879 14364879 + 29 29 0
## 752 chr21 14364882 14364882 + 29 29 0
## 753 chr21 14364894 14364894 + 28 28 0
## 754 chr21 14364896 14364896 + 29 29 0
## 755 chr21 14364917 14364917 + 28 28 0
## 756 chr21 14364979 14364979 - 18 18 0
## 757 chr21 14364981 14364981 - 18 18 0
## 758 chr21 14364995 14364995 - 18 18 0
## 759 chr21 14365015 14365015 - 18 18 0
## 760 chr21 14365019 14365019 - 18 18 0
## 761 chr21 14365025 14365025 - 18 18 0
## 762 chr21 14365260 14365260 - 25 22 3
## 763 chr21 14365278 14365278 + 11 11 0
## 764 chr21 14365296 14365296 + 11 11 0
## 765 chr21 14365648 14365648 + 39 35 4
## 766 chr21 14365665 14365665 + 40 37 3
## 767 chr21 14365671 14365671 + 40 37 3
## 768 chr21 14365683 14365683 + 40 37 3
## 769 chr21 14365688 14365688 + 39 36 3
## 770 chr21 14365696 14365696 + 34 29 5
## 771 chr21 14365748 14365748 - 30 30 0
## 772 chr21 14365755 14365755 - 32 31 1
## 773 chr21 14365777 14365777 - 31 31 0
## 774 chr21 14365790 14365790 + 10 6 4
## 775 chr21 14365791 14365791 - 32 32 0
## 776 chr21 14366063 14366063 + 31 31 0
## 777 chr21 14366096 14366096 + 30 30 0
## 778 chr21 14366108 14366108 + 30 28 2
## 779 chr21 14366211 14366211 - 21 20 1
## 780 chr21 14366216 14366216 - 22 21 1
## 781 chr21 14366245 14366245 - 22 18 4
## 782 chr21 14366254 14366254 - 20 19 1
## 783 chr21 14366580 14366580 + 29 27 2
## 784 chr21 14366581 14366581 - 15 10 5
## 785 chr21 14366595 14366595 + 29 29 0
## 786 chr21 14366605 14366605 + 36 36 0
## 787 chr21 14366612 14366612 + 36 36 0
## 788 chr21 14366625 14366625 + 35 35 0
## 789 chr21 14366626 14366626 - 15 15 0
## 790 chr21 14366655 14366655 - 15 14 1
## 791 chr21 14366658 14366658 - 15 15 0
## 792 chr21 14366663 14366663 - 15 15 0
## 793 chr21 14366922 14366922 + 10 10 0
## 794 chr21 14366923 14366923 - 20 18 2
## 795 chr21 14366939 14366939 + 10 8 2
## 796 chr21 14366940 14366940 - 19 16 3
## 797 chr21 14366954 14366954 + 10 10 0
## 798 chr21 14366955 14366955 - 20 17 3
## 799 chr21 14366972 14366972 + 36 36 0
## 800 chr21 14366973 14366973 - 20 20 0
## 801 chr21 14367008 14367008 + 37 37 0
## 802 chr21 14367009 14367009 - 31 27 4
## 803 chr21 14367010 14367010 + 37 35 2
## 804 chr21 14367011 14367011 - 32 28 4
## 805 chr21 14367014 14367014 + 37 37 0
## 806 chr21 14367015 14367015 - 32 28 4
## 807 chr21 14367045 14367045 - 33 28 5
## 808 chr21 14367054 14367054 - 32 26 6
## 809 chr21 14367059 14367059 + 13 13 0
## 810 chr21 14367062 14367062 + 13 13 0
## 811 chr21 14367092 14367092 + 13 12 1
## 812 chr21 14367099 14367099 + 12 12 0
## 813 chr21 14367103 14367103 + 12 12 0
## 814 chr21 14367105 14367105 + 11 11 0
## 815 chr21 14367175 14367175 + 19 13 6
## 816 chr21 14367190 14367190 + 19 13 6
## 817 chr21 14367207 14367207 + 17 9 8
## 818 chr21 14367213 14367213 + 17 9 8
## 819 chr21 14367220 14367220 + 17 10 7
## 820 chr21 14367242 14367242 - 15 14 1
## 821 chr21 14367268 14367268 - 15 15 0
## 822 chr21 14367278 14367278 + 11 11 0
## 823 chr21 14367279 14367279 - 25 17 8
## 824 chr21 14367292 14367292 + 11 11 0
## 825 chr21 14367293 14367293 - 10 10 0
## 826 chr21 14367296 14367296 + 11 11 0
## 827 chr21 14367297 14367297 - 10 9 1
## 828 chr21 14367310 14367310 + 11 11 0
## 829 chr21 14367311 14367311 - 10 10 0
## 830 chr21 14367325 14367325 + 18 11 7
## 831 chr21 14367326 14367326 - 10 10 0
## 832 chr21 14367414 14367414 - 12 12 0
## 833 chr21 14367456 14367456 - 14 13 1
## 834 chr21 14367459 14367459 - 14 14 0
## 835 chr21 14367464 14367464 + 26 23 3
## 836 chr21 14367497 14367497 + 25 25 0
## 837 chr21 14367498 14367498 - 22 22 0
## 838 chr21 14367515 14367515 - 23 22 1
## 839 chr21 14367517 14367517 - 23 23 0
## 840 chr21 14367534 14367534 - 21 21 0
## 841 chr21 14367621 14367621 - 13 9 4
## 842 chr21 14367628 14367628 - 14 13 1
## 843 chr21 14367665 14367665 + 15 12 3
## 844 chr21 14367666 14367666 - 14 13 1
## 845 chr21 14367701 14367701 + 15 15 0
## 846 chr21 14367702 14367702 - 16 14 2
## 847 chr21 14367705 14367705 + 15 14 1
## 848 chr21 14367706 14367706 - 16 15 1
## 849 chr21 14367716 14367716 - 16 12 4
## 850 chr21 14367739 14367739 + 11 8 3
## 851 chr21 14367740 14367740 - 17 16 1
## 852 chr21 14367882 14367882 - 38 38 0
## 853 chr21 14367894 14367894 - 37 37 0
## 854 chr21 14367961 14367961 - 14 14 0
## 855 chr21 14367979 14367979 - 14 14 0
## 856 chr21 14368400 14368400 + 17 14 3
## 857 chr21 14368415 14368415 + 16 7 9
## 858 chr21 14368438 14368438 + 17 13 4
## 859 chr21 14368445 14368445 + 12 10 2
## 860 chr21 14368446 14368446 - 14 7 7
## 861 chr21 14368475 14368475 - 12 11 1
## 862 chr21 14368481 14368481 - 12 11 1
## 863 chr21 14368490 14368490 - 12 11 1
## 864 chr21 14368539 14368539 + 13 13 0
## 865 chr21 14368552 14368552 + 13 10 3
## 866 chr21 14368553 14368553 - 10 10 0
## 867 chr21 14368590 14368590 - 10 10 0
## 868 chr21 14369784 14369784 + 10 6 4
## 869 chr21 14370025 14370025 + 28 25 3
## 870 chr21 14370143 14370143 - 23 21 2
## 871 chr21 14370153 14370153 - 23 21 2
## 872 chr21 14370169 14370169 - 238 205 33
## 873 chr21 14370181 14370181 - 238 213 25
## 874 chr21 14370193 14370193 + 11 11 0
## 875 chr21 14370195 14370195 + 11 11 0
## 876 chr21 14370310 14370310 + 18 17 1
## 877 chr21 14370318 14370318 + 17 15 2
## 878 chr21 14370352 14370352 + 14 14 0
## 879 chr21 14370353 14370353 - 24 24 0
## 880 chr21 14370359 14370359 + 14 13 1
## 881 chr21 14370360 14370360 - 25 25 0
## 882 chr21 14370376 14370376 - 25 25 0
## 883 chr21 14370415 14370415 + 11 11 0
## 884 chr21 14370521 14370521 + 15 13 2
## 885 chr21 14370528 14370528 + 15 15 0
## 886 chr21 14370529 14370529 - 13 9 4
## 887 chr21 14370535 14370535 + 15 15 0
## 888 chr21 14370536 14370536 - 14 13 1
## 889 chr21 14370559 14370559 + 13 13 0
## 890 chr21 14370560 14370560 - 14 13 1
## 891 chr21 14370576 14370576 - 14 13 1
## 892 chr21 14371662 14371662 + 18 10 8
## 893 chr21 14371668 14371668 + 18 10 8
## 894 chr21 14371678 14371678 + 18 18 0
## 895 chr21 14371698 14371698 + 18 18 0
## 896 chr21 14372027 14372027 + 10 9 1
## 897 chr21 14372031 14372031 + 10 9 1
## 898 chr21 14372039 14372039 + 10 9 1
## 899 chr21 14372041 14372041 + 10 9 1
## 900 chr21 14372043 14372043 + 10 9 1
## 901 chr21 14372845 14372845 - 11 0 11
## 902 chr21 14373909 14373909 + 18 18 0
## 903 chr21 14374515 14374515 + 20 20 0
## 904 chr21 14374542 14374542 + 19 19 0
## 905 chr21 14374557 14374557 + 20 20 0
## 906 chr21 14374686 14374686 + 24 24 0
## 907 chr21 14374689 14374689 + 24 23 1
## 908 chr21 14374698 14374698 + 24 21 3
## 909 chr21 14376417 14376417 - 75 0 75
## 910 chr21 14376427 14376427 + 10 9 1
## 911 chr21 14376428 14376428 - 60 0 60
## 912 chr21 14376512 14376512 - 19 18 1
## 913 chr21 14376555 14376555 + 16 16 0
## 914 chr21 14376556 14376556 - 22 21 1
## 915 chr21 14376609 14376609 - 35 25 10
## 916 chr21 14376611 14376611 - 32 30 2
## 917 chr21 14376644 14376644 - 64 7 57
## 918 chr21 14377513 14377513 + 19 17 2
## 919 chr21 14377532 14377532 + 19 16 3
## 920 chr21 14377533 14377533 - 14 12 2
## 921 chr21 14377537 14377537 + 19 15 4
## 922 chr21 14377538 14377538 - 14 14 0
## 923 chr21 14377576 14377576 - 14 14 0
## 924 chr21 14377582 14377582 - 14 14 0
## 925 chr21 14378318 14378318 + 18 17 1
## 926 chr21 14378325 14378325 + 18 15 3
## 927 chr21 14379889 14379889 - 27 25 2
## 928 chr21 14379907 14379907 - 27 26 1
## 929 chr21 14381228 14381228 + 13 12 1
## 930 chr21 14381276 14381276 + 11 11 0
## 931 chr21 14510325 14510325 + 32 1 31
## 932 chr21 14510336 14510336 + 31 0 31
## 933 chr21 14510337 14510337 - 14 0 14
## 934 chr21 14510363 14510363 + 33 0 33
## 935 chr21 14510364 14510364 - 14 0 14
## 936 chr21 14510368 14510368 - 14 0 14
## 937 chr21 14510373 14510373 + 176 13 163
## 938 chr21 14510374 14510374 - 14 2 12
## 939 chr21 14510402 14510402 + 176 8 168
## 940 chr21 14510411 14510411 + 174 4 170
## 941 chr21 14510480 14510480 - 53 0 53
## 942 chr21 14510496 14510496 - 55 3 52
## 943 chr21 14510499 14510499 - 54 6 48
## 944 chr21 14510507 14510507 - 54 0 54
## 945 chr21 14510514 14510514 - 55 1 54
## 946 chr21 14510517 14510517 - 54 10 44
## 947 chr21 14510522 14510522 - 55 7 48
## 948 chr21 14568171 14568171 + 18 2 16
## 949 chr21 14568176 14568176 + 16 2 14
## 950 chr21 14568183 14568183 + 18 2 16
## 951 chr21 14568196 14568196 + 18 1 17
## 952 chr21 14568204 14568204 + 17 1 16
## 953 chr21 14568220 14568220 + 18 0 18
## 954 chr21 14576853 14576853 + 13 13 0
## 955 chr21 14576856 14576856 + 13 13 0
## 956 chr21 14576858 14576858 + 13 13 0
## 957 chr21 14576864 14576864 + 13 13 0
## 958 chr21 14576901 14576901 + 12 12 0
## 959 chr21 14677252 14677252 - 11 0 11
## 960 chr21 14677256 14677256 - 10 0 10
## 961 chr21 14677277 14677277 - 11 0 11
## 962 chr21 14677504 14677504 - 22 0 22
## 963 chr21 14677516 14677516 - 34 0 34
## 964 chr21 14677535 14677535 - 35 0 35
## 965 chr21 14677540 14677540 - 35 0 35
## 966 chr21 14677544 14677544 - 35 0 35
## 967 chr21 14773120 14773120 + 15 11 4
## 968 chr21 14791643 14791643 + 13 12 1
## 969 chr21 14791677 14791677 + 13 12 1
## 970 chr21 15176925 15176925 - 11 11 0
## 971 chr21 15181507 15181507 + 22 19 3
## 972 chr21 15181541 15181541 + 21 18 3
## 973 chr21 15181542 15181542 - 18 14 4
## 974 chr21 15181553 15181553 + 21 17 4
## 975 chr21 15181554 15181554 - 23 20 3
## 976 chr21 15181558 15181558 - 24 22 2
## 977 chr21 15241981 15241981 + 43 41 2
## 978 chr21 15286353 15286353 + 14 14 0
## 979 chr21 15286362 15286362 + 15 15 0
## 980 chr21 15286377 15286377 + 15 4 11
## 981 chr21 15358331 15358331 + 116 2 114
## 982 chr21 15358368 15358368 + 123 0 123
## 983 chr21 15358374 15358374 + 122 0 122
## 984 chr21 15358375 15358375 - 67 0 67
## 985 chr21 15358384 15358384 - 62 0 62
## 986 chr21 15358390 15358390 - 68 0 68
## 987 chr21 15358401 15358401 - 67 0 67
## 988 chr21 15358403 15358403 - 68 0 68
## 989 chr21 15358409 15358409 - 68 0 68
## 990 chr21 15358418 15358418 - 68 0 68
## 991 chr21 15358423 15358423 - 68 1 67
## 992 chr21 15358444 15358444 - 10 0 10
## 993 chr21 15358476 15358476 + 155 0 155
## 994 chr21 15358477 15358477 - 11 0 11
## 995 chr21 15358487 15358487 + 155 0 155
## 996 chr21 15358488 15358488 - 120 2 118
## 997 chr21 15358514 15358514 + 138 0 138
## 998 chr21 15358515 15358515 - 124 2 122
## 999 chr21 15358521 15358521 + 131 0 131
## 1000 chr21 15358522 15358522 - 113 1 112
## 1001 chr21 15358524 15358524 + 137 0 137
## 1002 chr21 15358525 15358525 - 120 1 119
## 1003 chr21 15358536 15358536 - 121 0 121
## 1004 chr21 15358552 15358552 - 19 0 19
## 1005 chr21 15358559 15358559 - 20 0 20
## 1006 chr21 15358562 15358562 - 20 0 20
## 1007 chr21 15358567 15358567 - 20 0 20
## 1008 chr21 15358571 15358571 - 20 0 20
## 1009 chr21 15358577 15358577 + 29 0 29
## 1010 chr21 15358578 15358578 - 42 0 42
## 1011 chr21 15358584 15358584 + 29 0 29
## 1012 chr21 15358585 15358585 - 22 0 22
## 1013 chr21 15358590 15358590 + 29 0 29
## 1014 chr21 15358591 15358591 - 22 0 22
## 1015 chr21 15358602 15358602 + 28 0 28
## 1016 chr21 15358603 15358603 - 22 0 22
## 1017 chr21 15358606 15358606 + 24 0 24
## 1018 chr21 15358607 15358607 - 22 0 22
## 1019 chr21 15358618 15358618 + 28 0 28
## 1020 chr21 15358619 15358619 - 22 0 22
## 1021 chr21 15358623 15358623 + 28 0 28
## 1022 chr21 15358624 15358624 - 22 4 18
## 1023 chr21 15358637 15358637 + 149 0 149
## 1024 chr21 15358643 15358643 + 149 0 149
## 1025 chr21 15358656 15358656 + 148 0 148
## 1026 chr21 15358658 15358658 + 141 0 141
## 1027 chr21 15358659 15358659 - 121 0 121
## 1028 chr21 15358661 15358661 + 145 0 145
## 1029 chr21 15358662 15358662 - 123 0 123
## 1030 chr21 15358665 15358665 + 148 0 148
## 1031 chr21 15358666 15358666 - 95 0 95
## 1032 chr21 15358689 15358689 - 144 0 144
## 1033 chr21 15358693 15358693 - 145 0 145
## 1034 chr21 15358697 15358697 - 144 0 142
## 1035 chr21 15358702 15358702 - 145 0 144
## 1036 chr21 15358707 15358707 - 144 2 142
## 1037 chr21 15358767 15358767 - 21 0 21
## 1038 chr21 15358795 15358795 + 125 0 125
## 1039 chr21 15358796 15358796 - 21 0 21
## 1040 chr21 15358801 15358801 + 127 0 127
## 1041 chr21 15358803 15358803 + 123 0 123
## 1042 chr21 15358812 15358812 + 126 0 126
## 1043 chr21 15358815 15358815 + 126 0 126
## 1044 chr21 15358820 15358820 + 126 0 126
## 1045 chr21 15358823 15358823 + 126 0 126
## 1046 chr21 15358829 15358829 + 123 0 123
## 1047 chr21 15358830 15358830 - 68 0 68
## 1048 chr21 15358831 15358831 + 123 0 123
## 1049 chr21 15358832 15358832 - 79 0 79
## 1050 chr21 15358849 15358849 - 84 0 84
## 1051 chr21 15358856 15358856 - 84 0 84
## 1052 chr21 15358859 15358859 - 84 0 84
## 1053 chr21 15358861 15358861 - 84 0 84
## 1054 chr21 15358866 15358866 - 80 0 80
## 1055 chr21 15358868 15358868 - 80 0 80
## 1056 chr21 15358874 15358874 - 82 0 82
## 1057 chr21 15358893 15358893 + 97 1 96
## 1058 chr21 15358896 15358896 + 96 1 95
## 1059 chr21 15358900 15358900 + 97 0 97
## 1060 chr21 15358912 15358912 + 96 0 96
## 1061 chr21 15358921 15358921 + 98 0 98
## 1062 chr21 15358922 15358922 - 115 0 115
## 1063 chr21 15358923 15358923 + 94 0 94
## 1064 chr21 15358924 15358924 - 145 0 145
## 1065 chr21 15358927 15358927 + 95 0 95
## 1066 chr21 15358928 15358928 - 156 0 156
## 1067 chr21 15358945 15358945 - 166 0 166
## 1068 chr21 15358948 15358948 - 176 0 176
## 1069 chr21 15358950 15358950 - 178 0 178
## 1070 chr21 15358963 15358963 + 68 0 68
## 1071 chr21 15358964 15358964 - 177 0 177
## 1072 chr21 15358969 15358969 + 68 0 68
## 1073 chr21 15358970 15358970 - 45 0 45
## 1074 chr21 15358978 15358978 + 68 0 68
## 1075 chr21 15358979 15358979 - 53 0 53
## 1076 chr21 15358983 15358983 + 68 0 68
## 1077 chr21 15358984 15358984 - 54 0 54
## 1078 chr21 15358985 15358985 + 66 0 66
## 1079 chr21 15358986 15358986 - 57 0 57
## 1080 chr21 15358991 15358991 + 58 0 58
## 1081 chr21 15358992 15358992 - 49 0 49
## 1082 chr21 15358998 15358998 + 68 0 68
## 1083 chr21 15358999 15358999 - 55 0 55
## 1084 chr21 15359003 15359003 + 68 0 68
## 1085 chr21 15359004 15359004 - 57 0 57
## 1086 chr21 15359014 15359014 - 57 0 57
## 1087 chr21 15359021 15359021 + 65 1 64
## 1088 chr21 15359024 15359024 + 65 0 65
## 1089 chr21 15359038 15359038 + 65 0 65
## 1090 chr21 15359045 15359045 + 57 0 57
## 1091 chr21 15359048 15359048 + 48 0 48
## 1092 chr21 15359054 15359054 + 63 0 63
## 1093 chr21 15359057 15359057 + 48 0 48
## 1094 chr21 15359062 15359062 + 38 0 38
## 1095 chr21 15359065 15359065 + 23 0 23
## 1096 chr21 15359069 15359069 + 22 0 22
## 1097 chr21 15359072 15359072 - 135 0 135
## 1098 chr21 15359080 15359080 - 137 0 137
## 1099 chr21 15359084 15359084 - 136 0 136
## 1100 chr21 15359088 15359088 - 137 0 137
## 1101 chr21 15359092 15359092 - 135 0 135
## 1102 chr21 15359094 15359094 - 135 0 135
## 1103 chr21 15359096 15359096 - 131 0 131
## 1104 chr21 15359098 15359098 - 138 0 138
## 1105 chr21 15359105 15359105 - 138 0 138
## 1106 chr21 15359109 15359109 - 138 0 138
## 1107 chr21 15359113 15359113 - 138 0 138
## 1108 chr21 15359121 15359121 - 137 0 137
## 1109 chr21 15359358 15359358 - 101 0 101
## 1110 chr21 15359382 15359382 + 101 0 101
## 1111 chr21 15359386 15359386 + 98 0 98
## 1112 chr21 15359402 15359402 + 98 0 98
## 1113 chr21 15359403 15359403 - 147 0 147
## 1114 chr21 15359419 15359419 + 99 0 99
## 1115 chr21 15359420 15359420 - 151 0 151
## 1116 chr21 15359435 15359435 - 147 0 147
## 1117 chr21 15359446 15359446 - 151 2 149
## 1118 chr21 15359473 15359473 + 69 0 69
## 1119 chr21 15359482 15359482 + 67 0 67
## 1120 chr21 15359483 15359483 - 52 0 52
## 1121 chr21 15359491 15359491 + 68 0 68
## 1122 chr21 15359492 15359492 - 48 0 48
## 1123 chr21 15359493 15359493 + 64 0 64
## 1124 chr21 15359494 15359494 - 52 0 52
## 1125 chr21 15359499 15359499 + 68 0 68
## 1126 chr21 15359500 15359500 - 54 0 54
## 1127 chr21 15359515 15359515 + 63 0 63
## 1128 chr21 15359516 15359516 - 54 0 54
## 1129 chr21 15359518 15359518 + 50 2 48
## 1130 chr21 15359519 15359519 - 52 0 52
## 1131 chr21 15359520 15359520 + 45 0 45
## 1132 chr21 15359521 15359521 - 54 0 54
## 1133 chr21 15359525 15359525 - 54 0 54
## 1134 chr21 15359536 15359536 + 105 0 105
## 1135 chr21 15359558 15359558 + 106 0 106
## 1136 chr21 15359559 15359559 - 37 0 37
## 1137 chr21 15359562 15359562 + 105 0 105
## 1138 chr21 15359563 15359563 - 89 1 88
## 1139 chr21 15359564 15359564 + 103 0 103
## 1140 chr21 15359565 15359565 - 119 0 119
## 1141 chr21 15359568 15359568 + 105 0 105
## 1142 chr21 15359569 15359569 - 119 0 119
## 1143 chr21 15359577 15359577 + 86 0 86
## 1144 chr21 15359578 15359578 - 121 0 121
## 1145 chr21 15359602 15359602 + 122 0 122
## 1146 chr21 15359603 15359603 - 126 0 126
## 1147 chr21 15359611 15359611 + 119 3 116
## 1148 chr21 15359629 15359629 + 121 0 121
## 1149 chr21 15359630 15359630 - 152 0 152
## 1150 chr21 15359639 15359639 + 122 0 122
## 1151 chr21 15359640 15359640 - 150 0 150
## 1152 chr21 15359662 15359662 - 152 0 152
## 1153 chr21 15360037 15360037 + 57 0 57
## 1154 chr21 15360106 15360106 - 58 0 58
## 1155 chr21 15360111 15360111 - 58 0 58
## 1156 chr21 15360118 15360118 - 58 0 58
## 1157 chr21 15360143 15360143 - 58 0 58
## 1158 chr21 15360149 15360149 - 58 0 58
## 1159 chr21 15447185 15447185 - 28 11 17
## 1160 chr21 15582781 15582781 - 11 10 1
## 1161 chr21 15582804 15582804 - 12 9 3
## 1162 chr21 15583554 15583554 - 28 28 0
## 1163 chr21 15583576 15583576 - 28 23 5
## 1164 chr21 15583584 15583584 + 13 9 4
## 1165 chr21 15583585 15583585 - 28 28 0
## 1166 chr21 15583593 15583593 + 13 13 0
## 1167 chr21 15583594 15583594 - 16 16 0
## 1168 chr21 15583628 15583628 - 17 8 9
## 1169 chr21 15583637 15583637 - 17 7 10
## 1170 chr21 15583673 15583673 - 11 4 7
## 1171 chr21 15583716 15583716 - 14 6 8
## 1172 chr21 15584415 15584415 + 15 10 5
## 1173 chr21 15584444 15584444 + 15 10 5
## 1174 chr21 15584560 15584560 - 16 8 8
## 1175 chr21 15584596 15584596 - 16 8 8
## 1176 chr21 15584720 15584720 + 12 8 4
## 1177 chr21 15715029 15715029 + 20 19 1
## 1178 chr21 15788691 15788691 + 13 13 0
## 1179 chr21 15788694 15788694 + 13 13 0
## 1180 chr21 15788697 15788697 + 13 10 3
## 1181 chr21 15788732 15788732 + 13 13 0
## 1182 chr21 15788733 15788733 - 21 19 2
## 1183 chr21 15788752 15788752 - 21 21 0
## 1184 chr21 15788755 15788755 - 21 21 0
## 1185 chr21 16014722 16014722 - 12 12 0
## 1186 chr21 16014734 16014734 - 12 12 0
## 1187 chr21 16014738 16014738 - 12 12 0
## 1188 chr21 16023758 16023758 + 82 0 82
## 1189 chr21 16023765 16023765 + 81 0 81
## 1190 chr21 16023775 16023775 + 81 0 81
## 1191 chr21 16023788 16023788 + 80 1 79
## 1192 chr21 16023845 16023845 - 20 0 20
## 1193 chr21 16023859 16023859 - 18 0 18
## 1194 chr21 16023871 16023871 - 20 0 20
## 1195 chr21 16023880 16023880 + 87 0 87
## 1196 chr21 16023881 16023881 - 20 0 20
## 1197 chr21 16023891 16023891 + 89 0 89
## 1198 chr21 16023893 16023893 + 89 0 89
## 1199 chr21 16023900 16023900 + 87 0 87
## 1200 chr21 16023901 16023901 - 101 0 101
## 1201 chr21 16023910 16023910 + 85 0 85
## 1202 chr21 16023911 16023911 - 54 0 54
## 1203 chr21 16023916 16023916 + 80 0 80
## 1204 chr21 16023917 16023917 - 98 0 98
## 1205 chr21 16023919 16023919 + 86 0 86
## 1206 chr21 16023920 16023920 - 68 0 68
## 1207 chr21 16023921 16023921 + 86 0 86
## 1208 chr21 16023922 16023922 - 78 0 78
## 1209 chr21 16023932 16023932 - 119 0 118
## 1210 chr21 16023940 16023940 - 119 0 119
## 1211 chr21 16023947 16023947 - 118 3 115
## 1212 chr21 16023980 16023980 + 83 1 82
## 1213 chr21 16023986 16023986 + 83 0 83
## 1214 chr21 16023995 16023995 + 83 0 83
## 1215 chr21 16023996 16023996 - 73 0 73
## 1216 chr21 16023997 16023997 + 82 0 82
## 1217 chr21 16023998 16023998 - 78 0 78
## 1218 chr21 16024009 16024009 + 33 0 33
## 1219 chr21 16024010 16024010 - 85 0 85
## 1220 chr21 16024022 16024022 + 75 0 75
## 1221 chr21 16024023 16024023 - 84 1 83
## 1222 chr21 16024028 16024028 + 45 0 45
## 1223 chr21 16024029 16024029 - 62 0 62
## 1224 chr21 16024034 16024034 - 78 0 78
## 1225 chr21 16024037 16024037 - 82 0 82
## 1226 chr21 16024041 16024041 - 84 0 83
## 1227 chr21 16024102 16024102 + 58 0 58
## 1228 chr21 16024106 16024106 + 57 0 57
## 1229 chr21 16024112 16024112 + 46 0 46
## 1230 chr21 16024121 16024121 + 55 0 55
## 1231 chr21 16024123 16024123 + 58 0 58
## 1232 chr21 16024149 16024149 + 50 0 50
## 1233 chr21 16024150 16024150 - 12 0 12
## 1234 chr21 16024151 16024151 + 33 0 33
## 1235 chr21 16024152 16024152 - 25 0 25
## 1236 chr21 16024169 16024169 - 24 0 24
## 1237 chr21 16024176 16024176 + 12 0 12
## 1238 chr21 16024177 16024177 - 26 1 25
## 1239 chr21 16024180 16024180 + 12 0 12
## 1240 chr21 16024184 16024184 + 12 0 12
## 1241 chr21 16024185 16024185 - 32 0 32
## 1242 chr21 16024193 16024193 + 12 0 12
## 1243 chr21 16024194 16024194 - 108 0 108
## 1244 chr21 16024224 16024224 - 114 1 113
## 1245 chr21 16024226 16024226 - 115 0 115
## 1246 chr21 16024230 16024230 + 179 2 177
## 1247 chr21 16024231 16024231 - 117 1 116
## 1248 chr21 16024244 16024244 + 178 0 178
## 1249 chr21 16024245 16024245 - 150 0 150
## 1250 chr21 16024248 16024248 + 177 0 177
## 1251 chr21 16024249 16024249 - 153 0 153
## 1252 chr21 16024252 16024252 + 175 0 175
## 1253 chr21 16024253 16024253 - 150 0 150
## 1254 chr21 16024254 16024254 + 170 0 170
## 1255 chr21 16024255 16024255 - 153 0 153
## 1256 chr21 16024259 16024259 + 177 0 177
## 1257 chr21 16024260 16024260 - 108 0 108
## 1258 chr21 16024262 16024262 + 177 0 177
## 1259 chr21 16024263 16024263 - 137 0 137
## 1260 chr21 16024266 16024266 + 174 0 174
## 1261 chr21 16024267 16024267 - 156 0 156
## 1262 chr21 16024275 16024275 + 170 1 169
## 1263 chr21 16024276 16024276 - 150 0 150
## 1264 chr21 16024277 16024277 + 156 0 156
## 1265 chr21 16024278 16024278 - 144 0 144
## 1266 chr21 16024286 16024286 - 156 0 156
## 1267 chr21 16024293 16024293 - 155 0 155
## 1268 chr21 16024317 16024317 + 72 0 72
## 1269 chr21 16024322 16024322 + 72 0 72
## 1270 chr21 16024325 16024325 + 72 0 72
## 1271 chr21 16024327 16024327 + 72 0 72
## 1272 chr21 16024331 16024331 + 72 0 72
## 1273 chr21 16024335 16024335 + 67 0 67
## 1274 chr21 16024338 16024338 + 64 0 64
## 1275 chr21 16024346 16024346 + 61 0 61
## 1276 chr21 16024351 16024351 + 71 0 71
## 1277 chr21 16024353 16024353 + 67 0 67
## 1278 chr21 16024361 16024361 + 66 0 66
## 1279 chr21 16024391 16024391 - 48 0 48
## 1280 chr21 16024394 16024394 - 35 0 35
## 1281 chr21 16024399 16024399 - 22 0 22
## 1282 chr21 16024402 16024402 - 28 0 28
## 1283 chr21 16024410 16024410 - 54 0 54
## 1284 chr21 16024412 16024412 - 55 0 55
## 1285 chr21 16024419 16024419 - 55 1 54
## 1286 chr21 16024421 16024421 - 55 0 55
## 1287 chr21 16024425 16024425 - 55 1 54
## 1288 chr21 16024427 16024427 + 87 0 87
## 1289 chr21 16024428 16024428 - 55 2 53
## 1290 chr21 16024431 16024431 + 86 0 86
## 1291 chr21 16024437 16024437 + 87 0 87
## 1292 chr21 16024443 16024443 + 87 0 87
## 1293 chr21 16024445 16024445 + 86 0 86
## 1294 chr21 16024456 16024456 + 80 0 80
## 1295 chr21 16024462 16024462 + 42 0 42
## 1296 chr21 16024479 16024479 - 92 0 92
## 1297 chr21 16024483 16024483 - 91 0 91
## 1298 chr21 16024485 16024485 - 84 0 84
## 1299 chr21 16024488 16024488 - 84 0 84
## 1300 chr21 16024502 16024502 - 107 0 107
## 1301 chr21 16024505 16024505 - 108 0 108
## 1302 chr21 16024510 16024510 - 109 0 109
## 1303 chr21 16024512 16024512 - 109 0 109
## 1304 chr21 16024525 16024525 - 109 0 109
## 1305 chr21 16024543 16024543 + 109 0 109
## 1306 chr21 16024547 16024547 + 105 0 105
## 1307 chr21 16024553 16024553 + 109 0 109
## 1308 chr21 16024556 16024556 + 108 0 108
## 1309 chr21 16024559 16024559 + 108 0 108
## 1310 chr21 16024562 16024562 + 109 0 109
## 1311 chr21 16024565 16024565 + 107 0 107
## 1312 chr21 16024568 16024568 + 108 0 108
## 1313 chr21 16024571 16024571 + 107 0 107
## 1314 chr21 16024574 16024574 + 103 0 103
## 1315 chr21 16024576 16024576 + 87 0 87
## 1316 chr21 16024589 16024589 + 83 0 83
## 1317 chr21 16024602 16024602 - 50 0 50
## 1318 chr21 16024617 16024617 - 51 0 51
## 1319 chr21 16024619 16024619 - 53 0 53
## 1320 chr21 16024622 16024622 - 51 0 51
## 1321 chr21 16024636 16024636 - 54 0 54
## 1322 chr21 16024642 16024642 - 53 0 53
## 1323 chr21 16024648 16024648 + 16 0 16
## 1324 chr21 16024649 16024649 - 54 0 54
## 1325 chr21 16024651 16024651 + 16 0 16
## 1326 chr21 16024670 16024670 + 16 0 16
## 1327 chr21 16024673 16024673 + 16 0 16
## 1328 chr21 16024688 16024688 + 16 0 16
## 1329 chr21 16024689 16024689 - 11 0 11
## 1330 chr21 16024693 16024693 + 13 0 13
## 1331 chr21 16024694 16024694 - 11 0 11
## 1332 chr21 16024697 16024697 + 14 0 14
## 1333 chr21 16024698 16024698 - 11 0 11
## 1334 chr21 16033213 16033213 - 18 14 4
## 1335 chr21 16076341 16076341 + 11 11 0
## 1336 chr21 16076355 16076355 + 11 11 0
## 1337 chr21 16167784 16167784 + 18 18 0
## 1338 chr21 16167808 16167808 + 18 15 3
## 1339 chr21 16167809 16167809 - 12 12 0
## 1340 chr21 16167845 16167845 - 12 12 0
## 1341 chr21 16202432 16202432 + 33 21 12
## 1342 chr21 16291204 16291204 + 31 24 7
## 1343 chr21 16291208 16291208 + 31 27 4
## 1344 chr21 16291220 16291220 + 27 26 1
## 1345 chr21 16291265 16291265 - 22 17 5
## 1346 chr21 16291273 16291273 - 22 18 4
## 1347 chr21 16439628 16439628 - 22 16 6
## 1348 chr21 16454568 16454568 + 16 10 6
## 1349 chr21 16454615 16454615 + 10 8 2
## 1350 chr21 16454616 16454616 - 10 6 4
## 1351 chr21 16454628 16454628 - 10 7 3
## 1352 chr21 16454650 16454650 - 10 8 2
## 1353 chr21 16713387 16713387 + 15 9 6
## 1354 chr21 16713421 16713421 + 14 8 6
## 1355 chr21 16713423 16713423 + 14 9 5
## 1356 chr21 16751394 16751394 + 11 11 0
## 1357 chr21 16847597 16847597 + 13 5 8
## 1358 chr21 16849374 16849374 + 36 11 25
## 1359 chr21 16858301 16858301 - 11 0 11
## 1360 chr21 16858305 16858305 - 11 11 0
## 1361 chr21 16858315 16858315 - 11 9 2
## 1362 chr21 16858317 16858317 - 11 10 1
## 1363 chr21 16858350 16858350 - 11 11 0
## 1364 chr21 16869844 16869844 + 39 34 5
## 1365 chr21 16869860 16869860 + 39 26 13
## 1366 chr21 16869890 16869890 + 34 13 21
## 1367 chr21 16932700 16932700 - 10 4 6
## 1368 chr21 16932730 16932730 - 11 5 6
## 1369 chr21 17024371 17024371 - 11 7 4
## 1370 chr21 17024381 17024381 - 11 10 1
## 1371 chr21 17024391 17024391 - 11 10 1
## 1372 chr21 17024403 17024403 - 11 10 1
## 1373 chr21 17043261 17043261 + 12 12 0
## 1374 chr21 17043269 17043269 + 12 12 0
## 1375 chr21 17208072 17208072 + 11 2 9
## 1376 chr21 17385151 17385151 - 235 206 29
## 1377 chr21 17419460 17419460 - 14 7 7
## 1378 chr21 17419551 17419551 - 11 5 6
## 1379 chr21 17655329 17655329 - 10 9 1
## 1380 chr21 17655337 17655337 - 10 10 0
## 1381 chr21 17655352 17655352 - 10 10 0
## 1382 chr21 17655355 17655355 - 10 10 0
## 1383 chr21 17668415 17668415 - 16 4 12
## 1384 chr21 17668448 17668448 - 13 13 0
## 1385 chr21 17704512 17704512 + 11 10 1
## 1386 chr21 17704519 17704519 + 11 7 4
## 1387 chr21 17704561 17704561 + 10 0 10
## 1388 chr21 17719707 17719707 - 18 18 0
## 1389 chr21 17719743 17719743 - 18 16 2
## 1390 chr21 17731698 17731698 + 14 7 7
## 1391 chr21 17731701 17731701 + 14 8 6
## 1392 chr21 17731709 17731709 + 13 8 5
## 1393 chr21 17757042 17757042 + 11 10 1
## 1394 chr21 17757046 17757046 + 11 10 1
## 1395 chr21 17757048 17757048 + 11 10 1
## 1396 chr21 17757058 17757058 + 11 8 3
## 1397 chr21 17757083 17757083 + 11 9 2
## 1398 chr21 17771970 17771970 + 14 14 0
## 1399 chr21 17771990 17771990 + 14 14 0
## 1400 chr21 17772004 17772004 + 14 9 5
## 1401 chr21 17772015 17772015 + 14 14 0
## 1402 chr21 17803059 17803059 - 56 16 40
## 1403 chr21 17803088 17803088 - 55 15 40
## 1404 chr21 17804725 17804725 + 36 28 8
## 1405 chr21 17804726 17804726 - 23 20 3
## 1406 chr21 17804800 17804800 - 43 42 1
## 1407 chr21 17804826 17804826 - 16 16 0
## 1408 chr21 17804828 17804828 - 16 16 0
## 1409 chr21 17804839 17804839 + 47 44 3
## 1410 chr21 17804840 17804840 - 16 14 2
## 1411 chr21 17807150 17807150 + 54 0 54
## 1412 chr21 17807159 17807159 + 54 0 54
## 1413 chr21 17807161 17807161 + 52 0 52
## 1414 chr21 17807174 17807174 + 54 0 54
## 1415 chr21 17807178 17807178 + 51 0 51
## 1416 chr21 17807197 17807197 + 16 0 16
## 1417 chr21 17807401 17807401 + 36 1 35
## 1418 chr21 17807407 17807407 + 32 0 32
## 1419 chr21 17807411 17807411 + 36 0 36
## 1420 chr21 17807418 17807418 + 36 0 36
## 1421 chr21 17807420 17807420 + 36 0 36
## 1422 chr21 17807442 17807442 + 33 0 33
## 1423 chr21 17807488 17807488 - 40 0 40
## 1424 chr21 17807491 17807491 - 40 0 40
## 1425 chr21 17807496 17807496 - 40 0 40
## 1426 chr21 17807516 17807516 - 40 1 39
## 1427 chr21 17807529 17807529 + 72 1 71
## 1428 chr21 17807530 17807530 - 40 0 40
## 1429 chr21 17807533 17807533 + 68 2 66
## 1430 chr21 17807534 17807534 - 70 0 70
## 1431 chr21 17807545 17807545 + 71 1 70
## 1432 chr21 17807546 17807546 - 89 0 89
## 1433 chr21 17807550 17807550 + 72 0 72
## 1434 chr21 17807551 17807551 - 89 0 89
## 1435 chr21 17807560 17807560 + 71 0 71
## 1436 chr21 17807561 17807561 - 86 0 86
## 1437 chr21 17807562 17807562 + 68 0 68
## 1438 chr21 17807563 17807563 - 86 4 82
## 1439 chr21 17807569 17807569 + 62 0 62
## 1440 chr21 17807570 17807570 - 89 0 89
## 1441 chr21 17807571 17807571 + 64 0 64
## 1442 chr21 17807572 17807572 - 88 3 85
## 1443 chr21 17807573 17807573 + 63 1 62
## 1444 chr21 17807574 17807574 - 90 0 90
## 1445 chr21 17807582 17807582 + 68 0 68
## 1446 chr21 17807583 17807583 - 90 0 90
## 1447 chr21 17807588 17807588 + 68 0 68
## 1448 chr21 17807597 17807597 + 67 0 67
## 1449 chr21 17807600 17807600 + 61 0 61
## 1450 chr21 17807613 17807613 + 20 0 20
## 1451 chr21 17807623 17807623 + 60 0 60
## 1452 chr21 17807627 17807627 + 42 0 42
## 1453 chr21 17807628 17807628 - 79 0 79
## 1454 chr21 17807629 17807629 + 32 0 32
## 1455 chr21 17807630 17807630 - 85 0 85
## 1456 chr21 17807633 17807633 - 85 0 85
## 1457 chr21 17807638 17807638 - 84 0 84
## 1458 chr21 17807644 17807644 - 63 0 63
## 1459 chr21 17807649 17807649 - 78 0 78
## 1460 chr21 17807676 17807676 + 64 0 64
## 1461 chr21 17807677 17807677 - 88 2 86
## 1462 chr21 17807681 17807681 + 64 0 64
## 1463 chr21 17807703 17807703 + 64 0 64
## 1464 chr21 17807706 17807706 + 61 0 61
## 1465 chr21 17807725 17807725 + 62 1 61
## 1466 chr21 17807801 17807801 - 116 2 114
## 1467 chr21 17807821 17807821 - 114 0 114
## 1468 chr21 17807850 17807850 - 119 0 119
## 1469 chr21 17809554 17809554 + 43 43 0
## 1470 chr21 17809596 17809596 + 46 44 2
## 1471 chr21 17809723 17809723 - 49 49 0
## 1472 chr21 17809874 17809874 - 29 29 0
## 1473 chr21 17809912 17809912 - 30 27 3
## 1474 chr21 17809920 17809920 - 30 30 0
## 1475 chr21 17832330 17832330 - 13 12 1
## 1476 chr21 17898354 17898354 + 11 11 0
## 1477 chr21 17898358 17898358 + 11 9 2
## 1478 chr21 17898370 17898370 + 11 11 0
## 1479 chr21 17898384 17898384 + 11 11 0
## 1480 chr21 17898491 17898491 - 11 11 0
## 1481 chr21 17906151 17906151 + 66 0 66
## 1482 chr21 17906156 17906156 + 69 0 69
## 1483 chr21 17906165 17906165 + 72 0 72
## 1484 chr21 17906168 17906168 + 71 0 71
## 1485 chr21 17906353 17906353 + 51 0 51
## 1486 chr21 17906387 17906387 + 51 0 51
## 1487 chr21 17906408 17906408 - 57 0 57
## 1488 chr21 17906420 17906420 - 57 0 57
## 1489 chr21 17906430 17906430 - 56 0 56
## 1490 chr21 17906439 17906439 - 59 0 59
## 1491 chr21 17906441 17906441 - 61 0 61
## 1492 chr21 17906447 17906447 - 60 0 60
## 1493 chr21 17906456 17906456 - 61 0 61
## 1494 chr21 17906490 17906490 + 113 0 113
## 1495 chr21 17906496 17906496 + 109 0 109
## 1496 chr21 17906498 17906498 + 115 0 115
## 1497 chr21 17906517 17906517 + 109 0 109
## 1498 chr21 17906537 17906537 + 109 0 109
## 1499 chr21 17906541 17906541 - 34 0 34
## 1500 chr21 17906543 17906543 - 36 0 36
## 1501 chr21 17906571 17906571 - 37 0 37
## 1502 chr21 17906573 17906573 - 37 2 35
## 1503 chr21 17906582 17906582 - 32 0 32
## 1504 chr21 17906589 17906589 + 96 0 96
## 1505 chr21 17906590 17906590 - 37 0 37
## 1506 chr21 17906595 17906595 + 91 0 91
## 1507 chr21 17906605 17906605 + 106 0 106
## 1508 chr21 17906607 17906607 + 105 0 105
## 1509 chr21 17906619 17906619 + 105 0 105
## 1510 chr21 17906688 17906688 - 116 0 116
## 1511 chr21 17906693 17906693 - 79 0 79
## 1512 chr21 17906699 17906699 - 116 0 116
## 1513 chr21 17906701 17906701 - 121 0 121
## 1514 chr21 17906706 17906706 - 120 0 120
## 1515 chr21 17906714 17906714 - 111 0 111
## 1516 chr21 17906722 17906722 - 122 0 122
## 1517 chr21 17906725 17906725 - 122 0 122
## 1518 chr21 17906729 17906729 + 11 0 11
## 1519 chr21 17906730 17906730 - 139 0 139
## 1520 chr21 17906775 17906775 + 45 0 45
## 1521 chr21 17906776 17906776 - 18 0 18
## 1522 chr21 17906780 17906780 + 35 0 35
## 1523 chr21 17906783 17906783 + 36 0 36
## 1524 chr21 17906787 17906787 + 36 0 36
## 1525 chr21 17906793 17906793 + 36 0 36
## 1526 chr21 17906801 17906801 + 35 0 35
## 1527 chr21 17906804 17906804 + 36 0 36
## 1528 chr21 17906810 17906810 + 36 0 36
## 1529 chr21 17906818 17906818 + 26 0 26
## 1530 chr21 17906821 17906821 + 30 0 30
## 1531 chr21 17906838 17906838 - 43 0 43
## 1532 chr21 17906853 17906853 - 45 0 45
## 1533 chr21 17906872 17906872 - 47 0 47
## 1534 chr21 17906946 17906946 + 121 0 121
## 1535 chr21 17906962 17906962 + 120 0 120
## 1536 chr21 17906970 17906970 + 121 0 121
## 1537 chr21 17906971 17906971 - 127 0 127
## 1538 chr21 17906976 17906976 + 121 0 121
## 1539 chr21 17906977 17906977 - 127 0 127
## 1540 chr21 17906978 17906978 + 117 0 117
## 1541 chr21 17906979 17906979 - 125 0 125
## 1542 chr21 17906995 17906995 + 81 0 81
## 1543 chr21 17906996 17906996 - 128 2 126
## 1544 chr21 17907001 17907001 - 131 0 131
## 1545 chr21 17907008 17907008 - 131 0 131
## 1546 chr21 17907011 17907011 + 95 0 95
## 1547 chr21 17907012 17907012 - 130 1 129
## 1548 chr21 17907014 17907014 + 95 0 95
## 1549 chr21 17907016 17907016 + 95 0 95
## 1550 chr21 17907020 17907020 + 94 0 94
## 1551 chr21 17907027 17907027 + 94 0 94
## 1552 chr21 17907029 17907029 + 92 0 92
## 1553 chr21 17907032 17907032 + 94 0 94
## 1554 chr21 17907035 17907035 + 93 0 93
## 1555 chr21 17907039 17907039 + 87 0 87
## 1556 chr21 17907047 17907047 + 89 1 88
## 1557 chr21 17907049 17907049 + 90 0 90
## 1558 chr21 17907051 17907051 + 86 0 86
## 1559 chr21 17907087 17907087 - 74 0 74
## 1560 chr21 17907095 17907095 - 72 0 72
## 1561 chr21 17907108 17907108 - 71 0 71
## 1562 chr21 17907111 17907111 - 69 0 69
## 1563 chr21 17907113 17907113 - 70 0 70
## 1564 chr21 17907124 17907124 - 69 0 69
## 1565 chr21 17907126 17907126 - 70 0 70
## 1566 chr21 17907131 17907131 - 71 0 71
## 1567 chr21 17907154 17907154 + 71 0 71
## 1568 chr21 17907168 17907168 + 70 0 70
## 1569 chr21 17907169 17907169 - 79 0 79
## 1570 chr21 17907170 17907170 + 65 0 65
## 1571 chr21 17907171 17907171 - 80 0 80
## 1572 chr21 17907172 17907172 + 69 1 68
## 1573 chr21 17907173 17907173 - 78 0 78
## 1574 chr21 17907174 17907174 + 67 0 67
## 1575 chr21 17907175 17907175 - 81 0 81
## 1576 chr21 17907180 17907180 + 69 0 69
## 1577 chr21 17907181 17907181 - 82 0 82
## 1578 chr21 17907185 17907185 + 68 0 68
## 1579 chr21 17907186 17907186 - 80 0 80
## 1580 chr21 17907187 17907187 + 65 0 65
## 1581 chr21 17907188 17907188 - 81 0 81
## 1582 chr21 17907189 17907189 + 63 0 63
## 1583 chr21 17907190 17907190 - 81 0 81
## 1584 chr21 17907200 17907200 + 11 0 11
## 1585 chr21 17907201 17907201 - 81 0 81
## 1586 chr21 17907206 17907206 + 10 0 10
## 1587 chr21 17907207 17907207 - 82 0 82
## 1588 chr21 17907213 17907213 - 16 0 16
## 1589 chr21 17907220 17907220 + 10 0 10
## 1590 chr21 17907221 17907221 - 16 0 16
## 1591 chr21 17907223 17907223 + 10 0 10
## 1592 chr21 17907224 17907224 - 16 0 16
## 1593 chr21 17907225 17907225 + 10 0 10
## 1594 chr21 17907226 17907226 - 16 0 16
## 1595 chr21 17907237 17907237 + 10 0 10
## 1596 chr21 17907238 17907238 - 16 0 16
## 1597 chr21 17907244 17907244 + 62 0 62
## 1598 chr21 17907245 17907245 - 16 0 16
## 1599 chr21 17907254 17907254 + 62 0 62
## 1600 chr21 17907258 17907258 + 62 0 62
## 1601 chr21 17907265 17907265 + 37 0 37
## 1602 chr21 17907276 17907276 + 62 0 62
## 1603 chr21 17907306 17907306 - 58 0 58
## 1604 chr21 17907310 17907310 - 61 0 61
## 1605 chr21 17907312 17907312 - 76 0 76
## 1606 chr21 17907315 17907315 - 80 0 80
## 1607 chr21 17907318 17907318 - 72 0 72
## 1608 chr21 17907320 17907320 - 101 0 101
## 1609 chr21 17907325 17907325 - 98 0 98
## 1610 chr21 17907333 17907333 - 95 0 95
## 1611 chr21 17907335 17907335 - 103 0 103
## 1612 chr21 17907349 17907349 - 100 0 100
## 1613 chr21 17907387 17907387 + 76 0 76
## 1614 chr21 17907391 17907391 + 76 0 76
## 1615 chr21 17907398 17907398 + 74 2 72
## 1616 chr21 17907402 17907402 + 74 0 74
## 1617 chr21 17907408 17907408 + 72 7 65
## 1618 chr21 17907414 17907414 + 74 3 71
## 1619 chr21 17907421 17907421 + 74 2 72
## 1620 chr21 17907423 17907423 + 72 0 72
## 1621 chr21 17907425 17907425 + 73 4 69
## 1622 chr21 17907505 17907505 - 289 98 191
## 1623 chr21 17907513 17907513 - 290 85 205
## 1624 chr21 17907539 17907539 + 31 0 31
## 1625 chr21 17907540 17907540 - 283 6 277
## 1626 chr21 17907558 17907558 + 31 5 26
## 1627 chr21 17907564 17907564 + 31 5 26
## 1628 chr21 17907567 17907567 + 31 8 23
## 1629 chr21 17907732 17907732 - 86 16 70
## 1630 chr21 17907754 17907754 - 88 42 46
## 1631 chr21 17951753 17951753 - 11 11 0
## 1632 chr21 17983314 17983314 + 10 1 9
## 1633 chr21 17988383 17988383 + 12 12 0
## 1634 chr21 17988393 17988393 + 12 8 4
## 1635 chr21 17988405 17988405 + 12 11 1
## 1636 chr21 18012188 18012188 + 104 101 3
## 1637 chr21 18012222 18012222 + 88 85 3
## 1638 chr21 18014317 18014317 + 60 51 9
## 1639 chr21 18014321 18014321 + 59 55 4
## 1640 chr21 18014329 18014329 + 60 55 5
## 1641 chr21 18014358 18014358 + 54 54 0
## 1642 chr21 18014366 18014366 + 47 0 47
## 1643 chr21 18016826 18016826 + 17 16 1
## 1644 chr21 18016832 18016832 + 17 16 1
## 1645 chr21 18038986 18038986 + 80 77 3
## 1646 chr21 18038990 18038990 + 80 79 1
## 1647 chr21 18038995 18038995 + 79 74 5
## 1648 chr21 18039018 18039018 + 78 69 9
## 1649 chr21 18039019 18039019 - 38 34 4
## 1650 chr21 18039055 18039055 - 37 35 2
## 1651 chr21 18052003 18052003 - 18 15 3
## 1652 chr21 18052040 18052040 + 12 10 2
## 1653 chr21 18052041 18052041 - 18 17 1
## 1654 chr21 18052058 18052058 + 12 12 0
## 1655 chr21 18052081 18052081 + 11 11 0
## 1656 chr21 18052172 18052172 - 13 12 1
## 1657 chr21 18052188 18052188 - 13 9 4
## 1658 chr21 18052194 18052194 - 13 10 3
## 1659 chr21 18052433 18052433 - 10 0 10
## 1660 chr21 18113002 18113002 - 14 0 14
## 1661 chr21 18113006 18113006 - 14 0 14
## 1662 chr21 18113013 18113013 - 13 0 13
## 1663 chr21 18113016 18113016 - 13 0 13
## 1664 chr21 18113031 18113031 - 13 0 13
## 1665 chr21 18113040 18113040 - 13 0 13
## 1666 chr21 18113048 18113048 - 13 1 12
## 1667 chr21 18113070 18113070 + 15 0 15
## 1668 chr21 18113084 18113084 + 15 0 15
## 1669 chr21 18113085 18113085 - 12 0 12
## 1670 chr21 18113091 18113091 + 15 0 15
## 1671 chr21 18113092 18113092 - 12 0 12
## 1672 chr21 18113095 18113095 + 13 0 13
## 1673 chr21 18113096 18113096 - 12 0 12
## 1674 chr21 18113098 18113098 + 12 0 12
## 1675 chr21 18113099 18113099 - 12 0 12
## 1676 chr21 18113106 18113106 + 12 0 12
## 1677 chr21 18113107 18113107 - 12 0 12
## 1678 chr21 18113113 18113113 - 12 0 12
## 1679 chr21 18113173 18113173 + 10 0 10
## 1680 chr21 18113191 18113191 + 10 0 10
## 1681 chr21 18113199 18113199 + 10 0 10
## 1682 chr21 18113204 18113204 + 24 0 24
## 1683 chr21 18113230 18113230 + 23 0 23
## 1684 chr21 18113234 18113234 + 23 0 23
## 1685 chr21 18113243 18113243 + 23 0 23
## 1686 chr21 18113246 18113246 + 23 0 23
## 1687 chr21 18113250 18113250 + 15 0 15
## 1688 chr21 18113293 18113293 + 52 4 48
## 1689 chr21 18113300 18113300 + 52 0 52
## 1690 chr21 18113303 18113303 + 52 0 52
## 1691 chr21 18113312 18113312 + 49 0 49
## 1692 chr21 18113314 18113314 + 50 0 50
## 1693 chr21 18113316 18113316 + 52 0 52
## 1694 chr21 18113320 18113320 + 51 0 51
## 1695 chr21 18113334 18113334 + 52 0 52
## 1696 chr21 18113336 18113336 + 50 0 50
## 1697 chr21 18113450 18113450 - 19 0 19
## 1698 chr21 18113480 18113480 - 18 0 18
## 1699 chr21 18113492 18113492 - 18 0 18
## 1700 chr21 18113495 18113495 + 86 0 86
## 1701 chr21 18113496 18113496 - 18 1 17
## 1702 chr21 18113501 18113501 + 86 0 86
## 1703 chr21 18113503 18113503 + 84 0 84
## 1704 chr21 18113514 18113514 + 86 0 86
## 1705 chr21 18113515 18113515 - 135 0 135
## 1706 chr21 18113517 18113517 + 84 0 84
## 1707 chr21 18113518 18113518 - 154 0 154
## 1708 chr21 18113522 18113522 + 85 0 85
## 1709 chr21 18113523 18113523 - 155 0 155
## 1710 chr21 18113526 18113526 + 84 0 84
## 1711 chr21 18113527 18113527 - 155 0 155
## 1712 chr21 18113531 18113531 + 85 0 85
## 1713 chr21 18113532 18113532 - 157 0 157
## 1714 chr21 18113553 18113553 + 32 0 32
## 1715 chr21 18113554 18113554 - 209 2 207
## 1716 chr21 18113557 18113557 + 32 0 32
## 1717 chr21 18113558 18113558 - 52 0 52
## 1718 chr21 18113568 18113568 + 25 0 25
## 1719 chr21 18113569 18113569 - 62 0 62
## 1720 chr21 18113572 18113572 + 32 0 32
## 1721 chr21 18113573 18113573 - 58 0 58
## 1722 chr21 18113575 18113575 + 32 0 32
## 1723 chr21 18113576 18113576 - 61 0 61
## 1724 chr21 18113577 18113577 + 30 0 30
## 1725 chr21 18113578 18113578 - 60 0 60
## 1726 chr21 18113579 18113579 + 28 0 28
## 1727 chr21 18113580 18113580 - 54 0 54
## 1728 chr21 18113589 18113589 + 32 0 32
## 1729 chr21 18113590 18113590 - 55 0 55
## 1730 chr21 18113599 18113599 + 30 0 30
## 1731 chr21 18113600 18113600 - 63 0 63
## 1732 chr21 18113602 18113602 + 167 0 167
## 1733 chr21 18113603 18113603 - 64 0 64
## 1734 chr21 18113613 18113613 + 142 0 142
## 1735 chr21 18113624 18113624 + 138 0 138
## 1736 chr21 18113625 18113625 - 109 0 109
## 1737 chr21 18113636 18113636 + 94 0 94
## 1738 chr21 18113637 18113637 - 101 0 101
## 1739 chr21 18113638 18113638 + 126 0 126
## 1740 chr21 18113639 18113639 - 110 0 110
## 1741 chr21 18113642 18113642 + 134 0 134
## 1742 chr21 18113643 18113643 - 106 0 106
## 1743 chr21 18113644 18113644 + 126 0 126
## 1744 chr21 18113645 18113645 - 111 0 111
## 1745 chr21 18113657 18113657 - 110 0 110
## 1746 chr21 18113659 18113659 - 110 1 109
## 1747 chr21 18113663 18113663 - 111 0 111
## 1748 chr21 18113665 18113665 + 84 0 84
## 1749 chr21 18113666 18113666 - 111 0 111
## 1750 chr21 18113681 18113681 + 84 0 84
## 1751 chr21 18113689 18113689 + 84 0 84
## 1752 chr21 18113700 18113700 + 84 0 84
## 1753 chr21 18113703 18113703 + 84 0 84
## 1754 chr21 18113704 18113704 - 106 0 106
## 1755 chr21 18113706 18113706 + 83 0 83
## 1756 chr21 18113707 18113707 - 117 0 117
## 1757 chr21 18113714 18113714 + 81 0 81
## 1758 chr21 18113715 18113715 - 123 0 122
## 1759 chr21 18113720 18113720 - 138 0 138
## 1760 chr21 18113738 18113738 - 137 0 137
## 1761 chr21 18113750 18113750 + 28 0 28
## 1762 chr21 18113751 18113751 - 135 1 134
## 1763 chr21 18113763 18113763 + 28 0 28
## 1764 chr21 18113775 18113775 + 27 0 27
## 1765 chr21 18113817 18113817 - 93 0 93
## 1766 chr21 18113823 18113823 - 96 0 96
## 1767 chr21 18113863 18113863 - 97 0 97
## 1768 chr21 18127776 18127776 + 21 21 0
## 1769 chr21 18174355 18174355 + 11 8 3
## 1770 chr21 18178516 18178516 + 20 17 3
## 1771 chr21 18178528 18178528 + 20 14 6
## 1772 chr21 18178550 18178550 + 20 15 5
## 1773 chr21 18178551 18178551 - 11 9 2
## 1774 chr21 18178575 18178575 - 11 7 4
## 1775 chr21 18178579 18178579 - 11 4 7
## 1776 chr21 18178587 18178587 - 11 8 3
## 1777 chr21 18178591 18178591 - 11 8 3
## 1778 chr21 18292976 18292976 - 24 0 24
## 1779 chr21 18440282 18440282 - 14 13 1
## 1780 chr21 18440305 18440305 - 17 17 0
## 1781 chr21 18440313 18440313 - 17 14 3
## 1782 chr21 18440317 18440317 - 17 14 3
## 1783 chr21 18491657 18491657 + 14 11 3
## 1784 chr21 18491661 18491661 + 14 11 3
## 1785 chr21 18498731 18498731 - 13 6 7
## 1786 chr21 18513019 18513019 - 14 13 1
## 1787 chr21 18522637 18522637 + 10 10 0
## 1788 chr21 18522659 18522659 + 10 0 10
## 1789 chr21 18525466 18525466 - 10 10 0
## 1790 chr21 18525469 18525469 - 10 10 0
## 1791 chr21 18538816 18538816 + 12 1 11
## 1792 chr21 18538970 18538970 + 160 0 160
## 1793 chr21 18538978 18538978 + 160 1 159
## 1794 chr21 18538983 18538983 + 157 0 157
## 1795 chr21 18538998 18538998 + 158 12 146
## 1796 chr21 18539009 18539009 + 147 19 128
## 1797 chr21 18539010 18539010 - 58 7 51
## 1798 chr21 18539016 18539016 + 140 31 109
## 1799 chr21 18539017 18539017 - 62 8 54
## 1800 chr21 18539048 18539048 - 61 16 45
## 1801 chr21 18539053 18539053 - 62 11 51
## 1802 chr21 18539117 18539117 + 25 0 25
## 1803 chr21 18539131 18539131 + 25 1 24
## 1804 chr21 18539148 18539148 + 25 0 25
## 1805 chr21 18539208 18539208 - 75 45 30
## 1806 chr21 18539218 18539218 - 79 37 42
## 1807 chr21 18539220 18539220 - 79 32 47
## 1808 chr21 18539222 18539222 - 80 28 52
## 1809 chr21 18539231 18539231 - 79 38 41
## 1810 chr21 18539233 18539233 - 80 24 56
## 1811 chr21 18539240 18539240 - 77 14 63
## 1812 chr21 18539249 18539249 + 96 45 51
## 1813 chr21 18539250 18539250 - 77 39 38
## 1814 chr21 18539253 18539253 + 94 51 43
## 1815 chr21 18539259 18539259 + 90 19 71
## 1816 chr21 18539264 18539264 + 94 22 72
## 1817 chr21 18539265 18539265 - 71 15 56
## 1818 chr21 18539266 18539266 + 93 50 43
## 1819 chr21 18539267 18539267 - 71 39 32
## 1820 chr21 18539276 18539276 + 96 7 89
## 1821 chr21 18539277 18539277 - 72 3 69
## 1822 chr21 18539290 18539290 + 89 2 87
## 1823 chr21 18539291 18539291 - 73 3 70
## 1824 chr21 18539294 18539294 + 68 7 61
## 1825 chr21 18539295 18539295 - 73 6 67
## 1826 chr21 18539297 18539297 + 27 2 25
## 1827 chr21 18539298 18539298 - 73 1 72
## 1828 chr21 18539310 18539310 + 18 0 18
## 1829 chr21 18539311 18539311 - 73 2 71
## 1830 chr21 18539344 18539344 + 18 3 15
## 1831 chr21 18539346 18539346 + 16 0 16
## 1832 chr21 18539353 18539353 + 16 1 15
## 1833 chr21 18539540 18539540 - 11 0 11
## 1834 chr21 18539545 18539545 - 10 0 10
## 1835 chr21 18539547 18539547 - 11 0 11
## 1836 chr21 18539550 18539550 - 11 1 10
## 1837 chr21 18539565 18539565 - 12 0 12
## 1838 chr21 18539569 18539569 - 12 0 12
## 1839 chr21 18539580 18539580 - 12 0 12
## 1840 chr21 18539585 18539585 - 10 0 10
## 1841 chr21 18539612 18539612 + 17 0 17
## 1842 chr21 18539613 18539613 - 10 0 10
## 1843 chr21 18539622 18539622 + 16 0 16
## 1844 chr21 18539801 18539801 - 21 2 19
## 1845 chr21 18539804 18539804 + 124 14 110
## 1846 chr21 18539805 18539805 - 21 0 21
## 1847 chr21 18539838 18539838 + 80 4 76
## 1848 chr21 18539933 18539933 - 111 62 49
## 1849 chr21 18539960 18539960 - 116 37 79
## 1850 chr21 18539970 18539970 - 111 4 107
## 1851 chr21 18539976 18539976 - 117 17 100
## 1852 chr21 18550887 18550887 + 15 8 7
## 1853 chr21 18572834 18572834 + 10 8 2
## 1854 chr21 18572866 18572866 + 10 10 0
## 1855 chr21 18576141 18576141 + 15 15 0
## 1856 chr21 18576145 18576145 + 14 12 2
## 1857 chr21 18576157 18576157 + 14 14 0
## 1858 chr21 18576178 18576178 + 11 11 0
## 1859 chr21 18576179 18576179 - 11 11 0
## 1860 chr21 18576203 18576203 - 11 5 6
## 1861 chr21 18576211 18576211 - 11 11 0
## 1862 chr21 18582970 18582970 + 18 12 6
## 1863 chr21 18582975 18582975 + 18 11 7
## 1864 chr21 18615698 18615698 - 24 16 8
## 1865 chr21 18713624 18713624 + 12 11 1
## 1866 chr21 18713634 18713634 + 12 12 0
## 1867 chr21 18713648 18713648 + 12 8 4
## 1868 chr21 18737675 18737675 + 12 7 5
## 1869 chr21 18737730 18737730 - 13 8 5
## 1870 chr21 18763228 18763228 + 11 0 11
## 1871 chr21 18763284 18763284 - 10 0 10
## 1872 chr21 18957876 18957876 + 16 9 7
## 1873 chr21 19280741 19280741 + 10 8 2
## 1874 chr21 19282009 19282009 + 11 11 0
## 1875 chr21 19390691 19390691 - 14 8 6
## 1876 chr21 19390694 19390694 - 14 9 5
## 1877 chr21 19422671 19422671 - 17 15 2
## 1878 chr21 19539486 19539486 - 10 7 3
## 1879 chr21 19539511 19539511 - 12 8 4
## 1880 chr21 19539528 19539528 - 12 12 0
## 1881 chr21 19579289 19579289 + 23 10 13
## 1882 chr21 19579345 19579345 - 10 8 2
## 1883 chr21 19579352 19579352 - 12 6 6
## 1884 chr21 19579357 19579357 - 12 10 2
## 1885 chr21 19579361 19579361 - 12 10 2
## 1886 chr21 19626227 19626227 - 25 22 3
## 1887 chr21 19626384 19626384 - 22 22 0
## 1888 chr21 19711774 19711774 + 10 2 8
## 1889 chr21 19739452 19739452 + 10 9 1
## 1890 chr21 19788073 19788073 - 45 42 3
## 1891 chr21 19788076 19788076 - 44 42 2
## 1892 chr21 19788089 19788089 - 45 45 0
## 1893 chr21 19788093 19788093 - 45 45 0
## 1894 chr21 19788113 19788113 - 46 46 0
## 1895 chr21 19788117 19788117 - 46 46 0
## 1896 chr21 19788122 19788122 - 46 43 3
## 1897 chr21 19834041 19834041 + 10 10 0
## 1898 chr21 19834070 19834070 + 10 3 7
## 1899 chr21 19855659 19855659 + 10 10 0
## 1900 chr21 19855690 19855690 + 27 26 1
## 1901 chr21 19855706 19855706 + 27 27 0
## 1902 chr21 19855711 19855711 + 18 18 0
## 1903 chr21 19943653 19943653 + 12 12 0
## 1904 chr21 19943695 19943695 + 12 11 1
You can also convert methylDB objects to their in-memory equivalents. Since that requires an additional parameter (the directory where the files will be located), we have a different function, named makeMethylDB
to achieve this goal. Below, we convert a methylBase object to methylBaseDB
and saving it at “exMethylDB” directory.
data(methylKit)
objDB=makeMethylDB(methylBase.obj,"exMethylDB")
We can also select rows from methylRaw
, methylBase
and methylDiff
objects and methylDB pendants with select
function. An appropriate methylKit object will be returned as a result of select
function. Or you can use the '['
notation to subset the methylKit objects.
select(meth,1:5) # get first 10 rows of a methylBase object
## chr start end strand coverage1 numCs1 numTs1 coverage2 numCs2
## 1 chr21 9853296 9853296 + 17 10 7 333 268
## 2 chr21 9853326 9853326 + 17 12 5 329 249
## 3 chr21 9860126 9860126 + 39 38 1 83 78
## 4 chr21 9906604 9906604 + 68 42 26 111 97
## 5 chr21 9906616 9906616 + 68 52 16 111 104
## numTs2 coverage3 numCs3 numTs3 coverage4 numCs4 numTs4
## 1 65 18 16 2 395 341 54
## 2 79 16 14 2 379 284 95
## 3 5 83 83 0 41 40 1
## 4 14 23 18 5 37 33 4
## 5 7 23 14 9 37 27 10
myDiff[21:25,] # get 5 rows of a methylDiff object
## chr start end strand pvalue qvalue meth.diff
## 21 chr21 9913543 9913543 + 1.254379e-02 2.632641e-02 -13.343109
## 22 chr21 9913575 9913575 + 2.755448e-01 3.161628e-01 -5.442623
## 23 chr21 9927527 9927527 + 1.120126e-07 9.257475e-07 -46.109840
## 24 chr21 9944505 9944505 + 7.907909e-20 7.515975e-18 -51.017943
## 25 chr21 9944663 9944663 - 1.790779e-05 7.678302e-05 -28.099911
Important: Using select
or '['
on methylDB objects will return its normal methylKit
pendant, to avoid overhead of database operations.
We can select rows from any methylKit object, that lie inside the ranges of a GRanges
object from GenomicRanges
package with selectByOverlap
function. An appropriate methylKit object will be returned as a result of selectByOverlap
function.
library(GenomicRanges)
my.win=GRanges(seqnames="chr21",
ranges=IRanges(start=seq(from=9764513,by=10000,length.out=20),width=5000) )
# selects the records that lie inside the regions
selectByOverlap(myobj[[1]],my.win)
Important: Using selectByOverlap
on methylDB objects will return its normal methylKit
pendant, to avoid overhead of database operations.
The methylBase
and methylRawList
, as well as methylDB pendants can be reorganized by reorganize
function. The function can subset the objects based on provided sample ids, it also creates a new treatment vector determining which samples belong to which group. Order of sample ids should match the treatment vector order.
# creates a new methylRawList object
myobj2=reorganize(myobj,sample.ids=c("test1","ctrl2"),treatment=c(1,0) )
# creates a new methylBase object
meth2 =reorganize(meth,sample.ids=c("test1","ctrl2"),treatment=c(1,0) )
Percent methylation values can be extracted from methylBase
object by using percMethylation
function.
# creates a matrix containing percent methylation values
perc.meth=percMethylation(meth)
Methylation or differential methylation profiles can be segmented to sections that contain similar CpGs with respect to their methylation profiles. This kind of segmentation could help us find interesting regions. For example, segmentation analysis will usually reveal high or low methylated regions, where low methylated regions could be interesting for gene regulation. The algorithm first finds segments that have CpGs with similar methylation levels, then those segments are classified to segment groups based on their mean methylation levels. This enables us to group segments with similar methylation levels to the same class.
See more at http://zvfak.blogspot.de/2015/06/segmentation-of-methylation-profiles.html
download.file("https://dl.dropboxusercontent.com/u/1373164/H1.chr21.chr22.rds",
destfile="H1.chr21.chr22.rds",method="curl")
mbw=readRDS("H1.chr21.chr22.rds")
# it finds the optimal number of componets as 6
res=methSeg(mbw,diagnostic.plot=TRUE,maxInt=100,minSeg=10)
# however the BIC stabilizes after 4, we can also try 4 componets
res=methSeg(mbw,diagnostic.plot=TRUE,maxInt=100,minSeg=10,G=1:4)
# get segments to BED file
methSeg2bed(res,filename="H1.chr21.chr22.trial.seg.bed")
Detailed answers to some of the frequently asked questions and various HOW-TO’s can be found at http://zvfak.blogspot.com/search/label/methylKit. In addition, http://code.google.com/p/methylkit/ has online documentation and links to tutorials and other related material. You can also check methylKit Q&A forum for answers https://groups.google.com/forum/#!forum/methylkit_discussion.
Apart from those here are some of the frequently asked questions.
methylRaw
or methylBase
objects ?See ?select
or help("[", package = "methylKit")
Currently, we will be able to tell you if your regions/bases overlap with the genomic features or not. See ?getMembers
.
See ?genomation::getAssociationWithTSS
Promoters are defined by options at genomation::readTranscriptFeatures
function. The default option is to take -1000,+1000bp around the TSS and you can change that. Same goes for CpG islands when reading them in via genomation::readFeatureFlank
function. Default is to take 2000bp flanking regions on each side of the CpG island as shores. But you can change that as well.
Check the Bismark (Krueger and Andrews 2011) website and there are also example files that ship with the package. Look at their formats and try to run different variations of processBismarkAln()
command on the example files.
methylRawList
or methylBase
objects ?See ?reorganize
methylKit
comes with a simple normalizeCoverage()
function to normalize read coverage distributions between samples. Ideally, you should first filter bases with extreme coverage to account for PCR bias using filterByCoverage()
function, then run normalizeCoverage()
function to normalize coverage between samples. These two functions will help reduce the bias in the statistical tests that might occur due to systematic over-sampling of reads in certain samples.
methylKit
decides which test to use based on number of samples per group. In order to use Fisher’s exact there must be one sample in each of the test and control groups. So if you have multiple samples for group, the package will employ Logistic Regression based test. However, you can use pool()
function to pool samples in each group so that you have one representative sample per group. pool()
function will sum up number of Cs and Ts in each group. We recommend using filterByCoverage()
and normalizeCoverage()
functions prior to using pool()
. See ?pool
Yes, you can. methylKit can read any generic methylation percentage/ratio file as long as that text file contains columns for chromosome, start, end, strand, coverage and number of methylated cytosines. However, methylKit can only process SAM files from Bismark. For other aligners, you need to get a text file containing the minimal information described above. Some aligners will come with scripts or built-in tools to provide such files. See http://zvfak.blogspot.com/2012/10/how-to-read-bsmap-methylation-ratio.html for how to read methylation ratio files from BSMAP (Xi and Li 2009) aligner.
Yes, you can. Many functions of the analysis workflow provide an save.db
argument, which allows you to save the output as methylDB object. For example see ?unite
and also check the ...
argument section for further details.
This package is initially developed at Weill Cornell Medical College by Altuna Akalin with important code contributions from Matthias Kormaksson(mk375@cornell.edu) and Sheng Li (shl2018@med.cornell.edu). We wish to thank especially Maria E. Figueroa, Francine Garret-Bakelman, Christopher Mason and Ari Melnick for their contribution of ideas, data and support. Their support and discussions lead to development of methylKit.
sessionInfo()
## R version 3.3.1 (2016-06-21)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 16.04.1 LTS
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 LC_COLLATE=C
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## attached base packages:
## [1] grid parallel stats4 stats graphics grDevices utils
## [8] datasets methods base
##
## other attached packages:
## [1] genomation_1.6.0 methylKit_1.0.0 GenomicRanges_1.26.0
## [4] GenomeInfoDb_1.10.0 IRanges_2.8.0 S4Vectors_0.12.0
## [7] BiocGenerics_0.20.0
##
## loaded via a namespace (and not attached):
## [1] SummarizedExperiment_1.4.0 qvalue_2.6.0
## [3] gtools_3.5.0 reshape2_1.4.1
## [5] splines_3.3.1 lattice_0.20-34
## [7] seqPattern_1.6.0 colorspace_1.2-7
## [9] htmltools_0.3.5 rtracklayer_1.34.0
## [11] yaml_2.1.13 chron_2.3-47
## [13] XML_3.98-1.4 R.oo_1.20.0
## [15] R.utils_2.4.0 BiocParallel_1.8.0
## [17] fastseg_1.20.0 matrixStats_0.51.0
## [19] plyr_1.8.4 stringr_1.1.0
## [21] zlibbioc_1.20.0 Biostrings_2.42.0
## [23] munsell_0.4.3 gtable_0.2.0
## [25] R.methodsS3_1.7.1 coda_0.18-1
## [27] evaluate_0.10 Biobase_2.34.0
## [29] knitr_1.14 Rcpp_0.12.7
## [31] readr_1.0.0 KernSmooth_2.23-15
## [33] scales_0.4.0 BSgenome_1.42.0
## [35] formatR_1.4 plotrix_3.6-3
## [37] limma_3.30.0 XVector_0.14.0
## [39] Rsamtools_1.26.0 impute_1.48.0
## [41] ggplot2_2.1.0 digest_0.6.10
## [43] stringi_1.1.2 numDeriv_2016.8-1
## [45] tools_3.3.1 bitops_1.0-6
## [47] bbmle_1.0.18 magrittr_1.5
## [49] RCurl_1.95-4.8 tibble_1.2
## [51] MASS_7.3-45 Matrix_1.2-7.1
## [53] gridBase_0.4-7 data.table_1.9.6
## [55] assertthat_0.1 rmarkdown_1.1
## [57] emdbook_1.3.9 mclust_5.2
## [59] GenomicAlignments_1.10.0
Krueger, Felix, and Simon R Andrews. 2011. “Bismark: a flexible aligner and methylation caller for Bisulfite-Seq applications.” Bioinformatics (Oxford, England) 27 (11): 1571–2. doi:10.1093/bioinformatics/btr167.
Li, Heng. 2011. “Tabix: Fast Retrieval of Sequence Features from Generic TAB-Delimited Files.” Bioinformatics 27 (5): 718–19. doi:10.1093/bioinformatics/btq671.
McCullagh, P., and John A. Nelder. 1989. Generalized Linear Models, Second Edition (Chapman & Hall/CRC Monographs on Statistics & Applied Probability). Chapman; Hall/CRC. http://www.amazon.com/Generalized-Chapman-Monographs-Statistics-Probability/dp/0412317605%3FSubscriptionId%3D0JYN1NVW651KCA56C102%26tag%3Dtechkie-20%26linkCode%3Dxm2%26camp%3D2025%26creative%3D165953%26creativeASIN%3D0412317605.
Wang, Hong-Qiang, Lindsey K Tuominen, and Chung-Jui Tsai. 2011. “SLIM: a sliding linear model for estimating the proportion of true null hypotheses in datasets with dependence structures.” Bioinformatics (Oxford, England) 27 (2): 225–31. doi:10.1093/bioinformatics/btq650.
Xi, Yuanxin, and Wei Li. 2009. “BSMAP: whole genome bisulfite sequence MAPping program.” BMC Bioinformatics 10 (1): 232. doi:10.1186/1471-2105-10-232.
Author of the vignette. See Acknowledgements for a list of contributors.↩