3cc7ef1d470f2115f792a21a50ab4d2b201dd928
[alexxy/gromacs.git] / src / gmxlib / trnio.c
1 /*
2  * 
3  *                This source code is part of
4  * 
5  *                 G   R   O   M   A   C   S
6  * 
7  *          GROningen MAchine for Chemical Simulations
8  * 
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  * 
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  * 
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  * 
32  * And Hey:
33  * GROningen Mixture of Alchemy and Childrens' Stories
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <string.h>
40 #include "sysstuff.h"
41 #include "smalloc.h"
42 #include "gmx_fatal.h"
43 #include "txtdump.h"
44 #include "names.h"
45 #include "futil.h"
46 #include "trnio.h"
47 #include "gmxfio.h"
48
49 #define BUFSIZE         128
50 #define GROMACS_MAGIC   1993
51
52 static int nFloatSize(t_trnheader *sh)
53 {
54   int nflsize=0;
55   
56   if (sh->box_size)
57     nflsize = sh->box_size/(DIM*DIM);
58   else if (sh->x_size)
59     nflsize = sh->x_size/(sh->natoms*DIM);
60   else if (sh->v_size)
61     nflsize = sh->v_size/(sh->natoms*DIM);
62   else if (sh->f_size)
63     nflsize = sh->f_size/(sh->natoms*DIM);
64   else 
65     gmx_file("Can not determine precision of trn file");
66   
67   if (((nflsize != sizeof(float)) && (nflsize != sizeof(double))))
68     gmx_fatal(FARGS,"Float size %d. Maybe different CPU?",nflsize);
69       
70   return nflsize;
71 }
72
73 static gmx_bool do_trnheader(t_fileio *fio,gmx_bool bRead,t_trnheader *sh, gmx_bool *bOK)
74 {
75   int magic=GROMACS_MAGIC;
76   static gmx_bool bFirst=TRUE;
77   char buf[256];
78   
79   *bOK=TRUE;
80
81   gmx_fio_checktype(fio);
82
83   if (!gmx_fio_do_int(fio,magic) || magic!=GROMACS_MAGIC)
84     return FALSE;
85   
86   if (bRead) {
87     *bOK = *bOK && gmx_fio_do_string(fio,buf);
88     if (bFirst)
89       fprintf(stderr,"trn version: %s ",buf);
90   }
91   else {
92     sprintf(buf,"GMX_trn_file");
93     *bOK = *bOK && gmx_fio_do_string(fio,buf);
94   }
95   *bOK = *bOK && gmx_fio_do_int(fio,sh->ir_size);
96   *bOK = *bOK && gmx_fio_do_int(fio,sh->e_size);
97   *bOK = *bOK && gmx_fio_do_int(fio,sh->box_size);
98   *bOK = *bOK && gmx_fio_do_int(fio,sh->vir_size);
99   *bOK = *bOK && gmx_fio_do_int(fio,sh->pres_size);
100   *bOK = *bOK && gmx_fio_do_int(fio,sh->top_size); 
101   *bOK = *bOK && gmx_fio_do_int(fio,sh->sym_size); 
102   *bOK = *bOK && gmx_fio_do_int(fio,sh->x_size); 
103   *bOK = *bOK && gmx_fio_do_int(fio,sh->v_size); 
104   *bOK = *bOK && gmx_fio_do_int(fio,sh->f_size); 
105   *bOK = *bOK && gmx_fio_do_int(fio,sh->natoms);
106
107   if (!*bOK) return *bOK; 
108   sh->bDouble = (nFloatSize(sh) == sizeof(double));
109   gmx_fio_setprecision(fio,sh->bDouble);
110
111   if (bRead && bFirst) {
112     fprintf(stderr,"(%s precision)\n",sh->bDouble ? "double" : "single");
113     bFirst = FALSE;
114   }
115   
116   *bOK = *bOK && gmx_fio_do_int(fio,sh->step); 
117   *bOK = *bOK && gmx_fio_do_int(fio,sh->nre); 
118   *bOK = *bOK && gmx_fio_do_real(fio,sh->t); 
119   *bOK = *bOK && gmx_fio_do_real(fio,sh->lambda); 
120   
121   return *bOK;
122 }
123
124 void pr_trnheader(FILE *fp,int indent,char *title,t_trnheader *sh)
125 {
126   if (sh) {
127     indent=pr_title(fp,indent,title);
128     (void) pr_indent(fp,indent);
129     (void) fprintf(fp,"box_size    = %d\n",sh->box_size);
130     (void) pr_indent(fp,indent);
131     (void) fprintf(fp,"x_size      = %d\n",sh->x_size);
132     (void) pr_indent(fp,indent);
133     (void) fprintf(fp,"v_size      = %d\n",sh->v_size);
134     (void) pr_indent(fp,indent);
135     (void) fprintf(fp,"f_size      = %d\n",sh->f_size);
136     
137     (void) pr_indent(fp,indent);
138     (void) fprintf(fp,"natoms      = %d\n",sh->natoms);
139     (void) pr_indent(fp,indent);
140     (void) fprintf(fp,"step        = %d\n",sh->step);
141     (void) pr_indent(fp,indent);
142     (void) fprintf(fp,"t           = %e\n",sh->t);
143     (void) pr_indent(fp,indent);
144     (void) fprintf(fp,"lambda      = %e\n",sh->lambda);
145   }
146 }
147
148 static gmx_bool do_htrn(t_fileio *fio,gmx_bool bRead,t_trnheader *sh,
149                     rvec *box,rvec *x,rvec *v,rvec *f)
150 {
151   matrix pv;
152   gmx_bool bOK;
153
154   bOK = TRUE;
155   if (sh->box_size != 0) bOK = bOK && gmx_fio_ndo_rvec(fio,box,DIM);
156   if (sh->vir_size != 0) bOK = bOK && gmx_fio_ndo_rvec(fio,pv,DIM);
157   if (sh->pres_size!= 0) bOK = bOK && gmx_fio_ndo_rvec(fio,pv,DIM);
158   if (sh->x_size   != 0) bOK = bOK && gmx_fio_ndo_rvec(fio,x,sh->natoms);
159   if (sh->v_size   != 0) bOK = bOK && gmx_fio_ndo_rvec(fio,v,sh->natoms);
160   if (sh->f_size   != 0) bOK = bOK && gmx_fio_ndo_rvec(fio,f,sh->natoms);
161
162   return bOK;
163 }
164
165 static gmx_bool do_trn(t_fileio *fio,gmx_bool bRead,int *step,real *t,real *lambda,
166                    rvec *box,int *natoms,rvec *x,rvec *v,rvec *f)
167 {
168   t_trnheader *sh;
169   gmx_bool bOK;
170   
171   snew(sh,1);
172   if (!bRead) {
173     sh->box_size=(box)?sizeof(matrix):0;
174     sh->x_size=((x)?(*natoms*sizeof(x[0])):0);
175     sh->v_size=((v)?(*natoms*sizeof(v[0])):0);
176     sh->f_size=((f)?(*natoms*sizeof(f[0])):0);
177     sh->natoms = *natoms;
178     sh->step   = *step;
179     sh->nre    = 0;
180     sh->t      = *t;
181     sh->lambda = *lambda;
182   }
183   if (!do_trnheader(fio,bRead,sh,&bOK))
184     return FALSE;
185   if (bRead) {
186     *natoms = sh->natoms;
187     *step   = sh->step;
188     *t      = sh->t;
189     *lambda = sh->lambda;
190     if (sh->ir_size)
191       gmx_file("inputrec in trn file");
192     if (sh->e_size)
193       gmx_file("energies in trn file");
194     if (sh->top_size)
195       gmx_file("topology in trn file");
196     if (sh->sym_size)
197       gmx_file("symbol table in trn file");
198   }
199   bOK = do_htrn(fio,bRead,sh,box,x,v,f);
200
201   sfree(sh);
202   
203   return bOK;
204 }
205
206 /************************************************************
207  *
208  *  The following routines are the exported ones
209  *
210  ************************************************************/
211  
212 void read_trnheader(const char *fn,t_trnheader *trn)
213 {
214   t_fileio *fio;
215   gmx_bool bOK;
216   
217   fio = open_trn(fn,"r");
218   if (!do_trnheader(fio,TRUE,trn,&bOK))
219     gmx_fatal(FARGS,"Empty file %s",fn);
220   close_trn(fio);
221 }
222
223 gmx_bool fread_trnheader(t_fileio *fio,t_trnheader *trn, gmx_bool *bOK)
224 {
225   return do_trnheader(fio,TRUE,trn,bOK);
226 }
227
228 void write_trn(const char *fn,int step,real t,real lambda,
229                rvec *box,int natoms,rvec *x,rvec *v,rvec *f)
230 {
231   t_fileio *fio;
232   
233   fio = open_trn(fn,"w");
234   do_trn(fio,FALSE,&step,&t,&lambda,box,&natoms,x,v,f);
235   close_trn(fio);
236 }
237
238 void read_trn(const char *fn,int *step,real *t,real *lambda,
239               rvec *box,int *natoms,rvec *x,rvec *v,rvec *f)
240 {
241   t_fileio *fio;
242   
243   fio = open_trn(fn,"r");
244   (void) do_trn(fio,TRUE,step,t,lambda,box,natoms,x,v,f);
245   close_trn(fio);
246 }
247
248 void fwrite_trn(t_fileio *fio,int step,real t,real lambda,
249                 rvec *box,int natoms,rvec *x,rvec *v,rvec *f)
250 {
251   if( do_trn(fio,FALSE,&step,&t,&lambda,box,&natoms,x,v,f) == FALSE)
252   {
253       gmx_file("Cannot write trajectory frame; maybe you are out of disk space?");
254   }
255 }
256
257
258 gmx_bool fread_trn(t_fileio *fio,int *step,real *t,real *lambda,
259                rvec *box,int *natoms,rvec *x,rvec *v,rvec *f)
260 {
261   return do_trn(fio,TRUE,step,t,lambda,box,natoms,x,v,f);
262 }
263
264 gmx_bool fread_htrn(t_fileio *fio,t_trnheader *trn,rvec *box,rvec *x,rvec *v,
265                 rvec *f)
266 {
267   return do_htrn(fio,TRUE,trn,box,x,v,f);
268 }
269
270 t_fileio *open_trn(const char *fn,const char *mode)
271 {
272   return gmx_fio_open(fn,mode);
273 }
274
275 void close_trn(t_fileio *fio)
276 {
277   gmx_fio_close(fio);
278 }