Upped the version to 3.2.0
[gromacs.git] / src / gmxlib / trnio.c
blob363272ac8b7f5f125937b64f8c937beaa78290a0
1 /*
2 * $Id$
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * GROningen Mixture of Alchemy and Childrens' Stories
37 #include <string.h>
38 #include "sysstuff.h"
39 #include "smalloc.h"
40 #include "fatal.h"
41 #include "txtdump.h"
42 #include "names.h"
43 #include "futil.h"
44 #include "trnio.h"
45 #include "gmxfio.h"
47 #define BUFSIZE 128
48 #define GROMACS_MAGIC 1993
50 static int nFloatSize(t_trnheader *sh)
52 int nflsize=0;
54 if (sh->box_size)
55 nflsize = sh->box_size/(DIM*DIM);
56 else if (sh->x_size)
57 nflsize = sh->x_size/(sh->natoms*DIM);
58 else if (sh->v_size)
59 nflsize = sh->v_size/(sh->natoms*DIM);
60 else if (sh->f_size)
61 nflsize = sh->f_size/(sh->natoms*DIM);
62 else
63 fatal_error(0,"Can not determine precision of trn file, quit!");
65 if (((nflsize != sizeof(float)) && (nflsize != sizeof(double))))
66 fatal_error(0,"Float size %d. Maybe different CPU?",nflsize);
68 return nflsize;
71 static bool do_trnheader(int fp,bool bRead,t_trnheader *sh, bool *bOK)
73 const int magic=GROMACS_MAGIC;
74 const char *version = "GMX_trn_file";
75 static bool bFirst=TRUE;
76 char buf[256];
77 bool bDouble;
79 *bOK=TRUE;
81 fio_select(fp);
82 if (!do_int(magic))
83 return FALSE;
85 if (bRead) {
86 *bOK = *bOK && do_string(buf);
87 if (bFirst)
88 fprintf(stderr,"trn version: %s ",buf);
90 else
91 *bOK = *bOK && do_string(version);
92 *bOK = *bOK && do_int(sh->ir_size);
93 *bOK = *bOK && do_int(sh->e_size);
94 *bOK = *bOK && do_int(sh->box_size);
95 *bOK = *bOK && do_int(sh->vir_size);
96 *bOK = *bOK && do_int(sh->pres_size);
97 *bOK = *bOK && do_int(sh->top_size);
98 *bOK = *bOK && do_int(sh->sym_size);
99 *bOK = *bOK && do_int(sh->x_size);
100 *bOK = *bOK && do_int(sh->v_size);
101 *bOK = *bOK && do_int(sh->f_size);
103 if (!*bOK) return *bOK;
104 bDouble = (nFloatSize(sh) == sizeof(double));
105 fio_setprecision(fp,bDouble);
107 if (bRead && bFirst) {
108 fprintf(stderr,"(%s precision)\n",bDouble ? "double" : "single");
109 bFirst = FALSE;
112 *bOK = *bOK && do_int(sh->natoms);
113 *bOK = *bOK && do_int(sh->step);
114 *bOK = *bOK && do_int(sh->nre);
115 *bOK = *bOK && do_real(sh->t);
116 *bOK = *bOK && do_real(sh->lambda);
118 return *bOK;
121 void pr_trnheader(FILE *fp,int indent,char *title,t_trnheader *sh)
123 if (sh) {
124 indent=pr_title(fp,indent,title);
125 (void) pr_indent(fp,indent);
126 (void) fprintf(fp,"box_size = %d\n",sh->box_size);
127 (void) pr_indent(fp,indent);
128 (void) fprintf(fp,"x_size = %d\n",sh->x_size);
129 (void) pr_indent(fp,indent);
130 (void) fprintf(fp,"v_size = %d\n",sh->v_size);
131 (void) pr_indent(fp,indent);
132 (void) fprintf(fp,"f_size = %d\n",sh->f_size);
134 (void) pr_indent(fp,indent);
135 (void) fprintf(fp,"natoms = %d\n",sh->natoms);
136 (void) pr_indent(fp,indent);
137 (void) fprintf(fp,"step = %d\n",sh->step);
138 (void) pr_indent(fp,indent);
139 (void) fprintf(fp,"t = %e\n",sh->t);
140 (void) pr_indent(fp,indent);
141 (void) fprintf(fp,"lambda = %e\n",sh->lambda);
145 static bool do_htrn(int fp,bool bRead,t_trnheader *sh,
146 rvec *box,rvec *x,rvec *v,rvec *f)
148 matrix pv;
149 bool bOK;
151 bOK = TRUE;
152 if (sh->box_size != 0) bOK = bOK && ndo_rvec(box,DIM);
153 if (sh->vir_size != 0) bOK = bOK && ndo_rvec(pv,DIM);
154 if (sh->pres_size!= 0) bOK = bOK && ndo_rvec(pv,DIM);
155 if (sh->x_size != 0) bOK = bOK && ndo_rvec(x,sh->natoms);
156 if (sh->v_size != 0) bOK = bOK && ndo_rvec(v,sh->natoms);
157 if (sh->f_size != 0) bOK = bOK && ndo_rvec(f,sh->natoms);
159 return bOK;
162 static bool do_trn(int fp,bool bRead,int *step,real *t,real *lambda,
163 rvec *box,int *natoms,rvec *x,rvec *v,rvec *f)
165 t_trnheader *sh;
166 bool bOK;
168 snew(sh,1);
169 if (!bRead) {
170 sh->box_size=(box)?sizeof(matrix):0;
171 sh->x_size=((x)?(*natoms*sizeof(x[0])):0);
172 sh->v_size=((v)?(*natoms*sizeof(v[0])):0);
173 sh->f_size=((f)?(*natoms*sizeof(f[0])):0);
174 sh->natoms = *natoms;
175 sh->step = *step;
176 sh->nre = 0;
177 sh->t = *t;
178 sh->lambda = *lambda;
180 if (!do_trnheader(fp,bRead,sh,&bOK))
181 return FALSE;
182 if (bRead) {
183 *natoms = sh->natoms;
184 *step = sh->step;
185 *t = sh->t;
186 *lambda = sh->lambda;
187 if (sh->ir_size)
188 fatal_error(0,"inputrec in trn file");
189 if (sh->e_size)
190 fatal_error(0,"energies in trn file");
191 if (sh->top_size)
192 fatal_error(0,"topology in trn file");
193 if (sh->sym_size)
194 fatal_error(0,"symbol table in trn file");
196 bOK = do_htrn(fp,bRead,sh,box,x,v,f);
198 sfree(sh);
200 return bOK;
203 /************************************************************
205 * The following routines are the exported ones
207 ************************************************************/
209 void read_trnheader(char *fn,t_trnheader *trn)
211 int fp;
212 bool bOK;
214 fp = open_trn(fn,"r");
215 if (!do_trnheader(fp,TRUE,trn,&bOK))
216 fatal_error(0,"Empty file %s",fn);
217 close_trn(fp);
220 bool fread_trnheader(int fp,t_trnheader *trn, bool *bOK)
222 return do_trnheader(fp,TRUE,trn,bOK);
225 void write_trn(char *fn,int step,real t,real lambda,
226 rvec *box,int natoms,rvec *x,rvec *v,rvec *f)
228 int fp;
230 fp = open_trn(fn,"w");
231 do_trn(fp,FALSE,&step,&t,&lambda,box,&natoms,x,v,f);
232 close_trn(fp);
235 void read_trn(char *fn,int *step,real *t,real *lambda,
236 rvec *box,int *natoms,rvec *x,rvec *v,rvec *f)
238 int fp;
240 fp = open_trn(fn,"r");
241 (void) do_trn(fp,TRUE,step,t,lambda,box,natoms,x,v,f);
242 close_trn(fp);
245 void fwrite_trn(int fp,int step,real t,real lambda,
246 rvec *box,int natoms,rvec *x,rvec *v,rvec *f)
248 (void) do_trn(fp,FALSE,&step,&t,&lambda,box,&natoms,x,v,f);
252 bool fread_trn(int fp,int *step,real *t,real *lambda,
253 rvec *box,int *natoms,rvec *x,rvec *v,rvec *f)
255 return do_trn(fp,TRUE,step,t,lambda,box,natoms,x,v,f);
258 bool fread_htrn(int fp,t_trnheader *trn,rvec *box,rvec *x,rvec *v,rvec *f)
260 return do_htrn(fp,TRUE,trn,box,x,v,f);
263 int open_trn(char *fn,char *mode)
265 return fio_open(fn,mode);
268 void close_trn(int fp)
270 fio_close(fp);