Separate job script files for gmxapi package versions.
[gromacs.git] / src / gromacs / fileio / xtcio.cpp
blobc5ceb8338d1ff16430859bbe6023a54390d82c95
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 by the GROMACS development team.
7 * Copyright (c) 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.
38 #include "gmxpre.h"
40 #include "xtcio.h"
42 #include <cstring>
44 #include "gromacs/fileio/gmxfio.h"
45 #include "gromacs/fileio/gmxfio_xdr.h"
46 #include "gromacs/fileio/xdrf.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/utility/fatalerror.h"
49 #include "gromacs/utility/futil.h"
50 #include "gromacs/utility/smalloc.h"
52 #define XTC_MAGIC 1995
55 static int xdr_r2f(XDR* xdrs, real* r, gmx_bool gmx_unused bRead)
57 #if GMX_DOUBLE
58 float f;
59 int ret;
61 if (!bRead)
63 f = *r;
65 ret = xdr_float(xdrs, &f);
66 if (bRead)
68 *r = f;
71 return ret;
72 #else
73 return xdr_float(xdrs, static_cast<float*>(r));
74 #endif
78 t_fileio* open_xtc(const char* fn, const char* mode)
80 return gmx_fio_open(fn, mode);
83 void close_xtc(t_fileio* fio)
85 gmx_fio_close(fio);
88 static void check_xtc_magic(int magic)
90 if (magic != XTC_MAGIC)
92 gmx_fatal(FARGS, "Magic Number Error in XTC file (read %d, should be %d)", magic, XTC_MAGIC);
96 static int xtc_check(const char* str, gmx_bool bResult, const char* file, int line)
98 if (!bResult)
100 if (debug)
102 fprintf(debug,
103 "\nXTC error: read/write of %s failed, "
104 "source file %s, line %d\n",
105 str, file, line);
107 return 0;
109 return 1;
112 #define XTC_CHECK(s, b) xtc_check(s, b, __FILE__, __LINE__)
114 static int xtc_header(XDR* xd, int* magic, int* natoms, int64_t* step, real* time, gmx_bool bRead, gmx_bool* bOK)
116 int result;
118 if (xdr_int(xd, magic) == 0)
120 return 0;
122 result = XTC_CHECK("natoms", xdr_int(xd, natoms)); /* number of atoms */
123 if (result)
125 /* Note that XTC wasn't defined to be extensible, so we can't
126 * fix the fact that we used xdr_int for the step number,
127 * which is defined to be signed and 32 bit. */
128 int intStep = *step;
129 result = XTC_CHECK("step", xdr_int(xd, &intStep)); /* frame number */
130 *step = intStep;
132 if (result)
134 result = XTC_CHECK("time", xdr_r2f(xd, time, bRead)); /* time */
136 *bOK = (result != 0);
138 return result;
141 static int xtc_coord(XDR* xd, int* natoms, rvec* box, rvec* x, real* prec, gmx_bool bRead)
143 int i, j, result;
144 #if GMX_DOUBLE
145 float* ftmp;
146 float fprec;
147 #endif
149 /* box */
150 result = 1;
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));
159 if (!result)
161 return result;
164 #if GMX_DOUBLE
165 /* allocate temp. single-precision array */
166 snew(ftmp, (*natoms) * DIM);
168 /* Copy data to temp. array if writing */
169 if (!bRead)
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];
177 fprec = *prec;
179 result = XTC_CHECK("x", xdr3dfcoord(xd, ftmp, natoms, &fprec));
181 /* Copy from temp. array if reading */
182 if (bRead)
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];
190 *prec = fprec;
192 sfree(ftmp);
193 #else
194 result = XTC_CHECK("x", xdr3dfcoord(xd, x[0], natoms, prec));
195 #endif
197 return result;
201 int write_xtc(t_fileio* fio, int natoms, int64_t step, real time, const rvec* box, const rvec* x, real prec)
203 int magic_number = XTC_MAGIC;
204 XDR* xd;
205 gmx_bool bDum;
206 int bOK;
208 if (!fio)
210 /* This means the fio object is not being used, e.g. because
211 we are actually writing TNG output. We still have to return
212 a pseudo-success value, to keep some callers happy. */
213 return 1;
216 xd = gmx_fio_getxdr(fio);
217 /* write magic number and xtc identidier */
218 if (xtc_header(xd, &magic_number, &natoms, &step, &time, FALSE, &bDum) == 0)
220 return 0;
223 /* write data */
224 bOK = xtc_coord(xd, &natoms, const_cast<rvec*>(box), const_cast<rvec*>(x), &prec,
225 FALSE); /* bOK will be 1 if writing went well */
227 if (bOK)
229 if (gmx_fio_flush(fio) != 0)
231 bOK = 0;
234 return bOK; /* 0 if bad, 1 if writing went well */
237 int read_first_xtc(t_fileio* fio, int* natoms, int64_t* step, real* time, matrix box, rvec** x, real* prec, gmx_bool* bOK)
239 int magic;
240 XDR* xd;
242 *bOK = TRUE;
243 xd = gmx_fio_getxdr(fio);
245 /* read header and malloc x */
246 if (!xtc_header(xd, &magic, natoms, step, time, TRUE, bOK))
248 return 0;
251 /* Check magic number */
252 check_xtc_magic(magic);
254 snew(*x, *natoms);
256 *bOK = (xtc_coord(xd, natoms, box, *x, prec, TRUE) != 0);
258 return static_cast<int>(*bOK);
261 int read_next_xtc(t_fileio* fio, int natoms, int64_t* step, real* time, matrix box, rvec* x, real* prec, gmx_bool* bOK)
263 int magic;
264 int n;
265 XDR* xd;
267 *bOK = TRUE;
268 xd = gmx_fio_getxdr(fio);
270 /* read header */
271 if (!xtc_header(xd, &magic, &n, step, time, TRUE, bOK))
273 return 0;
276 /* Check magic number */
277 check_xtc_magic(magic);
279 if (n > natoms)
281 gmx_fatal(FARGS, "Frame contains more atoms (%d) than expected (%d)", n, natoms);
284 *bOK = (xtc_coord(xd, &natoms, box, x, prec, TRUE) != 0);
286 return static_cast<int>(*bOK);