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