forked from soedinglab/MMseqs2
-
Notifications
You must be signed in to change notification settings - Fork 0
/
convertalignments.cpp
705 lines (648 loc) · 34.5 KB
/
convertalignments.cpp
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
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
#include "Util.h"
#include "Parameters.h"
#include "Matcher.h"
#include "Debug.h"
#include "DBReader.h"
#include "DBWriter.h"
#include "IndexReader.h"
#include "FileUtil.h"
#include "TranslateNucl.h"
#include "Sequence.h"
#include "Orf.h"
#include "MemoryMapped.h"
#include "NcbiTaxonomy.h"
#include <map>
#ifdef OPENMP
#include <omp.h>
#endif
void printSeqBasedOnAln(std::string &out, const char *seq, unsigned int offset,
const std::string &bt, bool reverse, bool isReverseStrand,
bool translateSequence, const TranslateNucl &translateNucl) {
unsigned int seqPos = 0;
char codon[3];
for (uint32_t i = 0; i < bt.size(); ++i) {
char seqChar = (isReverseStrand == true) ? Orf::complement(seq[offset - seqPos]) : seq[offset + seqPos];
if (translateSequence) {
codon[0] = (isReverseStrand == true) ? Orf::complement(seq[offset - seqPos]) : seq[offset + seqPos];
codon[1] = (isReverseStrand == true) ? Orf::complement(seq[offset - (seqPos+1)]) : seq[offset + (seqPos+1)];
codon[2] = (isReverseStrand == true) ? Orf::complement(seq[offset - (seqPos+2)]) : seq[offset + (seqPos+2)];
seqChar = translateNucl.translateSingleCodon(codon);
}
switch (bt[i]) {
case 'M':
out.append(1, seqChar);
seqPos += (translateSequence) ? 3 : 1;
break;
case 'I':
if (reverse == true) {
out.append(1, '-');
} else {
out.append(1, seqChar);
seqPos += (translateSequence) ? 3 : 1;
}
break;
case 'D':
if (reverse == true) {
out.append(1, seqChar);
seqPos += (translateSequence) ? 3 : 1;
} else {
out.append(1, '-');
}
break;
}
}
}
/*
query Query sequence label
target Target sequenc label
evalue E-value
gapopen Number of gap opens
pident Percentage of identical matches
nident Number of identical matches
qstart 1-based start position of alignment in query sequence
qend 1-based end position of alignment in query sequence
qlen Query sequence length
tstart 1-based start position of alignment in target sequence
tend 1-based end position of alignment in target sequence
tlen Target sequence length
alnlen Number of alignment columns
raw Raw alignment score
bits Bit score
cigar Alignment as string M=letter pair, D=delete (gap in query), I=insert (gap in target)
qseq Full-length query sequence
tseq Full-length target sequence
qheader Header of Query sequence
theader Header of Target sequence
qaln Aligned query sequence with gaps
taln Aligned target sequence with gaps
qframe Query frame (-3 to +3)
tframe Target frame (-3 to +3)
mismatch Number of mismatches
qcov Fraction of query sequence covered by alignment
tcov Fraction of target sequence covered by alignment
qset Query set
tset Target set
*/
std::map<unsigned int, unsigned int> readKeyToSet(const std::string& file) {
std::map<unsigned int, unsigned int> mapping;
if (file.length() == 0) {
return mapping;
}
MemoryMapped lookup(file, MemoryMapped::WholeFile, MemoryMapped::SequentialScan);
char* data = (char *) lookup.getData();
const char* entry[255];
while (*data != '\0') {
const size_t columns = Util::getWordsOfLine(data, entry, 255);
if (columns < 3) {
Debug(Debug::WARNING) << "Not enough columns in lookup file " << file << "\n";
continue;
}
mapping.emplace(Util::fast_atoi<unsigned int>(entry[0]), Util::fast_atoi<unsigned int>(entry[2]));
data = Util::skipLine(data);
}
lookup.close();
return mapping;
}
std::map<unsigned int, std::string> readSetToSource(const std::string& file) {
std::map<unsigned int, std::string> mapping;
if (file.length() == 0) {
return mapping;
}
MemoryMapped source(file, MemoryMapped::WholeFile, MemoryMapped::SequentialScan);
char* data = (char *) source.getData();
const char* entry[255];
while (*data != '\0') {
const size_t columns = Util::getWordsOfLine(data, entry, 255);
if (columns < 2) {
Debug(Debug::WARNING) << "Not enough columns in lookup file " << file << "\n";
continue;
}
data = Util::skipLine(data);
std::string source(entry[1], data - entry[1] - 1);
mapping.emplace(Util::fast_atoi<unsigned int>(entry[0]), source);
}
source.close();
return mapping;
}
static bool compareToFirstInt(const std::pair<unsigned int, unsigned int>& lhs, const std::pair<unsigned int, unsigned int>& rhs){
return (lhs.first <= rhs.first);
}
int convertalignments(int argc, const char **argv, const Command &command) {
Parameters &par = Parameters::getInstance();
par.parseParameters(argc, argv, command, true, 0, 0);
const bool sameDB = par.db1.compare(par.db2) == 0 ? true : false;
const int format = par.formatAlignmentMode;
const bool touch = (par.preloadMode != Parameters::PRELOAD_MODE_MMAP);
bool needSequenceDB = false;
bool needBacktrace = false;
bool needFullHeaders = false;
bool needLookup = false;
bool needSource = false;
bool needTaxonomy = false;
bool needTaxonomyMapping = false;
const std::vector<int> outcodes = Parameters::getOutputFormat(par.outfmt, needSequenceDB, needBacktrace, needFullHeaders,
needLookup, needSource, needTaxonomyMapping, needTaxonomy);
if(format == Parameters::FORMAT_ALIGNMENT_SAM){
needSequenceDB = true;
needBacktrace = true;
}
NcbiTaxonomy * t = NULL;
std::vector<std::pair<unsigned int, unsigned int>> mapping;
if(needTaxonomy){
t = NcbiTaxonomy::openTaxonomy(par.db2);
}
if(needTaxonomy || needTaxonomyMapping){
if(FileUtil::fileExists(std::string(par.db2 + "_mapping").c_str()) == false){
Debug(Debug::ERROR) << par.db2 + "_mapping" << " does not exist. Please create the taxonomy mapping!\n";
EXIT(EXIT_FAILURE);
}
bool isSorted = Util::readMapping( par.db2 + "_mapping", mapping);
if(isSorted == false){
std::stable_sort(mapping.begin(), mapping.end(), compareToFirstInt);
}
}
bool isTranslatedSearch = false;
int dbaccessMode = needSequenceDB ? (DBReader<unsigned int>::USE_INDEX | DBReader<unsigned int>::USE_DATA) : (DBReader<unsigned int>::USE_INDEX);
std::map<unsigned int, unsigned int> qKeyToSet;
std::map<unsigned int, unsigned int> tKeyToSet;
if (needLookup) {
std::string file1 = par.db1 + ".lookup";
std::string file2 = par.db2 + ".lookup";
qKeyToSet = readKeyToSet(file1);
tKeyToSet = readKeyToSet(file2);
}
std::map<unsigned int, std::string> qSetToSource;
std::map<unsigned int, std::string> tSetToSource;
if (needSource) {
std::string file1 = par.db1 + ".source";
std::string file2 = par.db2 + ".source";
qSetToSource = readSetToSource(file1);
tSetToSource = readSetToSource(file2);
}
IndexReader qDbr(par.db1, par.threads, IndexReader::SRC_SEQUENCES, (touch) ? (IndexReader::PRELOAD_INDEX | IndexReader::PRELOAD_DATA) : 0, dbaccessMode);
IndexReader qDbrHeader(par.db1, par.threads, IndexReader::SRC_HEADERS , (touch) ? (IndexReader::PRELOAD_INDEX | IndexReader::PRELOAD_DATA) : 0);
IndexReader *tDbr;
IndexReader *tDbrHeader;
if (sameDB) {
tDbr = &qDbr;
tDbrHeader= &qDbrHeader;
} else {
tDbr = new IndexReader(par.db2, par.threads, IndexReader::SRC_SEQUENCES, (touch) ? (IndexReader::PRELOAD_INDEX | IndexReader::PRELOAD_DATA) : 0, dbaccessMode);
tDbrHeader = new IndexReader(par.db2, par.threads, IndexReader::SRC_HEADERS, (touch) ? (IndexReader::PRELOAD_INDEX | IndexReader::PRELOAD_DATA) : 0);
}
const bool queryNucs = Parameters::isEqualDbtype(qDbr.sequenceReader->getDbtype(), Parameters::DBTYPE_NUCLEOTIDES);
const bool targetNucs = Parameters::isEqualDbtype(tDbr->sequenceReader->getDbtype(), Parameters::DBTYPE_NUCLEOTIDES);
if (needSequenceDB) {
// try to figure out if search was translated. This is can not be solved perfectly.
bool seqtargetAA = false;
if(Parameters::isEqualDbtype(tDbr->getDbtype(), Parameters::DBTYPE_INDEX_DB)){
IndexReader tseqDbr(par.db2, par.threads, IndexReader::SEQUENCES, 0, IndexReader::PRELOAD_INDEX);
seqtargetAA = Parameters::isEqualDbtype(tseqDbr.sequenceReader->getDbtype(), Parameters::DBTYPE_AMINO_ACIDS);
} else if(targetNucs == true && queryNucs == true && par.searchType == Parameters::SEARCH_TYPE_AUTO){
Debug(Debug::WARNING) << "It is unclear from the input if a translated or nucleotide search was performed\n "
"Please provide the parameter --search-type 2 (translated) or 3 (nucleotide)\n";
EXIT(EXIT_FAILURE);
} else if(par.searchType == Parameters::SEARCH_TYPE_TRANSLATED){
seqtargetAA = true;
}
if((targetNucs == true && queryNucs == false ) || (targetNucs == false && queryNucs == true ) || (targetNucs == true && seqtargetAA == true && queryNucs == true ) ){
isTranslatedSearch = true;
}
}
SubstitutionMatrix * subMat= NULL;
if (targetNucs == true && queryNucs == true && isTranslatedSearch == false) {
subMat = new NucleotideMatrix(par.scoringMatrixFile.nucleotides, 1.0, 0.0);
if(par.PARAM_GAP_OPEN.wasSet==false){
par.gapOpen = 5;
}
if(par.PARAM_GAP_EXTEND.wasSet==false){
par.gapExtend = 2;
}
}else{
subMat = new SubstitutionMatrix(par.scoringMatrixFile.aminoacids, 2.0, 0.0);
}
EvalueComputation *evaluer = NULL;
bool queryProfile = false;
bool targetProfile = false;
if (needSequenceDB) {
queryProfile = Parameters::isEqualDbtype(qDbr.sequenceReader->getDbtype(), Parameters::DBTYPE_HMM_PROFILE);
targetProfile = Parameters::isEqualDbtype(tDbr->sequenceReader->getDbtype(), Parameters::DBTYPE_HMM_PROFILE);
evaluer = new EvalueComputation(tDbr->sequenceReader->getAminoAcidDBSize(), subMat, par.gapOpen, par.gapExtend);
}
DBReader<unsigned int> alnDbr(par.db3.c_str(), par.db3Index.c_str(), par.threads, DBReader<unsigned int>::USE_INDEX|DBReader<unsigned int>::USE_DATA);
alnDbr.open(DBReader<unsigned int>::LINEAR_ACCCESS);
unsigned int localThreads = 1;
#ifdef OPENMP
localThreads = std::min((unsigned int)par.threads, (unsigned int)alnDbr.getSize());
#endif
const bool shouldCompress = par.dbOut == true && par.compressed == true;
const int dbType = par.dbOut == true ? Parameters::DBTYPE_GENERIC_DB : Parameters::DBTYPE_OMIT_FILE;
DBWriter resultWriter(par.db4.c_str(), par.db4Index.c_str(), localThreads, shouldCompress, dbType);
resultWriter.open();
const bool isDb = par.dbOut;
TranslateNucl translateNucl(static_cast<TranslateNucl::GenCode>(par.translationTable));
if (format == Parameters::FORMAT_ALIGNMENT_SAM) {
char buffer[1024];
unsigned int lastKey = tDbr->sequenceReader->getLastKey();
bool *headerWritten = new bool[lastKey + 1];
memset(headerWritten, 0, sizeof(bool) * (lastKey + 1));
resultWriter.writeStart(0);
std::string header = "@HD\tVN:1.4\tSO:queryname\n";
resultWriter.writeAdd(header.c_str(), header.size(), 0);
for (size_t i = 0; i < alnDbr.getSize(); i++) {
char *data = alnDbr.getData(i, 0);
while (*data != '\0') {
char dbKeyBuffer[255 + 1];
Util::parseKey(data, dbKeyBuffer);
const unsigned int dbKey = (unsigned int) strtoul(dbKeyBuffer, NULL, 10);
if (headerWritten[dbKey] == false) {
headerWritten[dbKey] = true;
unsigned int tId = tDbr->sequenceReader->getId(dbKey);
unsigned int seqLen = tDbr->sequenceReader->getSeqLen(tId);
unsigned int tHeaderId = tDbrHeader->sequenceReader->getId(dbKey);
const char *tHeader = tDbrHeader->sequenceReader->getData(tHeaderId, 0);
std::string targetId = Util::parseFastaHeader(tHeader);
int count = snprintf(buffer, sizeof(buffer), "@SQ\tSN:%s\tLN:%d\n", targetId.c_str(),
(int32_t) seqLen);
if (count < 0 || static_cast<size_t>(count) >= sizeof(buffer)) {
Debug(Debug::WARNING) << "Truncated line in header " << i << "!\n";
continue;
}
resultWriter.writeAdd(buffer, count, 0);
}
resultWriter.writeEnd(0, 0, false, 0);
data = Util::skipLine(data);
}
}
delete[] headerWritten;
}
Debug::Progress progress(alnDbr.getSize());
#pragma omp parallel num_threads(localThreads)
{
unsigned int thread_idx = 0;
#ifdef OPENMP
thread_idx = static_cast<unsigned int>(omp_get_thread_num());
#endif
char buffer[1024];
std::string result;
result.reserve(1024*1024);
std::string queryProfData;
queryProfData.reserve(1024);
std::string queryBuffer;
queryProfData.reserve(1024);
std::string queryHeaderBuffer;
queryHeaderBuffer.reserve(1024);
std::string targetProfData;
targetProfData.reserve(1024);
std::string newBacktrace;
newBacktrace.reserve(1024);
const TaxonNode * taxonNode = NULL;
#pragma omp for schedule(dynamic, 10)
for (size_t i = 0; i < alnDbr.getSize(); i++) {
progress.updateProgress();
const unsigned int queryKey = alnDbr.getDbKey(i);
char *querySeqData = NULL;
queryProfData.clear();
if (needSequenceDB) {
size_t qId = qDbr.sequenceReader->getId(queryKey);
querySeqData = qDbr.sequenceReader->getData(qId, thread_idx);
if(sameDB && qDbr.sequenceReader->isCompressed()){
size_t querySeqDataLen = qDbr.sequenceReader->getSeqLen(qId);
queryBuffer.assign(querySeqData, querySeqDataLen);
querySeqData = (char*) queryBuffer.c_str();
}
if (queryProfile) {
Sequence::extractProfileConsensus(querySeqData, *subMat, queryProfData);
}
}
size_t qHeaderId = qDbrHeader.sequenceReader->getId(queryKey);
const char *qHeader = qDbrHeader.sequenceReader->getData(qHeaderId, thread_idx);
size_t qHeaderLen = qDbrHeader.sequenceReader->getSeqLen(qHeaderId);
std::string queryId = Util::parseFastaHeader(qHeader);
if (sameDB && needFullHeaders) {
queryHeaderBuffer.assign(qHeader, qHeaderLen);
qHeader = (char*) queryHeaderBuffer.c_str();
}
char *data = alnDbr.getData(i, thread_idx);
while (*data != '\0') {
Matcher::result_t res = Matcher::parseAlignmentRecord(data, true);
data = Util::skipLine(data);
if (res.backtrace.empty() && needBacktrace == true) {
Debug(Debug::ERROR) << "Backtrace cigar is missing in the alignment result. Please recompute the alignment with the -a flag.\n"
"Command: mmseqs align " << par.db1 << " " << par.db2 << " " << par.db3 << " " << "alnNew -a\n";
EXIT(EXIT_FAILURE);
}
size_t tHeaderId = tDbrHeader->sequenceReader->getId(res.dbKey);
const char *tHeader = tDbrHeader->sequenceReader->getData(tHeaderId, thread_idx);
size_t tHeaderLen = tDbrHeader->sequenceReader->getSeqLen(tHeaderId);
std::string targetId = Util::parseFastaHeader(tHeader);
unsigned int gapOpenCount = 0;
unsigned int alnLen = res.alnLength;
unsigned int missMatchCount = 0;
unsigned int identical = 0;
if (res.backtrace.empty() == false) {
size_t matchCount = 0;
alnLen = 0;
for (size_t pos = 0; pos < res.backtrace.size(); pos++) {
int cnt = 0;
if (isdigit(res.backtrace[pos])) {
cnt += Util::fast_atoi<int>(res.backtrace.c_str() + pos);
while (isdigit(res.backtrace[pos])) {
pos++;
}
}
alnLen += cnt;
switch (res.backtrace[pos]) {
case 'M':
matchCount += cnt;
break;
case 'D':
case 'I':
gapOpenCount += 1;
break;
}
}
// res.seqId = X / alnLen;
identical = static_cast<unsigned int>(res.seqId * static_cast<float>(alnLen) + 0.5);
//res.alnLength = alnLen;
missMatchCount = static_cast<unsigned int>( matchCount - identical);
} else {
const int adjustQstart = (res.qStartPos == -1) ? 0 : res.qStartPos;
const int adjustDBstart = (res.dbStartPos == -1) ? 0 : res.dbStartPos;
const float bestMatchEstimate = static_cast<float>(std::min(abs(res.qEndPos - adjustQstart), abs(res.dbEndPos - adjustDBstart)));
missMatchCount = static_cast<unsigned int>(bestMatchEstimate * (1.0f - res.seqId) + 0.5);
}
switch (format) {
case Parameters::FORMAT_ALIGNMENT_BLAST_TAB: {
if (outcodes.empty()) {
int count = snprintf(buffer, sizeof(buffer),
"%s\t%s\t%1.3f\t%d\t%d\t%d\t%d\t%d\t%d\t%d\t%.2E\t%d\n",
queryId.c_str(), targetId.c_str(), res.seqId, alnLen,
missMatchCount, gapOpenCount,
res.qStartPos + 1, res.qEndPos + 1,
res.dbStartPos + 1, res.dbEndPos + 1,
res.eval, res.score);
if (count < 0 || static_cast<size_t>(count) >= sizeof(buffer)) {
Debug(Debug::WARNING) << "Truncated line in entry" << i << "!\n";
continue;
}
result.append(buffer, count);
} else {
char *targetSeqData = NULL;
targetProfData.clear();
unsigned int taxon = 0;
if(needTaxonomy || needTaxonomyMapping) {
std::pair<unsigned int, unsigned int> val;
val.first = res.dbKey;
std::vector<std::pair<unsigned int, unsigned int>>::iterator mappingIt;
mappingIt = std::upper_bound(mapping.begin(), mapping.end(), val, compareToFirstInt);
if (mappingIt == mapping.end() || mappingIt->first != val.first) {
taxon = 0;
}else{
taxon = mappingIt->second;
if(needTaxonomy){
taxonNode = t->taxonNode(taxon, false);
}
}
}
if (needSequenceDB) {
size_t tId = tDbr->sequenceReader->getId(res.dbKey);
targetSeqData = tDbr->sequenceReader->getData(tId, thread_idx);
if (targetProfile) {
Sequence::extractProfileConsensus(targetSeqData, *subMat, targetProfData);
}
}
for(size_t i = 0; i < outcodes.size(); i++) {
switch (outcodes[i]) {
case Parameters::OUTFMT_QUERY:
result.append(queryId);
break;
case Parameters::OUTFMT_TARGET:
result.append(targetId);
break;
case Parameters::OUTFMT_EVALUE:
result.append(SSTR(res.eval));
break;
case Parameters::OUTFMT_GAPOPEN:
result.append(SSTR(gapOpenCount));
break;
case Parameters::OUTFMT_PIDENT:
result.append(SSTR(res.seqId));
break;
case Parameters::OUTFMT_NIDENT:
result.append(SSTR(identical));
break;
case Parameters::OUTFMT_QSTART:
result.append(SSTR(res.qStartPos + 1));
break;
case Parameters::OUTFMT_QEND:
result.append(SSTR(res.qEndPos + 1));
break;
case Parameters::OUTFMT_QLEN:
result.append(SSTR(res.qLen));
break;
case Parameters::OUTFMT_TSTART:
result.append(SSTR(res.dbStartPos + 1));
break;
case Parameters::OUTFMT_TEND:
result.append(SSTR(res.dbEndPos + 1));
break;
case Parameters::OUTFMT_TLEN:
result.append(SSTR(res.dbLen));
break;
case Parameters::OUTFMT_ALNLEN:
result.append(SSTR(alnLen));
break;
case Parameters::OUTFMT_RAW:
result.append(SSTR(static_cast<int>(evaluer->computeRawScoreFromBitScore(res.score) + 0.5)));
break;
case Parameters::OUTFMT_BITS:
result.append(SSTR(res.score));
break;
case Parameters::OUTFMT_CIGAR:
result.append(SSTR(res.backtrace));
newBacktrace.clear();
break;
case Parameters::OUTFMT_QSEQ:
if (queryProfile) {
result.append(queryProfData.c_str(), res.qLen);
} else {
result.append(querySeqData, res.qLen);
}
break;
case Parameters::OUTFMT_TSEQ:
if (targetProfile) {
result.append(targetProfData.c_str(), res.dbLen);
} else {
result.append(targetSeqData, res.dbLen);
}
break;
case Parameters::OUTFMT_QHEADER:
result.append(qHeader, qHeaderLen);
break;
case Parameters::OUTFMT_THEADER:
result.append(tHeader, tHeaderLen);
break;
case Parameters::OUTFMT_QALN:
if (queryProfile) {
printSeqBasedOnAln(result, queryProfData.c_str(), res.qStartPos,
Matcher::uncompressAlignment(res.backtrace), false, (res.qStartPos > res.qEndPos),
(isTranslatedSearch == true && queryNucs == true), translateNucl);
} else {
printSeqBasedOnAln(result, querySeqData, res.qStartPos,
Matcher::uncompressAlignment(res.backtrace), false, (res.qStartPos > res.qEndPos),
(isTranslatedSearch == true && queryNucs == true), translateNucl);
}
break;
case Parameters::OUTFMT_TALN: {
if (targetProfile) {
printSeqBasedOnAln(result, targetProfData.c_str(), res.dbStartPos,
Matcher::uncompressAlignment(res.backtrace), true,
(res.dbStartPos > res.dbEndPos),
(isTranslatedSearch == true && targetNucs == true), translateNucl);
} else {
printSeqBasedOnAln(result, targetSeqData, res.dbStartPos,
Matcher::uncompressAlignment(res.backtrace), true,
(res.dbStartPos > res.dbEndPos),
(isTranslatedSearch == true && targetNucs == true), translateNucl);
}
break;
}
case Parameters::OUTFMT_MISMATCH:
result.append(SSTR(missMatchCount));
break;
case Parameters::OUTFMT_QCOV:
result.append(SSTR(res.qcov));
break;
case Parameters::OUTFMT_TCOV:
result.append(SSTR(res.dbcov));
break;
case Parameters::OUTFMT_QSET:
result.append(SSTR(qSetToSource[qKeyToSet[queryKey]]));
break;
case Parameters::OUTFMT_QSETID:
result.append(SSTR(qKeyToSet[queryKey]));
break;
case Parameters::OUTFMT_TSET:
result.append(SSTR(tSetToSource[tKeyToSet[res.dbKey]]));
break;
case Parameters::OUTFMT_TSETID:
result.append(SSTR(tKeyToSet[res.dbKey]));
break;
case Parameters::OUTFMT_TAXID:
result.append(SSTR(taxon));
break;
case Parameters::OUTFMT_TAXNAME:
result.append((taxonNode != NULL) ? taxonNode->name : "unclassified");
break;
case Parameters::OUTFMT_TAXLIN:
result.append((taxonNode != NULL) ? t->taxLineage(taxonNode) : "unclassified");
break;
case Parameters::OUTFMT_EMPTY:
result.push_back('-');
break;
}
if (i < outcodes.size() - 1) {
result.push_back('\t');
}
}
result.push_back('\n');
}
break;
}
case Parameters::FORMAT_ALIGNMENT_BLAST_WITH_LEN: {
int count = snprintf(buffer, sizeof(buffer),
"%s\t%s\t%1.3f\t%d\t%d\t%d\t%d\t%d\t%d\t%d\t%.2E\t%d\t%d\t%d\n",
queryId.c_str(), targetId.c_str(), res.seqId, alnLen,
missMatchCount, gapOpenCount,
res.qStartPos + 1, res.qEndPos + 1,
res.dbStartPos + 1, res.dbEndPos + 1,
res.eval, res.score,
res.qLen, res.dbLen);
if (count < 0 || static_cast<size_t>(count) >= sizeof(buffer)) {
Debug(Debug::WARNING) << "Truncated line in entry" << i << "!\n";
continue;
}
result.append(buffer, count);
break;
}
case Parameters::FORMAT_ALIGNMENT_SAM: {
bool strand = res.qEndPos > res.qStartPos;
int rawScore = static_cast<int>(evaluer->computeRawScoreFromBitScore(res.score) + 0.5);
uint32_t mapq = -4.343 * log(exp(static_cast<double>(-rawScore)));
mapq = (uint32_t) (mapq + 4.99);
mapq = mapq < 254 ? mapq : 254;
int count = snprintf(buffer, sizeof(buffer), "%s\t%d\t%s\t%d\t%d\t", queryId.c_str(), (strand) ? 16: 0, targetId.c_str(), res.dbStartPos + 1, mapq);
if (count < 0 || static_cast<size_t>(count) >= sizeof(buffer)) {
Debug(Debug::WARNING) << "Truncated line in entry" << i << "!\n";
continue;
}
result.append(buffer, count);
result.append(res.backtrace);
result.append("\t*\t0\t0\t");
int start = std::min(res.qStartPos, res.qEndPos);
int end = std::max(res.qStartPos, res.qEndPos);
if (queryProfile) {
result.append(queryProfData.c_str() + start, (end + 1) - start);
} else {
result.append(querySeqData + start, (end + 1) - start);
}
count = snprintf(buffer, sizeof(buffer), "\t*\tAS:i:%d\tNM:i:%d\n", rawScore, missMatchCount);
if (count < 0 || static_cast<size_t>(count) >= sizeof(buffer)) {
Debug(Debug::WARNING) << "Truncated line in entry" << i << "!\n";
continue;
}
result.append(buffer, count);
break;
}
// case Parameters::FORMAT_ALIGNMENT_GFF:{
// // for TBLASTX
// bool strand = res.qEndPos > res.qStartPos;
// int currStart = std::min(res.qStartPos, res.qEndPos);
// int currEnd = std::max(res.qStartPos, res.qEndPos);
// int currLen = currEnd - currStart;
// result.append(queryId);
// result.append("\tconserve\tprotein_match\t");
// result.append(SSTR(currStart+1));
// result.push_back('\t');
// result.append(SSTR(currEnd+1));
// result.push_back('\t');
// result.append(SSTR(currLen));
// result.push_back('\t');
// result.push_back((strand) ? '-' : '+');
// result.append("\t.\t");
// result.append("ID=");
// result.append(queryId);
// result.append(":hsp:");
// result.append(SSTR(counter));
// result.append(";");
// break;
// }
default:
Debug(Debug::ERROR) << "Not implemented yet";
EXIT(EXIT_FAILURE);
}
}
resultWriter.writeData(result.c_str(), result.size(), queryKey, thread_idx, isDb);
result.clear();
}
}
// tsv output
resultWriter.close(true);
if (isDb == false) {
FileUtil::remove(par.db4Index.c_str());
}
if(needTaxonomy){
delete t;
}
alnDbr.close();
if (sameDB == false) {
delete tDbr;
delete tDbrHeader;
}
if (needSequenceDB) {
delete evaluer;
}
delete subMat;
return EXIT_SUCCESS;
}