File: options.hpp

package info (click to toggle)
lambda-align 1.0.3-6
  • links: PTS, VCS
  • area: main
  • in suites: bookworm, bullseye, trixie
  • size: 788 kB
  • sloc: cpp: 4,653; sh: 70; makefile: 27
file content (1451 lines) | stat: -rw-r--r-- 57,645 bytes parent folder | download | duplicates (3)
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
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
// ==========================================================================
//                                  lambda
// ==========================================================================
// Copyright (c) 2013-2017, Hannes Hauswedell <h2 @ fsfe.org>
// Copyright (c) 2016-2017, Knut Reinert and Freie Universität Berlin
// All rights reserved.
//
// This file is part of Lambda.
//
// Lambda is Free Software: you can redistribute it and/or modify it
// under the terms found in the LICENSE[.md|.rst] file distributed
// together with this file.
//
// Lambda is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
//
// ==========================================================================
// options.h: contains the options and argument parser
// ==========================================================================


#ifndef SEQAN_LAMBDA_OPTIONS_H_
#define SEQAN_LAMBDA_OPTIONS_H_

#include <cstdio>
#include <unistd.h>
#include <bitset>

#include <seqan/basic.h>
#include <seqan/translation.h>
#include <seqan/arg_parse.h>
#include <seqan/index.h>
#include <seqan/blast.h>

// ==========================================================================
// Forwards
// ==========================================================================

template <typename T>
struct SamBamExtraTags;

// ==========================================================================
// Metafunctions
// ==========================================================================

// SIZE TYPES

// Expected Number of Sequences
template <typename TAlph>
using SizeTypeNum_ = uint32_t;

// Expected Lengths of Sequences
template <typename T>
struct SizeTypePosMeta_
{
    using Type = uint16_t;
};

template <>
struct SizeTypePosMeta_<Dna5>
{
    // DNA sequences are expected to be longer
    using Type = uint32_t;
};

template <typename TAlph>
using SizeTypePos_ = typename SizeTypePosMeta_<TAlph>::Type;


// suffix array overloads
namespace seqan
{

template<typename TSpec1, typename TSpec2, typename TSpec3>
struct SAValue<StringSet<String<ReducedAminoAcid<TSpec1>, TSpec2>, TSpec3> >
{
    typedef Pair<SizeTypeNum_<TSpec1>, SizeTypePos_<TSpec1>, Pack> Type;
};

template<typename TSpec1, typename TSpec2, typename TSpec3, typename TFunctor>
struct SAValue<StringSet<ModifiedString<String<TSpec1, TSpec2>, TFunctor>, TSpec3> >
{
    typedef Pair<SizeTypeNum_<TSpec1>, SizeTypePos_<TSpec1>, Pack> Type;
};

template<typename TSpec1, typename TSpec2, typename TSpec3>
struct SAValue<StringSet<String<TSpec1, TSpec2>, TSpec3> >
{
    typedef Pair<SizeTypeNum_<TSpec1>, SizeTypePos_<TSpec1>, Pack> Type;
};

template <typename TString, typename TSpec>
struct DefaultIndexStringSpec<StringSet<TString, TSpec>>
{
#if !defined(LAMBDA_INDEXER) && defined(LAMBDA_MMAPPED_DB)
    using Type    = MMap<>;
#else
    using Type    = Alloc<>;
#endif
};

// our custom Bam Overload
template <typename TDirection, typename TStorageSpec>
struct FormattedFileContext<FormattedFile<Bam, TDirection, BlastTabular>, TStorageSpec>
{
    typedef StringSet<Segment<String<char, MMap<> >, InfixSegment> >   TNameStore;
    typedef NameStoreCache<TNameStore>                                  TNameStoreCache;
    typedef BamIOContext<TNameStore, TNameStoreCache, TStorageSpec>     Type;
};

}

using namespace seqan;

// Index Specs
struct LambdaFMIndexConfig
{
    using LengthSum = size_t;
#if !defined(LAMBDA_INDEXER) && defined(LAMBDA_MMAPPED_DB)
    using TAlloc    = MMap<>;
#else
    using TAlloc    = Alloc<>;
#endif
    using Bwt       = WaveletTree<void, WTRDConfig<LengthSum, TAlloc> >;
    using Sentinels = Levels<void, LevelsRDConfig<LengthSum, TAlloc> >;

    static const unsigned SAMPLING = 10;
};

template <typename TSpec = void>
using TFMIndex = FMIndex<TSpec, LambdaFMIndexConfig>;

// lazy...
template <typename TString>
using TCDStringSet = StringSet<TString, Owner<ConcatDirect<> > >;

template <BlastProgram p>
using OrigQryAlph = typename std::conditional<
                                           (p == BlastProgram::BLASTN) ||
                                           (p == BlastProgram::BLASTX) ||
                                           (p == BlastProgram::TBLASTX),
                                           Dna5,
                                           AminoAcid>::type;

template <BlastProgram p>
using OrigSubjAlph = typename std::conditional<
                                           (p == BlastProgram::BLASTN) ||
                                           (p == BlastProgram::TBLASTN) ||
                                           (p == BlastProgram::TBLASTX),
                                           Dna5,
                                           AminoAcid>::type;

template <BlastProgram p>
using TransAlph = typename std::conditional<(p == BlastProgram::BLASTN),
                                            Dna5,
                                            AminoAcid>::type;


template <BlastProgram p, typename TRedAlph_>
using RedAlph = typename std::conditional<(p == BlastProgram::BLASTN),
                                          Dna5,
                                          TRedAlph_>::type;

template <typename TString>
void getCwd(TString & string)
{
    char cwd[1000];

#ifdef PLATFORM_WINDOWS
    _getcwd(cwd, 1000);
#else
    getcwd(cwd, 1000);
#endif

    assign(string, cwd);
}

template <typename TString, typename TValue>
bool setEnv(TString const & key, TValue & value)
{
#ifdef PLATFORM_WINDOWS
    return !_putenv_s(toCString(key), toCString(value));
#else
    return !setenv(toCString(key), toCString(value), true);
#endif
}

// ==========================================================================
// Classes
// ==========================================================================

// --------------------------------------------------------------------------
// Class LambdaOptions
// --------------------------------------------------------------------------

// This struct stores the options from the command line.

struct SharedOptions
{
    // Verbosity level.  0 -- quiet, 1 -- normal, 2 -- verbose, 3 -- very verbose.
    int verbosity = 1;

    std::string commandLine;

    std::string dbFile;

    int      dbIndexType = 0;
    // for indexer, the file format of database sequences
    // for main app, the file format of query sequences
    // 0 -- fasta, 1 -- fastq
//     int      fileFormat = 0;

    int      alphReduction = 0;

    GeneticCodeSpec geneticCode = CANONICAL;

    BlastProgram blastProgram = BlastProgram::BLASTX;

    bool        isTerm = true;
    unsigned    terminalCols = 80;

    unsigned        threads     = 1;

    SharedOptions()
    {
        isTerm = isTerminal();
        if (isTerm)
        {
            unsigned _rows;
            getTerminalSize(terminalCols, _rows);
        }
    }
};


struct LambdaOptions : public SharedOptions
{

    std::string     queryFile;
    bool            revComp     = true;

    int             outFileFormat; // 0 = BLAST, 1 = SAM, 2 = BAM
    std::string     output;
    std::vector<BlastMatchField<>::Enum> columns;
    std::string     outputBam;
    std::bitset<64> samBamTags;
    bool            samWithRefHeader;
    unsigned        samBamSeq;
    bool            samBamHardClip;
    bool            versionInformationToOutputFile;

    unsigned        queryPart = 0;

//     bool            semiGlobal;

    bool            doubleIndexing = true;

    unsigned        seedLength  = 0;
    unsigned        maxSeedDist = 1;
    bool            hammingOnly = true;

    int             seedGravity     = 0;
    unsigned        seedOffset      = 0;
    unsigned        minSeedLength   = 0;

//     unsigned int    minSeedEVal     = 0;
//     double          minSeedBitS     = -1;

    // 0 = manual, positive X = blosumX, negative Y = pamY
    int             scoringMethod   = 62;
    // scores
    int             gapOpen         = -11;
    int             gapExtend       = -1;
    int             match           = 0; // only for manual
    int             misMatch        = 0; // only for manual

    int             xDropOff    = 0;
    int             band        = -1;
    double          eCutOff     = 0;
    int             idCutOff    = 0;
    unsigned long   maxMatches  = 500;

    bool            filterPutativeDuplicates = true;
    bool            filterPutativeAbundant = true;

    int             preScoring = 0; // 0 = off, 1 = seed, 2 = region (
    double          preScoringThresh    = 0.0;

    LambdaOptions() :
        SharedOptions()
    {
    }
};

struct LambdaIndexerOptions : public SharedOptions
{
    std::string     segFile = "";
    std::string     algo = "";

    bool            truncateIDs;

    LambdaIndexerOptions()
        : SharedOptions()
    {}
};

// ==========================================================================
// Functions
// ==========================================================================

// --------------------------------------------------------------------------
// Function displayCopyright()
// --------------------------------------------------------------------------

void
sharedSetup(ArgumentParser & parser)
{
    // Set short description, version, and date.
    std::string versionString = std::string(SEQAN_APP_VERSION) + " (Git commit " +
                                std::string(SEQAN_REVISION) + ")";
    setVersion(parser, versionString);
    setDate(parser, __DATE__);
    setShortCopyright(parser, "2013-2017 Hannes Hauswedell, released under the GNU GPL v3 (or later); "
                              "2016-2017 Knut Reinert and Freie Universität Berlin, released under the 3-clause-BSDL");

    setCitation(parser, "Hauswedell et al (2014); doi: 10.1093/bioinformatics/btu439");

    setLongCopyright(parser,
        " Copyright (c) 2013-2017, Hannes Hauswedell\n"
        " All rights reserved.\n"
        "\n"
        " Lambda is free software: you can redistribute it and/or modify\n"
        " it under the terms of the GNU General Public License as published by\n"
        " the Free Software Foundation, either version 3 of the License, or\n"
        " (at your option) any later version.\n"
        "\n"
        " Lambda is distributed in the hope that it will be useful,\n"
        " but WITHOUT ANY WARRANTY; without even the implied warranty of\n"
        " MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the\n"
        " GNU General Public License for more details.\n"
        "\n"
        " You should have received a copy of the GNU General Public License\n"
        " along with Lambda.  If not, see <http://www.gnu.org/licenses/>.\n"
        "\n"
        " Copyright (c) 2016-2017 Knut Reinert and Freie Universität Berlin\n"
        " All rights reserved.\n"
        "\n"
        " Redistribution and use in source and binary forms, with or without\n"
        " modification, are permitted provided that the following conditions are met:\n"
        "\n"
        " * Redistributions of source code must retain the above copyright\n"
        "   notice, this list of conditions and the following disclaimer.\n"
        " * Redistributions in binary form must reproduce the above copyright\n"
        "   notice, this list of conditions and the following disclaimer in the\n"
        "   documentation and/or other materials provided with the distribution.\n"
        " * Neither the name of Knut Reinert or the FU Berlin nor the names of\n"
        "   its contributors may be used to endorse or promote products derived\n"
        "   from this software without specific prior written permission.\n"
        "\n"
        " THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS \"AS IS\"\n"
        " AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE\n"
        " IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE\n"
        " ARE DISCLAIMED. IN NO EVENT SHALL KNUT REINERT OR THE FU BERLIN BE LIABLE\n"
        " FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL\n"
        " DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR\n"
        " SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER\n"
        " CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT\n"
        " LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY\n"
        " OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH\n"
        " DAMAGE.\n");

    addOption(parser, ArgParseOption("v", "verbosity",
        "Display more/less diagnostic output during operation: 0 [only errors]; 1 [default]; 2 "
        "[+run-time, options and statistics].",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "verbosity", "1");
    setMinValue(parser, "verbosity", "0");
    setMaxValue(parser, "verbosity", "2");

    addDescription(parser, "Lambda is a local aligner optimized for many query "
        "sequences and searches in protein space. It is compatible to BLAST, but "
        "much faster than BLAST and many other comparable tools.");

    addDescription(parser, "Detailed information is available in the wiki: "
        "<https://github.com/seqan/lambda/wiki>");
}

// --------------------------------------------------------------------------
// Function parseCommandLine()
// --------------------------------------------------------------------------

// SHARED
ArgumentParser::ParseResult
parseCommandLineShared(SharedOptions & options, ArgumentParser & parser);

ArgumentParser::ParseResult
parseCommandLine(LambdaOptions & options, int argc, char const ** argv)
{
    // save commandLine
    for (int i = 0; i < argc; ++i)
        options.commandLine += std::string(argv[i]) + " ";
    eraseBack(options.commandLine);

    // Setup ArgumentParser.
    ArgumentParser parser("lambda");
    // Set short description, version, and date.
    setShortDescription(parser, "the Local Aligner for Massive Biological "
    "DatA");

    // Define usage line and long description.
    addUsageLine(parser, "[\\fIOPTIONS\\fP] \\fI-q QUERY.fasta\\fP "
                         "\\fI-d DATABASE.fasta\\fP "
                         "[\\fI-o output.m8\\fP]");

    sharedSetup(parser);

    addSection(parser, "Input Options");
    addOption(parser, ArgParseOption("q", "query",
        "Query sequences.",
        ArgParseArgument::INPUT_FILE,
        "IN"));
    setValidValues(parser, "query", getFileExtensions(SeqFileIn()));
    setRequired(parser, "q");

    addOption(parser, ArgParseOption("d", "database",
        "Path to original database sequences (a precomputed index with .sa or .fm needs to exist!).",
        ArgParseArgument::INPUT_FILE,
        "IN"));
    setValidValues(parser, "database", getFileExtensions(SeqFileIn()));
    setRequired(parser, "d");

    addOption(parser, ArgParseOption("di", "db-index-type",
        "database index is in this format.",
//         "(auto means \"try sa first then fm\").",
        ArgParseArgument::STRING,
        "STR"));
    setValidValues(parser, "db-index-type", "sa fm");
    setDefaultValue(parser, "db-index-type", "fm");
    setAdvanced(parser, "db-index-type");

    addSection(parser, "Output Options");
    addOption(parser, ArgParseOption("o", "output",
        "File to hold reports on hits (.m* are blastall -m* formats; .m8 is tab-seperated, .m9 is tab-seperated with "
        "with comments, .m0 is pairwise format).",
        ArgParseArgument::OUTPUT_FILE,
        "OUT"));
    auto exts = getFileExtensions(BlastTabularFileOut<>());
    append(exts, getFileExtensions(BlastReportFileOut<>()));
    append(exts, getFileExtensions(BamFileOut()));
    CharString extsConcat;
    // remove .sam.bam, .sam.vcf.gz, .sam.tbi
    for (auto const & ext : exts)
    {
        if ((!endsWith(ext, ".bam") || startsWith(ext, ".bam")) &&
            (!endsWith(ext, ".vcf.gz")) &&
            (!endsWith(ext, ".sam.tbi")))
        {
            append(extsConcat, ext);
            appendValue(extsConcat, ' ');
        }
    }
    setValidValues(parser, "output", toCString(extsConcat));
    setDefaultValue(parser, "output", "output.m8");

    addOption(parser, ArgParseOption("oc", "output-columns",
        "Print specified column combination and/or order (.m8 and .m9 outputs only); call -oc help for more details.",
        ArgParseArgument::STRING,
        "STR"));
    setDefaultValue(parser, "output-columns", "std");
    setAdvanced(parser, "output-columns");

    addOption(parser, ArgParseOption("id", "percent-identity",
        "Output only matches above this threshold (checked before e-value "
        "check).",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "percent-identity", "0");
    setMinValue(parser, "percent-identity", "0");
    setMaxValue(parser, "percent-identity", "100");

    addOption(parser, ArgParseOption("e", "e-value",
        "Output only matches that score below this threshold.",
        ArgParseArgument::DOUBLE));
    setDefaultValue(parser, "e-value", "0.1");
    setMinValue(parser, "e-value", "0");

    addOption(parser, ArgParseOption("nm", "num-matches",
        "Print at most this number of matches per query.",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "num-matches", "500");
    setMinValue(parser, "num-matches", "1");

    addOption(parser, ArgParseOption("", "sam-with-refheader",
        "BAM files require all subject names to be written to the header. For SAM this is not required, so Lambda does "
        "not automatically do it to save space (especially for protein database this is a lot!). If you still want "
        "them with SAM, e.g. for better BAM compatibility, use this option.",
        ArgParseArgument::STRING,
        "STR"));
    setValidValues(parser, "sam-with-refheader", "on off");
    setDefaultValue(parser, "sam-with-refheader", "off");
    setAdvanced(parser, "sam-with-refheader");

    addOption(parser, ArgParseOption("", "sam-bam-seq",
        "Write matching DNA subsequence into SAM/BAM file (BLASTN). For BLASTX and TBLASTX the matching protein "
        "sequence is \"untranslated\" and positions retransformed to the original sequence. For BLASTP and TBLASTN "
        "there is no DNA sequence so a \"*\" is written to the SEQ column. The matching protein sequence can be "
        "written as an optional tag, see --sam-bam-tags. If set to uniq than "
        "the sequence is omitted iff it is identical to the previous match's subsequence.",
        ArgParseArgument::STRING,
        "STR"));
    setValidValues(parser, "sam-bam-seq", "always uniq never");
    setDefaultValue(parser, "sam-bam-seq", "uniq");
    setAdvanced(parser, "sam-bam-seq");

    addOption(parser, ArgParseOption("", "sam-bam-tags",
        "Write the specified optional columns to the SAM/BAM file. Call --sam-bam-tags help for more details.",
        ArgParseArgument::STRING,
        "STR"));
    setDefaultValue(parser, "sam-bam-tags", "AS NM ZE ZI ZF");
    setAdvanced(parser, "sam-bam-tags");

    addOption(parser, ArgParseOption("", "sam-bam-clip",
        "Whether to hard-clip or soft-clip the regions beyond the local match. Soft-clipping retains the full sequence "
        "in the output file, but obviously uses more space.",
        ArgParseArgument::STRING,
        "STR"));
    setValidValues(parser, "sam-bam-clip", "hard soft");
    setDefaultValue(parser, "sam-bam-clip", "hard");
    setAdvanced(parser, "sam-bam-clip");

    addOption(parser, ArgParseOption("", "version-to-outputfile",
        "Write the Lambda program tag and version number to the output file.",
        ArgParseArgument::STRING,
        "STR"));
    setValidValues(parser, "version-to-outputfile", "on off");
    setDefaultValue(parser, "version-to-outputfile", "on");
    hideOption(parser, "version-to-outputfile");

    addSection(parser, "General Options");
#ifdef _OPENMP
    addOption(parser, ArgParseOption("t", "threads",
        "number of threads to run concurrently.",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "threads", omp_get_max_threads());
#else
    addOption(parser, ArgParseOption("t", "threads",
        "LAMBDA BUILT WITHOUT OPENMP; setting this option has no effect.",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "threads", 1);
#endif
    setAdvanced(parser, "threads");

    addOption(parser, ArgParseOption("qi", "query-index-type",
        "controls double-indexing.",
        ArgParseArgument::STRING));
    setValidValues(parser, "query-index-type", "radix none");
    setDefaultValue(parser, "query-index-type", "none");
    setAdvanced(parser, "query-index-type");

    addOption(parser, ArgParseOption("qp", "query-partitions",
        "Divide the query into qp number of blocks before processing; should be"
        " a multiple of the number of threads, defaults to one per thread. "
        "Only used with double-indexing; strong influence on memory, see below.",
        ArgParseArgument::INTEGER));
#ifdef _OPENMP
    setDefaultValue(parser, "query-partitions", omp_get_max_threads());
#else
    setDefaultValue(parser, "query-partitions", 1);
#endif
    hideOption(parser, "query-partitions"); // HIDDEN

    addSection(parser, "Alphabets and Translation");
    addOption(parser, ArgParseOption("p", "program",
        "Blast Operation Mode.",
        ArgParseArgument::STRING,
        "STR"));
#ifdef FASTBUILD
    setValidValues(parser, "program", "blastp blastx");
#else
    setValidValues(parser, "program", "blastn blastp blastx tblastn tblastx");
#endif
    setDefaultValue(parser, "program", "blastx");

//     addOption(parser, ArgParseOption("qa", "query-alphabet",
//         "original alphabet of the query sequences",
//         ArgParseArgument::STRING,
//         "STR"));
//     setValidValues(parser, "query-alphabet", "dna5 aminoacid");
//     setDefaultValue(parser, "query-alphabet", "dna5");
//
//     addOption(parser, ArgParseOption("da", "db-alphabet",
//         "original alphabet of the subject sequences",
//         ArgParseArgument::STRING,
//         "STR"));
//     setValidValues(parser, "db-alphabet", "dna5 aminoacid");
//     setDefaultValue(parser, "db-alphabet", "aminoacid");
//
//     addOption(parser, ArgParseOption("sa", "seeding-alphabet",
//         "alphabet to use during seeding (reduction possible)",
//         ArgParseArgument::STRING,
//         "STR"));
//     setValidValues(parser, "seeding-alphabet", "dna5 aminoacid");
//     setDefaultValue(parser, "seeding-alphabet", "murphy10");

    addOption(parser, ArgParseOption("g", "genetic-code",
        "The translation table to use for nucl -> amino acid translation"
        "(not for BlastN, BlastP). See "
        "https://www.ncbi.nlm.nih.gov/Taxonomy/Utils/wprintgc.cgi?mode=c"
        " for ids (default is generic). Six frames are generated.",
        ArgParseArgument::INTEGER));
//     setValidValues(parser, "alph", "0 10");
    setDefaultValue(parser, "genetic-code", "1");
    setAdvanced(parser, "genetic-code");

    addOption(parser, ArgParseOption("ar", "alphabet-reduction",
        "Alphabet Reduction for seeding phase (ignored for BLASTN).",
        ArgParseArgument::STRING,
        "STR"));
    setValidValues(parser, "alphabet-reduction", "none murphy10");
    setDefaultValue(parser, "alphabet-reduction", "murphy10");
    setAdvanced(parser, "alphabet-reduction");

    addSection(parser, "Seeding / Filtration");
//     addOption(parser, ArgParseOption("su",
//                                             "ungapped-seeds",
//                                             "allow only mismatches in seeds.",
//                                             ArgParseArgument::INTEGER));
//     setDefaultValue(parser, "ungapped-seeds", "1");

    addOption(parser, ArgParseOption("sl", "seed-length",
        "Length of the seeds (default = 14 for BLASTN).",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "seed-length", "10");
    setAdvanced(parser, "seed-length");

    addOption(parser, ArgParseOption("so", "seed-offset",
        "Offset for seeding (if unset = seed-length, non-overlapping; "
        "default = 5 for BLASTN).",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "seed-offset", "10");
    setAdvanced(parser, "seed-offset");

    addOption(parser, ArgParseOption("sd", "seed-delta",
        "maximum seed distance.",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "seed-delta", "1");
    setAdvanced(parser, "seed-delta");

    addOption(parser, ArgParseOption("sg", "seed-gravity",
        "Seeds closer than this are merged into region (if unset = "
        "seed-length).",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "seed-gravity", "10");
    hideOption(parser, "seed-gravity"); // HIDDEN

    addOption(parser, ArgParseOption("sm", "seed-min-length",
        "after postproc shorter seeds are discarded (if unset = seed-length).",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "seed-min-length", "10");
    hideOption(parser, "seed-min-length"); // HIDDEN

    addSection(parser, "Miscellaneous Heuristics");

    addOption(parser, ArgParseOption("ps", "pre-scoring",
        "evaluate score of a region NUM times the size of the seed "
        "before extension (0 -> no pre-scoring, 1 -> evaluate seed, n-> area "
        "around seed, as well; default = 1 if no reduction is used).",
        ArgParseArgument::INTEGER));
    setMinValue(parser, "pre-scoring", "1");
    setDefaultValue(parser, "pre-scoring", "2");
    setAdvanced(parser, "pre-scoring");

    addOption(parser, ArgParseOption("pt", "pre-scoring-threshold",
        "minimum average score per position in pre-scoring region.",
        ArgParseArgument::DOUBLE));
    setDefaultValue(parser, "pre-scoring-threshold", "2");
    setAdvanced(parser, "pre-scoring-threshold");

    addOption(parser, ArgParseOption("pd", "filter-putative-duplicates",
        "filter hits that will likely duplicate a match already found.",
        ArgParseArgument::STRING));
    setValidValues(parser, "filter-putative-duplicates", "on off");
    setDefaultValue(parser, "filter-putative-duplicates", "on");
    setAdvanced(parser, "filter-putative-duplicates");

    addOption(parser, ArgParseOption("pa", "filter-putative-abundant",
        "If the maximum number of matches per query are found already, "
        "stop searching if the remaining realm looks unfeasable.",
        ArgParseArgument::STRING));
    setValidValues(parser, "filter-putative-abundant", "on off");
    setDefaultValue(parser, "filter-putative-abundant", "on");
    setAdvanced(parser, "filter-putative-abundant");

//     addOption(parser, ArgParseOption("se",
//                                             "seedminevalue",
//                                             "after postproc worse seeds are "
//                                             "discarded"
//                                             "(0 -> off).",
//                                             ArgParseArgument::INTEGER));
//     setDefaultValue(parser, "seedminevalue", "100000");

//     addOption(parser, ArgParseOption("sb",
//                                             "seedminbits",
//                                             "after postproc worse seeds are "
//                                             "discarded"
//                                             "(-1 -> off).",
//                                             ArgParseArgument::DOUBLE));
//     setDefaultValue(parser, "seedminbits", "-1");

    addSection(parser, "Scoring");

    addOption(parser, ArgParseOption("sc", "scoring-scheme",
        "use '45' for Blosum45; '62' for Blosum62 (default); '80' for Blosum80; "
        "[ignored for BlastN]",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "scoring-scheme", "62");
    setAdvanced(parser, "scoring-scheme");

    addOption(parser, ArgParseOption("ge", "score-gap",
        "Score per gap character (default = -2 for BLASTN).",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "score-gap", "-1");
    setAdvanced(parser, "score-gap");

    addOption(parser, ArgParseOption("go", "score-gap-open",
        "Additional cost for opening gap (default = -5 for BLASTN).",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "score-gap-open", "-11");
    setAdvanced(parser, "score-gap-open");

    addOption(parser, ArgParseOption("ma", "score-match",
        "Match score [only BLASTN])",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "score-match", "2");
    setAdvanced(parser, "score-match");

    addOption(parser, ArgParseOption("mi", "score-mismatch",
        "Mismatch score [only BLASTN]",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "score-mismatch", "-3");
    setAdvanced(parser, "score-mismatch");

    addSection(parser, "Extension");

    addOption(parser, ArgParseOption("x", "x-drop",
        "Stop Banded extension if score x below the maximum seen (-1 means no "
        "xdrop).",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "x-drop", "30");
    setMinValue(parser, "x-drop", "-1");
    setAdvanced(parser, "x-drop");

    addOption(parser, ArgParseOption("b", "band",
        "Size of the DP-band used in extension (-3 means log2 of query length; "
        "-2 means sqrt of query length; -1 means full dp; n means band of size "
        "2n+1)",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "band", "-3");
    setMinValue(parser, "band", "-3");
    setAdvanced(parser, "band");

    addTextSection(parser, "Tuning");
    addText(parser, "Tuning the seeding parameters and (de)activating alphabet "
                    "reduction has a strong "
                    "influence on both speed and sensitivity. We recommend the "
                    "following alternative profiles for protein searches:");
    addText(parser, "fast (high similarity):       -ar none -sl 7 -sd 0");
    addText(parser, "sensitive (lower similarity): -so 5");

    addText(parser, "For further information see the wiki: <https://github.com/seqan/lambda/wiki>");
//         addTextSection(parser, "Speed VS memory requirements");
//         addText(parser, "Lambda requires approximately the following amount of RAM:"
//                         " \033[1msize(queryFile) + size(dbIDs) + 2 * size(dbSeqs)\033[0m. "
//                         "If you have more RAM, use double indexing and SA:\n"
//                         "\033[1m-di sa -qi radix\033[0m "
//                         "which will result in an additional speed-up of up to 30% "
//                         "compared to the published version (you need to run the "
//                         "indexer with \033[1m-di sa \033[0m, as well). The amount "
//                         "of RAM required will be: "
//                         "\033[1msize(queryFile) + size(dbIDs) + 7 * size(dbSeqs) + n\033[0m "
//                         "where n grows slowly but linearly with input size. "
//                         "Note that size(dbSeqs) refers to the total "
//                         "sequence length and does not include IDs (so it is less "
//                         "than the size of the file).");
//         addText(parser, "To save more RAM, you can define "
//                         "LAMBDA_BITCOPMRESSED_STRINGS while compiling lambda. "
//                         "This will reduce memory usage by about:"
//                         " \033[1m0.3 * ( size(queryFile) + size(dbSeqs) )\033[0m,"
//                         " but slow down lambda by about 10%.");

    // Parse command line.
    ArgumentParser::ParseResult res = parse(parser, argc, argv);

    // Only extract  options if the program will continue after parseCommandLine()
    if (res != ArgumentParser::PARSE_OK)
        return res;

    // Options shared by lambda and its indexer
    res = parseCommandLineShared(options, parser);
    if (res != ArgumentParser::PARSE_OK)
        return res;

    std::string buffer;

    // Extract option values.
    getOptionValue(options.queryFile, parser, "query");
//     if (endsWith(options.queryFile, ".fastq") ||
//         endsWith(options.queryFile, ".fq"))
//         options.fileFormat = 1;
//     else
//         options.fileFormat = 0;

    getOptionValue(options.output, parser, "output");
    buffer = options.output;
    if (endsWith(buffer, ".gz"))
        buffer.resize(length(buffer) - 3);
    else if (endsWith(buffer, ".bz2"))
        buffer.resize(length(buffer) - 4);

    if (endsWith(buffer, ".sam"))
        options.outFileFormat = 1;
    else if (endsWith(buffer, ".bam"))
        options.outFileFormat = 2;
    else
        options.outFileFormat = 0;

    clear(buffer);
    getOptionValue(buffer, parser, "sam-with-refheader");
    options.samWithRefHeader = (buffer == "on");

    clear(buffer);
    getOptionValue(buffer, parser, "sam-bam-seq");
    if (buffer == "never")
        options.samBamSeq = 0;
    else if (buffer == "uniq")
        options.samBamSeq = 1;
    else
        options.samBamSeq = 2;

    clear(buffer);
    getOptionValue(buffer, parser, "sam-bam-clip");
    options.samBamHardClip = (buffer == "hard");

    clear(buffer);
    getOptionValue(buffer, parser, "output-columns");
    if (buffer == "help")
    {
        std::cout << "Please specify the columns in this format -oc 'column1 column2', i.e. space-seperated and "
                  << "enclosed in single quotes.\nThe specifiers are the same as in NCBI Blast, currently "
                  << "the following are supported:\n";
        for (unsigned i = 0; i < length(BlastMatchField<>::implemented); ++i)
        {
            if (BlastMatchField<>::implemented[i])
            {
                std::cout << "\t" << BlastMatchField<>::optionLabels[i]
                          << (length(BlastMatchField<>::optionLabels[i]) >= 8 ? "\t" : "\t\t")
                          << BlastMatchField<>::descriptions[i] << "\n";
            }
        }
        return ArgumentParser::PARSE_HELP;
    }
    else
    {
        StringSet<CharString> fields;
        strSplit(fields, buffer, IsSpace(), false);
        for (auto str : fields)
        {
            bool resolved = false;
            for (unsigned i = 0; i < length(BlastMatchField<>::optionLabels); ++i)
            {
                if (BlastMatchField<>::optionLabels[i] == str)
                {
                    appendValue(options.columns, static_cast<BlastMatchField<>::Enum>(i));
                    resolved = true;
                    break;
                }
            }
            if (!resolved)
            {
                std::cerr << "Unknown column specifier \"" << str << "\". Please see -oc help for valid options.\n";
                return ArgumentParser::PARSE_ERROR;
            }
        }
    }
    clear(buffer);

    getOptionValue(buffer, parser, "sam-bam-tags");
    if (buffer == "help")
    {
        std::cout << "Please specify the tags in this format -oc 'tag1 tag2', i.e. space-seperated and "
                  << "enclosed in quotes. The order of tags is not preserved.\nThe following specifiers are "
                  << "supported:\n";

        for (auto const & c : SamBamExtraTags<>::keyDescPairs)
            std::cout << "\t" << std::get<0>(c) << "\t" << std::get<1>(c) << "\n";

        return ArgumentParser::PARSE_HELP;
    }
    else
    {
        StringSet<CharString> fields;
        strSplit(fields, buffer, IsSpace(), false);
        for (auto str : fields)
        {
            bool resolved = false;
            for (unsigned i = 0; i < length(SamBamExtraTags<>::keyDescPairs); ++i)
            {
                if (std::get<0>(SamBamExtraTags<>::keyDescPairs[i]) == str)
                {
                    options.samBamTags[i] = true;
                    resolved = true;
                    break;
                }
            }
            if (!resolved)
            {
                std::cerr << "Unknown column specifier \"" << str
                          << "\". Please see \"--sam-bam-tags help\" for valid options.\n";
                return ArgumentParser::PARSE_ERROR;
            }
        }
    }
    // if original is protein, than only write if explicitly asked for
    if (((options.blastProgram == BlastProgram::BLASTP) || (options.blastProgram == BlastProgram::TBLASTN)) &&
        (!options.samBamTags[SamBamExtraTags<>::Q_AA_CIGAR]))
        options.samBamSeq = 0;

    getOptionValue(buffer, parser, "version-to-outputfile");
    options.versionInformationToOutputFile = (buffer == "on");

    clear(buffer);
    getOptionValue(options.seedLength, parser, "seed-length");
    if ((!isSet(parser, "seed-length")) &&
        (options.blastProgram == BlastProgram::BLASTN))
        options.seedLength = 14;

    if (isSet(parser, "seed-offset"))
        getOptionValue(options.seedOffset, parser, "seed-offset");
    else
        options.seedOffset = options.seedLength;

    if (isSet(parser, "seed-gravity"))
        getOptionValue(options.seedGravity, parser, "seed-gravity");
    else
        options.seedGravity = options.seedLength;

    if (isSet(parser, "seed-min-length"))
        getOptionValue(options.minSeedLength, parser, "seed-min-length");
    else
        options.minSeedLength = options.seedLength;

    getOptionValue(options.maxSeedDist, parser, "seed-delta");


    getOptionValue(buffer, parser, "query-index-type");
    options.doubleIndexing = (buffer == "radix");

    getOptionValue(options.eCutOff, parser, "e-value");
    getOptionValue(options.idCutOff, parser, "percent-identity");

    getOptionValue(options.xDropOff, parser, "x-drop");
//     if ((!isSet(parser, "x-drop")) &&
//         (options.blastProgram == BlastProgram::BLASTN))
//         options.xDropOff = 16;

    getOptionValue(options.band, parser, "band");

    if (options.doubleIndexing)
    {
        if (isSet(parser, "query-partitions"))
            getOptionValue(options.queryPart, parser, "query-partitions");
        else
            options.queryPart = options.threads;
        if ((options.queryPart % options.threads) != 0)
            std::cout << "-qp not a multiple of -t; expect suboptimal performance.\n";
    } else
    {
        options.queryPart = 1;
    }

    getOptionValue(options.scoringMethod, parser, "scoring-scheme");
    if (options.blastProgram == BlastProgram::BLASTN)
        options.scoringMethod = 0;
    switch (options.scoringMethod)
    {
        case 0:
            getOptionValue(options.misMatch, parser, "score-mismatch");
            getOptionValue(options.match, parser, "score-match");
            break;
        case 45: case 62: case 80: break;
        default:
            std::cerr << "Unsupported Scoring Scheme selected.\n";
            return ArgumentParser::PARSE_ERROR;
    }

    getOptionValue(options.gapExtend, parser, "score-gap");
    if ((!isSet(parser, "score-gap")) &&
        (options.blastProgram == BlastProgram::BLASTN))
        options.gapExtend = -2;

    getOptionValue(options.gapOpen, parser, "score-gap-open");
    if ((!isSet(parser, "score-gap-open")) &&
        (options.blastProgram == BlastProgram::BLASTN))
        options.gapOpen = -5;

    getOptionValue(buffer, parser, "filter-putative-duplicates");
    options.filterPutativeDuplicates = (buffer == "on");

    getOptionValue(buffer, parser, "filter-putative-abundant");
    options.filterPutativeAbundant = (buffer == "on");

    // TODO always prescore 1
    getOptionValue(options.preScoring, parser, "pre-scoring");
    if ((!isSet(parser, "pre-scoring")) &&
        (options.alphReduction == 0))
        options.preScoring = 1;

    getOptionValue(options.preScoringThresh, parser, "pre-scoring-threshold");
//     if (options.preScoring == 0)
//         options.preScoringThresh = 4;

    int numbuf;
    getOptionValue(numbuf, parser, "num-matches");
    options.maxMatches = static_cast<unsigned long>(numbuf);

    return ArgumentParser::PARSE_OK;
}

// INDEXER
ArgumentParser::ParseResult
parseCommandLine(LambdaIndexerOptions & options, int argc, char const ** argv)
{
    // Setup ArgumentParser.
    ArgumentParser parser("lambda_indexer");

    setShortDescription(parser, "indexer for creating lambda-compatible databases");

    // Define usage line and long description.
    addUsageLine(parser, "[\\fIOPTIONS\\fP] \\-d DATABASE.fasta\\fP");

    sharedSetup(parser);

    addDescription(parser, "This is the indexer_binary for creating lambda-compatible databases.");

    addSection(parser, "Input Options");
    addOption(parser, ArgParseOption("d", "database",
        "Database sequences.",
        ArgParseArgument::INPUT_FILE,
        "IN"));
    setRequired(parser, "database");
    setValidValues(parser, "database", getFileExtensions(SeqFileIn()));

    addOption(parser, ArgParseOption("s",
        "segfile",
        "SEG intervals for database"
        "(optional).",
        ArgParseArgument::INPUT_FILE));

    setValidValues(parser, "segfile", "seg");

    addSection(parser, "Output Options");
//     addOption(parser, ArgParseOption("o",
//                                             "output",
//                                             "Index of database sequences",
//                                             ArgParseArgument::OUTPUT_FILE,
//                                             "OUT"));
//     setValidValues(parser, "output", "sa fm");

    addOption(parser, ArgParseOption("di", "db-index-type",
        "Suffix array or full-text minute space.",
        ArgParseArgument::STRING,
        "type"));
    setValidValues(parser, "db-index-type", "sa fm");
    setDefaultValue(parser, "db-index-type", "fm");
    setAdvanced(parser, "db-index-type");

    addOption(parser, ArgParseOption("", "truncate-ids",
        "Truncate IDs at first whitespace. This saves a lot of space and is irrelevant for all LAMBDA output formats "
        "other than BLAST Pairwise (.m0).",
        ArgParseArgument::STRING,
        "STR"));
    setValidValues(parser, "truncate-ids", "on off");
    setDefaultValue(parser, "truncate-ids", "on");

    addSection(parser, "Alphabets and Translation");
    addOption(parser, ArgParseOption("p", "program",
        "Blast Operation Mode.",
        ArgParseArgument::STRING,
        "program"));
    setValidValues(parser, "program", "blastn blastp blastx tblastn tblastx");
    setDefaultValue(parser, "program", "blastx");
    addOption(parser, ArgParseOption("g", "genetic-code",
        "The translation table to use (not for BlastN, BlastP). See "
        "https://www.ncbi.nlm.nih.gov/Taxonomy/Utils/wprintgc.cgi?mode=c"
        " for ids (default is generic).",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "genetic-code", "1");
    setAdvanced(parser, "genetic-code");

    addOption(parser, ArgParseOption("ar", "alphabet-reduction",
        "Alphabet Reduction for seeding phase (ignored for BLASTN).",
        ArgParseArgument::STRING,
        "STR"));
    setValidValues(parser, "alphabet-reduction", "none murphy10");
    setDefaultValue(parser, "alphabet-reduction", "murphy10");
    setAdvanced(parser, "alphabet-reduction");

    addSection(parser, "Algorithm");
    addOption(parser, ArgParseOption("a", "algorithm",
        "Algorithm for SA construction (also used for FM; see Memory "
        " Requirements below!).",
        ArgParseArgument::STRING,
        "STR"));
    setValidValues(parser, "algorithm", "mergesort quicksortbuckets quicksort radixsort skew7ext");
    setDefaultValue(parser, "algorithm", "radixsort");
    setAdvanced(parser, "algorithm");

#ifdef _OPENMP
    addOption(parser, ArgParseOption("t", "threads",
        "number of threads to run concurrently (ignored if a == skew7ext).",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "threads", omp_get_max_threads());
#else
    addOption(parser, ArgParseOption("t", "threads",
        "LAMBDA BUILT WITHOUT OPENMP; setting this option has no effect.",
        ArgParseArgument::INTEGER));
    setDefaultValue(parser, "threads", 1);
#endif
    setAdvanced(parser, "threads");

    std::string tmpdir;
    getCwd(tmpdir);
    addOption(parser, ArgParseOption("td", "tmp-dir",
        "temporary directory used by skew, defaults to working directory.",
        ArgParseArgument::STRING,
        "STR"));
    setDefaultValue(parser, "tmp-dir", tmpdir);
    setAdvanced(parser, "tmp-dir");

    //TODO move manual / auto-detect
//     addTextSection(parser, "Memory requirements and Speed");
//     addText(parser, "\033[1mmergesort [RAM]:\033[0m"
//                     "\t14 * size(dbSeqs)");
//     addText(parser, "\033[1mmergesort [speed]:\033[0m"
//                     "\tup to t threads");
//     addText(parser, "\033[1mquicksort and quicksortbuckets [RAM]:\033[0m"
//                     "\t7 * size(dbSeqs)");
//     addText(parser, "\033[1mquicksort [speed]:\033[0m"
//                     "\t1-2 threads");
//     addText(parser, "\033[1mquicksortbuckets [speed]:\033[0m"
//                     "\t1-2 threads for initial sort, up to t for buckets");
//     addText(parser, "\033[1mskew7ext [RAM]:\033[0m"
//                     "\t2 * size(dbSeqs)");
//     addText(parser, "\033[1mskew7ext [DISK]:\033[0m"
//                     "\t30 * size(dbSeqs)");
//     addText(parser, "\033[1mskew7ext [speed]:\033[0m"
//                     "\tnot parallelized");
//     addText(parser, "size(dbSeqs) refers to the total "
//                     "sequence length and does not include IDs (which can "
//                     "account for >50% of the file size for protein databases). "
//                     "The space is the maximum obseverved factor, for many "
//                     "databases the factor is smaller." );
//     addText(parser, "Use mergesort if you have enough memory! If not, you will "
//                     "probably want to use skew. For small databases and only a "
//                     "few cores the quicksorts might be a good tradeoff. "
//                     "mergesort and quicksortbuckets provide a rough progress "
//                     "estimate.");
// //     addText(parser, "Disk space required is in TMPDIR which you can set as "
// //                     "an environment variable.");

    addTextSection(parser, "Remarks");
    addText(parser, "Please see the wiki (<https://github.com/seqan/lambda/wiki>) for more information on which indexes"
        " to chose and which algorithms to pick.");

    addText(parser, "Note that the indexes created are binary and not compatible between different CPU endiannesses. "
        "Also the on-disk format is still subject to change between Lambda versions.");

    // Parse command line.
    ArgumentParser::ParseResult res = parse(parser, argc, argv);

    // Only extract  options if the program will continue after parseCommandLine()
    if (res != ArgumentParser::PARSE_OK)
        return res;

    // Options shared by lambda and its indexer
    res = parseCommandLineShared(options, parser);
    if (res != ArgumentParser::PARSE_OK)
        return res;

    // Extract option values
    getOptionValue(options.segFile, parser, "segfile");
    getOptionValue(options.algo, parser, "algorithm");
    if ((options.algo == "mergesort") || (options.algo == "quicksort") || (options.algo == "quicksortbuckets"))
    {
        std::cerr << "WARNING: " << options.algo << " tag is deprecated and superseded by \"radixsort\", please "
                  << "adapt your program calls.\n";
        options.algo = "radixsort";
    }

    getOptionValue(tmpdir, parser, "tmp-dir");
    setEnv("TMPDIR", tmpdir);

    std::string buffer;
    getOptionValue(buffer, parser, "truncate-ids");
    options.truncateIDs = (buffer == "on");

    return ArgumentParser::PARSE_OK;
}

// SHARED
ArgumentParser::ParseResult
parseCommandLineShared(SharedOptions & options, ArgumentParser & parser)
{
    int buf = 0;
    std::string buffer;

    getOptionValue(options.dbFile, parser, "database");

    getOptionValue(buffer, parser, "db-index-type");
    if (buffer == "sa")
        options.dbIndexType = 0;
    else // if fm
        options.dbIndexType = 1;

    getOptionValue(buffer, parser, "program");
    if (buffer == "blastn")
        options.blastProgram = BlastProgram::BLASTN;
    else if (buffer == "blastp")
        options.blastProgram = BlastProgram::BLASTP;
    else if (buffer == "blastx")
        options.blastProgram = BlastProgram::BLASTX;
    else if (buffer == "tblastn")
        options.blastProgram = BlastProgram::TBLASTN;
    else if (buffer == "tblastx")
        options.blastProgram = BlastProgram::TBLASTX;
    else
        return ArgumentParser::PARSE_ERROR;

    getOptionValue(buffer, parser, "alphabet-reduction");
    if ((buffer == "murphy10") &&
        (options.blastProgram != BlastProgram::BLASTN))
        options.alphReduction = 2;
    else
        options.alphReduction = 0;

    getOptionValue(buf, parser, "genetic-code");
    switch (buf)
    {
        case 1: case 2: case 3: case 4: case 5: case 6:
        case 9: case 10: case 11: case 12: case 13: case 14: case 15: case 16:
        case 21: case 22: case 23: case 24 : case 25:
            options.geneticCode = static_cast<GeneticCodeSpec>(buf);
            break;
        default:
            std::cerr << "Invalid genetic code. See trans_table vars at "
                      << "https://www.ncbi.nlm.nih.gov/Taxonomy/Utils/wprintgc.cgi?mode=c"
                      << std::endl;
            return ArgumentParser::PARSE_ERROR;
    }

#ifdef _OPENMP
    getOptionValue(options.threads, parser, "threads");
    omp_set_num_threads(options.threads);
#else
    options.threads = 1;
#endif

    getOptionValue(buf, parser, "verbosity");
    switch(buf)
    {
        case 0: options.verbosity = 0; break;
        case 2: options.verbosity = 2; break;
        default: options.verbosity = 1; break;
    }

    return ArgumentParser::PARSE_OK;
}

constexpr const char *
_alphName(AminoAcid const & /**/)
{
    return "aminoacid";
}

constexpr const char *
_alphName(ReducedAminoAcid<Murphy10> const & /**/)
{
    return "murphy10";
}

// constexpr const char *
// _alphName(ReducedAminoAcid<ClusterReduction<8>> const & /**/)
// {
//     return "lambda08";
// }
//
// constexpr const char *
// _alphName(ReducedAminoAcid<ClusterReduction<10>> const & /**/)
// {
//     return "lambda10";
// }
//
// constexpr const char *
// _alphName(ReducedAminoAcid<ClusterReduction<12>> const & /**/)
// {
//     return "lambda12";
// }

constexpr const char *
_alphName(Dna const & /**/)
{
    return "dna4";
}

constexpr const char *
_alphName(Dna5 const & /**/)
{
    return "dna5";
}

template <typename TLH>
inline void
printOptions(LambdaOptions const & options)
{
    using TGH = typename TLH::TGlobalHolder;

    std::string bandStr;
    switch(options.band)
    {
        case -3: bandStr = "2 * log(queryLength) + 1"; break;
        case -2: bandStr = "2 * sqrt(queryLength) + 1"; break;
        case -1: bandStr = "no band"; break;
        default: bandStr = std::to_string(2 * options.band + 1); break;
    }

    std::cout << "OPTIONS\n"
              << " INPUT\n"
              << "  query file:               " << options.queryFile << "\n"
              << "  db file:                  " << options.dbFile << "\n"
              << "  db index type:            " << (TGH::indexIsFM
                                                    ? "FM-Index\n"
                                                    : "SA-Index\n")
              << " OUTPUT (file)\n"
              << "  output file:              " << options.output << "\n"
              << "  minimum % identity:       " << options.idCutOff << "\n"
              << "  maximum e-value:          " << options.eCutOff << "\n"
              << "  max #matches per query:   " << options.maxMatches << "\n"
              << "  include subj names in sam:" << options.samWithRefHeader << "\n"
              << "  include seq in sam/bam:   " << options.samBamSeq << "\n"
              << " OUTPUT (stdout)\n"
              << "  stdout is terminal:       " << options.isTerm << "\n"
              << "  terminal width:           " << options.terminalCols << "\n"
              << "  verbosity:                " << options.verbosity << "\n"
              << " GENERAL\n"
              << "  double indexing:          " << options.doubleIndexing << "\n"
              << "  threads:                  " << uint(options.threads) << "\n"
              << "  query partitions:         " << (options.doubleIndexing
                                                    ? std::to_string(options.queryPart)
                                                    : std::string("n/a")) << "\n"
              << " TRANSLATION AND ALPHABETS\n"
              << "  genetic code:             "
              << ((TGH::blastProgram != BlastProgram::BLASTN) &&
                  (TGH::blastProgram != BlastProgram::BLASTP)
                 ? std::to_string(options.geneticCode)
                 : std::string("n/a")) << "\n"
              << "  blast mode:               " << _programTagToString(TGH::blastProgram)
              << "\n"
              << "  original alphabet (query):" << _alphName(OrigQryAlph<TGH::blastProgram>())
              << "\n"
              << "  original alphabet (subj): " << _alphName(OrigSubjAlph<TGH::blastProgram>())
              << "\n"
              << "  translated alphabet:      " << _alphName(TransAlph<TGH::blastProgram>())
              << "\n"
              << "  reduced alphabet:         " << _alphName(typename TGH::TRedAlph())
              << "\n"
              << " SEEDING\n"
              << "  seed length:              " << uint(options.seedLength) << "\n"
              << "  seed offset:              " << uint(options.seedOffset) << "\n"
              << "  seed delta:               " << uint(options.maxSeedDist) << "\n"
              << "  seeds ungapped:           " << uint(options.hammingOnly) << "\n"
              << "  seed gravity:             " << uint(options.seedGravity) << "\n"
              << "  min seed length:          " << uint(options.minSeedLength) << "\n"
              << " MISCELLANEOUS HEURISTICS\n"
              << "  pre-scoring:              " << (options.preScoring
                                                    ? std::string("on")
                                                    : std::string("off")) << "\n"
              << "  pre-scoring-region:       " << (options.preScoring
                                                    ? std::to_string(
                                                        options.preScoring *
                                                        options.seedLength)
                                                    : std::string("n/a")) << "\n"
              << "  pre-scoring-threshold:    " << (options.preScoring
                                                    ? std::to_string(
                                                       options.preScoringThresh)
                                                    : std::string("n/a")) << "\n"
              << "  putative-abundancy:       " << (options.filterPutativeAbundant
                                                    ? std::string("on")
                                                    : std::string("off")) << "\n"
              << "  putative-duplicates:      " << (options.filterPutativeDuplicates
                                                    ? std::string("on")
                                                    : std::string("off")) << "\n"
              << " SCORING\n"
              << "  scoring scheme:           " << options.scoringMethod << "\n"
              << "  score-match:              " << (options.scoringMethod
                                                    ? std::string("n/a")
                                                    : std::to_string(options.match)) << "\n"
              << "  score-mismatch:           " << (options.scoringMethod
                                                    ? std::string("n/a")
                                                    : std::to_string(options.misMatch)) << "\n"
              << "  score-gap:                " << options.gapExtend << "\n"
              << "  score-gap-open:           " << options.gapOpen << "\n"
              << " EXTENSION\n"
              << "  x-drop:                   " << options.xDropOff << "\n"
              << "  band:                     " << bandStr << "\n"
              << " BUILD OPTIONS:\n"
              << "  cmake_build_type:         " << std::string(CMAKE_BUILD_TYPE) << "\n"
              << "  fastbuild:                "
    #if defined(FASTBUILD)
              << "on\n"
    #else
              << "off\n"
    #endif
              << "  native_build:             "
    #if defined(LAMBDA_NATIVE_BUILD)
              << "on\n"
    #else
              << "off\n"
    #endif
              << "  static_build:             "
    #if defined(LAMBDA_STATIC_BUILD)
              << "on\n"
    #else
              << "off\n"
    #endif
              << "  mmapped_db:               "
    #if defined(LAMBDA_MMAPPED_DB)
              << "on\n"
    #else
              << "off\n"
    #endif
              << "  lingaps_opt:              "
    #if defined(LAMBDA_LINGAPS_OPT)
              << "on\n"
    #else
              << "off\n"
    #endif
              << "\n";
}


#endif // header guard