Skip to content

Commit

Permalink
Merge pull request #138 from yixuan-chen-elisa/read.modkit
Browse files Browse the repository at this point in the history
Update strand information in read.modkit.R
  • Loading branch information
kasperdanielhansen authored May 9, 2024
2 parents d32c2c6 + b8886bc commit 7c309d3
Showing 1 changed file with 14 additions and 7 deletions.
21 changes: 14 additions & 7 deletions R/read.modkit.R
Original file line number Diff line number Diff line change
Expand Up @@ -9,21 +9,24 @@ read.modkit <- function(files,
}

for (i in seq_along(files)){
data <- read.table(files[i], header = FALSE, sep="\t",
data <- fread(files[i], header = FALSE, sep="\t",
stringsAsFactors=FALSE, quote="")
data$V6[data$V6 == "."] <- "*"

if (length(unique(data$V4)) == 2){
gr <- GRanges(seqnames = data[data$V4 == "m", ]$V1,
ranges = IRanges(start = data[data$V4 == "m", ]$V2+1,
end = data[data$V4 == "m", ]$V3))
end = data[data$V4 == "m", ]$V3),
strand = data[data$V4 == "m", ]$V6)

mcols(gr)$m <- data[data$V4 == "m", ]$V12
mcols(gr)$h <- data[data$V4 != "m", ]$V12
mcols(gr)$u <- data[data$V4 == "m", ]$V13
mcols(gr)$filter <- data[data$V4 == "m", ]$V16
}else{
gr <- GRanges(seqnames = data$V1,
ranges = IRanges(start = data$V2+1, end = data$V3))
ranges = IRanges(start = data$V2+1, end = data$V3),
strand = data$V6)

mcols(gr)$m <- data$V12
mcols(gr)$u <- data$V13
Expand Down Expand Up @@ -58,19 +61,23 @@ read.modkit <- function(files,
Filtered = as.matrix(filter),
coef = NULL, se.coef = NULL,
pos = start(overlap_gr), trans = NULL,
parameters = NULL, pData = colData, gr = NULL,
parameters = NULL, pData = colData, gr = overlap_gr,
chr = as.vector(seqnames(overlap_gr)),
sampleNames = sampleNames, rmZeroCov = rmZeroCov)
if (strandCollapse) {strandCollapse(bsseq_obj)}
if (strandCollapse) {
bsseq_obj <- strandCollapse(bsseq_obj)
}
}else{
bsseq_obj <- BSseq(M = as.matrix(m), Cov = as.matrix(u + m),
Filtered = as.matrix(filter),
coef = NULL, se.coef = NULL,
pos = start(overlap_gr), trans = NULL,
parameters = NULL, pData = colData, gr = NULL,
parameters = NULL, pData = colData, gr = overlap_gr,
chr = as.vector(seqnames(overlap_gr)),
sampleNames = sampleNames, rmZeroCov = rmZeroCov)
if (strandCollapse) {strandCollapse(bsseq_obj)}
if (strandCollapse) {
bsseq_obj <- strandCollapse(bsseq_obj)
}
}

return(bsseq_obj)
Expand Down

0 comments on commit 7c309d3

Please sign in to comment.