Split espresso I/O routines from confio.cpp
[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(t_trnheader *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 trn 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 do_trnheader(t_fileio *fio, gmx_bool bRead, t_trnheader *sh, gmx_bool *bOK)
88 {
89     int             magic  = GROMACS_MAGIC;
90     static gmx_bool bFirst = TRUE;
91     char            buf[256];
92
93     *bOK = TRUE;
94
95     if (!gmx_fio_do_int(fio, magic) || magic != GROMACS_MAGIC)
96     {
97         return FALSE;
98     }
99
100     if (bRead)
101     {
102         *bOK = *bOK && gmx_fio_do_string(fio, buf);
103         if (bFirst)
104         {
105             fprintf(stderr, "trn version: %s ", buf);
106         }
107     }
108     else
109     {
110         sprintf(buf, "GMX_trn_file");
111         *bOK = *bOK && gmx_fio_do_string(fio, buf);
112     }
113     *bOK = *bOK && gmx_fio_do_int(fio, sh->ir_size);
114     *bOK = *bOK && gmx_fio_do_int(fio, sh->e_size);
115     *bOK = *bOK && gmx_fio_do_int(fio, sh->box_size);
116     *bOK = *bOK && gmx_fio_do_int(fio, sh->vir_size);
117     *bOK = *bOK && gmx_fio_do_int(fio, sh->pres_size);
118     *bOK = *bOK && gmx_fio_do_int(fio, sh->top_size);
119     *bOK = *bOK && gmx_fio_do_int(fio, sh->sym_size);
120     *bOK = *bOK && gmx_fio_do_int(fio, sh->x_size);
121     *bOK = *bOK && gmx_fio_do_int(fio, sh->v_size);
122     *bOK = *bOK && gmx_fio_do_int(fio, sh->f_size);
123     *bOK = *bOK && gmx_fio_do_int(fio, sh->natoms);
124
125     if (!*bOK)
126     {
127         return *bOK;
128     }
129     sh->bDouble = (nFloatSize(sh) == sizeof(double));
130     gmx_fio_setprecision(fio, sh->bDouble);
131
132     if (bRead && bFirst)
133     {
134         fprintf(stderr, "(%s precision)\n", sh->bDouble ? "double" : "single");
135         bFirst = FALSE;
136     }
137
138     *bOK = *bOK && gmx_fio_do_int(fio, sh->step);
139     *bOK = *bOK && gmx_fio_do_int(fio, sh->nre);
140     *bOK = *bOK && gmx_fio_do_real(fio, sh->t);
141     *bOK = *bOK && gmx_fio_do_real(fio, sh->lambda);
142
143     return *bOK;
144 }
145
146 static gmx_bool do_htrn(t_fileio *fio, t_trnheader *sh,
147                         rvec *box, rvec *x, rvec *v, rvec *f)
148 {
149     matrix   pv;
150     gmx_bool bOK;
151
152     bOK = TRUE;
153     if (sh->box_size != 0)
154     {
155         bOK = bOK && gmx_fio_ndo_rvec(fio, box, DIM);
156     }
157     if (sh->vir_size != 0)
158     {
159         bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
160     }
161     if (sh->pres_size != 0)
162     {
163         bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
164     }
165     if (sh->x_size   != 0)
166     {
167         bOK = bOK && gmx_fio_ndo_rvec(fio, x, sh->natoms);
168     }
169     if (sh->v_size   != 0)
170     {
171         bOK = bOK && gmx_fio_ndo_rvec(fio, v, sh->natoms);
172     }
173     if (sh->f_size   != 0)
174     {
175         bOK = bOK && gmx_fio_ndo_rvec(fio, f, sh->natoms);
176     }
177
178     return bOK;
179 }
180
181 static gmx_bool do_trn(t_fileio *fio, gmx_bool bRead, int *step, real *t, real *lambda,
182                        rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
183 {
184     t_trnheader *sh;
185     gmx_bool     bOK;
186
187     snew(sh, 1);
188     if (!bRead)
189     {
190         sh->box_size = (box) ? sizeof(matrix) : 0;
191         sh->x_size   = ((x) ? (*natoms*sizeof(x[0])) : 0);
192         sh->v_size   = ((v) ? (*natoms*sizeof(v[0])) : 0);
193         sh->f_size   = ((f) ? (*natoms*sizeof(f[0])) : 0);
194         sh->natoms   = *natoms;
195         sh->step     = *step;
196         sh->nre      = 0;
197         sh->t        = *t;
198         sh->lambda   = *lambda;
199     }
200     if (!do_trnheader(fio, bRead, sh, &bOK))
201     {
202         return FALSE;
203     }
204     if (bRead)
205     {
206         *natoms = sh->natoms;
207         *step   = sh->step;
208         *t      = sh->t;
209         *lambda = sh->lambda;
210         if (sh->ir_size)
211         {
212             gmx_file("inputrec in trn file");
213         }
214         if (sh->e_size)
215         {
216             gmx_file("energies in trn file");
217         }
218         if (sh->top_size)
219         {
220             gmx_file("topology in trn file");
221         }
222         if (sh->sym_size)
223         {
224             gmx_file("symbol table in trn file");
225         }
226     }
227     bOK = do_htrn(fio, sh, box, x, v, f);
228
229     sfree(sh);
230
231     return bOK;
232 }
233
234 /************************************************************
235  *
236  *  The following routines are the exported ones
237  *
238  ************************************************************/
239
240 void read_trnheader(const char *fn, t_trnheader *trn)
241 {
242     t_fileio *fio;
243     gmx_bool  bOK;
244
245     fio = open_trn(fn, "r");
246     if (!do_trnheader(fio, TRUE, trn, &bOK))
247     {
248         gmx_fatal(FARGS, "Empty file %s", fn);
249     }
250     close_trn(fio);
251 }
252
253 gmx_bool fread_trnheader(t_fileio *fio, t_trnheader *trn, gmx_bool *bOK)
254 {
255     return do_trnheader(fio, TRUE, trn, bOK);
256 }
257
258 void write_trn(const char *fn, int step, real t, real lambda,
259                rvec *box, int natoms, rvec *x, rvec *v, rvec *f)
260 {
261     t_fileio *fio;
262
263     fio = open_trn(fn, "w");
264     do_trn(fio, FALSE, &step, &t, &lambda, box, &natoms, x, v, f);
265     close_trn(fio);
266 }
267
268 void read_trn(const char *fn, int *step, real *t, real *lambda,
269               rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
270 {
271     t_fileio *fio;
272
273     fio = open_trn(fn, "r");
274     (void) do_trn(fio, TRUE, step, t, lambda, box, natoms, x, v, f);
275     close_trn(fio);
276 }
277
278 void fwrite_trn(t_fileio *fio, int step, real t, real lambda,
279                 rvec *box, int natoms, rvec *x, rvec *v, rvec *f)
280 {
281     if (do_trn(fio, FALSE, &step, &t, &lambda, box, &natoms, x, v, f) == FALSE)
282     {
283         gmx_file("Cannot write trajectory frame; maybe you are out of disk space?");
284     }
285 }
286
287
288 gmx_bool fread_trn(t_fileio *fio, int *step, real *t, real *lambda,
289                    rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
290 {
291     return do_trn(fio, TRUE, step, t, lambda, box, natoms, x, v, f);
292 }
293
294 gmx_bool fread_htrn(t_fileio *fio, t_trnheader *trn, rvec *box, rvec *x, rvec *v,
295                     rvec *f)
296 {
297     return do_htrn(fio, trn, box, x, v, f);
298 }
299
300 t_fileio *open_trn(const char *fn, const char *mode)
301 {
302     return gmx_fio_open(fn, mode);
303 }
304
305 void close_trn(t_fileio *fio)
306 {
307     gmx_fio_close(fio);
308 }