-
Notifications
You must be signed in to change notification settings - Fork 20
/
_demux.py
483 lines (402 loc) · 19.6 KB
/
_demux.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
# ----------------------------------------------------------------------------
# Copyright (c) 2016-2020, QIIME 2 development team.
#
# Distributed under the terms of the Modified BSD License.
#
# The full license is in the file LICENSE, distributed with this software.
# ----------------------------------------------------------------------------
import gzip
import yaml
import itertools
import collections
import collections.abc
import random
import resource
import skbio
import psutil
import qiime2
from q2_types.per_sample_sequences import (
SingleLanePerSampleSingleEndFastqDirFmt,
SingleLanePerSamplePairedEndFastqDirFmt,
FastqManifestFormat, YamlFormat)
from ._ecc import GolayDecoder
from ._format import ErrorCorrectionDetailsFmt
FastqHeader = collections.namedtuple('FastqHeader', ['id', 'description'])
class ECDetails:
COLUMNS = ['id',
'sample',
'barcode-sequence-id',
'barcode-uncorrected',
'barcode-corrected',
'barcode-errors']
def __init__(self, fmt):
self._fp = open(str(fmt), 'w')
self._write_header()
def write(self, parts):
self._fp.write('\t'.join([str(part) for part in parts]))
self._fp.write('\n')
def _write_header(self):
self.write(self.COLUMNS)
def _read_fastq_seqs(filepath):
# This function is adapted from @jairideout's SO post:
# http://stackoverflow.com/a/39302117/3424666
fh = gzip.open(filepath, 'rt')
for seq_header, seq, qual_header, qual in itertools.zip_longest(*[fh] * 4):
yield (seq_header.strip(), seq.strip(), qual_header.strip(),
qual.strip())
def _trim_id(id):
return id.rsplit('/', 1)[0]
def _trim_description(desc):
# The first number of ':' seperated description is the read number
if ':' in desc:
desc = desc.split(':', 1)[1]
return desc.rsplit('/', 1)[0]
def _record_to_fastq_header(record):
tokens = record[0][1:].split(' ', maxsplit=1)
if len(tokens) == 1:
id, = tokens
description = None
else:
id, description = tokens
return FastqHeader(id=id, description=description)
# This is global so that it can be tested without changing the actual ulimits.
# NOTE: UNIX only
OPEN_FH_LIMIT, _ = resource.getrlimit(resource.RLIMIT_NOFILE)
def _maintain_open_fh_count(per_sample_fastqs, paired=False):
files_to_open = 1 if not paired else 2
# NOTE: UNIX only
if psutil.Process().num_fds() + files_to_open < OPEN_FH_LIMIT:
return
# currently open per-sample files
if not paired:
open_fhs = [fh for fh in per_sample_fastqs.values()
if not fh.closed]
else:
open_fhs = [fh for fh in per_sample_fastqs.values()
if not fh[0].closed]
# If the number of open files reaches the allotted resources limit,
# close around 15% of the open files. 15% was chosen because if you
# only close a single file it will start to have to do it on every loop
# and on a 35 file benchmark using a hard coded limit of 10 files,
# only closing one file added an increased runtime of 160%
n_to_close = round(0.15 * len(open_fhs))
if paired:
n_to_close //= 2
# Never close more than files than are open, also if closing,
# close at least the number of files that will need to be opened.
n_to_close = min(len(open_fhs), max(n_to_close, files_to_open))
for rand_fh in random.sample(open_fhs, n_to_close):
if paired:
fwd, rev = rand_fh
fwd.close()
rev.close()
else:
rand_fh.close()
class BarcodeSequenceFastqIterator(collections.abc.Iterable):
def __init__(self, barcode_generator, sequence_generator,
ignore_description_mismatch=False):
self.barcode_generator = barcode_generator
self.sequence_generator = sequence_generator
self.ignore_description_mismatch = ignore_description_mismatch
def __iter__(self):
# Adapted from q2-types
for barcode_record, sequence_record in itertools.zip_longest(
self.barcode_generator, self.sequence_generator):
if barcode_record is None:
raise ValueError('More sequences were provided than barcodes.')
if sequence_record is None:
raise ValueError('More barcodes were provided than sequences.')
# The id or description fields may end with "/read-number", which
# will differ between the sequence and barcode reads. Confirm that
# they are identical up until the last /
barcode_header = _record_to_fastq_header(barcode_record)
sequence_header = _record_to_fastq_header(sequence_record)
# confirm that the id fields are equal
if _trim_id(barcode_header.id) != \
_trim_id(sequence_header.id):
raise ValueError(
'Mismatched sequence ids: %s and %s' %
(_trim_id(barcode_header.id),
_trim_id(sequence_header.id)))
if not self.ignore_description_mismatch:
# if a description field is present, confirm that they're equal
if barcode_header.description is None and \
sequence_header.description is None:
pass
elif barcode_header.description is None:
raise ValueError(
'Barcode header lines do not contain description '
'fields but sequence header lines do.')
elif sequence_header.description is None:
raise ValueError(
'Sequence header lines do not contain description '
'fields but barcode header lines do.')
elif _trim_description(barcode_header.description) != \
_trim_description(sequence_header.description):
raise ValueError(
'Mismatched sequence descriptions: %s and %s' %
(_trim_description(barcode_header.description),
_trim_description(sequence_header.description)))
yield barcode_record, sequence_record
class BarcodePairedSequenceFastqIterator(collections.abc.Iterable):
def __init__(self, barcode_generator, forward_generator,
reverse_generator, ignore_description_mismatch=False):
self.barcode_generator = barcode_generator
self.forward_generator = forward_generator
self.reverse_generator = reverse_generator
self.ignore_description_mismatch = ignore_description_mismatch
def __iter__(self):
# Adapted from q2-types
for barcode_record, forward_record, reverse_record \
in itertools.zip_longest(self.barcode_generator,
self.forward_generator,
self.reverse_generator):
if barcode_record is None:
raise ValueError('More sequences were provided than barcodes.')
if forward_record is None:
raise ValueError('More barcodes were provided than '
'forward-sequences.')
elif reverse_record is None:
raise ValueError('More barcodes were provided than '
'reverse-sequences.')
# The id or description fields may end with "/read-number", which
# will differ between the sequence and barcode reads. Confirm that
# they are identical up until the last /
barcode_header = _record_to_fastq_header(barcode_record)
forward_header = _record_to_fastq_header(forward_record)
reverse_header = _record_to_fastq_header(reverse_record)
# confirm that the id fields are equal
if not (_trim_id(barcode_header.id) ==
_trim_id(forward_header.id) ==
_trim_id(reverse_header.id)):
raise ValueError(
'Mismatched sequence ids: %s, %s, and %s' %
(_trim_id(barcode_header.id),
_trim_id(forward_header.id),
_trim_id(reverse_header.id)))
if not self.ignore_description_mismatch:
# if a description field is present, confirm that they're equal
if barcode_header.description is None and \
forward_header.description is None and \
reverse_header.description is None:
pass
elif barcode_header.description is None:
raise ValueError(
'Barcode header lines do not contain description '
'fields but sequence header lines do.')
elif forward_header.description is None:
raise ValueError(
'Forward-read header lines do not contain description '
'fields but barcode header lines do.')
elif reverse_header.description is None:
raise ValueError(
'Reverse-read header lines do not contain description '
'fields but barcode header lines do.')
elif not (_trim_description(barcode_header.description) ==
_trim_description(forward_header.description) ==
_trim_description(reverse_header.description)):
raise ValueError(
'Mismatched sequence descriptions: %s, %s, and %s' %
(_trim_description(barcode_header.description),
_trim_description(forward_header.description),
_trim_description(reverse_header.description)))
yield barcode_record, forward_record, reverse_record
def _make_barcode_map(barcodes, rev_comp_mapping_barcodes):
barcode_map = {}
barcode_len = None
for sample_id, barcode in barcodes.to_series().iteritems():
if barcode_len is None:
barcode_len = len(barcode)
elif len(barcode) != barcode_len:
raise ValueError('Barcodes of different lengths were detected: '
'%d != %d. Variable length barcodes are not '
'supported.' % (len(barcode), barcode_len))
if rev_comp_mapping_barcodes:
barcode = str(skbio.DNA(barcode).reverse_complement())
if barcode in barcode_map:
raise ValueError('A duplicate barcode was detected. The barcode '
'%s was observed for samples %s and %s.'
% (barcode, sample_id, barcode_map[barcode]))
barcode_map[barcode] = sample_id
return barcode_map, barcode_len
def _write_metadata_yaml(dir_fmt):
metadata = YamlFormat()
metadata.path.write_text(yaml.dump({'phred-offset': 33}))
dir_fmt.metadata.write_data(metadata, YamlFormat)
def emp_single(seqs: BarcodeSequenceFastqIterator,
barcodes: qiime2.CategoricalMetadataColumn,
golay_error_correction: bool = True,
rev_comp_barcodes: bool = False,
rev_comp_mapping_barcodes: bool = False,
ignore_description_mismatch: bool = False
) -> (SingleLanePerSampleSingleEndFastqDirFmt,
ErrorCorrectionDetailsFmt):
seqs.ignore_description_mismatch = ignore_description_mismatch
result = SingleLanePerSampleSingleEndFastqDirFmt()
barcode_map, barcode_len = _make_barcode_map(
barcodes, rev_comp_mapping_barcodes)
if golay_error_correction:
decoder = GolayDecoder()
manifest = FastqManifestFormat()
manifest_fh = manifest.open()
manifest_fh.write('sample-id,filename,direction\n')
manifest_fh.write('# direction is not meaningful in this file as these\n')
manifest_fh.write('# data may be derived from forward, reverse, or \n')
manifest_fh.write('# joined reads\n')
per_sample_fastqs = {}
ec_details_fmt = ErrorCorrectionDetailsFmt()
ec_details = ECDetails(ec_details_fmt)
for i, (barcode_record, sequence_record) in enumerate(seqs, start=1):
barcode_read = barcode_record[1]
if rev_comp_barcodes:
barcode_read = str(skbio.DNA(barcode_read).reverse_complement())
raw_barcode_read = barcode_read[:barcode_len]
if golay_error_correction:
# A three bit filter is implicitly used by the decoder. See Hamady
# and Knight 2009 Genome Research for the justification:
#
# https://genome.cshlp.org/content/19/7/1141.full
#
# Specifically that "...Golay codes of 12 bases can correct all
# triple-bit errors and detect all quadruple-bit errors."
barcode_read, ecc_errors = decoder.decode(raw_barcode_read)
golay_stats = [barcode_read, ecc_errors]
else:
barcode_read = raw_barcode_read
golay_stats = [None, None]
sample_id = barcode_map.get(barcode_read)
record = [
f'record-{i}',
sample_id,
barcode_record[0],
raw_barcode_read,
]
ec_details.write(record + golay_stats)
if sample_id is None:
continue
if sample_id not in per_sample_fastqs:
# The barcode id, lane number and read number are not relevant
# here. We might ultimately want to use a dir format other than
# SingleLanePerSampleSingleEndFastqDirFmt which doesn't care
# about this information. Similarly, the direction of the read
# isn't relevant here anymore.
barcode_id = len(per_sample_fastqs) + 1
path = result.sequences.path_maker(sample_id=sample_id,
barcode_id=barcode_id,
lane_number=1,
read_number=1)
_maintain_open_fh_count(per_sample_fastqs)
per_sample_fastqs[sample_id] = gzip.open(str(path), mode='a')
manifest_fh.write('%s,%s,%s\n' % (sample_id, path.name, 'forward'))
if per_sample_fastqs[sample_id].closed:
_maintain_open_fh_count(per_sample_fastqs)
per_sample_fastqs[sample_id] = gzip.open(
per_sample_fastqs[sample_id].name, mode='a')
fastq_lines = '\n'.join(sequence_record) + '\n'
fastq_lines = fastq_lines.encode('utf-8')
per_sample_fastqs[sample_id].write(fastq_lines)
if len(per_sample_fastqs) == 0:
raise ValueError('No sequences were mapped to samples. Check that '
'your barcodes are in the correct orientation (see '
'the rev_comp_barcodes and/or '
'rev_comp_mapping_barcodes options). If barcodes are '
'NOT Golay format set golay_error_correction '
'to False.')
for fh in per_sample_fastqs.values():
fh.close()
manifest_fh.close()
result.manifest.write_data(manifest, FastqManifestFormat)
_write_metadata_yaml(result)
return result, ec_details_fmt
def emp_paired(seqs: BarcodePairedSequenceFastqIterator,
barcodes: qiime2.CategoricalMetadataColumn,
golay_error_correction: bool = True,
rev_comp_barcodes: bool = False,
rev_comp_mapping_barcodes: bool = False,
ignore_description_mismatch: bool = False
) -> (SingleLanePerSamplePairedEndFastqDirFmt,
ErrorCorrectionDetailsFmt):
seqs.ignore_description_mismatch = ignore_description_mismatch
result = SingleLanePerSamplePairedEndFastqDirFmt()
barcode_map, barcode_len = _make_barcode_map(
barcodes, rev_comp_mapping_barcodes)
if golay_error_correction:
decoder = GolayDecoder()
manifest = FastqManifestFormat()
manifest_fh = manifest.open()
manifest_fh.write('sample-id,filename,direction\n')
per_sample_fastqs = {}
ec_details_fmt = ErrorCorrectionDetailsFmt()
ec_details = ECDetails(ec_details_fmt)
for i, record in enumerate(seqs, start=1):
barcode_record, forward_record, reverse_record = record
barcode_read = barcode_record[1]
if rev_comp_barcodes:
barcode_read = str(skbio.DNA(barcode_read).reverse_complement())
raw_barcode_read = barcode_read[:barcode_len]
if golay_error_correction:
# A three bit filter is implicitly used by the decoder. See Hamady
# and Knight 2009 Genome Research for the justification:
#
# https://genome.cshlp.org/content/19/7/1141.full
#
# Specifically that "...Golay codes of 12 bases can correct all
# triple-bit errors and detect all quadruple-bit errors."
barcode_read, ecc_errors = decoder.decode(raw_barcode_read)
golay_stats = [barcode_read, ecc_errors]
else:
barcode_read = raw_barcode_read
golay_stats = [None, None]
sample_id = barcode_map.get(barcode_read)
record = [
f'record-{i}',
sample_id,
barcode_record[0],
raw_barcode_read,
]
ec_details.write(record + golay_stats)
if sample_id is None:
continue
if sample_id not in per_sample_fastqs:
barcode_id = len(per_sample_fastqs) + 1
fwd_path = result.sequences.path_maker(sample_id=sample_id,
barcode_id=barcode_id,
lane_number=1,
read_number=1)
rev_path = result.sequences.path_maker(sample_id=sample_id,
barcode_id=barcode_id,
lane_number=1,
read_number=2)
_maintain_open_fh_count(per_sample_fastqs, paired=True)
per_sample_fastqs[sample_id] = (
gzip.open(str(fwd_path), mode='a'),
gzip.open(str(rev_path), mode='a')
)
manifest_fh.write('%s,%s,%s\n' % (sample_id, fwd_path.name,
'forward'))
manifest_fh.write('%s,%s,%s\n' % (sample_id, rev_path.name,
'reverse'))
if per_sample_fastqs[sample_id][0].closed:
_maintain_open_fh_count(per_sample_fastqs, paired=True)
fwd, rev = per_sample_fastqs[sample_id]
per_sample_fastqs[sample_id] = (
gzip.open(fwd.name, mode='a'),
gzip.open(rev.name, mode='a')
)
fwd, rev = per_sample_fastqs[sample_id]
fwd.write(('\n'.join(forward_record) + '\n').encode('utf-8'))
rev.write(('\n'.join(reverse_record) + '\n').encode('utf-8'))
if len(per_sample_fastqs) == 0:
raise ValueError('No sequences were mapped to samples. Check that '
'your barcodes are in the correct orientation (see '
'the rev_comp_barcodes and/or '
'rev_comp_mapping_barcodes options). If barcodes are '
'NOT Golay format set golay_error_correction '
'to False.')
for fwd, rev in per_sample_fastqs.values():
fwd.close()
rev.close()
manifest_fh.close()
result.manifest.write_data(manifest, FastqManifestFormat)
_write_metadata_yaml(result)
return result, ec_details_fmt