e9f227e4caa654c0f30b603de59628505c8e8f95
[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.
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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #include "gmxpre.h"
39
40 #include "trrio.h"
41
42 #include <cstring>
43
44 #include "gromacs/fileio/gmxfio.h"
45 #include "gromacs/fileio/gmxfio_xdr.h"
46 #include "gromacs/mdtypes/md_enums.h"
47 #include "gromacs/utility/fatalerror.h"
48 #include "gromacs/utility/futil.h"
49 #include "gromacs/utility/smalloc.h"
50
51 #define BUFSIZE 128
52
53 static int nFloatSize(gmx_trr_header_t* sh)
54 {
55     int nflsize = 0;
56
57     if (sh->box_size)
58     {
59         nflsize = sh->box_size / (DIM * DIM);
60     }
61     else if (sh->x_size)
62     {
63         nflsize = sh->x_size / (sh->natoms * DIM);
64     }
65     else if (sh->v_size)
66     {
67         nflsize = sh->v_size / (sh->natoms * DIM);
68     }
69     else if (sh->f_size)
70     {
71         nflsize = sh->f_size / (sh->natoms * DIM);
72     }
73     else
74     {
75         gmx_file("Can not determine precision of trr file");
76     }
77
78     if (((nflsize != sizeof(float)) && (nflsize != sizeof(double))))
79     {
80         gmx_fatal(FARGS, "Float size %d. Maybe different CPU?", nflsize);
81     }
82
83     return nflsize;
84 }
85
86 /* Returns whether a valid frame header was read. Upon exit, *bOK is
87    TRUE if a normal outcome resulted. Usually that is the same thing,
88    but encountering the end of the file before reading the magic
89    integer is a normal outcome for TRR reading, and it does not
90    produce a valid frame header, so the values differ in that case.
91    That does not exclude the possibility of a reading error between
92    frames, but the trajectory-handling infrastructure needs an
93    overhaul before we can handle that. */
94 static gmx_bool 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,
115                   "Failed to find GROMACS magic number in trr frame header, so this is not a trr "
116                   "file!\n");
117     }
118
119     if (bRead)
120     {
121         *bOK = *bOK && gmx_fio_do_string(fio, buf);
122         if (bFirst)
123         {
124             fprintf(stderr, "trr version: %s ", buf);
125         }
126     }
127     else
128     {
129         sprintf(buf, "GMX_trn_file");
130         *bOK = *bOK && gmx_fio_do_string(fio, buf);
131     }
132     *bOK = *bOK && gmx_fio_do_int(fio, sh->ir_size);
133     *bOK = *bOK && gmx_fio_do_int(fio, sh->e_size);
134     *bOK = *bOK && gmx_fio_do_int(fio, sh->box_size);
135     *bOK = *bOK && gmx_fio_do_int(fio, sh->vir_size);
136     *bOK = *bOK && gmx_fio_do_int(fio, sh->pres_size);
137     *bOK = *bOK && gmx_fio_do_int(fio, sh->top_size);
138     *bOK = *bOK && gmx_fio_do_int(fio, sh->sym_size);
139     *bOK = *bOK && gmx_fio_do_int(fio, sh->x_size);
140     *bOK = *bOK && gmx_fio_do_int(fio, sh->v_size);
141     *bOK = *bOK && gmx_fio_do_int(fio, sh->f_size);
142     *bOK = *bOK && gmx_fio_do_int(fio, sh->natoms);
143
144     if (!*bOK)
145     {
146         return *bOK;
147     }
148     sh->bDouble = (nFloatSize(sh) == sizeof(double));
149     gmx_fio_setprecision(fio, sh->bDouble);
150
151     if (bRead && bFirst)
152     {
153         fprintf(stderr, "(%s precision)\n", sh->bDouble ? "double" : "single");
154         bFirst = FALSE;
155     }
156
157     /* Note that TRR wasn't defined to be extensible, so we can't fix
158      * the fact that we used a default int for the step number, which
159      * is typically defined to be signed and 32 bit. */
160     int intStep = sh->step;
161     *bOK        = *bOK && gmx_fio_do_int(fio, intStep);
162     sh->step    = intStep;
163     *bOK        = *bOK && gmx_fio_do_int(fio, sh->nre);
164     *bOK        = *bOK && gmx_fio_do_real(fio, sh->t);
165     *bOK        = *bOK && gmx_fio_do_real(fio, sh->lambda);
166
167     return *bOK;
168 }
169
170 static gmx_bool do_trr_frame_data(t_fileio* fio, gmx_trr_header_t* sh, 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 do_trr_frame(t_fileio* fio,
205                              bool      bRead,
206                              int64_t*  step,
207                              real*     t,
208                              real*     lambda,
209                              rvec*     box,
210                              int*      natoms,
211                              rvec*     x,
212                              rvec*     v,
213                              rvec*     f)
214 {
215     gmx_trr_header_t* sh;
216     gmx_bool          bOK;
217
218     snew(sh, 1);
219     if (!bRead)
220     {
221         sh->box_size = (box) ? sizeof(matrix) : 0;
222         sh->x_size   = ((x) ? (*natoms * sizeof(x[0])) : 0);
223         sh->v_size   = ((v) ? (*natoms * sizeof(v[0])) : 0);
224         sh->f_size   = ((f) ? (*natoms * sizeof(f[0])) : 0);
225         sh->natoms   = *natoms;
226         sh->step     = *step;
227         sh->nre      = 0;
228         sh->t        = *t;
229         sh->lambda   = *lambda;
230     }
231     if (!do_trr_frame_header(fio, bRead, sh, &bOK))
232     {
233         return FALSE;
234     }
235     if (bRead)
236     {
237         *natoms = sh->natoms;
238         *step   = sh->step;
239         *t      = sh->t;
240         *lambda = sh->lambda;
241         if (sh->ir_size)
242         {
243             gmx_file("inputrec in trr file");
244         }
245         if (sh->e_size)
246         {
247             gmx_file("energies in trr file");
248         }
249         if (sh->top_size)
250         {
251             gmx_file("topology in trr file");
252         }
253         if (sh->sym_size)
254         {
255             gmx_file("symbol table in trr file");
256         }
257     }
258     bOK = do_trr_frame_data(fio, sh, box, x, v, f);
259
260     sfree(sh);
261
262     return bOK;
263 }
264
265 /************************************************************
266  *
267  *  The following routines are the exported ones
268  *
269  ************************************************************/
270
271 void gmx_trr_read_single_header(const char* fn, gmx_trr_header_t* header)
272 {
273     t_fileio* fio = gmx_trr_open(fn, "r");
274     gmx_bool  bOK;
275     if (!do_trr_frame_header(fio, true, header, &bOK))
276     {
277         gmx_fatal(FARGS, "Empty file %s", fn);
278     }
279     gmx_trr_close(fio);
280 }
281
282 gmx_bool gmx_trr_read_frame_header(t_fileio* fio, gmx_trr_header_t* header, gmx_bool* bOK)
283 {
284     return do_trr_frame_header(fio, true, header, bOK);
285 }
286
287 void gmx_trr_write_single_frame(const char* fn,
288                                 int64_t     step,
289                                 real        t,
290                                 real        lambda,
291                                 const rvec* box,
292                                 int         natoms,
293                                 const rvec* x,
294                                 const rvec* v,
295                                 const rvec* f)
296 {
297     t_fileio* fio = gmx_trr_open(fn, "w");
298     do_trr_frame(fio, false, &step, &t, &lambda, const_cast<rvec*>(box), &natoms,
299                  const_cast<rvec*>(x), const_cast<rvec*>(v), const_cast<rvec*>(f));
300     gmx_trr_close(fio);
301 }
302
303 void gmx_trr_read_single_frame(const char* fn,
304                                int64_t*    step,
305                                real*       t,
306                                real*       lambda,
307                                rvec*       box,
308                                int*        natoms,
309                                rvec*       x,
310                                rvec*       v,
311                                rvec*       f)
312 {
313     t_fileio* fio = gmx_trr_open(fn, "r");
314     do_trr_frame(fio, true, step, t, lambda, box, natoms, x, v, f);
315     gmx_trr_close(fio);
316 }
317
318 void gmx_trr_write_frame(t_fileio*   fio,
319                          int64_t     step,
320                          real        t,
321                          real        lambda,
322                          const rvec* box,
323                          int         natoms,
324                          const rvec* x,
325                          const rvec* v,
326                          const rvec* f)
327 {
328     if (!do_trr_frame(fio, false, &step, &t, &lambda, const_cast<rvec*>(box), &natoms,
329                       const_cast<rvec*>(x), const_cast<rvec*>(v), const_cast<rvec*>(f)))
330     {
331         gmx_file("Cannot write trajectory frame; maybe you are out of disk space?");
332     }
333 }
334
335
336 gmx_bool gmx_trr_read_frame(t_fileio* fio,
337                             int64_t*  step,
338                             real*     t,
339                             real*     lambda,
340                             rvec*     box,
341                             int*      natoms,
342                             rvec*     x,
343                             rvec*     v,
344                             rvec*     f)
345 {
346     return do_trr_frame(fio, true, step, t, lambda, box, natoms, x, v, f);
347 }
348
349 gmx_bool gmx_trr_read_frame_data(t_fileio* fio, gmx_trr_header_t* header, rvec* box, rvec* x, rvec* v, rvec* f)
350 {
351     return do_trr_frame_data(fio, header, box, x, v, f);
352 }
353
354 t_fileio* gmx_trr_open(const char* fn, const char* mode)
355 {
356     return gmx_fio_open(fn, mode);
357 }
358
359 void gmx_trr_close(t_fileio* fio)
360 {
361     gmx_fio_close(fio);
362 }