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, 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.
54 #include "gromacs/legacyheaders/statutil.h"
62 #include "checkpoint.h"
66 /* defines for frame counter output */
71 /* Globals for gromos-87 input */
73 effXYZ, effXYZBox, effG87, effG87Box, effNR
86 char *persistent_line; /* Persistent line for reading g96 trajectories */
89 /* utility functions */
91 int check_times2(real t, real t0, gmx_bool bDouble)
96 /* since t is float, we can not use double precision for bRmod */
101 if ((!bTimeSet(TBEGIN) || (t >= rTimeValue(TBEGIN))) &&
102 (!bTimeSet(TEND) || (t <= rTimeValue(TEND))))
104 if (bTimeSet(TDELTA) && !bRmod_fd(t, t0, rTimeValue(TDELTA), bDouble))
113 else if (bTimeSet(TEND) && (t >= rTimeValue(TEND)))
119 fprintf(debug, "t=%g, t0=%g, b=%g, e=%g, dt=%g: r=%d\n",
120 t, t0, rTimeValue(TBEGIN), rTimeValue(TEND), rTimeValue(TDELTA), r);
125 int check_times(real t)
127 return check_times2(t, t, FALSE);
130 static void initcount(t_trxstatus *status)
132 status->__frame = -1;
135 static void status_init(t_trxstatus *status)
138 status->xframe = NULL;
140 status->__frame = -1;
141 status->persistent_line = NULL;
145 int nframes_read(t_trxstatus *status)
147 return status->__frame;
150 static void printcount_(t_trxstatus *status, const output_env_t oenv,
151 const char *l, real t)
153 if ((status->__frame < 2*SKIP1 || status->__frame % SKIP1 == 0) &&
154 (status->__frame < 2*SKIP2 || status->__frame % SKIP2 == 0) &&
155 (status->__frame < 2*SKIP3 || status->__frame % SKIP3 == 0))
157 fprintf(stderr, "\r%-14s %6d time %8.3f ", l, status->__frame,
158 output_env_conv_time(oenv, t));
162 static void printcount(t_trxstatus *status, const output_env_t oenv, real t,
166 printcount_(status, oenv, bSkip ? "Skipping frame" : "Reading frame", t);
169 static void printlast(t_trxstatus *status, const output_env_t oenv, real t)
171 printcount_(status, oenv, "Last frame", t);
172 fprintf(stderr, "\n");
175 static void printincomp(t_trxstatus *status, t_trxframe *fr)
177 if (fr->not_ok & HEADER_NOT_OK)
179 fprintf(stderr, "WARNING: Incomplete header: nr %d time %g\n",
180 status->__frame+1, fr->time);
184 fprintf(stderr, "WARNING: Incomplete frame: nr %d time %g\n",
185 status->__frame+1, fr->time);
189 int prec2ndec(real prec)
193 gmx_fatal(FARGS, "DEATH HORROR prec (%g) <= 0 in prec2ndec", prec);
196 return (int)(log(prec)/log(10.0)+0.5);
200 t_fileio *trx_get_fileio(t_trxstatus *status)
207 void clear_trxframe(t_trxframe *fr, gmx_bool bFirst)
214 fr->bFepState = FALSE;
245 void set_trxframe_ePBC(t_trxframe *fr, int ePBC)
247 fr->bPBC = (ePBC == -1);
251 int write_trxframe_indexed(t_trxstatus *status, t_trxframe *fr, int nind,
252 const atom_id *ind, gmx_conect gc)
255 rvec *xout = NULL, *vout = NULL, *fout = NULL;
268 switch (gmx_fio_getftp(status->fio))
276 gmx_fatal(FARGS, "Need coordinates to write a %s trajectory",
277 ftp2ext(gmx_fio_getftp(status->fio)));
282 switch (gmx_fio_getftp(status->fio))
289 for (i = 0; i < nind; i++)
291 copy_rvec(fr->v[ind[i]], vout[i]);
297 for (i = 0; i < nind; i++)
299 copy_rvec(fr->f[ind[i]], fout[i]);
308 for (i = 0; i < nind; i++)
310 copy_rvec(fr->x[ind[i]], xout[i]);
318 switch (gmx_fio_getftp(status->fio))
321 write_xtc(status->fio, nind, fr->step, fr->time, fr->box, xout, prec);
325 fwrite_trn(status->fio, nframes_read(status),
326 fr->time, fr->step, fr->box, nind, xout, vout, fout);
334 gmx_fatal(FARGS, "Can not write a %s file without atom names",
335 ftp2ext(gmx_fio_getftp(status->fio)));
337 sprintf(title, "frame t= %.3f", fr->time);
338 if (gmx_fio_getftp(status->fio) == efGRO)
340 write_hconf_indexed_p(gmx_fio_getfp(status->fio), title, fr->atoms, nind, ind,
342 fr->x, fr->bV ? fr->v : NULL, fr->box);
346 write_pdbfile_indexed(gmx_fio_getfp(status->fio), title, fr->atoms,
347 fr->x, -1, fr->box, ' ', fr->step, nind, ind, gc, TRUE);
351 write_gms(gmx_fio_getfp(status->fio), nind, xout, fr->box);
354 write_g96_conf(gmx_fio_getfp(status->fio), fr, nind, ind);
357 gmx_fatal(FARGS, "Sorry, write_trxframe_indexed can not write %s",
358 ftp2ext(gmx_fio_getftp(status->fio)));
362 switch (gmx_fio_getftp(status->fio))
387 int write_trxframe(t_trxstatus *status, t_trxframe *fr, gmx_conect gc)
401 switch (gmx_fio_getftp(status->fio))
409 gmx_fatal(FARGS, "Need coordinates to write a %s trajectory",
410 ftp2ext(gmx_fio_getftp(status->fio)));
415 switch (gmx_fio_getftp(status->fio))
418 write_xtc(status->fio, fr->natoms, fr->step, fr->time, fr->box, fr->x, prec);
422 fwrite_trn(status->fio, fr->step, fr->time, fr->lambda, fr->box, fr->natoms,
423 fr->bX ? fr->x : NULL, fr->bV ? fr->v : NULL, fr->bF ? fr->f : NULL);
431 gmx_fatal(FARGS, "Can not write a %s file without atom names",
432 ftp2ext(gmx_fio_getftp(status->fio)));
434 sprintf(title, "frame t= %.3f", fr->time);
435 if (gmx_fio_getftp(status->fio) == efGRO)
437 write_hconf_p(gmx_fio_getfp(status->fio), title, fr->atoms,
438 prec2ndec(prec), fr->x, fr->bV ? fr->v : NULL, fr->box);
442 write_pdbfile(gmx_fio_getfp(status->fio), title,
443 fr->atoms, fr->x, fr->bPBC ? fr->ePBC : -1, fr->box,
444 ' ', fr->step, gc, TRUE);
448 write_gms(gmx_fio_getfp(status->fio), fr->natoms, fr->x, fr->box);
451 write_g96_conf(gmx_fio_getfp(status->fio), fr, -1, NULL);
454 gmx_fatal(FARGS, "Sorry, write_trxframe can not write %s",
455 ftp2ext(gmx_fio_getftp(status->fio)));
462 int write_trx(t_trxstatus *status, int nind, const atom_id *ind, t_atoms *atoms,
463 int step, real time, matrix box, rvec x[], rvec *v,
468 clear_trxframe(&fr, TRUE);
473 fr.bAtoms = atoms != NULL;
480 copy_mat(box, fr.box);
482 return write_trxframe_indexed(status, &fr, nind, ind, gc);
485 void close_trx(t_trxstatus *status)
487 gmx_fio_close(status->fio);
491 t_trxstatus *open_trx(const char *outfile, const char *filemode)
494 if (filemode[0] != 'w' && filemode[0] != 'a' && filemode[1] != '+')
496 gmx_fatal(FARGS, "Sorry, write_trx can only write");
502 stat->fio = gmx_fio_open(outfile, filemode);
506 static gmx_bool gmx_next_frame(t_trxstatus *status, t_trxframe *fr)
513 if (fread_trnheader(status->fio, &sh, &bOK))
515 fr->bDouble = sh.bDouble;
516 fr->natoms = sh.natoms;
522 fr->bFepState = TRUE;
523 fr->lambda = sh.lambda;
524 fr->bBox = sh.box_size > 0;
525 if (fr->flags & (TRX_READ_X | TRX_NEED_X))
529 snew(fr->x, sh.natoms);
531 fr->bX = sh.x_size > 0;
533 if (fr->flags & (TRX_READ_V | TRX_NEED_V))
537 snew(fr->v, sh.natoms);
539 fr->bV = sh.v_size > 0;
541 if (fr->flags & (TRX_READ_F | TRX_NEED_F))
545 snew(fr->f, sh.natoms);
547 fr->bF = sh.f_size > 0;
549 if (fread_htrn(status->fio, &sh, fr->box, fr->x, fr->v, fr->f))
555 fr->not_ok = DATA_NOT_OK;
561 fr->not_ok = HEADER_NOT_OK;
567 static void choose_file_format(FILE *fp)
575 printf(" Select File Format\n");
576 printf("---------------------------\n");
577 printf("1. XYZ File\n");
578 printf("2. XYZ File with Box\n");
579 printf("3. Gromos-87 Ascii Trajectory\n");
580 printf("4. Gromos-87 Ascii Trajectory with Box\n");
587 printf("\nChoice: ");
591 rc = scanf("%d", &i);
596 while ((i < 0) || (i >= effNR));
599 stat->eFF = (eFileFormat) i;
601 for (m = 0; (m < DIM); m++)
606 stat->bReadBox = (stat->eFF == effG87Box) || (stat->eFF == effXYZBox);
612 if (5 != fscanf(fp, "%d%lf%lf%lf%lf", &stat->NATOMS, &stat->BOX[XX], &stat->BOX[YY], &stat->BOX[ZZ], &stat->DT))
614 gmx_fatal(FARGS, "Error reading natoms/box in file");
619 printf("GROMOS! OH DEAR...\n\n");
620 printf("Number of atoms ? ");
622 if (1 != scanf("%d", &stat->NATOMS))
624 gmx_fatal(FARGS, "Error reading natoms in file");
627 printf("Time between timeframes ? ");
629 if (1 != scanf("%lf", &stat->DT))
631 gmx_fatal(FARGS, "Error reading dt from file");
634 if (stat->eFF == effG87)
636 printf("Box X Y Z ? ");
638 if (3 != scanf("%lf%lf%lf", &stat->BOX[XX], &stat->BOX[YY], &stat->BOX[ZZ]))
640 gmx_fatal(FARGS, "Error reading box in file");
653 printf("Hellow World\n");
657 static gmx_bool do_read_xyz(t_trxstatus *status, FILE *fp, int natoms,
658 rvec x[], matrix box)
663 for (i = 0; (i < natoms); i++)
665 for (m = 0; (m < DIM); m++)
667 if (fscanf(fp, "%lf", &x0) != 1)
671 fprintf(stderr, "error reading statusfile: x[%d][%d]\n", i, m);
679 if (status->bReadBox)
681 for (m = 0; (m < DIM); m++)
683 if (fscanf(fp, "%lf", &x0) != 1)
693 static gmx_bool xyz_next_x(t_trxstatus *status, FILE *fp, const output_env_t oenv,
694 real *t, int natoms, rvec x[], matrix box)
695 /* Reads until a new x can be found (return TRUE)
696 * or eof (return FALSE)
702 while (!bTimeSet(TBEGIN) || (*t < rTimeValue(TBEGIN)))
704 if (!do_read_xyz(status, fp, natoms, x, box))
708 printcount(status, oenv, *t, FALSE);
712 if (!bTimeSet(TEND) || (*t <= rTimeValue(TEND)))
714 if (!do_read_xyz(status, fp, natoms, x, box))
716 printlast(status, oenv, *t);
719 printcount(status, oenv, *t, FALSE);
724 printlast(status, oenv, pt);
728 static int xyz_first_x(t_trxstatus *status, FILE *fp, const output_env_t oenv,
729 real *t, rvec **x, matrix box)
730 /* Reads fp, mallocs x, and returns x and box
731 * Returns natoms when successful, FALSE otherwise
739 choose_file_format(fp);
741 for (m = 0; (m < DIM); m++)
743 box[m][m] = status->BOX[m];
746 snew(*x, status->NATOMS);
748 if (!xyz_next_x(status, fp, oenv, t, status->NATOMS, *x, box))
754 return status->NATOMS;
757 static gmx_bool pdb_next_x(t_trxstatus *status, FILE *fp, t_trxframe *fr)
761 int ePBC, model_nr, na;
762 char title[STRLEN], *time;
765 atoms.nr = fr->natoms;
767 atoms.pdbinfo = NULL;
768 /* the other pointers in atoms should not be accessed if these are NULL */
770 na = read_pdbfile(fp, title, &model_nr, &atoms, fr->x, &ePBC, boxpdb, TRUE, NULL);
771 set_trxframe_ePBC(fr, ePBC);
772 if (nframes_read(status) == 0)
774 fprintf(stderr, " '%s', %d atoms\n", title, fr->natoms);
779 fr->bBox = (boxpdb[XX][XX] != 0.0);
782 copy_mat(boxpdb, fr->box);
785 if (model_nr != NOTSET)
790 time = strstr(title, " t= ");
794 sscanf(time+4, "%lf", &dbl);
795 fr->time = (real)dbl;
800 /* this is a bit dirty, but it will work: if no time is read from
801 comment line in pdb file, set time to current frame number */
804 fr->time = (real)fr->step;
808 fr->time = (real)nframes_read(status);
817 if (na != fr->natoms)
819 gmx_fatal(FARGS, "Number of atoms in pdb frame %d is %d instead of %d",
820 nframes_read(status), na, fr->natoms);
826 static int pdb_first_x(t_trxstatus *status, FILE *fp, t_trxframe *fr)
830 fprintf(stderr, "Reading frames from pdb file");
832 get_pdb_coordnum(fp, &fr->natoms);
835 gmx_fatal(FARGS, "\nNo coordinates in pdb file\n");
838 snew(fr->x, fr->natoms);
839 pdb_next_x(status, fp, fr);
844 gmx_bool read_next_frame(const output_env_t oenv, t_trxstatus *status, t_trxframe *fr)
848 gmx_bool bOK, bRet, bMissingData = FALSE, bSkip = FALSE;
856 clear_trxframe(fr, FALSE);
860 switch (gmx_fio_getftp(status->fio))
864 bRet = gmx_next_frame(status, fr);
867 /* Checkpoint files can not contain mulitple frames */
870 read_g96_conf(gmx_fio_getfp(status->fio), NULL, fr,
871 status->persistent_line);
872 bRet = (fr->natoms > 0);
875 bRet = xyz_next_x(status, gmx_fio_getfp(status->fio), oenv, &fr->time,
876 fr->natoms, fr->x, fr->box);
883 * Sometimes is off by one frame
884 * and sometimes reports frame not present/file not seekable
886 /* DvdS 2005-05-31: this has been fixed along with the increased
887 * accuracy of the control over -b and -e options.
889 if (bTimeSet(TBEGIN) && (fr->time < rTimeValue(TBEGIN)))
891 if (xtc_seek_time(status->fio, rTimeValue(TBEGIN), fr->natoms, TRUE))
893 gmx_fatal(FARGS, "Specified frame (time %f) doesn't exist or file corrupt/inconsistent.",
898 bRet = read_next_xtc(status->fio, fr->natoms, &fr->step, &fr->time, fr->box,
899 fr->x, &fr->prec, &bOK);
900 fr->bPrec = (bRet && fr->prec > 0);
907 /* Actually the header could also be not ok,
908 but from bOK from read_next_xtc this can't be distinguished */
909 fr->not_ok = DATA_NOT_OK;
913 bRet = pdb_next_x(status, gmx_fio_getfp(status->fio), fr);
916 bRet = gro_next_x_or_v(gmx_fio_getfp(status->fio), fr);
919 #ifdef GMX_USE_PLUGINS
920 bRet = read_next_vmd_frame(fr);
922 gmx_fatal(FARGS, "DEATH HORROR in read_next_frame ftp=%s,status=%s",
923 ftp2ext(gmx_fio_getftp(status->fio)),
924 gmx_fio_getname(status->fio));
930 bMissingData = (((fr->flags & TRX_NEED_X) && !fr->bX) ||
931 ((fr->flags & TRX_NEED_V) && !fr->bV) ||
932 ((fr->flags & TRX_NEED_F) && !fr->bF));
936 ct = check_times2(fr->time, fr->t0, fr->bDouble);
937 if (ct == 0 || ((fr->flags & TRX_DONT_SKIP) && ct < 0))
939 printcount(status, oenv, fr->time, FALSE);
947 printcount(status, oenv, fr->time, TRUE);
954 while (bRet && (bMissingData || bSkip));
958 printlast(status, oenv, pt);
961 printincomp(status, fr);
968 int read_first_frame(const output_env_t oenv, t_trxstatus **status,
969 const char *fn, t_trxframe *fr, int flags)
972 gmx_bool bFirst, bOK;
975 clear_trxframe(fr, TRUE);
982 status_init( *status );
983 (*status)->nxframe = 1;
986 fio = (*status)->fio = gmx_fio_open(fn, "r");
987 switch (gmx_fio_getftp(fio))
993 read_checkpoint_trxframe(fio, fr);
997 /* Can not rewind a compressed file, so open it twice */
998 if (!(*status)->persistent_line)
1000 /* allocate the persistent line */
1001 snew((*status)->persistent_line, STRLEN+1);
1003 read_g96_conf(gmx_fio_getfp(fio), fn, fr, (*status)->persistent_line);
1005 clear_trxframe(fr, FALSE);
1006 if (flags & (TRX_READ_X | TRX_NEED_X))
1008 snew(fr->x, fr->natoms);
1010 if (flags & (TRX_READ_V | TRX_NEED_V))
1012 snew(fr->v, fr->natoms);
1014 fio = (*status)->fio = gmx_fio_open(fn, "r");
1017 fr->natoms = xyz_first_x(*status, gmx_fio_getfp(fio), oenv, &fr->time,
1024 printcount(*status, oenv, fr->time, FALSE);
1029 if (read_first_xtc(fio, &fr->natoms, &fr->step, &fr->time, fr->box, &fr->x,
1030 &fr->prec, &bOK) == 0)
1034 gmx_fatal(FARGS, "No XTC!\n");
1038 fr->not_ok = DATA_NOT_OK;
1044 printincomp(*status, fr);
1048 fr->bPrec = (fr->prec > 0);
1053 printcount(*status, oenv, fr->time, FALSE);
1058 pdb_first_x(*status, gmx_fio_getfp(fio), fr);
1061 printcount(*status, oenv, fr->time, FALSE);
1066 if (gro_first_x_or_v(gmx_fio_getfp(fio), fr))
1068 printcount(*status, oenv, fr->time, FALSE);
1073 #ifdef GMX_USE_PLUGINS
1074 fprintf(stderr, "The file format of %s is not a known trajectory format to GROMACS.\n"
1075 "Please make sure that the file is a trajectory!\n"
1076 "GROMACS will now assume it to be a trajectory and will try to open it using the VMD plug-ins.\n"
1077 "This will only work in case the VMD plugins are found and it is a trajectory format supported by VMD.\n", fn);
1078 gmx_fio_fp_close(fio); /*only close the file without removing FIO entry*/
1079 if (!read_first_vmd_frame(fn, fr))
1081 gmx_fatal(FARGS, "Not supported in read_first_frame: %s", fn);
1084 gmx_fatal(FARGS, "Not supported in read_first_frame: %s. Please make sure that the file is a trajectory.\n"
1085 "GROMACS is not compiled with plug-in support. Thus it cannot read non-GROMACS trajectory formats using the VMD plug-ins.\n"
1086 "Please compile with plug-in support if you want to read non-GROMACS trajectory formats.\n", fn);
1091 /* Return FALSE if we read a frame that's past the set ending time. */
1092 if (!bFirst && (!(fr->flags & TRX_DONT_SKIP) && check_times(fr->time) > 0))
1099 (!(fr->flags & TRX_DONT_SKIP) && check_times(fr->time) < 0))
1101 /* Read a frame when no frame was read or the first was skipped */
1102 if (!read_next_frame(oenv, *status, fr))
1109 return (fr->natoms > 0);
1112 /***** C O O R D I N A T E S T U F F *****/
1114 int read_first_x(const output_env_t oenv, t_trxstatus **status, const char *fn,
1115 real *t, rvec **x, matrix box)
1119 read_first_frame(oenv, status, fn, &fr, TRX_NEED_X);
1121 snew((*status)->xframe, 1);
1122 (*status)->nxframe = 1;
1123 (*(*status)->xframe) = fr;
1124 *t = (*status)->xframe->time;
1125 *x = (*status)->xframe->x;
1126 copy_mat((*status)->xframe->box, box);
1128 return (*status)->xframe->natoms;
1131 gmx_bool read_next_x(const output_env_t oenv, t_trxstatus *status, real *t,
1132 rvec x[], matrix box)
1136 status->xframe->x = x;
1137 /*xframe[status].x = x;*/
1138 bRet = read_next_frame(oenv, status, status->xframe);
1139 *t = status->xframe->time;
1140 copy_mat(status->xframe->box, box);
1145 void close_trj(t_trxstatus *status)
1147 gmx_fio_close(status->fio);
1148 /* The memory in status->xframe is lost here,
1149 * but the read_first_x/read_next_x functions are deprecated anyhow.
1150 * read_first_frame/read_next_frame and close_trx should be used.
1155 void rewind_trj(t_trxstatus *status)
1159 gmx_fio_rewind(status->fio);