1 /* -*- mode: c; tab-width: 4; c-basic-offset: 4; indent-tabs-mode: nil -*- */
3 /*********************************************************************
4 * Clustal Omega - Multiple sequence alignment
6 * Copyright (C) 2010 University College Dublin
8 * Clustal-Omega is free software; you can redistribute it and/or
9 * modify it under the terms of the GNU General Public License as
10 * published by the Free Software Foundation; either version 2 of the
11 * License, or (at your option) any later version.
13 * This file is part of Clustal-Omega.
15 ********************************************************************/
18 * RCS $Id: hhalign.cpp 250 2011-06-17 13:06:20Z fabian $
22 * Align a multiple alignment to an alignment or HMM
23 * Print out aligned input sequences in a3m format
24 * Compile: g++ hhalign.C -o hhalign -I/usr/include/ -L/usr/lib -lpng -lz -O3 -fno-strict-aliasing
25 * Compile with efence: g++ hhalign.C -o hhalign -I/usr/include/ -lefence -L/usr/lib -lpng -lz -O -g
27 * Error codes: 0: ok 1: file format error 2: file access error
28 * 3: memory error 4: internal numeric error
29 * 5: command line error
31 #undef PNG /* include options for making png files?
32 (will need the png library) */
34 #include <iostream> // cin, cout, cerr
35 #include <fstream> // ofstream, ifstream
36 #include <string.h> // strcmp, strstr
37 #include <stdio.h> // printf
38 #include <stdlib.h> // exit
39 #include <math.h> // sqrt, pow
40 #include <limits.h> // INT_MIN
41 #include <float.h> // FLT_MIN
42 #include <time.h> // clock
43 #include <ctype.h> // islower, isdigit etc
55 #include "util-C.h" /* imax, fmax, iround, iceil, ifloor, strint, strscn,
56 strcut, substr, uprstr, uprchr, Basename etc. */
57 #include "list-C.h" // list data structure
58 #include "hash-C.h" // hash data structure
60 #include "hhdecl-C.h" // Constants, global variables, struct Parameters
61 #include "hhutil-C.h" /* MatchChr, InsertChr, aa2i, i2aa, log2,
62 fast_log2, WriteToScreen, */
63 #include "hhmatrices-C.h" // BLOSUM50, GONNET, HSDM
64 #include "hhhmm.h" // class HMM
65 #include "hhhit.h" // class Hit
66 #include "hhalignment.h" // class Alignment
67 #include "hhhalfalignment.h" // class HalfAlignment
68 #include "hhfullalignment.h" // class FullAlignment
69 #include "hhhitlist.h" // class Hit
70 #include "hhhmm-C.h" // class HMM
71 #include "hhalignment-C.h" // class Alignment
72 #include "hhhit-C.h" // class Hit
73 #include "hhhalfalignment-C.h" // class HalfAlignment
74 #include "hhfullalignment-C.h" // class FullAlignment
75 #include "hhhitlist-C.h" // class HitList
76 #include "hhfunc-C.h" // some functions common to hh programs
79 #include "pngwriter.h" //PNGWriter (http://pngwriter.sourceforge.net/)
80 #include "pngwriter.cc" //PNGWriter (http://pngwriter.sourceforge.net/)
83 ////////////////////////////////////////////////////////////////////////
85 ////////////////////////////////////////////////////////////////////////
86 HMM *q; // Create query HMM with maximum of MAXRES match states
87 HMM *t; /* Create template HMM with maximum of MAXRES
89 Alignment *qali; /* (query alignment might be needed outside of hhfunc.C
91 Hit hit; // Ceate new hit object pointed at by hit
92 HitList hitlist; /* list of hits with one Hit object for each
93 pairwise comparison done */
94 char aliindices[256]; /* hash containing indices of all alignments
95 which to show in dot plot */
96 char* dmapfile=NULL; /* where to write the coordinates for the HTML map file
97 (to click the alignments) */
98 char* alitabfile=NULL; // where to write pairs of aligned residues
99 char* strucfile=NULL; // where to read structure scores
100 char* pngfile=NULL; // pointer to pngfile
101 char* tcfile=NULL; // TCoffee output file name
102 float probmin_tc=0.05; /* 5% minimum posterior probability for printing
103 pairs of residues for TCoffee */
105 int dotW=10; // average score of dot plot over window [i-W..i+W]
106 float dotthr=0.5; // score threshold for dot plot
107 int dotscale=600; // size scale of dotplot
108 char dotali=0; // show no alignments in dotplot
109 float dotsat=0.3; // saturation of grid and alignments in dot plot
110 float pself=0.001; // maximum p-value of 2nd and following self-alignments
111 int Nstochali=0; // # of stochastically traced alignments in dot plot
112 float** Pstruc=NULL; /* structure matrix which can be multiplied to prob
113 ratios from aa column comparisons in Forward() */
114 float** Sstruc=NULL; /* structure matrix which can be added to log odds
115 from aa column comparisons in Forward() */
121 ////////////////////////////////////////////////////////////////////////////
123 ////////////////////////////////////////////////////////////////////////////
125 * @brief general help function, not accessible in Clustal-Omega
132 printf("HHalign %s\n",VERSION_AND_DATE);
133 printf("Align a query alignment/HMM to a template alignment/HMM by HMM-HMM alignment\n");
134 printf("If only one alignment/HMM is given it is compared to itself and the best\n");
135 printf("off-diagonal alignment plus all further non-overlapping alignments above \n");
136 printf("significance threshold are shown.\n");
137 printf("%s",REFERENCE);
138 printf("%s",COPYRIGHT);
140 printf("Usage: %s -i query [-t template] [options] \n",program_name);
141 printf(" -i <file> input query alignment (fasta/a2m/a3m) or HMM file (.hhm)\n");
142 printf(" -t <file> input template alignment (fasta/a2m/a3m) or HMM file (.hhm)\n");
144 printf(" -png <file> write dotplot into PNG-file (default=none) \n");
147 printf("Output options: \n");
148 printf(" -o <file> write output alignment to file\n");
149 printf(" -ofas <file> write alignments in FASTA, A2M (-oa2m) or A3M (-oa3m) format \n");
150 printf(" -a <file> write query alignment in a3m format to file (default=none)\n");
151 printf(" -aa <file> append query alignment in a3m format to file (default=none)\n");
152 printf(" -atab <file> write alignment as a table (with posteriors) to file (default=none)\n");
153 printf(" -v <int> verbose mode: 0:no screen output 1:only warings 2: verbose\n");
154 printf(" -seq [1,inf[ max. number of query/template sequences displayed (def=%i) \n",par.nseqdis);
155 printf(" -nocons don't show consensus sequence in alignments (default=show) \n");
156 printf(" -nopred don't show predicted 2ndary structure in alignments (default=show) \n");
157 printf(" -nodssp don't show DSSP 2ndary structure in alignments (default=show) \n");
158 printf(" -aliw int number of columns per line in alignment list (def=%i)\n",par.aliwidth);
159 printf(" -P <float> for self-comparison: max p-value of alignments (def=%.2g\n",pself);
160 printf(" -p <float> minimum probability in summary and alignment list (def=%G) \n",par.p);
161 printf(" -E <float> maximum E-value in summary and alignment list (def=%G) \n",par.E);
162 printf(" -Z <int> maximum number of lines in summary hit list (def=%i) \n",par.Z);
163 printf(" -z <int> minimum number of lines in summary hit list (def=%i) \n",par.z);
164 printf(" -B <int> maximum number of alignments in alignment list (def=%i) \n",par.B);
165 printf(" -b <int> minimum number of alignments in alignment list (def=%i) \n",par.b);
166 printf(" -rank int specify rank of alignment to write with -a or -aa option (default=1)\n");
169 printf("Dotplot options:\n");
170 printf(" -dwin <int> average score in dotplot over window [i-W..i+W] (def=%i) \n",dotW);
171 printf(" -dthr <float> score threshold for dotplot (default=%.2f) \n",dotthr);
172 printf(" -dsca <int> if value <= 20: size of dot plot unit box in pixels \n");
173 printf(" if value > 20: maximum dot plot size in pixels (default=%i) \n",dotscale);
174 printf(" -dali <list> show alignments with indices in <list> in dot plot \n");
175 printf(" <list> = <index1> ... <indexN> or <list> = all \n");
178 printf("Filter input alignment (options can be combined): \n");
179 printf(" -id [0,100] maximum pairwise sequence identity (%%) (def=%i) \n",par.max_seqid);
180 printf(" -diff [0,inf[ filter most diverse set of sequences, keeping at least this \n");
181 printf(" many sequences in each block of >50 columns (def=%i)\n",par.Ndiff);
182 printf(" -cov [0,100] minimum coverage with query (%%) (def=%i) \n",par.coverage);
183 printf(" -qid [0,100] minimum sequence identity with query (%%) (def=%i) \n",par.qid);
184 printf(" -qsc [0,100] minimum score per column with query (def=%.1f)\n",par.qsc);
185 // printf(" -csc [0,100] minimum score per column with core alignment (def=%-.2f)\n",par.coresc);
186 // printf(" -qscc [0,100] minimum score per column of core sequence with query (def=%-.2f)\n",par.qsc_core);
188 printf("Input alignment format: \n");
189 printf(" -M a2m use A2M/A3M (default): upper case = Match; lower case = Insert;\n");
190 printf(" '-' = Delete; '.' = gaps aligned to inserts (may be omitted) \n");
191 printf(" -M first use FASTA: columns with residue in 1st sequence are match states\n");
192 printf(" -M [0,100] use FASTA: columns with fewer than X%% gaps are match states \n");
194 printf("HMM-HMM alignment options: \n");
195 printf(" -glob/-loc global or local alignment mode (def=local) \n");
196 printf(" -alt <int> show up to this number of alternative alignments (def=%i) \n",par.altali);
197 printf(" -vit use Viterbi algorithm for alignment instead of MAC algorithm \n");
198 printf(" -mac use Maximum Accuracy (MAC) alignment (default) \n");
199 printf(" -mact [0,1[ posterior probability threshold for MAC alignment (def=%.3f) \n",par.mact);
200 printf(" A threshold value of 0.0 yields global alignments.\n");
201 printf(" -sto <int> use global stochastic sampling algorithm to sample this many alignments\n");
202 printf(" -excl <range> exclude query positions from the alignment, e.g. '1-33,97-168'\n");
203 printf(" -shift [-1,1] score offset (def=%-.3f) \n",par.shift);
204 printf(" -corr [0,1] weight of term for pair correlations (def=%.2f) \n",par.corr);
205 printf(" -ssm 0-4 0:no ss scoring [default=%i] \n",par.ssm);
206 printf(" 1:ss scoring after alignment \n");
207 printf(" 2:ss scoring during alignment \n");
208 printf(" -ssw [0,1] weight of ss score (def=%-.2f) \n",par.ssw);
210 printf(" -def read default options from ./.hhdefaults or <home>/.hhdefault. \n");
212 printf("Example: %s -i T0187.a3m -t d1hz4a_.hhm -png T0187pdb.png \n",program_name);
214 // printf("More help: \n");
215 // printf(" -h out output options \n");
216 // printf(" -h hmm options for building HMM from multiple alignment \n");
217 // printf(" -h gap options for setting gap penalties \n");
218 // printf(" -h ali options for HMM-HMM alignment \n");
219 // printf(" -h all all options \n");
224 * @brief helpt for output, not accessible in Clustal-Omega
231 printf("Output options: \n");
232 printf(" -v verbose mode (default: show only warnings) \n");
233 printf(" -v 0 suppress all screen output \n");
234 printf(" -nocons don't show consensus sequence in alignments (default=show) \n");
235 printf(" -nopred don't show predicted 2ndary structure in alignments (default=show) \n");
236 printf(" -nodssp don't show DSSP SS 2ndary structure in alignments (default=show) \n");
237 printf(" -seq [1,inf[ max. number of query/template sequences displayed (def=%i) \n",par.nseqdis);
238 printf(" -aliw [40,..[ number of columns per line in alignment list (def=%i)\n",par.aliwidth);
239 printf(" -P <float> for self-comparison: max p-value of alignments (def=%.2g\n",pself);
240 printf(" -p <float> minimum probability in summary and alignment list (def=%G) \n",par.p);
241 printf(" -E <float> maximum E-value in summary and alignment list (def=%G) \n",par.E);
242 printf(" -Z <int> maximum number of lines in summary hit list (def=%i) \n",par.Z);
243 printf(" -z <int> minimum number of lines in summary hit list (def=%i) \n",par.z);
244 printf(" -B <int> maximum number of alignments in alignment list (def=%i) \n",par.B);
245 printf(" -b <int> minimum number of alignments in alignment list (def=%i) \n",par.b);
246 printf(" -rank int specify rank of alignment to write with -a or -aa option (def=1)\n");
247 printf(" -tc <file> write a TCoffee library file for the pairwise comparison \n");
248 printf(" -tct [0,100] min. probobability of residue pairs for TCoffee (def=%i%%)\n",iround(100*probmin_tc));
251 printf("Dotplot options:\n");
252 printf(" -dwin int average score in dotplot over window [i-W..i+W] (def=%i) \n",dotW);
253 printf(" -dthr float score threshold for dotplot (default=%.2f) \n",dotthr);
254 printf(" -dsca int size of dot plot box in pixels (default=%i) \n",dotscale);
255 printf(" -dali <list> show alignments with indices in <list> in dot plot\n");
256 printf(" <list> = <index1> ... <indexN> or <list> = all \n");
257 printf(" -dmap <file> print list of coordinates in png plot \n");
263 * @brief help hit HMM options, not accessible in Clustal-Omega
270 printf("Options to filter input alignment (options can be combined): \n");
271 printf(" -id [0,100] maximum pairwise sequence identity (%%) (def=%i) \n",par.max_seqid);
272 printf(" -diff [0,inf[ filter most diverse set of sequences, keeping at least this \n");
273 printf(" many sequences in each block of >50 columns (def=%i)\n",par.Ndiff);
274 printf(" -cov [0,100] minimum coverage with query (%%) (def=%i) \n",par.coverage);
275 printf(" -qid [0,100] minimum sequence identity with query (%%) (def=%i) \n",par.qid);
276 printf(" -qsc [0,100] minimum score per column with query (def=%.1f)\n",par.qsc);
277 // printf(" -csc [0,100] minimum score per column with core alignment (def=%-.2f)\n",par.coresc);
278 // printf(" -qscc [0,100] minimum score per column of core sequence with query (def=%-.2f)\n",par.qsc_core);
280 printf("HMM-building options: \n");
281 printf(" -M a2m use A2M/A3M (default): upper case = Match; lower case = Insert;\n");
282 printf(" '-' = Delete; '.' = gaps aligned to inserts (may be omitted) \n");
283 printf(" -M first use FASTA: columns with residue in 1st sequence are match states\n");
284 printf(" -M [0,100] use FASTA: columns with fewer than X%% gaps are match states \n");
285 printf(" -tags do NOT neutralize His-, C-myc-, FLAG-tags, and \n");
286 printf(" trypsin recognition sequence to background distribution \n");
288 printf("Pseudocount options: \n");
289 printf(" -Gonnet use the Gonnet substitution matrix (default) \n");
290 printf(" -Blosum50 use the Blosum50 substitution matrix \n");
291 printf(" -Blosum62 use the Blosum62 substitution matrix \n");
292 printf(" -HSDM use the structure-derived HSDM substitution matrix \n");
293 printf(" -pcm 0-3 Pseudocount mode (default=%-i) \n",par.pcm);
294 printf(" tau = substitution matrix pseudocount admixture \n");
295 printf(" 0: no pseudo counts: tau = 0 \n");
296 printf(" 1: constant tau = a \n");
297 printf(" 2: divergence-dependent: tau = a/(1 + ((Neff-1)/b)^c) \n");
298 printf(" Neff=( (Neff_q^d+Neff_t^d)/2 )^(1/d) \n");
299 printf(" Neff_q = av number of different AAs per column in query \n");
300 printf(" 3: column-specific: tau = \'2\' * (Neff(i)/Neff)^(w*Neff/20)\n");
301 printf(" -pca [0,1] set a (overall admixture) (def=%-.1f) \n",par.pca);
302 printf(" -pcb [1,inf[ set b (threshold for Neff) (def=%-.1f) \n",par.pcb);
303 printf(" -pcc [0,3] set c (extinction exponent for tau(Neff)) (def=%-.1f) \n",par.pcc);
304 printf(" -pcw [0,3] set w (weight of pos-specificity for pcs) (def=%-.1f) \n",par.pcw);
309 * @brief help with gap handling, not accessible in Clustal-Omega
316 printf("Gap cost options: \n");
317 printf(" -gapb [0,inf[ transition pseudocount admixture (def=%-.2f) \n",par.gapb);
318 printf(" -gapd [0,inf[ Transition pseudocount admixture for opening gap (default=%-.2f)\n",par.gapd);
319 printf(" -gape [0,1.5] Transition pseudocount admixture for extending gap (def=%-.1f)\n",par.gape);
320 printf(" -gapf ]0,inf] factor for increasing/reducing the gap open penalty for deletes (def=%-.2f)\n",par.gapf);
321 printf(" -gapg ]0,inf] factor for increasing/reducing the gap open penalty for deletes (def=%-.2f)\n",par.gapg);
322 printf(" -gaph ]0,inf] factor for increasing/reducing the gap extension penalty for deletes(def=%-.2f)\n",par.gaph);
323 printf(" -gapi ]0,inf] factor for increasing/reducing the gap extension penalty for inserts(def=%-.2f)\n",par.gapi);
324 printf(" -egq [0,inf[ penalty (bits) for end gaps aligned to query residues (def=%-.2f)\n",par.egq);
325 printf(" -egt [0,inf[ penalty (bits) for end gaps aligned to template residues (def=%-.2f)\n",par.egt);
330 * @brief help with alignment options, not accessible in Clustal-Omega
337 printf("Alignment options: \n");
338 printf(" -glob/-loc global or local alignment mode (def=global) \n");
339 printf(" -mac use Maximum Accuracy (MAC) alignment instead of Viterbi\n");
340 printf(" -mact [0,1] posterior prob threshold for MAC alignment (def=%.3f) \n",par.mact);
341 printf(" -sto <int> use global stochastic sampling algorithm to sample this many alignments\n");
342 printf(" -sc <int> amino acid score (tja: template HMM at column j) (def=%i)\n",par.columnscore);
343 printf(" 0 = log2 Sum(tja*qia/pa) (pa: aa background frequencies) \n");
344 printf(" 1 = log2 Sum(tja*qia/pqa) (pqa = 1/2*(pa+ta) ) \n");
345 printf(" 2 = log2 Sum(tja*qia/ta) (ta: av. aa freqs in template) \n");
346 printf(" 3 = log2 Sum(tja*qia/qa) (qa: av. aa freqs in query) \n");
347 printf(" -corr [0,1] weight of term for pair correlations (def=%.2f) \n",par.corr);
348 printf(" -shift [-1,1] score offset (def=%-.3f) \n",par.shift);
349 printf(" -r repeat identification: multiple hits not treated as independent\n");
350 printf(" -ssm 0-2 0:no ss scoring [default=%i] \n",par.ssm);
351 printf(" 1:ss scoring after alignment \n");
352 printf(" 2:ss scoring during alignment \n");
353 printf(" -ssw [0,1] weight of ss score compared to column score (def=%-.2f) \n",par.ssw);
354 printf(" -ssa [0,1] ss confusion matrix = (1-ssa)*I + ssa*psipred-confusion-matrix [def=%-.2f)\n",par.ssa);
359 * @brief general help menu, not accessible in Clustal-Omega
371 printf("Default options can be specified in './.hhdefaults' or '~/.hhdefaults'\n");
375 /////////////////////////////////////////////////////////////////////////////
376 //// Processing input options from command line and .hhdefaults file
377 /////////////////////////////////////////////////////////////////////////////
379 * @brief process arguments from commandline, not accessible from Clustal-Omega
383 ProcessArguments(int argc, char** argv)
385 //Processing command line input
386 for (int i=1; i<argc; i++)
388 if (v>=4) cout<<i<<" "<<argv[i]<<endl; //PRINT
389 if (!strcmp(argv[i],"-i"))
391 if (++i>=argc || argv[i][0]=='-')
392 {help(); cerr<<endl<<"Error in "<<program_name<<": no query file following -i\n"; exit(4);}
393 else strcpy(par.infile,argv[i]);
395 else if (!strcmp(argv[i],"-t"))
397 if (++i>=argc || argv[i][0]=='-')
398 {help(); cerr<<endl<<"Error in "<<program_name<<": no template file following -d\n"; exit(4);}
399 else strcpy(par.tfile,argv[i]); /* FS, ProcessArguments */
401 else if (!strcmp(argv[i],"-o"))
404 {help(); cerr<<endl<<"Error in "<<program_name<<": no filename following -o\n"; exit(4);}
405 else strcpy(par.outfile,argv[i]);
407 else if (!strcmp(argv[i],"-ofas"))
410 if (++i>=argc || argv[i][0]=='-')
411 {help() ; cerr<<endl<<"Error in "<<program_name<<": no output file following -o\n"; exit(4);}
412 else strcpy(par.pairwisealisfile,argv[i]);
414 else if (!strcmp(argv[i],"-oa2m"))
417 if (++i>=argc || argv[i][0]=='-')
418 {help() ; cerr<<endl<<"Error in "<<program_name<<": no output file following -o\n"; exit(4);}
419 else strcpy(par.pairwisealisfile,argv[i]);
421 else if (!strcmp(argv[i],"-oa3m"))
424 if (++i>=argc || argv[i][0]=='-')
425 {help() ; cerr<<endl<<"Error in "<<program_name<<": no output file following -o\n"; exit(4);}
426 else strcpy(par.pairwisealisfile,argv[i]);
428 else if (!strcmp(argv[i],"-rank") && (i<argc-1)) par.hitrank=atoi(argv[++i]);
429 else if (!strcmp(argv[i],"-Oa3m"))
432 if (++i>=argc || argv[i][0]=='-')
433 {help() ; cerr<<endl<<"Error in "<<program_name<<": no output file following -Oa3m\n"; exit(4);}
434 else strcpy(par.alnfile,argv[i]);
436 else if (!strcmp(argv[i],"-Aa3m"))
439 if (++i>=argc || argv[i][0]=='-')
440 {help() ; cerr<<endl<<"Error in "<<program_name<<": no output file following -Aa3m\n"; exit(4);}
441 else strcpy(par.alnfile,argv[i]);
443 else if (!strcmp(argv[i],"-Ohhm"))
446 if (++i>=argc || argv[i][0]=='-')
447 {help() ; cerr<<endl<<"Error in "<<program_name<<": no output file following -Ohhm\n"; exit(4);}
448 else strcpy(par.hhmfile,argv[i]);
450 else if (!strcmp(argv[i],"-Ahhm"))
453 if (++i>=argc || argv[i][0]=='-')
454 {help() ; cerr<<endl<<"Error in "<<program_name<<": no output file following -Ahhm\n"; exit(4);}
455 else strcpy(par.hhmfile,argv[i]);
457 else if (!strcmp(argv[i],"-Opsi"))
460 if (++i>=argc || argv[i][0]=='-')
461 {help() ; cerr<<endl<<"Error in "<<program_name<<": no output file following -Opsi\n"; exit(4);}
462 else strcpy(par.psifile,argv[i]);
464 else if (!strcmp(argv[i],"-Apsi"))
467 if (++i>=argc || argv[i][0]=='-')
468 {help() ; cerr<<endl<<"Error in "<<program_name<<": no output file following -Apsi\n"; exit(4);}
469 else strcpy(par.psifile,argv[i]);
471 else if (!strcmp(argv[i],"-png"))
474 {help(); cerr<<endl<<"Error in "<<program_name<<": no filename following -png\n"; exit(4);}
477 pngfile = new(char[strlen(argv[i])+1]);
478 strcpy(pngfile,argv[i]);
481 else if (!strcmp(argv[i],"-Struc"))
483 if (++i>=argc || argv[i][0]=='-')
484 {help(); cerr<<endl<<"Error in "<<program_name<<": no query file following -Struc\n"; exit(4);}
487 strucfile = new(char[strlen(argv[i])+1]);
488 strcpy(strucfile,argv[i]);
491 else if (!strcmp(argv[i],"-atab") || !strcmp(argv[i],"-Aliout"))
493 if (++i>=argc || argv[i][0]=='-')
494 {help(); cerr<<endl<<"Error in "<<program_name<<": no query file following -Struc\n"; exit(4);}
497 alitabfile = new(char[strlen(argv[i])+1]);
498 strcpy(alitabfile,argv[i]);
501 else if (!strcmp(argv[i],"-tc"))
503 if (++i>=argc || argv[i][0]=='-')
504 {help() ; cerr<<endl<<"Error in "<<program_name<<": no output file following -Opsi\n"; exit(4);}
507 tcfile = new(char[strlen(argv[i])+1]);
508 strcpy(tcfile,argv[i]);
511 else if (!strcmp(argv[i],"-h")|| !strcmp(argv[i],"--help"))
513 if (++i>=argc) {help(); exit(0);}
514 if (!strcmp(argv[i],"out")) {help_out(); exit(0);}
515 if (!strcmp(argv[i],"hmm")) {help_hmm(); exit(0);}
516 if (!strcmp(argv[i],"gap")) {help_gap(); exit(0);}
517 if (!strcmp(argv[i],"ali")) {help_ali(); exit(0);}
518 if (!strcmp(argv[i],"all")) {help_all(); exit(0);}
519 else {help(); exit(0);}
521 else if (!strcmp(argv[i],"-excl"))
523 if (++i>=argc) {help(); exit(4);}
524 par.exclstr = new(char[strlen(argv[i])+1]);
525 strcpy(par.exclstr,argv[i]);
527 else if (!strcmp(argv[i],"-v") && (i<argc-1) && argv[i+1][0]!='-' ) v=atoi(argv[++i]);
528 else if (!strcmp(argv[i],"-v")) v=2;
529 else if (!strcmp(argv[i],"-v0")) v=0;
530 else if (!strcmp(argv[i],"-v1")) v=1;
531 else if (!strcmp(argv[i],"-v2")) v=2;
532 else if (!strcmp(argv[i],"-v3")) v=3;
533 else if (!strcmp(argv[i],"-v4")) v=4;
534 else if (!strcmp(argv[i],"-v5")) v=5;
535 else if (!strcmp(argv[i],"-P") && (i<argc-1)) pself=atof(argv[++i]);
536 else if (!strcmp(argv[i],"-p") && (i<argc-1)) par.p = atof(argv[++i]);
537 else if (!strcmp(argv[i],"-e") && (i<argc-1)) par.E = atof(argv[++i]);
538 else if (!strcmp(argv[i],"-E") && (i<argc-1)) par.E = atof(argv[++i]);
539 else if (!strcmp(argv[i],"-b") && (i<argc-1)) par.b = atoi(argv[++i]);
540 else if (!strcmp(argv[i],"-B") && (i<argc-1)) par.B = atoi(argv[++i]);
541 else if (!strcmp(argv[i],"-z") && (i<argc-1)) par.z = atoi(argv[++i]);
542 else if (!strcmp(argv[i],"-Z") && (i<argc-1)) par.Z = atoi(argv[++i]);
543 else if (!strncmp(argv[i],"-nocons",7)) par.showcons=0;
544 else if (!strncmp(argv[i],"-nopred",7)) par.showpred=0;
545 else if (!strncmp(argv[i],"-nodssp",7)) par.showdssp=0;
546 else if (!strncmp(argv[i],"-mark",7)) par.mark=1;
547 else if (!strcmp(argv[i],"-seq") && (i<argc-1)) par.nseqdis=atoi(argv[++i]);
548 else if (!strcmp(argv[i],"-aliw") && (i<argc-1)) par.aliwidth=atoi(argv[++i]);
549 else if (!strcmp(argv[i],"-id") && (i<argc-1)) par.max_seqid=atoi(argv[++i]);
550 else if (!strcmp(argv[i],"-tct") && (i<argc-1)) probmin_tc=atoi(argv[++i]);
551 else if (!strcmp(argv[i],"-dwin") && (i<argc-1)) dotW=atoi(argv[++i]);
552 else if (!strcmp(argv[i],"-dsca") && (i<argc-1)) dotscale=atoi(argv[++i]);
553 else if (!strcmp(argv[i],"-dthr") && (i<argc-1)) dotthr=atof(argv[++i]);
554 else if (!strcmp(argv[i],"-dali") && (i<argc-1))
557 for (int index=0; index<256; index++) aliindices[index]=0;
558 while (i+1<argc && argv[i+1][0]!='-') // adds index to hash aliindices
561 if (strcmp(argv[i],"all")) aliindices[atoi(argv[i])]=1;
565 else if (!strcmp(argv[i],"-dmap"))
568 {help(); cerr<<endl<<"Error in "<<program_name<<": no filename following -o\n"; exit(4);}
571 dmapfile = new(char[strlen(argv[i])+1]);
572 strcpy(dmapfile,argv[i]);
575 else if (!strcmp(argv[i],"-dsat") && (i<argc-1)) dotsat=atof(argv[++i]);
576 else if (!strcmp(argv[i],"-qid") && (i<argc-1)) par.qid=atoi(argv[++i]);
577 else if (!strcmp(argv[i],"-qsc") && (i<argc-1)) par.qsc=atof(argv[++i]);
578 else if (!strcmp(argv[i],"-cov") && (i<argc-1)) par.coverage=atoi(argv[++i]);
579 else if (!strcmp(argv[i],"-diff") && (i<argc-1)) par.Ndiff=atoi(argv[++i]);
580 else if (!strcmp(argv[i],"-qscc") && (i<argc-1)) par.qsc_core=atof(argv[++i]);
581 else if (!strcmp(argv[i],"-csc") && (i<argc-1)) par.coresc=atof(argv[++i]);
582 else if (!strcmp(argv[i],"-Gonnet")) par.matrix=0;
583 else if (!strcmp(argv[i],"-HSDM")) par.matrix=1;
584 else if (!strcmp(argv[i],"-BLOSUM50")) par.matrix=2;
585 else if (!strcmp(argv[i],"-Blosum50")) par.matrix=2;
586 else if (!strcmp(argv[i],"-B50")) par.matrix=2;
587 else if (!strcmp(argv[i],"-BLOSUM62")) par.matrix=3;
588 else if (!strcmp(argv[i],"-Blosum62")) par.matrix=3;
589 else if (!strcmp(argv[i],"-B62")) par.matrix=3;
590 else if (!strcmp(argv[i],"-pcm") && (i<argc-1)) par.pcm=atoi(argv[++i]);
591 else if (!strcmp(argv[i],"-pca") && (i<argc-1)) par.pca=atof(argv[++i]);
592 else if (!strcmp(argv[i],"-pcb") && (i<argc-1)) par.pcb=atof(argv[++i]);
593 else if (!strcmp(argv[i],"-pcc") && (i<argc-1)) par.pcc=atof(argv[++i]);
594 else if (!strcmp(argv[i],"-pcw") && (i<argc-1)) par.pcw=atof(argv[++i]);
595 else if (!strcmp(argv[i],"-gapb") && (i<argc-1)) { par.gapb=atof(argv[++i]); if (par.gapb<=0.01) par.gapb=0.01;}
596 else if (!strcmp(argv[i],"-gapd") && (i<argc-1)) par.gapd=atof(argv[++i]);
597 else if (!strcmp(argv[i],"-gape") && (i<argc-1)) par.gape=atof(argv[++i]);
598 else if (!strcmp(argv[i],"-gapf") && (i<argc-1)) par.gapf=atof(argv[++i]);
599 else if (!strcmp(argv[i],"-gapg") && (i<argc-1)) par.gapg=atof(argv[++i]);
600 else if (!strcmp(argv[i],"-gaph") && (i<argc-1)) par.gaph=atof(argv[++i]);
601 else if (!strcmp(argv[i],"-gapi") && (i<argc-1)) par.gapi=atof(argv[++i]);
602 else if (!strcmp(argv[i],"-egq") && (i<argc-1)) par.egq=atof(argv[++i]);
603 else if (!strcmp(argv[i],"-egt") && (i<argc-1)) par.egt=atof(argv[++i]);
604 else if (!strcmp(argv[i],"-ssgap")) par.ssgap=1;
605 else if (!strcmp(argv[i],"-ssgapd") && (i<argc-1)) par.ssgapd=atof(argv[++i]);
606 else if (!strcmp(argv[i],"-ssgape") && (i<argc-1)) par.ssgape=atof(argv[++i]);
607 else if (!strcmp(argv[i],"-ssgapi") && (i<argc-1)) par.ssgapi=atoi(argv[++i]);
608 else if (!strcmp(argv[i],"-ssm") && (i<argc-1)) par.ssm=atoi(argv[++i]);
609 else if (!strcmp(argv[i],"-ssw") && (i<argc-1)) par.ssw=atof(argv[++i]);
610 else if (!strcmp(argv[i],"-ssa") && (i<argc-1)) par.ssa=atof(argv[++i]);
611 else if (!strncmp(argv[i],"-glo",3)) {par.loc=0; if (par.mact>0.3 && par.mact<0.301) {par.mact=0;} }
612 else if (!strncmp(argv[i],"-loc",3)) par.loc=1;
613 else if (!strncmp(argv[i],"-alt",4) && (i<argc-1)) par.altali=atoi(argv[++i]);
614 else if (!strcmp(argv[i],"-map") || !strcmp(argv[i],"-MAP")) par.forward=2;
615 else if (!strcmp(argv[i],"-mac") || !strcmp(argv[i],"-MAC")) par.forward=2;
616 else if (!strcmp(argv[i],"-vit")) par.forward=0;
617 else if (!strcmp(argv[i],"-sto") && (i<argc-1)) {Nstochali=atoi(argv[++i]); par.forward=1;}
618 else if (!strcmp(argv[i],"-r")) par.repmode=1;
619 else if (!strcmp(argv[i],"-M") && (i<argc-1))
620 if (!strcmp(argv[++i],"a2m") || !strcmp(argv[i],"a3m")) par.M=1;
621 else if(!strcmp(argv[i],"first")) par.M=3;
622 else if (argv[i][0]>='0' && argv[i][0]<='9') {par.Mgaps=atoi(argv[i]); par.M=2;}
623 else cerr<<endl<<"WARNING: Ignoring unknown argument: -M "<<argv[i]<<"\n";
624 else if (!strcmp(argv[i],"-shift") && (i<argc-1)) par.shift=atof(argv[++i]);
625 else if (!strcmp(argv[i],"-mact") && (i<argc-1)) {par.mact=atof(argv[++i]); par.forward=2;}
626 else if (!strcmp(argv[i],"-wstruc") && (i<argc-1)) par.wstruc=atof(argv[++i]);
627 else if (!strcmp(argv[i],"-opt") && (i<argc-1)) par.opt=atoi(argv[++i]);
628 else if (!strcmp(argv[i],"-sc") && (i<argc-1)) par.columnscore=atoi(argv[++i]);
629 else if (!strcmp(argv[i],"-corr") && (i<argc-1)) par.corr=atof(argv[++i]);
630 else if (!strcmp(argv[i],"-def")) par.readdefaultsfile=1;
631 else if (!strcmp(argv[i],"-ovlp") && (i<argc-1)) par.min_overlap=atoi(argv[++i]);
632 else if (!strcmp(argv[i],"-tags")) par.notags=0;
633 else if (!strcmp(argv[i],"-notags")) par.notags=1;
634 else cerr<<endl<<"WARNING: Ignoring unknown option "<<argv[i]<<" ...\n";
635 if (v>=4) cout<<i<<" "<<argv[i]<<endl; //PRINT
636 } // end of for-loop for command line input
643 /////////////////////////////////////////////////////////////////////////////
645 /////////////////////////////////////////////////////////////////////////////
650 * @param[in,out] ppcFirstProf
651 * first profile to be aligned
652 * @param[in] iFirstCnt
653 * number of sequences in 1st profile
654 * @param[in,out] ppcSecndProf
655 * second profile to be aligned
656 * @param[in] iSecndCnt
657 * number of sequences in 2nd profile
658 * @param[out] dScore_p
661 * HMM info of external HMM (background)
662 * @param[in] pcPrealigned1
663 * (gapped) 1st sequence aligned to background HMM
664 * @param[in] pcRepresent1
665 * (gapped) sequence representing background HMM aligned to 1st sequence
666 * @param[in] pcPrealigned2
667 * (gapped) 2nd sequence aligned to background HMM
668 * @param[in] pcRepresent2
669 * (gapped) sequence representing background HMM aligned to 2nd sequence
670 * @param[in] rHhalignPara,
671 * various parameters passed down from commandline
674 * iFlag,zcAux,zcError are debugging arguments
676 * @return Non-zero on error
679 hhalign(char **ppcFirstProf, int iFirstCnt, double *pdWeightsL,
680 char **ppcSecndProf, int iSecndCnt, double *pdWeightsR,
681 double *dScore_p, hmm_light *prHMM,
682 char *pcPrealigned1, char *pcRepresent1,
683 char *pcPrealigned2, char *pcRepresent2,
684 hhalign_para rHhalignPara,
686 char zcAux[], char zcError[]) {
690 int iRetVal = RETURN_OK;
692 #ifndef CLUSTALO_NOFILE
695 argv = (char **)malloc(argc*sizeof(char *));
696 for (int i = 0; i < argc; i++){
697 argv[i] = (char *)malloc(100);
699 strcpy(argv[0], "./hhalign");
700 strcpy(argv[1], "-t");
701 strcpy(argv[2], "hhalign.C");
702 strcpy(argv[3], "-i");
703 strcpy(argv[4], "hhalign.C");
704 strcpy(argv[5], "-ofas");
705 strcpy(argv[6], "out");
708 /*char* argv_conf[MAXOPT];*/ /* Input arguments from .hhdefaults file
709 (first=1: argv_conf[0] is not used) */
710 /*int argc_conf;*/ // Number of arguments in argv_conf
711 /*char inext[IDLEN]="";*/ // Extension of query input file (hhm or a3m)
712 /*char text[IDLEN]="";*/ // Extension of template input file (hhm or a3m)
713 /*int** ali=NULL;*/ // ali[i][j]=1 if (i,j) is part of an alignment
714 /*int** alisto=NULL;*/ // ali[i][j]=1 if (i,j) is part of an alignment
715 int Nali; /* number of normally backtraced alignments
719 strcpy(par.tfile,""); /* FS, Initialise Argument */
720 strcpy(par.alnfile,"");
721 par.p=0.0 ; /* minimum threshold for inclusion in hit list
722 and alignment listing */
723 par.E=1e6; /* maximum threshold for inclusion in hit list
724 and alignment listing */
725 par.b=1; // min number of alignments
726 par.B=100; // max number of alignments
727 par.z=1; // min number of lines in hit list
728 par.Z=100; // max number of lines in hit list
729 par.append=0; /* append alignment to output file
731 par.altali=1; /* find only ONE (possibly overlapping)
733 par.hitrank=0; /* rank of hit to be printed as a3m alignment
735 par.outformat=3; // default output format for alignment is a3m
736 hit.self=0; // no self-alignment
737 par.forward=2; /* 0: Viterbi algorithm;
738 1: Viterbi+stochastic sampling;
739 2:Maximum Accuracy (MAC) algorithm */
741 // Make command line input globally available
742 #ifndef CLUSTALO_NOFILE
745 RemovePathAndExtension(program_name,argv[0]);
748 #ifndef CLUSTALO_NOFILE
749 /* Enable changing verbose mode before defaults file
750 and command line are processed */
751 for (int i=1; i<argc; i++)
753 if (!strcmp(argv[i],"-def")) par.readdefaultsfile=1;
754 else if (argc>1 && !strcmp(argv[i],"-v0")) v=0;
755 else if (argc>1 && !strcmp(argv[i],"-v1")) v=1;
756 else if (argc>2 && !strcmp(argv[i],"-v")) v=atoi(argv[i+1]);
759 // Read .hhdefaults file?
760 if (par.readdefaultsfile)
762 // Process default otpions from .hhconfig file
763 ReadDefaultsFile(argc_conf,argv_conf);
764 ProcessArguments(argc_conf,argv_conf);
767 /* Process command line options
768 (they override defaults from .hhdefaults file) */
769 #ifndef CLUSTALO_NOFILE
770 ProcessArguments(argc,argv);
775 int iAuxLen1 = strlen(ppcFirstProf[0]);
776 int iAuxLen2 = strlen(ppcSecndProf[0]);
777 if ( (0 == iAuxLen1) || (0 == iAuxLen2) ){ /* problem with empty profiles, FS, r249 -> r250 */
778 sprintf(zcError, "%s:%s:%d: strlen(prof1)=%d, strlen(prof2)=%d -- Nothing to align!\n",
779 __FUNCTION__, __FILE__, __LINE__, iAuxLen1, iAuxLen2);
780 iRetVal = RETURN_UNKNOWN;
781 /* Note: at this stage cannot do 'goto this_is_the_end;'
782 because would cross initialisation of several variables */
785 par.maxResLen = iAuxLen1 > iAuxLen2 ? iAuxLen1 : iAuxLen2;
786 const int ciGoodMeasureSeq = 10;
787 par.maxResLen += ciGoodMeasureSeq;
788 par.maxColCnt = iAuxLen1 + iAuxLen2 + ciGoodMeasureSeq;
789 //par.max_seqid=iFirstCnt+iSecndCnt+3; /* -id */
790 par.max_seqid=DEFAULT_FILTER; /* -id */
791 par.loc=0; par.mact=0; /* -glob */
792 par.nseqdis=iFirstCnt+iSecndCnt; /* -seq */
793 par.showcons=0; /* -nocons */
794 par.showdssp=0; /* -nodssp */
795 par.Mgaps=100; /* -M */
797 par.pdWg1=pdWeightsL; /* tree wg */
798 par.pdWg2=pdWeightsR; /* tree wg */
800 /* NOTE: *qali declared globally but only created here,
801 pass in number of sequences to get rid of statically
804 Alignment qali(iFirstCnt+iSecndCnt);
805 HMM q(iFirstCnt+iSecndCnt);
806 HMM t(iFirstCnt+iSecndCnt);
810 #ifndef CLUSTALO_NOFILE
811 // Check command line input and default values
813 {help(); cerr<<endl<<"Error in "<<program_name<<": no query alignment file given (-i file)\n"; exit(4);}
814 if (par.forward==2 && par.loc==0)
816 if (par.mact<0.301 || par.mact>0.300)
817 if (v>=1) fprintf(stderr,"REMARK: in -mac -global mode -mact is forced to 0\n");
819 // par.loc=1; // global forward-backward algorithm seems to work fine! use it in place of local version.
822 // Get rootname (no directory path, no extension) and extension of infile
823 RemoveExtension(q.file,par.infile);
824 RemoveExtension(t.file,par.tfile); /* FS, NOFILE2 (commented out) */
825 Extension(inext,par.infile);
826 Extension(text,par.tfile); /*FS, NOFILE2 (commented out) */
828 // Check option compatibilities
829 if (par.nseqdis>MAXSEQDIS-3-par.showcons) par.nseqdis=MAXSEQDIS-3-par.showcons; //3 reserved for secondary structure
830 if (par.aliwidth<20) par.aliwidth=20;
831 if (par.pca<0.001) par.pca=0.001; // to avoid log(0)
832 if (par.b>par.B) par.B=par.b;
833 if (par.z>par.Z) par.Z=par.z;
834 if (par.hitrank>0) par.altali=0;
839 cout<<"query file : "<<par.infile<<"\n";
840 cout<<"template file: "<<par.tfile<<"\n";/*FS, NOFILE2 (commented out) */
841 cout<<"Output file: "<<par.outfile<<"\n";
842 cout<<"Alignment file: "<<par.alnfile<<"\n";
847 // Set (global variable) substitution matrix and derived matrices
848 SetSubstitutionMatrix();
851 // Set secondary structure substitution matrix
852 SetSecStrucSubstitutionMatrix();
855 /* moved Viterbi switch from after RnP() to here,
856 switch after RnP() ineffectual as RnP decides log/lin of transition,
857 however log/lin of transitions depends on MAC/Viterbi,
861 qL = strlen(ppcFirstProf[0]);
867 tL = strlen(ppcSecndProf[0]);
874 // const float MEMSPACE_DYNPROG = 512*1024*1024;
875 /* determine amount of memory available for MAC on command-line; FS, r240 -> r241 */
876 const float MEMSPACE_DYNPROG = (double)1024*1024*rHhalignPara.iMacRamMB;
877 // longest allowable length of database HMM
878 int Lmaxmem=(int)((float)MEMSPACE_DYNPROG/qL/6/8);
879 if (par.forward==2 && tL+2>=Lmaxmem) {
881 cerr<<"WARNING: Not sufficient memory to realign with MAC algorithm. Using Viterbi algorithm."<<endl;
884 /* use different 'fudge' parameters for Viterbi, FS, r228 -> r229 */
890 par.gapb = par.gapbV;
891 par.gapd = par.gapdV;
892 par.gape = par.gapeV;
893 par.gapf = par.gapfV;
894 par.gapg = par.gapgV;
895 par.gaph = par.gaphV;
896 par.gapi = par.gapiV;
901 // Read input file (HMM, HHM, or alignment format), and add pseudocounts etc.
903 if (OK != ReadAndPrepare(INTERN_ALN_2_HMM,
904 ppcFirstProf, iFirstCnt, prHMM,
905 pcPrealigned1, pcRepresent1, pdWeightsL,
906 (char*)(""), q, &qali)) {
907 sprintf(zcError, "%s:%s:%d: Problem Reading/Preparing q-profile (len=%d)\n",
908 __FUNCTION__, __FILE__, __LINE__, qL);
909 iRetVal = RETURN_FROM_RNP;
910 goto this_is_the_end;
912 // Set query columns in His-tags etc to Null model distribution
913 if (par.notags) q.NeutralizeTags();
915 // Do self-comparison?
916 if (0 /* !*par.tfile // FS, 2010-03-10 */)
918 // Deep-copy q into t
921 // Find overlapping alternative alignments
924 // Read template alignment/HMM t and add pseudocounts
928 /* Read input file (HMM, HHM, or alignment format),
929 and add pseudocounts etc. */
931 if (OK != ReadAndPrepare(INTERN_ALN_2_HMM,
932 ppcSecndProf, iSecndCnt, prHMM,
933 pcPrealigned2, pcRepresent2, pdWeightsR,
935 sprintf(zcError, "%s:%s:%d: Problem Reading/Preparing t-profile (len=%d)\n",
936 __FUNCTION__, __FILE__, __LINE__, tL);
937 iRetVal = RETURN_FROM_RNP;
938 goto this_is_the_end;
942 // Factor Null model into HMM t
943 t.IncludeNullModelInHMM(q,t);
945 /* switch at this stage is ineffectual as log/lin already decided in RnP().
947 /*const float MEMSPACE_DYNPROG = 512*1024*1024;
948 // longest allowable length of database HMM
949 int Lmaxmem=(int)((float)MEMSPACE_DYNPROG/q.L/6/8);
950 if (par.forward==2 && t.L+2>=Lmaxmem)
953 cerr<<"WARNING: Not sufficient memory to realign with MAC algorithm. Using Viterbi algorithm."<<endl;
957 // Allocate memory for dynamic programming matrix
958 hit.AllocateBacktraceMatrix(q.L+2,t.L+2); // ...with a separate dynamic programming matrix (memory!!)
959 if (par.forward>=1 || Nstochali)
960 hit.AllocateForwardMatrix(q.L+2,t.L+2);
962 hit.AllocateBackwardMatrix(q.L+2,t.L+2);
964 #ifndef CLUSTALO_NOFILE
965 // Read structure file for Forward() function?
966 if (strucfile && par.wstruc>0)
969 Pstruc = new(float*[q.L+2]);
970 for (int i=0; i<q.L+2; i++) Pstruc[i] = new(float[t.L+2]);
971 Sstruc = new(float*[q.L+2]);
972 for (int i=0; i<q.L+2; i++) Sstruc[i] = new(float[t.L+2]);
974 if (strcmp(strucfile,"stdin"))
976 strucf = fopen(strucfile, "r");
977 if (!strucf) OpenFileError(strucfile);
982 if (v>=2) printf("Reading structure matrix from standard input ... (for UNIX use ^D for 'end-of-file')\n");
984 for (int i=1; i<=q.L; i++)
986 for (int j=1; j<=t.L; j++)
989 if (fscanf(strucf,"%f",&f) <=0 )
991 fprintf(stderr,"Error: too few numbers in file %s while reading line %i, column %i\n",strucfile,i,j);
995 Pstruc[i][j]=fmax(f,PMIN);
997 Pstruc[i][j]=fmax(pow(f,par.wstruc),PMIN);
998 // printf("%10.2E ",f);
999 Sstruc[i][j] = par.wstruc * log2(f);
1004 } /* (strucfile && par.wstruc>0) */
1007 /* Do (self-)comparison, store results if score>SMIN,
1008 and try next best alignment */
1009 /* FIXME very ambigous and possibly faulty if-else */
1011 if (par.forward==2) {
1012 printf("Using maximum accuracy (MAC) alignment algorithm ...\n");
1014 else if (par.forward==0) {
1015 printf("Using Viterbi algorithm ...\n");
1017 else if (par.forward==1) {
1018 printf("Using stochastic sampling algorithm ...\n");
1021 printf("\nWhat alignment algorithm are we using??\n");
1026 for (hit.irep=1; hit.irep<=imax(par.hitrank,par.altali); hit.irep++)
1030 // generate Viterbi alignment
1031 hit.Viterbi(q,t,Sstruc);
1034 else if (par.forward==1)
1036 // generate a single stochastically sampled alignment
1037 hit.Forward(q,t,Pstruc);
1038 srand( time(NULL) ); // initialize random generator
1039 hit.StochasticBacktrace(q,t);
1040 hitlist.Push(hit); /* insert hit at beginning of list
1041 (last repeats first!) */
1045 else if (par.forward==2)
1047 // generate forward alignment
1048 if (OK != hit.Forward(q,t,Pstruc)){
1049 fprintf(stderr, "%s:%s:%d: cannot complete hit.Forward\n",
1050 __FUNCTION__, __FILE__, __LINE__);
1051 iRetVal = RETURN_FROM_MAC; /* spot double overflow in Forward(). FS, r241 -> r243 */
1052 goto this_is_the_end;
1054 if (OK != hit.Backward(q,t)){
1055 fprintf(stderr, "%s:%s:%d: cannot complete hit.backward\n",
1056 __FUNCTION__, __FILE__, __LINE__);
1057 iRetVal = RETURN_FROM_MAC; /* spot double overflow in hit.Backward(). FS, r241 -> r243 */
1058 goto this_is_the_end;
1060 hit.MACAlignment(q,t);
1061 if ((isnan(hit.score)) || (isnan(hit.Pforward))){
1062 printf("nan after MAC\n");
1064 hit.BacktraceMAC(q,t);
1065 if ((isnan(hit.score)) || (isnan(hit.Pforward))){
1066 printf("nan after backtrace\n");
1068 } /* use MAC algorithm */
1069 // printf ("%-12.12s %-12.12s irep=%-2i score=%6.2f hit.Pvalt=%.2g\n",hit.name,hit.fam,hit.irep,hit.score,hit.Pvalt);
1070 *dScore_p = hit.score;
1072 if (hit.irep<=par.hitrank || hit.score>SMIN || (hit.Pvalt<pself && hit.score>0 ))
1073 hitlist.Push(hit); // insert hit at beginning of list (last repeats first!) and do next alignment
1076 if (hit.irep==1) hitlist.Push(hit); // first hit will be inserted into hitlist anyway, even if not significant
1077 hit = hitlist.ReadLast(); // last alignment was not significant => read last (significant) hit from list
1079 /* FIXME: memory leak during normal alignment: both push'es above are not freed:
1080 * valgrind --leak-check=full --show-reachable=yes ./src/clustalo -i ~/void/protein.fa -o /dev/null -v
1081 ==9506== 1,456 bytes in 1 blocks are still reachable in loss record 4 of 5
1082 ==9506== at 0x4C2726C: operator new(unsigned long) (vg_replace_malloc.c:230)
1083 ==9506== by 0x4618C5: List<Hit>::Push(Hit) (list-C.h:134)
1084 ==9506== by 0x45E87A: hhalign (hhalign.cpp:914)
1085 ==9506== by 0x413F8C: HhalignWrapper (hhalign_wrapper.c:405)
1086 ==9506== by 0x41A7B5: MyMain (mymain.c:676)
1087 ==9506== by 0x4031A6: main (main.cpp:37)
1090 ==9506== 4,442 (4,368 direct, 74 indirect) bytes in 3 blocks are definitely lost in loss record 5 of 5
1091 ==9506== at 0x4C2726C: operator new(unsigned long) (vg_replace_malloc.c:230)
1092 ==9506== by 0x4618C5: List<Hit>::Push(Hit) (list-C.h:134)
1093 ==9506== by 0x45E7C0: hhalign (hhalign.cpp:911)
1094 ==9506== by 0x413F8C: HhalignWrapper (hhalign_wrapper.c:405)
1095 ==9506== by 0x41A7B5: MyMain (mymain.c:676)
1096 ==9506== by 0x4031A6: main (main.cpp:37)
1099 } /* 1 <= hit.irep <= imax(par.hitrank,par.altali) */
1103 #ifndef CLUSTALO_NOFILE
1104 // Write posterior probability matrix as TCoffee library file
1107 if (v>=2) printf("Writing TCoffee library file to %s\n",tcfile);
1110 if (strcmp(tcfile,"stdout")) tcf = fopen(tcfile, "w"); else tcf = stdout;
1111 if (!tcf) OpenFileError(tcfile);
1112 fprintf(tcf,"! TC_LIB_FORMAT_01\n");
1113 fprintf(tcf,"%i\n",2); // two sequences in library file
1114 fprintf(tcf,"%s %i %s\n",q.name,q.L,q.seq[q.nfirst]+1);
1115 fprintf(tcf,"%s %i %s\n",hit.name,hit.L,hit.seq[hit.nfirst]+1);
1116 fprintf(tcf,"#1 2\n");
1117 for (i=1; i<=q.L; i++) // print all pairs (i,j) with probability above PROBTCMIN
1118 for (j=1; j<=t.L; j++)
1119 if (hit.B_MM[i][j]>probmin_tc)
1120 fprintf(tcf,"%5i %5i %5i\n",i,j,iround(100.0*hit.B_MM[i][j]));
1121 for (int step=hit.nsteps; step>=1; step--) // print all pairs on MAC alignment which were not yet printed
1123 i=hit.i[step]; j=hit.j[step];
1124 // printf("%5i %5i %5i %i\n",i,j,iround(100.0*hit.B_MM[i][j]),hit.states[step]);
1125 if (hit.states[step]>=MM && hit.B_MM[i][j]<=probmin_tc)
1126 fprintf(tcf,"%5i %5i %5i\n",i,j,iround(100.0*hit.B_MM[i][j]));
1130 fprintf(tcf,"! SEQ_1_TO_N\n");
1132 // for (i=1; i<=q.L; i++)
1135 // for (j=1; j<=t.L; j++) sum+=hit.B_MM[i][j];
1136 // printf("i=%-3i sum=%7.4f\n",i,sum);
1141 // Write last alignment into alitabfile
1145 if (strcmp(alitabfile,"stdout")) alitabf = fopen(alitabfile, "w"); else alitabf = stdout;
1146 if (!alitabf) OpenFileError(alitabfile);
1149 fprintf(alitabf," i j score SS probab\n");
1150 for (int step=hit.nsteps; step>=1; step--)
1151 if (hit.states[step]>=MM)
1152 fprintf(alitabf,"%5i %5i %6.2f %6.2f %7.4f\n",hit.i[step],hit.j[step],hit.S[step],hit.S_ss[step],hit.P_posterior[step]);
1156 fprintf(alitabf," i j score SS\n");
1157 for (int step=hit.nsteps; step>=1; step--)
1158 if (hit.states[step]>=MM)
1159 fprintf(alitabf,"%5i %5i %6.2f %6.2f\n",hit.i[step],hit.j[step],hit.S[step],hit.S_ss[step]);
1162 } /* if (alitabfile) */
1165 // Do Stochastic backtracing?
1167 for (int i=1; i<Nstochali; i++)
1169 hit.StochasticBacktrace(q,t);
1170 hitlist.Push(hit); //insert hit at beginning of list (last repeats first!)
1173 else // Set P-value, E-value and probability
1175 if (q.mu) hitlist.GetPvalsFromCalibration(q);
1176 else if (t.mu) hitlist.GetPvalsFromCalibration(t);
1179 // Print FASTA or A2M alignments?
1180 #ifndef CLUSTALO_NOFILE
1181 if (*par.pairwisealisfile)
1185 cout<<"Printing alignments in " <<
1186 (par.outformat==1? "FASTA" : par.outformat==2?"A2M" :"A3M") <<
1187 " format to "<<par.pairwisealisfile<<"\n";
1189 int iPrAliRtn = hitlist.PrintAlignments(
1191 ppcFirstProf, ppcSecndProf,
1193 q, par.pairwisealisfile, par.outformat);
1194 if (OK != iPrAliRtn){
1195 fprintf(stderr, "%s:%s:%d: Could not print alignments\n",
1196 __FUNCTION__, __FILE__, __LINE__);
1197 iRetVal = RETURN_FROM_PRINT_ALI; /* this is where mis-alignment was originally spotted,
1198 hope to trap it now earlier. FS, r241 -> r243 */
1199 goto this_is_the_end;
1201 } /* if (*par.pairwisealisfile) */
1203 #ifndef CLUSTALO_NOFILE
1204 // Print hit list and alignments
1207 hitlist.PrintHitList(q, par.outfile);
1208 hitlist.PrintAlignments(
1210 ppcFirstProf, ppcSecndProf,
1214 WriteToScreen(par.outfile,1000); //write only hit list to screen
1220 ///////////////////////////////////////////////////////////////////////
1221 #ifndef CLUSTALO_NOFILE
1222 // Show results for hit with rank par.hitrank
1223 if (par.hitrank==0) hit=hitlist.Read(1); else hit=hitlist.Read(par.hitrank);
1225 // Generate output alignment or HMM file?
1226 if (*par.alnfile || *par.psifile || *par.hhmfile)
1230 if (v>=2 && *par.alnfile) printf("Merging template to query alignment and writing resulting alignment in A3M format to %s...\n",par.alnfile);
1231 if (v>=2 && *par.psifile) printf("Merging template to query alignment and writing resulting alignment in PSI format to %s...\n",par.psifile);
1235 if (v>=2 && *par.alnfile) printf("Merging template to query alignment and appending template alignment in A3M format to %s...\n",par.alnfile);
1236 if (v>=2 && *par.psifile) printf("Merging template to query alignment and appending template alignment in PSI format to %s...\n",par.psifile);
1239 // Read query alignment into Qali
1240 Alignment Qali; // output A3M generated by merging A3M alignments for significant hits to the query alignment
1241 char qa3mfile[NAMELEN];
1242 RemoveExtension(qa3mfile,par.infile); // directory??
1243 strcat(qa3mfile,".a3m");
1244 FILE* qa3mf=fopen(qa3mfile,"r");
1245 if (!qa3mf) OpenFileError(qa3mfile);
1246 Qali.Read(qa3mf,qa3mfile);
1249 // Align query with template in master-slave mode
1250 Qali.MergeMasterSlave(hit,par.tfile); /*FS, NOFILE2 (commented out) */
1252 // Write output A3M alignment?
1254 Qali.WriteToFile(par.alnfile,"a3m");
1259 /* Convert ASCII to int (0-20),throw out all insert states,
1260 record their number in I[k][i] */
1261 Qali.Compress("merged A3M file");
1263 // Write output PSI-BLAST-formatted alignment?
1264 Qali.WriteToFile(par.psifile,"psi");
1266 } /* if (*par.alnfile || *par.psifile || *par.hhmfile) */
1267 #endif /* NOFILE2 */
1268 //////////////////////////////////////////////////////////////////////////
1273 // double log2Pvalue;
1274 // if (par.ssm && (par.ssm1 || par.ssm2))
1276 // log2Pvalue=hit.logPval/0.693147181+0.45*(4.0*par.ssw/0.15-hit.score_ss);
1278 // printf("Aligned %s with %s:\nApproximate P-value INCLUDING SS SCORE = %7.2g\n",q.name,t.name,pow(2.0,log2Pvalue));
1281 // printf("Aligned %s with %s:\nApproximate P-value (without SS score) = %7.2g\n",q.name,t.name,hit.Pval);
1287 printf("Aligned %s with %s: Score = %-7.2f P-value = %-7.2g\n",
1288 q.name,t.name,hit.score,hit.Pval);
1290 printf("Aligned %s with %s (rank %i): Score = %-7.2f P-value = %-7.2g\n",
1291 q.name,t.name,par.hitrank,hit.score,hit.Pval);
1297 // Delete memory for dynamic programming matrix
1298 hit.DeleteBacktraceMatrix(q.L+2);
1299 if (par.forward>=1 || Nstochali)
1300 hit.DeleteForwardMatrix(q.L+2);
1302 hit.DeleteBackwardMatrix(q.L+2);
1304 for (int i=0; i<q.L+2; i++) delete[](Pstruc[i]); delete[](Pstruc);}
1307 // Delete content of hits in hitlist
1309 while (!hitlist.End())
1310 hitlist.ReadNext().Delete(); // Delete content of hit object
1312 if (strucfile && par.wstruc>0)
1314 for (int i=0; i<q.L+2; i++){
1315 delete[] Pstruc[i]; Pstruc[i] = NULL;
1317 delete[] Pstruc; Pstruc = NULL;
1318 for (int i=0; i<q.L+2; i++){
1319 delete[] Sstruc[i]; Sstruc[i] = NULL;
1321 delete[] Sstruc; Sstruc = NULL;
1322 delete[] strucfile; strucfile = NULL;
1326 delete[] pngfile; pngfile = NULL;
1329 delete[] alitabfile; alitabfile = NULL;
1332 delete[] tcfile; tcfile = NULL;
1335 delete[] par.exclstr; par.exclstr = NULL;
1338 #ifndef CLUSTALO_NOFILE
1341 if (!strcmp(par.outfile,"stdout")) printf("Done!\n");
1346 outf=fopen(par.outfile,"a"); //open for append
1347 fprintf(outf,"Done!\n");
1350 if (v>=2) printf("Done\n");
1354 //qali.ClobberGlobal();
1355 hit.ClobberGlobal();
1362 hitlist.ClobberGlobal();
1366 } /* this is the end of hhalign() //end main */
1368 //////////////////////////////////////////////////////////////////////////////
1370 //////////////////////////////////////////////////////////////////////////////