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.
49 #ifdef GMX_USE_PLUGINS
60 #include "tngio_for_tools.h"
67 #include "checkpoint.h"
70 #include "gromacs/fileio/timecontrol.h"
72 /* defines for frame counter output */
87 char *persistent_line; /* Persistent line for reading g96 trajectories */
90 /* utility functions */
92 gmx_bool bRmod_fd(double a, double b, double c, gmx_bool bDouble)
97 tol = 2*(bDouble ? GMX_DOUBLE_EPS : GMX_FLOAT_EPS);
99 iq = (int)((a - b + tol*a)/c);
101 if (fabs(a - b - c*iq) <= tol*fabs(a))
113 int check_times2(real t, real t0, gmx_bool bDouble)
118 /* since t is float, we can not use double precision for bRmod */
123 if ((!bTimeSet(TBEGIN) || (t >= rTimeValue(TBEGIN))) &&
124 (!bTimeSet(TEND) || (t <= rTimeValue(TEND))))
126 if (bTimeSet(TDELTA) && !bRmod_fd(t, t0, rTimeValue(TDELTA), bDouble))
135 else if (bTimeSet(TEND) && (t >= rTimeValue(TEND)))
141 fprintf(debug, "t=%g, t0=%g, b=%g, e=%g, dt=%g: r=%d\n",
142 t, t0, rTimeValue(TBEGIN), rTimeValue(TEND), rTimeValue(TDELTA), r);
147 int check_times(real t)
149 return check_times2(t, t, FALSE);
152 static void initcount(t_trxstatus *status)
154 status->__frame = -1;
157 static void status_init(t_trxstatus *status)
160 status->xframe = NULL;
162 status->__frame = -1;
163 status->persistent_line = NULL;
168 int nframes_read(t_trxstatus *status)
170 return status->__frame;
173 static void printcount_(t_trxstatus *status, const output_env_t oenv,
174 const char *l, real t)
176 if ((status->__frame < 2*SKIP1 || status->__frame % SKIP1 == 0) &&
177 (status->__frame < 2*SKIP2 || status->__frame % SKIP2 == 0) &&
178 (status->__frame < 2*SKIP3 || status->__frame % SKIP3 == 0))
180 fprintf(stderr, "\r%-14s %6d time %8.3f ", l, status->__frame,
181 output_env_conv_time(oenv, t));
185 static void printcount(t_trxstatus *status, const output_env_t oenv, real t,
189 printcount_(status, oenv, bSkip ? "Skipping frame" : "Reading frame", t);
192 static void printlast(t_trxstatus *status, const output_env_t oenv, real t)
194 printcount_(status, oenv, "Last frame", t);
195 fprintf(stderr, "\n");
198 static void printincomp(t_trxstatus *status, t_trxframe *fr)
200 if (fr->not_ok & HEADER_NOT_OK)
202 fprintf(stderr, "WARNING: Incomplete header: nr %d time %g\n",
203 status->__frame+1, fr->time);
207 fprintf(stderr, "WARNING: Incomplete frame: nr %d time %g\n",
208 status->__frame+1, fr->time);
212 int prec2ndec(real prec)
216 gmx_fatal(FARGS, "DEATH HORROR prec (%g) <= 0 in prec2ndec", prec);
219 return (int)(log(prec)/log(10.0)+0.5);
222 real ndec2prec(int ndec)
224 return pow(10.0, ndec);
227 t_fileio *trx_get_fileio(t_trxstatus *status)
232 float trx_get_time_of_final_frame(t_trxstatus *status)
234 t_fileio *stfio = trx_get_fileio(status);
235 int filetype = gmx_fio_getftp(stfio);
239 if (filetype == efXTC)
242 xdr_xtc_get_last_frame_time(gmx_fio_getfp(stfio),
243 gmx_fio_getxdr(stfio),
244 status->xframe->natoms, &bOK);
247 gmx_fatal(FARGS, "Error reading last frame. Maybe seek not supported." );
250 else if (filetype == efTNG)
252 tng_trajectory_t tng = status->tng;
255 gmx_fatal(FARGS, "Error opening TNG file.");
257 lasttime = gmx_tng_get_time_of_final_frame(tng);
261 gmx_incons("Only supported for TNG and XTC");
266 void clear_trxframe(t_trxframe *fr, gmx_bool bFirst)
273 fr->bFepState = FALSE;
304 void set_trxframe_ePBC(t_trxframe *fr, int ePBC)
306 fr->bPBC = (ePBC == -1);
310 int write_trxframe_indexed(t_trxstatus *status, t_trxframe *fr, int nind,
311 const atom_id *ind, gmx_conect gc)
314 rvec *xout = NULL, *vout = NULL, *fout = NULL;
331 else if (status->fio)
333 ftp = gmx_fio_getftp(status->fio);
337 gmx_incons("No input file available");
349 gmx_fatal(FARGS, "Need coordinates to write a %s trajectory",
363 for (i = 0; i < nind; i++)
365 copy_rvec(fr->v[ind[i]], vout[i]);
371 for (i = 0; i < nind; i++)
373 copy_rvec(fr->f[ind[i]], fout[i]);
381 for (i = 0; i < nind; i++)
383 copy_rvec(fr->x[ind[i]], xout[i]);
394 gmx_write_tng_from_trxframe(status->tng, fr, nind);
397 write_xtc(status->fio, nind, fr->step, fr->time, fr->box, xout, prec);
401 fwrite_trn(status->fio, nframes_read(status),
402 fr->time, fr->step, fr->box, nind, xout, vout, fout);
410 gmx_fatal(FARGS, "Can not write a %s file without atom names",
413 sprintf(title, "frame t= %.3f", fr->time);
416 write_hconf_indexed_p(gmx_fio_getfp(status->fio), title, fr->atoms, nind, ind,
418 fr->x, fr->bV ? fr->v : NULL, fr->box);
422 write_pdbfile_indexed(gmx_fio_getfp(status->fio), title, fr->atoms,
423 fr->x, -1, fr->box, ' ', fr->step, nind, ind, gc, TRUE);
427 write_g96_conf(gmx_fio_getfp(status->fio), fr, nind, ind);
430 gmx_fatal(FARGS, "Sorry, write_trxframe_indexed can not write %s",
460 void trjtools_gmx_prepare_tng_writing(const char *filename,
466 const gmx_mtop_t *mtop,
467 const atom_id *index,
468 const char *index_group_name)
470 if (filemode != 'w' && filemode != 'a')
472 gmx_incons("Sorry, can only prepare for TNG output.");
483 gmx_prepare_tng_writing(filename,
492 else if (efTNG == fn2ftp(infile))
494 tng_trajectory_t tng_in;
495 gmx_tng_open(infile, 'r', &tng_in);
497 gmx_prepare_tng_writing(filename,
508 void write_tng_frame(t_trxstatus *status,
511 gmx_write_tng_from_trxframe(status->tng, frame, -1);
514 int write_trxframe(t_trxstatus *status, t_trxframe *fr, gmx_conect gc)
530 gmx_tng_set_compression_precision(status->tng, prec);
531 write_tng_frame(status, fr);
536 switch (gmx_fio_getftp(status->fio))
544 gmx_fatal(FARGS, "Need coordinates to write a %s trajectory",
545 ftp2ext(gmx_fio_getftp(status->fio)));
550 switch (gmx_fio_getftp(status->fio))
553 write_xtc(status->fio, fr->natoms, fr->step, fr->time, fr->box, fr->x, prec);
557 fwrite_trn(status->fio, fr->step, fr->time, fr->lambda, fr->box, fr->natoms,
558 fr->bX ? fr->x : NULL, fr->bV ? fr->v : NULL, fr->bF ? fr->f : NULL);
566 gmx_fatal(FARGS, "Can not write a %s file without atom names",
567 ftp2ext(gmx_fio_getftp(status->fio)));
569 sprintf(title, "frame t= %.3f", fr->time);
570 if (gmx_fio_getftp(status->fio) == efGRO)
572 write_hconf_p(gmx_fio_getfp(status->fio), title, fr->atoms,
573 prec2ndec(prec), fr->x, fr->bV ? fr->v : NULL, fr->box);
577 write_pdbfile(gmx_fio_getfp(status->fio), title,
578 fr->atoms, fr->x, fr->bPBC ? fr->ePBC : -1, fr->box,
579 ' ', fr->step, gc, TRUE);
583 write_g96_conf(gmx_fio_getfp(status->fio), fr, -1, NULL);
586 gmx_fatal(FARGS, "Sorry, write_trxframe can not write %s",
587 ftp2ext(gmx_fio_getftp(status->fio)));
594 int write_trx(t_trxstatus *status, int nind, const atom_id *ind, t_atoms *atoms,
595 int step, real time, matrix box, rvec x[], rvec *v,
600 clear_trxframe(&fr, TRUE);
605 fr.bAtoms = atoms != NULL;
612 copy_mat(box, fr.box);
614 return write_trxframe_indexed(status, &fr, nind, ind, gc);
617 void close_trx(t_trxstatus *status)
619 gmx_tng_close(&status->tng);
622 gmx_fio_close(status->fio);
627 t_trxstatus *open_trx(const char *outfile, const char *filemode)
630 if (filemode[0] != 'w' && filemode[0] != 'a' && filemode[1] != '+')
632 gmx_fatal(FARGS, "Sorry, write_trx can only write");
638 stat->fio = gmx_fio_open(outfile, filemode);
642 static gmx_bool gmx_next_frame(t_trxstatus *status, t_trxframe *fr)
649 if (fread_trnheader(status->fio, &sh, &bOK))
651 fr->bDouble = sh.bDouble;
652 fr->natoms = sh.natoms;
658 fr->bFepState = TRUE;
659 fr->lambda = sh.lambda;
660 fr->bBox = sh.box_size > 0;
661 if (fr->flags & (TRX_READ_X | TRX_NEED_X))
665 snew(fr->x, sh.natoms);
667 fr->bX = sh.x_size > 0;
669 if (fr->flags & (TRX_READ_V | TRX_NEED_V))
673 snew(fr->v, sh.natoms);
675 fr->bV = sh.v_size > 0;
677 if (fr->flags & (TRX_READ_F | TRX_NEED_F))
681 snew(fr->f, sh.natoms);
683 fr->bF = sh.f_size > 0;
685 if (fread_htrn(status->fio, &sh, fr->box, fr->x, fr->v, fr->f))
691 fr->not_ok = DATA_NOT_OK;
697 fr->not_ok = HEADER_NOT_OK;
703 static gmx_bool pdb_next_x(t_trxstatus *status, FILE *fp, t_trxframe *fr)
707 int ePBC, model_nr, na;
708 char title[STRLEN], *time;
711 atoms.nr = fr->natoms;
713 atoms.pdbinfo = NULL;
714 /* the other pointers in atoms should not be accessed if these are NULL */
716 na = read_pdbfile(fp, title, &model_nr, &atoms, fr->x, &ePBC, boxpdb, TRUE, NULL);
717 set_trxframe_ePBC(fr, ePBC);
718 if (nframes_read(status) == 0)
720 fprintf(stderr, " '%s', %d atoms\n", title, fr->natoms);
725 fr->bBox = (boxpdb[XX][XX] != 0.0);
728 copy_mat(boxpdb, fr->box);
731 if (model_nr != NOTSET)
736 time = strstr(title, " t= ");
740 sscanf(time+4, "%lf", &dbl);
741 fr->time = (real)dbl;
746 /* this is a bit dirty, but it will work: if no time is read from
747 comment line in pdb file, set time to current frame number */
750 fr->time = (real)fr->step;
754 fr->time = (real)nframes_read(status);
763 if (na != fr->natoms)
765 gmx_fatal(FARGS, "Number of atoms in pdb frame %d is %d instead of %d",
766 nframes_read(status), na, fr->natoms);
772 static int pdb_first_x(t_trxstatus *status, FILE *fp, t_trxframe *fr)
776 fprintf(stderr, "Reading frames from pdb file");
778 get_pdb_coordnum(fp, &fr->natoms);
781 gmx_fatal(FARGS, "\nNo coordinates in pdb file\n");
784 snew(fr->x, fr->natoms);
785 pdb_next_x(status, fp, fr);
790 gmx_bool read_next_frame(const output_env_t oenv, t_trxstatus *status, t_trxframe *fr)
794 gmx_bool bOK, bRet, bMissingData = FALSE, bSkip = FALSE;
803 clear_trxframe(fr, FALSE);
809 /* Special treatment for TNG files */
814 ftp = gmx_fio_getftp(status->fio);
820 bRet = gmx_next_frame(status, fr);
823 /* Checkpoint files can not contain mulitple frames */
826 read_g96_conf(gmx_fio_getfp(status->fio), NULL, fr,
827 status->persistent_line);
828 bRet = (fr->natoms > 0);
832 * Sometimes is off by one frame
833 * and sometimes reports frame not present/file not seekable
835 /* DvdS 2005-05-31: this has been fixed along with the increased
836 * accuracy of the control over -b and -e options.
838 if (bTimeSet(TBEGIN) && (fr->time < rTimeValue(TBEGIN)))
840 if (xtc_seek_time(status->fio, rTimeValue(TBEGIN), fr->natoms, TRUE))
842 gmx_fatal(FARGS, "Specified frame (time %f) doesn't exist or file corrupt/inconsistent.",
847 bRet = read_next_xtc(status->fio, fr->natoms, &fr->step, &fr->time, fr->box,
848 fr->x, &fr->prec, &bOK);
849 fr->bPrec = (bRet && fr->prec > 0);
856 /* Actually the header could also be not ok,
857 but from bOK from read_next_xtc this can't be distinguished */
858 fr->not_ok = DATA_NOT_OK;
862 bRet = gmx_read_next_tng_frame(status->tng, fr, NULL, 0);
865 bRet = pdb_next_x(status, gmx_fio_getfp(status->fio), fr);
868 bRet = gro_next_x_or_v(gmx_fio_getfp(status->fio), fr);
871 #ifdef GMX_USE_PLUGINS
872 bRet = read_next_vmd_frame(fr);
874 gmx_fatal(FARGS, "DEATH HORROR in read_next_frame ftp=%s,status=%s",
875 ftp2ext(gmx_fio_getftp(status->fio)),
876 gmx_fio_getname(status->fio));
882 bMissingData = (((fr->flags & TRX_NEED_X) && !fr->bX) ||
883 ((fr->flags & TRX_NEED_V) && !fr->bV) ||
884 ((fr->flags & TRX_NEED_F) && !fr->bF));
888 ct = check_times2(fr->time, fr->t0, fr->bDouble);
889 if (ct == 0 || ((fr->flags & TRX_DONT_SKIP) && ct < 0))
891 printcount(status, oenv, fr->time, FALSE);
899 printcount(status, oenv, fr->time, TRUE);
906 while (bRet && (bMissingData || bSkip));
910 printlast(status, oenv, pt);
913 printincomp(status, fr);
920 int read_first_frame(const output_env_t oenv, t_trxstatus **status,
921 const char *fn, t_trxframe *fr, int flags)
924 gmx_bool bFirst, bOK;
926 int ftp = fn2ftp(fn);
927 gmx_int64_t *tng_ids;
929 clear_trxframe(fr, TRUE);
936 status_init( *status );
937 (*status)->nxframe = 1;
942 /* Special treatment for TNG files */
943 gmx_tng_open(fn, 'r', &(*status)->tng);
947 fio = (*status)->fio = gmx_fio_open(fn, "r");
955 read_checkpoint_trxframe(fio, fr);
959 /* Can not rewind a compressed file, so open it twice */
960 if (!(*status)->persistent_line)
962 /* allocate the persistent line */
963 snew((*status)->persistent_line, STRLEN+1);
965 read_g96_conf(gmx_fio_getfp(fio), fn, fr, (*status)->persistent_line);
967 clear_trxframe(fr, FALSE);
968 if (flags & (TRX_READ_X | TRX_NEED_X))
970 snew(fr->x, fr->natoms);
972 if (flags & (TRX_READ_V | TRX_NEED_V))
974 snew(fr->v, fr->natoms);
976 fio = (*status)->fio = gmx_fio_open(fn, "r");
979 if (read_first_xtc(fio, &fr->natoms, &fr->step, &fr->time, fr->box, &fr->x,
980 &fr->prec, &bOK) == 0)
983 fr->not_ok = DATA_NOT_OK;
988 printincomp(*status, fr);
992 fr->bPrec = (fr->prec > 0);
997 printcount(*status, oenv, fr->time, FALSE);
1003 if (!gmx_read_next_tng_frame((*status)->tng, fr, NULL, 0))
1005 fr->not_ok = DATA_NOT_OK;
1007 printincomp(*status, fr);
1011 printcount(*status, oenv, fr->time, FALSE);
1016 pdb_first_x(*status, gmx_fio_getfp(fio), fr);
1019 printcount(*status, oenv, fr->time, FALSE);
1024 if (gro_first_x_or_v(gmx_fio_getfp(fio), fr))
1026 printcount(*status, oenv, fr->time, FALSE);
1031 #ifdef GMX_USE_PLUGINS
1032 fprintf(stderr, "The file format of %s is not a known trajectory format to GROMACS.\n"
1033 "Please make sure that the file is a trajectory!\n"
1034 "GROMACS will now assume it to be a trajectory and will try to open it using the VMD plug-ins.\n"
1035 "This will only work in case the VMD plugins are found and it is a trajectory format supported by VMD.\n", fn);
1036 gmx_fio_fp_close(fio); /*only close the file without removing FIO entry*/
1037 if (!read_first_vmd_frame(fn, fr))
1039 gmx_fatal(FARGS, "Not supported in read_first_frame: %s", fn);
1042 gmx_fatal(FARGS, "Not supported in read_first_frame: %s. Please make sure that the file is a trajectory.\n"
1043 "GROMACS is not compiled with plug-in support. Thus it cannot read non-GROMACS trajectory formats using the VMD plug-ins.\n"
1044 "Please compile with plug-in support if you want to read non-GROMACS trajectory formats.\n", fn);
1049 /* Return FALSE if we read a frame that's past the set ending time. */
1050 if (!bFirst && (!(fr->flags & TRX_DONT_SKIP) && check_times(fr->time) > 0))
1057 (!(fr->flags & TRX_DONT_SKIP) && check_times(fr->time) < 0))
1059 /* Read a frame when no frame was read or the first was skipped */
1060 if (!read_next_frame(oenv, *status, fr))
1067 return (fr->natoms > 0);
1070 /***** C O O R D I N A T E S T U F F *****/
1072 int read_first_x(const output_env_t oenv, t_trxstatus **status, const char *fn,
1073 real *t, rvec **x, matrix box)
1077 read_first_frame(oenv, status, fn, &fr, TRX_NEED_X);
1079 snew((*status)->xframe, 1);
1080 (*status)->nxframe = 1;
1081 (*(*status)->xframe) = fr;
1082 *t = (*status)->xframe->time;
1083 *x = (*status)->xframe->x;
1084 copy_mat((*status)->xframe->box, box);
1086 return (*status)->xframe->natoms;
1089 gmx_bool read_next_x(const output_env_t oenv, t_trxstatus *status, real *t,
1090 rvec x[], matrix box)
1094 status->xframe->x = x;
1095 /*xframe[status].x = x;*/
1096 bRet = read_next_frame(oenv, status, status->xframe);
1097 *t = status->xframe->time;
1098 copy_mat(status->xframe->box, box);
1103 void close_trj(t_trxstatus *status)
1105 gmx_tng_close(&status->tng);
1108 gmx_fio_close(status->fio);
1111 /* The memory in status->xframe is lost here,
1112 * but the read_first_x/read_next_x functions are deprecated anyhow.
1113 * read_first_frame/read_next_frame and close_trx should be used.
1118 void rewind_trj(t_trxstatus *status)
1122 gmx_fio_rewind(status->fio);
1125 /***** T O P O L O G Y S T U F F ******/
1127 t_topology *read_top(const char *fn, int *ePBC)
1133 epbc = read_tpx_top(fn, NULL, NULL, &natoms, NULL, NULL, NULL, top);