Merge branch release-4-6
[alexxy/gromacs.git] / src / gromacs / fileio / trnio.c
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, 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 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <string.h>
42 #include "sysstuff.h"
43 #include "smalloc.h"
44 #include "gmx_fatal.h"
45 #include "txtdump.h"
46 #include "names.h"
47 #include "futil.h"
48 #include "trnio.h"
49 #include "gmxfio.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     gmx_fio_checktype(fio);
96
97     if (!gmx_fio_do_int(fio, magic) || magic != GROMACS_MAGIC)
98     {
99         return FALSE;
100     }
101
102     if (bRead)
103     {
104         *bOK = *bOK && gmx_fio_do_string(fio, buf);
105         if (bFirst)
106         {
107             fprintf(stderr, "trn version: %s ", buf);
108         }
109     }
110     else
111     {
112         sprintf(buf, "GMX_trn_file");
113         *bOK = *bOK && gmx_fio_do_string(fio, buf);
114     }
115     *bOK = *bOK && gmx_fio_do_int(fio, sh->ir_size);
116     *bOK = *bOK && gmx_fio_do_int(fio, sh->e_size);
117     *bOK = *bOK && gmx_fio_do_int(fio, sh->box_size);
118     *bOK = *bOK && gmx_fio_do_int(fio, sh->vir_size);
119     *bOK = *bOK && gmx_fio_do_int(fio, sh->pres_size);
120     *bOK = *bOK && gmx_fio_do_int(fio, sh->top_size);
121     *bOK = *bOK && gmx_fio_do_int(fio, sh->sym_size);
122     *bOK = *bOK && gmx_fio_do_int(fio, sh->x_size);
123     *bOK = *bOK && gmx_fio_do_int(fio, sh->v_size);
124     *bOK = *bOK && gmx_fio_do_int(fio, sh->f_size);
125     *bOK = *bOK && gmx_fio_do_int(fio, sh->natoms);
126
127     if (!*bOK)
128     {
129         return *bOK;
130     }
131     sh->bDouble = (nFloatSize(sh) == sizeof(double));
132     gmx_fio_setprecision(fio, sh->bDouble);
133
134     if (bRead && bFirst)
135     {
136         fprintf(stderr, "(%s precision)\n", sh->bDouble ? "double" : "single");
137         bFirst = FALSE;
138     }
139
140     *bOK = *bOK && gmx_fio_do_int(fio, sh->step);
141     *bOK = *bOK && gmx_fio_do_int(fio, sh->nre);
142     *bOK = *bOK && gmx_fio_do_real(fio, sh->t);
143     *bOK = *bOK && gmx_fio_do_real(fio, sh->lambda);
144
145     return *bOK;
146 }
147
148 void pr_trnheader(FILE *fp, int indent, char *title, t_trnheader *sh)
149 {
150     if (sh)
151     {
152         indent = pr_title(fp, indent, title);
153         (void) pr_indent(fp, indent);
154         (void) fprintf(fp, "box_size    = %d\n", sh->box_size);
155         (void) pr_indent(fp, indent);
156         (void) fprintf(fp, "x_size      = %d\n", sh->x_size);
157         (void) pr_indent(fp, indent);
158         (void) fprintf(fp, "v_size      = %d\n", sh->v_size);
159         (void) pr_indent(fp, indent);
160         (void) fprintf(fp, "f_size      = %d\n", sh->f_size);
161
162         (void) pr_indent(fp, indent);
163         (void) fprintf(fp, "natoms      = %d\n", sh->natoms);
164         (void) pr_indent(fp, indent);
165         (void) fprintf(fp, "step        = %d\n", sh->step);
166         (void) pr_indent(fp, indent);
167         (void) fprintf(fp, "t           = %e\n", sh->t);
168         (void) pr_indent(fp, indent);
169         (void) fprintf(fp, "lambda      = %e\n", sh->lambda);
170     }
171 }
172
173 static gmx_bool do_htrn(t_fileio *fio, t_trnheader *sh,
174                         rvec *box, rvec *x, rvec *v, rvec *f)
175 {
176     matrix   pv;
177     gmx_bool bOK;
178
179     bOK = TRUE;
180     if (sh->box_size != 0)
181     {
182         bOK = bOK && gmx_fio_ndo_rvec(fio, box, DIM);
183     }
184     if (sh->vir_size != 0)
185     {
186         bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
187     }
188     if (sh->pres_size != 0)
189     {
190         bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
191     }
192     if (sh->x_size   != 0)
193     {
194         bOK = bOK && gmx_fio_ndo_rvec(fio, x, sh->natoms);
195     }
196     if (sh->v_size   != 0)
197     {
198         bOK = bOK && gmx_fio_ndo_rvec(fio, v, sh->natoms);
199     }
200     if (sh->f_size   != 0)
201     {
202         bOK = bOK && gmx_fio_ndo_rvec(fio, f, sh->natoms);
203     }
204
205     return bOK;
206 }
207
208 static gmx_bool do_trn(t_fileio *fio, gmx_bool bRead, int *step, real *t, real *lambda,
209                        rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
210 {
211     t_trnheader *sh;
212     gmx_bool     bOK;
213
214     snew(sh, 1);
215     if (!bRead)
216     {
217         sh->box_size = (box) ? sizeof(matrix) : 0;
218         sh->x_size   = ((x) ? (*natoms*sizeof(x[0])) : 0);
219         sh->v_size   = ((v) ? (*natoms*sizeof(v[0])) : 0);
220         sh->f_size   = ((f) ? (*natoms*sizeof(f[0])) : 0);
221         sh->natoms   = *natoms;
222         sh->step     = *step;
223         sh->nre      = 0;
224         sh->t        = *t;
225         sh->lambda   = *lambda;
226     }
227     if (!do_trnheader(fio, bRead, sh, &bOK))
228     {
229         return FALSE;
230     }
231     if (bRead)
232     {
233         *natoms = sh->natoms;
234         *step   = sh->step;
235         *t      = sh->t;
236         *lambda = sh->lambda;
237         if (sh->ir_size)
238         {
239             gmx_file("inputrec in trn file");
240         }
241         if (sh->e_size)
242         {
243             gmx_file("energies in trn file");
244         }
245         if (sh->top_size)
246         {
247             gmx_file("topology in trn file");
248         }
249         if (sh->sym_size)
250         {
251             gmx_file("symbol table in trn file");
252         }
253     }
254     bOK = do_htrn(fio, sh, box, x, v, f);
255
256     sfree(sh);
257
258     return bOK;
259 }
260
261 /************************************************************
262  *
263  *  The following routines are the exported ones
264  *
265  ************************************************************/
266
267 void read_trnheader(const char *fn, t_trnheader *trn)
268 {
269     t_fileio *fio;
270     gmx_bool  bOK;
271
272     fio = open_trn(fn, "r");
273     if (!do_trnheader(fio, TRUE, trn, &bOK))
274     {
275         gmx_fatal(FARGS, "Empty file %s", fn);
276     }
277     close_trn(fio);
278 }
279
280 gmx_bool fread_trnheader(t_fileio *fio, t_trnheader *trn, gmx_bool *bOK)
281 {
282     return do_trnheader(fio, TRUE, trn, bOK);
283 }
284
285 void write_trn(const char *fn, int step, real t, real lambda,
286                rvec *box, int natoms, rvec *x, rvec *v, rvec *f)
287 {
288     t_fileio *fio;
289
290     fio = open_trn(fn, "w");
291     do_trn(fio, FALSE, &step, &t, &lambda, box, &natoms, x, v, f);
292     close_trn(fio);
293 }
294
295 void read_trn(const char *fn, int *step, real *t, real *lambda,
296               rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
297 {
298     t_fileio *fio;
299
300     fio = open_trn(fn, "r");
301     (void) do_trn(fio, TRUE, step, t, lambda, box, natoms, x, v, f);
302     close_trn(fio);
303 }
304
305 void fwrite_trn(t_fileio *fio, int step, real t, real lambda,
306                 rvec *box, int natoms, rvec *x, rvec *v, rvec *f)
307 {
308     if (do_trn(fio, FALSE, &step, &t, &lambda, box, &natoms, x, v, f) == FALSE)
309     {
310         gmx_file("Cannot write trajectory frame; maybe you are out of disk space?");
311     }
312 }
313
314
315 gmx_bool fread_trn(t_fileio *fio, int *step, real *t, real *lambda,
316                    rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
317 {
318     return do_trn(fio, TRUE, step, t, lambda, box, natoms, x, v, f);
319 }
320
321 gmx_bool fread_htrn(t_fileio *fio, t_trnheader *trn, rvec *box, rvec *x, rvec *v,
322                     rvec *f)
323 {
324     return do_htrn(fio, trn, box, x, v, f);
325 }
326
327 t_fileio *open_trn(const char *fn, const char *mode)
328 {
329     return gmx_fio_open(fn, mode);
330 }
331
332 void close_trn(t_fileio *fio)
333 {
334     gmx_fio_close(fio);
335 }