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