2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
43 #include "gromacs/utility/smalloc.h"
44 #include "gmx_fatal.h"
52 #define GROMACS_MAGIC 1993
54 static int nFloatSize(t_trnheader *sh)
60 nflsize = sh->box_size/(DIM*DIM);
64 nflsize = sh->x_size/(sh->natoms*DIM);
68 nflsize = sh->v_size/(sh->natoms*DIM);
72 nflsize = sh->f_size/(sh->natoms*DIM);
76 gmx_file("Can not determine precision of trn file");
79 if (((nflsize != sizeof(float)) && (nflsize != sizeof(double))))
81 gmx_fatal(FARGS, "Float size %d. Maybe different CPU?", nflsize);
87 static gmx_bool do_trnheader(t_fileio *fio, gmx_bool bRead, t_trnheader *sh, gmx_bool *bOK)
89 int magic = GROMACS_MAGIC;
90 static gmx_bool bFirst = TRUE;
95 gmx_fio_checktype(fio);
97 if (!gmx_fio_do_int(fio, magic) || magic != GROMACS_MAGIC)
104 *bOK = *bOK && gmx_fio_do_string(fio, buf);
107 fprintf(stderr, "trn version: %s ", buf);
112 sprintf(buf, "GMX_trn_file");
113 *bOK = *bOK && gmx_fio_do_string(fio, buf);
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);
131 sh->bDouble = (nFloatSize(sh) == sizeof(double));
132 gmx_fio_setprecision(fio, sh->bDouble);
136 fprintf(stderr, "(%s precision)\n", sh->bDouble ? "double" : "single");
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);
148 void pr_trnheader(FILE *fp, int indent, char *title, t_trnheader *sh)
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);
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);
173 static gmx_bool do_htrn(t_fileio *fio, t_trnheader *sh,
174 rvec *box, rvec *x, rvec *v, rvec *f)
180 if (sh->box_size != 0)
182 bOK = bOK && gmx_fio_ndo_rvec(fio, box, DIM);
184 if (sh->vir_size != 0)
186 bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
188 if (sh->pres_size != 0)
190 bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
194 bOK = bOK && gmx_fio_ndo_rvec(fio, x, sh->natoms);
198 bOK = bOK && gmx_fio_ndo_rvec(fio, v, sh->natoms);
202 bOK = bOK && gmx_fio_ndo_rvec(fio, f, sh->natoms);
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)
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;
225 sh->lambda = *lambda;
227 if (!do_trnheader(fio, bRead, sh, &bOK))
233 *natoms = sh->natoms;
236 *lambda = sh->lambda;
239 gmx_file("inputrec in trn file");
243 gmx_file("energies in trn file");
247 gmx_file("topology in trn file");
251 gmx_file("symbol table in trn file");
254 bOK = do_htrn(fio, sh, box, x, v, f);
261 /************************************************************
263 * The following routines are the exported ones
265 ************************************************************/
267 void read_trnheader(const char *fn, t_trnheader *trn)
272 fio = open_trn(fn, "r");
273 if (!do_trnheader(fio, TRUE, trn, &bOK))
275 gmx_fatal(FARGS, "Empty file %s", fn);
280 gmx_bool fread_trnheader(t_fileio *fio, t_trnheader *trn, gmx_bool *bOK)
282 return do_trnheader(fio, TRUE, trn, bOK);
285 void write_trn(const char *fn, int step, real t, real lambda,
286 rvec *box, int natoms, rvec *x, rvec *v, rvec *f)
290 fio = open_trn(fn, "w");
291 do_trn(fio, FALSE, &step, &t, &lambda, box, &natoms, x, v, f);
295 void read_trn(const char *fn, int *step, real *t, real *lambda,
296 rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
300 fio = open_trn(fn, "r");
301 (void) do_trn(fio, TRUE, step, t, lambda, box, natoms, x, v, f);
305 void fwrite_trn(t_fileio *fio, int step, real t, real lambda,
306 rvec *box, int natoms, rvec *x, rvec *v, rvec *f)
308 if (do_trn(fio, FALSE, &step, &t, &lambda, box, &natoms, x, v, f) == FALSE)
310 gmx_file("Cannot write trajectory frame; maybe you are out of disk space?");
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)
318 return do_trn(fio, TRUE, step, t, lambda, box, natoms, x, v, f);
321 gmx_bool fread_htrn(t_fileio *fio, t_trnheader *trn, rvec *box, rvec *x, rvec *v,
324 return do_htrn(fio, trn, box, x, v, f);
327 t_fileio *open_trn(const char *fn, const char *mode)
329 return gmx_fio_open(fn, mode);
332 void close_trn(t_fileio *fio)