Code beautification with uncrustify
[alexxy/gromacs.git] / src / gromacs / gmxlib / trnio.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * GROningen Mixture of Alchemy and Childrens' Stories
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <string.h>
40 #include "sysstuff.h"
41 #include "smalloc.h"
42 #include "gmx_fatal.h"
43 #include "txtdump.h"
44 #include "names.h"
45 #include "futil.h"
46 #include "trnio.h"
47 #include "gmxfio.h"
48
49 #define BUFSIZE     128
50 #define GROMACS_MAGIC   1993
51
52 static int nFloatSize(t_trnheader *sh)
53 {
54     int nflsize = 0;
55
56     if (sh->box_size)
57     {
58         nflsize = sh->box_size/(DIM*DIM);
59     }
60     else if (sh->x_size)
61     {
62         nflsize = sh->x_size/(sh->natoms*DIM);
63     }
64     else if (sh->v_size)
65     {
66         nflsize = sh->v_size/(sh->natoms*DIM);
67     }
68     else if (sh->f_size)
69     {
70         nflsize = sh->f_size/(sh->natoms*DIM);
71     }
72     else
73     {
74         gmx_file("Can not determine precision of trn file");
75     }
76
77     if (((nflsize != sizeof(float)) && (nflsize != sizeof(double))))
78     {
79         gmx_fatal(FARGS, "Float size %d. Maybe different CPU?", nflsize);
80     }
81
82     return nflsize;
83 }
84
85 static gmx_bool do_trnheader(t_fileio *fio, gmx_bool bRead, t_trnheader *sh, gmx_bool *bOK)
86 {
87     int             magic  = GROMACS_MAGIC;
88     static gmx_bool bFirst = TRUE;
89     char            buf[256];
90
91     *bOK = TRUE;
92
93     gmx_fio_checktype(fio);
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 void pr_trnheader(FILE *fp, int indent, char *title, t_trnheader *sh)
147 {
148     if (sh)
149     {
150         indent = pr_title(fp, indent, title);
151         (void) pr_indent(fp, indent);
152         (void) fprintf(fp, "box_size    = %d\n", sh->box_size);
153         (void) pr_indent(fp, indent);
154         (void) fprintf(fp, "x_size      = %d\n", sh->x_size);
155         (void) pr_indent(fp, indent);
156         (void) fprintf(fp, "v_size      = %d\n", sh->v_size);
157         (void) pr_indent(fp, indent);
158         (void) fprintf(fp, "f_size      = %d\n", sh->f_size);
159
160         (void) pr_indent(fp, indent);
161         (void) fprintf(fp, "natoms      = %d\n", sh->natoms);
162         (void) pr_indent(fp, indent);
163         (void) fprintf(fp, "step        = %d\n", sh->step);
164         (void) pr_indent(fp, indent);
165         (void) fprintf(fp, "t           = %e\n", sh->t);
166         (void) pr_indent(fp, indent);
167         (void) fprintf(fp, "lambda      = %e\n", sh->lambda);
168     }
169 }
170
171 static gmx_bool do_htrn(t_fileio *fio, gmx_bool bRead, t_trnheader *sh,
172                         rvec *box, rvec *x, rvec *v, rvec *f)
173 {
174     matrix   pv;
175     gmx_bool bOK;
176
177     bOK = TRUE;
178     if (sh->box_size != 0)
179     {
180         bOK = bOK && gmx_fio_ndo_rvec(fio, box, DIM);
181     }
182     if (sh->vir_size != 0)
183     {
184         bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
185     }
186     if (sh->pres_size != 0)
187     {
188         bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
189     }
190     if (sh->x_size   != 0)
191     {
192         bOK = bOK && gmx_fio_ndo_rvec(fio, x, sh->natoms);
193     }
194     if (sh->v_size   != 0)
195     {
196         bOK = bOK && gmx_fio_ndo_rvec(fio, v, sh->natoms);
197     }
198     if (sh->f_size   != 0)
199     {
200         bOK = bOK && gmx_fio_ndo_rvec(fio, f, sh->natoms);
201     }
202
203     return bOK;
204 }
205
206 static gmx_bool do_trn(t_fileio *fio, gmx_bool bRead, int *step, real *t, real *lambda,
207                        rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
208 {
209     t_trnheader *sh;
210     gmx_bool     bOK;
211
212     snew(sh, 1);
213     if (!bRead)
214     {
215         sh->box_size = (box) ? sizeof(matrix) : 0;
216         sh->x_size   = ((x) ? (*natoms*sizeof(x[0])) : 0);
217         sh->v_size   = ((v) ? (*natoms*sizeof(v[0])) : 0);
218         sh->f_size   = ((f) ? (*natoms*sizeof(f[0])) : 0);
219         sh->natoms   = *natoms;
220         sh->step     = *step;
221         sh->nre      = 0;
222         sh->t        = *t;
223         sh->lambda   = *lambda;
224     }
225     if (!do_trnheader(fio, bRead, sh, &bOK))
226     {
227         return FALSE;
228     }
229     if (bRead)
230     {
231         *natoms = sh->natoms;
232         *step   = sh->step;
233         *t      = sh->t;
234         *lambda = sh->lambda;
235         if (sh->ir_size)
236         {
237             gmx_file("inputrec in trn file");
238         }
239         if (sh->e_size)
240         {
241             gmx_file("energies in trn file");
242         }
243         if (sh->top_size)
244         {
245             gmx_file("topology in trn file");
246         }
247         if (sh->sym_size)
248         {
249             gmx_file("symbol table in trn file");
250         }
251     }
252     bOK = do_htrn(fio, bRead, sh, box, x, v, f);
253
254     sfree(sh);
255
256     return bOK;
257 }
258
259 /************************************************************
260  *
261  *  The following routines are the exported ones
262  *
263  ************************************************************/
264
265 void read_trnheader(const char *fn, t_trnheader *trn)
266 {
267     t_fileio *fio;
268     gmx_bool  bOK;
269
270     fio = open_trn(fn, "r");
271     if (!do_trnheader(fio, TRUE, trn, &bOK))
272     {
273         gmx_fatal(FARGS, "Empty file %s", fn);
274     }
275     close_trn(fio);
276 }
277
278 gmx_bool fread_trnheader(t_fileio *fio, t_trnheader *trn, gmx_bool *bOK)
279 {
280     return do_trnheader(fio, TRUE, trn, bOK);
281 }
282
283 void write_trn(const char *fn, int step, real t, real lambda,
284                rvec *box, int natoms, rvec *x, rvec *v, rvec *f)
285 {
286     t_fileio *fio;
287
288     fio = open_trn(fn, "w");
289     do_trn(fio, FALSE, &step, &t, &lambda, box, &natoms, x, v, f);
290     close_trn(fio);
291 }
292
293 void read_trn(const char *fn, int *step, real *t, real *lambda,
294               rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
295 {
296     t_fileio *fio;
297
298     fio = open_trn(fn, "r");
299     (void) do_trn(fio, TRUE, step, t, lambda, box, natoms, x, v, f);
300     close_trn(fio);
301 }
302
303 void fwrite_trn(t_fileio *fio, int step, real t, real lambda,
304                 rvec *box, int natoms, rvec *x, rvec *v, rvec *f)
305 {
306     if (do_trn(fio, FALSE, &step, &t, &lambda, box, &natoms, x, v, f) == FALSE)
307     {
308         gmx_file("Cannot write trajectory frame; maybe you are out of disk space?");
309     }
310 }
311
312
313 gmx_bool fread_trn(t_fileio *fio, int *step, real *t, real *lambda,
314                    rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
315 {
316     return do_trn(fio, TRUE, step, t, lambda, box, natoms, x, v, f);
317 }
318
319 gmx_bool fread_htrn(t_fileio *fio, t_trnheader *trn, rvec *box, rvec *x, rvec *v,
320                     rvec *f)
321 {
322     return do_htrn(fio, TRUE, trn, box, x, v, f);
323 }
324
325 t_fileio *open_trn(const char *fn, const char *mode)
326 {
327     return gmx_fio_open(fn, mode);
328 }
329
330 void close_trn(t_fileio *fio)
331 {
332     gmx_fio_close(fio);
333 }