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