2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
45 #include "gromacs/fileio/filenm.h"
46 #include "gromacs/fileio/gmxfio.h"
47 #include "gromacs/fileio/pdbio.h"
48 #include "gromacs/fileio/tpxio.h"
49 #include "gromacs/fileio/trx.h"
50 #include "gromacs/fileio/trxio.h"
51 #include "gromacs/fileio/xdrf.h"
52 #include "gromacs/legacyheaders/copyrite.h"
53 #include "gromacs/legacyheaders/macros.h"
54 #include "gromacs/legacyheaders/typedefs.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/pbcutil/pbc.h"
57 #include "gromacs/topology/mtop_util.h"
58 #include "gromacs/topology/symtab.h"
59 #include "gromacs/topology/topology.h"
60 #include "gromacs/utility/cstringutil.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/futil.h"
63 #include "gromacs/utility/smalloc.h"
67 static int read_g96_pos(char line[], t_symtab *symtab,
68 FILE *fp, const char *infile,
73 int nwanted, natoms, atnr, resnr = 0, oldres, newres, shift;
74 char anm[STRLEN], resnm[STRLEN];
95 oldres = -666; /* Unlikely number for the first residue! */
97 while (!bEnd && fgets2(line, STRLEN, fp))
99 bEnd = (strncmp(line, "END", 3) == 0);
100 if (!bEnd && (line[0] != '#'))
102 if (sscanf(line+shift, "%15lf%15lf%15lf", &db1, &db2, &db3) != 3)
104 gmx_fatal(FARGS, "Did not find 3 coordinates for atom %d in %s\n",
107 if ((nwanted != -1) && (natoms >= nwanted))
110 "Found more coordinates (%d) in %s than expected %d\n",
111 natoms, infile, nwanted);
116 (sscanf(line, "%5d%c%5s%c%5s%7d", &resnr, &c1, resnm, &c2, anm, &atnr)
126 strncpy(resnm, "???", sizeof(resnm)-1);
128 strncpy(anm, "???", sizeof(anm)-1);
130 atoms->atomname[natoms] = put_symtab(symtab, anm);
135 if (newres >= atoms->nr)
137 gmx_fatal(FARGS, "More residues than atoms in %s (natoms = %d)",
140 atoms->atom[natoms].resind = newres;
141 if (newres+1 > atoms->nres)
143 atoms->nres = newres+1;
145 t_atoms_set_resinfo(atoms, natoms, symtab, resnm, resnr, ' ', 0, ' ');
149 atoms->atom[natoms].resind = newres;
154 fr->x[natoms][0] = db1;
155 fr->x[natoms][1] = db2;
156 fr->x[natoms][2] = db3;
161 if ((nwanted != -1) && natoms != nwanted)
164 "Warning: found less coordinates (%d) in %s than expected %d\n",
165 natoms, infile, nwanted);
174 static int read_g96_vel(char line[], FILE *fp, const char *infile,
178 int nwanted, natoms = -1, shift;
179 double db1, db2, db3;
181 nwanted = fr->natoms;
185 if (strcmp(line, "VELOCITYRED") == 0)
195 while (!bEnd && fgets2(line, STRLEN, fp))
197 bEnd = (strncmp(line, "END", 3) == 0);
198 if (!bEnd && (line[0] != '#'))
200 if (sscanf(line+shift, "%15lf%15lf%15lf", &db1, &db2, &db3) != 3)
202 gmx_fatal(FARGS, "Did not find 3 velocities for atom %d in %s",
205 if ((nwanted != -1) && (natoms >= nwanted))
207 gmx_fatal(FARGS, "Found more velocities (%d) in %s than expected %d\n",
208 natoms, infile, nwanted);
212 fr->v[natoms][0] = db1;
213 fr->v[natoms][1] = db2;
214 fr->v[natoms][2] = db3;
219 if ((nwanted != -1) && (natoms != nwanted))
222 "Warning: found less velocities (%d) in %s than expected %d\n",
223 natoms, infile, nwanted);
230 int read_g96_conf(FILE *fp, const char *infile, t_trxframe *fr, char *line)
232 t_symtab *symtab = NULL;
233 gmx_bool bAtStart, bTime, bAtoms, bPos, bVel, bBox, bEnd, bFinished;
235 double db1, db2, db3, db4, db5, db6, db7, db8, db9;
237 bAtStart = (ftell(fp) == 0);
239 clear_trxframe(fr, FALSE);
251 while (!fr->bTitle && fgets2(line, STRLEN, fp))
253 fr->bTitle = (strcmp(line, "TITLE") == 0);
255 if (fr->title == NULL)
257 fgets2(line, STRLEN, fp);
258 fr->title = gmx_strdup(line);
261 while (!bEnd && fgets2(line, STRLEN, fp))
263 bEnd = (strcmp(line, "END") == 0);
265 fgets2(line, STRLEN, fp);
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). */
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);
282 if (!fr->bTime && !fr->bX)
288 bFinished = (fgets2(line, STRLEN, fp) == NULL);
290 while (!bFinished && (line[0] == '#'));
291 sscanf(line, "%15d%15lf", &(fr->step), &db1);
305 natoms = read_g96_pos(line, symtab, fp, infile, fr);
315 natoms = read_g96_vel(line, fp, infile, fr);
322 while (!bEnd && fgets2(line, STRLEN, fp))
324 bEnd = (strncmp(line, "END", 3) == 0);
325 if (!bEnd && (line[0] != '#'))
327 nbp = sscanf(line, "%15lf%15lf%15lf%15lf%15lf%15lf%15lf%15lf%15lf",
328 &db1, &db2, &db3, &db4, &db5, &db6, &db7, &db8, &db9);
331 gmx_fatal(FARGS, "Found a BOX line, but no box in %s", infile);
333 fr->box[XX][XX] = db1;
334 fr->box[YY][YY] = db2;
335 fr->box[ZZ][ZZ] = db3;
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;
350 while (!bFinished && fgets2(line, STRLEN, fp));
359 void write_g96_conf(FILE *out, t_trxframe *fr,
360 int nindex, const atom_id *index)
378 fprintf(out, "TITLE\n%s\nEND\n", fr->title);
380 if (fr->bStep || fr->bTime)
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);
389 fprintf(out, "POSITION\n");
390 for (i = 0; i < nout; i++)
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]);
409 fprintf(out, "POSITIONRED\n");
410 for (i = 0; i < nout; i++)
420 fprintf(out, "%15.9f%15.9f%15.9f\n",
421 fr->x[a][XX], fr->x[a][YY], fr->x[a][ZZ]);
424 fprintf(out, "END\n");
430 fprintf(out, "VELOCITY\n");
431 for (i = 0; i < nout; i++)
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]);
450 fprintf(out, "VELOCITYRED\n");
451 for (i = 0; i < nout; i++)
461 fprintf(out, "%15.9f%15.9f%15.9f\n",
462 fr->v[a][XX], fr->v[a][YY], fr->v[a][ZZ]);
465 fprintf(out, "END\n");
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])
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]);
480 fprintf(out, "END\n");
484 static int get_espresso_word(FILE *fp, char word[])
496 if (i == ' ' || i == '\n' || i == '\t')
521 word[nc++] = (char)i;
525 while (i != EOF && ret == 0);
529 /* printf("'%s'\n",word); */
534 static int check_open_parenthesis(FILE *fp, int r,
535 const char *infile, const char *keyword)
547 r = get_espresso_word(fp, word);
554 gmx_fatal(FARGS, "Expected '{' after '%s' in file '%s'",
562 static int check_close_parenthesis(FILE *fp, int r,
563 const char *infile, const char *keyword)
575 r = get_espresso_word(fp, word);
582 gmx_fatal(FARGS, "Expected '}' after section '%s' in file '%s'",
591 espID, espPOS, espTYPE, espQ, espV, espF, espMOLECULE, espNR
593 const char *esp_prop[espNR] = {
594 "id", "pos", "type", "q", "v", "f",
598 static void read_espresso_conf(const char *infile,
599 t_atoms *atoms, rvec x[], rvec *v, matrix box)
601 t_symtab *symtab = NULL;
603 char word[STRLEN], buf[STRLEN];
604 int natoms, level, npar, r, nprop, p, i, m, molnr;
607 gmx_bool bFoundParticles, bFoundProp, bFoundVariable, bMol;
617 fp = gmx_fio_fopen(infile, "r");
619 bFoundParticles = FALSE;
620 bFoundVariable = FALSE;
623 while ((r = get_espresso_word(fp, word)))
625 if (level == 1 && strcmp(word, "particles") == 0 && !bFoundParticles)
627 bFoundParticles = TRUE;
628 level += check_open_parenthesis(fp, r, infile, "particles");
630 while (level == 2 && (r = get_espresso_word(fp, word)))
633 for (p = 0; p < espNR; p++)
635 if (strcmp(word, esp_prop[p]) == 0)
639 /* printf(" prop[%d] = %s\n",nprop-1,esp_prop[prop[nprop-1]]); */
642 if (!bFoundProp && word[0] != '}')
644 gmx_fatal(FARGS, "Can not read Espresso files with particle property '%s'", word);
646 if (bFoundProp && p == espMOLECULE)
657 while (level > 0 && (r = get_espresso_word(fp, word)))
669 for (p = 0; p < nprop; p++)
674 r = get_espresso_word(fp, word);
678 for (m = 0; m < 3; m++)
680 r = get_espresso_word(fp, word);
681 sscanf(word, "%lf", &d);
686 r = get_espresso_word(fp, word);
687 atoms->atom[i].type = strtol(word, NULL, 10);
690 r = get_espresso_word(fp, word);
691 sscanf(word, "%lf", &d);
692 atoms->atom[i].q = d;
695 for (m = 0; m < 3; m++)
697 r = get_espresso_word(fp, word);
698 sscanf(word, "%lf", &d);
703 for (m = 0; m < 3; m++)
705 r = get_espresso_word(fp, word);
710 r = get_espresso_word(fp, word);
711 molnr = strtol(word, NULL, 10);
713 atoms->resinfo[atoms->atom[i-1].resind].nr != molnr)
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? */
724 atoms->atom[i].resind = atoms->atom[i-1].resind;
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);
734 if (i == 0 || atoms->atom[i].resind != atoms->atom[i-1].resind)
736 atoms->resinfo[atoms->atom[i].resind].name =
737 put_symtab(symtab, "MOL");
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)
747 sprintf(buf, "T%c", 'A'+atoms->atom[i].type);
751 sprintf(buf, "T%c%c",
752 'A'+atoms->atom[i].type/26, 'A'+atoms->atom[i].type%26);
754 t_atoms_set_resinfo(atoms, i, symtab, buf, i, ' ', 0, ' ');
764 atoms->nres = atoms->nr;
768 gmx_fatal(FARGS, "Internal inconsistency in Espresso routines, read %d atoms, expected %d atoms", i, atoms->nr);
771 else if (level == 1 && strcmp(word, "variable") == 0 && !bFoundVariable)
773 bFoundVariable = TRUE;
774 level += check_open_parenthesis(fp, r, infile, "variable");
775 while (level == 2 && (r = get_espresso_word(fp, word)))
777 if (level == 2 && strcmp(word, "box_l") == 0)
779 for (m = 0; m < 3; m++)
781 r = get_espresso_word(fp, word);
782 sscanf(word, "%lf", &d);
785 level += check_close_parenthesis(fp, r, infile, "box_l");
799 if (!bFoundParticles)
801 fprintf(stderr, "Did not find a particles section in Espresso file '%s'\n",
808 static int get_espresso_coordnum(const char *infile)
812 int natoms, level, r;
813 gmx_bool bFoundParticles;
817 fp = gmx_fio_fopen(infile, "r");
819 bFoundParticles = FALSE;
821 while ((r = get_espresso_word(fp, word)) && !bFoundParticles)
823 if (level == 1 && strcmp(word, "particles") == 0 && !bFoundParticles)
825 bFoundParticles = TRUE;
826 level += check_open_parenthesis(fp, r, infile, "particles");
827 while (level > 0 && (r = get_espresso_word(fp, word)))
852 if (!bFoundParticles)
854 fprintf(stderr, "Did not find a particles section in Espresso file '%s'\n",
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)
869 fprintf(out, "# %s\n", title);
872 gmx_warning("The Espresso format does not support triclinic unit-cells");
874 fprintf(out, "{variable {box_l %f %f %f}}\n", box[0][0], box[1][1], box[2][2]);
876 fprintf(out, "{particles {id pos type q%s}\n", v ? " v" : "");
877 for (i = 0; i < nx; i++)
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);
892 fprintf(out, " %f %f %f", v[j][XX], v[j][YY], v[j][ZZ]);
899 static void get_coordnum_fp (FILE *in, char *title, int *natoms)
903 fgets2 (title, STRLEN, in);
904 fgets2 (line, STRLEN, in);
905 if (sscanf (line, "%d", natoms) != 1)
907 gmx_fatal(FARGS, "gro file does not have the number of atoms on the second line");
911 static void get_coordnum (const char *infile, int *natoms)
916 in = gmx_fio_fopen(infile, "r");
917 get_coordnum_fp(in, title, natoms);
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)
926 char line[STRLEN+1], *ptr;
928 double x1, y1, z1, x2, y2, z2;
930 int natoms, i, m, resnr, newres, oldres, ddist, c;
931 gmx_bool bFirst, bVel;
935 oldres = NOTSET; /* Unlikely number for the first residue! */
938 /* Read the title and number of atoms */
939 get_coordnum_fp(in, title, &natoms);
941 if (natoms > atoms->nr)
943 gmx_fatal(FARGS, "gro file contains more atoms (%d) than expected (%d)",
946 else if (natoms < atoms->nr)
948 fprintf(stderr, "Warning: gro file contains less atoms (%d) than expected"
949 " (%d)\n", natoms, atoms->nr);
956 /* just pray the arrays are big enough */
957 for (i = 0; (i < natoms); i++)
959 if ((fgets2 (line, STRLEN, in)) == NULL)
961 gmx_fatal(FARGS, "Unexpected end of file in file %s at line %d",
964 if (strlen(line) < 39)
966 gmx_fatal(FARGS, "Invalid line in %s for atom %d:\n%s", infile, i+1, line);
969 /* determine read precision from distance between periods
974 p1 = strchr(line, '.');
977 gmx_fatal(FARGS, "A coordinate in file %s does not contain a '.'", infile);
979 p2 = strchr(&p1[1], '.');
982 gmx_fatal(FARGS, "A coordinate in file %s does not contain a '.'", infile);
987 p3 = strchr(&p2[1], '.');
990 gmx_fatal(FARGS, "A coordinate in file %s does not contain a '.'", infile);
993 if (p3 - p2 != ddist)
995 gmx_fatal(FARGS, "The spacing of the decimal points in file %s is not consistent for x, y and z", infile);
1000 memcpy(name, line, 5);
1002 sscanf(name, "%d", &resnr);
1003 memcpy(name, line+5, 5);
1005 if (resnr != oldres)
1009 if (newres >= natoms)
1011 gmx_fatal(FARGS, "More residues than atoms in %s (natoms = %d)",
1014 atoms->atom[i].resind = newres;
1015 t_atoms_set_resinfo(atoms, i, symtab, name, resnr, ' ', 0, ' ');
1019 atoms->atom[i].resind = newres;
1023 memcpy(name, line+10, 5);
1024 atoms->atomname[i] = put_symtab(symtab, name);
1026 /* eventueel controle atomnumber met i+1 */
1028 /* coordinates (start after residue data) */
1030 /* Read fixed format */
1031 for (m = 0; m < DIM; m++)
1033 for (c = 0; (c < ddist && ptr[0]); c++)
1039 if (sscanf (buf, "%lf %lf", &x1, &x2) != 1)
1041 gmx_fatal(FARGS, "Something is wrong in the coordinate formatting of file %s. Note that gro is fixed format (see the manual)", infile);
1049 /* velocities (start after residues and coordinates) */
1052 /* Read fixed format */
1053 for (m = 0; m < DIM; m++)
1055 for (c = 0; (c < ddist && ptr[0]); c++)
1061 if (sscanf (buf, "%lf", &x1) != 1)
1073 atoms->nres = newres + 1;
1076 fgets2 (line, STRLEN, in);
1077 if (sscanf (line, "%lf%lf%lf", &x1, &y1, &z1) != 3)
1079 gmx_warning("Bad box in file %s", infile);
1081 /* Generate a cubic box */
1082 for (m = 0; (m < DIM); m++)
1084 xmin[m] = xmax[m] = x[0][m];
1086 for (i = 1; (i < atoms->nr); i++)
1088 for (m = 0; (m < DIM); m++)
1090 xmin[m] = min(xmin[m], x[i][m]);
1091 xmax[m] = max(xmax[m], x[i][m]);
1094 for (i = 0; i < DIM; i++)
1096 for (m = 0; m < DIM; m++)
1101 for (m = 0; (m < DIM); m++)
1103 box[m][m] = (xmax[m]-xmin[m]);
1105 fprintf(stderr, "Generated a cubic box %8.3f x %8.3f x %8.3f\n",
1106 box[XX][XX], box[YY][YY], box[ZZ][ZZ]);
1110 /* We found the first three values, the diagonal elements */
1114 if (sscanf (line, "%*f%*f%*f%lf%lf%lf%lf%lf%lf",
1115 &x1, &y1, &z1, &x2, &y2, &z2) != 6)
1117 x1 = y1 = z1 = x2 = y2 = z2 = 0.0;
1130 static void read_whole_conf(const char *infile, char *title,
1131 t_atoms *atoms, rvec x[], rvec *v, matrix box)
1138 in = gmx_fio_fopen(infile, "r");
1140 open_symtab(&symtab);
1141 get_w_conf(in, infile, title, &symtab, atoms, &ndec, x, v, box);
1142 /* We can't free the symbols, as they are still used in atoms, so
1143 * the only choice is to leak them. */
1144 free_symtab(&symtab);
1149 static gmx_bool gmx_one_before_eof(FILE *fp)
1154 if ((beof = fread(data, 1, 1, fp)) == 1)
1156 gmx_fseek(fp, -1, SEEK_CUR);
1161 gmx_bool gro_next_x_or_v(FILE *status, t_trxframe *fr)
1165 char title[STRLEN], *p;
1169 if (gmx_one_before_eof(status))
1174 open_symtab(&symtab);
1175 atoms.nr = fr->natoms;
1176 snew(atoms.atom, fr->natoms);
1177 atoms.nres = fr->natoms;
1178 snew(atoms.resinfo, fr->natoms);
1179 snew(atoms.atomname, fr->natoms);
1181 fr->bV = get_w_conf(status, title, title, &symtab, &atoms, &ndec, fr->x, fr->v, fr->box);
1184 /* prec = 10^ndec: */
1185 for (i = 0; i < ndec; i++)
1195 sfree(atoms.resinfo);
1196 sfree(atoms.atomname);
1197 done_symtab(&symtab);
1199 if ((p = strstr(title, "t=")) != NULL)
1202 if (sscanf(p, "%lf", &tt) == 1)
1214 if (atoms.nr != fr->natoms)
1216 gmx_fatal(FARGS, "Number of atoms in gro frame (%d) doesn't match the number in the previous frame (%d)", atoms.nr, fr->natoms);
1222 int gro_first_x_or_v(FILE *status, t_trxframe *fr)
1227 fprintf(stderr, "Reading frames from gro file");
1228 get_coordnum_fp(status, title, &fr->natoms);
1230 fprintf(stderr, " '%s', %d atoms.\n", title, fr->natoms);
1233 if (fr->natoms == 0)
1235 gmx_file("No coordinates in gro file");
1238 snew(fr->x, fr->natoms);
1239 snew(fr->v, fr->natoms);
1240 gro_next_x_or_v(status, fr);
1245 static void make_hconf_format(int pr, gmx_bool bVel, char format[])
1249 /* build format string for printing,
1250 something like "%8.3f" for x and "%8.4f" for v */
1263 sprintf(format, "%%%d.%df%%%d.%df%%%d.%df%%%d.%df%%%d.%df%%%d.%df\n",
1264 l, pr, l, pr, l, pr, l, vpr, l, vpr, l, vpr);
1268 sprintf(format, "%%%d.%df%%%d.%df%%%d.%df\n", l, pr, l, pr, l, pr);
1273 static void write_hconf_box(FILE *out, int pr, matrix box)
1284 if (box[XX][YY] || box[XX][ZZ] || box[YY][XX] || box[YY][ZZ] ||
1285 box[ZZ][XX] || box[ZZ][YY])
1287 sprintf(format, "%%%d.%df%%%d.%df%%%d.%df"
1288 "%%%d.%df%%%d.%df%%%d.%df%%%d.%df%%%d.%df%%%d.%df\n",
1289 l, pr, l, pr, l, pr, l, pr, l, pr, l, pr, l, pr, l, pr, l, pr);
1290 fprintf(out, format,
1291 box[XX][XX], box[YY][YY], box[ZZ][ZZ],
1292 box[XX][YY], box[XX][ZZ], box[YY][XX],
1293 box[YY][ZZ], box[ZZ][XX], box[ZZ][YY]);
1297 sprintf(format, "%%%d.%df%%%d.%df%%%d.%df\n", l, pr, l, pr, l, pr);
1298 fprintf(out, format,
1299 box[XX][XX], box[YY][YY], box[ZZ][ZZ]);
1303 void write_hconf_indexed_p(FILE *out, const char *title, t_atoms *atoms,
1304 int nx, const atom_id index[], int pr,
1305 rvec *x, rvec *v, matrix box)
1307 char resnm[6], nm[6], format[100];
1308 int ai, i, resind, resnr;
1310 bromacs(format, 99);
1311 fprintf (out, "%s\n", (title && title[0]) ? title : format);
1312 fprintf (out, "%5d\n", nx);
1314 make_hconf_format(pr, v != NULL, format);
1316 for (i = 0; (i < nx); i++)
1320 resind = atoms->atom[ai].resind;
1321 strncpy(resnm, " ??? ", sizeof(resnm)-1);
1322 if (resind < atoms->nres)
1324 strncpy(resnm, *atoms->resinfo[resind].name, sizeof(resnm)-1);
1325 resnr = atoms->resinfo[resind].nr;
1329 strncpy(resnm, " ??? ", sizeof(resnm)-1);
1335 strncpy(nm, *atoms->atomname[ai], sizeof(nm)-1);
1339 strncpy(nm, " ??? ", sizeof(nm)-1);
1342 fprintf(out, "%5d%-5.5s%5.5s%5d", resnr%100000, resnm, nm, (ai+1)%100000);
1343 /* next fprintf uses built format string */
1346 fprintf(out, format,
1347 x[ai][XX], x[ai][YY], x[ai][ZZ], v[ai][XX], v[ai][YY], v[ai][ZZ]);
1351 fprintf(out, format,
1352 x[ai][XX], x[ai][YY], x[ai][ZZ]);
1356 write_hconf_box(out, pr, box);
1361 static void write_hconf_mtop(FILE *out, const char *title, gmx_mtop_t *mtop,
1363 rvec *x, rvec *v, matrix box)
1367 gmx_mtop_atomloop_all_t aloop;
1369 char *atomname, *resname;
1371 bromacs(format, 99);
1372 fprintf (out, "%s\n", (title && title[0]) ? title : format);
1373 fprintf (out, "%5d\n", mtop->natoms);
1375 make_hconf_format(pr, v != NULL, format);
1377 aloop = gmx_mtop_atomloop_all_init(mtop);
1378 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
1380 gmx_mtop_atomloop_all_names(aloop, &atomname, &resnr, &resname);
1382 fprintf(out, "%5d%-5.5s%5.5s%5d",
1383 resnr%100000, resname, atomname, (i+1)%100000);
1384 /* next fprintf uses built format string */
1387 fprintf(out, format,
1388 x[i][XX], x[i][YY], x[i][ZZ], v[i][XX], v[i][YY], v[i][ZZ]);
1392 fprintf(out, format,
1393 x[i][XX], x[i][YY], x[i][ZZ]);
1397 write_hconf_box(out, pr, box);
1402 void write_hconf_p(FILE *out, const char *title, t_atoms *atoms, int pr,
1403 rvec *x, rvec *v, matrix box)
1408 snew(aa, atoms->nr);
1409 for (i = 0; (i < atoms->nr); i++)
1413 write_hconf_indexed_p(out, title, atoms, atoms->nr, aa, pr, x, v, box);
1417 void write_conf_p(const char *outfile, const char *title,
1418 t_atoms *atoms, int pr,
1419 rvec *x, rvec *v, matrix box)
1423 out = gmx_fio_fopen(outfile, "w");
1424 write_hconf_p(out, title, atoms, pr, x, v, box);
1426 gmx_fio_fclose (out);
1429 static void write_conf(const char *outfile, const char *title, t_atoms *atoms,
1430 rvec *x, rvec *v, matrix box)
1432 write_conf_p(outfile, title, atoms, 3, x, v, box);
1435 void write_sto_conf_indexed(const char *outfile, const char *title,
1437 rvec x[], rvec *v, int ePBC, matrix box,
1438 atom_id nindex, atom_id index[])
1444 ftp = fn2ftp(outfile);
1448 out = gmx_fio_fopen(outfile, "w");
1449 write_hconf_indexed_p(out, title, atoms, nindex, index, 3, x, v, box);
1450 gmx_fio_fclose(out);
1453 clear_trxframe(&fr, TRUE);
1456 fr.natoms = atoms->nr;
1467 copy_mat(box, fr.box);
1468 out = gmx_fio_fopen(outfile, "w");
1469 write_g96_conf(out, &fr, nindex, index);
1470 gmx_fio_fclose(out);
1476 out = gmx_fio_fopen(outfile, "w");
1477 write_pdbfile_indexed(out, title, atoms, x, ePBC, box, ' ', -1, nindex, index, NULL, TRUE);
1478 gmx_fio_fclose(out);
1481 out = gmx_fio_fopen(outfile, "w");
1482 write_espresso_conf_indexed(out, title, atoms, nindex, index, x, v, box);
1483 gmx_fio_fclose(out);
1486 gmx_fatal(FARGS, "Sorry, can not write a topology to %s", outfile);
1489 gmx_incons("Not supported in write_sto_conf_indexed");
1493 void write_sto_conf(const char *outfile, const char *title, t_atoms *atoms,
1494 rvec x[], rvec *v, int ePBC, matrix box)
1500 ftp = fn2ftp(outfile);
1504 write_conf(outfile, title, atoms, x, v, box);
1507 clear_trxframe(&fr, TRUE);
1510 fr.natoms = atoms->nr;
1521 copy_mat(box, fr.box);
1522 out = gmx_fio_fopen(outfile, "w");
1523 write_g96_conf(out, &fr, -1, NULL);
1524 gmx_fio_fclose(out);
1529 out = gmx_fio_fopen(outfile, "w");
1530 write_pdbfile(out, title, atoms, x, ePBC, box, ' ', -1, NULL, TRUE);
1531 gmx_fio_fclose(out);
1534 out = gmx_fio_fopen(outfile, "w");
1535 write_espresso_conf_indexed(out, title, atoms, atoms->nr, NULL, x, v, box);
1536 gmx_fio_fclose(out);
1539 gmx_fatal(FARGS, "Sorry, can not write a topology to %s", outfile);
1542 gmx_incons("Not supported in write_sto_conf");
1546 void write_sto_conf_mtop(const char *outfile, const char *title,
1548 rvec x[], rvec *v, int ePBC, matrix box)
1554 ftp = fn2ftp(outfile);
1558 out = gmx_fio_fopen(outfile, "w");
1559 write_hconf_mtop(out, title, mtop, 3, x, v, box);
1560 gmx_fio_fclose(out);
1563 /* This is a brute force approach which requires a lot of memory.
1564 * We should implement mtop versions of all writing routines.
1566 atoms = gmx_mtop_global_atoms(mtop);
1568 write_sto_conf(outfile, title, &atoms, x, v, ePBC, box);
1575 void get_stx_coordnum(const char *infile, int *natoms)
1578 int ftp, tpxver, tpxgen;
1580 char g96_line[STRLEN+1];
1582 ftp = fn2ftp(infile);
1583 range_check(ftp, 0, efNR);
1587 get_coordnum(infile, natoms);
1590 in = gmx_fio_fopen(infile, "r");
1597 *natoms = read_g96_conf(in, infile, &fr, g96_line);
1603 in = gmx_fio_fopen(infile, "r");
1604 get_pdb_coordnum(in, natoms);
1608 *natoms = get_espresso_coordnum(infile);
1614 read_tpxheader(infile, &tpx, TRUE, &tpxver, &tpxgen);
1615 *natoms = tpx.natoms;
1619 gmx_fatal(FARGS, "File type %s not supported in get_stx_coordnum",
1624 void read_stx_conf(const char *infile, char *title, t_atoms *atoms,
1625 rvec x[], rvec *v, int *ePBC, matrix box)
1634 char g96_line[STRLEN+1];
1638 fprintf(stderr, "Warning: Number of atoms in %s is 0\n", infile);
1640 else if (atoms->atom == NULL)
1642 gmx_mem("Uninitialized array atom");
1650 ftp = fn2ftp(infile);
1654 read_whole_conf(infile, title, atoms, x, v, box);
1658 fr.natoms = atoms->nr;
1663 in = gmx_fio_fopen(infile, "r");
1664 read_g96_conf(in, infile, &fr, g96_line);
1666 copy_mat(fr.box, box);
1667 strncpy(title, fr.title, STRLEN);
1672 read_pdb_conf(infile, title, atoms, x, ePBC, box, TRUE, NULL);
1675 read_espresso_conf(infile, atoms, x, v, box);
1679 i = read_tpx(infile, NULL, box, &natoms, x, v, NULL, mtop);
1685 strcpy(title, *(mtop->name));
1687 /* Free possibly allocated memory */
1690 *atoms = gmx_mtop_global_atoms(mtop);
1691 top = gmx_mtop_t_to_t_topology(mtop);
1692 tpx_make_chain_identifiers(atoms, &top.mols);
1695 /* The strings in the symtab are still in use in the returned t_atoms
1696 * structure, so we should not free them. But there is no place to put the
1697 * symbols; the only choice is to leak the memory...
1698 * So we clear the symbol table before freeing the topology structure. */
1699 free_symtab(&top.symtab);
1704 gmx_incons("Not supported in read_stx_conf");