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.
47 #include "gromacs/utility/futil.h"
54 #include "gromacs/topology/mtop_util.h"
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"
68 static int read_g96_pos(char line[], t_symtab *symtab,
69 FILE *fp, const char *infile,
74 int nwanted, natoms, atnr, resnr = 0, oldres, newres, shift;
75 char anm[STRLEN], resnm[STRLEN];
96 oldres = -666; /* Unlikely number for the first residue! */
98 while (!bEnd && fgets2(line, STRLEN, fp))
100 bEnd = (strncmp(line, "END", 3) == 0);
101 if (!bEnd && (line[0] != '#'))
103 if (sscanf(line+shift, "%15lf%15lf%15lf", &db1, &db2, &db3) != 3)
105 gmx_fatal(FARGS, "Did not find 3 coordinates for atom %d in %s\n",
108 if ((nwanted != -1) && (natoms >= nwanted))
111 "Found more coordinates (%d) in %s than expected %d\n",
112 natoms, infile, nwanted);
117 (sscanf(line, "%5d%c%5s%c%5s%7d", &resnr, &c1, resnm, &c2, anm, &atnr)
127 strncpy(resnm, "???", sizeof(resnm)-1);
129 strncpy(anm, "???", sizeof(anm)-1);
131 atoms->atomname[natoms] = put_symtab(symtab, anm);
136 if (newres >= atoms->nr)
138 gmx_fatal(FARGS, "More residues than atoms in %s (natoms = %d)",
141 atoms->atom[natoms].resind = newres;
142 if (newres+1 > atoms->nres)
144 atoms->nres = newres+1;
146 t_atoms_set_resinfo(atoms, natoms, symtab, resnm, resnr, ' ', 0, ' ');
150 atoms->atom[natoms].resind = newres;
155 fr->x[natoms][0] = db1;
156 fr->x[natoms][1] = db2;
157 fr->x[natoms][2] = db3;
162 if ((nwanted != -1) && natoms != nwanted)
165 "Warning: found less coordinates (%d) in %s than expected %d\n",
166 natoms, infile, nwanted);
175 static int read_g96_vel(char line[], FILE *fp, const char *infile,
179 int nwanted, natoms = -1, shift;
180 double db1, db2, db3;
182 nwanted = fr->natoms;
186 if (strcmp(line, "VELOCITYRED") == 0)
196 while (!bEnd && fgets2(line, STRLEN, fp))
198 bEnd = (strncmp(line, "END", 3) == 0);
199 if (!bEnd && (line[0] != '#'))
201 if (sscanf(line+shift, "%15lf%15lf%15lf", &db1, &db2, &db3) != 3)
203 gmx_fatal(FARGS, "Did not find 3 velocities for atom %d in %s",
206 if ((nwanted != -1) && (natoms >= nwanted))
208 gmx_fatal(FARGS, "Found more velocities (%d) in %s than expected %d\n",
209 natoms, infile, nwanted);
213 fr->v[natoms][0] = db1;
214 fr->v[natoms][1] = db2;
215 fr->v[natoms][2] = db3;
220 if ((nwanted != -1) && (natoms != nwanted))
223 "Warning: found less velocities (%d) in %s than expected %d\n",
224 natoms, infile, nwanted);
231 int read_g96_conf(FILE *fp, const char *infile, t_trxframe *fr, char *line)
233 t_symtab *symtab = NULL;
234 gmx_bool bAtStart, bTime, bAtoms, bPos, bVel, bBox, bEnd, bFinished;
236 double db1, db2, db3, db4, db5, db6, db7, db8, db9;
238 bAtStart = (ftell(fp) == 0);
240 clear_trxframe(fr, FALSE);
252 while (!fr->bTitle && fgets2(line, STRLEN, fp))
254 fr->bTitle = (strcmp(line, "TITLE") == 0);
256 if (fr->title == NULL)
258 fgets2(line, STRLEN, fp);
259 fr->title = gmx_strdup(line);
262 while (!bEnd && fgets2(line, STRLEN, fp))
264 bEnd = (strcmp(line, "END") == 0);
266 fgets2(line, STRLEN, fp);
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). */
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);
283 if (!fr->bTime && !fr->bX)
289 bFinished = (fgets2(line, STRLEN, fp) == NULL);
291 while (!bFinished && (line[0] == '#'));
292 sscanf(line, "%15d%15lf", &(fr->step), &db1);
306 natoms = read_g96_pos(line, symtab, fp, infile, fr);
316 natoms = read_g96_vel(line, fp, infile, fr);
323 while (!bEnd && fgets2(line, STRLEN, fp))
325 bEnd = (strncmp(line, "END", 3) == 0);
326 if (!bEnd && (line[0] != '#'))
328 nbp = sscanf(line, "%15lf%15lf%15lf%15lf%15lf%15lf%15lf%15lf%15lf",
329 &db1, &db2, &db3, &db4, &db5, &db6, &db7, &db8, &db9);
332 gmx_fatal(FARGS, "Found a BOX line, but no box in %s", infile);
334 fr->box[XX][XX] = db1;
335 fr->box[YY][YY] = db2;
336 fr->box[ZZ][ZZ] = db3;
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;
351 while (!bFinished && fgets2(line, STRLEN, fp));
360 void write_g96_conf(FILE *out, t_trxframe *fr,
361 int nindex, const atom_id *index)
379 fprintf(out, "TITLE\n%s\nEND\n", fr->title);
381 if (fr->bStep || fr->bTime)
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);
390 fprintf(out, "POSITION\n");
391 for (i = 0; i < nout; i++)
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]);
410 fprintf(out, "POSITIONRED\n");
411 for (i = 0; i < nout; i++)
421 fprintf(out, "%15.9f%15.9f%15.9f\n",
422 fr->x[a][XX], fr->x[a][YY], fr->x[a][ZZ]);
425 fprintf(out, "END\n");
431 fprintf(out, "VELOCITY\n");
432 for (i = 0; i < nout; i++)
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]);
451 fprintf(out, "VELOCITYRED\n");
452 for (i = 0; i < nout; i++)
462 fprintf(out, "%15.9f%15.9f%15.9f\n",
463 fr->v[a][XX], fr->v[a][YY], fr->v[a][ZZ]);
466 fprintf(out, "END\n");
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])
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]);
481 fprintf(out, "END\n");
485 static int get_espresso_word(FILE *fp, char word[])
497 if (i == ' ' || i == '\n' || i == '\t')
522 word[nc++] = (char)i;
526 while (i != EOF && ret == 0);
530 /* printf("'%s'\n",word); */
535 static int check_open_parenthesis(FILE *fp, int r,
536 const char *infile, const char *keyword)
548 r = get_espresso_word(fp, word);
555 gmx_fatal(FARGS, "Expected '{' after '%s' in file '%s'",
563 static int check_close_parenthesis(FILE *fp, int r,
564 const char *infile, const char *keyword)
576 r = get_espresso_word(fp, word);
583 gmx_fatal(FARGS, "Expected '}' after section '%s' in file '%s'",
592 espID, espPOS, espTYPE, espQ, espV, espF, espMOLECULE, espNR
594 const char *esp_prop[espNR] = {
595 "id", "pos", "type", "q", "v", "f",
599 static void read_espresso_conf(const char *infile,
600 t_atoms *atoms, rvec x[], rvec *v, matrix box)
602 t_symtab *symtab = NULL;
604 char word[STRLEN], buf[STRLEN];
605 int natoms, level, npar, r, nprop, p, i, m, molnr;
608 gmx_bool bFoundParticles, bFoundProp, bFoundVariable, bMol;
618 fp = gmx_fio_fopen(infile, "r");
620 bFoundParticles = FALSE;
621 bFoundVariable = FALSE;
624 while ((r = get_espresso_word(fp, word)))
626 if (level == 1 && strcmp(word, "particles") == 0 && !bFoundParticles)
628 bFoundParticles = TRUE;
629 level += check_open_parenthesis(fp, r, infile, "particles");
631 while (level == 2 && (r = get_espresso_word(fp, word)))
634 for (p = 0; p < espNR; p++)
636 if (strcmp(word, esp_prop[p]) == 0)
640 /* printf(" prop[%d] = %s\n",nprop-1,esp_prop[prop[nprop-1]]); */
643 if (!bFoundProp && word[0] != '}')
645 gmx_fatal(FARGS, "Can not read Espresso files with particle property '%s'", word);
647 if (bFoundProp && p == espMOLECULE)
658 while (level > 0 && (r = get_espresso_word(fp, word)))
670 for (p = 0; p < nprop; p++)
675 r = get_espresso_word(fp, word);
679 for (m = 0; m < 3; m++)
681 r = get_espresso_word(fp, word);
682 sscanf(word, "%lf", &d);
687 r = get_espresso_word(fp, word);
688 atoms->atom[i].type = strtol(word, NULL, 10);
691 r = get_espresso_word(fp, word);
692 sscanf(word, "%lf", &d);
693 atoms->atom[i].q = d;
696 for (m = 0; m < 3; m++)
698 r = get_espresso_word(fp, word);
699 sscanf(word, "%lf", &d);
704 for (m = 0; m < 3; m++)
706 r = get_espresso_word(fp, word);
711 r = get_espresso_word(fp, word);
712 molnr = strtol(word, NULL, 10);
714 atoms->resinfo[atoms->atom[i-1].resind].nr != molnr)
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? */
725 atoms->atom[i].resind = atoms->atom[i-1].resind;
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);
735 if (i == 0 || atoms->atom[i].resind != atoms->atom[i-1].resind)
737 atoms->resinfo[atoms->atom[i].resind].name =
738 put_symtab(symtab, "MOL");
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)
748 sprintf(buf, "T%c", 'A'+atoms->atom[i].type);
752 sprintf(buf, "T%c%c",
753 'A'+atoms->atom[i].type/26, 'A'+atoms->atom[i].type%26);
755 t_atoms_set_resinfo(atoms, i, symtab, buf, i, ' ', 0, ' ');
765 atoms->nres = atoms->nr;
769 gmx_fatal(FARGS, "Internal inconsistency in Espresso routines, read %d atoms, expected %d atoms", i, atoms->nr);
772 else if (level == 1 && strcmp(word, "variable") == 0 && !bFoundVariable)
774 bFoundVariable = TRUE;
775 level += check_open_parenthesis(fp, r, infile, "variable");
776 while (level == 2 && (r = get_espresso_word(fp, word)))
778 if (level == 2 && strcmp(word, "box_l") == 0)
780 for (m = 0; m < 3; m++)
782 r = get_espresso_word(fp, word);
783 sscanf(word, "%lf", &d);
786 level += check_close_parenthesis(fp, r, infile, "box_l");
800 if (!bFoundParticles)
802 fprintf(stderr, "Did not find a particles section in Espresso file '%s'\n",
809 static int get_espresso_coordnum(const char *infile)
813 int natoms, level, r;
814 gmx_bool bFoundParticles;
818 fp = gmx_fio_fopen(infile, "r");
820 bFoundParticles = FALSE;
822 while ((r = get_espresso_word(fp, word)) && !bFoundParticles)
824 if (level == 1 && strcmp(word, "particles") == 0 && !bFoundParticles)
826 bFoundParticles = TRUE;
827 level += check_open_parenthesis(fp, r, infile, "particles");
828 while (level > 0 && (r = get_espresso_word(fp, word)))
853 if (!bFoundParticles)
855 fprintf(stderr, "Did not find a particles section in Espresso file '%s'\n",
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)
870 fprintf(out, "# %s\n", title);
873 gmx_warning("The Espresso format does not support triclinic unit-cells");
875 fprintf(out, "{variable {box_l %f %f %f}}\n", box[0][0], box[1][1], box[2][2]);
877 fprintf(out, "{particles {id pos type q%s}\n", v ? " v" : "");
878 for (i = 0; i < nx; i++)
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);
893 fprintf(out, " %f %f %f", v[j][XX], v[j][YY], v[j][ZZ]);
900 static void get_coordnum_fp (FILE *in, char *title, int *natoms)
904 fgets2 (title, STRLEN, in);
905 fgets2 (line, STRLEN, in);
906 if (sscanf (line, "%d", natoms) != 1)
908 gmx_fatal(FARGS, "gro file does not have the number of atoms on the second line");
912 static void get_coordnum (const char *infile, int *natoms)
917 in = gmx_fio_fopen(infile, "r");
918 get_coordnum_fp(in, title, natoms);
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)
927 char line[STRLEN+1], *ptr;
929 double x1, y1, z1, x2, y2, z2;
931 int natoms, i, m, resnr, newres, oldres, ddist, c;
932 gmx_bool bFirst, bVel;
936 oldres = NOTSET; /* Unlikely number for the first residue! */
939 /* Read the title and number of atoms */
940 get_coordnum_fp(in, title, &natoms);
942 if (natoms > atoms->nr)
944 gmx_fatal(FARGS, "gro file contains more atoms (%d) than expected (%d)",
947 else if (natoms < atoms->nr)
949 fprintf(stderr, "Warning: gro file contains less atoms (%d) than expected"
950 " (%d)\n", natoms, atoms->nr);
957 /* just pray the arrays are big enough */
958 for (i = 0; (i < natoms); i++)
960 if ((fgets2 (line, STRLEN, in)) == NULL)
962 gmx_fatal(FARGS, "Unexpected end of file in file %s at line %d",
965 if (strlen(line) < 39)
967 gmx_fatal(FARGS, "Invalid line in %s for atom %d:\n%s", infile, i+1, line);
970 /* determine read precision from distance between periods
975 p1 = strchr(line, '.');
978 gmx_fatal(FARGS, "A coordinate in file %s does not contain a '.'", infile);
980 p2 = strchr(&p1[1], '.');
983 gmx_fatal(FARGS, "A coordinate in file %s does not contain a '.'", infile);
988 p3 = strchr(&p2[1], '.');
991 gmx_fatal(FARGS, "A coordinate in file %s does not contain a '.'", infile);
994 if (p3 - p2 != ddist)
996 gmx_fatal(FARGS, "The spacing of the decimal points in file %s is not consistent for x, y and z", infile);
1001 memcpy(name, line, 5);
1003 sscanf(name, "%d", &resnr);
1004 memcpy(name, line+5, 5);
1006 if (resnr != oldres)
1010 if (newres >= natoms)
1012 gmx_fatal(FARGS, "More residues than atoms in %s (natoms = %d)",
1015 atoms->atom[i].resind = newres;
1016 t_atoms_set_resinfo(atoms, i, symtab, name, resnr, ' ', 0, ' ');
1020 atoms->atom[i].resind = newres;
1024 memcpy(name, line+10, 5);
1025 atoms->atomname[i] = put_symtab(symtab, name);
1027 /* eventueel controle atomnumber met i+1 */
1029 /* coordinates (start after residue data) */
1031 /* Read fixed format */
1032 for (m = 0; m < DIM; m++)
1034 for (c = 0; (c < ddist && ptr[0]); c++)
1040 if (sscanf (buf, "%lf %lf", &x1, &x2) != 1)
1042 gmx_fatal(FARGS, "Something is wrong in the coordinate formatting of file %s. Note that gro is fixed format (see the manual)", infile);
1050 /* velocities (start after residues and coordinates) */
1053 /* Read fixed format */
1054 for (m = 0; m < DIM; m++)
1056 for (c = 0; (c < ddist && ptr[0]); c++)
1062 if (sscanf (buf, "%lf", &x1) != 1)
1074 atoms->nres = newres + 1;
1077 fgets2 (line, STRLEN, in);
1078 if (sscanf (line, "%lf%lf%lf", &x1, &y1, &z1) != 3)
1080 gmx_warning("Bad box in file %s", infile);
1082 /* Generate a cubic box */
1083 for (m = 0; (m < DIM); m++)
1085 xmin[m] = xmax[m] = x[0][m];
1087 for (i = 1; (i < atoms->nr); i++)
1089 for (m = 0; (m < DIM); m++)
1091 xmin[m] = min(xmin[m], x[i][m]);
1092 xmax[m] = max(xmax[m], x[i][m]);
1095 for (i = 0; i < DIM; i++)
1097 for (m = 0; m < DIM; m++)
1102 for (m = 0; (m < DIM); m++)
1104 box[m][m] = (xmax[m]-xmin[m]);
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]);
1111 /* We found the first three values, the diagonal elements */
1115 if (sscanf (line, "%*f%*f%*f%lf%lf%lf%lf%lf%lf",
1116 &x1, &y1, &z1, &x2, &y2, &z2) != 6)
1118 x1 = y1 = z1 = x2 = y2 = z2 = 0.0;
1131 static void read_whole_conf(const char *infile, char *title,
1132 t_atoms *atoms, rvec x[], rvec *v, matrix box)
1139 in = gmx_fio_fopen(infile, "r");
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);
1150 static gmx_bool gmx_one_before_eof(FILE *fp)
1155 if ((beof = fread(data, 1, 1, fp)) == 1)
1157 gmx_fseek(fp, -1, SEEK_CUR);
1162 gmx_bool gro_next_x_or_v(FILE *status, t_trxframe *fr)
1166 char title[STRLEN], *p;
1170 if (gmx_one_before_eof(status))
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);
1182 fr->bV = get_w_conf(status, title, title, &symtab, &atoms, &ndec, fr->x, fr->v, fr->box);
1185 /* prec = 10^ndec: */
1186 for (i = 0; i < ndec; i++)
1196 sfree(atoms.resinfo);
1197 sfree(atoms.atomname);
1198 done_symtab(&symtab);
1200 if ((p = strstr(title, "t=")) != NULL)
1203 if (sscanf(p, "%lf", &tt) == 1)
1215 if (atoms.nr != fr->natoms)
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);
1223 int gro_first_x_or_v(FILE *status, t_trxframe *fr)
1228 fprintf(stderr, "Reading frames from gro file");
1229 get_coordnum_fp(status, title, &fr->natoms);
1231 fprintf(stderr, " '%s', %d atoms.\n", title, fr->natoms);
1234 if (fr->natoms == 0)
1236 gmx_file("No coordinates in gro file");
1239 snew(fr->x, fr->natoms);
1240 snew(fr->v, fr->natoms);
1241 gro_next_x_or_v(status, fr);
1246 static void make_hconf_format(int pr, gmx_bool bVel, char format[])
1250 /* build format string for printing,
1251 something like "%8.3f" for x and "%8.4f" for v */
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);
1269 sprintf(format, "%%%d.%df%%%d.%df%%%d.%df\n", l, pr, l, pr, l, pr);
1274 static void write_hconf_box(FILE *out, int pr, matrix box)
1285 if (box[XX][YY] || box[XX][ZZ] || box[YY][XX] || box[YY][ZZ] ||
1286 box[ZZ][XX] || box[ZZ][YY])
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]);
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]);
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)
1308 char resnm[6], nm[6], format[100];
1309 int ai, i, resind, resnr;
1311 bromacs(format, 99);
1312 fprintf (out, "%s\n", (title && title[0]) ? title : format);
1313 fprintf (out, "%5d\n", nx);
1315 make_hconf_format(pr, v != NULL, format);
1317 for (i = 0; (i < nx); i++)
1321 resind = atoms->atom[ai].resind;
1322 strncpy(resnm, " ??? ", sizeof(resnm)-1);
1323 if (resind < atoms->nres)
1325 strncpy(resnm, *atoms->resinfo[resind].name, sizeof(resnm)-1);
1326 resnr = atoms->resinfo[resind].nr;
1330 strncpy(resnm, " ??? ", sizeof(resnm)-1);
1336 strncpy(nm, *atoms->atomname[ai], sizeof(nm)-1);
1340 strncpy(nm, " ??? ", sizeof(nm)-1);
1343 fprintf(out, "%5d%-5.5s%5.5s%5d", resnr%100000, resnm, nm, (ai+1)%100000);
1344 /* next fprintf uses built format string */
1347 fprintf(out, format,
1348 x[ai][XX], x[ai][YY], x[ai][ZZ], v[ai][XX], v[ai][YY], v[ai][ZZ]);
1352 fprintf(out, format,
1353 x[ai][XX], x[ai][YY], x[ai][ZZ]);
1357 write_hconf_box(out, pr, box);
1362 static void write_hconf_mtop(FILE *out, const char *title, gmx_mtop_t *mtop,
1364 rvec *x, rvec *v, matrix box)
1368 gmx_mtop_atomloop_all_t aloop;
1370 char *atomname, *resname;
1372 bromacs(format, 99);
1373 fprintf (out, "%s\n", (title && title[0]) ? title : format);
1374 fprintf (out, "%5d\n", mtop->natoms);
1376 make_hconf_format(pr, v != NULL, format);
1378 aloop = gmx_mtop_atomloop_all_init(mtop);
1379 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
1381 gmx_mtop_atomloop_all_names(aloop, &atomname, &resnr, &resname);
1383 fprintf(out, "%5d%-5.5s%5.5s%5d",
1384 resnr%100000, resname, atomname, (i+1)%100000);
1385 /* next fprintf uses built format string */
1388 fprintf(out, format,
1389 x[i][XX], x[i][YY], x[i][ZZ], v[i][XX], v[i][YY], v[i][ZZ]);
1393 fprintf(out, format,
1394 x[i][XX], x[i][YY], x[i][ZZ]);
1398 write_hconf_box(out, pr, box);
1403 void write_hconf_p(FILE *out, const char *title, t_atoms *atoms, int pr,
1404 rvec *x, rvec *v, matrix box)
1409 snew(aa, atoms->nr);
1410 for (i = 0; (i < atoms->nr); i++)
1414 write_hconf_indexed_p(out, title, atoms, atoms->nr, aa, pr, x, v, box);
1418 void write_conf_p(const char *outfile, const char *title,
1419 t_atoms *atoms, int pr,
1420 rvec *x, rvec *v, matrix box)
1424 out = gmx_fio_fopen(outfile, "w");
1425 write_hconf_p(out, title, atoms, pr, x, v, box);
1427 gmx_fio_fclose (out);
1430 static void write_conf(const char *outfile, const char *title, t_atoms *atoms,
1431 rvec *x, rvec *v, matrix box)
1433 write_conf_p(outfile, title, atoms, 3, x, v, box);
1436 void write_sto_conf_indexed(const char *outfile, const char *title,
1438 rvec x[], rvec *v, int ePBC, matrix box,
1439 atom_id nindex, atom_id index[])
1445 ftp = fn2ftp(outfile);
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);
1454 clear_trxframe(&fr, TRUE);
1457 fr.natoms = atoms->nr;
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);
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);
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);
1489 gmx_fatal(FARGS, "Sorry, can not write a topology to %s", outfile);
1492 gmx_incons("Not supported in write_sto_conf_indexed");
1496 void write_sto_conf(const char *outfile, const char *title, t_atoms *atoms,
1497 rvec x[], rvec *v, int ePBC, matrix box)
1503 ftp = fn2ftp(outfile);
1507 write_conf(outfile, title, atoms, x, v, box);
1510 clear_trxframe(&fr, TRUE);
1513 fr.natoms = atoms->nr;
1524 copy_mat(box, fr.box);
1525 out = gmx_fio_fopen(outfile, "w");
1526 write_g96_conf(out, &fr, -1, NULL);
1527 gmx_fio_fclose(out);
1532 out = gmx_fio_fopen(outfile, "w");
1533 write_pdbfile(out, title, atoms, x, ePBC, box, ' ', -1, NULL, TRUE);
1534 gmx_fio_fclose(out);
1537 out = gmx_fio_fopen(outfile, "w");
1538 write_espresso_conf_indexed(out, title, atoms, atoms->nr, NULL, x, v, box);
1539 gmx_fio_fclose(out);
1544 gmx_fatal(FARGS, "Sorry, can not write a topology to %s", outfile);
1547 gmx_incons("Not supported in write_sto_conf");
1551 void write_sto_conf_mtop(const char *outfile, const char *title,
1553 rvec x[], rvec *v, int ePBC, matrix box)
1559 ftp = fn2ftp(outfile);
1563 out = gmx_fio_fopen(outfile, "w");
1564 write_hconf_mtop(out, title, mtop, 3, x, v, box);
1565 gmx_fio_fclose(out);
1568 /* This is a brute force approach which requires a lot of memory.
1569 * We should implement mtop versions of all writing routines.
1571 atoms = gmx_mtop_global_atoms(mtop);
1573 write_sto_conf(outfile, title, &atoms, x, v, ePBC, box);
1580 void get_stx_coordnum(const char *infile, int *natoms)
1583 int ftp, tpxver, tpxgen;
1585 char g96_line[STRLEN+1];
1587 ftp = fn2ftp(infile);
1588 range_check(ftp, 0, efNR);
1592 get_coordnum(infile, natoms);
1595 in = gmx_fio_fopen(infile, "r");
1602 *natoms = read_g96_conf(in, infile, &fr, g96_line);
1608 in = gmx_fio_fopen(infile, "r");
1609 get_pdb_coordnum(in, natoms);
1613 *natoms = get_espresso_coordnum(infile);
1621 read_tpxheader(infile, &tpx, TRUE, &tpxver, &tpxgen);
1622 *natoms = tpx.natoms;
1626 gmx_fatal(FARGS, "File type %s not supported in get_stx_coordnum",
1631 void read_stx_conf(const char *infile, char *title, t_atoms *atoms,
1632 rvec x[], rvec *v, int *ePBC, matrix box)
1641 char g96_line[STRLEN+1];
1645 fprintf(stderr, "Warning: Number of atoms in %s is 0\n", infile);
1647 else if (atoms->atom == NULL)
1649 gmx_mem("Uninitialized array atom");
1657 ftp = fn2ftp(infile);
1661 read_whole_conf(infile, title, atoms, x, v, box);
1665 fr.natoms = atoms->nr;
1670 in = gmx_fio_fopen(infile, "r");
1671 read_g96_conf(in, infile, &fr, g96_line);
1673 copy_mat(fr.box, box);
1674 strncpy(title, fr.title, STRLEN);
1679 read_pdb_conf(infile, title, atoms, x, ePBC, box, TRUE, NULL);
1682 read_espresso_conf(infile, atoms, x, v, box);
1688 i = read_tpx(infile, NULL, box, &natoms, x, v, NULL, mtop);
1694 strcpy(title, *(mtop->name));
1696 /* Free possibly allocated memory */
1699 *atoms = gmx_mtop_global_atoms(mtop);
1700 top = gmx_mtop_t_to_t_topology(mtop);
1701 tpx_make_chain_identifiers(atoms, &top.mols);
1704 /* The strings in the symtab are still in use in the returned t_atoms
1705 * structure, so we should not free them. But there is no place to put the
1706 * symbols; the only choice is to leak the memory...
1707 * So we clear the symbol table before freeing the topology structure. */
1708 free_symtab(&top.symtab);
1713 gmx_incons("Not supported in read_stx_conf");