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.
46 #ifdef GMX_USE_PLUGINS
54 #include "tngio_for_tools.h"
55 #include "gromacs/legacyheaders/names.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/utility/futil.h"
61 #include "gromacs/legacyheaders/checkpoint.h"
64 #include "gromacs/fileio/timecontrol.h"
65 #include "gromacs/fileio/trx.h"
66 #include "gromacs/topology/atoms.h"
67 #include "gromacs/utility/fatalerror.h"
68 #include "gromacs/utility/smalloc.h"
70 /* defines for frame counter output */
85 char *persistent_line; /* Persistent line for reading g96 trajectories */
88 /* utility functions */
90 gmx_bool bRmod_fd(double a, double b, double c, gmx_bool bDouble)
95 tol = 2*(bDouble ? GMX_DOUBLE_EPS : GMX_FLOAT_EPS);
97 iq = (int)((a - b + tol*a)/c);
99 if (fabs(a - b - c*iq) <= tol*fabs(a))
111 int check_times2(real t, real t0, gmx_bool bDouble)
116 /* since t is float, we can not use double precision for bRmod */
121 if ((!bTimeSet(TBEGIN) || (t >= rTimeValue(TBEGIN))) &&
122 (!bTimeSet(TEND) || (t <= rTimeValue(TEND))))
124 if (bTimeSet(TDELTA) && !bRmod_fd(t, t0, rTimeValue(TDELTA), bDouble))
133 else if (bTimeSet(TEND) && (t >= rTimeValue(TEND)))
139 fprintf(debug, "t=%g, t0=%g, b=%g, e=%g, dt=%g: r=%d\n",
140 t, t0, rTimeValue(TBEGIN), rTimeValue(TEND), rTimeValue(TDELTA), r);
145 int check_times(real t)
147 return check_times2(t, t, FALSE);
150 static void initcount(t_trxstatus *status)
152 status->__frame = -1;
155 static void status_init(t_trxstatus *status)
158 status->xframe = NULL;
160 status->__frame = -1;
161 status->persistent_line = NULL;
166 int nframes_read(t_trxstatus *status)
168 return status->__frame;
171 static void printcount_(t_trxstatus *status, const output_env_t oenv,
172 const char *l, real t)
174 if ((status->__frame < 2*SKIP1 || status->__frame % SKIP1 == 0) &&
175 (status->__frame < 2*SKIP2 || status->__frame % SKIP2 == 0) &&
176 (status->__frame < 2*SKIP3 || status->__frame % SKIP3 == 0))
178 fprintf(stderr, "\r%-14s %6d time %8.3f ", l, status->__frame,
179 output_env_conv_time(oenv, t));
183 static void printcount(t_trxstatus *status, const output_env_t oenv, real t,
187 printcount_(status, oenv, bSkip ? "Skipping frame" : "Reading frame", t);
190 static void printlast(t_trxstatus *status, const output_env_t oenv, real t)
192 printcount_(status, oenv, "Last frame", t);
193 fprintf(stderr, "\n");
196 static void printincomp(t_trxstatus *status, t_trxframe *fr)
198 if (fr->not_ok & HEADER_NOT_OK)
200 fprintf(stderr, "WARNING: Incomplete header: nr %d time %g\n",
201 status->__frame+1, fr->time);
205 fprintf(stderr, "WARNING: Incomplete frame: nr %d time %g\n",
206 status->__frame+1, fr->time);
210 int prec2ndec(real prec)
214 gmx_fatal(FARGS, "DEATH HORROR prec (%g) <= 0 in prec2ndec", prec);
217 return (int)(log(prec)/log(10.0)+0.5);
220 real ndec2prec(int ndec)
222 return pow(10.0, ndec);
225 t_fileio *trx_get_fileio(t_trxstatus *status)
230 float trx_get_time_of_final_frame(t_trxstatus *status)
232 t_fileio *stfio = trx_get_fileio(status);
233 int filetype = gmx_fio_getftp(stfio);
237 if (filetype == efXTC)
240 xdr_xtc_get_last_frame_time(gmx_fio_getfp(stfio),
241 gmx_fio_getxdr(stfio),
242 status->xframe->natoms, &bOK);
245 gmx_fatal(FARGS, "Error reading last frame. Maybe seek not supported." );
248 else if (filetype == efTNG)
250 tng_trajectory_t tng = status->tng;
253 gmx_fatal(FARGS, "Error opening TNG file.");
255 lasttime = gmx_tng_get_time_of_final_frame(tng);
259 gmx_incons("Only supported for TNG and XTC");
264 void clear_trxframe(t_trxframe *fr, gmx_bool bFirst)
271 fr->bFepState = FALSE;
303 void set_trxframe_ePBC(t_trxframe *fr, int ePBC)
305 fr->bPBC = (ePBC == -1);
309 int write_trxframe_indexed(t_trxstatus *status, t_trxframe *fr, int nind,
310 const atom_id *ind, gmx_conect gc)
313 rvec *xout = NULL, *vout = NULL, *fout = NULL;
330 else if (status->fio)
332 ftp = gmx_fio_getftp(status->fio);
336 gmx_incons("No input file available");
347 gmx_fatal(FARGS, "Need coordinates to write a %s trajectory",
360 for (i = 0; i < nind; i++)
362 copy_rvec(fr->v[ind[i]], vout[i]);
368 for (i = 0; i < nind; i++)
370 copy_rvec(fr->f[ind[i]], fout[i]);
378 for (i = 0; i < nind; i++)
380 copy_rvec(fr->x[ind[i]], xout[i]);
391 gmx_write_tng_from_trxframe(status->tng, fr, nind);
394 write_xtc(status->fio, nind, fr->step, fr->time, fr->box, xout, prec);
397 fwrite_trn(status->fio, nframes_read(status),
398 fr->time, fr->step, fr->box, nind, xout, vout, fout);
406 gmx_fatal(FARGS, "Can not write a %s file without atom names",
409 sprintf(title, "frame t= %.3f", fr->time);
412 write_hconf_indexed_p(gmx_fio_getfp(status->fio), title, fr->atoms, nind, ind,
414 fr->x, fr->bV ? fr->v : NULL, fr->box);
418 write_pdbfile_indexed(gmx_fio_getfp(status->fio), title, fr->atoms,
419 fr->x, -1, fr->box, ' ', fr->step, nind, ind, gc, TRUE);
423 write_g96_conf(gmx_fio_getfp(status->fio), fr, nind, ind);
426 gmx_fatal(FARGS, "Sorry, write_trxframe_indexed can not write %s",
454 void trjtools_gmx_prepare_tng_writing(const char *filename,
460 const gmx_mtop_t *mtop,
461 const atom_id *index,
462 const char *index_group_name)
464 if (filemode != 'w' && filemode != 'a')
466 gmx_incons("Sorry, can only prepare for TNG output.");
477 gmx_prepare_tng_writing(filename,
486 else if (efTNG == fn2ftp(infile))
488 tng_trajectory_t tng_in;
489 gmx_tng_open(infile, 'r', &tng_in);
491 gmx_prepare_tng_writing(filename,
502 void write_tng_frame(t_trxstatus *status,
505 gmx_write_tng_from_trxframe(status->tng, frame, -1);
508 int write_trxframe(t_trxstatus *status, t_trxframe *fr, gmx_conect gc)
524 gmx_tng_set_compression_precision(status->tng, prec);
525 write_tng_frame(status, fr);
530 switch (gmx_fio_getftp(status->fio))
537 gmx_fatal(FARGS, "Need coordinates to write a %s trajectory",
538 ftp2ext(gmx_fio_getftp(status->fio)));
543 switch (gmx_fio_getftp(status->fio))
546 write_xtc(status->fio, fr->natoms, fr->step, fr->time, fr->box, fr->x, prec);
549 fwrite_trn(status->fio, fr->step, fr->time, fr->lambda, fr->box, fr->natoms,
550 fr->bX ? fr->x : NULL, fr->bV ? fr->v : NULL, fr->bF ? fr->f : NULL);
558 gmx_fatal(FARGS, "Can not write a %s file without atom names",
559 ftp2ext(gmx_fio_getftp(status->fio)));
561 sprintf(title, "frame t= %.3f", fr->time);
562 if (gmx_fio_getftp(status->fio) == efGRO)
564 write_hconf_p(gmx_fio_getfp(status->fio), title, fr->atoms,
565 prec2ndec(prec), fr->x, fr->bV ? fr->v : NULL, fr->box);
569 write_pdbfile(gmx_fio_getfp(status->fio), title,
570 fr->atoms, fr->x, fr->bPBC ? fr->ePBC : -1, fr->box,
571 ' ', fr->step, gc, TRUE);
575 write_g96_conf(gmx_fio_getfp(status->fio), fr, -1, NULL);
578 gmx_fatal(FARGS, "Sorry, write_trxframe can not write %s",
579 ftp2ext(gmx_fio_getftp(status->fio)));
586 int write_trx(t_trxstatus *status, int nind, const atom_id *ind, t_atoms *atoms,
587 int step, real time, matrix box, rvec x[], rvec *v,
592 clear_trxframe(&fr, TRUE);
597 fr.bAtoms = atoms != NULL;
604 copy_mat(box, fr.box);
606 return write_trxframe_indexed(status, &fr, nind, ind, gc);
609 void close_trx(t_trxstatus *status)
611 gmx_tng_close(&status->tng);
614 gmx_fio_close(status->fio);
619 t_trxstatus *open_trx(const char *outfile, const char *filemode)
622 if (filemode[0] != 'w' && filemode[0] != 'a' && filemode[1] != '+')
624 gmx_fatal(FARGS, "Sorry, write_trx can only write");
630 stat->fio = gmx_fio_open(outfile, filemode);
634 static gmx_bool gmx_next_frame(t_trxstatus *status, t_trxframe *fr)
641 if (fread_trnheader(status->fio, &sh, &bOK))
643 fr->bDouble = sh.bDouble;
644 fr->natoms = sh.natoms;
650 fr->bFepState = TRUE;
651 fr->lambda = sh.lambda;
652 fr->bBox = sh.box_size > 0;
653 if (fr->flags & (TRX_READ_X | TRX_NEED_X))
657 snew(fr->x, sh.natoms);
659 fr->bX = sh.x_size > 0;
661 if (fr->flags & (TRX_READ_V | TRX_NEED_V))
665 snew(fr->v, sh.natoms);
667 fr->bV = sh.v_size > 0;
669 if (fr->flags & (TRX_READ_F | TRX_NEED_F))
673 snew(fr->f, sh.natoms);
675 fr->bF = sh.f_size > 0;
677 if (fread_htrn(status->fio, &sh, fr->box, fr->x, fr->v, fr->f))
683 fr->not_ok = DATA_NOT_OK;
689 fr->not_ok = HEADER_NOT_OK;
695 static gmx_bool pdb_next_x(t_trxstatus *status, FILE *fp, t_trxframe *fr)
699 int ePBC, model_nr, na;
700 char title[STRLEN], *time;
703 atoms.nr = fr->natoms;
705 atoms.pdbinfo = NULL;
706 /* the other pointers in atoms should not be accessed if these are NULL */
708 na = read_pdbfile(fp, title, &model_nr, &atoms, fr->x, &ePBC, boxpdb, TRUE, NULL);
709 set_trxframe_ePBC(fr, ePBC);
710 if (nframes_read(status) == 0)
712 fprintf(stderr, " '%s', %d atoms\n", title, fr->natoms);
717 fr->bBox = (boxpdb[XX][XX] != 0.0);
720 copy_mat(boxpdb, fr->box);
723 if (model_nr != NOTSET)
728 time = strstr(title, " t= ");
732 sscanf(time+4, "%lf", &dbl);
733 fr->time = (real)dbl;
738 /* this is a bit dirty, but it will work: if no time is read from
739 comment line in pdb file, set time to current frame number */
742 fr->time = (real)fr->step;
746 fr->time = (real)nframes_read(status);
755 if (na != fr->natoms)
757 gmx_fatal(FARGS, "Number of atoms in pdb frame %d is %d instead of %d",
758 nframes_read(status), na, fr->natoms);
764 static int pdb_first_x(t_trxstatus *status, FILE *fp, t_trxframe *fr)
768 fprintf(stderr, "Reading frames from pdb file");
770 get_pdb_coordnum(fp, &fr->natoms);
773 gmx_fatal(FARGS, "\nNo coordinates in pdb file\n");
776 snew(fr->x, fr->natoms);
777 pdb_next_x(status, fp, fr);
782 gmx_bool read_next_frame(const output_env_t oenv, t_trxstatus *status, t_trxframe *fr)
786 gmx_bool bOK, bRet, bMissingData = FALSE, bSkip = FALSE;
795 clear_trxframe(fr, FALSE);
801 /* Special treatment for TNG files */
806 ftp = gmx_fio_getftp(status->fio);
811 bRet = gmx_next_frame(status, fr);
814 /* Checkpoint files can not contain mulitple frames */
817 read_g96_conf(gmx_fio_getfp(status->fio), NULL, fr,
818 status->persistent_line);
819 bRet = (fr->natoms > 0);
823 * Sometimes is off by one frame
824 * and sometimes reports frame not present/file not seekable
826 /* DvdS 2005-05-31: this has been fixed along with the increased
827 * accuracy of the control over -b and -e options.
829 if (bTimeSet(TBEGIN) && (fr->tf < rTimeValue(TBEGIN)))
831 if (xtc_seek_time(status->fio, rTimeValue(TBEGIN), fr->natoms, TRUE))
833 gmx_fatal(FARGS, "Specified frame (time %f) doesn't exist or file corrupt/inconsistent.",
838 bRet = read_next_xtc(status->fio, fr->natoms, &fr->step, &fr->time, fr->box,
839 fr->x, &fr->prec, &bOK);
840 fr->bPrec = (bRet && fr->prec > 0);
847 /* Actually the header could also be not ok,
848 but from bOK from read_next_xtc this can't be distinguished */
849 fr->not_ok = DATA_NOT_OK;
853 bRet = gmx_read_next_tng_frame(status->tng, fr, NULL, 0);
856 bRet = pdb_next_x(status, gmx_fio_getfp(status->fio), fr);
859 bRet = gro_next_x_or_v(gmx_fio_getfp(status->fio), fr);
862 #ifdef GMX_USE_PLUGINS
863 bRet = read_next_vmd_frame(fr);
865 gmx_fatal(FARGS, "DEATH HORROR in read_next_frame ftp=%s,status=%s",
866 ftp2ext(gmx_fio_getftp(status->fio)),
867 gmx_fio_getname(status->fio));
874 bMissingData = (((fr->flags & TRX_NEED_X) && !fr->bX) ||
875 ((fr->flags & TRX_NEED_V) && !fr->bV) ||
876 ((fr->flags & TRX_NEED_F) && !fr->bF));
880 ct = check_times2(fr->time, fr->t0, fr->bDouble);
881 if (ct == 0 || ((fr->flags & TRX_DONT_SKIP) && ct < 0))
883 printcount(status, oenv, fr->time, FALSE);
891 printcount(status, oenv, fr->time, TRUE);
898 while (bRet && (bMissingData || bSkip));
902 printlast(status, oenv, pt);
905 printincomp(status, fr);
912 int read_first_frame(const output_env_t oenv, t_trxstatus **status,
913 const char *fn, t_trxframe *fr, int flags)
916 gmx_bool bFirst, bOK;
918 int ftp = fn2ftp(fn);
919 gmx_int64_t *tng_ids;
921 clear_trxframe(fr, TRUE);
928 status_init( *status );
929 (*status)->nxframe = 1;
934 /* Special treatment for TNG files */
935 gmx_tng_open(fn, 'r', &(*status)->tng);
939 fio = (*status)->fio = gmx_fio_open(fn, "r");
946 read_checkpoint_trxframe(fio, fr);
950 /* Can not rewind a compressed file, so open it twice */
951 if (!(*status)->persistent_line)
953 /* allocate the persistent line */
954 snew((*status)->persistent_line, STRLEN+1);
956 read_g96_conf(gmx_fio_getfp(fio), fn, fr, (*status)->persistent_line);
958 clear_trxframe(fr, FALSE);
959 if (flags & (TRX_READ_X | TRX_NEED_X))
961 snew(fr->x, fr->natoms);
963 if (flags & (TRX_READ_V | TRX_NEED_V))
965 snew(fr->v, fr->natoms);
967 fio = (*status)->fio = gmx_fio_open(fn, "r");
970 if (read_first_xtc(fio, &fr->natoms, &fr->step, &fr->time, fr->box, &fr->x,
971 &fr->prec, &bOK) == 0)
974 fr->not_ok = DATA_NOT_OK;
979 printincomp(*status, fr);
983 fr->bPrec = (fr->prec > 0);
988 printcount(*status, oenv, fr->time, FALSE);
994 if (!gmx_read_next_tng_frame((*status)->tng, fr, NULL, 0))
996 fr->not_ok = DATA_NOT_OK;
998 printincomp(*status, fr);
1002 printcount(*status, oenv, fr->time, FALSE);
1007 pdb_first_x(*status, gmx_fio_getfp(fio), fr);
1010 printcount(*status, oenv, fr->time, FALSE);
1015 if (gro_first_x_or_v(gmx_fio_getfp(fio), fr))
1017 printcount(*status, oenv, fr->time, FALSE);
1022 #ifdef GMX_USE_PLUGINS
1023 fprintf(stderr, "The file format of %s is not a known trajectory format to GROMACS.\n"
1024 "Please make sure that the file is a trajectory!\n"
1025 "GROMACS will now assume it to be a trajectory and will try to open it using the VMD plug-ins.\n"
1026 "This will only work in case the VMD plugins are found and it is a trajectory format supported by VMD.\n", fn);
1027 gmx_fio_fp_close(fio); /*only close the file without removing FIO entry*/
1028 if (!read_first_vmd_frame(fn, fr))
1030 gmx_fatal(FARGS, "Not supported in read_first_frame: %s", fn);
1033 gmx_fatal(FARGS, "Not supported in read_first_frame: %s. Please make sure that the file is a trajectory.\n"
1034 "GROMACS is not compiled with plug-in support. Thus it cannot read non-GROMACS trajectory formats using the VMD plug-ins.\n"
1035 "Please compile with plug-in support if you want to read non-GROMACS trajectory formats.\n", fn);
1041 /* Return FALSE if we read a frame that's past the set ending time. */
1042 if (!bFirst && (!(fr->flags & TRX_DONT_SKIP) && check_times(fr->time) > 0))
1049 (!(fr->flags & TRX_DONT_SKIP) && check_times(fr->time) < 0))
1051 /* Read a frame when no frame was read or the first was skipped */
1052 if (!read_next_frame(oenv, *status, fr))
1059 return (fr->natoms > 0);
1062 /***** C O O R D I N A T E S T U F F *****/
1064 int read_first_x(const output_env_t oenv, t_trxstatus **status, const char *fn,
1065 real *t, rvec **x, matrix box)
1069 read_first_frame(oenv, status, fn, &fr, TRX_NEED_X);
1071 snew((*status)->xframe, 1);
1072 (*status)->nxframe = 1;
1073 (*(*status)->xframe) = fr;
1074 *t = (*status)->xframe->time;
1075 *x = (*status)->xframe->x;
1076 copy_mat((*status)->xframe->box, box);
1078 return (*status)->xframe->natoms;
1081 gmx_bool read_next_x(const output_env_t oenv, t_trxstatus *status, real *t,
1082 rvec x[], matrix box)
1086 status->xframe->x = x;
1087 /*xframe[status].x = x;*/
1088 bRet = read_next_frame(oenv, status, status->xframe);
1089 *t = status->xframe->time;
1090 copy_mat(status->xframe->box, box);
1095 void close_trj(t_trxstatus *status)
1097 gmx_tng_close(&status->tng);
1100 gmx_fio_close(status->fio);
1103 /* The memory in status->xframe is lost here,
1104 * but the read_first_x/read_next_x functions are deprecated anyhow.
1105 * read_first_frame/read_next_frame and close_trx should be used.
1110 void rewind_trj(t_trxstatus *status)
1114 gmx_fio_rewind(status->fio);
1117 /***** T O P O L O G Y S T U F F ******/
1119 t_topology *read_top(const char *fn, int *ePBC)
1125 epbc = read_tpx_top(fn, NULL, NULL, &natoms, NULL, NULL, NULL, top);