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