| countFragmentOverlaps {FourCSeq} | R Documentation |
countFragmentOverlaps counts the number of reads mapping to each
fragment end in rowRanges of the FourC object.
countFragmentOverlaps(object, trim=0, minMapq=0, shift=0)
object |
A |
trim |
Number of bases that should be trimmed at the start of a read.
This is necessary if the read still contains the restriction enzyme sequence.
Default is |
minMapq |
Minimum mapping quality required for counting the read.
Default is |
shift |
Maximum difference in starts or ends between read and fragment
positions. Default is |
Updated FourC object that contains two new assays
countsLeftFragmentEnd and countsRightFragmentEnd with the
count values at the respective fragment end.
Felix A. Klein, felix.klein@embl.de
FourC, findViewpointFragments,
countFragmentOverlapsSecondCutter
metadata <- list(projectPath=tempdir(),
fragmentDir="re_fragments",
referenceGenomeFile=system.file("extdata/dm3_chr2L_1-6900.fa",
package="FourCSeq"),
reSequence1="GATC",
reSequence2="CATG",
primerFile=system.file("extdata/primer.fa",
package="FourCSeq"),
bamFilePath=system.file("extdata/bam", package="FourCSeq"))
colData <- DataFrame(viewpoint = "testdata",
condition = factor(rep(c("WE_68h", "MESO_68h", "WE_34h"),
each=2),
levels = c("WE_68h", "MESO_68h", "WE_34h")),
replicate = rep(c(1, 2),
3),
bamFile = c("CRM_ap_ApME680_WE_6-8h_1_testdata.bam",
"CRM_ap_ApME680_WE_6-8h_2_testdata.bam",
"CRM_ap_ApME680_MESO_6-8h_1_testdata.bam",
"CRM_ap_ApME680_MESO_6-8h_2_testdata.bam",
"CRM_ap_ApME680_WE_3-4h_1_testdata.bam",
"CRM_ap_ApME680_WE_3-4h_2_testdata.bam"),
sequencingPrimer="first")
fc <- FourC(colData, metadata)
fc
fc <- addFragments(fc)
findViewpointFragments(fc)
fc <- addViewpointFrags(fc)
fc
fc <- countFragmentOverlaps(fc, trim=4, minMapq=30)
fc