3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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
33 * GROningen Mixture of Alchemy and Childrens' Stories
42 #include "gmx_fatal.h"
50 #define GROMACS_MAGIC 1993
52 static int nFloatSize(t_trnheader
*sh
)
57 nflsize
= sh
->box_size
/(DIM
*DIM
);
59 nflsize
= sh
->x_size
/(sh
->natoms
*DIM
);
61 nflsize
= sh
->v_size
/(sh
->natoms
*DIM
);
63 nflsize
= sh
->f_size
/(sh
->natoms
*DIM
);
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
);
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
;
81 gmx_fio_checktype(fio
);
83 if (!gmx_fio_do_int(fio
,magic
) || magic
!=GROMACS_MAGIC
)
87 *bOK
= *bOK
&& gmx_fio_do_string(fio
,buf
);
89 fprintf(stderr
,"trn version: %s ",buf
);
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");
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
);
124 void pr_trnheader(FILE *fp
,int indent
,char *title
,t_trnheader
*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
)
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
);
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
)
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
;
181 sh
->lambda
= *lambda
;
183 if (!do_trnheader(fio
,bRead
,sh
,&bOK
))
186 *natoms
= sh
->natoms
;
189 *lambda
= sh
->lambda
;
191 gmx_file("inputrec in trn file");
193 gmx_file("energies in trn file");
195 gmx_file("topology in trn file");
197 gmx_file("symbol table in trn file");
199 bOK
= do_htrn(fio
,bRead
,sh
,box
,x
,v
,f
);
206 /************************************************************
208 * The following routines are the exported ones
210 ************************************************************/
212 void read_trnheader(const char *fn
,t_trnheader
*trn
)
217 fio
= open_trn(fn
,"r");
218 if (!do_trnheader(fio
,TRUE
,trn
,&bOK
))
219 gmx_fatal(FARGS
,"Empty file %s",fn
);
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
)
233 fio
= open_trn(fn
,"w");
234 do_trn(fio
,FALSE
,&step
,&t
,&lambda
,box
,&natoms
,x
,v
,f
);
238 void read_trn(const char *fn
,int *step
,real
*t
,real
*lambda
,
239 rvec
*box
,int *natoms
,rvec
*x
,rvec
*v
,rvec
*f
)
243 fio
= open_trn(fn
,"r");
244 (void) do_trn(fio
,TRUE
,step
,t
,lambda
,box
,natoms
,x
,v
,f
);
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
,
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
)