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.
43 #include "gromacs/fileio/gmxfio.h"
44 #include "gromacs/fileio/xdrf.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/utility/fatalerror.h"
47 #include "gromacs/utility/futil.h"
48 #include "gromacs/utility/smalloc.h"
50 #define XTC_MAGIC 1995
53 static int xdr_r2f(XDR *xdrs, real *r, gmx_bool gmx_unused bRead)
63 ret = xdr_float(xdrs, &f);
71 return xdr_float(xdrs, (float *)r);
76 t_fileio *open_xtc(const char *fn, const char *mode)
78 return gmx_fio_open(fn, mode);
81 void close_xtc(t_fileio *fio)
86 static void check_xtc_magic(int magic)
88 if (magic != XTC_MAGIC)
90 gmx_fatal(FARGS, "Magic Number Error in XTC file (read %d, should be %d)",
95 int xtc_check(const char *str, gmx_bool bResult, const char *file, int line)
101 fprintf(debug, "\nXTC error: read/write of %s failed, "
102 "source file %s, line %d\n", str, file, line);
109 void xtc_check_fat_err(const char *str, gmx_bool bResult, const char *file, int line)
113 gmx_fatal(FARGS, "XTC read/write of %s failed, "
114 "source file %s, line %d\n", str, file, line);
118 static int xtc_header(XDR *xd, int *magic, int *natoms, int *step, real *time,
119 gmx_bool bRead, gmx_bool *bOK)
123 if (xdr_int(xd, magic) == 0)
127 result = XTC_CHECK("natoms", xdr_int(xd, natoms)); /* number of atoms */
130 result = XTC_CHECK("step", xdr_int(xd, step)); /* frame number */
134 result = XTC_CHECK("time", xdr_r2f(xd, time, bRead)); /* time */
136 *bOK = (result != 0);
141 static int xtc_coord(XDR *xd, int *natoms, matrix box, rvec *x, real *prec, gmx_bool bRead)
151 for (i = 0; ((i < DIM) && result); i++)
153 for (j = 0; ((j < DIM) && result); j++)
155 result = XTC_CHECK("box", xdr_r2f(xd, &(box[i][j]), bRead));
165 /* allocate temp. single-precision array */
166 snew(ftmp, (*natoms)*DIM);
168 /* Copy data to temp. array if writing */
171 for (i = 0; (i < *natoms); i++)
173 ftmp[DIM*i+XX] = x[i][XX];
174 ftmp[DIM*i+YY] = x[i][YY];
175 ftmp[DIM*i+ZZ] = x[i][ZZ];
179 result = XTC_CHECK("x", xdr3dfcoord(xd, ftmp, natoms, &fprec));
181 /* Copy from temp. array if reading */
184 for (i = 0; (i < *natoms); i++)
186 x[i][XX] = ftmp[DIM*i+XX];
187 x[i][YY] = ftmp[DIM*i+YY];
188 x[i][ZZ] = ftmp[DIM*i+ZZ];
194 result = XTC_CHECK("x", xdr3dfcoord(xd, x[0], natoms, prec));
202 int write_xtc(t_fileio *fio,
203 int natoms, int step, real time,
204 matrix box, rvec *x, real prec)
206 int magic_number = XTC_MAGIC;
213 /* This means the fio object is not being used, e.g. because
214 we are actually writing TNG output. We still have to return
215 a pseudo-success value, to keep some callers happy. */
219 xd = gmx_fio_getxdr(fio);
220 /* write magic number and xtc identidier */
221 if (xtc_header(xd, &magic_number, &natoms, &step, &time, FALSE, &bDum) == 0)
227 bOK = xtc_coord(xd, &natoms, box, x, &prec, FALSE); /* bOK will be 1 if writing went well */
231 if (gmx_fio_flush(fio) != 0)
236 return bOK; /* 0 if bad, 1 if writing went well */
239 int read_first_xtc(t_fileio *fio, int *natoms, int *step, real *time,
240 matrix box, rvec **x, real *prec, gmx_bool *bOK)
246 xd = gmx_fio_getxdr(fio);
248 /* read header and malloc x */
249 if (!xtc_header(xd, &magic, natoms, step, time, TRUE, bOK))
254 /* Check magic number */
255 check_xtc_magic(magic);
259 *bOK = xtc_coord(xd, natoms, box, *x, prec, TRUE);
264 int read_next_xtc(t_fileio* fio,
265 int natoms, int *step, real *time,
266 matrix box, rvec *x, real *prec, gmx_bool *bOK)
273 xd = gmx_fio_getxdr(fio);
276 if (!xtc_header(xd, &magic, &n, step, time, TRUE, bOK))
281 /* Check magic number */
282 check_xtc_magic(magic);
286 gmx_fatal(FARGS, "Frame contains more atoms (%d) than expected (%d)",
290 *bOK = xtc_coord(xd, &natoms, box, x, prec, TRUE);