-
Notifications
You must be signed in to change notification settings - Fork 1
/
blaster.py
676 lines (566 loc) · 23.9 KB
/
blaster.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
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
#!/usr/bin/python
# -*- coding: utf-8 -*-
'''BLASTing batch script.
Dependencies:
- Python 2.7.1
- Biopython 1.56
- BLAST+ 2.2.25
Configuration (Linux):
- Regular packages for Python and Biopython.
- Add BLAST+ commands to path putting this code to .bashrc:
PATH=$PATH:/home/nelas/Downloads/ncbi-blast-2.2.25+/bin
export PATH
See README file for details or execute the command for usage and arguments:
python blaster.py -h
'''
import getopt
import logging
import os
import sys
from datetime import datetime
from shutil import move
from Bio import Entrez, SeqIO
from Bio.Seq import Seq
from Bio.SeqRecord import SeqRecord
from Bio.Blast import NCBIWWW, NCBIXML
from Bio.Blast.Applications import NcbiblastnCommandline, NcbiblastpCommandline, NcbiblastxCommandline, NcbitblastnCommandline, NcbitblastxCommandline
class Sequence(object):
'''Represents a known sequence.
filepath: path to fasta file.
ref: sequence reference number.
description: long description.
sequence: the sequence itself.
gene_id: NCBI gene identification.
gene_name: NCBI gene name.
blast_output: path to .xml file.
loci: list of loci with recurrency?.
'''
EVALUE_THRESH = 0.001
# Cache folder in user home.
cache_folder = os.path.join(os.environ['HOME'], '.blaster/seqs')
# Check if cache folder exists.
try:
os.makedirs(cache_folder)
except:
pass
#TODO add __str__ to classes.
def __init__(self, filepath=None, ref=None, database=None):
self.limit = 3
self.loci = {}
self.gene_name = ''
self.database = database
if filepath:
# Candidate gene sequences are instantiated with a filepath.
self.filepath = filepath
self.parse_genbank(self.filepath)
# Fallback to candidate gene filename, if needed.
if not self.gene_name:
self.gene_name = os.path.splitext(os.path.basename(self.filepath))[0]
else:
# Resulting genes from reverse blasts are queried by ref.
cache_filename = ref + '.gb'
cache_path = os.path.join(self.cache_folder, cache_filename)
try:
# Try to get cached file.
self.parse_genbank(cache_path)
except:
# Fetch data, create cache, and parse data.
handle_string = efetcher(ref)
cache_file = open(cache_path, 'wb')
cache_file.write(handle_string)
cache_file.close()
self.parse_genbank(cache_path)
def parse_genbank(self, filepath):
'''Parse data from GENBANK file.'''
record = SeqIO.read(filepath, 'genbank')
# Define class attributes.
self.description = record.description
self.ref = record.id
self.sequence = str(record.seq)
self.gene_name, self.gene_id = self.get_gene_name_id(record)
self.organism = record.description.split('[')[1].split(']')[0]
def parse_fasta(self, filepath):
'''Parse data from FASTA file.'''
record = SeqIO.read(filepath, 'fasta')
# Define attributes.
self.description = record.description.split('|')[-1]
self.ref = record.id.split('|')[3]
self.sequence = str(record.seq)
self.gene_name = filepath.split('/')[-1][:-3]
self.organism = record.description.split('[')[1].split(']')[0]
def get_gene_name_id(self, record):
'''Return the gene name and id from a SeqRecord.'''
# Extract CDS feature. Check if CDS exists?
cds = [feature for feature in record.features if feature.type == 'CDS'][0]
# Get gene name.
try:
gene_name = cds.qualifiers.get('gene')[0]
gene_name = parse_filename(gene_name)
except:
logger.critical('No gene name for %s', record.id)
gene_name = None
# Find gene id.
try:
gene_db_xref = cds.qualifiers.get('db_xref')
for xref in gene_db_xref:
if xref.startswith('GeneID'):
gene_id = xref.split(':')[-1]
break
else:
gene_id = '[empty]'
logger.critical('No gene id for %s', record.id)
except:
gene_id = '[empty]'
logger.critical('No gene id for %s', record.id)
return gene_name, gene_id
def set_reciprocal_db(self):
'''Set the reciprocal database according to candidate gene organism.'''
# Locally define organism.
organism = self.organism
base_path = '/home/nelas/Biologia/Doutorado/databases/'
reciprocals = {
'Homo sapiens': 'human_protein.fa',
'Drosophila melanogaster': 'drosophila.fa',
'Danio rerio': 'zebrafish.fa',
'Mus musculus': 'mouse.fa',
'Strongylocentrotus purpuratus': 'urchin.fa',
'Ciona intestinalis': 'ciona.fa',
'Prostheceraeus vittatus': 'prostheceraeus.fa',
}
try:
self.reciprocal_db = os.path.join(base_path, reciprocals[organism])
except:
logger.warning('Database missing for %s!', organism)
self.reciprocal_db = None
def parse_blast(self):
'''Parse BLAST output file.'''
blast_records = NCBIXML.parse(open(self.blast_output))
n = 0
# Iterate over BLAST results from database.
for blast_record in blast_records:
for alignment in blast_record.alignments:
for hsp in alignment.hsps:
if n < self.limit:
if hsp.expect < self.EVALUE_THRESH:
print 'Instantiating >> %s' % alignment.title
locus_id = alignment.title.split()[0]
# Save locus to dictionary if key does not exists.
# Only first hsp found for locus is added, other are skipped.
if not locus_id in self.loci.keys():
self.loci[locus_id] = {
'score': hsp.score,
'evalue': hsp.expect,
'frame': hsp.frame[1], # hit frame
}
n += 1
class Locus(object):
'''Represents a new sequence.
candidate: instance of Candidate.
description: long description.
sequence: the sequence itself.
score: from Candidate BLAST.
evalue: from Candidate BLAST.
reverse_blast_output: reverse BLAST to human protein db.
'''
EVALUE_THRESH = 0.001
# Reverse BLAST folder.
reverse_folder = 'reverse'
def __init__(self, id, frame, candidate, database):
# Linkback to candidates.
self.candidates = {}
# Store BLASTs by reciprocal databases (ie, species name).
self.reciprocal_blasts = {}
self.update_candidates(candidate)
self.reciprocal = False
self.rank = 0
self.id = id
if frame > 0:
self.frame = '+%d' % frame
else:
self.frame = str(frame)
self.database = database
# Create filename before the rest.
self.filename = parse_filename(self.id)
self.reverse_results_folder = os.path.join(self.reverse_folder, 'results')
# Check if reverse folder exists.
if not os.path.isdir(self.reverse_folder):
os.mkdir(self.reverse_folder)
# Check if reverse results folder exists.
if not os.path.isdir(self.reverse_results_folder):
os.mkdir(self.reverse_results_folder)
self.set_filepath()
self.set_sequence()
self.write_fasta()
def update_candidates(self, candidate):
'''Update list of linked candidate genes.'''
if not candidate.gene_name in self.candidates.keys():
self.candidates[candidate.gene_name] = candidate
def set_sequence(self):
'''Get and set sequence from database.'''
parsed = SeqIO.parse(self.database, 'fasta')
for locus in parsed:
if locus.id == self.id:
self.sequence = str(locus.seq)
try:
dummy_seq = self.sequence
except:
print 'SEQUENCE NOT SET! Check the IDs.'
sys.exit(2)
def set_blast_output(self, organism):
'''Build and return filepath for reverse blast output with organism name.'''
if not organism in self.reciprocal_blasts.keys():
o = organism.replace(' ', '_')
blast_output = os.path.join(
self.reverse_results_folder, '%s-%s.xml' % (self.filename, o)
)
self.reciprocal_blasts[organism] = {
'path': blast_output,
'sequences': []
}
return self.reciprocal_blasts[organism]['path']
def set_filepath(self):
'''Build filepath from folder and filename.'''
locus_filepath = os.path.join(self.reverse_folder, '%s.fa' % self.filename)
self.filepath = locus_filepath
def write_fasta(self):
'''Write FASTA file to filepath.'''
SeqIO.write(SeqRecord(Seq(self.sequence), id=self.id, description=''), self.filepath, 'fasta')
def parse_blast(self, organism):
'''Parse reciprocal BLAST.'''
print 'Parsing reciprocal BLAST: %s' % self.filename
reverse_blast_output = self.reciprocal_blasts[organism]['path']
blast_records = NCBIXML.parse(open(reverse_blast_output))
for blast_record in blast_records:
for alignment in blast_record.alignments:
for hsp in alignment.hsps:
if hsp.expect < self.EVALUE_THRESH:
try:
ref = alignment.title.split('|')[3]
except:
ref = alignment.title
sequence = Sequence(ref=ref)
sequence.evalue = hsp.expect
sequence.score = hsp.score
self.reciprocal_blasts[organism]['sequences'].append(sequence)
def parse_filename(string):
'''Parse string to a valid filepath.'''
# Remove spaces and slashes to avoid filepath issues.
parsed = string.replace(' ', '_').replace('/', '-')
return parsed
def local_blast(blast_type, arguments):
'''Execute BLAST command.'''
# Instantiate the BLAST command.
if blast_type == 'blastn':
cline = NcbiblastnCommandline(**arguments)
elif blast_type == 'blastp':
cline = NcbiblastpCommandline(**arguments)
elif blast_type == 'blastx':
cline = NcbiblastxCommandline(**arguments)
elif blast_type == 'tblastn':
cline = NcbitblastnCommandline(**arguments)
elif blast_type == 'tblastx':
cline = NcbitblastxCommandline(**arguments)
# Execute BLAST.
stdout, stderr = cline()
print '%s BLASTed!' % arguments['query']
def prepare(candidates, candidates_folder):
'''Check candidate gene files (convert to FASTA, if needed).'''
for gene in candidates:
gene_filepath = os.path.join(candidates_folder, gene)
gene_gb = os.path.splitext(gene_filepath)[0] + '.gb'
if gene.endswith(('.gb', '.genbank')):
continue
elif gene.endswith(('.fa', '.txt')):
# Only fetch if file does not exist.
try:
open(gene_gb, 'r')
pass
except IOError:
logger.debug('Found FASTA: %s; converting to GENBANK...', gene_filepath)
# 1. Parse ref.
record = SeqIO.read(gene_filepath, 'fasta')
#XXX Candidate genes id should be formatted as NCBI (with reference number).
gene_ref = record.id.split('|')[3]
# 2. Fetch entry with ref and save handle to string.
handle_string = efetcher(ref=gene_ref)
# 3. Write genbank file.
if handle_string:
f = open(gene_gb, 'w')
f.write(handle_string)
f.close()
else:
logger.debug('File type not supported: %s', gene)
def efetcher(ref, db='protein', rettype='gp', retmode='txt'):
'''Container function for Entrez.efetch.
Tries to access Entrez, retry in case of failure.
Return handle string ready to be written into file.
'''
logger.debug('Initiating Entrez via efetch...')
logger.debug('ref=%s, db=%s, rettype=%s, retmode=%s', ref, db, rettype, retmode)
# Trying to contact Entrez.
try:
handle = Entrez.efetch(db=db, id=ref, rettype=rettype, retmode=retmode)
handle_string = handle.read()
logger.debug('Fetched %s!', ref)
return handle_string
except:
logger.debug('Entrez connection failed. Retrying...')
handle = Entrez.efetch(db=db, id=ref, rettype=rettype, retmode=retmode)
handle_string = handle.read()
logger.debug('Fetched %s!', ref)
return handle_string
else:
logger.critical('Entrez fetching failed for %s, giving up.', ref)
return None
def backup(filepath):
'''Store previous results file.'''
# Check if results.txt exists.
# If yes, rename it with timestamp and mv to hidden.
raw_timestamp = os.path.getmtime(filepath)
timestamp = datetime.fromtimestamp(raw_timestamp)
timestamp_string = timestamp.strftime('%Y-%m-%d_%Hh%Mm%Ss')
new_name = '.results_%s.txt' % timestamp_string
move(filepath, new_name)
return new_name
def usage():
'''Explanation for arguments.'''
print
print 'USAGE: ./blaster.py -d Membranipora.fasta -b tblastn'
print
print ' -c, --candidates \n\tFolder with `.fa` or `.gb` files of candidate genes. One gene per file.'
print ' -d, --database \n\tLocal database with new data (eg, transcriptome).'
print ' -b, --blast \n\tBLAST command (blastn, blastp, blastx, tblastn, tblastx).'
print ' -e, --email \n\tYour email is required for Entrez.'
print
#TODO Limit number of loci from candidate blast results.
#TODO Specify threshold evalue.
#TODO Limit the number of genes looked up during reciprocal blast.
def main(argv):
# Available BLASTs.
blast_types = ['blastn', 'blastp', 'blastx', 'tblastn', 'tblastx']
# Database name.
database = None
# Default BLAST command.
blast_type = 'tblastn'
# Folder with candidate-genes.
candidates_folder = 'candidates'
# Folder with candidate-genes' BLASTs.
candidates_results_folder = 'blasts'
# Default email.
Entrez.email = 'your@email.com'
# Parse arguments.
try:
opts, args = getopt.getopt(argv, 'hc:d:b:e:', [
'help',
'candidates=',
'database=',
'blast=',
'email=',
])
except getopt.GetoptError:
usage()
logger.critical('Argument error: "%s". Aborting...',
' '.join(argv))
sys.exit(2)
# Argument handler.
for opt, arg in opts:
if opt in ('-h', '--help'):
usage()
sys.exit()
elif opt in ('-c', '--candidates'):
candidates_folder = arg
elif opt in ('-d', '--database'):
database = os.path.abspath(arg)
elif opt in ('-b', '--blast'):
blast_type = arg
elif opt in ('-e', '--email'):
Entrez.email = arg
# Print summary of arguments.
logger.debug('Arguments: candidates=%s, database=%s, blast=%s, email=%s', candidates_folder, database, blast_type, Entrez.email)
# Check if BLAST command was specified.
if not blast_type:
logger.critical('BLAST command was not specified (use "-b"). Aborting...')
sys.exit(2)
else:
if not blast_type in blast_types:
logger.critical('Unknown BLAST command: %s', blast_type)
logger.info('Available BLAST commands: %s', ', '.join(blast_types))
sys.exit(2)
# Get candidate genes.
#TODO Make it recursively search sub-directories.
try:
candidates = os.listdir(candidates_folder)
except OSError:
logger.critical('Folder "%s" does not exist! Aborting...', candidates_folder)
sys.exit(2)
# Issue error if there are no candidate genes.
if not candidates:
logger.critical('No candidate genes at %s folder! Aborting...',
candidates_folder)
sys.exit(2)
# Check if results folder exists.
if not os.path.isdir(candidates_results_folder):
os.mkdir(candidates_results_folder)
# Process candidate genes.
prepare(candidates, candidates_folder)
# Get proper genes, now.
candidates = os.listdir(candidates_folder)
# Only read GENBANK files.
candidates = [file for file in candidates if file.endswith(('.gb', '.genbank'))]
# Print info before starting.
logger.info('%d genes to be BLASTed against %s database!',
len(candidates), database)
# Main object to store genes and loci.
genes = {}
loci = {}
# BLAST each gene against local database.
for genefile in candidates:
print
logger.info('BLASTing %s against %s...', genefile, database)
# Instantiate Sequence.
gene_filepath = os.path.join(candidates_folder, genefile)
candidate = Sequence(filepath=gene_filepath, database=database)
# Define variables.
output_filepath = os.path.join(candidates_results_folder, '%s.xml' % candidate.gene_name)
candidate.blast_output = output_filepath
# Only BLAST if needed.
try:
blastfile = open(candidate.blast_output)
blastfile.close()
except:
arguments = {
'query': candidate.filepath,
'db': database,
'out': candidate.blast_output,
'outfmt': 5, # Export in XML
}
# Execute BLAST command.
local_blast(blast_type, arguments)
logger.info('Parsing XML...')
print '\nCandidate >> Gene: %s, ID:%s, Ref: %s, Description: %s' % (
candidate.gene_name, candidate.gene_id,
candidate.ref, candidate.description
)
candidate.parse_blast()
# Work on loci.
for locus_id in candidate.loci.keys():
# Instantiate Locus object.
if locus_id in loci.keys():
locus = loci[locus_id]
# Update candidate list.
locus.update_candidates(candidate)
else:
frame = candidate.loci[locus_id]['frame']
locus = Locus(locus_id, frame, candidate, database)
# Generate file for BLAST output.
reverse_blast_output = locus.set_blast_output(candidate.organism)
# Reciprocal BLAST querying locus sequence against human database.
try:
blastfile = open(reverse_blast_output)
blastfile.close()
except:
# 1. Check available databases using set_reciprocal_db.
candidate.set_reciprocal_db()
if candidate.reciprocal_db:
# Sanitize filepaths for shell (ie handling "|" in locus IDs).
reverse_args = {
'query': '"%s"' % locus.filepath,
'db': candidate.reciprocal_db,
'out': '"%s"' % reverse_blast_output,
'outfmt': 5,
}
#XXX What to do when recirocal database is not protein?
if blast_type == 'tblastn':
local_blast('blastx', reverse_args)
elif blast_type == 'blastp':
local_blast('blastp', reverse_args)
else:
# 2. If None, prepare for BLAST online based on the organism.
#XXX Handle blastp as above, also.
print 'BLASTing over NCBI... (may take a while).'
logger.info('BLASTing %s over NCBI (restricted to %s).',
locus.filepath, candidate.organism)
handle = NCBIWWW.qblast('blastx', 'refseq_protein', open(locus.filepath).read(), entrez_query='%s[Organism]' % candidate.organism)
reverse_file = open(reverse_blast_output, 'w')
reverse_file.write(handle.read())
reverse_file.close()
handle.close()
# Parse reverse BLAST output.
locus.parse_blast(candidate.organism)
# Add locus to main list.
loci[locus_id] = locus
# Add to main dictionary.
genes[candidate.gene_name] = candidate
# Print loci equivalent to candidate genes.
logger.info('Creating result files...')
fasta_loci = []
results_filepath = 'results.txt'
try:
previous = open(results_filepath)
backed_up = backup(results_filepath)
if backed_up:
logger.info('Old results backed up: %s', backed_up)
else:
logger.warning('No backup for %s', results_filepath)
except:
pass
output_file = open(results_filepath, 'w')
#TODO Print date and variables used to results.txt.
for locus_id, locus in loci.iteritems():
locus_description = '| frame: %s | candidates: %s' % (locus.frame, ', '.join(locus.candidates.keys()))
fasta_loci.append(SeqRecord(Seq(locus.sequence), id=locus_id, description=locus_description))
output_file.write('%s %s\n\n' % (locus_id, locus_description))
output_file.write('\t\t\t\t{0:20}{1:20}{2:20}{3:1}\n'.format('gene', 'id', 'accession', 'e-value'))
#output_file.write('\t\t\t\tgene\t\t\tid\t\taccession\t\te-value\n')
for organism, reciprocal_blast in locus.reciprocal_blasts.iteritems():
output_file.write(' {0:>2}\n'.format(organism))
sequences = reciprocal_blast['sequences']
# Instantiate first sequence to get gene id.
try:
first = sequences[0]
except:
output_file.write('\t\t\t\t{0}'.format('No reciprocal sequences found.\n'))
continue
current_gene_id = first.gene_id
n_genes = 0
for sequence in sequences:
# Limits to 5 most similar genes.
if sequence.gene_id == current_gene_id or n_genes < 5:
output_file.write('\t\t\t\t{0:20}{1:20}{2:20}{3:<.2e}'.format(sequence.gene_name, sequence.gene_id, sequence.ref, sequence.evalue))
accession_numbers = [seq.ref for seq in locus.candidates.values()]
if sequence.ref in accession_numbers:
output_file.write(' {0:>}'.format('<<'))
output_file.write('\n')
else:
break
current_gene_id = sequence.gene_id
n_genes += 1
output_file.write('\n\n')
output_file.close()
wrote_fasta = SeqIO.write(fasta_loci, 'results.fa', 'fasta')
logger.info('%d records written to results.fa', wrote_fasta)
logger.info('Done, bye!')
if __name__ == '__main__':
# Create logger.
logger = logging.getLogger('blaster')
logger.setLevel(logging.DEBUG)
logger.propagate = False
# Define message format.
formatter = logging.Formatter('[%(levelname)s] %(asctime)s @ %(module)s %(funcName)s (l%(lineno)d): %(message)s')
# Create logger console handler.
console_handler = logging.StreamHandler()
console_handler.setLevel(logging.DEBUG)
# Format console output.
console_handler.setFormatter(formatter)
# Add console to logger.
logger.addHandler(console_handler)
# Create logger file handler.
file_handler = logging.FileHandler('log_blaster.log')
file_handler.setLevel(logging.DEBUG)
# Format file output.
file_handler.setFormatter(formatter)
# Add file to logger.
logger.addHandler(file_handler)
# Initiates main function.
logger.info('BLASTer is running...')
main(sys.argv[1:])