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