4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
12 * Copyright (c) 1991-1999
13 * BIOSON Research Institute, Dept. of Biophysical Chemistry
14 * University of Groningen, The Netherlands
17 * GROMACS: A message-passing parallel molecular dynamics implementation
18 * H.J.C. Berendsen, D. van der Spoel and R. van Drunen
19 * Comp. Phys. Comm. 91, 43-56 (1995)
21 * Also check out our WWW page:
22 * http://md.chem.rug.nl/~gmx
27 * Green Red Orange Magenta Azure Cyan Skyblue
29 static char *SRCID_trnio_c
= "$Id$";
43 #define GROMACS_MAGIC 1993
45 static int nFloatSize(t_trnheader
*sh
)
50 nflsize
= sh
->box_size
/(DIM
*DIM
);
52 nflsize
= sh
->x_size
/(sh
->natoms
*DIM
);
54 nflsize
= sh
->v_size
/(sh
->natoms
*DIM
);
56 nflsize
= sh
->f_size
/(sh
->natoms
*DIM
);
58 fatal_error(0,"Can not determine precision of trn file, quit!");
60 if (((nflsize
!= sizeof(float)) && (nflsize
!= sizeof(double))))
61 fatal_error(0,"Float size %d. Maybe different CPU?",nflsize
);
66 static bool do_trnheader(int fp
,bool bRead
,t_trnheader
*sh
, bool *bOK
)
68 static int magic
=GROMACS_MAGIC
;
69 static char *version
= "GMX_trn_file";
70 static bool bFirst
=TRUE
;
80 *bOK
= *bOK
&& do_string(buf
);
82 fprintf(stderr
,"trn version: %s\n",buf
);
87 *bOK
= *bOK
&& do_string(version
);
88 *bOK
= *bOK
&& do_int(sh
->ir_size
);
89 *bOK
= *bOK
&& do_int(sh
->e_size
);
90 *bOK
= *bOK
&& do_int(sh
->box_size
);
91 *bOK
= *bOK
&& do_int(sh
->vir_size
);
92 *bOK
= *bOK
&& do_int(sh
->pres_size
);
93 *bOK
= *bOK
&& do_int(sh
->top_size
);
94 *bOK
= *bOK
&& do_int(sh
->sym_size
);
95 *bOK
= *bOK
&& do_int(sh
->x_size
);
96 *bOK
= *bOK
&& do_int(sh
->v_size
);
97 *bOK
= *bOK
&& do_int(sh
->f_size
);
99 if (!*bOK
) return *bOK
;
100 fio_setprecision(fp
,(nFloatSize(sh
) == sizeof(double)));
102 *bOK
= *bOK
&& do_int(sh
->natoms
);
103 *bOK
= *bOK
&& do_int(sh
->step
);
104 *bOK
= *bOK
&& do_int(sh
->nre
);
105 *bOK
= *bOK
&& do_real(sh
->t
);
106 *bOK
= *bOK
&& do_real(sh
->lambda
);
111 void pr_trnheader(FILE *fp
,int indent
,char *title
,t_trnheader
*sh
)
114 indent
=pr_title(fp
,indent
,title
);
115 (void) pr_indent(fp
,indent
);
116 (void) fprintf(fp
,"box_size = %d\n",sh
->box_size
);
117 (void) pr_indent(fp
,indent
);
118 (void) fprintf(fp
,"x_size = %d\n",sh
->x_size
);
119 (void) pr_indent(fp
,indent
);
120 (void) fprintf(fp
,"v_size = %d\n",sh
->v_size
);
121 (void) pr_indent(fp
,indent
);
122 (void) fprintf(fp
,"f_size = %d\n",sh
->f_size
);
124 (void) pr_indent(fp
,indent
);
125 (void) fprintf(fp
,"natoms = %d\n",sh
->natoms
);
126 (void) pr_indent(fp
,indent
);
127 (void) fprintf(fp
,"step = %d\n",sh
->step
);
128 (void) pr_indent(fp
,indent
);
129 (void) fprintf(fp
,"t = %e\n",sh
->t
);
130 (void) pr_indent(fp
,indent
);
131 (void) fprintf(fp
,"lambda = %e\n",sh
->lambda
);
135 static bool do_htrn(int fp
,bool bRead
,t_trnheader
*sh
,
136 rvec
*box
,rvec
*x
,rvec
*v
,rvec
*f
)
142 if (sh
->box_size
!= 0) bOK
= bOK
&& ndo_rvec(box
,DIM
);
143 if (sh
->vir_size
!= 0) bOK
= bOK
&& ndo_rvec(pv
,DIM
);
144 if (sh
->pres_size
!= 0) bOK
= bOK
&& ndo_rvec(pv
,DIM
);
145 if (sh
->x_size
!= 0) bOK
= bOK
&& ndo_rvec(x
,sh
->natoms
);
146 if (sh
->v_size
!= 0) bOK
= bOK
&& ndo_rvec(v
,sh
->natoms
);
147 if (sh
->f_size
!= 0) bOK
= bOK
&& ndo_rvec(f
,sh
->natoms
);
152 static bool do_trn(int fp
,bool bRead
,int *step
,real
*t
,real
*lambda
,
153 rvec
*box
,int *natoms
,rvec
*x
,rvec
*v
,rvec
*f
)
160 sh
->box_size
=(box
)?sizeof(matrix
):0;
161 sh
->x_size
=((x
)?(*natoms
*sizeof(x
[0])):0);
162 sh
->v_size
=((v
)?(*natoms
*sizeof(v
[0])):0);
163 sh
->f_size
=((f
)?(*natoms
*sizeof(f
[0])):0);
164 sh
->natoms
= *natoms
;
168 sh
->lambda
= *lambda
;
170 if (!do_trnheader(fp
,bRead
,sh
,&bOK
))
173 *natoms
= sh
->natoms
;
176 *lambda
= sh
->lambda
;
178 fatal_error(0,"inputrec in trn file");
180 fatal_error(0,"energies in trn file");
182 fatal_error(0,"topology in trn file");
184 fatal_error(0,"symbol table in trn file");
186 bOK
= do_htrn(fp
,bRead
,sh
,box
,x
,v
,f
);
193 /************************************************************
195 * The following routines are the exported ones
197 ************************************************************/
199 void read_trnheader(char *fn
,t_trnheader
*trn
)
204 fp
= open_trn(fn
,"r");
205 if (!do_trnheader(fp
,TRUE
,trn
,&bOK
))
206 fatal_error(0,"Empty file %s",fn
);
210 bool fread_trnheader(int fp
,t_trnheader
*trn
, bool *bOK
)
212 return do_trnheader(fp
,TRUE
,trn
,bOK
);
215 void write_trn(char *fn
,int step
,real t
,real lambda
,
216 rvec
*box
,int natoms
,rvec
*x
,rvec
*v
,rvec
*f
)
220 fp
= open_trn(fn
,"w");
221 do_trn(fp
,FALSE
,&step
,&t
,&lambda
,box
,&natoms
,x
,v
,f
);
225 void read_trn(char *fn
,int *step
,real
*t
,real
*lambda
,
226 rvec
*box
,int *natoms
,rvec
*x
,rvec
*v
,rvec
*f
)
230 fp
= open_trn(fn
,"r");
231 (void) do_trn(fp
,TRUE
,step
,t
,lambda
,box
,natoms
,x
,v
,f
);
235 void fwrite_trn(int fp
,int step
,real t
,real lambda
,
236 rvec
*box
,int natoms
,rvec
*x
,rvec
*v
,rvec
*f
)
238 (void) do_trn(fp
,FALSE
,&step
,&t
,&lambda
,box
,&natoms
,x
,v
,f
);
242 bool fread_trn(int fp
,int *step
,real
*t
,real
*lambda
,
243 rvec
*box
,int *natoms
,rvec
*x
,rvec
*v
,rvec
*f
)
245 return do_trn(fp
,TRUE
,step
,t
,lambda
,box
,natoms
,x
,v
,f
);
248 bool fread_htrn(int fp
,t_trnheader
*trn
,rvec
*box
,rvec
*x
,rvec
*v
,rvec
*f
)
250 return do_htrn(fp
,TRUE
,trn
,box
,x
,v
,f
);
253 int open_trn(char *fn
,char *mode
)
255 return fio_open(fn
,mode
);
258 void close_trn(int fp
)