3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * GROningen Mixture of Alchemy and Childrens' Stories
55 #include "gmx_fatal.h"
60 #include "mtop_util.h"
65 static int read_g96_pos(char line[], t_symtab *symtab,
66 FILE *fp, const char *infile,
71 int nwanted, natoms, atnr, resnr = 0, oldres, newres, shift;
72 char anm[STRLEN], resnm[STRLEN];
93 oldres = -666; /* Unlikely number for the first residue! */
95 while (!bEnd && fgets2(line, STRLEN, fp))
97 bEnd = (strncmp(line, "END", 3) == 0);
98 if (!bEnd && (line[0] != '#'))
100 if (sscanf(line+shift, "%15lf%15lf%15lf", &db1, &db2, &db3) != 3)
102 gmx_fatal(FARGS, "Did not find 3 coordinates for atom %d in %s\n",
105 if ((nwanted != -1) && (natoms >= nwanted))
108 "Found more coordinates (%d) in %s than expected %d\n",
109 natoms, infile, nwanted);
114 (sscanf(line, "%5d%c%5s%c%5s%7d", &resnr, &c1, resnm, &c2, anm, &atnr)
124 strncpy(resnm, "???", sizeof(resnm)-1);
126 strncpy(anm, "???", sizeof(anm)-1);
128 atoms->atomname[natoms] = put_symtab(symtab, anm);
133 if (newres >= atoms->nr)
135 gmx_fatal(FARGS, "More residues than atoms in %s (natoms = %d)",
138 atoms->atom[natoms].resind = newres;
139 if (newres+1 > atoms->nres)
141 atoms->nres = newres+1;
143 t_atoms_set_resinfo(atoms, natoms, symtab, resnm, resnr, ' ', 0, ' ');
147 atoms->atom[natoms].resind = newres;
152 fr->x[natoms][0] = db1;
153 fr->x[natoms][1] = db2;
154 fr->x[natoms][2] = db3;
159 if ((nwanted != -1) && natoms != nwanted)
162 "Warning: found less coordinates (%d) in %s than expected %d\n",
163 natoms, infile, nwanted);
172 static int read_g96_vel(char line[], FILE *fp, const char *infile,
176 int nwanted, natoms = -1, shift;
177 double db1, db2, db3;
179 nwanted = fr->natoms;
183 if (strcmp(line, "VELOCITYRED") == 0)
193 while (!bEnd && fgets2(line, STRLEN, fp))
195 bEnd = (strncmp(line, "END", 3) == 0);
196 if (!bEnd && (line[0] != '#'))
198 if (sscanf(line+shift, "%15lf%15lf%15lf", &db1, &db2, &db3) != 3)
200 gmx_fatal(FARGS, "Did not find 3 velocities for atom %d in %s",
203 if ((nwanted != -1) && (natoms >= nwanted))
205 gmx_fatal(FARGS, "Found more velocities (%d) in %s than expected %d\n",
206 natoms, infile, nwanted);
210 fr->v[natoms][0] = db1;
211 fr->v[natoms][1] = db2;
212 fr->v[natoms][2] = db3;
217 if ((nwanted != -1) && (natoms != nwanted))
220 "Warning: found less velocities (%d) in %s than expected %d\n",
221 natoms, infile, nwanted);
228 int read_g96_conf(FILE *fp, const char *infile, t_trxframe *fr, char *line)
230 t_symtab *symtab = NULL;
231 gmx_bool bAtStart, bTime, bAtoms, bPos, bVel, bBox, bEnd, bFinished;
233 double db1, db2, db3, db4, db5, db6, db7, db8, db9;
235 bAtStart = (ftell(fp) == 0);
237 clear_trxframe(fr, FALSE);
249 while (!fr->bTitle && fgets2(line, STRLEN, fp))
251 fr->bTitle = (strcmp(line, "TITLE") == 0);
253 if (fr->title == NULL)
255 fgets2(line, STRLEN, fp);
256 fr->title = strdup(line);
259 while (!bEnd && fgets2(line, STRLEN, fp))
261 bEnd = (strcmp(line, "END") == 0);
263 fgets2(line, STRLEN, fp);
266 /* Do not get a line if we are not at the start of the file, *
267 * because without a parameter file we don't know what is in *
268 * the trajectory and we have already read the line in the *
269 * previous call (VERY DIRTY). */
273 bTime = (strcmp(line, "TIMESTEP") == 0);
274 bAtoms = (strcmp(line, "POSITION") == 0);
275 bPos = (bAtoms || (strcmp(line, "POSITIONRED") == 0));
276 bVel = (strncmp(line, "VELOCITY", 8) == 0);
277 bBox = (strcmp(line, "BOX") == 0);
280 if (!fr->bTime && !fr->bX)
286 bFinished = (fgets2(line, STRLEN, fp) == NULL);
288 while (!bFinished && (line[0] == '#'));
289 sscanf(line, "%15d%15lf", &(fr->step), &db1);
303 natoms = read_g96_pos(line, symtab, fp, infile, fr);
313 natoms = read_g96_vel(line, fp, infile, fr);
320 while (!bEnd && fgets2(line, STRLEN, fp))
322 bEnd = (strncmp(line, "END", 3) == 0);
323 if (!bEnd && (line[0] != '#'))
325 nbp = sscanf(line, "%15lf%15lf%15lf%15lf%15lf%15lf%15lf%15lf%15lf",
326 &db1, &db2, &db3, &db4, &db5, &db6, &db7, &db8, &db9);
329 gmx_fatal(FARGS, "Found a BOX line, but no box in %s", infile);
331 fr->box[XX][XX] = db1;
332 fr->box[YY][YY] = db2;
333 fr->box[ZZ][ZZ] = db3;
336 fr->box[XX][YY] = db4;
337 fr->box[XX][ZZ] = db5;
338 fr->box[YY][XX] = db6;
339 fr->box[YY][ZZ] = db7;
340 fr->box[ZZ][XX] = db8;
341 fr->box[ZZ][YY] = db9;
348 while (!bFinished && fgets2(line, STRLEN, fp));
357 void write_g96_conf(FILE *out, t_trxframe *fr,
358 int nindex, const atom_id *index)
376 fprintf(out, "TITLE\n%s\nEND\n", fr->title);
378 if (fr->bStep || fr->bTime)
380 /* Officially the time format is %15.9, which is not enough for 10 ns */
381 fprintf(out, "TIMESTEP\n%15d%15.6f\nEND\n", fr->step, fr->time);
387 fprintf(out, "POSITION\n");
388 for (i = 0; i < nout; i++)
398 fprintf(out, "%5d %-5s %-5s%7d%15.9f%15.9f%15.9f\n",
399 (atoms->resinfo[atoms->atom[a].resind].nr) % 100000,
400 *atoms->resinfo[atoms->atom[a].resind].name,
401 *atoms->atomname[a], (i+1) % 10000000,
402 fr->x[a][XX], fr->x[a][YY], fr->x[a][ZZ]);
407 fprintf(out, "POSITIONRED\n");
408 for (i = 0; i < nout; i++)
418 fprintf(out, "%15.9f%15.9f%15.9f\n",
419 fr->x[a][XX], fr->x[a][YY], fr->x[a][ZZ]);
422 fprintf(out, "END\n");
428 fprintf(out, "VELOCITY\n");
429 for (i = 0; i < nout; i++)
439 fprintf(out, "%5d %-5s %-5s%7d%15.9f%15.9f%15.9f\n",
440 (atoms->resinfo[atoms->atom[a].resind].nr) % 100000,
441 *atoms->resinfo[atoms->atom[a].resind].name,
442 *atoms->atomname[a], (i+1) % 10000000,
443 fr->v[a][XX], fr->v[a][YY], fr->v[a][ZZ]);
448 fprintf(out, "VELOCITYRED\n");
449 for (i = 0; i < nout; i++)
459 fprintf(out, "%15.9f%15.9f%15.9f\n",
460 fr->v[a][XX], fr->v[a][YY], fr->v[a][ZZ]);
463 fprintf(out, "END\n");
467 fprintf(out, "BOX\n");
468 fprintf(out, "%15.9f%15.9f%15.9f",
469 fr->box[XX][XX], fr->box[YY][YY], fr->box[ZZ][ZZ]);
470 if (fr->box[XX][YY] || fr->box[XX][ZZ] || fr->box[YY][XX] ||
471 fr->box[YY][ZZ] || fr->box[ZZ][XX] || fr->box[ZZ][YY])
473 fprintf(out, "%15.9f%15.9f%15.9f%15.9f%15.9f%15.9f",
474 fr->box[XX][YY], fr->box[XX][ZZ], fr->box[YY][XX],
475 fr->box[YY][ZZ], fr->box[ZZ][XX], fr->box[ZZ][YY]);
478 fprintf(out, "END\n");
482 static int get_espresso_word(FILE *fp, char word[])
494 if (i == ' ' || i == '\n' || i == '\t')
519 word[nc++] = (char)i;
523 while (i != EOF && ret == 0);
527 /* printf("'%s'\n",word); */
532 static int check_open_parenthesis(FILE *fp, int r,
533 const char *infile, const char *keyword)
545 r = get_espresso_word(fp, word);
552 gmx_fatal(FARGS, "Expected '{' after '%s' in file '%s'",
560 static int check_close_parenthesis(FILE *fp, int r,
561 const char *infile, const char *keyword)
573 r = get_espresso_word(fp, word);
580 gmx_fatal(FARGS, "Expected '}' after section '%s' in file '%s'",
589 espID, espPOS, espTYPE, espQ, espV, espF, espMOLECULE, espNR
591 const char *esp_prop[espNR] = {
592 "id", "pos", "type", "q", "v", "f",
596 static void read_espresso_conf(const char *infile,
597 t_atoms *atoms, rvec x[], rvec *v, matrix box)
599 t_symtab *symtab = NULL;
601 char word[STRLEN], buf[STRLEN];
602 int natoms, level, npar, r, nprop, p, i, m, molnr;
605 gmx_bool bFoundParticles, bFoundProp, bFoundVariable, bMol;
615 fp = gmx_fio_fopen(infile, "r");
617 bFoundParticles = FALSE;
618 bFoundVariable = FALSE;
621 while ((r = get_espresso_word(fp, word)))
623 if (level == 1 && strcmp(word, "particles") == 0 && !bFoundParticles)
625 bFoundParticles = TRUE;
626 level += check_open_parenthesis(fp, r, infile, "particles");
628 while (level == 2 && (r = get_espresso_word(fp, word)))
631 for (p = 0; p < espNR; p++)
633 if (strcmp(word, esp_prop[p]) == 0)
637 /* printf(" prop[%d] = %s\n",nprop-1,esp_prop[prop[nprop-1]]); */
640 if (!bFoundProp && word[0] != '}')
642 gmx_fatal(FARGS, "Can not read Espresso files with particle property '%s'", word);
644 if (bFoundProp && p == espMOLECULE)
655 while (level > 0 && (r = get_espresso_word(fp, word)))
667 for (p = 0; p < nprop; p++)
672 r = get_espresso_word(fp, word);
676 for (m = 0; m < 3; m++)
678 r = get_espresso_word(fp, word);
679 sscanf(word, "%lf", &d);
684 r = get_espresso_word(fp, word);
685 atoms->atom[i].type = strtol(word, NULL, 10);
688 r = get_espresso_word(fp, word);
689 sscanf(word, "%lf", &d);
690 atoms->atom[i].q = d;
693 for (m = 0; m < 3; m++)
695 r = get_espresso_word(fp, word);
696 sscanf(word, "%lf", &d);
701 for (m = 0; m < 3; m++)
703 r = get_espresso_word(fp, word);
708 r = get_espresso_word(fp, word);
709 molnr = strtol(word, NULL, 10);
711 atoms->resinfo[atoms->atom[i-1].resind].nr != molnr)
713 atoms->atom[i].resind =
714 (i == 0 ? 0 : atoms->atom[i-1].resind+1);
715 atoms->resinfo[atoms->atom[i].resind].nr = molnr;
716 atoms->resinfo[atoms->atom[i].resind].ic = ' ';
717 atoms->resinfo[atoms->atom[i].resind].chainid = ' ';
718 atoms->resinfo[atoms->atom[i].resind].chainnum = molnr; /* Not sure if this is right? */
722 atoms->atom[i].resind = atoms->atom[i-1].resind;
727 /* Generate an atom name from the particle type */
728 sprintf(buf, "T%d", atoms->atom[i].type);
729 atoms->atomname[i] = put_symtab(symtab, buf);
732 if (i == 0 || atoms->atom[i].resind != atoms->atom[i-1].resind)
734 atoms->resinfo[atoms->atom[i].resind].name =
735 put_symtab(symtab, "MOL");
740 /* Residue number is the atom number */
741 atoms->atom[i].resind = i;
742 /* Generate an residue name from the particle type */
743 if (atoms->atom[i].type < 26)
745 sprintf(buf, "T%c", 'A'+atoms->atom[i].type);
749 sprintf(buf, "T%c%c",
750 'A'+atoms->atom[i].type/26, 'A'+atoms->atom[i].type%26);
752 t_atoms_set_resinfo(atoms, i, symtab, buf, i, ' ', 0, ' ');
762 atoms->nres = atoms->nr;
766 gmx_fatal(FARGS, "Internal inconsistency in Espresso routines, read %d atoms, expected %d atoms", i, atoms->nr);
769 else if (level == 1 && strcmp(word, "variable") == 0 && !bFoundVariable)
771 bFoundVariable = TRUE;
772 level += check_open_parenthesis(fp, r, infile, "variable");
773 while (level == 2 && (r = get_espresso_word(fp, word)))
775 if (level == 2 && strcmp(word, "box_l") == 0)
777 for (m = 0; m < 3; m++)
779 r = get_espresso_word(fp, word);
780 sscanf(word, "%lf", &d);
783 level += check_close_parenthesis(fp, r, infile, "box_l");
797 if (!bFoundParticles)
799 fprintf(stderr, "Did not find a particles section in Espresso file '%s'\n",
806 static int get_espresso_coordnum(const char *infile)
810 int natoms, level, r;
811 gmx_bool bFoundParticles;
815 fp = gmx_fio_fopen(infile, "r");
817 bFoundParticles = FALSE;
819 while ((r = get_espresso_word(fp, word)) && !bFoundParticles)
821 if (level == 1 && strcmp(word, "particles") == 0 && !bFoundParticles)
823 bFoundParticles = TRUE;
824 level += check_open_parenthesis(fp, r, infile, "particles");
825 while (level > 0 && (r = get_espresso_word(fp, word)))
850 if (!bFoundParticles)
852 fprintf(stderr, "Did not find a particles section in Espresso file '%s'\n",
861 static void write_espresso_conf_indexed(FILE *out, const char *title,
862 t_atoms *atoms, int nx, atom_id *index,
863 rvec *x, rvec *v, matrix box)
867 fprintf(out, "# %s\n", title);
870 gmx_warning("The Espresso format does not support triclinic unit-cells");
872 fprintf(out, "{variable {box_l %f %f %f}}\n", box[0][0], box[1][1], box[2][2]);
874 fprintf(out, "{particles {id pos type q%s}\n", v ? " v" : "");
875 for (i = 0; i < nx; i++)
885 fprintf(out, "\t{%d %f %f %f %d %g",
886 j, x[j][XX], x[j][YY], x[j][ZZ],
887 atoms->atom[j].type, atoms->atom[j].q);
890 fprintf(out, " %f %f %f", v[j][XX], v[j][YY], v[j][ZZ]);
897 static void get_coordnum_fp (FILE *in, char *title, int *natoms)
901 fgets2 (title, STRLEN, in);
902 fgets2 (line, STRLEN, in);
903 if (sscanf (line, "%d", natoms) != 1)
905 gmx_fatal(FARGS, "gro file does not have the number of atoms on the second line");
909 static void get_coordnum (const char *infile, int *natoms)
914 in = gmx_fio_fopen(infile, "r");
915 get_coordnum_fp(in, title, natoms);
919 static gmx_bool get_w_conf(FILE *in, const char *infile, char *title,
920 t_symtab *symtab, t_atoms *atoms, int *ndec,
921 rvec x[], rvec *v, matrix box)
924 char line[STRLEN+1], *ptr;
926 double x1, y1, z1, x2, y2, z2;
928 int natoms, i, m, resnr, newres, oldres, ddist, c;
929 gmx_bool bFirst, bVel;
933 oldres = NOTSET; /* Unlikely number for the first residue! */
936 /* Read the title and number of atoms */
937 get_coordnum_fp(in, title, &natoms);
939 if (natoms > atoms->nr)
941 gmx_fatal(FARGS, "gro file contains more atoms (%d) than expected (%d)",
944 else if (natoms < atoms->nr)
946 fprintf(stderr, "Warning: gro file contains less atoms (%d) than expected"
947 " (%d)\n", natoms, atoms->nr);
954 /* just pray the arrays are big enough */
955 for (i = 0; (i < natoms); i++)
957 if ((fgets2 (line, STRLEN, in)) == NULL)
959 unexpected_eof(infile, i+2);
961 if (strlen(line) < 39)
963 gmx_fatal(FARGS, "Invalid line in %s for atom %d:\n%s", infile, i+1, line);
966 /* determine read precision from distance between periods
971 p1 = strchr(line, '.');
974 gmx_fatal(FARGS, "A coordinate in file %s does not contain a '.'", infile);
976 p2 = strchr(&p1[1], '.');
979 gmx_fatal(FARGS, "A coordinate in file %s does not contain a '.'", infile);
984 p3 = strchr(&p2[1], '.');
987 gmx_fatal(FARGS, "A coordinate in file %s does not contain a '.'", infile);
990 if (p3 - p2 != ddist)
992 gmx_fatal(FARGS, "The spacing of the decimal points in file %s is not consistent for x, y and z", infile);
997 memcpy(name, line, 5);
999 sscanf(name, "%d", &resnr);
1000 memcpy(name, line+5, 5);
1002 if (resnr != oldres)
1006 if (newres >= natoms)
1008 gmx_fatal(FARGS, "More residues than atoms in %s (natoms = %d)",
1011 atoms->atom[i].resind = newres;
1012 t_atoms_set_resinfo(atoms, i, symtab, name, resnr, ' ', 0, ' ');
1016 atoms->atom[i].resind = newres;
1020 memcpy(name, line+10, 5);
1021 atoms->atomname[i] = put_symtab(symtab, name);
1023 /* eventueel controle atomnumber met i+1 */
1025 /* coordinates (start after residue shit) */
1027 /* Read fixed format */
1028 for (m = 0; m < DIM; m++)
1030 for (c = 0; (c < ddist && ptr[0]); c++)
1036 if (sscanf (buf, "%lf %lf", &x1, &x2) != 1)
1038 gmx_fatal(FARGS, "Something is wrong in the coordinate formatting of file %s. Note that gro is fixed format (see the manual)", infile);
1046 /* velocities (start after residues and coordinates) */
1049 /* Read fixed format */
1050 for (m = 0; m < DIM; m++)
1052 for (c = 0; (c < ddist && ptr[0]); c++)
1058 if (sscanf (buf, "%lf", &x1) != 1)
1070 atoms->nres = newres + 1;
1073 fgets2 (line, STRLEN, in);
1074 if (sscanf (line, "%lf%lf%lf", &x1, &y1, &z1) != 3)
1076 gmx_warning("Bad box in file %s", infile);
1078 /* Generate a cubic box */
1079 for (m = 0; (m < DIM); m++)
1081 xmin[m] = xmax[m] = x[0][m];
1083 for (i = 1; (i < atoms->nr); i++)
1085 for (m = 0; (m < DIM); m++)
1087 xmin[m] = min(xmin[m], x[i][m]);
1088 xmax[m] = max(xmax[m], x[i][m]);
1091 for (i = 0; i < DIM; i++)
1093 for (m = 0; m < DIM; m++)
1098 for (m = 0; (m < DIM); m++)
1100 box[m][m] = (xmax[m]-xmin[m]);
1102 fprintf(stderr, "Generated a cubic box %8.3f x %8.3f x %8.3f\n",
1103 box[XX][XX], box[YY][YY], box[ZZ][ZZ]);
1107 /* We found the first three values, the diagonal elements */
1111 if (sscanf (line, "%*f%*f%*f%lf%lf%lf%lf%lf%lf",
1112 &x1, &y1, &z1, &x2, &y2, &z2) != 6)
1114 x1 = y1 = z1 = x2 = y2 = z2 = 0.0;
1127 static void read_whole_conf(const char *infile, char *title,
1128 t_atoms *atoms, rvec x[], rvec *v, matrix box)
1135 in = gmx_fio_fopen(infile, "r");
1137 open_symtab(&symtab);
1138 get_w_conf(in, infile, title, &symtab, atoms, &ndec, x, v, box);
1139 /* We can't free the symbols, as they are still used in atoms, so
1140 * the only choice is to leak them. */
1141 free_symtab(&symtab);
1146 gmx_bool gro_next_x_or_v(FILE *status, t_trxframe *fr)
1150 char title[STRLEN], *p;
1154 if (gmx_eof(status))
1159 open_symtab(&symtab);
1160 atoms.nr = fr->natoms;
1161 snew(atoms.atom, fr->natoms);
1162 atoms.nres = fr->natoms;
1163 snew(atoms.resinfo, fr->natoms);
1164 snew(atoms.atomname, fr->natoms);
1166 fr->bV = get_w_conf(status, title, title, &symtab, &atoms, &ndec, fr->x, fr->v, fr->box);
1169 /* prec = 10^ndec: */
1170 for (i = 0; i < ndec; i++)
1180 sfree(atoms.resinfo);
1181 sfree(atoms.atomname);
1182 done_symtab(&symtab);
1184 if ((p = strstr(title, "t=")) != NULL)
1187 if (sscanf(p, "%lf", &tt) == 1)
1199 if (atoms.nr != fr->natoms)
1201 gmx_fatal(FARGS, "Number of atoms in gro frame (%d) doesn't match the number in the previous frame (%d)", atoms.nr, fr->natoms);
1207 int gro_first_x_or_v(FILE *status, t_trxframe *fr)
1212 fprintf(stderr, "Reading frames from gro file");
1213 get_coordnum_fp(status, title, &fr->natoms);
1215 fprintf(stderr, " '%s', %d atoms.\n", title, fr->natoms);
1218 if (fr->natoms == 0)
1220 gmx_file("No coordinates in gro file");
1223 snew(fr->x, fr->natoms);
1224 snew(fr->v, fr->natoms);
1225 gro_next_x_or_v(status, fr);
1230 static void make_hconf_format(int pr, gmx_bool bVel, char format[])
1234 /* build format string for printing,
1235 something like "%8.3f" for x and "%8.4f" for v */
1248 sprintf(format, "%%%d.%df%%%d.%df%%%d.%df%%%d.%df%%%d.%df%%%d.%df\n",
1249 l, pr, l, pr, l, pr, l, vpr, l, vpr, l, vpr);
1253 sprintf(format, "%%%d.%df%%%d.%df%%%d.%df\n", l, pr, l, pr, l, pr);
1258 static void write_hconf_box(FILE *out, int pr, matrix box)
1269 if (box[XX][YY] || box[XX][ZZ] || box[YY][XX] || box[YY][ZZ] ||
1270 box[ZZ][XX] || box[ZZ][YY])
1272 sprintf(format, "%%%d.%df%%%d.%df%%%d.%df"
1273 "%%%d.%df%%%d.%df%%%d.%df%%%d.%df%%%d.%df%%%d.%df\n",
1274 l, pr, l, pr, l, pr, l, pr, l, pr, l, pr, l, pr, l, pr, l, pr);
1275 fprintf(out, format,
1276 box[XX][XX], box[YY][YY], box[ZZ][ZZ],
1277 box[XX][YY], box[XX][ZZ], box[YY][XX],
1278 box[YY][ZZ], box[ZZ][XX], box[ZZ][YY]);
1282 sprintf(format, "%%%d.%df%%%d.%df%%%d.%df\n", l, pr, l, pr, l, pr);
1283 fprintf(out, format,
1284 box[XX][XX], box[YY][YY], box[ZZ][ZZ]);
1288 void write_hconf_indexed_p(FILE *out, const char *title, t_atoms *atoms,
1289 int nx, const atom_id index[], int pr,
1290 rvec *x, rvec *v, matrix box)
1292 char resnm[6], nm[6], format[100];
1293 int ai, i, resind, resnr;
1295 bromacs(format, 99);
1296 fprintf (out, "%s\n", (title && title[0]) ? title : format);
1297 fprintf (out, "%5d\n", nx);
1299 make_hconf_format(pr, v != NULL, format);
1301 for (i = 0; (i < nx); i++)
1305 resind = atoms->atom[ai].resind;
1306 strncpy(resnm, " ??? ", sizeof(resnm)-1);
1307 if (resind < atoms->nres)
1309 strncpy(resnm, *atoms->resinfo[resind].name, sizeof(resnm)-1);
1310 resnr = atoms->resinfo[resind].nr;
1314 strncpy(resnm, " ??? ", sizeof(resnm)-1);
1320 strncpy(nm, *atoms->atomname[ai], sizeof(nm)-1);
1324 strncpy(nm, " ??? ", sizeof(nm)-1);
1327 fprintf(out, "%5d%-5.5s%5.5s%5d", resnr%100000, resnm, nm, (ai+1)%100000);
1328 /* next fprintf uses built format string */
1331 fprintf(out, format,
1332 x[ai][XX], x[ai][YY], x[ai][ZZ], v[ai][XX], v[ai][YY], v[ai][ZZ]);
1336 fprintf(out, format,
1337 x[ai][XX], x[ai][YY], x[ai][ZZ]);
1341 write_hconf_box(out, pr, box);
1346 static void write_hconf_mtop(FILE *out, const char *title, gmx_mtop_t *mtop,
1348 rvec *x, rvec *v, matrix box)
1352 gmx_mtop_atomloop_all_t aloop;
1354 char *atomname, *resname;
1356 bromacs(format, 99);
1357 fprintf (out, "%s\n", (title && title[0]) ? title : format);
1358 fprintf (out, "%5d\n", mtop->natoms);
1360 make_hconf_format(pr, v != NULL, format);
1362 aloop = gmx_mtop_atomloop_all_init(mtop);
1363 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
1365 gmx_mtop_atomloop_all_names(aloop, &atomname, &resnr, &resname);
1367 fprintf(out, "%5d%-5.5s%5.5s%5d",
1368 resnr%100000, resname, atomname, (i+1)%100000);
1369 /* next fprintf uses built format string */
1372 fprintf(out, format,
1373 x[i][XX], x[i][YY], x[i][ZZ], v[i][XX], v[i][YY], v[i][ZZ]);
1377 fprintf(out, format,
1378 x[i][XX], x[i][YY], x[i][ZZ]);
1382 write_hconf_box(out, pr, box);
1387 void write_hconf_p(FILE *out, const char *title, t_atoms *atoms, int pr,
1388 rvec *x, rvec *v, matrix box)
1393 snew(aa, atoms->nr);
1394 for (i = 0; (i < atoms->nr); i++)
1398 write_hconf_indexed_p(out, title, atoms, atoms->nr, aa, pr, x, v, box);
1402 void write_conf_p(const char *outfile, const char *title,
1403 t_atoms *atoms, int pr,
1404 rvec *x, rvec *v, matrix box)
1408 out = gmx_fio_fopen(outfile, "w");
1409 write_hconf_p(out, title, atoms, pr, x, v, box);
1411 gmx_fio_fclose (out);
1414 static void write_conf(const char *outfile, const char *title, t_atoms *atoms,
1415 rvec *x, rvec *v, matrix box)
1417 write_conf_p(outfile, title, atoms, 3, x, v, box);
1420 void write_sto_conf_indexed(const char *outfile, const char *title,
1422 rvec x[], rvec *v, int ePBC, matrix box,
1423 atom_id nindex, atom_id index[])
1429 ftp = fn2ftp(outfile);
1433 out = gmx_fio_fopen(outfile, "w");
1434 write_hconf_indexed_p(out, title, atoms, nindex, index, 3, x, v, box);
1435 gmx_fio_fclose(out);
1438 clear_trxframe(&fr, TRUE);
1441 fr.natoms = atoms->nr;
1452 copy_mat(box, fr.box);
1453 out = gmx_fio_fopen(outfile, "w");
1454 write_g96_conf(out, &fr, nindex, index);
1455 gmx_fio_fclose(out);
1461 out = gmx_fio_fopen(outfile, "w");
1462 write_pdbfile_indexed(out, title, atoms, x, ePBC, box, ' ', -1, nindex, index, NULL, TRUE);
1463 gmx_fio_fclose(out);
1466 out = gmx_fio_fopen(outfile, "w");
1467 write_espresso_conf_indexed(out, title, atoms, nindex, index, x, v, box);
1468 gmx_fio_fclose(out);
1473 gmx_fatal(FARGS, "Sorry, can not write a topology to %s", outfile);
1476 gmx_incons("Not supported in write_sto_conf_indexed");
1480 static void write_xyz_conf(const char *outfile, const char *title,
1481 t_atoms *atoms, rvec *x)
1487 gmx_atomprop_t aps = gmx_atomprop_init();
1489 fp = gmx_fio_fopen(outfile, "w");
1490 fprintf(fp, "%3d\n", atoms->nr);
1491 fprintf(fp, "%s\n", title);
1492 for (i = 0; (i < atoms->nr); i++)
1494 anr = atoms->atom[i].atomnumber;
1495 name = *atoms->atomname[i];
1498 if (gmx_atomprop_query(aps, epropElement, "???", name, &value))
1500 anr = gmx_nint(value);
1503 if ((ptr = gmx_atomprop_element(aps, anr)) == NULL)
1507 fprintf(fp, "%3s%10.5f%10.5f%10.5f\n", ptr,
1508 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ]);
1511 gmx_atomprop_destroy(aps);
1514 void write_sto_conf(const char *outfile, const char *title, t_atoms *atoms,
1515 rvec x[], rvec *v, int ePBC, matrix box)
1521 ftp = fn2ftp(outfile);
1525 write_conf(outfile, title, atoms, x, v, box);
1528 clear_trxframe(&fr, TRUE);
1531 fr.natoms = atoms->nr;
1542 copy_mat(box, fr.box);
1543 out = gmx_fio_fopen(outfile, "w");
1544 write_g96_conf(out, &fr, -1, NULL);
1545 gmx_fio_fclose(out);
1548 write_xyz_conf(outfile, (strlen(title) > 0) ? title : outfile, atoms, x);
1553 out = gmx_fio_fopen(outfile, "w");
1554 write_pdbfile(out, title, atoms, x, ePBC, box, ' ', -1, NULL, TRUE);
1555 gmx_fio_fclose(out);
1558 out = gmx_fio_fopen(outfile, "w");
1559 write_espresso_conf_indexed(out, title, atoms, atoms->nr, NULL, x, v, box);
1560 gmx_fio_fclose(out);
1565 gmx_fatal(FARGS, "Sorry, can not write a topology to %s", outfile);
1568 gmx_incons("Not supported in write_sto_conf");
1572 void write_sto_conf_mtop(const char *outfile, const char *title,
1574 rvec x[], rvec *v, int ePBC, matrix box)
1580 ftp = fn2ftp(outfile);
1584 out = gmx_fio_fopen(outfile, "w");
1585 write_hconf_mtop(out, title, mtop, 3, x, v, box);
1586 gmx_fio_fclose(out);
1589 /* This is a brute force approach which requires a lot of memory.
1590 * We should implement mtop versions of all writing routines.
1592 atoms = gmx_mtop_global_atoms(mtop);
1594 write_sto_conf(outfile, title, &atoms, x, v, ePBC, box);
1601 static int get_xyz_coordnum(const char *infile)
1606 fp = gmx_fio_fopen(infile, "r");
1607 if (fscanf(fp, "%d", &n) != 1)
1609 gmx_fatal(FARGS, "Can not read number of atoms from %s", infile);
1616 static void read_xyz_conf(const char *infile, char *title,
1617 t_atoms *atoms, rvec *x)
1623 char atomnm[32], buf[STRLEN];
1626 fp = gmx_fio_fopen(infile, "r");
1627 fgets2(buf, STRLEN-1, fp);
1628 if (sscanf(buf, "%d", &n) != 1)
1630 gmx_fatal(FARGS, "Can not read number of atoms from %s", infile);
1632 fgets2(buf, STRLEN-1, fp);
1634 for (i = 0; (i < n); i++)
1636 fgets2(buf, STRLEN-1, fp);
1637 if (sscanf(buf, "%s%lf%lf%lf", atomnm, &xx, &yy, &zz) != 4)
1639 gmx_fatal(FARGS, "Can not read coordinates from %s", infile);
1641 atoms->atomname[i] = put_symtab(tab, atomnm);
1649 void get_stx_coordnum(const char *infile, int *natoms)
1652 int ftp, tpxver, tpxgen;
1654 char g96_line[STRLEN+1];
1656 ftp = fn2ftp(infile);
1657 range_check(ftp, 0, efNR);
1661 get_coordnum(infile, natoms);
1664 in = gmx_fio_fopen(infile, "r");
1671 *natoms = read_g96_conf(in, infile, &fr, g96_line);
1675 *natoms = get_xyz_coordnum(infile);
1680 in = gmx_fio_fopen(infile, "r");
1681 get_pdb_coordnum(in, natoms);
1685 *natoms = get_espresso_coordnum(infile);
1693 read_tpxheader(infile, &tpx, TRUE, &tpxver, &tpxgen);
1694 *natoms = tpx.natoms;
1698 gmx_fatal(FARGS, "File type %s not supported in get_stx_coordnum",
1703 void read_stx_conf(const char *infile, char *title, t_atoms *atoms,
1704 rvec x[], rvec *v, int *ePBC, matrix box)
1713 char g96_line[STRLEN+1];
1717 fprintf(stderr, "Warning: Number of atoms in %s is 0\n", infile);
1719 else if (atoms->atom == NULL)
1721 gmx_mem("Uninitialized array atom");
1729 ftp = fn2ftp(infile);
1733 read_whole_conf(infile, title, atoms, x, v, box);
1736 read_xyz_conf(infile, title, atoms, x);
1740 fr.natoms = atoms->nr;
1745 in = gmx_fio_fopen(infile, "r");
1746 read_g96_conf(in, infile, &fr, g96_line);
1748 copy_mat(fr.box, box);
1749 strncpy(title, fr.title, STRLEN);
1754 read_pdb_conf(infile, title, atoms, x, ePBC, box, TRUE, NULL);
1757 read_espresso_conf(infile, atoms, x, v, box);
1763 i = read_tpx(infile, NULL, box, &natoms, x, v, NULL, mtop);
1769 strcpy(title, *(mtop->name));
1771 /* Free possibly allocated memory */
1774 *atoms = gmx_mtop_global_atoms(mtop);
1775 top = gmx_mtop_t_to_t_topology(mtop);
1776 tpx_make_chain_identifiers(atoms, &top.mols);
1779 /* The strings in the symtab are still in use in the returned t_atoms
1780 * structure, so we should not free them. But there is no place to put the
1781 * symbols; the only choice is to leak the memory...
1782 * So we clear the symbol table before freeing the topology structure. */
1783 free_symtab(&top.symtab);
1788 gmx_incons("Not supported in read_stx_conf");