changed reading hint
[gromacs/adressmacs.git] / src / kernel / gmxdump.c
blobbdd56864c9cc1e16ab1c254e363fe61e98cb927c
1 /*
2 * $Id$
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 2.0
12 * Copyright (c) 1991-1999
13 * BIOSON Research Institute, Dept. of Biophysical Chemistry
14 * University of Groningen, The Netherlands
16 * Please refer to:
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
23 * or e-mail to:
24 * gromacs@chem.rug.nl
26 * And Hey:
27 * GRowing Old MAkes el Chrono Sweat
29 static char *SRCID_gmxdump_c = "$Id$";
31 #include <stdio.h>
32 #include <string.h>
33 #include <math.h>
34 #include "main.h"
35 #include "macros.h"
36 #include "futil.h"
37 #include "statutil.h"
38 #include "copyrite.h"
39 #include "sysstuff.h"
40 #include "txtdump.h"
41 #include "fatal.h"
42 #include "xtcio.h"
43 #include "enxio.h"
44 #include "assert.h"
45 #include "smalloc.h"
46 #include "names.h"
47 #include "gmxfio.h"
48 #include "tpxio.h"
49 #include "trnio.h"
50 #include "txtdump.h"
52 static void list_tpx(char *fn, bool bAltLayout)
54 int step,natoms,fp,indent,i,j,**gcount,atot;
55 real t,lambda;
56 rvec *x,*v,*f;
57 matrix box;
58 t_inputrec ir;
59 t_tpxheader tpx;
60 t_topology top;
62 read_tpxheader(fn,&tpx);
63 snew(x,tpx.natoms);
64 snew(v,tpx.natoms);
65 snew(f,tpx.natoms);
66 read_tpx(fn,&step,&t,&lambda,
67 tpx.bIr ? &ir : NULL,
68 tpx.bBox ? box : NULL,
69 &natoms,
70 tpx.bX ? x : NULL,
71 tpx.bV ? v : NULL,
72 tpx.bF ? f : NULL,
73 tpx.bTop ? &top: NULL);
75 if (bAltLayout) {
76 printf("----------- begin of %s ------------\n",fn);
77 fp = open_tpx(NULL,"w");
78 fio_setdebug(fp,TRUE);
79 fwrite_tpx(fp,step,t,lambda,
80 tpx.bIr ? &ir : NULL,
81 tpx.bBox ? box : NULL,
82 natoms,
83 tpx.bX ? x : NULL,
84 tpx.bV ? v : NULL,
85 tpx.bF ? f : NULL,
86 tpx.bTop ? &top : NULL);
87 close_tpx(fp);
88 printf("----------- end of %s ------------\n",fn);
89 } else {
90 if (available(stdout,&tpx,fn))
92 indent=0;
93 indent=pr_title(stdout,indent,fn);
94 pr_header(stdout,indent,"header",&(tpx));
95 pr_inputrec(stdout,indent,"ir",&(ir));
96 pr_rvecs(stdout,indent,"box",box,DIM);
97 pr_rvecs(stdout,indent,"x",x,natoms);
98 pr_rvecs(stdout,indent,"v",v,natoms);
99 pr_rvecs(stdout,indent,"f",f,natoms);
100 pr_top(stdout,indent,"topology",&(top));
103 snew(gcount,egcNR);
104 for(i=0; (i<egcNR); i++)
105 snew(gcount[i],top.atoms.grps[i].nr);
107 for(i=0; (i<top.atoms.nr); i++) {
108 for(j=0; (j<egcNR); j++)
109 gcount[j][top.atoms.atom[i].grpnr[j]]++;
111 printf("Group statistics\n");
112 for(i=0; (i<egcNR); i++) {
113 atot=0;
114 printf("%-12s: ",gtypes[i]);
115 for(j=0; (j<top.atoms.grps[i].nr); j++) {
116 printf(" %5d",gcount[i][j]);
117 atot+=gcount[i][j];
119 printf(" (total %d atoms)\n",atot);
120 sfree(gcount[i]);
122 sfree(gcount);
124 sfree(x);
125 sfree(v);
126 sfree(f);
129 static void list_trn(char *fn,bool bAltLayout)
131 int fpread,fpwrite,nframe,indent;
132 char buf[256];
133 rvec *x,*v,*f;
134 matrix box;
135 t_trnheader trn;
136 bool bOK;
138 fpread = open_trn(fn,"r");
139 fpwrite = open_tpx(NULL,"w");
140 fio_setdebug(fpwrite,TRUE);
142 nframe = 0;
143 while (fread_trnheader(fpread,&trn,&bOK)) {
144 snew(x,trn.natoms);
145 snew(v,trn.natoms);
146 snew(f,trn.natoms);
147 if (fread_htrn(fpread,&trn,
148 trn.box_size ? box : NULL,
149 trn.x_size ? x : NULL,
150 trn.v_size ? v : NULL,
151 trn.f_size ? f : NULL)) {
152 if (bAltLayout) {
153 printf("----------- begin of %s frame %d ------------\n",fn,nframe);
154 fwrite_tpx(fpwrite,trn.step,trn.t,trn.lambda,NULL,
155 trn.box_size ? box : NULL,
156 trn.natoms,
157 trn.x_size ? x : NULL,
158 trn.v_size ? v : NULL,
159 trn.f_size ? f : NULL,
160 NULL);
161 printf("----------- end of %s frame %d ------------\n",fn,nframe);
162 } else {
163 sprintf(buf,"%s frame %d",fn,nframe);
164 indent=0;
165 indent=pr_title(stdout,indent,buf);
166 pr_indent(stdout,indent);
167 fprintf(stdout,"natoms=%10d step=%10d time=%10g lambda=%10g\n",
168 trn.natoms,trn.step,trn.t,trn.lambda);
169 if (trn.box_size)
170 pr_rvecs(stdout,indent,"box",box,DIM);
171 if (trn.x_size)
172 pr_rvecs(stdout,indent,"x",x,trn.natoms);
173 if (trn.v_size)
174 pr_rvecs(stdout,indent,"v",v,trn.natoms);
175 if (trn.f_size)
176 pr_rvecs(stdout,indent,"f",f,trn.natoms);
178 } else
179 fprintf(stderr,"\nWARNING: Incomplete frame: nr %d, t=%g\n",
180 nframe,trn.t);
182 sfree(x);
183 sfree(v);
184 sfree(f);
185 nframe++;
187 if (!bOK)
188 fprintf(stderr,"\nWARNING: Incomplete frame header: nr %d, t=%g\n",
189 nframe,trn.t);
190 close_tpx(fpwrite);
191 close_trn(fpread);
194 void list_xtc(char *fn, bool bXVG)
196 int xd,indent;
197 char buf[256];
198 rvec *x;
199 matrix box;
200 int nframe,natoms,step;
201 real prec,time;
202 bool bOK;
204 xd = open_xtc(fn,"r");
205 read_first_xtc(xd,&natoms,&step,&time,box,&x,&prec,&bOK);
207 nframe=0;
208 do {
209 if (bXVG) {
210 int i,d;
212 fprintf(stdout,"%g",time);
213 for(i=0; i<natoms; i++)
214 for(d=0; d<DIM; d++)
215 fprintf(stdout," %g",x[i][d]);
216 fprintf(stdout,"\n");
217 } else {
218 sprintf(buf,"%s frame %d",fn,nframe);
219 indent=0;
220 indent=pr_title(stdout,indent,buf);
221 pr_indent(stdout,indent);
222 fprintf(stdout,"natoms=%10d step=%10d time=%10g prec=%10g\n",
223 natoms,step,time,prec);
224 pr_rvecs(stdout,indent,"box",box,DIM);
225 pr_rvecs(stdout,indent,"x",x,natoms);
227 nframe++;
228 } while (read_next_xtc(xd,&natoms,&step,&time,box,x,&prec,&bOK));
229 if (!bOK)
230 fprintf(stderr,"\nWARNING: Incomplete frame at time %g\n",time);
231 close_xtc(xd);
234 void list_trx(char *fn,bool bAltLayout,bool bXVG)
236 int ftp;
238 ftp = fn2ftp(fn);
239 if (ftp == efXTC)
240 list_xtc(fn,bXVG);
241 else if ((ftp == efTRR) || (ftp == efTRJ))
242 list_trn(fn,bAltLayout);
243 else
244 fprintf(stderr,"File %s not supported. Try using more %s\n",
245 fn,fn);
248 void list_ene(char *fn,bool bEDR)
250 int in,ndr;
251 bool bCont;
252 t_energy *ener;
253 t_drblock drblock;
254 int i,nre,nre2,step;
255 real t,rav,minthird;
256 char **enm=NULL;
258 printf("gmxdump: %s\n",fn);
259 in = open_enx(fn,"r");
260 do_enxnms(in,&nre,&enm);
262 printf("energy components:\n");
263 for(i=0; (i<nre); i++)
264 printf("%5d %s\n",i,enm[i]);
266 drblock.ndr=0;
267 minthird=-1.0/3.0;
268 snew(ener,nre);
269 do {
270 bCont=do_enx(in,&t,&step,&nre2,ener,&ndr,&drblock);
271 assert(nre2==nre);
273 if (bCont) {
274 printf("\n%24s %12.5e %12s %12d\n","time:",t,"step:",step);
275 printf("%24s %12s %12s %12s\n",
276 "Component","Energy","Av. Energy","Sum Energy");
277 for(i=0; (i<nre); i++)
278 printf("%24s %12.5e %12.5e %12.5e\n",
279 enm[i],ener[i].e,ener[i].eav,ener[i].esum);
280 if (ndr > 0) {
281 printf("Restraint %8s %8s\n","r(t)","< r >");
282 for(i=0; (i<drblock.ndr); i++) {
283 rav=pow(drblock.rav[i],minthird);
284 printf("%8d %8.4f %8.4f\n",i,drblock.rt[i],rav);
288 } while (bCont);
290 close_enx(in);
292 sfree(ener);
293 sfree(enm);
294 if (drblock.ndr) {
295 sfree(drblock.rav);
296 sfree(drblock.rt);
300 int main(int argc,char *argv[])
302 static char *desc[] = {
303 "gmxdump reads a run input file ([TT].tpa[tt]/[TT].tpr[tt]/[TT].tpb[tt]),",
304 "a trajectory ([TT].trj[tt]/[TT].trr[tt]/[TT].xtc[tt]) or an energy",
305 "file ([TT].ene[tt]/[TT].edr[tt]) and prints that to standard",
306 "output in a readable format. This program is essential for",
307 "checking your run input file in case of problems.[PAR]"
309 t_filenm fnm[] = {
310 { efTPX, "-s", NULL, ffOPTRD },
311 { efTRX, "-f", NULL, ffOPTRD },
312 { efENX, "-e", NULL, ffOPTRD },
314 #define NFILE asize(fnm)
316 /* Command line options */
317 static bool bAltLayout=FALSE,bXVG=FALSE;
318 static bool bShowNumbers=TRUE;
319 t_pargs pa[] = {
320 { "-a", FALSE, etBOOL, {&bAltLayout}, "HIDDENAlternative layout for run startup files" },
321 { "-xvg", FALSE, etBOOL, {&bXVG}, "HIDDENXVG layout for xtc" },
322 { "-nr",FALSE, etBOOL, {&bShowNumbers},"Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" }
325 CopyRight(stderr,argv[0]);
326 parse_common_args(&argc,argv,0,FALSE,NFILE,fnm,asize(pa),pa,
327 asize(desc),desc,0,NULL);
329 pr_shownumbers(bShowNumbers);
330 if (ftp2bSet(efTPX,NFILE,fnm))
331 list_tpx(ftp2fn(efTPX,NFILE,fnm), bAltLayout);
333 if (ftp2bSet(efTRX,NFILE,fnm))
334 list_trx(ftp2fn(efTRX,NFILE,fnm), bAltLayout, bXVG);
336 if (ftp2bSet(efENX,NFILE,fnm))
337 list_ene(ftp2fn(efENX,NFILE,fnm), FALSE);
339 thanx(stderr);
341 return 0;