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,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
48 #include "gromacs/fileio/checkpoint.h"
49 #include "gromacs/fileio/confio.h"
50 #include "gromacs/fileio/filetypes.h"
51 #include "gromacs/fileio/g96io.h"
52 #include "gromacs/fileio/gmxfio.h"
53 #include "gromacs/fileio/gmxfio_xdr.h"
54 #include "gromacs/fileio/groio.h"
55 #include "gromacs/fileio/oenv.h"
56 #include "gromacs/fileio/pdbio.h"
57 #include "gromacs/fileio/timecontrol.h"
58 #include "gromacs/fileio/tngio.h"
59 #include "gromacs/fileio/tpxio.h"
60 #include "gromacs/fileio/trrio.h"
61 #include "gromacs/fileio/xdrf.h"
62 #include "gromacs/fileio/xtcio.h"
63 #include "gromacs/math/vec.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/topology/atoms.h"
66 #include "gromacs/topology/symtab.h"
67 #include "gromacs/topology/topology.h"
68 #include "gromacs/trajectory/trajectoryframe.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/futil.h"
71 #include "gromacs/utility/gmxassert.h"
72 #include "gromacs/utility/smalloc.h"
75 # include "gromacs/fileio/vmdio.h"
78 /* defines for frame counter output */
85 int flags; /* flags for read_first/next_frame */
87 real t0; /* time of the first frame, needed *
88 * for skipping frames with -dt */
89 real tf; /* internal frame time */
92 gmx_tng_trajectory_t tng;
96 char* persistent_line; /* Persistent line for reading g96 trajectories */
98 gmx_vmdplugin_t* vmdplugin;
102 /* utility functions */
104 gmx_bool bRmod_fd(double a, double b, double c, gmx_bool bDouble)
109 tol = 2 * (bDouble ? GMX_DOUBLE_EPS : GMX_FLOAT_EPS);
111 iq = static_cast<int>((a - b + tol * a) / c);
113 return fabs(a - b - c * iq) <= tol * fabs(a);
117 int check_times2(real t, real t0, gmx_bool bDouble)
122 /* since t is float, we can not use double precision for bRmod */
127 if ((!bTimeSet(TBEGIN) || (t >= rTimeValue(TBEGIN))) && (!bTimeSet(TEND) || (t <= rTimeValue(TEND))))
129 if (bTimeSet(TDELTA) && !bRmod_fd(t, t0, rTimeValue(TDELTA), bDouble))
138 else if (bTimeSet(TEND) && (t >= rTimeValue(TEND)))
144 fprintf(debug, "t=%g, t0=%g, b=%g, e=%g, dt=%g: r=%d\n", t, t0, rTimeValue(TBEGIN),
145 rTimeValue(TEND), rTimeValue(TDELTA), r);
150 int check_times(real t)
152 return check_times2(t, t, FALSE);
155 static void initcount(t_trxstatus* status)
157 status->__frame = -1;
160 static void status_init(t_trxstatus* status)
163 status->xframe = nullptr;
164 status->fio = nullptr;
165 status->__frame = -1;
168 status->persistent_line = nullptr;
169 status->tng = nullptr;
173 int nframes_read(t_trxstatus* status)
175 return status->__frame;
178 static void printcount_(t_trxstatus* status, const gmx_output_env_t* oenv, const char* l, real t)
180 if ((status->__frame < 2 * SKIP1 || status->__frame % SKIP1 == 0)
181 && (status->__frame < 2 * SKIP2 || status->__frame % SKIP2 == 0)
182 && (status->__frame < 2 * SKIP3 || status->__frame % SKIP3 == 0)
183 && output_env_get_trajectory_io_verbosity(oenv) != 0)
185 fprintf(stderr, "\r%-14s %6d time %8.3f ", l, status->__frame, output_env_conv_time(oenv, t));
190 static void printcount(t_trxstatus* status, const gmx_output_env_t* oenv, real t, gmx_bool bSkip)
193 printcount_(status, oenv, bSkip ? "Skipping frame" : "Reading frame", t);
196 static void printlast(t_trxstatus* status, const gmx_output_env_t* oenv, real t)
198 printcount_(status, oenv, "Last frame", t);
199 fprintf(stderr, "\n");
203 static void printincomp(t_trxstatus* status, t_trxframe* fr)
205 if (fr->not_ok & HEADER_NOT_OK)
207 fprintf(stderr, "WARNING: Incomplete header: nr %d time %g\n", status->__frame + 1, fr->time);
211 fprintf(stderr, "WARNING: Incomplete frame: nr %d time %g\n", status->__frame + 1, fr->time);
216 int prec2ndec(real prec)
220 gmx_fatal(FARGS, "DEATH HORROR prec (%g) <= 0 in prec2ndec", prec);
223 return gmx::roundToInt(log(prec) / log(10.0));
226 real ndec2prec(int ndec)
228 return pow(10.0, ndec);
231 t_fileio* trx_get_fileio(t_trxstatus* status)
236 float trx_get_time_of_final_frame(t_trxstatus* status)
238 t_fileio* stfio = trx_get_fileio(status);
239 int filetype = gmx_fio_getftp(stfio);
243 if (filetype == efXTC)
245 lasttime = xdr_xtc_get_last_frame_time(gmx_fio_getfp(stfio), gmx_fio_getxdr(stfio),
246 status->natoms, &bOK);
249 gmx_fatal(FARGS, "Error reading last frame. Maybe seek not supported.");
252 else if (filetype == efTNG)
254 gmx_tng_trajectory_t tng = status->tng;
257 gmx_fatal(FARGS, "Error opening TNG file.");
259 lasttime = gmx_tng_get_time_of_final_frame(tng);
263 gmx_incons("Only supported for TNG and XTC");
268 void clear_trxframe(t_trxframe* fr, gmx_bool bFirst)
274 fr->bFepState = FALSE;
300 void set_trxframe_ePBC(t_trxframe* fr, int ePBC)
302 fr->bPBC = (ePBC == -1);
306 int write_trxframe_indexed(t_trxstatus* status, const t_trxframe* fr, int nind, const int* ind, gmx_conect gc)
309 rvec *xout = nullptr, *vout = nullptr, *fout = nullptr;
326 else if (status->fio)
328 ftp = gmx_fio_getftp(status->fio);
332 gmx_incons("No input file available");
342 gmx_fatal(FARGS, "Need coordinates to write a %s trajectory", ftp2ext(ftp));
354 for (i = 0; i < nind; i++)
356 copy_rvec(fr->v[ind[i]], vout[i]);
362 for (i = 0; i < nind; i++)
364 copy_rvec(fr->f[ind[i]], fout[i]);
370 for (i = 0; i < nind; i++)
372 copy_rvec(fr->x[ind[i]], xout[i]);
380 for (i = 0; i < nind; i++)
382 copy_rvec(fr->x[ind[i]], xout[i]);
391 case efTNG: gmx_write_tng_from_trxframe(status->tng, fr, nind); break;
392 case efXTC: write_xtc(status->fio, nind, fr->step, fr->time, fr->box, xout, prec); break;
394 gmx_trr_write_frame(status->fio, nframes_read(status), fr->time, fr->step, fr->box,
395 nind, xout, vout, fout);
403 gmx_fatal(FARGS, "Can not write a %s file without atom names", ftp2ext(ftp));
405 sprintf(title, "frame t= %.3f", fr->time);
408 write_hconf_indexed_p(gmx_fio_getfp(status->fio), title, fr->atoms, nind, ind,
409 fr->x, fr->bV ? fr->v : nullptr, fr->box);
413 write_pdbfile_indexed(gmx_fio_getfp(status->fio), title, fr->atoms, fr->x, -1,
414 fr->box, ' ', fr->step, nind, ind, gc, FALSE);
418 sprintf(title, "frame t= %.3f", fr->time);
419 write_g96_conf(gmx_fio_getfp(status->fio), title, fr, nind, ind);
421 default: gmx_fatal(FARGS, "Sorry, write_trxframe_indexed can not write %s", ftp2ext(ftp));
438 case efXTC: sfree(xout); break;
445 t_trxstatus* trjtools_gmx_prepare_tng_writing(const char* filename,
450 const gmx_mtop_t* mtop,
451 gmx::ArrayRef<const int> index,
452 const char* index_group_name)
454 if (filemode != 'w' && filemode != 'a')
456 gmx_incons("Sorry, can only prepare for TNG output.");
464 gmx_prepare_tng_writing(filename, filemode, &in->tng, &out->tng, natoms, mtop, index,
467 else if ((infile) && (efTNG == fn2ftp(infile)))
469 gmx_tng_trajectory_t tng_in;
470 gmx_tng_open(infile, 'r', &tng_in);
472 gmx_prepare_tng_writing(filename, filemode, &tng_in, &out->tng, natoms, mtop, index,
477 // we start from a file that is not a tng file or have been unable to load the
478 // input file, so we need to populate the fields independently of it
479 gmx_prepare_tng_writing(filename, filemode, nullptr, &out->tng, natoms, mtop, index,
485 void write_tng_frame(t_trxstatus* status, t_trxframe* frame)
487 gmx_write_tng_from_trxframe(status->tng, frame, -1);
490 int write_trxframe(t_trxstatus* status, t_trxframe* fr, gmx_conect gc)
507 gmx_tng_set_compression_precision(status->tng, prec);
508 write_tng_frame(status, fr);
513 switch (gmx_fio_getftp(status->fio))
519 gmx_fatal(FARGS, "Need coordinates to write a %s trajectory",
520 ftp2ext(gmx_fio_getftp(status->fio)));
525 switch (gmx_fio_getftp(status->fio))
528 write_xtc(status->fio, fr->natoms, fr->step, fr->time, fr->box, fr->x, prec);
531 gmx_trr_write_frame(status->fio, fr->step, fr->time, fr->lambda, fr->box, fr->natoms,
532 fr->bX ? fr->x : nullptr, fr->bV ? fr->v : nullptr,
533 fr->bF ? fr->f : nullptr);
541 gmx_fatal(FARGS, "Can not write a %s file without atom names",
542 ftp2ext(gmx_fio_getftp(status->fio)));
544 sprintf(title, "frame t= %.3f", fr->time);
545 if (gmx_fio_getftp(status->fio) == efGRO)
547 write_hconf_p(gmx_fio_getfp(status->fio), title, fr->atoms, fr->x,
548 fr->bV ? fr->v : nullptr, fr->box);
552 write_pdbfile(gmx_fio_getfp(status->fio), title, fr->atoms, fr->x,
553 fr->bPBC ? fr->ePBC : -1, fr->box, ' ', fr->step, gc);
556 case efG96: write_g96_conf(gmx_fio_getfp(status->fio), title, fr, -1, nullptr); break;
558 gmx_fatal(FARGS, "Sorry, write_trxframe can not write %s",
559 ftp2ext(gmx_fio_getftp(status->fio)));
565 int write_trx(t_trxstatus* status,
568 const t_atoms* atoms,
578 clear_trxframe(&fr, TRUE);
583 fr.bAtoms = atoms != nullptr;
584 fr.atoms = const_cast<t_atoms*>(atoms);
587 fr.bV = v != nullptr;
590 copy_mat(box, fr.box);
592 return write_trxframe_indexed(status, &fr, nind, ind, gc);
595 void close_trx(t_trxstatus* status)
597 if (status == nullptr)
601 gmx_tng_close(&status->tng);
604 gmx_fio_close(status->fio);
606 sfree(status->persistent_line);
608 sfree(status->vmdplugin);
610 /* The memory in status->xframe is lost here,
611 * but the read_first_x/read_next_x functions are deprecated anyhow.
612 * read_first_frame/read_next_frame and close_trx should be used.
617 t_trxstatus* open_trx(const char* outfile, const char* filemode)
620 if (filemode[0] != 'w' && filemode[0] != 'a' && filemode[1] != '+')
622 gmx_fatal(FARGS, "Sorry, write_trx can only write");
628 stat->fio = gmx_fio_open(outfile, filemode);
632 static gmx_bool gmx_next_frame(t_trxstatus* status, t_trxframe* fr)
639 if (gmx_trr_read_frame_header(status->fio, &sh, &bOK))
641 fr->bDouble = sh.bDouble;
642 fr->natoms = sh.natoms;
648 fr->bFepState = TRUE;
649 fr->lambda = sh.lambda;
650 fr->bBox = sh.box_size > 0;
651 if (status->flags & (TRX_READ_X | TRX_NEED_X))
653 if (fr->x == nullptr)
655 snew(fr->x, sh.natoms);
657 fr->bX = sh.x_size > 0;
659 if (status->flags & (TRX_READ_V | TRX_NEED_V))
661 if (fr->v == nullptr)
663 snew(fr->v, sh.natoms);
665 fr->bV = sh.v_size > 0;
667 if (status->flags & (TRX_READ_F | TRX_NEED_F))
669 if (fr->f == nullptr)
671 snew(fr->f, sh.natoms);
673 fr->bF = sh.f_size > 0;
675 if (gmx_trr_read_frame_data(status->fio, &sh, fr->box, fr->x, fr->v, fr->f))
681 fr->not_ok = DATA_NOT_OK;
686 fr->not_ok = HEADER_NOT_OK;
692 static gmx_bool pdb_next_x(t_trxstatus* status, FILE* fp, t_trxframe* fr)
697 // Initiate model_nr to -1 rather than NOTSET.
698 // It is not worthwhile introducing extra variables in the
699 // read_pdbfile call to verify that a model_nr was read.
700 int ePBC, model_nr = -1, na;
701 char title[STRLEN], *time, *step;
704 atoms.nr = fr->natoms;
705 atoms.atom = nullptr;
706 atoms.pdbinfo = nullptr;
707 /* the other pointers in atoms should not be accessed if these are NULL */
710 na = read_pdbfile(fp, title, &model_nr, &atoms, symtab, fr->x, &ePBC, boxpdb, TRUE, nullptr);
713 set_trxframe_ePBC(fr, ePBC);
714 if (nframes_read(status) == 0)
716 fprintf(stderr, " '%s', %d atoms\n", title, fr->natoms);
721 fr->bBox = (boxpdb[XX][XX] != 0.0);
724 copy_mat(boxpdb, fr->box);
728 step = std::strstr(title, " step= ");
729 fr->bStep = ((step != nullptr) && sscanf(step + 7, "%" SCNd64, &fr->step) == 1);
732 time = std::strstr(title, " t= ");
733 fr->bTime = ((time != nullptr) && sscanf(time + 4, "%lf", &dbl) == 1);
742 if (na != fr->natoms)
744 gmx_fatal(FARGS, "Number of atoms in pdb frame %d is %d instead of %d",
745 nframes_read(status), na, fr->natoms);
751 static int pdb_first_x(t_trxstatus* status, FILE* fp, t_trxframe* fr)
755 fprintf(stderr, "Reading frames from pdb file");
757 get_pdb_coordnum(fp, &fr->natoms);
760 gmx_fatal(FARGS, "\nNo coordinates in pdb file\n");
763 snew(fr->x, fr->natoms);
764 pdb_next_x(status, fp, fr);
769 bool read_next_frame(const gmx_output_env_t* oenv, t_trxstatus* status, t_trxframe* fr)
773 gmx_bool bOK, bMissingData = FALSE, bSkip = FALSE;
781 clear_trxframe(fr, FALSE);
785 /* Special treatment for TNG files */
790 ftp = gmx_fio_getftp(status->fio);
794 case efTRR: bRet = gmx_next_frame(status, fr); break;
796 /* Checkpoint files can not contain mulitple frames */
800 t_symtab* symtab = nullptr;
801 read_g96_conf(gmx_fio_getfp(status->fio), nullptr, nullptr, fr, symtab,
802 status->persistent_line);
803 bRet = (fr->natoms > 0);
807 if (bTimeSet(TBEGIN) && (status->tf < rTimeValue(TBEGIN)))
809 if (xtc_seek_time(status->fio, rTimeValue(TBEGIN), fr->natoms, TRUE))
812 "Specified frame (time %f) doesn't exist or file "
813 "corrupt/inconsistent.",
818 bRet = (read_next_xtc(status->fio, fr->natoms, &fr->step, &fr->time, fr->box, fr->x,
821 fr->bPrec = (bRet && fr->prec > 0);
828 /* Actually the header could also be not ok,
829 but from bOK from read_next_xtc this can't be distinguished */
830 fr->not_ok = DATA_NOT_OK;
833 case efTNG: bRet = gmx_read_next_tng_frame(status->tng, fr, nullptr, 0); break;
834 case efPDB: bRet = pdb_next_x(status, gmx_fio_getfp(status->fio), fr); break;
835 case efGRO: bRet = gro_next_x_or_v(gmx_fio_getfp(status->fio), fr); break;
838 bRet = read_next_vmd_frame(status->vmdplugin, fr);
840 gmx_fatal(FARGS, "DEATH HORROR in read_next_frame ftp=%s,status=%s",
841 ftp2ext(gmx_fio_getftp(status->fio)), gmx_fio_getname(status->fio));
844 status->tf = fr->time;
848 bMissingData = ((((status->flags & TRX_NEED_X) != 0) && !fr->bX)
849 || (((status->flags & TRX_NEED_V) != 0) && !fr->bV)
850 || (((status->flags & TRX_NEED_F) != 0) && !fr->bF));
854 ct = check_times2(fr->time, status->t0, fr->bDouble);
855 if (ct == 0 || ((status->flags & TRX_DONT_SKIP) && ct < 0))
857 printcount(status, oenv, fr->time, FALSE);
865 printcount(status, oenv, fr->time, TRUE);
871 } while (bRet && (bMissingData || bSkip));
875 printlast(status, oenv, pt);
878 printincomp(status, fr);
885 bool read_first_frame(const gmx_output_env_t* oenv, t_trxstatus** status, const char* fn, t_trxframe* fr, int flags)
887 t_fileio* fio = nullptr;
888 gmx_bool bFirst, bOK;
889 int ftp = fn2ftp(fn);
891 clear_trxframe(fr, TRUE);
897 status_init(*status);
899 (*status)->flags = flags;
903 /* Special treatment for TNG files */
904 gmx_tng_open(fn, 'r', &(*status)->tng);
908 fio = (*status)->fio = gmx_fio_open(fn, "r");
914 read_checkpoint_trxframe(fio, fr);
919 /* Can not rewind a compressed file, so open it twice */
920 if (!(*status)->persistent_line)
922 /* allocate the persistent line */
923 snew((*status)->persistent_line, STRLEN + 1);
925 t_symtab* symtab = nullptr;
926 read_g96_conf(gmx_fio_getfp(fio), fn, nullptr, fr, symtab, (*status)->persistent_line);
928 clear_trxframe(fr, FALSE);
929 if (flags & (TRX_READ_X | TRX_NEED_X))
931 snew(fr->x, fr->natoms);
933 if (flags & (TRX_READ_V | TRX_NEED_V))
935 snew(fr->v, fr->natoms);
937 (*status)->fio = gmx_fio_open(fn, "r");
941 if (read_first_xtc(fio, &fr->natoms, &fr->step, &fr->time, fr->box, &fr->x, &fr->prec, &bOK) == 0)
943 GMX_RELEASE_ASSERT(!bOK,
944 "Inconsistent results - OK status from read_first_xtc, but 0 "
946 fr->not_ok = DATA_NOT_OK;
951 printincomp(*status, fr);
955 fr->bPrec = (fr->prec > 0);
960 printcount(*status, oenv, fr->time, FALSE);
966 if (!gmx_read_next_tng_frame((*status)->tng, fr, nullptr, 0))
968 fr->not_ok = DATA_NOT_OK;
970 printincomp(*status, fr);
974 printcount(*status, oenv, fr->time, FALSE);
979 pdb_first_x(*status, gmx_fio_getfp(fio), fr);
982 printcount(*status, oenv, fr->time, FALSE);
987 if (gro_first_x_or_v(gmx_fio_getfp(fio), fr))
989 printcount(*status, oenv, fr->time, FALSE);
996 "The file format of %s is not a known trajectory format to GROMACS.\n"
997 "Please make sure that the file is a trajectory!\n"
998 "GROMACS will now assume it to be a trajectory and will try to open it using "
999 "the VMD plug-ins.\n"
1000 "This will only work in case the VMD plugins are found and it is a trajectory "
1001 "format supported by VMD.\n",
1003 gmx_fio_fp_close(fio); /*only close the file without removing FIO entry*/
1004 if (!read_first_vmd_frame(fn, &(*status)->vmdplugin, fr))
1006 gmx_fatal(FARGS, "Not supported in read_first_frame: %s", fn);
1010 "Not supported in read_first_frame: %s. Please make sure that the file is a "
1012 "GROMACS is not compiled with plug-in support. Thus it cannot read "
1013 "non-GROMACS trajectory formats using the VMD plug-ins.\n"
1014 "Please compile with plug-in support if you want to read non-GROMACS "
1015 "trajectory formats.\n",
1019 (*status)->tf = fr->time;
1021 /* Return FALSE if we read a frame that's past the set ending time. */
1022 if (!bFirst && (!(flags & TRX_DONT_SKIP) && check_times(fr->time) > 0))
1024 (*status)->t0 = fr->time;
1028 if (bFirst || (!(flags & TRX_DONT_SKIP) && check_times(fr->time) < 0))
1030 /* Read a frame when no frame was read or the first was skipped */
1031 if (!read_next_frame(oenv, *status, fr))
1036 (*status)->t0 = fr->time;
1038 /* We need the number of atoms for random-access XTC searching, even when
1039 * we don't have access to the actual frame data.
1041 (*status)->natoms = fr->natoms;
1043 return (fr->natoms > 0);
1046 /***** C O O R D I N A T E S T U F F *****/
1048 int read_first_x(const gmx_output_env_t* oenv, t_trxstatus** status, const char* fn, real* t, rvec** x, matrix box)
1052 read_first_frame(oenv, status, fn, &fr, TRX_NEED_X);
1054 snew((*status)->xframe, 1);
1055 (*(*status)->xframe) = fr;
1056 *t = (*status)->xframe->time;
1057 *x = (*status)->xframe->x;
1058 copy_mat((*status)->xframe->box, box);
1060 return (*status)->xframe->natoms;
1063 gmx_bool read_next_x(const gmx_output_env_t* oenv, t_trxstatus* status, real* t, rvec x[], matrix box)
1067 status->xframe->x = x;
1068 /*xframe[status].x = x;*/
1069 bRet = read_next_frame(oenv, status, status->xframe);
1070 *t = status->xframe->time;
1071 copy_mat(status->xframe->box, box);
1076 void rewind_trj(t_trxstatus* status)
1080 gmx_fio_rewind(status->fio);
1083 /***** T O P O L O G Y S T U F F ******/
1085 t_topology* read_top(const char* fn, int* ePBC)
1091 epbc = read_tpx_top(fn, nullptr, nullptr, &natoms, nullptr, nullptr, top);