Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / fileio / confio.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  * Copyright (c) 2013, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42 #include "sysstuff.h"
43 #include "typedefs.h"
44 #include "smalloc.h"
45 #include "sysstuff.h"
46 #include <errno.h>
47 #include "macros.h"
48 #include "string2.h"
49 #include "confio.h"
50 #include "vec.h"
51 #include "symtab.h"
52 #include "futil.h"
53 #include "xdrf.h"
54 #include "filenm.h"
55 #include "pdbio.h"
56 #include "tpxio.h"
57 #include "trxio.h"
58 #include "gmx_fatal.h"
59 #include "copyrite.h"
60 #include "statutil.h"
61 #include "pbc.h"
62 #include "mtop_util.h"
63 #include "gmxfio.h"
64
65 #define CHAR_SHIFT 24
66
67 static int read_g96_pos(char line[], t_symtab *symtab,
68                         FILE *fp, const char *infile,
69                         t_trxframe *fr)
70 {
71     t_atoms   *atoms;
72     gmx_bool   bEnd;
73     int        nwanted, natoms, atnr, resnr = 0, oldres, newres, shift;
74     char       anm[STRLEN], resnm[STRLEN];
75     char       c1, c2;
76     double     db1, db2, db3;
77
78     nwanted = fr->natoms;
79
80     atoms = fr->atoms;
81
82     natoms = 0;
83
84     if (fr->bX)
85     {
86         if (fr->bAtoms)
87         {
88             shift = CHAR_SHIFT;
89         }
90         else
91         {
92             shift = 0;
93         }
94         newres  = -1;
95         oldres  = -666; /* Unlikely number for the first residue! */
96         bEnd    = FALSE;
97         while (!bEnd && fgets2(line, STRLEN, fp))
98         {
99             bEnd = (strncmp(line, "END", 3) == 0);
100             if (!bEnd  && (line[0] != '#'))
101             {
102                 if (sscanf(line+shift, "%15lf%15lf%15lf", &db1, &db2, &db3) != 3)
103                 {
104                     gmx_fatal(FARGS, "Did not find 3 coordinates for atom %d in %s\n",
105                               natoms+1, infile);
106                 }
107                 if ((nwanted != -1) && (natoms >= nwanted))
108                 {
109                     gmx_fatal(FARGS,
110                               "Found more coordinates (%d) in %s than expected %d\n",
111                               natoms, infile, nwanted);
112                 }
113                 if (atoms)
114                 {
115                     if (fr->bAtoms &&
116                         (sscanf(line, "%5d%c%5s%c%5s%7d", &resnr, &c1, resnm, &c2, anm, &atnr)
117                          != 6))
118                     {
119                         if (oldres >= 0)
120                         {
121                             resnr = oldres;
122                         }
123                         else
124                         {
125                             resnr    = 1;
126                             strncpy(resnm, "???", sizeof(resnm)-1);
127                         }
128                         strncpy(anm, "???", sizeof(anm)-1);
129                     }
130                     atoms->atomname[natoms] = put_symtab(symtab, anm);
131                     if (resnr != oldres)
132                     {
133                         oldres = resnr;
134                         newres++;
135                         if (newres >= atoms->nr)
136                         {
137                             gmx_fatal(FARGS, "More residues than atoms in %s (natoms = %d)",
138                                       infile, atoms->nr);
139                         }
140                         atoms->atom[natoms].resind = newres;
141                         if (newres+1 > atoms->nres)
142                         {
143                             atoms->nres = newres+1;
144                         }
145                         t_atoms_set_resinfo(atoms, natoms, symtab, resnm, resnr, ' ', 0, ' ');
146                     }
147                     else
148                     {
149                         atoms->atom[natoms].resind = newres;
150                     }
151                 }
152                 if (fr->x)
153                 {
154                     fr->x[natoms][0] = db1;
155                     fr->x[natoms][1] = db2;
156                     fr->x[natoms][2] = db3;
157                 }
158                 natoms++;
159             }
160         }
161         if ((nwanted != -1) && natoms != nwanted)
162         {
163             fprintf(stderr,
164                     "Warning: found less coordinates (%d) in %s than expected %d\n",
165                     natoms, infile, nwanted);
166         }
167     }
168
169     fr->natoms = natoms;
170
171     return natoms;
172 }
173
174 static int read_g96_vel(char line[], FILE *fp, const char *infile,
175                         t_trxframe *fr)
176 {
177     gmx_bool   bEnd;
178     int        nwanted, natoms = -1, shift;
179     double     db1, db2, db3;
180
181     nwanted = fr->natoms;
182
183     if (fr->v && fr->bV)
184     {
185         if (strcmp(line, "VELOCITYRED") == 0)
186         {
187             shift = 0;
188         }
189         else
190         {
191             shift = CHAR_SHIFT;
192         }
193         natoms = 0;
194         bEnd   = FALSE;
195         while (!bEnd && fgets2(line, STRLEN, fp))
196         {
197             bEnd = (strncmp(line, "END", 3) == 0);
198             if (!bEnd && (line[0] != '#'))
199             {
200                 if (sscanf(line+shift, "%15lf%15lf%15lf", &db1, &db2, &db3) != 3)
201                 {
202                     gmx_fatal(FARGS, "Did not find 3 velocities for atom %d in %s",
203                               natoms+1, infile);
204                 }
205                 if ((nwanted != -1) && (natoms >= nwanted))
206                 {
207                     gmx_fatal(FARGS, "Found more velocities (%d) in %s than expected %d\n",
208                               natoms, infile, nwanted);
209                 }
210                 if (fr->v)
211                 {
212                     fr->v[natoms][0] = db1;
213                     fr->v[natoms][1] = db2;
214                     fr->v[natoms][2] = db3;
215                 }
216                 natoms++;
217             }
218         }
219         if ((nwanted != -1) && (natoms != nwanted))
220         {
221             fprintf(stderr,
222                     "Warning: found less velocities (%d) in %s than expected %d\n",
223                     natoms, infile, nwanted);
224         }
225     }
226
227     return natoms;
228 }
229
230 int read_g96_conf(FILE *fp, const char *infile, t_trxframe *fr, char *line)
231 {
232     t_symtab  *symtab = NULL;
233     gmx_bool   bAtStart, bTime, bAtoms, bPos, bVel, bBox, bEnd, bFinished;
234     int        natoms, nbp;
235     double     db1, db2, db3, db4, db5, db6, db7, db8, db9;
236
237     bAtStart = (ftell(fp) == 0);
238
239     clear_trxframe(fr, FALSE);
240
241     if (!symtab)
242     {
243         snew(symtab, 1);
244         open_symtab(symtab);
245     }
246
247     natoms = 0;
248
249     if (bAtStart)
250     {
251         while (!fr->bTitle && fgets2(line, STRLEN, fp))
252         {
253             fr->bTitle = (strcmp(line, "TITLE") == 0);
254         }
255         if (fr->title == NULL)
256         {
257             fgets2(line, STRLEN, fp);
258             fr->title = strdup(line);
259         }
260         bEnd = FALSE;
261         while (!bEnd && fgets2(line, STRLEN, fp))
262         {
263             bEnd = (strcmp(line, "END") == 0);
264         }
265         fgets2(line, STRLEN, fp);
266     }
267
268     /* Do not get a line if we are not at the start of the file, *
269      * because without a parameter file we don't know what is in *
270      * the trajectory and we have already read the line in the   *
271      * previous call (VERY DIRTY).                               */
272     bFinished = FALSE;
273     do
274     {
275         bTime  = (strcmp(line, "TIMESTEP") == 0);
276         bAtoms = (strcmp(line, "POSITION") == 0);
277         bPos   = (bAtoms || (strcmp(line, "POSITIONRED") == 0));
278         bVel   = (strncmp(line, "VELOCITY", 8) == 0);
279         bBox   = (strcmp(line, "BOX") == 0);
280         if (bTime)
281         {
282             if (!fr->bTime && !fr->bX)
283             {
284                 fr->bStep = bTime;
285                 fr->bTime = bTime;
286                 do
287                 {
288                     bFinished = (fgets2(line, STRLEN, fp) == NULL);
289                 }
290                 while (!bFinished && (line[0] == '#'));
291                 sscanf(line, "%15d%15lf", &(fr->step), &db1);
292                 fr->time = db1;
293             }
294             else
295             {
296                 bFinished = TRUE;
297             }
298         }
299         if (bPos)
300         {
301             if (!fr->bX)
302             {
303                 fr->bAtoms = bAtoms;
304                 fr->bX     = bPos;
305                 natoms     = read_g96_pos(line, symtab, fp, infile, fr);
306             }
307             else
308             {
309                 bFinished = TRUE;
310             }
311         }
312         if (fr->v && bVel)
313         {
314             fr->bV = bVel;
315             natoms = read_g96_vel(line, fp, infile, fr);
316         }
317         if (bBox)
318         {
319             fr->bBox = bBox;
320             clear_mat(fr->box);
321             bEnd = FALSE;
322             while (!bEnd && fgets2(line, STRLEN, fp))
323             {
324                 bEnd = (strncmp(line, "END", 3) == 0);
325                 if (!bEnd && (line[0] != '#'))
326                 {
327                     nbp = sscanf(line, "%15lf%15lf%15lf%15lf%15lf%15lf%15lf%15lf%15lf",
328                                  &db1, &db2, &db3, &db4, &db5, &db6, &db7, &db8, &db9);
329                     if (nbp < 3)
330                     {
331                         gmx_fatal(FARGS, "Found a BOX line, but no box in %s", infile);
332                     }
333                     fr->box[XX][XX] = db1;
334                     fr->box[YY][YY] = db2;
335                     fr->box[ZZ][ZZ] = db3;
336                     if (nbp == 9)
337                     {
338                         fr->box[XX][YY] = db4;
339                         fr->box[XX][ZZ] = db5;
340                         fr->box[YY][XX] = db6;
341                         fr->box[YY][ZZ] = db7;
342                         fr->box[ZZ][XX] = db8;
343                         fr->box[ZZ][YY] = db9;
344                     }
345                 }
346             }
347             bFinished = TRUE;
348         }
349     }
350     while (!bFinished && fgets2(line, STRLEN, fp));
351
352     free_symtab(symtab);
353
354     fr->natoms = natoms;
355
356     return natoms;
357 }
358
359 void write_g96_conf(FILE *out, t_trxframe *fr,
360                     int nindex, const atom_id *index)
361 {
362     t_atoms *atoms;
363     int      nout, i, a;
364
365     atoms = fr->atoms;
366
367     if (index)
368     {
369         nout = nindex;
370     }
371     else
372     {
373         nout = fr->natoms;
374     }
375
376     if (fr->bTitle)
377     {
378         fprintf(out, "TITLE\n%s\nEND\n", fr->title);
379     }
380     if (fr->bStep || fr->bTime)
381     {
382         /* Officially the time format is %15.9, which is not enough for 10 ns */
383         fprintf(out, "TIMESTEP\n%15d%15.6f\nEND\n", fr->step, fr->time);
384     }
385     if (fr->bX)
386     {
387         if (fr->bAtoms)
388         {
389             fprintf(out, "POSITION\n");
390             for (i = 0; i < nout; i++)
391             {
392                 if (index)
393                 {
394                     a = index[i];
395                 }
396                 else
397                 {
398                     a = i;
399                 }
400                 fprintf(out, "%5d %-5s %-5s%7d%15.9f%15.9f%15.9f\n",
401                         (atoms->resinfo[atoms->atom[a].resind].nr) % 100000,
402                         *atoms->resinfo[atoms->atom[a].resind].name,
403                         *atoms->atomname[a], (i+1) % 10000000,
404                         fr->x[a][XX], fr->x[a][YY], fr->x[a][ZZ]);
405             }
406         }
407         else
408         {
409             fprintf(out, "POSITIONRED\n");
410             for (i = 0; i < nout; i++)
411             {
412                 if (index)
413                 {
414                     a = index[i];
415                 }
416                 else
417                 {
418                     a = i;
419                 }
420                 fprintf(out, "%15.9f%15.9f%15.9f\n",
421                         fr->x[a][XX], fr->x[a][YY], fr->x[a][ZZ]);
422             }
423         }
424         fprintf(out, "END\n");
425     }
426     if (fr->bV)
427     {
428         if (fr->bAtoms)
429         {
430             fprintf(out, "VELOCITY\n");
431             for (i = 0; i < nout; i++)
432             {
433                 if (index)
434                 {
435                     a = index[i];
436                 }
437                 else
438                 {
439                     a = i;
440                 }
441                 fprintf(out, "%5d %-5s %-5s%7d%15.9f%15.9f%15.9f\n",
442                         (atoms->resinfo[atoms->atom[a].resind].nr) % 100000,
443                         *atoms->resinfo[atoms->atom[a].resind].name,
444                         *atoms->atomname[a], (i+1) % 10000000,
445                         fr->v[a][XX], fr->v[a][YY], fr->v[a][ZZ]);
446             }
447         }
448         else
449         {
450             fprintf(out, "VELOCITYRED\n");
451             for (i = 0; i < nout; i++)
452             {
453                 if (index)
454                 {
455                     a = index[i];
456                 }
457                 else
458                 {
459                     a = i;
460                 }
461                 fprintf(out, "%15.9f%15.9f%15.9f\n",
462                         fr->v[a][XX], fr->v[a][YY], fr->v[a][ZZ]);
463             }
464         }
465         fprintf(out, "END\n");
466     }
467     if (fr->bBox)
468     {
469         fprintf(out, "BOX\n");
470         fprintf(out, "%15.9f%15.9f%15.9f",
471                 fr->box[XX][XX], fr->box[YY][YY], fr->box[ZZ][ZZ]);
472         if (fr->box[XX][YY] || fr->box[XX][ZZ] || fr->box[YY][XX] ||
473             fr->box[YY][ZZ] || fr->box[ZZ][XX] || fr->box[ZZ][YY])
474         {
475             fprintf(out, "%15.9f%15.9f%15.9f%15.9f%15.9f%15.9f",
476                     fr->box[XX][YY], fr->box[XX][ZZ], fr->box[YY][XX],
477                     fr->box[YY][ZZ], fr->box[ZZ][XX], fr->box[ZZ][YY]);
478         }
479         fprintf(out, "\n");
480         fprintf(out, "END\n");
481     }
482 }
483
484 static int get_espresso_word(FILE *fp, char word[])
485 {
486     int  ret, nc, i;
487
488     ret = 0;
489     nc  = 0;
490
491     do
492     {
493         i = fgetc(fp);
494         if (i != EOF)
495         {
496             if (i == ' ' || i == '\n' || i == '\t')
497             {
498                 if (nc > 0)
499                 {
500                     ret = 1;
501                 }
502             }
503             else if (i == '{')
504             {
505                 if (nc == 0)
506                 {
507                     word[nc++] = '{';
508                 }
509                 ret = 2;
510             }
511             else if (i == '}')
512             {
513                 if (nc == 0)
514                 {
515                     word[nc++] = '}';
516                 }
517                 ret = 3;
518             }
519             else
520             {
521                 word[nc++] = (char)i;
522             }
523         }
524     }
525     while (i != EOF && ret == 0);
526
527     word[nc] = '\0';
528
529     /*  printf("'%s'\n",word); */
530
531     return ret;
532 }
533
534 static int check_open_parenthesis(FILE *fp, int r,
535                                   const char *infile, const char *keyword)
536 {
537     int  level_inc;
538     char word[STRLEN];
539
540     level_inc = 0;
541     if (r == 2)
542     {
543         level_inc++;
544     }
545     else
546     {
547         r = get_espresso_word(fp, word);
548         if (r == 2)
549         {
550             level_inc++;
551         }
552         else
553         {
554             gmx_fatal(FARGS, "Expected '{' after '%s' in file '%s'",
555                       keyword, infile);
556         }
557     }
558
559     return level_inc;
560 }
561
562 static int check_close_parenthesis(FILE *fp, int r,
563                                    const char *infile, const char *keyword)
564 {
565     int  level_inc;
566     char word[STRLEN];
567
568     level_inc = 0;
569     if (r == 3)
570     {
571         level_inc--;
572     }
573     else
574     {
575         r = get_espresso_word(fp, word);
576         if (r == 3)
577         {
578             level_inc--;
579         }
580         else
581         {
582             gmx_fatal(FARGS, "Expected '}' after section '%s' in file '%s'",
583                       keyword, infile);
584         }
585     }
586
587     return level_inc;
588 }
589
590 enum {
591     espID, espPOS, espTYPE, espQ, espV, espF, espMOLECULE, espNR
592 };
593 const char *esp_prop[espNR] = {
594     "id", "pos", "type", "q", "v", "f",
595     "molecule"
596 };
597
598 static void read_espresso_conf(const char *infile,
599                                t_atoms *atoms, rvec x[], rvec *v, matrix box)
600 {
601     t_symtab *symtab = NULL;
602     FILE     *fp;
603     char      word[STRLEN], buf[STRLEN];
604     int       natoms, level, npar, r, nprop, p, i, m, molnr;
605     int       prop[32];
606     double    d;
607     gmx_bool  bFoundParticles, bFoundProp, bFoundVariable, bMol;
608
609     if (!symtab)
610     {
611         snew(symtab, 1);
612         open_symtab(symtab);
613     }
614
615     clear_mat(box);
616
617     fp = gmx_fio_fopen(infile, "r");
618
619     bFoundParticles = FALSE;
620     bFoundVariable  = FALSE;
621     bMol            = FALSE;
622     level           = 0;
623     while ((r = get_espresso_word(fp, word)))
624     {
625         if (level == 1 && strcmp(word, "particles") == 0 && !bFoundParticles)
626         {
627             bFoundParticles = TRUE;
628             level          += check_open_parenthesis(fp, r, infile, "particles");
629             nprop           = 0;
630             while (level == 2 && (r = get_espresso_word(fp, word)))
631             {
632                 bFoundProp = FALSE;
633                 for (p = 0; p < espNR; p++)
634                 {
635                     if (strcmp(word, esp_prop[p]) == 0)
636                     {
637                         bFoundProp    = TRUE;
638                         prop[nprop++] = p;
639                         /* printf("  prop[%d] = %s\n",nprop-1,esp_prop[prop[nprop-1]]); */
640                     }
641                 }
642                 if (!bFoundProp && word[0] != '}')
643                 {
644                     gmx_fatal(FARGS, "Can not read Espresso files with particle property '%s'", word);
645                 }
646                 if (bFoundProp && p == espMOLECULE)
647                 {
648                     bMol = TRUE;
649                 }
650                 if (r == 3)
651                 {
652                     level--;
653                 }
654             }
655
656             i = 0;
657             while (level > 0 && (r = get_espresso_word(fp, word)))
658             {
659                 if (r == 2)
660                 {
661                     level++;
662                 }
663                 else if (r == 3)
664                 {
665                     level--;
666                 }
667                 if (level == 2)
668                 {
669                     for (p = 0; p < nprop; p++)
670                     {
671                         switch (prop[p])
672                         {
673                             case espID:
674                                 r = get_espresso_word(fp, word);
675                                 /* Not used */
676                                 break;
677                             case espPOS:
678                                 for (m = 0; m < 3; m++)
679                                 {
680                                     r = get_espresso_word(fp, word);
681                                     sscanf(word, "%lf", &d);
682                                     x[i][m] = d;
683                                 }
684                                 break;
685                             case espTYPE:
686                                 r                   = get_espresso_word(fp, word);
687                                 atoms->atom[i].type = strtol(word, NULL, 10);
688                                 break;
689                             case espQ:
690                                 r = get_espresso_word(fp, word);
691                                 sscanf(word, "%lf", &d);
692                                 atoms->atom[i].q = d;
693                                 break;
694                             case espV:
695                                 for (m = 0; m < 3; m++)
696                                 {
697                                     r = get_espresso_word(fp, word);
698                                     sscanf(word, "%lf", &d);
699                                     v[i][m] = d;
700                                 }
701                                 break;
702                             case espF:
703                                 for (m = 0; m < 3; m++)
704                                 {
705                                     r = get_espresso_word(fp, word);
706                                     /* not used */
707                                 }
708                                 break;
709                             case espMOLECULE:
710                                 r     = get_espresso_word(fp, word);
711                                 molnr = strtol(word, NULL, 10);
712                                 if (i == 0 ||
713                                     atoms->resinfo[atoms->atom[i-1].resind].nr != molnr)
714                                 {
715                                     atoms->atom[i].resind =
716                                         (i == 0 ? 0 : atoms->atom[i-1].resind+1);
717                                     atoms->resinfo[atoms->atom[i].resind].nr       = molnr;
718                                     atoms->resinfo[atoms->atom[i].resind].ic       = ' ';
719                                     atoms->resinfo[atoms->atom[i].resind].chainid  = ' ';
720                                     atoms->resinfo[atoms->atom[i].resind].chainnum = molnr; /* Not sure if this is right? */
721                                 }
722                                 else
723                                 {
724                                     atoms->atom[i].resind = atoms->atom[i-1].resind;
725                                 }
726                                 break;
727                         }
728                     }
729                     /* Generate an atom name from the particle type */
730                     sprintf(buf, "T%d", atoms->atom[i].type);
731                     atoms->atomname[i] = put_symtab(symtab, buf);
732                     if (bMol)
733                     {
734                         if (i == 0 || atoms->atom[i].resind != atoms->atom[i-1].resind)
735                         {
736                             atoms->resinfo[atoms->atom[i].resind].name =
737                                 put_symtab(symtab, "MOL");
738                         }
739                     }
740                     else
741                     {
742                         /* Residue number is the atom number */
743                         atoms->atom[i].resind = i;
744                         /* Generate an residue name from the particle type */
745                         if (atoms->atom[i].type < 26)
746                         {
747                             sprintf(buf, "T%c", 'A'+atoms->atom[i].type);
748                         }
749                         else
750                         {
751                             sprintf(buf, "T%c%c",
752                                     'A'+atoms->atom[i].type/26, 'A'+atoms->atom[i].type%26);
753                         }
754                         t_atoms_set_resinfo(atoms, i, symtab, buf, i, ' ', 0, ' ');
755                     }
756
757                     if (r == 3)
758                     {
759                         level--;
760                     }
761                     i++;
762                 }
763             }
764             atoms->nres = atoms->nr;
765
766             if (i != atoms->nr)
767             {
768                 gmx_fatal(FARGS, "Internal inconsistency in Espresso routines, read %d atoms, expected %d atoms", i, atoms->nr);
769             }
770         }
771         else if (level == 1 && strcmp(word, "variable") == 0 && !bFoundVariable)
772         {
773             bFoundVariable = TRUE;
774             level         += check_open_parenthesis(fp, r, infile, "variable");
775             while (level == 2 && (r = get_espresso_word(fp, word)))
776             {
777                 if (level == 2 && strcmp(word, "box_l") == 0)
778                 {
779                     for (m = 0; m < 3; m++)
780                     {
781                         r = get_espresso_word(fp, word);
782                         sscanf(word, "%lf", &d);
783                         box[m][m] = d;
784                     }
785                     level += check_close_parenthesis(fp, r, infile, "box_l");
786                 }
787             }
788         }
789         else if (r == 2)
790         {
791             level++;
792         }
793         else if (r == 3)
794         {
795             level--;
796         }
797     }
798
799     if (!bFoundParticles)
800     {
801         fprintf(stderr, "Did not find a particles section in Espresso file '%s'\n",
802                 infile);
803     }
804
805     gmx_fio_fclose(fp);
806 }
807
808 static int get_espresso_coordnum(const char *infile)
809 {
810     FILE    *fp;
811     char     word[STRLEN];
812     int      natoms, level, r;
813     gmx_bool bFoundParticles;
814
815     natoms = 0;
816
817     fp = gmx_fio_fopen(infile, "r");
818
819     bFoundParticles = FALSE;
820     level           = 0;
821     while ((r = get_espresso_word(fp, word)) && !bFoundParticles)
822     {
823         if (level == 1 && strcmp(word, "particles") == 0 && !bFoundParticles)
824         {
825             bFoundParticles = TRUE;
826             level          += check_open_parenthesis(fp, r, infile, "particles");
827             while (level > 0 && (r = get_espresso_word(fp, word)))
828             {
829                 if (r == 2)
830                 {
831                     level++;
832                     if (level == 2)
833                     {
834                         natoms++;
835                     }
836                 }
837                 else if (r == 3)
838                 {
839                     level--;
840                 }
841             }
842         }
843         else if (r == 2)
844         {
845             level++;
846         }
847         else if (r == 3)
848         {
849             level--;
850         }
851     }
852     if (!bFoundParticles)
853     {
854         fprintf(stderr, "Did not find a particles section in Espresso file '%s'\n",
855                 infile);
856     }
857
858     gmx_fio_fclose(fp);
859
860     return natoms;
861 }
862
863 static void write_espresso_conf_indexed(FILE *out, const char *title,
864                                         t_atoms *atoms, int nx, atom_id *index,
865                                         rvec *x, rvec *v, matrix box)
866 {
867     int i, j;
868
869     fprintf(out, "# %s\n", title);
870     if (TRICLINIC(box))
871     {
872         gmx_warning("The Espresso format does not support triclinic unit-cells");
873     }
874     fprintf(out, "{variable {box_l %f %f %f}}\n", box[0][0], box[1][1], box[2][2]);
875
876     fprintf(out, "{particles {id pos type q%s}\n", v ? " v" : "");
877     for (i = 0; i < nx; i++)
878     {
879         if (index)
880         {
881             j = index[i];
882         }
883         else
884         {
885             j = i;
886         }
887         fprintf(out, "\t{%d %f %f %f %d %g",
888                 j, x[j][XX], x[j][YY], x[j][ZZ],
889                 atoms->atom[j].type, atoms->atom[j].q);
890         if (v)
891         {
892             fprintf(out, " %f %f %f", v[j][XX], v[j][YY], v[j][ZZ]);
893         }
894         fprintf(out, "}\n");
895     }
896     fprintf(out, "}\n");
897 }
898
899 static void get_coordnum_fp (FILE *in, char *title, int *natoms)
900 {
901     char line[STRLEN+1];
902
903     fgets2 (title, STRLEN, in);
904     fgets2 (line, STRLEN, in);
905     if (sscanf (line, "%d", natoms) != 1)
906     {
907         gmx_fatal(FARGS, "gro file does not have the number of atoms on the second line");
908     }
909 }
910
911 static void get_coordnum (const char *infile, int *natoms)
912 {
913     FILE *in;
914     char  title[STRLEN];
915
916     in = gmx_fio_fopen(infile, "r");
917     get_coordnum_fp(in, title, natoms);
918     gmx_fio_fclose (in);
919 }
920
921 static gmx_bool get_w_conf(FILE *in, const char *infile, char *title,
922                            t_symtab *symtab, t_atoms *atoms, int *ndec,
923                            rvec x[], rvec *v, matrix box)
924 {
925     char       name[6];
926     char       line[STRLEN+1], *ptr;
927     char       buf[256];
928     double     x1, y1, z1, x2, y2, z2;
929     rvec       xmin, xmax;
930     int        natoms, i, m, resnr, newres, oldres, ddist, c;
931     gmx_bool   bFirst, bVel;
932     char      *p1, *p2, *p3;
933
934     newres  = -1;
935     oldres  = NOTSET; /* Unlikely number for the first residue! */
936     ddist   = 0;
937
938     /* Read the title and number of atoms */
939     get_coordnum_fp(in, title, &natoms);
940
941     if (natoms > atoms->nr)
942     {
943         gmx_fatal(FARGS, "gro file contains more atoms (%d) than expected (%d)",
944                   natoms, atoms->nr);
945     }
946     else if (natoms <  atoms->nr)
947     {
948         fprintf(stderr, "Warning: gro file contains less atoms (%d) than expected"
949                 " (%d)\n", natoms, atoms->nr);
950     }
951
952     bFirst = TRUE;
953
954     bVel = FALSE;
955
956     /* just pray the arrays are big enough */
957     for (i = 0; (i < natoms); i++)
958     {
959         if ((fgets2 (line, STRLEN, in)) == NULL)
960         {
961             unexpected_eof(infile, i+2);
962         }
963         if (strlen(line) < 39)
964         {
965             gmx_fatal(FARGS, "Invalid line in %s for atom %d:\n%s", infile, i+1, line);
966         }
967
968         /* determine read precision from distance between periods
969            (decimal points) */
970         if (bFirst)
971         {
972             bFirst = FALSE;
973             p1     = strchr(line, '.');
974             if (p1 == NULL)
975             {
976                 gmx_fatal(FARGS, "A coordinate in file %s does not contain a '.'", infile);
977             }
978             p2 = strchr(&p1[1], '.');
979             if (p2 == NULL)
980             {
981                 gmx_fatal(FARGS, "A coordinate in file %s does not contain a '.'", infile);
982             }
983             ddist = p2 - p1;
984             *ndec = ddist - 5;
985
986             p3 = strchr(&p2[1], '.');
987             if (p3 == NULL)
988             {
989                 gmx_fatal(FARGS, "A coordinate in file %s does not contain a '.'", infile);
990             }
991
992             if (p3 - p2 != ddist)
993             {
994                 gmx_fatal(FARGS, "The spacing of the decimal points in file %s is not consistent for x, y and z", infile);
995             }
996         }
997
998         /* residue number*/
999         memcpy(name, line, 5);
1000         name[5] = '\0';
1001         sscanf(name, "%d", &resnr);
1002         memcpy(name, line+5, 5);
1003         name[5] = '\0';
1004         if (resnr != oldres)
1005         {
1006             oldres = resnr;
1007             newres++;
1008             if (newres >= natoms)
1009             {
1010                 gmx_fatal(FARGS, "More residues than atoms in %s (natoms = %d)",
1011                           infile, natoms);
1012             }
1013             atoms->atom[i].resind = newres;
1014             t_atoms_set_resinfo(atoms, i, symtab, name, resnr, ' ', 0, ' ');
1015         }
1016         else
1017         {
1018             atoms->atom[i].resind = newres;
1019         }
1020
1021         /* atomname */
1022         memcpy(name, line+10, 5);
1023         atoms->atomname[i] = put_symtab(symtab, name);
1024
1025         /* eventueel controle atomnumber met i+1 */
1026
1027         /* coordinates (start after residue shit) */
1028         ptr = line + 20;
1029         /* Read fixed format */
1030         for (m = 0; m < DIM; m++)
1031         {
1032             for (c = 0; (c < ddist && ptr[0]); c++)
1033             {
1034                 buf[c] = ptr[0];
1035                 ptr++;
1036             }
1037             buf[c] = '\0';
1038             if (sscanf (buf, "%lf %lf", &x1, &x2) != 1)
1039             {
1040                 gmx_fatal(FARGS, "Something is wrong in the coordinate formatting of file %s. Note that gro is fixed format (see the manual)", infile);
1041             }
1042             else
1043             {
1044                 x[i][m] = x1;
1045             }
1046         }
1047
1048         /* velocities (start after residues and coordinates) */
1049         if (v)
1050         {
1051             /* Read fixed format */
1052             for (m = 0; m < DIM; m++)
1053             {
1054                 for (c = 0; (c < ddist && ptr[0]); c++)
1055                 {
1056                     buf[c] = ptr[0];
1057                     ptr++;
1058                 }
1059                 buf[c] = '\0';
1060                 if (sscanf (buf, "%lf", &x1) != 1)
1061                 {
1062                     v[i][m] = 0;
1063                 }
1064                 else
1065                 {
1066                     v[i][m] = x1;
1067                     bVel    = TRUE;
1068                 }
1069             }
1070         }
1071     }
1072     atoms->nres = newres + 1;
1073
1074     /* box */
1075     fgets2 (line, STRLEN, in);
1076     if (sscanf (line, "%lf%lf%lf", &x1, &y1, &z1) != 3)
1077     {
1078         gmx_warning("Bad box in file %s", infile);
1079
1080         /* Generate a cubic box */
1081         for (m = 0; (m < DIM); m++)
1082         {
1083             xmin[m] = xmax[m] = x[0][m];
1084         }
1085         for (i = 1; (i < atoms->nr); i++)
1086         {
1087             for (m = 0; (m < DIM); m++)
1088             {
1089                 xmin[m] = min(xmin[m], x[i][m]);
1090                 xmax[m] = max(xmax[m], x[i][m]);
1091             }
1092         }
1093         for (i = 0; i < DIM; i++)
1094         {
1095             for (m = 0; m < DIM; m++)
1096             {
1097                 box[i][m] = 0.0;
1098             }
1099         }
1100         for (m = 0; (m < DIM); m++)
1101         {
1102             box[m][m] = (xmax[m]-xmin[m]);
1103         }
1104         fprintf(stderr, "Generated a cubic box %8.3f x %8.3f x %8.3f\n",
1105                 box[XX][XX], box[YY][YY], box[ZZ][ZZ]);
1106     }
1107     else
1108     {
1109         /* We found the first three values, the diagonal elements */
1110         box[XX][XX] = x1;
1111         box[YY][YY] = y1;
1112         box[ZZ][ZZ] = z1;
1113         if (sscanf (line, "%*f%*f%*f%lf%lf%lf%lf%lf%lf",
1114                     &x1, &y1, &z1, &x2, &y2, &z2) != 6)
1115         {
1116             x1 = y1 = z1 = x2 = y2 = z2 = 0.0;
1117         }
1118         box[XX][YY] = x1;
1119         box[XX][ZZ] = y1;
1120         box[YY][XX] = z1;
1121         box[YY][ZZ] = x2;
1122         box[ZZ][XX] = y2;
1123         box[ZZ][YY] = z2;
1124     }
1125
1126     return bVel;
1127 }
1128
1129 static void read_whole_conf(const char *infile, char *title,
1130                             t_atoms *atoms, rvec x[], rvec *v, matrix box)
1131 {
1132     FILE    *in;
1133     int      ndec;
1134     t_symtab symtab;
1135
1136     /* open file */
1137     in = gmx_fio_fopen(infile, "r");
1138
1139     open_symtab(&symtab);
1140     get_w_conf(in, infile, title, &symtab, atoms, &ndec, x, v, box);
1141     /* We can't free the symbols, as they are still used in atoms, so
1142      * the only choice is to leak them. */
1143     free_symtab(&symtab);
1144
1145     gmx_fio_fclose(in);
1146 }
1147
1148 gmx_bool gro_next_x_or_v(FILE *status, t_trxframe *fr)
1149 {
1150     t_atoms  atoms;
1151     t_symtab symtab;
1152     char     title[STRLEN], *p;
1153     double   tt;
1154     int      ndec = 0, i;
1155
1156     if (gmx_eof(status))
1157     {
1158         return FALSE;
1159     }
1160
1161     open_symtab(&symtab);
1162     atoms.nr = fr->natoms;
1163     snew(atoms.atom, fr->natoms);
1164     atoms.nres = fr->natoms;
1165     snew(atoms.resinfo, fr->natoms);
1166     snew(atoms.atomname, fr->natoms);
1167
1168     fr->bV    = get_w_conf(status, title, title, &symtab, &atoms, &ndec, fr->x, fr->v, fr->box);
1169     fr->bPrec = TRUE;
1170     fr->prec  = 1;
1171     /* prec = 10^ndec: */
1172     for (i = 0; i < ndec; i++)
1173     {
1174         fr->prec *= 10;
1175     }
1176     fr->title  = title;
1177     fr->bTitle = TRUE;
1178     fr->bX     = TRUE;
1179     fr->bBox   = TRUE;
1180
1181     sfree(atoms.atom);
1182     sfree(atoms.resinfo);
1183     sfree(atoms.atomname);
1184     done_symtab(&symtab);
1185
1186     if ((p = strstr(title, "t=")) != NULL)
1187     {
1188         p += 2;
1189         if (sscanf(p, "%lf", &tt) == 1)
1190         {
1191             fr->time  = tt;
1192             fr->bTime = TRUE;
1193         }
1194         else
1195         {
1196             fr->time  = 0;
1197             fr->bTime = FALSE;
1198         }
1199     }
1200
1201     if (atoms.nr != fr->natoms)
1202     {
1203         gmx_fatal(FARGS, "Number of atoms in gro frame (%d) doesn't match the number in the previous frame (%d)", atoms.nr, fr->natoms);
1204     }
1205
1206     return TRUE;
1207 }
1208
1209 int gro_first_x_or_v(FILE *status, t_trxframe *fr)
1210 {
1211     char title[STRLEN];
1212
1213     frewind(status);
1214     fprintf(stderr, "Reading frames from gro file");
1215     get_coordnum_fp(status, title, &fr->natoms);
1216     frewind(status);
1217     fprintf(stderr, " '%s', %d atoms.\n", title, fr->natoms);
1218     fr->bTitle = TRUE;
1219     fr->title  = title;
1220     if (fr->natoms == 0)
1221     {
1222         gmx_file("No coordinates in gro file");
1223     }
1224
1225     snew(fr->x, fr->natoms);
1226     snew(fr->v, fr->natoms);
1227     gro_next_x_or_v(status, fr);
1228
1229     return fr->natoms;
1230 }
1231
1232 static void make_hconf_format(int pr, gmx_bool bVel, char format[])
1233 {
1234     int l, vpr;
1235
1236     /* build format string for printing,
1237        something like "%8.3f" for x and "%8.4f" for v */
1238     if (pr < 0)
1239     {
1240         pr = 0;
1241     }
1242     if (pr > 30)
1243     {
1244         pr = 30;
1245     }
1246     l   = pr+5;
1247     vpr = pr+1;
1248     if (bVel)
1249     {
1250         sprintf(format, "%%%d.%df%%%d.%df%%%d.%df%%%d.%df%%%d.%df%%%d.%df\n",
1251                 l, pr, l, pr, l, pr, l, vpr, l, vpr, l, vpr);
1252     }
1253     else
1254     {
1255         sprintf(format, "%%%d.%df%%%d.%df%%%d.%df\n", l, pr, l, pr, l, pr);
1256     }
1257
1258 }
1259
1260 static void write_hconf_box(FILE *out, int pr, matrix box)
1261 {
1262     char format[100];
1263     int  l;
1264
1265     if (pr < 5)
1266     {
1267         pr = 5;
1268     }
1269     l = pr+5;
1270
1271     if (box[XX][YY] || box[XX][ZZ] || box[YY][XX] || box[YY][ZZ] ||
1272         box[ZZ][XX] || box[ZZ][YY])
1273     {
1274         sprintf(format, "%%%d.%df%%%d.%df%%%d.%df"
1275                 "%%%d.%df%%%d.%df%%%d.%df%%%d.%df%%%d.%df%%%d.%df\n",
1276                 l, pr, l, pr, l, pr, l, pr, l, pr, l, pr, l, pr, l, pr, l, pr);
1277         fprintf(out, format,
1278                 box[XX][XX], box[YY][YY], box[ZZ][ZZ],
1279                 box[XX][YY], box[XX][ZZ], box[YY][XX],
1280                 box[YY][ZZ], box[ZZ][XX], box[ZZ][YY]);
1281     }
1282     else
1283     {
1284         sprintf(format, "%%%d.%df%%%d.%df%%%d.%df\n", l, pr, l, pr, l, pr);
1285         fprintf(out, format,
1286                 box[XX][XX], box[YY][YY], box[ZZ][ZZ]);
1287     }
1288 }
1289
1290 void write_hconf_indexed_p(FILE *out, const char *title, t_atoms *atoms,
1291                            int nx, const atom_id index[], int pr,
1292                            rvec *x, rvec *v, matrix box)
1293 {
1294     char resnm[6], nm[6], format[100];
1295     int  ai, i, resind, resnr;
1296
1297     bromacs(format, 99);
1298     fprintf (out, "%s\n", (title && title[0]) ? title : format);
1299     fprintf (out, "%5d\n", nx);
1300
1301     make_hconf_format(pr, v != NULL, format);
1302
1303     for (i = 0; (i < nx); i++)
1304     {
1305         ai = index[i];
1306
1307         resind = atoms->atom[ai].resind;
1308         strncpy(resnm, " ??? ", sizeof(resnm)-1);
1309         if (resind < atoms->nres)
1310         {
1311             strncpy(resnm, *atoms->resinfo[resind].name, sizeof(resnm)-1);
1312             resnr = atoms->resinfo[resind].nr;
1313         }
1314         else
1315         {
1316             strncpy(resnm, " ??? ", sizeof(resnm)-1);
1317             resnr = resind + 1;
1318         }
1319
1320         if (atoms->atom)
1321         {
1322             strncpy(nm, *atoms->atomname[ai], sizeof(nm)-1);
1323         }
1324         else
1325         {
1326             strncpy(nm, " ??? ", sizeof(nm)-1);
1327         }
1328
1329         fprintf(out, "%5d%-5.5s%5.5s%5d", resnr%100000, resnm, nm, (ai+1)%100000);
1330         /* next fprintf uses built format string */
1331         if (v)
1332         {
1333             fprintf(out, format,
1334                     x[ai][XX], x[ai][YY], x[ai][ZZ], v[ai][XX], v[ai][YY], v[ai][ZZ]);
1335         }
1336         else
1337         {
1338             fprintf(out, format,
1339                     x[ai][XX], x[ai][YY], x[ai][ZZ]);
1340         }
1341     }
1342
1343     write_hconf_box(out, pr, box);
1344
1345     fflush(out);
1346 }
1347
1348 static void write_hconf_mtop(FILE *out, const char *title, gmx_mtop_t *mtop,
1349                              int pr,
1350                              rvec *x, rvec *v, matrix box)
1351 {
1352     char                    format[100];
1353     int                     i, resnr;
1354     gmx_mtop_atomloop_all_t aloop;
1355     t_atom                 *atom;
1356     char                   *atomname, *resname;
1357
1358     bromacs(format, 99);
1359     fprintf (out, "%s\n", (title && title[0]) ? title : format);
1360     fprintf (out, "%5d\n", mtop->natoms);
1361
1362     make_hconf_format(pr, v != NULL, format);
1363
1364     aloop = gmx_mtop_atomloop_all_init(mtop);
1365     while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
1366     {
1367         gmx_mtop_atomloop_all_names(aloop, &atomname, &resnr, &resname);
1368
1369         fprintf(out, "%5d%-5.5s%5.5s%5d",
1370                 resnr%100000, resname, atomname, (i+1)%100000);
1371         /* next fprintf uses built format string */
1372         if (v)
1373         {
1374             fprintf(out, format,
1375                     x[i][XX], x[i][YY], x[i][ZZ], v[i][XX], v[i][YY], v[i][ZZ]);
1376         }
1377         else
1378         {
1379             fprintf(out, format,
1380                     x[i][XX], x[i][YY], x[i][ZZ]);
1381         }
1382     }
1383
1384     write_hconf_box(out, pr, box);
1385
1386     fflush(out);
1387 }
1388
1389 void write_hconf_p(FILE *out, const char *title, t_atoms *atoms, int pr,
1390                    rvec *x, rvec *v, matrix box)
1391 {
1392     atom_id *aa;
1393     int      i;
1394
1395     snew(aa, atoms->nr);
1396     for (i = 0; (i < atoms->nr); i++)
1397     {
1398         aa[i] = i;
1399     }
1400     write_hconf_indexed_p(out, title, atoms, atoms->nr, aa, pr, x, v, box);
1401     sfree(aa);
1402 }
1403
1404 void write_conf_p(const char *outfile, const char *title,
1405                   t_atoms *atoms, int pr,
1406                   rvec *x, rvec *v, matrix box)
1407 {
1408     FILE *out;
1409
1410     out = gmx_fio_fopen(outfile, "w");
1411     write_hconf_p(out, title, atoms, pr, x, v, box);
1412
1413     gmx_fio_fclose (out);
1414 }
1415
1416 static void write_conf(const char *outfile, const char *title, t_atoms *atoms,
1417                        rvec *x, rvec *v, matrix box)
1418 {
1419     write_conf_p(outfile, title, atoms, 3, x, v, box);
1420 }
1421
1422 void write_sto_conf_indexed(const char *outfile, const char *title,
1423                             t_atoms *atoms,
1424                             rvec x[], rvec *v, int ePBC, matrix box,
1425                             atom_id nindex, atom_id index[])
1426 {
1427     FILE       *out;
1428     int         ftp;
1429     t_trxframe  fr;
1430
1431     ftp = fn2ftp(outfile);
1432     switch (ftp)
1433     {
1434         case efGRO:
1435             out = gmx_fio_fopen(outfile, "w");
1436             write_hconf_indexed_p(out, title, atoms, nindex, index, 3, x, v, box);
1437             gmx_fio_fclose(out);
1438             break;
1439         case efG96:
1440             clear_trxframe(&fr, TRUE);
1441             fr.bTitle = TRUE;
1442             fr.title  = title;
1443             fr.natoms = atoms->nr;
1444             fr.bAtoms = TRUE;
1445             fr.atoms  = atoms;
1446             fr.bX     = TRUE;
1447             fr.x      = x;
1448             if (v)
1449             {
1450                 fr.bV = TRUE;
1451                 fr.v  = v;
1452             }
1453             fr.bBox = TRUE;
1454             copy_mat(box, fr.box);
1455             out = gmx_fio_fopen(outfile, "w");
1456             write_g96_conf(out, &fr, nindex, index);
1457             gmx_fio_fclose(out);
1458             break;
1459         case efPDB:
1460         case efBRK:
1461         case efENT:
1462         case efPQR:
1463             out = gmx_fio_fopen(outfile, "w");
1464             write_pdbfile_indexed(out, title, atoms, x, ePBC, box, ' ', -1, nindex, index, NULL, TRUE);
1465             gmx_fio_fclose(out);
1466             break;
1467         case efESP:
1468             out = gmx_fio_fopen(outfile, "w");
1469             write_espresso_conf_indexed(out, title, atoms, nindex, index, x, v, box);
1470             gmx_fio_fclose(out);
1471             break;
1472         case efTPR:
1473         case efTPB:
1474         case efTPA:
1475             gmx_fatal(FARGS, "Sorry, can not write a topology to %s", outfile);
1476             break;
1477         default:
1478             gmx_incons("Not supported in write_sto_conf_indexed");
1479     }
1480 }
1481
1482 static void write_xyz_conf(const char *outfile, const char *title,
1483                            t_atoms *atoms, rvec *x)
1484 {
1485     FILE          *fp;
1486     int            i, anr;
1487     real           value;
1488     char          *ptr, *name;
1489     gmx_atomprop_t aps = gmx_atomprop_init();
1490
1491     fp = gmx_fio_fopen(outfile, "w");
1492     fprintf(fp, "%3d\n", atoms->nr);
1493     fprintf(fp, "%s\n", title);
1494     for (i = 0; (i < atoms->nr); i++)
1495     {
1496         anr  = atoms->atom[i].atomnumber;
1497         name = *atoms->atomname[i];
1498         if (anr == NOTSET)
1499         {
1500             if (gmx_atomprop_query(aps, epropElement, "???", name, &value))
1501             {
1502                 anr = gmx_nint(value);
1503             }
1504         }
1505         if ((ptr = gmx_atomprop_element(aps, anr)) == NULL)
1506         {
1507             ptr = name;
1508         }
1509         fprintf(fp, "%3s%10.5f%10.5f%10.5f\n", ptr,
1510                 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ]);
1511     }
1512     gmx_fio_fclose(fp);
1513     gmx_atomprop_destroy(aps);
1514 }
1515
1516 void write_sto_conf(const char *outfile, const char *title, t_atoms *atoms,
1517                     rvec x[], rvec *v, int ePBC, matrix box)
1518 {
1519     FILE       *out;
1520     int         ftp;
1521     t_trxframe  fr;
1522
1523     ftp = fn2ftp(outfile);
1524     switch (ftp)
1525     {
1526         case efGRO:
1527             write_conf(outfile, title, atoms, x, v, box);
1528             break;
1529         case efG96:
1530             clear_trxframe(&fr, TRUE);
1531             fr.bTitle = TRUE;
1532             fr.title  = title;
1533             fr.natoms = atoms->nr;
1534             fr.bAtoms = TRUE;
1535             fr.atoms  = atoms;
1536             fr.bX     = TRUE;
1537             fr.x      = x;
1538             if (v)
1539             {
1540                 fr.bV = TRUE;
1541                 fr.v  = v;
1542             }
1543             fr.bBox = TRUE;
1544             copy_mat(box, fr.box);
1545             out = gmx_fio_fopen(outfile, "w");
1546             write_g96_conf(out, &fr, -1, NULL);
1547             gmx_fio_fclose(out);
1548             break;
1549         case efXYZ:
1550             write_xyz_conf(outfile, (strlen(title) > 0) ? title : outfile, atoms, x);
1551             break;
1552         case efPDB:
1553         case efBRK:
1554         case efENT:
1555             out = gmx_fio_fopen(outfile, "w");
1556             write_pdbfile(out, title, atoms, x, ePBC, box, ' ', -1, NULL, TRUE);
1557             gmx_fio_fclose(out);
1558             break;
1559         case efESP:
1560             out = gmx_fio_fopen(outfile, "w");
1561             write_espresso_conf_indexed(out, title, atoms, atoms->nr, NULL, x, v, box);
1562             gmx_fio_fclose(out);
1563             break;
1564         case efTPR:
1565         case efTPB:
1566         case efTPA:
1567             gmx_fatal(FARGS, "Sorry, can not write a topology to %s", outfile);
1568             break;
1569         default:
1570             gmx_incons("Not supported in write_sto_conf");
1571     }
1572 }
1573
1574 void write_sto_conf_mtop(const char *outfile, const char *title,
1575                          gmx_mtop_t *mtop,
1576                          rvec x[], rvec *v, int ePBC, matrix box)
1577 {
1578     int     ftp;
1579     FILE   *out;
1580     t_atoms atoms;
1581
1582     ftp = fn2ftp(outfile);
1583     switch (ftp)
1584     {
1585         case efGRO:
1586             out = gmx_fio_fopen(outfile, "w");
1587             write_hconf_mtop(out, title, mtop, 3, x, v, box);
1588             gmx_fio_fclose(out);
1589             break;
1590         default:
1591             /* This is a brute force approach which requires a lot of memory.
1592              * We should implement mtop versions of all writing routines.
1593              */
1594             atoms = gmx_mtop_global_atoms(mtop);
1595
1596             write_sto_conf(outfile, title, &atoms, x, v, ePBC, box);
1597
1598             done_atom(&atoms);
1599             break;
1600     }
1601 }
1602
1603 static int get_xyz_coordnum(const char *infile)
1604 {
1605     FILE *fp;
1606     int   n;
1607
1608     fp = gmx_fio_fopen(infile, "r");
1609     if (fscanf(fp, "%d", &n) != 1)
1610     {
1611         gmx_fatal(FARGS, "Can not read number of atoms from %s", infile);
1612     }
1613     gmx_fio_fclose(fp);
1614
1615     return n;
1616 }
1617
1618 static void read_xyz_conf(const char *infile, char *title,
1619                           t_atoms *atoms, rvec *x)
1620 {
1621     FILE     *fp;
1622     int       i, n;
1623     double    xx, yy, zz;
1624     t_symtab *tab;
1625     char      atomnm[32], buf[STRLEN];
1626
1627     snew(tab, 1);
1628     fp = gmx_fio_fopen(infile, "r");
1629     fgets2(buf, STRLEN-1, fp);
1630     if (sscanf(buf, "%d", &n) != 1)
1631     {
1632         gmx_fatal(FARGS, "Can not read number of atoms from %s", infile);
1633     }
1634     fgets2(buf, STRLEN-1, fp);
1635     strcpy(title, buf);
1636     for (i = 0; (i < n); i++)
1637     {
1638         fgets2(buf, STRLEN-1, fp);
1639         if (sscanf(buf, "%s%lf%lf%lf", atomnm, &xx, &yy, &zz) != 4)
1640         {
1641             gmx_fatal(FARGS, "Can not read coordinates from %s", infile);
1642         }
1643         atoms->atomname[i] = put_symtab(tab, atomnm);
1644         x[i][XX]           = xx*0.1;
1645         x[i][YY]           = yy*0.1;
1646         x[i][ZZ]           = zz*0.1;
1647     }
1648     gmx_fio_fclose(fp);
1649 }
1650
1651 void get_stx_coordnum(const char *infile, int *natoms)
1652 {
1653     FILE      *in;
1654     int        ftp, tpxver, tpxgen;
1655     t_trxframe fr;
1656     char       g96_line[STRLEN+1];
1657
1658     ftp = fn2ftp(infile);
1659     range_check(ftp, 0, efNR);
1660     switch (ftp)
1661     {
1662         case efGRO:
1663             get_coordnum(infile, natoms);
1664             break;
1665         case efG96:
1666             in        = gmx_fio_fopen(infile, "r");
1667             fr.title  = NULL;
1668             fr.natoms = -1;
1669             fr.atoms  = NULL;
1670             fr.x      = NULL;
1671             fr.v      = NULL;
1672             fr.f      = NULL;
1673             *natoms   = read_g96_conf(in, infile, &fr, g96_line);
1674             gmx_fio_fclose(in);
1675             break;
1676         case efXYZ:
1677             *natoms = get_xyz_coordnum(infile);
1678             break;
1679         case efPDB:
1680         case efBRK:
1681         case efENT:
1682             in = gmx_fio_fopen(infile, "r");
1683             get_pdb_coordnum(in, natoms);
1684             gmx_fio_fclose(in);
1685             break;
1686         case efESP:
1687             *natoms = get_espresso_coordnum(infile);
1688             break;
1689         case efTPA:
1690         case efTPB:
1691         case efTPR:
1692         {
1693             t_tpxheader tpx;
1694
1695             read_tpxheader(infile, &tpx, TRUE, &tpxver, &tpxgen);
1696             *natoms = tpx.natoms;
1697             break;
1698         }
1699         default:
1700             gmx_fatal(FARGS, "File type %s not supported in get_stx_coordnum",
1701                       ftp2ext(ftp));
1702     }
1703 }
1704
1705 void read_stx_conf(const char *infile, char *title, t_atoms *atoms,
1706                    rvec x[], rvec *v, int *ePBC, matrix box)
1707 {
1708     FILE       *in;
1709     char        buf[256];
1710     gmx_mtop_t *mtop;
1711     t_topology  top;
1712     t_trxframe  fr;
1713     int         i, ftp, natoms;
1714     real        d;
1715     char        g96_line[STRLEN+1];
1716
1717     if (atoms->nr == 0)
1718     {
1719         fprintf(stderr, "Warning: Number of atoms in %s is 0\n", infile);
1720     }
1721     else if (atoms->atom == NULL)
1722     {
1723         gmx_mem("Uninitialized array atom");
1724     }
1725
1726     if (ePBC)
1727     {
1728         *ePBC = -1;
1729     }
1730
1731     ftp = fn2ftp(infile);
1732     switch (ftp)
1733     {
1734         case efGRO:
1735             read_whole_conf(infile, title, atoms, x, v, box);
1736             break;
1737         case efXYZ:
1738             read_xyz_conf(infile, title, atoms, x);
1739             break;
1740         case efG96:
1741             fr.title  = NULL;
1742             fr.natoms = atoms->nr;
1743             fr.atoms  = atoms;
1744             fr.x      = x;
1745             fr.v      = v;
1746             fr.f      = NULL;
1747             in        = gmx_fio_fopen(infile, "r");
1748             read_g96_conf(in, infile, &fr, g96_line);
1749             gmx_fio_fclose(in);
1750             copy_mat(fr.box, box);
1751             strncpy(title, fr.title, STRLEN);
1752             break;
1753         case efPDB:
1754         case efBRK:
1755         case efENT:
1756             read_pdb_conf(infile, title, atoms, x, ePBC, box, TRUE, NULL);
1757             break;
1758         case efESP:
1759             read_espresso_conf(infile, atoms, x, v, box);
1760             break;
1761         case efTPR:
1762         case efTPB:
1763         case efTPA:
1764             snew(mtop, 1);
1765             i = read_tpx(infile, NULL, box, &natoms, x, v, NULL, mtop);
1766             if (ePBC)
1767             {
1768                 *ePBC = i;
1769             }
1770
1771             strcpy(title, *(mtop->name));
1772
1773             /* Free possibly allocated memory */
1774             done_atom(atoms);
1775
1776             *atoms = gmx_mtop_global_atoms(mtop);
1777             top    = gmx_mtop_t_to_t_topology(mtop);
1778             tpx_make_chain_identifiers(atoms, &top.mols);
1779
1780             sfree(mtop);
1781             /* The strings in the symtab are still in use in the returned t_atoms
1782              * structure, so we should not free them. But there is no place to put the
1783              * symbols; the only choice is to leak the memory...
1784              * So we clear the symbol table before freeing the topology structure. */
1785             free_symtab(&top.symtab);
1786             done_top(&top);
1787
1788             break;
1789         default:
1790             gmx_incons("Not supported in read_stx_conf");
1791     }
1792 }