gmxlib: clean up -Wunused-parameter warnings
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_do_dssp.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team,
6  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include "sysstuff.h"
43 #include "typedefs.h"
44 #include "string2.h"
45 #include "strdb.h"
46 #include "macros.h"
47 #include "smalloc.h"
48 #include "mshift.h"
49 #include "statutil.h"
50 #include "pdbio.h"
51 #include "gmx_fatal.h"
52 #include "xvgr.h"
53 #include "matio.h"
54 #include "index.h"
55 #include "gstat.h"
56 #include "tpxio.h"
57 #include "viewit.h"
58
59
60 static int strip_dssp(char *dsspfile, int nres,
61                       gmx_bool bPhobres[], real t,
62                       real *acc, FILE *fTArea,
63                       t_matrix *mat, int average_area[],
64                       const output_env_t oenv)
65 {
66     static gmx_bool bFirst = TRUE;
67     static char    *ssbuf;
68     FILE           *tapeout;
69     static int      xsize, frame;
70     char            buf[STRLEN+1];
71     char            SSTP;
72     int             i, nr, iacc, nresidues;
73     int             naccf, naccb; /* Count hydrophobic and hydrophilic residues */
74     real            iaccf, iaccb;
75     t_xpmelmt       c;
76
77     tapeout = ffopen(dsspfile, "r");
78
79     /* Skip header */
80     do
81     {
82         fgets2(buf, STRLEN, tapeout);
83     }
84     while (strstr(buf, "KAPPA") == NULL);
85     if (bFirst)
86     {
87         /* Since we also have empty lines in the dssp output (temp) file,
88          * and some line content is saved to the ssbuf variable,
89          * we need more memory than just nres elements. To be shure,
90          * we allocate 2*nres-1, since for each chain there is a
91          * separating line in the temp file. (At most each residue
92          * could have been defined as a separate chain.) */
93         snew(ssbuf, 2*nres-1);
94     }
95
96     iaccb     = iaccf = 0;
97     nresidues = 0;
98     naccf     = 0;
99     naccb     = 0;
100     for (nr = 0; (fgets2(buf, STRLEN, tapeout) != NULL); nr++)
101     {
102         if (buf[13] == '!') /* Chain separator line has '!' at pos. 13 */
103         {
104             SSTP = '=';     /* Chain separator sign '=' */
105         }
106         else
107         {
108             SSTP = buf[16] == ' ' ? '~' : buf[16];
109         }
110         ssbuf[nr] = SSTP;
111
112         buf[39] = '\0';
113
114         /* Only calculate solvent accessible area if needed */
115         if ((NULL != acc) && (buf[13] != '!'))
116         {
117             sscanf(&(buf[34]), "%d", &iacc);
118             acc[nr] = iacc;
119             /* average_area and bPhobres are counted from 0...nres-1 */
120             average_area[nresidues] += iacc;
121             if (bPhobres[nresidues])
122             {
123                 naccb++;
124                 iaccb += iacc;
125             }
126             else
127             {
128                 naccf++;
129                 iaccf += iacc;
130             }
131             /* Keep track of the residue number (do not count chain separator lines '!') */
132             nresidues++;
133         }
134     }
135     ssbuf[nr] = '\0';
136
137     if (bFirst)
138     {
139         if (0 != acc)
140         {
141             fprintf(stderr, "%d residues were classified as hydrophobic and %d as hydrophilic.\n", naccb, naccf);
142         }
143
144         sprintf(mat->title, "Secondary structure");
145         mat->legend[0] = 0;
146         sprintf(mat->label_x, "%s", output_env_get_time_label(oenv));
147         sprintf(mat->label_y, "Residue");
148         mat->bDiscrete = TRUE;
149         mat->ny        = nr;
150         snew(mat->axis_y, nr);
151         for (i = 0; i < nr; i++)
152         {
153             mat->axis_y[i] = i+1;
154         }
155         mat->axis_x = NULL;
156         mat->matrix = NULL;
157         xsize       = 0;
158         frame       = 0;
159         bFirst      = FALSE;
160     }
161     if (frame >= xsize)
162     {
163         xsize += 10;
164         srenew(mat->axis_x, xsize);
165         srenew(mat->matrix, xsize);
166     }
167     mat->axis_x[frame] = t;
168     snew(mat->matrix[frame], nr);
169     c.c2 = 0;
170     for (i = 0; i < nr; i++)
171     {
172         c.c1                  = ssbuf[i];
173         mat->matrix[frame][i] = max(0, searchcmap(mat->nmap, mat->map, c));
174     }
175     frame++;
176     mat->nx = frame;
177
178     if (fTArea)
179     {
180         fprintf(fTArea, "%10g  %10g  %10g\n", t, 0.01*iaccb, 0.01*iaccf);
181     }
182     ffclose(tapeout);
183
184     /* Return the number of lines found in the dssp file (i.e. number
185      * of redidues plus chain separator lines).
186      * This is the number of y elements needed for the area xpm file */
187     return nr;
188 }
189
190 static gmx_bool *bPhobics(t_atoms *atoms)
191 {
192     int         i, nb;
193     char      **cb;
194     gmx_bool   *bb;
195
196
197     nb = get_strings("phbres.dat", &cb);
198     snew(bb, atoms->nres);
199
200     for (i = 0; (i < atoms->nres); i++)
201     {
202         if (-1 != search_str(nb, cb, *atoms->resinfo[i].name) )
203         {
204             bb[i] = TRUE;
205         }
206     }
207     return bb;
208 }
209
210 static void check_oo(t_atoms *atoms)
211 {
212     char *OOO;
213
214     int   i;
215
216     OOO = strdup("O");
217
218     for (i = 0; (i < atoms->nr); i++)
219     {
220         if (strcmp(*(atoms->atomname[i]), "OXT") == 0)
221         {
222             *atoms->atomname[i] = OOO;
223         }
224         else if (strcmp(*(atoms->atomname[i]), "O1") == 0)
225         {
226             *atoms->atomname[i] = OOO;
227         }
228         else if (strcmp(*(atoms->atomname[i]), "OC1") == 0)
229         {
230             *atoms->atomname[i] = OOO;
231         }
232     }
233 }
234
235 static void norm_acc(t_atoms *atoms, int nres,
236                      real av_area[], real norm_av_area[])
237 {
238     int     i, n, n_surf;
239
240     char    surffn[] = "surface.dat";
241     char  **surf_res, **surf_lines;
242     double *surf;
243
244     n_surf = get_lines(surffn, &surf_lines);
245     snew(surf, n_surf);
246     snew(surf_res, n_surf);
247     for (i = 0; (i < n_surf); i++)
248     {
249         snew(surf_res[i], 5);
250         sscanf(surf_lines[i], "%s %lf", surf_res[i], &surf[i]);
251     }
252
253     for (i = 0; (i < nres); i++)
254     {
255         n = search_str(n_surf, surf_res, *atoms->resinfo[i].name);
256         if (n != -1)
257         {
258             norm_av_area[i] = av_area[i] / surf[n];
259         }
260         else
261         {
262             fprintf(stderr, "Residue %s not found in surface database (%s)\n",
263                     *atoms->resinfo[i].name, surffn);
264         }
265     }
266 }
267
268 void prune_ss_legend(t_matrix *mat)
269 {
270     gmx_bool  *present;
271     int       *newnum;
272     int        i, r, f, newnmap;
273     t_mapping *newmap;
274
275     snew(present, mat->nmap);
276     snew(newnum, mat->nmap);
277
278     for (f = 0; f < mat->nx; f++)
279     {
280         for (r = 0; r < mat->ny; r++)
281         {
282             present[mat->matrix[f][r]] = TRUE;
283         }
284     }
285
286     newnmap = 0;
287     newmap  = NULL;
288     for (i = 0; i < mat->nmap; i++)
289     {
290         newnum[i] = -1;
291         if (present[i])
292         {
293             newnum[i] = newnmap;
294             newnmap++;
295             srenew(newmap, newnmap);
296             newmap[newnmap-1] = mat->map[i];
297         }
298     }
299     if (newnmap != mat->nmap)
300     {
301         mat->nmap = newnmap;
302         mat->map  = newmap;
303         for (f = 0; f < mat->nx; f++)
304         {
305             for (r = 0; r < mat->ny; r++)
306             {
307                 mat->matrix[f][r] = newnum[mat->matrix[f][r]];
308             }
309         }
310     }
311 }
312
313 void write_sas_mat(const char *fn, real **accr, int nframe, int nres, t_matrix *mat)
314 {
315     real  lo, hi;
316     int   i, j, nlev;
317     t_rgb rlo = {1, 1, 1}, rhi = {0, 0, 0};
318     FILE *fp;
319
320     if (fn)
321     {
322         hi = lo = accr[0][0];
323         for (i = 0; i < nframe; i++)
324         {
325             for (j = 0; j < nres; j++)
326             {
327                 lo = min(lo, accr[i][j]);
328                 hi = max(hi, accr[i][j]);
329             }
330         }
331         fp   = ffopen(fn, "w");
332         nlev = hi-lo+1;
333         write_xpm(fp, 0, "Solvent Accessible Surface", "Surface (A^2)",
334                   "Time", "Residue Index", nframe, nres,
335                   mat->axis_x, mat->axis_y, accr, lo, hi, rlo, rhi, &nlev);
336         ffclose(fp);
337     }
338 }
339
340 void analyse_ss(const char *outfile, t_matrix *mat, const char *ss_string,
341                 const output_env_t oenv)
342 {
343     FILE        *fp;
344     t_mapping   *map;
345     int          f, r, *count, *total, ss_count, total_count;
346     size_t       s;
347     const char** leg;
348
349     map = mat->map;
350     snew(count, mat->nmap);
351     snew(total, mat->nmap);
352     snew(leg, mat->nmap+1);
353     leg[0] = "Structure";
354     for (s = 0; s < (size_t)mat->nmap; s++)
355     {
356         leg[s+1] = strdup(map[s].desc);
357     }
358
359     fp = xvgropen(outfile, "Secondary Structure",
360                   output_env_get_xvgr_tlabel(oenv), "Number of Residues", oenv);
361     if (output_env_get_print_xvgr_codes(oenv))
362     {
363         fprintf(fp, "@ subtitle \"Structure = ");
364     }
365     for (s = 0; s < strlen(ss_string); s++)
366     {
367         if (s > 0)
368         {
369             fprintf(fp, " + ");
370         }
371         for (f = 0; f < mat->nmap; f++)
372         {
373             if (ss_string[s] == map[f].code.c1)
374             {
375                 fprintf(fp, "%s", map[f].desc);
376             }
377         }
378     }
379     fprintf(fp, "\"\n");
380     xvgr_legend(fp, mat->nmap+1, leg, oenv);
381
382     total_count = 0;
383     for (s = 0; s < (size_t)mat->nmap; s++)
384     {
385         total[s] = 0;
386     }
387     for (f = 0; f < mat->nx; f++)
388     {
389         ss_count = 0;
390         for (s = 0; s < (size_t)mat->nmap; s++)
391         {
392             count[s] = 0;
393         }
394         for (r = 0; r < mat->ny; r++)
395         {
396             count[mat->matrix[f][r]]++;
397             total[mat->matrix[f][r]]++;
398         }
399         for (s = 0; s < (size_t)mat->nmap; s++)
400         {
401             if (strchr(ss_string, map[s].code.c1))
402             {
403                 ss_count    += count[s];
404                 total_count += count[s];
405             }
406         }
407         fprintf(fp, "%8g %5d", mat->axis_x[f], ss_count);
408         for (s = 0; s < (size_t)mat->nmap; s++)
409         {
410             fprintf(fp, " %5d", count[s]);
411         }
412         fprintf(fp, "\n");
413     }
414     /* now print column totals */
415     fprintf(fp, "%-8s %5d", "# Totals", total_count);
416     for (s = 0; s < (size_t)mat->nmap; s++)
417     {
418         fprintf(fp, " %5d", total[s]);
419     }
420     fprintf(fp, "\n");
421
422     /* now print percentages */
423     fprintf(fp, "%-8s %5.2f", "# SS %", total_count / (real) (mat->nx * mat->ny));
424     for (s = 0; s < (size_t)mat->nmap; s++)
425     {
426         fprintf(fp, " %5.2f", total[s] / (real) (mat->nx * mat->ny));
427     }
428     fprintf(fp, "\n");
429
430     ffclose(fp);
431     sfree(leg);
432     sfree(count);
433 }
434
435 int gmx_do_dssp(int argc, char *argv[])
436 {
437     const char        *desc[] = {
438         "[TT]do_dssp[tt] ",
439         "reads a trajectory file and computes the secondary structure for",
440         "each time frame ",
441         "calling the dssp program. If you do not have the dssp program,",
442         "get it from http://swift.cmbi.ru.nl/gv/dssp. [TT]do_dssp[tt] assumes ",
443         "that the dssp executable is located in ",
444         "[TT]/usr/local/bin/dssp[tt]. If this is not the case, then you should",
445         "set an environment variable [TT]DSSP[tt] pointing to the dssp",
446         "executable, e.g.: [PAR]",
447         "[TT]setenv DSSP /opt/dssp/bin/dssp[tt][PAR]",
448         "Since version 2.0.0, dssp is invoked with a syntax that differs",
449         "from earlier versions. If you have an older version of dssp,",
450         "use the [TT]-ver[tt] option to direct do_dssp to use the older syntax.",
451         "By default, do_dssp uses the syntax introduced with version 2.0.0.",
452         "Even newer versions (which at the time of writing are not yet released)",
453         "are assumed to have the same syntax as 2.0.0.[PAR]",
454         "The structure assignment for each residue and time is written to an",
455         "[TT].xpm[tt] matrix file. This file can be visualized with for instance",
456         "[TT]xv[tt] and can be converted to postscript with [TT]xpm2ps[tt].",
457         "Individual chains are separated by light grey lines in the [TT].xpm[tt] and",
458         "postscript files.",
459         "The number of residues with each secondary structure type and the",
460         "total secondary structure ([TT]-sss[tt]) count as a function of",
461         "time are also written to file ([TT]-sc[tt]).[PAR]",
462         "Solvent accessible surface (SAS) per residue can be calculated, both in",
463         "absolute values (A^2) and in fractions of the maximal accessible",
464         "surface of a residue. The maximal accessible surface is defined as",
465         "the accessible surface of a residue in a chain of glycines.",
466         "[BB]Note[bb] that the program [TT]g_sas[tt] can also compute SAS",
467         "and that is more efficient.[PAR]",
468         "Finally, this program can dump the secondary structure in a special file",
469         "[TT]ssdump.dat[tt] for usage in the program [TT]g_chi[tt]. Together",
470         "these two programs can be used to analyze dihedral properties as a",
471         "function of secondary structure type."
472     };
473     static gmx_bool    bVerbose;
474     static const char *ss_string   = "HEBT";
475     static int         dsspVersion = 2;
476     t_pargs            pa[]        = {
477         { "-v",  FALSE, etBOOL, {&bVerbose},
478           "HIDDENGenerate miles of useless information" },
479         { "-sss", FALSE, etSTR, {&ss_string},
480           "Secondary structures for structure count"},
481         { "-ver", FALSE, etINT, {&dsspVersion},
482           "DSSP major version. Syntax changed with version 2"}
483     };
484
485     t_trxstatus       *status;
486     FILE              *tapein;
487     FILE              *ss, *acc, *fTArea, *tmpf;
488     const char        *fnSCount, *fnArea, *fnTArea, *fnAArea;
489     const char        *leg[] = { "Phobic", "Phylic" };
490     t_topology         top;
491     int                ePBC;
492     t_atoms           *atoms;
493     t_matrix           mat;
494     int                nres, nr0, naccr, nres_plus_separators;
495     gmx_bool          *bPhbres, bDoAccSurf;
496     real               t;
497     int                i, j, natoms, nframe = 0;
498     matrix             box;
499     int                gnx;
500     char              *grpnm, *ss_str;
501     atom_id           *index;
502     rvec              *xp, *x;
503     int               *average_area;
504     real             **accr, *accr_ptr = NULL, *av_area, *norm_av_area;
505     char               pdbfile[32], tmpfile[32], title[256];
506     char               dssp[256];
507     const char        *dptr;
508     output_env_t       oenv;
509     gmx_rmpbc_t        gpbc = NULL;
510
511     t_filenm           fnm[] = {
512         { efTRX, "-f",   NULL,      ffREAD },
513         { efTPS, NULL,   NULL,      ffREAD },
514         { efNDX, NULL,   NULL,      ffOPTRD },
515         { efDAT, "-ssdump", "ssdump", ffOPTWR },
516         { efMAP, "-map", "ss",      ffLIBRD },
517         { efXPM, "-o",   "ss",      ffWRITE },
518         { efXVG, "-sc",  "scount",  ffWRITE },
519         { efXPM, "-a",   "area",    ffOPTWR },
520         { efXVG, "-ta",  "totarea", ffOPTWR },
521         { efXVG, "-aa",  "averarea", ffOPTWR }
522     };
523 #define NFILE asize(fnm)
524
525     parse_common_args(&argc, argv,
526                       PCA_CAN_TIME | PCA_CAN_VIEW | PCA_TIME_UNIT | PCA_BE_NICE,
527                       NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
528     fnSCount   = opt2fn("-sc", NFILE, fnm);
529     fnArea     = opt2fn_null("-a", NFILE, fnm);
530     fnTArea    = opt2fn_null("-ta", NFILE, fnm);
531     fnAArea    = opt2fn_null("-aa", NFILE, fnm);
532     bDoAccSurf = (fnArea || fnTArea || fnAArea);
533
534     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xp, NULL, box, FALSE);
535     atoms = &(top.atoms);
536     check_oo(atoms);
537     bPhbres = bPhobics(atoms);
538
539     get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpnm);
540     nres = 0;
541     nr0  = -1;
542     for (i = 0; (i < gnx); i++)
543     {
544         if (atoms->atom[index[i]].resind != nr0)
545         {
546             nr0 = atoms->atom[index[i]].resind;
547             nres++;
548         }
549     }
550     fprintf(stderr, "There are %d residues in your selected group\n", nres);
551
552     strcpy(pdbfile, "ddXXXXXX");
553     gmx_tmpnam(pdbfile);
554     if ((tmpf = fopen(pdbfile, "w")) == NULL)
555     {
556         sprintf(pdbfile, "%ctmp%cfilterXXXXXX", DIR_SEPARATOR, DIR_SEPARATOR);
557         gmx_tmpnam(pdbfile);
558         if ((tmpf = fopen(pdbfile, "w")) == NULL)
559         {
560             gmx_fatal(FARGS, "Can not open tmp file %s", pdbfile);
561         }
562     }
563     else
564     {
565         fclose(tmpf);
566     }
567
568     strcpy(tmpfile, "ddXXXXXX");
569     gmx_tmpnam(tmpfile);
570     if ((tmpf = fopen(tmpfile, "w")) == NULL)
571     {
572         sprintf(tmpfile, "%ctmp%cfilterXXXXXX", DIR_SEPARATOR, DIR_SEPARATOR);
573         gmx_tmpnam(tmpfile);
574         if ((tmpf = fopen(tmpfile, "w")) == NULL)
575         {
576             gmx_fatal(FARGS, "Can not open tmp file %s", tmpfile);
577         }
578     }
579     else
580     {
581         fclose(tmpf);
582     }
583
584     if ((dptr = getenv("DSSP")) == NULL)
585     {
586         dptr = "/usr/local/bin/dssp";
587     }
588     if (!gmx_fexist(dptr))
589     {
590         gmx_fatal(FARGS, "DSSP executable (%s) does not exist (use setenv DSSP)",
591                   dptr);
592     }
593     if (dsspVersion >= 2)
594     {
595         if (dsspVersion > 2)
596         {
597             printf("\nWARNING: You use DSSP version %d, which is not explicitly\nsupported by do_dssp. Assuming version 2 syntax.\n\n", dsspVersion);
598         }
599
600         sprintf(dssp, "%s -i %s -o %s > /dev/null %s",
601                 dptr, pdbfile, tmpfile, bVerbose ? "" : "2> /dev/null");
602     }
603     else
604     {
605         sprintf(dssp, "%s %s %s %s > /dev/null %s",
606                 dptr, bDoAccSurf ? "" : "-na", pdbfile, tmpfile, bVerbose ? "" : "2> /dev/null");
607
608     }
609     fprintf(stderr, "dssp cmd='%s'\n", dssp);
610
611     if (fnTArea)
612     {
613         fTArea = xvgropen(fnTArea, "Solvent Accessible Surface Area",
614                           output_env_get_xvgr_tlabel(oenv), "Area (nm\\S2\\N)", oenv);
615         xvgr_legend(fTArea, 2, leg, oenv);
616     }
617     else
618     {
619         fTArea = NULL;
620     }
621
622     mat.map  = NULL;
623     mat.nmap = getcmap(libopen(opt2fn("-map", NFILE, fnm)),
624                        opt2fn("-map", NFILE, fnm), &(mat.map));
625
626     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
627     if (natoms > atoms->nr)
628     {
629         gmx_fatal(FARGS, "\nTrajectory does not match topology!");
630     }
631     if (gnx > natoms)
632     {
633         gmx_fatal(FARGS, "\nTrajectory does not match selected group!");
634     }
635
636     snew(average_area, atoms->nres);
637     snew(av_area, atoms->nres);
638     snew(norm_av_area, atoms->nres);
639     accr  = NULL;
640     naccr = 0;
641
642     gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
643     do
644     {
645         t = output_env_conv_time(oenv, t);
646         if (bDoAccSurf && nframe >= naccr)
647         {
648             naccr += 10;
649             srenew(accr, naccr);
650             for (i = naccr-10; i < naccr; i++)
651             {
652                 snew(accr[i], 2*atoms->nres-1);
653             }
654         }
655         gmx_rmpbc(gpbc, natoms, box, x);
656         tapein = ffopen(pdbfile, "w");
657         write_pdbfile_indexed(tapein, NULL, atoms, x, ePBC, box, ' ', -1, gnx, index, NULL, TRUE);
658         ffclose(tapein);
659
660 #ifdef GMX_NO_SYSTEM
661         printf("Warning-- No calls to system(3) supported on this platform.");
662         printf("Warning-- Skipping execution of 'system(\"%s\")'.", dssp);
663         exit(1);
664 #else
665         if (0 != system(dssp))
666         {
667             gmx_fatal(FARGS, "Failed to execute command: %s\n",
668                       "Try specifying your dssp version with the -ver option.", dssp);
669         }
670 #endif
671
672         /* strip_dssp returns the number of lines found in the dssp file, i.e.
673          * the number of residues plus the separator lines */
674
675         if (bDoAccSurf)
676         {
677             accr_ptr = accr[nframe];
678         }
679
680         nres_plus_separators = strip_dssp(tmpfile, nres, bPhbres, t,
681                                           accr_ptr, fTArea, &mat, average_area, oenv);
682         remove(tmpfile);
683         remove(pdbfile);
684         nframe++;
685     }
686     while (read_next_x(oenv, status, &t, x, box));
687     fprintf(stderr, "\n");
688     close_trj(status);
689     if (fTArea)
690     {
691         ffclose(fTArea);
692     }
693     gmx_rmpbc_done(gpbc);
694
695     prune_ss_legend(&mat);
696
697     ss        = opt2FILE("-o", NFILE, fnm, "w");
698     mat.flags = 0;
699     write_xpm_m(ss, mat);
700     ffclose(ss);
701
702     if (opt2bSet("-ssdump", NFILE, fnm))
703     {
704         ss = opt2FILE("-ssdump", NFILE, fnm, "w");
705         snew(ss_str, nres+1);
706         fprintf(ss, "%d\n", nres);
707         for (j = 0; j < mat.nx; j++)
708         {
709             for (i = 0; (i < mat.ny); i++)
710             {
711                 ss_str[i] = mat.map[mat.matrix[j][i]].code.c1;
712             }
713             ss_str[i] = '\0';
714             fprintf(ss, "%s\n", ss_str);
715         }
716         ffclose(ss);
717         sfree(ss_str);
718     }
719     analyse_ss(fnSCount, &mat, ss_string, oenv);
720
721     if (bDoAccSurf)
722     {
723         write_sas_mat(fnArea, accr, nframe, nres_plus_separators, &mat);
724
725         for (i = 0; i < atoms->nres; i++)
726         {
727             av_area[i] = (average_area[i] / (real)nframe);
728         }
729
730         norm_acc(atoms, nres, av_area, norm_av_area);
731
732         if (fnAArea)
733         {
734             acc = xvgropen(fnAArea, "Average Accessible Area",
735                            "Residue", "A\\S2", oenv);
736             for (i = 0; (i < nres); i++)
737             {
738                 fprintf(acc, "%5d  %10g %10g\n", i+1, av_area[i], norm_av_area[i]);
739             }
740             ffclose(acc);
741         }
742     }
743
744     view_all(oenv, NFILE, fnm);
745
746     return 0;
747 }