Rename trnio routines
[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, 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/legacyheaders/names.h"
46 #include "gromacs/legacyheaders/txtdump.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 #define GROMACS_MAGIC   1993
53
54 static int nFloatSize(gmx_trr_header_t *sh)
55 {
56     int nflsize = 0;
57
58     if (sh->box_size)
59     {
60         nflsize = sh->box_size/(DIM*DIM);
61     }
62     else if (sh->x_size)
63     {
64         nflsize = sh->x_size/(sh->natoms*DIM);
65     }
66     else if (sh->v_size)
67     {
68         nflsize = sh->v_size/(sh->natoms*DIM);
69     }
70     else if (sh->f_size)
71     {
72         nflsize = sh->f_size/(sh->natoms*DIM);
73     }
74     else
75     {
76         gmx_file("Can not determine precision of trr file");
77     }
78
79     if (((nflsize != sizeof(float)) && (nflsize != sizeof(double))))
80     {
81         gmx_fatal(FARGS, "Float size %d. Maybe different CPU?", nflsize);
82     }
83
84     return nflsize;
85 }
86
87 static gmx_bool
88 do_trr_frame_header(t_fileio *fio, bool bRead, gmx_trr_header_t *sh, gmx_bool *bOK)
89 {
90     int             magic  = GROMACS_MAGIC;
91     static gmx_bool bFirst = TRUE;
92     char            buf[256];
93
94     *bOK = TRUE;
95
96     if (!gmx_fio_do_int(fio, magic) || magic != GROMACS_MAGIC)
97     {
98         return FALSE;
99     }
100
101     if (bRead)
102     {
103         *bOK = *bOK && gmx_fio_do_string(fio, buf);
104         if (bFirst)
105         {
106             fprintf(stderr, "trr version: %s ", buf);
107         }
108     }
109     else
110     {
111         sprintf(buf, "GMX_trn_file");
112         *bOK = *bOK && gmx_fio_do_string(fio, buf);
113     }
114     *bOK = *bOK && gmx_fio_do_int(fio, sh->ir_size);
115     *bOK = *bOK && gmx_fio_do_int(fio, sh->e_size);
116     *bOK = *bOK && gmx_fio_do_int(fio, sh->box_size);
117     *bOK = *bOK && gmx_fio_do_int(fio, sh->vir_size);
118     *bOK = *bOK && gmx_fio_do_int(fio, sh->pres_size);
119     *bOK = *bOK && gmx_fio_do_int(fio, sh->top_size);
120     *bOK = *bOK && gmx_fio_do_int(fio, sh->sym_size);
121     *bOK = *bOK && gmx_fio_do_int(fio, sh->x_size);
122     *bOK = *bOK && gmx_fio_do_int(fio, sh->v_size);
123     *bOK = *bOK && gmx_fio_do_int(fio, sh->f_size);
124     *bOK = *bOK && gmx_fio_do_int(fio, sh->natoms);
125
126     if (!*bOK)
127     {
128         return *bOK;
129     }
130     sh->bDouble = (nFloatSize(sh) == sizeof(double));
131     gmx_fio_setprecision(fio, sh->bDouble);
132
133     if (bRead && bFirst)
134     {
135         fprintf(stderr, "(%s precision)\n", sh->bDouble ? "double" : "single");
136         bFirst = FALSE;
137     }
138
139     *bOK = *bOK && gmx_fio_do_int(fio, sh->step);
140     *bOK = *bOK && gmx_fio_do_int(fio, sh->nre);
141     *bOK = *bOK && gmx_fio_do_real(fio, sh->t);
142     *bOK = *bOK && gmx_fio_do_real(fio, sh->lambda);
143
144     return *bOK;
145 }
146
147 static gmx_bool
148 do_trr_frame_data(t_fileio *fio, gmx_trr_header_t *sh,
149                   rvec *box, rvec *x, rvec *v, rvec *f)
150 {
151     matrix   pv;
152     gmx_bool bOK;
153
154     bOK = TRUE;
155     if (sh->box_size != 0)
156     {
157         bOK = bOK && gmx_fio_ndo_rvec(fio, box, DIM);
158     }
159     if (sh->vir_size != 0)
160     {
161         bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
162     }
163     if (sh->pres_size != 0)
164     {
165         bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
166     }
167     if (sh->x_size   != 0)
168     {
169         bOK = bOK && gmx_fio_ndo_rvec(fio, x, sh->natoms);
170     }
171     if (sh->v_size   != 0)
172     {
173         bOK = bOK && gmx_fio_ndo_rvec(fio, v, sh->natoms);
174     }
175     if (sh->f_size   != 0)
176     {
177         bOK = bOK && gmx_fio_ndo_rvec(fio, f, sh->natoms);
178     }
179
180     return bOK;
181 }
182
183 static gmx_bool
184 do_trr_frame(t_fileio *fio, bool bRead, int *step, real *t, real *lambda,
185              rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
186 {
187     gmx_trr_header_t *sh;
188     gmx_bool          bOK;
189
190     snew(sh, 1);
191     if (!bRead)
192     {
193         sh->box_size = (box) ? sizeof(matrix) : 0;
194         sh->x_size   = ((x) ? (*natoms*sizeof(x[0])) : 0);
195         sh->v_size   = ((v) ? (*natoms*sizeof(v[0])) : 0);
196         sh->f_size   = ((f) ? (*natoms*sizeof(f[0])) : 0);
197         sh->natoms   = *natoms;
198         sh->step     = *step;
199         sh->nre      = 0;
200         sh->t        = *t;
201         sh->lambda   = *lambda;
202     }
203     if (!do_trr_frame_header(fio, bRead, sh, &bOK))
204     {
205         return FALSE;
206     }
207     if (bRead)
208     {
209         *natoms = sh->natoms;
210         *step   = sh->step;
211         *t      = sh->t;
212         *lambda = sh->lambda;
213         if (sh->ir_size)
214         {
215             gmx_file("inputrec in trr file");
216         }
217         if (sh->e_size)
218         {
219             gmx_file("energies in trr file");
220         }
221         if (sh->top_size)
222         {
223             gmx_file("topology in trr file");
224         }
225         if (sh->sym_size)
226         {
227             gmx_file("symbol table in trr file");
228         }
229     }
230     bOK = do_trr_frame_data(fio, sh, box, x, v, f);
231
232     sfree(sh);
233
234     return bOK;
235 }
236
237 /************************************************************
238  *
239  *  The following routines are the exported ones
240  *
241  ************************************************************/
242
243 void gmx_trr_read_single_header(const char *fn, gmx_trr_header_t *header)
244 {
245     t_fileio *fio = gmx_trr_open(fn, "r");
246     gmx_bool  bOK;
247     if (!do_trr_frame_header(fio, true, header, &bOK))
248     {
249         gmx_fatal(FARGS, "Empty file %s", fn);
250     }
251     gmx_trr_close(fio);
252 }
253
254 gmx_bool gmx_trr_read_frame_header(t_fileio *fio, gmx_trr_header_t *header, gmx_bool *bOK)
255 {
256     return do_trr_frame_header(fio, true, header, bOK);
257 }
258
259 void gmx_trr_write_single_frame(const char *fn, int step, real t, real lambda,
260                                 rvec *box, int natoms, rvec *x, rvec *v, rvec *f)
261 {
262     t_fileio *fio = gmx_trr_open(fn, "w");
263     do_trr_frame(fio, false, &step, &t, &lambda, box, &natoms, x, v, f);
264     gmx_trr_close(fio);
265 }
266
267 void gmx_trr_read_single_frame(const char *fn, int *step, real *t, real *lambda,
268                                rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
269 {
270     t_fileio *fio = gmx_trr_open(fn, "r");
271     do_trr_frame(fio, true, step, t, lambda, box, natoms, x, v, f);
272     gmx_trr_close(fio);
273 }
274
275 void gmx_trr_write_frame(t_fileio *fio, int step, real t, real lambda,
276                          rvec *box, int natoms, rvec *x, rvec *v, rvec *f)
277 {
278     if (!do_trr_frame(fio, false, &step, &t, &lambda, box, &natoms, x, v, f))
279     {
280         gmx_file("Cannot write trajectory frame; maybe you are out of disk space?");
281     }
282 }
283
284
285 gmx_bool gmx_trr_read_frame(t_fileio *fio, int *step, real *t, real *lambda,
286                             rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
287 {
288     return do_trr_frame(fio, true, step, t, lambda, box, natoms, x, v, f);
289 }
290
291 gmx_bool gmx_trr_read_frame_data(t_fileio *fio, gmx_trr_header_t *header,
292                                  rvec *box, rvec *x, rvec *v, rvec *f)
293 {
294     return do_trr_frame_data(fio, header, box, x, v, f);
295 }
296
297 t_fileio *gmx_trr_open(const char *fn, const char *mode)
298 {
299     return gmx_fio_open(fn, mode);
300 }
301
302 void gmx_trr_close(t_fileio *fio)
303 {
304     gmx_fio_close(fio);
305 }