Add TNG writing and reading support
[gromacs.git] / src / gromacs / fileio / xtcio.c
blobd7a866597227015529072bb5c7cf59ba2ae29abb
1 /*
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.
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
41 #include <string.h>
42 #include "typedefs.h"
43 #include "xdrf.h"
44 #include "gmxfio.h"
45 #include "xtcio.h"
46 #include "smalloc.h"
47 #include "vec.h"
48 #include "futil.h"
49 #include "gmx_fatal.h"
51 #define XTC_MAGIC 1995
54 static int xdr_r2f(XDR *xdrs, real *r, gmx_bool gmx_unused bRead)
56 #ifdef GMX_DOUBLE
57 float f;
58 int ret;
60 if (!bRead)
62 f = *r;
64 ret = xdr_float(xdrs, &f);
65 if (bRead)
67 *r = f;
70 return ret;
71 #else
72 return xdr_float(xdrs, (float *)r);
73 #endif
77 t_fileio *open_xtc(const char *fn, const char *mode)
79 return gmx_fio_open(fn, mode);
82 void close_xtc(t_fileio *fio)
84 gmx_fio_close(fio);
87 static void check_xtc_magic(int magic)
89 if (magic != XTC_MAGIC)
91 gmx_fatal(FARGS, "Magic Number Error in XTC file (read %d, should be %d)",
92 magic, XTC_MAGIC);
96 int xtc_check(const char *str, gmx_bool bResult, const char *file, int line)
98 if (!bResult)
100 if (debug)
102 fprintf(debug, "\nXTC error: read/write of %s failed, "
103 "source file %s, line %d\n", str, file, line);
105 return 0;
107 return 1;
110 void xtc_check_fat_err(const char *str, gmx_bool bResult, const char *file, int line)
112 if (!bResult)
114 gmx_fatal(FARGS, "XTC read/write of %s failed, "
115 "source file %s, line %d\n", str, file, line);
119 static int xtc_header(XDR *xd, int *magic, int *natoms, int *step, real *time,
120 gmx_bool bRead, gmx_bool *bOK)
122 int result;
124 if (xdr_int(xd, magic) == 0)
126 return 0;
128 result = XTC_CHECK("natoms", xdr_int(xd, natoms)); /* number of atoms */
129 if (result)
131 result = XTC_CHECK("step", xdr_int(xd, step)); /* frame number */
133 if (result)
135 result = XTC_CHECK("time", xdr_r2f(xd, time, bRead)); /* time */
137 *bOK = (result != 0);
139 return result;
142 static int xtc_coord(XDR *xd, int *natoms, matrix box, rvec *x, real *prec, gmx_bool bRead)
144 int i, j, result;
145 #ifdef GMX_DOUBLE
146 float *ftmp;
147 float fprec;
148 #endif
150 /* box */
151 result = 1;
152 for (i = 0; ((i < DIM) && result); i++)
154 for (j = 0; ((j < DIM) && result); j++)
156 result = XTC_CHECK("box", xdr_r2f(xd, &(box[i][j]), bRead));
160 if (!result)
162 return result;
165 #ifdef GMX_DOUBLE
166 /* allocate temp. single-precision array */
167 snew(ftmp, (*natoms)*DIM);
169 /* Copy data to temp. array if writing */
170 if (!bRead)
172 for (i = 0; (i < *natoms); i++)
174 ftmp[DIM*i+XX] = x[i][XX];
175 ftmp[DIM*i+YY] = x[i][YY];
176 ftmp[DIM*i+ZZ] = x[i][ZZ];
178 fprec = *prec;
180 result = XTC_CHECK("x", xdr3dfcoord(xd, ftmp, natoms, &fprec));
182 /* Copy from temp. array if reading */
183 if (bRead)
185 for (i = 0; (i < *natoms); i++)
187 x[i][XX] = ftmp[DIM*i+XX];
188 x[i][YY] = ftmp[DIM*i+YY];
189 x[i][ZZ] = ftmp[DIM*i+ZZ];
191 *prec = fprec;
193 sfree(ftmp);
194 #else
195 result = XTC_CHECK("x", xdr3dfcoord(xd, x[0], natoms, prec));
196 #endif
198 return result;
203 int write_xtc(t_fileio *fio,
204 int natoms, int step, real time,
205 matrix box, rvec *x, real prec)
207 int magic_number = XTC_MAGIC;
208 XDR *xd;
209 gmx_bool bDum;
210 int bOK;
212 if (!fio)
214 /* This means the fio object is not being used, e.g. because
215 we are actually writing TNG output. We still have to return
216 a pseudo-success value, to keep some callers happy. */
217 return 1;
220 xd = gmx_fio_getxdr(fio);
221 /* write magic number and xtc identidier */
222 if (xtc_header(xd, &magic_number, &natoms, &step, &time, FALSE, &bDum) == 0)
224 return 0;
227 /* write data */
228 bOK = xtc_coord(xd, &natoms, box, x, &prec, FALSE); /* bOK will be 1 if writing went well */
230 if (bOK)
232 if (gmx_fio_flush(fio) != 0)
234 bOK = 0;
237 return bOK; /* 0 if bad, 1 if writing went well */
240 int read_first_xtc(t_fileio *fio, int *natoms, int *step, real *time,
241 matrix box, rvec **x, real *prec, gmx_bool *bOK)
243 int magic;
244 XDR *xd;
246 *bOK = TRUE;
247 xd = gmx_fio_getxdr(fio);
249 /* read header and malloc x */
250 if (!xtc_header(xd, &magic, natoms, step, time, TRUE, bOK))
252 return 0;
255 /* Check magic number */
256 check_xtc_magic(magic);
258 snew(*x, *natoms);
260 *bOK = xtc_coord(xd, natoms, box, *x, prec, TRUE);
262 return *bOK;
265 int read_next_xtc(t_fileio* fio,
266 int natoms, int *step, real *time,
267 matrix box, rvec *x, real *prec, gmx_bool *bOK)
269 int magic;
270 int n;
271 XDR *xd;
273 *bOK = TRUE;
274 xd = gmx_fio_getxdr(fio);
276 /* read header */
277 if (!xtc_header(xd, &magic, &n, step, time, TRUE, bOK))
279 return 0;
282 /* Check magic number */
283 check_xtc_magic(magic);
285 if (n > natoms)
287 gmx_fatal(FARGS, "Frame contains more atoms (%d) than expected (%d)",
288 n, natoms);
291 *bOK = xtc_coord(xd, &natoms, box, x, prec, TRUE);
293 return *bOK;