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