-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathsimIntSiteReads.R
217 lines (181 loc) · 9.21 KB
/
simIntSiteReads.R
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
#### simulate R1, R2, I1 fastq files for intSiteCaller, and save the truth file ####
#' Approach:
#' 1. Get 1000 sites, either random or from intsites_miseq database;
#' 2. Simulate template lengths ~ N(300, 100) within [30, 1000], 100 templates for each site;
#' 3. Get sequences for the templates, attach primer and linkers according to strands;
#' 4. Get R1 and R2 according to strand, simulate I1 based on barcode;
#' 5' Read names assigned to
#' R1: @M03249:site(#id):template(#id):(#chr):(#strand):(#pos) 1:N:0:0
#' R2: @M03249:site(#id):template(#id):(#chr):(#strand):(#pos) 2:N:0:0
#' I1: @M03249:site(#id):template(#id):(#chr):(#strand):(#pos) 1:N:0:0
#' 6' Write to fastq.gz files
#'
#' set default and/or get arguments from command line input
#' @param commandline arguments or none for interactive session
#' @return list of default values
#' freeze reference genome
#' sites number to integration sites to simulate
#' sonicLength number of templates to simulate for a site
#' group group in the ~/.my.cnf file
options(stringsAsFactors=FALSE)
get_args <- function() {
suppressMessages(library(argparse))
codeDir <- dirname(sub("--file=", "", grep("--file=", commandArgs(trailingOnly=FALSE), value=T)))
if( length(codeDir)!=1 ) codeDir <- list.files(path="~", pattern="simIntSiteReads$", recursive=TRUE, include.dirs=TRUE, full.names=TRUE)
stopifnot(file.exists(file.path(codeDir, "simIntSiteReads.R")))
stopifnot(file.exists(file.path(codeDir, "processingParams.tsv")))
refGenome <- read.table(file.path(codeDir, "processingParams.tsv"),
header=TRUE)$refGenome
stopifnot(length(refGenome)==1)
parser <- ArgumentParser(formatter_class='argparse.RawTextHelpFormatter')
parser$add_argument("-f", "--freeze", type="character", nargs=1,
default=refGenome,
help="hg18, etc, default read from processingParams.tsv")
parser$add_argument("-o", "--outFolder", type="character", nargs=1,
default="intSiteSimulation",
help="output folder")
parser$add_argument("-s", "--sites", type="integer", nargs=1,
default=5000,
help="number of integration sites")
parser$add_argument("-l", "--sonicLength", type="integer", nargs=1,
default=100,
help="number of sonic lengths for each site")
parser$add_argument("-1", "--R1L", type="integer", nargs=1,
default=179,
help="R1 read length")
parser$add_argument("-2", "--R2L", type="integer", nargs=1,
default=144,
help="R2 read length")
parser$add_argument("-m", "--info", type="character", nargs=1,
default="sampleInfo.tsv",
help="metadata file, default sampleInfo.tsv")
parser$add_argument("-g", "--group", type="character", nargs=1,
default="intsites_miseq.read",
help="number of sonic lengths for each site")
parser$add_argument("-c", "--codeDir", type="character", nargs=1,
default=codeDir,
help="Directory of code")
parser$add_argument("-e", "--errRate", type="character", nargs=1,
default=0,
help="error rate, [0,1]")
args <- parser$parse_args(commandArgs(trailingOnly=TRUE))
args$errRate <- as.numeric(args$errRate)
args$seed <- 123457
return(args)
}
args <- get_args()
print(t(as.data.frame(args)), quote=FALSE)
libs <- c("stringr",
"RMySQL",
"dplyr",
"data.table",
"ggplot2",
"GenomicRanges",
"ShortRead",
"BiocParallel",
"BSgenome",
sprintf("BSgenome.Hsapiens.UCSC.%s", args$freeze))
null <- suppressMessages(sapply(libs, library, character.only=TRUE))
source(file.path(args$codeDir, "simIntSiteReads_func.R"))
source(file.path(args$codeDir, "sequencing_error.R"))
set.seed(args$seed)
options(dplyr.width = Inf)
#' @note this is from one line of the sample information file
#' alias,linkerSequence,bcSeq,gender,primer,ltrBit,largeLTRFrag,vectorSeq
#' GTSP0308-1,GAACGAGCACTAGTAAGCCCNNNNNNNNNNNNCTCCGCTTAAGGGACT,GTATTCGACTTG,m,GAAAATC,TCTAGCA,TGCTAGAGATTTTCCACACTGACTAAAAGGGTCT,vector_WasLenti.fa
sampleInfo <- read.table(file.path(args$codeDir, args$info), header=TRUE)
sampleInfo$linkerSequence <- gsub("N", "T", sampleInfo$linkerSequence)
#' @note this is specific to the integration protocol
oligo <- data.frame(P5="AATGATACGGCGACCACCGA",
P7="CAAGCAGAAGACGGCATACGAGAT",
Spacer="AGTCAGTCAGCC",
SP1="ATCTACACCAGGACTGACGCTATGGTAATTGT",
SP2="AGACCCTTTTAGTCAGTGTG",
IP="CACACTGACTAAAAGGGTCTGGCTGACTGACT",
Linker=sampleInfo$linkerSequence,
BC=sampleInfo$bcSeq,
Primer=sampleInfo$primer,
LTRBit=sampleInfo$ltrBit)
oligo$R2Start <- with(oligo, 1+nchar(paste0(P7, BC, Spacer, SP2))) #! 1 based
oligo$R1Start <- with(oligo, 1+nchar(paste0(P5, SP1))) #! 1 based
## validated known sites from a prep
##clone1 Clonal 293T cells with known integration at chr1+52699700
##clone2 Clonal 293T cells with known integration at chr17+77440127
##clone3 Clonal 293T cells with known integration at chr19+1330529
##clone4 Clonal 293T cells with known integration at chr1-153461600
##clone7 Clonal 293T cells with known integration at chr1+148889088
##site <- data.frame(chr=c("chr1", "chr17", "chr19", "chr1", "chr1"),
## position=c(52699700, 77440127, 1330529, 153461600, 148889088),
## strand=c("+", "+", "+", "-", "+") )
message("\nGenerate random sites")
site <- get_random_loci(sp=Hsapiens, n=as.integer(args$sites*1.2))
message("\nRemoving sites close to N regions")
checkNbase <- function(site, width=3000) {
seq.plus <- get_sequence_downstream(Hsapiens,
site$chr,
site$position,
"+",
width)
seq.minus <- get_sequence_downstream(Hsapiens,
site$chr,
site$position,
"-",
width)
Nclose <- (grepl('N', seq.plus$seq, ignore.case=TRUE) |
grepl('N', seq.minus$seq, ignore.case=TRUE) )
return( Nclose )
}
isNClose <- checkNbase(site)
site <- site[!isNClose,]
stopifnot(!any(checkNbase(site)))
site <- dplyr::sample_n(site, args$sites, replace=FALSE)
site <- as.data.table(site)
pos <- site[,
data.frame(chr,
position,
strand,
width=sort(get_random_width_gaussian(n=args$sonicLength))),
1:nrow(site)]
message("\nGenerate human sequences for sites")
intseq <- get_sequence_downstream(Hsapiens,
pos$chr,
pos$position,
pos$strand,
pos$width)
intseq <- dplyr::mutate(intseq,
siteid=as.integer(factor(paste0(chr, strand, position))),
sampleid=siteid%%nrow(oligo)+1)
message("\nGenerate machine sequences for sites")
intseq.list <- split(intseq, intseq$sampleid)
I1R1R2qName.list <- bplapply(seq(intseq.list), function(i)
{message(i, "\tof\t", length(intseq.list))
df <- make_miseq_reads(oligo[i,],
intseq.list[[i]],
R1L=args$R1L,
R2L=args$R2L)
return(df) }
,BPPARAM=MulticoreParam(5))
I1R1R2qNamedf <- dplyr::rbind_all(I1R1R2qName.list)
message("\nPlanting base errors with rate ", args$errRate)
I1R1R2qNamedf$I1 <- plant_base_error(I1R1R2qNamedf$I1, args$errRate)
I1R1R2qNamedf$R1 <- plant_base_error(I1R1R2qNamedf$R1, args$errRate)
I1R1R2qNamedf$R2 <- plant_base_error(I1R1R2qNamedf$R2, args$errRate)
## fix qname qid
I1R1R2qNamedf <- (I1R1R2qNamedf %>%
dplyr::mutate(qname=sub(":\\d+$", "", qname),
qname=paste0(qname, ":", 1:n())))
message("\nDump sequences to fastq files")
makeInputFolder(I1R1R2qNamedf, args$outFolder)
message("\nDump truth to bed file")
truth.bed <- (intseq %>%
dplyr::mutate(
breakpoint=ifelse(strand=="+", position+width, position-width),
start=pmin(position, breakpoint),
end=pmax(position, breakpoint),
note="sim",
score=500) %>%
arrange(chr, position, strand, breakpoint) %>%
select(chr, start, end, note, score, strand) )
write.table(truth.bed, file=file.path(args$outFolder, "truth.bed"),
row.names=FALSE, col.names=FALSE, sep="\t", quote=FALSE)
message(file.path(args$outFolder, "truth.bed"))