a477f568d5658f2769fd52b2fa4a40e498eca232
[alexxy/gromacs.git] / src / gromacs / fileio / xtcio.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "xtcio.h"
40
41 #include <string.h>
42
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"
50
51 #define XTC_MAGIC 1995
52
53
54 static int xdr_r2f(XDR *xdrs, real *r, gmx_bool gmx_unused bRead)
55 {
56 #ifdef GMX_DOUBLE
57     float f;
58     int   ret;
59
60     if (!bRead)
61     {
62         f = *r;
63     }
64     ret = xdr_float(xdrs, &f);
65     if (bRead)
66     {
67         *r = f;
68     }
69
70     return ret;
71 #else
72     return xdr_float(xdrs, (float *)r);
73 #endif
74 }
75
76
77 t_fileio *open_xtc(const char *fn, const char *mode)
78 {
79     return gmx_fio_open(fn, mode);
80 }
81
82 void close_xtc(t_fileio *fio)
83 {
84     gmx_fio_close(fio);
85 }
86
87 static void check_xtc_magic(int magic)
88 {
89     if (magic != XTC_MAGIC)
90     {
91         gmx_fatal(FARGS, "Magic Number Error in XTC file (read %d, should be %d)",
92                   magic, XTC_MAGIC);
93     }
94 }
95
96 static int xtc_check(const char *str, gmx_bool bResult, const char *file, int line)
97 {
98     if (!bResult)
99     {
100         if (debug)
101         {
102             fprintf(debug, "\nXTC error: read/write of %s failed, "
103                     "source file %s, line %d\n", str, file, line);
104         }
105         return 0;
106     }
107     return 1;
108 }
109
110 #define XTC_CHECK(s, b) xtc_check(s, b, __FILE__, __LINE__)
111
112 static int xtc_header(XDR *xd, int *magic, int *natoms, int *step, real *time,
113                       gmx_bool bRead, gmx_bool *bOK)
114 {
115     int result;
116
117     if (xdr_int(xd, magic) == 0)
118     {
119         return 0;
120     }
121     result = XTC_CHECK("natoms", xdr_int(xd, natoms)); /* number of atoms */
122     if (result)
123     {
124         result = XTC_CHECK("step",   xdr_int(xd, step)); /* frame number    */
125     }
126     if (result)
127     {
128         result = XTC_CHECK("time",   xdr_r2f(xd, time, bRead)); /* time */
129     }
130     *bOK = (result != 0);
131
132     return result;
133 }
134
135 static int xtc_coord(XDR *xd, int *natoms, matrix box, rvec *x, real *prec, gmx_bool bRead)
136 {
137     int    i, j, result;
138 #ifdef GMX_DOUBLE
139     float *ftmp;
140     float  fprec;
141 #endif
142
143     /* box */
144     result = 1;
145     for (i = 0; ((i < DIM) && result); i++)
146     {
147         for (j = 0; ((j < DIM) && result); j++)
148         {
149             result = XTC_CHECK("box", xdr_r2f(xd, &(box[i][j]), bRead));
150         }
151     }
152
153     if (!result)
154     {
155         return result;
156     }
157
158 #ifdef GMX_DOUBLE
159     /* allocate temp. single-precision array */
160     snew(ftmp, (*natoms)*DIM);
161
162     /* Copy data to temp. array if writing */
163     if (!bRead)
164     {
165         for (i = 0; (i < *natoms); i++)
166         {
167             ftmp[DIM*i+XX] = x[i][XX];
168             ftmp[DIM*i+YY] = x[i][YY];
169             ftmp[DIM*i+ZZ] = x[i][ZZ];
170         }
171         fprec = *prec;
172     }
173     result = XTC_CHECK("x", xdr3dfcoord(xd, ftmp, natoms, &fprec));
174
175     /* Copy from temp. array if reading */
176     if (bRead)
177     {
178         for (i = 0; (i < *natoms); i++)
179         {
180             x[i][XX] = ftmp[DIM*i+XX];
181             x[i][YY] = ftmp[DIM*i+YY];
182             x[i][ZZ] = ftmp[DIM*i+ZZ];
183         }
184         *prec = fprec;
185     }
186     sfree(ftmp);
187 #else
188     result = XTC_CHECK("x", xdr3dfcoord(xd, x[0], natoms, prec));
189 #endif
190
191     return result;
192 }
193
194
195
196 int write_xtc(t_fileio *fio,
197               int natoms, int step, real time,
198               matrix box, rvec *x, real prec)
199 {
200     int      magic_number = XTC_MAGIC;
201     XDR     *xd;
202     gmx_bool bDum;
203     int      bOK;
204
205     if (!fio)
206     {
207         /* This means the fio object is not being used, e.g. because
208            we are actually writing TNG output. We still have to return
209            a pseudo-success value, to keep some callers happy. */
210         return 1;
211     }
212
213     xd = gmx_fio_getxdr(fio);
214     /* write magic number and xtc identidier */
215     if (xtc_header(xd, &magic_number, &natoms, &step, &time, FALSE, &bDum) == 0)
216     {
217         return 0;
218     }
219
220     /* write data */
221     bOK = xtc_coord(xd, &natoms, box, x, &prec, FALSE); /* bOK will be 1 if writing went well */
222
223     if (bOK)
224     {
225         if (gmx_fio_flush(fio) != 0)
226         {
227             bOK = 0;
228         }
229     }
230     return bOK; /* 0 if bad, 1 if writing went well */
231 }
232
233 int read_first_xtc(t_fileio *fio, int *natoms, int *step, real *time,
234                    matrix box, rvec **x, real *prec, gmx_bool *bOK)
235 {
236     int  magic;
237     XDR *xd;
238
239     *bOK = TRUE;
240     xd   = gmx_fio_getxdr(fio);
241
242     /* read header and malloc x */
243     if (!xtc_header(xd, &magic, natoms, step, time, TRUE, bOK))
244     {
245         return 0;
246     }
247
248     /* Check magic number */
249     check_xtc_magic(magic);
250
251     snew(*x, *natoms);
252
253     *bOK = xtc_coord(xd, natoms, box, *x, prec, TRUE);
254
255     return *bOK;
256 }
257
258 int read_next_xtc(t_fileio* fio,
259                   int natoms, int *step, real *time,
260                   matrix box, rvec *x, real *prec, gmx_bool *bOK)
261 {
262     int  magic;
263     int  n;
264     XDR *xd;
265
266     *bOK = TRUE;
267     xd   = gmx_fio_getxdr(fio);
268
269     /* read header */
270     if (!xtc_header(xd, &magic, &n, step, time, TRUE, bOK))
271     {
272         return 0;
273     }
274
275     /* Check magic number */
276     check_xtc_magic(magic);
277
278     if (n > natoms)
279     {
280         gmx_fatal(FARGS, "Frame contains more atoms (%d) than expected (%d)",
281                   n, natoms);
282     }
283
284     *bOK = xtc_coord(xd, &natoms, box, x, prec, TRUE);
285
286     return *bOK;
287 }