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.
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"
49 #define GROMACS_MAGIC 1993
51 static int nFloatSize(t_trnheader *sh)
57 nflsize = sh->box_size/(DIM*DIM);
61 nflsize = sh->x_size/(sh->natoms*DIM);
65 nflsize = sh->v_size/(sh->natoms*DIM);
69 nflsize = sh->f_size/(sh->natoms*DIM);
73 gmx_file("Can not determine precision of trn file");
76 if (((nflsize != sizeof(float)) && (nflsize != sizeof(double))))
78 gmx_fatal(FARGS, "Float size %d. Maybe different CPU?", nflsize);
84 static gmx_bool do_trnheader(t_fileio *fio, gmx_bool bRead, t_trnheader *sh, gmx_bool *bOK)
86 int magic = GROMACS_MAGIC;
87 static gmx_bool bFirst = TRUE;
92 gmx_fio_checktype(fio);
94 if (!gmx_fio_do_int(fio, magic) || magic != GROMACS_MAGIC)
101 *bOK = *bOK && gmx_fio_do_string(fio, buf);
104 fprintf(stderr, "trn version: %s ", buf);
109 sprintf(buf, "GMX_trn_file");
110 *bOK = *bOK && gmx_fio_do_string(fio, buf);
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);
128 sh->bDouble = (nFloatSize(sh) == sizeof(double));
129 gmx_fio_setprecision(fio, sh->bDouble);
133 fprintf(stderr, "(%s precision)\n", sh->bDouble ? "double" : "single");
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);
145 void pr_trnheader(FILE *fp, int indent, char *title, t_trnheader *sh)
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);
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);
170 static gmx_bool do_htrn(t_fileio *fio, t_trnheader *sh,
171 rvec *box, rvec *x, rvec *v, rvec *f)
177 if (sh->box_size != 0)
179 bOK = bOK && gmx_fio_ndo_rvec(fio, box, DIM);
181 if (sh->vir_size != 0)
183 bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
185 if (sh->pres_size != 0)
187 bOK = bOK && gmx_fio_ndo_rvec(fio, pv, DIM);
191 bOK = bOK && gmx_fio_ndo_rvec(fio, x, sh->natoms);
195 bOK = bOK && gmx_fio_ndo_rvec(fio, v, sh->natoms);
199 bOK = bOK && gmx_fio_ndo_rvec(fio, f, sh->natoms);
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)
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;
222 sh->lambda = *lambda;
224 if (!do_trnheader(fio, bRead, sh, &bOK))
230 *natoms = sh->natoms;
233 *lambda = sh->lambda;
236 gmx_file("inputrec in trn file");
240 gmx_file("energies in trn file");
244 gmx_file("topology in trn file");
248 gmx_file("symbol table in trn file");
251 bOK = do_htrn(fio, sh, box, x, v, f);
258 /************************************************************
260 * The following routines are the exported ones
262 ************************************************************/
264 void read_trnheader(const char *fn, t_trnheader *trn)
269 fio = open_trn(fn, "r");
270 if (!do_trnheader(fio, TRUE, trn, &bOK))
272 gmx_fatal(FARGS, "Empty file %s", fn);
277 gmx_bool fread_trnheader(t_fileio *fio, t_trnheader *trn, gmx_bool *bOK)
279 return do_trnheader(fio, TRUE, trn, bOK);
282 void write_trn(const char *fn, int step, real t, real lambda,
283 rvec *box, int natoms, rvec *x, rvec *v, rvec *f)
287 fio = open_trn(fn, "w");
288 do_trn(fio, FALSE, &step, &t, &lambda, box, &natoms, x, v, f);
292 void read_trn(const char *fn, int *step, real *t, real *lambda,
293 rvec *box, int *natoms, rvec *x, rvec *v, rvec *f)
297 fio = open_trn(fn, "r");
298 (void) do_trn(fio, TRUE, step, t, lambda, box, natoms, x, v, f);
302 void fwrite_trn(t_fileio *fio, int step, real t, real lambda,
303 rvec *box, int natoms, rvec *x, rvec *v, rvec *f)
305 if (do_trn(fio, FALSE, &step, &t, &lambda, box, &natoms, x, v, f) == FALSE)
307 gmx_file("Cannot write trajectory frame; maybe you are out of disk space?");
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)
315 return do_trn(fio, TRUE, step, t, lambda, box, natoms, x, v, f);
318 gmx_bool fread_htrn(t_fileio *fio, t_trnheader *trn, rvec *box, rvec *x, rvec *v,
321 return do_htrn(fio, trn, box, x, v, f);
324 t_fileio *open_trn(const char *fn, const char *mode)
326 return gmx_fio_open(fn, mode);
329 void close_trn(t_fileio *fio)