a01cd2da6b0f959cbab05baa3da17ac0b682338f
[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,2014, 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 "config.h"
40
41 #include <string.h>
42 #include "gromacs/utility/smalloc.h"
43 #include "gromacs/utility/fatalerror.h"
44 #include "gromacs/legacyheaders/txtdump.h"
45 #include "gromacs/legacyheaders/names.h"
46 #include "gromacs/utility/futil.h"
47 #include "trnio.h"
48 #include "gmxfio.h"
49
50 #define BUFSIZE     128
51 #define GROMACS_MAGIC   1993
52
53 static int nFloatSize(t_trnheader *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 trn 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 static gmx_bool do_trnheader(t_fileio *fio, gmx_bool bRead, t_trnheader *sh, gmx_bool *bOK)
87 {
88     int             magic  = GROMACS_MAGIC;
89     static gmx_bool bFirst = TRUE;
90     char            buf[256];
91
92     *bOK = TRUE;
93
94     gmx_fio_checktype(fio);
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, "trn 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 void pr_trnheader(FILE *fp, int indent, char *title, t_trnheader *sh)
148 {
149     if (sh)
150     {
151         indent = pr_title(fp, indent, title);
152         (void) pr_indent(fp, indent);
153         (void) fprintf(fp, "box_size    = %d\n", sh->box_size);
154         (void) pr_indent(fp, indent);
155         (void) fprintf(fp, "x_size      = %d\n", sh->x_size);
156         (void) pr_indent(fp, indent);
157         (void) fprintf(fp, "v_size      = %d\n", sh->v_size);
158         (void) pr_indent(fp, indent);
159         (void) fprintf(fp, "f_size      = %d\n", sh->f_size);
160
161         (void) pr_indent(fp, indent);
162         (void) fprintf(fp, "natoms      = %d\n", sh->natoms);
163         (void) pr_indent(fp, indent);
164         (void) fprintf(fp, "step        = %d\n", sh->step);
165         (void) pr_indent(fp, indent);
166         (void) fprintf(fp, "t           = %e\n", sh->t);
167         (void) pr_indent(fp, indent);
168         (void) fprintf(fp, "lambda      = %e\n", sh->lambda);
169     }
170 }
171
172 static gmx_bool do_htrn(t_fileio *fio, t_trnheader *sh,
173                         rvec *box, rvec *x, rvec *v, rvec *f)
174 {
175     matrix   pv;
176     gmx_bool bOK;
177
178     bOK = TRUE;
179     if (sh->box_size != 0)
180     {
181         bOK = bOK && gmx_fio_ndo_rvec(fio, box, DIM);
182     }
183     if (sh->vir_size != 0)
184     {
185         bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
186     }
187     if (sh->pres_size != 0)
188     {
189         bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
190     }
191     if (sh->x_size   != 0)
192     {
193         bOK = bOK && gmx_fio_ndo_rvec(fio, x, sh->natoms);
194     }
195     if (sh->v_size   != 0)
196     {
197         bOK = bOK && gmx_fio_ndo_rvec(fio, v, sh->natoms);
198     }
199     if (sh->f_size   != 0)
200     {
201         bOK = bOK && gmx_fio_ndo_rvec(fio, f, sh->natoms);
202     }
203
204     return bOK;
205 }
206
207 static gmx_bool do_trn(t_fileio *fio, gmx_bool bRead, int *step, real *t, real *lambda,
208                        rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
209 {
210     t_trnheader *sh;
211     gmx_bool     bOK;
212
213     snew(sh, 1);
214     if (!bRead)
215     {
216         sh->box_size = (box) ? sizeof(matrix) : 0;
217         sh->x_size   = ((x) ? (*natoms*sizeof(x[0])) : 0);
218         sh->v_size   = ((v) ? (*natoms*sizeof(v[0])) : 0);
219         sh->f_size   = ((f) ? (*natoms*sizeof(f[0])) : 0);
220         sh->natoms   = *natoms;
221         sh->step     = *step;
222         sh->nre      = 0;
223         sh->t        = *t;
224         sh->lambda   = *lambda;
225     }
226     if (!do_trnheader(fio, bRead, sh, &bOK))
227     {
228         return FALSE;
229     }
230     if (bRead)
231     {
232         *natoms = sh->natoms;
233         *step   = sh->step;
234         *t      = sh->t;
235         *lambda = sh->lambda;
236         if (sh->ir_size)
237         {
238             gmx_file("inputrec in trn file");
239         }
240         if (sh->e_size)
241         {
242             gmx_file("energies in trn file");
243         }
244         if (sh->top_size)
245         {
246             gmx_file("topology in trn file");
247         }
248         if (sh->sym_size)
249         {
250             gmx_file("symbol table in trn file");
251         }
252     }
253     bOK = do_htrn(fio, sh, box, x, v, f);
254
255     sfree(sh);
256
257     return bOK;
258 }
259
260 /************************************************************
261  *
262  *  The following routines are the exported ones
263  *
264  ************************************************************/
265
266 void read_trnheader(const char *fn, t_trnheader *trn)
267 {
268     t_fileio *fio;
269     gmx_bool  bOK;
270
271     fio = open_trn(fn, "r");
272     if (!do_trnheader(fio, TRUE, trn, &bOK))
273     {
274         gmx_fatal(FARGS, "Empty file %s", fn);
275     }
276     close_trn(fio);
277 }
278
279 gmx_bool fread_trnheader(t_fileio *fio, t_trnheader *trn, gmx_bool *bOK)
280 {
281     return do_trnheader(fio, TRUE, trn, bOK);
282 }
283
284 void write_trn(const char *fn, int step, real t, real lambda,
285                rvec *box, int natoms, rvec *x, rvec *v, rvec *f)
286 {
287     t_fileio *fio;
288
289     fio = open_trn(fn, "w");
290     do_trn(fio, FALSE, &step, &t, &lambda, box, &natoms, x, v, f);
291     close_trn(fio);
292 }
293
294 void read_trn(const char *fn, int *step, real *t, real *lambda,
295               rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
296 {
297     t_fileio *fio;
298
299     fio = open_trn(fn, "r");
300     (void) do_trn(fio, TRUE, step, t, lambda, box, natoms, x, v, f);
301     close_trn(fio);
302 }
303
304 void fwrite_trn(t_fileio *fio, int step, real t, real lambda,
305                 rvec *box, int natoms, rvec *x, rvec *v, rvec *f)
306 {
307     if (do_trn(fio, FALSE, &step, &t, &lambda, box, &natoms, x, v, f) == FALSE)
308     {
309         gmx_file("Cannot write trajectory frame; maybe you are out of disk space?");
310     }
311 }
312
313
314 gmx_bool fread_trn(t_fileio *fio, int *step, real *t, real *lambda,
315                    rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
316 {
317     return do_trn(fio, TRUE, step, t, lambda, box, natoms, x, v, f);
318 }
319
320 gmx_bool fread_htrn(t_fileio *fio, t_trnheader *trn, rvec *box, rvec *x, rvec *v,
321                     rvec *f)
322 {
323     return do_htrn(fio, trn, box, x, v, f);
324 }
325
326 t_fileio *open_trn(const char *fn, const char *mode)
327 {
328     return gmx_fio_open(fn, mode);
329 }
330
331 void close_trn(t_fileio *fio)
332 {
333     gmx_fio_close(fio);
334 }