SYCL: Avoid using no_init read accessor in rocFFT
[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,2021, 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 static int nFloatSize(gmx_trr_header_t* sh)
52 {
53     int nflsize = 0;
54
55     if (sh->box_size)
56     {
57         nflsize = sh->box_size / (DIM * DIM);
58     }
59     else if (sh->x_size)
60     {
61         nflsize = sh->x_size / (sh->natoms * DIM);
62     }
63     else if (sh->v_size)
64     {
65         nflsize = sh->v_size / (sh->natoms * DIM);
66     }
67     else if (sh->f_size)
68     {
69         nflsize = sh->f_size / (sh->natoms * DIM);
70     }
71     else
72     {
73         gmx_file("Can not determine precision of trr file");
74     }
75
76     if (((nflsize != sizeof(float)) && (nflsize != sizeof(double))))
77     {
78         gmx_fatal(FARGS, "Float size %d. Maybe different CPU?", nflsize);
79     }
80
81     return nflsize;
82 }
83
84 /* Returns whether a valid frame header was read. Upon exit, *bOK is
85    TRUE if a normal outcome resulted. Usually that is the same thing,
86    but encountering the end of the file before reading the magic
87    integer is a normal outcome for TRR reading, and it does not
88    produce a valid frame header, so the values differ in that case.
89    That does not exclude the possibility of a reading error between
90    frames, but the trajectory-handling infrastructure needs an
91    overhaul before we can handle that. */
92 static gmx_bool do_trr_frame_header(t_fileio* fio, bool bRead, gmx_trr_header_t* sh, gmx_bool* bOK)
93 {
94     const int       magicValue = 1993;
95     int             magic      = magicValue;
96     static gmx_bool bFirst     = TRUE;
97     char            buf[256];
98
99     *bOK = TRUE;
100
101     if (!gmx_fio_do_int(fio, magic))
102     {
103         /* Failed to read an integer, which should be the magic
104            number, which usually means we've reached the end
105            of the file (but could be an I/O error that we now
106            might mishandle). */
107         return FALSE;
108     }
109     if (magic != magicValue)
110     {
111         *bOK = FALSE;
112         gmx_fatal(FARGS,
113                   "Failed to find GROMACS magic number in trr frame header, so this is not a trr "
114                   "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 do_trr_frame_data(t_fileio* fio, gmx_trr_header_t* sh, rvec* box, rvec* x, rvec* v, rvec* f)
169 {
170     matrix   pv;
171     gmx_bool bOK;
172
173     bOK = TRUE;
174     if (sh->box_size != 0)
175     {
176         bOK = bOK && gmx_fio_ndo_rvec(fio, box, DIM);
177     }
178     if (sh->vir_size != 0)
179     {
180         bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
181     }
182     if (sh->pres_size != 0)
183     {
184         bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
185     }
186     if (sh->x_size != 0)
187     {
188         bOK = bOK && gmx_fio_ndo_rvec(fio, x, sh->natoms);
189     }
190     if (sh->v_size != 0)
191     {
192         bOK = bOK && gmx_fio_ndo_rvec(fio, v, sh->natoms);
193     }
194     if (sh->f_size != 0)
195     {
196         bOK = bOK && gmx_fio_ndo_rvec(fio, f, sh->natoms);
197     }
198
199     return bOK;
200 }
201
202 static gmx_bool do_trr_frame(t_fileio* fio,
203                              bool      bRead,
204                              int64_t*  step,
205                              real*     t,
206                              real*     lambda,
207                              rvec*     box,
208                              int*      natoms,
209                              rvec*     x,
210                              rvec*     v,
211                              rvec*     f)
212 {
213     gmx_trr_header_t* sh;
214     gmx_bool          bOK;
215
216     snew(sh, 1);
217     if (!bRead)
218     {
219         sh->box_size = (box) ? sizeof(matrix) : 0;
220         sh->x_size   = ((x) ? (*natoms * sizeof(x[0])) : 0);
221         sh->v_size   = ((v) ? (*natoms * sizeof(v[0])) : 0);
222         sh->f_size   = ((f) ? (*natoms * sizeof(f[0])) : 0);
223         sh->natoms   = *natoms;
224         sh->step     = *step;
225         sh->nre      = 0;
226         sh->t        = *t;
227         sh->lambda   = *lambda;
228     }
229     if (!do_trr_frame_header(fio, bRead, sh, &bOK))
230     {
231         return FALSE;
232     }
233     if (bRead)
234     {
235         *natoms = sh->natoms;
236         *step   = sh->step;
237         *t      = sh->t;
238         *lambda = sh->lambda;
239         if (sh->ir_size)
240         {
241             gmx_file("inputrec in trr file");
242         }
243         if (sh->e_size)
244         {
245             gmx_file("energies in trr file");
246         }
247         if (sh->top_size)
248         {
249             gmx_file("topology in trr file");
250         }
251         if (sh->sym_size)
252         {
253             gmx_file("symbol table in trr file");
254         }
255     }
256     bOK = do_trr_frame_data(fio, sh, box, x, v, f);
257
258     sfree(sh);
259
260     return bOK;
261 }
262
263 /************************************************************
264  *
265  *  The following routines are the exported ones
266  *
267  ************************************************************/
268
269 void gmx_trr_read_single_header(const char* fn, gmx_trr_header_t* header)
270 {
271     t_fileio* fio = gmx_trr_open(fn, "r");
272     gmx_bool  bOK;
273     if (!do_trr_frame_header(fio, true, header, &bOK))
274     {
275         gmx_fatal(FARGS, "Empty file %s", fn);
276     }
277     gmx_trr_close(fio);
278 }
279
280 gmx_bool gmx_trr_read_frame_header(t_fileio* fio, gmx_trr_header_t* header, gmx_bool* bOK)
281 {
282     return do_trr_frame_header(fio, true, header, bOK);
283 }
284
285 void gmx_trr_write_single_frame(const char* fn,
286                                 int64_t     step,
287                                 real        t,
288                                 real        lambda,
289                                 const rvec* box,
290                                 int         natoms,
291                                 const rvec* x,
292                                 const rvec* v,
293                                 const rvec* f)
294 {
295     t_fileio* fio = gmx_trr_open(fn, "w");
296     do_trr_frame(fio,
297                  false,
298                  &step,
299                  &t,
300                  &lambda,
301                  const_cast<rvec*>(box),
302                  &natoms,
303                  const_cast<rvec*>(x),
304                  const_cast<rvec*>(v),
305                  const_cast<rvec*>(f));
306     gmx_trr_close(fio);
307 }
308
309 void gmx_trr_read_single_frame(const char* fn,
310                                int64_t*    step,
311                                real*       t,
312                                real*       lambda,
313                                rvec*       box,
314                                int*        natoms,
315                                rvec*       x,
316                                rvec*       v,
317                                rvec*       f)
318 {
319     t_fileio* fio = gmx_trr_open(fn, "r");
320     do_trr_frame(fio, true, step, t, lambda, box, natoms, x, v, f);
321     gmx_trr_close(fio);
322 }
323
324 void gmx_trr_write_frame(t_fileio*   fio,
325                          int64_t     step,
326                          real        t,
327                          real        lambda,
328                          const rvec* box,
329                          int         natoms,
330                          const rvec* x,
331                          const rvec* v,
332                          const rvec* f)
333 {
334     if (!do_trr_frame(fio,
335                       false,
336                       &step,
337                       &t,
338                       &lambda,
339                       const_cast<rvec*>(box),
340                       &natoms,
341                       const_cast<rvec*>(x),
342                       const_cast<rvec*>(v),
343                       const_cast<rvec*>(f)))
344     {
345         gmx_file("Cannot write trajectory frame; maybe you are out of disk space?");
346     }
347 }
348
349
350 gmx_bool gmx_trr_read_frame(t_fileio* fio,
351                             int64_t*  step,
352                             real*     t,
353                             real*     lambda,
354                             rvec*     box,
355                             int*      natoms,
356                             rvec*     x,
357                             rvec*     v,
358                             rvec*     f)
359 {
360     return do_trr_frame(fio, true, step, t, lambda, box, natoms, x, v, f);
361 }
362
363 gmx_bool gmx_trr_read_frame_data(t_fileio* fio, gmx_trr_header_t* header, rvec* box, rvec* x, rvec* v, rvec* f)
364 {
365     return do_trr_frame_data(fio, header, box, x, v, f);
366 }
367
368 t_fileio* gmx_trr_open(const char* fn, const char* mode)
369 {
370     return gmx_fio_open(fn, mode);
371 }
372
373 void gmx_trr_close(t_fileio* fio)
374 {
375     gmx_fio_close(fio);
376 }