50b3282af2e8d17fd2195aab1564251b27e338fc
[alexxy/gromacs.git] / src / gromacs / fileio / trrio.cpp
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,2016,2018, 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 "trrio.h"
40
41 #include <cstring>
42
43 #include "gromacs/fileio/gmxfio.h"
44 #include "gromacs/fileio/gmxfio-xdr.h"
45 #include "gromacs/mdtypes/md_enums.h"
46 #include "gromacs/utility/fatalerror.h"
47 #include "gromacs/utility/futil.h"
48 #include "gromacs/utility/smalloc.h"
49
50 #define BUFSIZE     128
51
52 static int nFloatSize(gmx_trr_header_t *sh)
53 {
54     int nflsize = 0;
55
56     if (sh->box_size)
57     {
58         nflsize = sh->box_size/(DIM*DIM);
59     }
60     else if (sh->x_size)
61     {
62         nflsize = sh->x_size/(sh->natoms*DIM);
63     }
64     else if (sh->v_size)
65     {
66         nflsize = sh->v_size/(sh->natoms*DIM);
67     }
68     else if (sh->f_size)
69     {
70         nflsize = sh->f_size/(sh->natoms*DIM);
71     }
72     else
73     {
74         gmx_file("Can not determine precision of trr file");
75     }
76
77     if (((nflsize != sizeof(float)) && (nflsize != sizeof(double))))
78     {
79         gmx_fatal(FARGS, "Float size %d. Maybe different CPU?", nflsize);
80     }
81
82     return nflsize;
83 }
84
85 /* Returns whether a valid frame header was read. Upon exit, *bOK is
86    TRUE if a normal outcome resulted. Usually that is the same thing,
87    but encountering the end of the file before reading the magic
88    integer is a normal outcome for TRR reading, and it does not
89    produce a valid frame header, so the values differ in that case.
90    That does not exclude the possibility of a reading error between
91    frames, but the trajectory-handling infrastructure needs an
92    overhaul before we can handle that. */
93 static gmx_bool
94 do_trr_frame_header(t_fileio *fio, bool bRead, gmx_trr_header_t *sh, gmx_bool *bOK)
95 {
96     const int       magicValue = 1993;
97     int             magic      = magicValue;
98     static gmx_bool bFirst     = TRUE;
99     char            buf[256];
100
101     *bOK = TRUE;
102
103     if (!gmx_fio_do_int(fio, magic))
104     {
105         /* Failed to read an integer, which should be the magic
106            number, which usually means we've reached the end
107            of the file (but could be an I/O error that we now
108            might mishandle). */
109         return FALSE;
110     }
111     if (magic != magicValue)
112     {
113         *bOK = FALSE;
114         gmx_fatal(FARGS, "Failed to find GROMACS magic number in trr frame header, so this is not a trr file!\n");
115     }
116
117     if (bRead)
118     {
119         *bOK = *bOK && gmx_fio_do_string(fio, buf);
120         if (bFirst)
121         {
122             fprintf(stderr, "trr version: %s ", buf);
123         }
124     }
125     else
126     {
127         sprintf(buf, "GMX_trn_file");
128         *bOK = *bOK && gmx_fio_do_string(fio, buf);
129     }
130     *bOK = *bOK && gmx_fio_do_int(fio, sh->ir_size);
131     *bOK = *bOK && gmx_fio_do_int(fio, sh->e_size);
132     *bOK = *bOK && gmx_fio_do_int(fio, sh->box_size);
133     *bOK = *bOK && gmx_fio_do_int(fio, sh->vir_size);
134     *bOK = *bOK && gmx_fio_do_int(fio, sh->pres_size);
135     *bOK = *bOK && gmx_fio_do_int(fio, sh->top_size);
136     *bOK = *bOK && gmx_fio_do_int(fio, sh->sym_size);
137     *bOK = *bOK && gmx_fio_do_int(fio, sh->x_size);
138     *bOK = *bOK && gmx_fio_do_int(fio, sh->v_size);
139     *bOK = *bOK && gmx_fio_do_int(fio, sh->f_size);
140     *bOK = *bOK && gmx_fio_do_int(fio, sh->natoms);
141
142     if (!*bOK)
143     {
144         return *bOK;
145     }
146     sh->bDouble = (nFloatSize(sh) == sizeof(double));
147     gmx_fio_setprecision(fio, sh->bDouble);
148
149     if (bRead && bFirst)
150     {
151         fprintf(stderr, "(%s precision)\n", sh->bDouble ? "double" : "single");
152         bFirst = FALSE;
153     }
154
155     /* Note that TRR wasn't defined to be extensible, so we can't fix
156      * the fact that we used a default int for the step number, which
157      * is typically defined to be signed and 32 bit. */
158     int intStep = sh->step;
159     *bOK     = *bOK && gmx_fio_do_int(fio, intStep);
160     sh->step = intStep;
161     *bOK     = *bOK && gmx_fio_do_int(fio, sh->nre);
162     *bOK     = *bOK && gmx_fio_do_real(fio, sh->t);
163     *bOK     = *bOK && gmx_fio_do_real(fio, sh->lambda);
164
165     return *bOK;
166 }
167
168 static gmx_bool
169 do_trr_frame_data(t_fileio *fio, gmx_trr_header_t *sh,
170                   rvec *box, rvec *x, rvec *v, rvec *f)
171 {
172     matrix   pv;
173     gmx_bool bOK;
174
175     bOK = TRUE;
176     if (sh->box_size != 0)
177     {
178         bOK = bOK && gmx_fio_ndo_rvec(fio, box, DIM);
179     }
180     if (sh->vir_size != 0)
181     {
182         bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
183     }
184     if (sh->pres_size != 0)
185     {
186         bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
187     }
188     if (sh->x_size   != 0)
189     {
190         bOK = bOK && gmx_fio_ndo_rvec(fio, x, sh->natoms);
191     }
192     if (sh->v_size   != 0)
193     {
194         bOK = bOK && gmx_fio_ndo_rvec(fio, v, sh->natoms);
195     }
196     if (sh->f_size   != 0)
197     {
198         bOK = bOK && gmx_fio_ndo_rvec(fio, f, sh->natoms);
199     }
200
201     return bOK;
202 }
203
204 static gmx_bool
205 do_trr_frame(t_fileio *fio, bool bRead, gmx_int64_t *step, real *t, real *lambda,
206              rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
207 {
208     gmx_trr_header_t *sh;
209     gmx_bool          bOK;
210
211     snew(sh, 1);
212     if (!bRead)
213     {
214         sh->box_size = (box) ? sizeof(matrix) : 0;
215         sh->x_size   = ((x) ? (*natoms*sizeof(x[0])) : 0);
216         sh->v_size   = ((v) ? (*natoms*sizeof(v[0])) : 0);
217         sh->f_size   = ((f) ? (*natoms*sizeof(f[0])) : 0);
218         sh->natoms   = *natoms;
219         sh->step     = *step;
220         sh->nre      = 0;
221         sh->t        = *t;
222         sh->lambda   = *lambda;
223     }
224     if (!do_trr_frame_header(fio, bRead, sh, &bOK))
225     {
226         return FALSE;
227     }
228     if (bRead)
229     {
230         *natoms = sh->natoms;
231         *step   = sh->step;
232         *t      = sh->t;
233         *lambda = sh->lambda;
234         if (sh->ir_size)
235         {
236             gmx_file("inputrec in trr file");
237         }
238         if (sh->e_size)
239         {
240             gmx_file("energies in trr file");
241         }
242         if (sh->top_size)
243         {
244             gmx_file("topology in trr file");
245         }
246         if (sh->sym_size)
247         {
248             gmx_file("symbol table in trr file");
249         }
250     }
251     bOK = do_trr_frame_data(fio, sh, box, x, v, f);
252
253     sfree(sh);
254
255     return bOK;
256 }
257
258 /************************************************************
259  *
260  *  The following routines are the exported ones
261  *
262  ************************************************************/
263
264 void gmx_trr_read_single_header(const char *fn, gmx_trr_header_t *header)
265 {
266     t_fileio *fio = gmx_trr_open(fn, "r");
267     gmx_bool  bOK;
268     if (!do_trr_frame_header(fio, true, header, &bOK))
269     {
270         gmx_fatal(FARGS, "Empty file %s", fn);
271     }
272     gmx_trr_close(fio);
273 }
274
275 gmx_bool gmx_trr_read_frame_header(t_fileio *fio, gmx_trr_header_t *header, gmx_bool *bOK)
276 {
277     return do_trr_frame_header(fio, true, header, bOK);
278 }
279
280 void gmx_trr_write_single_frame(const char *fn, gmx_int64_t step, real t, real lambda,
281                                 const rvec *box, int natoms, const rvec *x, const rvec *v, const rvec *f)
282 {
283     t_fileio *fio = gmx_trr_open(fn, "w");
284     do_trr_frame(fio, false, &step, &t, &lambda, const_cast<rvec *>(box), &natoms, const_cast<rvec *>(x), const_cast<rvec *>(v), const_cast<rvec *>(f));
285     gmx_trr_close(fio);
286 }
287
288 void gmx_trr_read_single_frame(const char *fn, gmx_int64_t *step, real *t, real *lambda,
289                                rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
290 {
291     t_fileio *fio = gmx_trr_open(fn, "r");
292     do_trr_frame(fio, true, step, t, lambda, box, natoms, x, v, f);
293     gmx_trr_close(fio);
294 }
295
296 void gmx_trr_write_frame(t_fileio *fio, gmx_int64_t step, real t, real lambda,
297                          const rvec *box, int natoms, const rvec *x, const rvec *v, const rvec *f)
298 {
299     if (!do_trr_frame(fio, false, &step, &t, &lambda, const_cast<rvec *>(box), &natoms, const_cast<rvec *>(x), const_cast<rvec *>(v), const_cast<rvec *>(f)))
300     {
301         gmx_file("Cannot write trajectory frame; maybe you are out of disk space?");
302     }
303 }
304
305
306 gmx_bool gmx_trr_read_frame(t_fileio *fio, gmx_int64_t *step, real *t, real *lambda,
307                             rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
308 {
309     return do_trr_frame(fio, true, step, t, lambda, box, natoms, x, v, f);
310 }
311
312 gmx_bool gmx_trr_read_frame_data(t_fileio *fio, gmx_trr_header_t *header,
313                                  rvec *box, rvec *x, rvec *v, rvec *f)
314 {
315     return do_trr_frame_data(fio, header, box, x, v, f);
316 }
317
318 t_fileio *gmx_trr_open(const char *fn, const char *mode)
319 {
320     return gmx_fio_open(fn, mode);
321 }
322
323 void gmx_trr_close(t_fileio *fio)
324 {
325     gmx_fio_close(fio);
326 }