Revised wording in pdb2gmx.c, hopefully clearer now.
[gromacs/rigid-bodies.git] / src / gmxlib / trnio.c
blobea72d36f706241b90a6ad5c76b7ed6c0b3574b7a
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * GROningen Mixture of Alchemy and Childrens' Stories
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
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"
49 #define BUFSIZE 128
50 #define GROMACS_MAGIC 1993
52 static int nFloatSize(t_trnheader *sh)
54 int nflsize=0;
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");
67 if (((nflsize != sizeof(float)) && (nflsize != sizeof(double))))
68 gmx_fatal(FARGS,"Float size %d. Maybe different CPU?",nflsize);
70 return nflsize;
73 static gmx_bool do_trnheader(t_fileio *fio,gmx_bool bRead,t_trnheader *sh, gmx_bool *bOK)
75 int magic=GROMACS_MAGIC;
76 static gmx_bool bFirst=TRUE;
77 char buf[256];
79 *bOK=TRUE;
81 gmx_fio_checktype(fio);
83 if (!gmx_fio_do_int(fio,magic) || magic!=GROMACS_MAGIC)
84 return FALSE;
86 if (bRead) {
87 *bOK = *bOK && gmx_fio_do_string(fio,buf);
88 if (bFirst)
89 fprintf(stderr,"trn version: %s ",buf);
91 else {
92 sprintf(buf,"GMX_trn_file");
93 *bOK = *bOK && gmx_fio_do_string(fio,buf);
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);
107 if (!*bOK) return *bOK;
108 sh->bDouble = (nFloatSize(sh) == sizeof(double));
109 gmx_fio_setprecision(fio,sh->bDouble);
111 if (bRead && bFirst) {
112 fprintf(stderr,"(%s precision)\n",sh->bDouble ? "double" : "single");
113 bFirst = FALSE;
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);
121 return *bOK;
124 void pr_trnheader(FILE *fp,int indent,char *title,t_trnheader *sh)
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);
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);
148 static gmx_bool do_htrn(t_fileio *fio,gmx_bool bRead,t_trnheader *sh,
149 rvec *box,rvec *x,rvec *v,rvec *f)
151 matrix pv;
152 gmx_bool bOK;
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);
162 return bOK;
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)
168 t_trnheader *sh;
169 gmx_bool bOK;
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;
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");
199 bOK = do_htrn(fio,bRead,sh,box,x,v,f);
201 sfree(sh);
203 return bOK;
206 /************************************************************
208 * The following routines are the exported ones
210 ************************************************************/
212 void read_trnheader(const char *fn,t_trnheader *trn)
214 t_fileio *fio;
215 gmx_bool bOK;
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);
223 gmx_bool fread_trnheader(t_fileio *fio,t_trnheader *trn, gmx_bool *bOK)
225 return do_trnheader(fio,TRUE,trn,bOK);
228 void write_trn(const char *fn,int step,real t,real lambda,
229 rvec *box,int natoms,rvec *x,rvec *v,rvec *f)
231 t_fileio *fio;
233 fio = open_trn(fn,"w");
234 do_trn(fio,FALSE,&step,&t,&lambda,box,&natoms,x,v,f);
235 close_trn(fio);
238 void read_trn(const char *fn,int *step,real *t,real *lambda,
239 rvec *box,int *natoms,rvec *x,rvec *v,rvec *f)
241 t_fileio *fio;
243 fio = open_trn(fn,"r");
244 (void) do_trn(fio,TRUE,step,t,lambda,box,natoms,x,v,f);
245 close_trn(fio);
248 void fwrite_trn(t_fileio *fio,int step,real t,real lambda,
249 rvec *box,int natoms,rvec *x,rvec *v,rvec *f)
251 if( do_trn(fio,FALSE,&step,&t,&lambda,box,&natoms,x,v,f) == FALSE)
253 gmx_file("Cannot write trajectory frame; maybe you are out of quota?");
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)
261 return do_trn(fio,TRUE,step,t,lambda,box,natoms,x,v,f);
264 gmx_bool fread_htrn(t_fileio *fio,t_trnheader *trn,rvec *box,rvec *x,rvec *v,
265 rvec *f)
267 return do_htrn(fio,TRUE,trn,box,x,v,f);
270 t_fileio *open_trn(const char *fn,const char *mode)
272 return gmx_fio_open(fn,mode);
275 void close_trn(t_fileio *fio)
277 gmx_fio_close(fio);