Copy of CI from master to 2020
[gromacs.git] / src / gromacs / fileio / xtcio.cpp
blob35c77deba6b0266719caf3c71e263d644d93e719
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,2014,2015,2016,2018,2019, 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 #include "gmxpre.h"
39 #include "xtcio.h"
41 #include <cstring>
43 #include "gromacs/fileio/gmxfio.h"
44 #include "gromacs/fileio/gmxfio_xdr.h"
45 #include "gromacs/fileio/xdrf.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/utility/fatalerror.h"
48 #include "gromacs/utility/futil.h"
49 #include "gromacs/utility/smalloc.h"
51 #define XTC_MAGIC 1995
54 static int xdr_r2f(XDR* xdrs, real* r, gmx_bool gmx_unused bRead)
56 #if 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, static_cast<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)", magic, XTC_MAGIC);
95 static int xtc_check(const char* str, gmx_bool bResult, const char* file, int line)
97 if (!bResult)
99 if (debug)
101 fprintf(debug,
102 "\nXTC error: read/write of %s failed, "
103 "source file %s, line %d\n",
104 str, file, line);
106 return 0;
108 return 1;
111 #define XTC_CHECK(s, b) xtc_check(s, b, __FILE__, __LINE__)
113 static int xtc_header(XDR* xd, int* magic, int* natoms, int64_t* step, real* time, gmx_bool bRead, gmx_bool* bOK)
115 int result;
117 if (xdr_int(xd, magic) == 0)
119 return 0;
121 result = XTC_CHECK("natoms", xdr_int(xd, natoms)); /* number of atoms */
122 if (result)
124 /* Note that XTC wasn't defined to be extensible, so we can't
125 * fix the fact that we used xdr_int for the step number,
126 * which is defined to be signed and 32 bit. */
127 int intStep = *step;
128 result = XTC_CHECK("step", xdr_int(xd, &intStep)); /* frame number */
129 *step = intStep;
131 if (result)
133 result = XTC_CHECK("time", xdr_r2f(xd, time, bRead)); /* time */
135 *bOK = (result != 0);
137 return result;
140 static int xtc_coord(XDR* xd, int* natoms, rvec* box, rvec* x, real* prec, gmx_bool bRead)
142 int i, j, result;
143 #if GMX_DOUBLE
144 float* ftmp;
145 float fprec;
146 #endif
148 /* box */
149 result = 1;
150 for (i = 0; ((i < DIM) && result); i++)
152 for (j = 0; ((j < DIM) && result); j++)
154 result = XTC_CHECK("box", xdr_r2f(xd, &(box[i][j]), bRead));
158 if (!result)
160 return result;
163 #if GMX_DOUBLE
164 /* allocate temp. single-precision array */
165 snew(ftmp, (*natoms) * DIM);
167 /* Copy data to temp. array if writing */
168 if (!bRead)
170 for (i = 0; (i < *natoms); i++)
172 ftmp[DIM * i + XX] = x[i][XX];
173 ftmp[DIM * i + YY] = x[i][YY];
174 ftmp[DIM * i + ZZ] = x[i][ZZ];
176 fprec = *prec;
178 result = XTC_CHECK("x", xdr3dfcoord(xd, ftmp, natoms, &fprec));
180 /* Copy from temp. array if reading */
181 if (bRead)
183 for (i = 0; (i < *natoms); i++)
185 x[i][XX] = ftmp[DIM * i + XX];
186 x[i][YY] = ftmp[DIM * i + YY];
187 x[i][ZZ] = ftmp[DIM * i + ZZ];
189 *prec = fprec;
191 sfree(ftmp);
192 #else
193 result = XTC_CHECK("x", xdr3dfcoord(xd, x[0], natoms, prec));
194 #endif
196 return result;
200 int write_xtc(t_fileio* fio, int natoms, int64_t step, real time, const rvec* box, const rvec* x, real prec)
202 int magic_number = XTC_MAGIC;
203 XDR* xd;
204 gmx_bool bDum;
205 int bOK;
207 if (!fio)
209 /* This means the fio object is not being used, e.g. because
210 we are actually writing TNG output. We still have to return
211 a pseudo-success value, to keep some callers happy. */
212 return 1;
215 xd = gmx_fio_getxdr(fio);
216 /* write magic number and xtc identidier */
217 if (xtc_header(xd, &magic_number, &natoms, &step, &time, FALSE, &bDum) == 0)
219 return 0;
222 /* write data */
223 bOK = xtc_coord(xd, &natoms, const_cast<rvec*>(box), const_cast<rvec*>(x), &prec,
224 FALSE); /* bOK will be 1 if writing went well */
226 if (bOK)
228 if (gmx_fio_flush(fio) != 0)
230 bOK = 0;
233 return bOK; /* 0 if bad, 1 if writing went well */
236 int read_first_xtc(t_fileio* fio, int* natoms, int64_t* step, real* time, matrix box, rvec** x, real* prec, gmx_bool* bOK)
238 int magic;
239 XDR* xd;
241 *bOK = TRUE;
242 xd = gmx_fio_getxdr(fio);
244 /* read header and malloc x */
245 if (!xtc_header(xd, &magic, natoms, step, time, TRUE, bOK))
247 return 0;
250 /* Check magic number */
251 check_xtc_magic(magic);
253 snew(*x, *natoms);
255 *bOK = (xtc_coord(xd, natoms, box, *x, prec, TRUE) != 0);
257 return static_cast<int>(*bOK);
260 int read_next_xtc(t_fileio* fio, int natoms, int64_t* step, real* time, matrix box, rvec* x, real* prec, gmx_bool* bOK)
262 int magic;
263 int n;
264 XDR* xd;
266 *bOK = TRUE;
267 xd = gmx_fio_getxdr(fio);
269 /* read header */
270 if (!xtc_header(xd, &magic, &n, step, time, TRUE, bOK))
272 return 0;
275 /* Check magic number */
276 check_xtc_magic(magic);
278 if (n > natoms)
280 gmx_fatal(FARGS, "Frame contains more atoms (%d) than expected (%d)", n, natoms);
283 *bOK = (xtc_coord(xd, &natoms, box, x, prec, TRUE) != 0);
285 return static_cast<int>(*bOK);