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 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
50 #include "gmx_fatal.h"
60 #include "checkpoint.h"
61 #include "mtop_util.h"
63 #include "gromacs/linearalgebra/mtxio.h"
64 #include "gromacs/linearalgebra/sparsematrix.h"
67 static void list_tpx(const char *fn, gmx_bool bShowNumbers,const char *mdpfn,
71 int fp,indent,i,j,**gcount,atot;
80 read_tpxheader(fn,&tpx,TRUE,NULL,NULL);
84 &state,tpx.bF ? f : NULL,
85 tpx.bTop ? &mtop: NULL);
87 if (mdpfn && tpx.bIr) {
88 gp = gmx_fio_fopen(mdpfn,"w");
89 pr_inputrec(gp,0,NULL,&(ir),TRUE);
95 top = gmx_mtop_t_to_t_topology(&mtop);
97 if (available(stdout,&tpx,0,fn)) {
99 indent=pr_title(stdout,indent,fn);
100 pr_inputrec(stdout,0,"inputrec",tpx.bIr ? &(ir) : NULL,FALSE);
103 pr_header(stdout,indent,"header",&(tpx));
106 pr_mtop(stdout,indent,"topology",&(mtop),bShowNumbers);
108 pr_top(stdout,indent,"topology",&(top),bShowNumbers);
110 pr_rvecs(stdout,indent,"box",tpx.bBox ? state.box : NULL,DIM);
111 pr_rvecs(stdout,indent,"box_rel",tpx.bBox ? state.box_rel : NULL,DIM);
112 pr_rvecs(stdout,indent,"boxv",tpx.bBox ? state.boxv : NULL,DIM);
113 pr_rvecs(stdout,indent,"pres_prev",tpx.bBox ? state.pres_prev : NULL,DIM);
114 pr_rvecs(stdout,indent,"svir_prev",tpx.bBox ? state.svir_prev : NULL,DIM);
115 pr_rvecs(stdout,indent,"fvir_prev",tpx.bBox ? state.fvir_prev : NULL,DIM);
116 /* leave nosehoover_xi in for now to match the tpr version */
117 pr_doubles(stdout,indent,"nosehoover_xi",state.nosehoover_xi,state.ngtc);
118 /*pr_doubles(stdout,indent,"nosehoover_vxi",state.nosehoover_vxi,state.ngtc);*/
119 /*pr_doubles(stdout,indent,"therm_integral",state.therm_integral,state.ngtc);*/
120 pr_rvecs(stdout,indent,"x",tpx.bX ? state.x : NULL,state.natoms);
121 pr_rvecs(stdout,indent,"v",tpx.bV ? state.v : NULL,state.natoms);
123 pr_rvecs(stdout,indent,"f",f,state.natoms);
127 groups = &mtop.groups;
130 for(i=0; (i<egcNR); i++)
131 snew(gcount[i],groups->grps[i].nr);
133 for(i=0; (i<mtop.natoms); i++) {
134 for(j=0; (j<egcNR); j++)
135 gcount[j][ggrpnr(groups,j,i)]++;
137 printf("Group statistics\n");
138 for(i=0; (i<egcNR); i++) {
140 printf("%-12s: ",gtypes[i]);
141 for(j=0; (j<groups->grps[i].nr); j++) {
142 printf(" %5d",gcount[i][j]);
145 printf(" (total %d atoms)\n",atot);
154 static void list_top(const char *fn)
160 char *cppopts[] = { NULL };
162 status = cpp_open_file(fn,&handle,cppopts);
164 gmx_fatal(FARGS,cpp_error(&handle,status));
166 status = cpp_read_line(&handle,BUFLEN,buf);
167 done = (status == eCPP_EOF);
169 if (status != eCPP_OK)
170 gmx_fatal(FARGS,cpp_error(&handle,status));
175 status = cpp_close_file(&handle);
176 if (status != eCPP_OK)
177 gmx_fatal(FARGS,cpp_error(&handle,status));
180 static void list_trn(const char *fn)
182 t_fileio *fpread, *fpwrite;
190 fpread = open_trn(fn,"r");
191 fpwrite = open_tpx(NULL,"w");
192 gmx_fio_setdebug(fpwrite,TRUE);
195 while (fread_trnheader(fpread,&trn,&bOK)) {
199 if (fread_htrn(fpread,&trn,
200 trn.box_size ? box : NULL,
201 trn.x_size ? x : NULL,
202 trn.v_size ? v : NULL,
203 trn.f_size ? f : NULL)) {
204 sprintf(buf,"%s frame %d",fn,nframe);
206 indent=pr_title(stdout,indent,buf);
207 pr_indent(stdout,indent);
208 fprintf(stdout,"natoms=%10d step=%10d time=%12.7e lambda=%10g\n",
209 trn.natoms,trn.step,trn.t,trn.lambda);
211 pr_rvecs(stdout,indent,"box",box,DIM);
213 pr_rvecs(stdout,indent,"x",x,trn.natoms);
215 pr_rvecs(stdout,indent,"v",v,trn.natoms);
217 pr_rvecs(stdout,indent,"f",f,trn.natoms);
220 fprintf(stderr,"\nWARNING: Incomplete frame: nr %d, t=%g\n",
229 fprintf(stderr,"\nWARNING: Incomplete frame header: nr %d, t=%g\n",
235 void list_xtc(const char *fn, gmx_bool bXVG)
242 int nframe,natoms,step;
246 xd = open_xtc(fn,"r");
247 read_first_xtc(xd,&natoms,&step,&time,box,&x,&prec,&bOK);
254 fprintf(stdout,"%g",time);
255 for(i=0; i<natoms; i++)
257 fprintf(stdout," %g",x[i][d]);
258 fprintf(stdout,"\n");
260 sprintf(buf,"%s frame %d",fn,nframe);
262 indent=pr_title(stdout,indent,buf);
263 pr_indent(stdout,indent);
264 fprintf(stdout,"natoms=%10d step=%10d time=%12.7e prec=%10g\n",
265 natoms,step,time,prec);
266 pr_rvecs(stdout,indent,"box",box,DIM);
267 pr_rvecs(stdout,indent,"x",x,natoms);
270 } while (read_next_xtc(xd,natoms,&step,&time,box,x,&prec,&bOK));
272 fprintf(stderr,"\nWARNING: Incomplete frame at time %g\n",time);
276 void list_trx(const char *fn,gmx_bool bXVG)
283 else if ((ftp == efTRR) || (ftp == efTRJ))
286 fprintf(stderr,"File %s is of an unsupported type. Try using the command\n 'less %s'\n",
290 void list_ene(const char *fn)
295 gmx_enxnm_t *enm=NULL;
301 printf("gmxdump: %s\n",fn);
302 in = open_enx(fn,"r");
303 do_enxnms(in,&nre,&enm);
306 printf("energy components:\n");
307 for(i=0; (i<nre); i++)
308 printf("%5d %-24s (%s)\n",i,enm[i].name,enm[i].unit);
317 printf("\n%24s %12.5e %12s %12s\n","time:",
318 fr->t,"step:",gmx_step_str(fr->step,buf));
319 printf("%24s %12s %12s %12s\n",
320 "","","nsteps:",gmx_step_str(fr->nsteps,buf));
321 printf("%24s %12.5e %12s %12s\n",
322 "delta_t:",fr->dt,"sum steps:",gmx_step_str(fr->nsum,buf));
323 if (fr->nre == nre) {
324 printf("%24s %12s %12s %12s\n",
325 "Component","Energy","Av. Energy","Sum Energy");
327 for(i=0; (i<nre); i++)
328 printf("%24s %12.5e %12.5e %12.5e\n",
329 enm[i].name,fr->ener[i].e,fr->ener[i].eav,
332 for(i=0; (i<nre); i++)
333 printf("%24s %12.5e\n",
334 enm[i].name,fr->ener[i].e);
337 for(b=0; b<fr->nblock; b++)
339 const char *typestr="";
341 t_enxblock *eb=&(fr->block[b]);
342 printf("Block data %2d (%3d subblocks, id=%d)\n",
343 b, eb->nsub, eb->id);
346 typestr=enx_block_id_name[eb->id];
347 printf(" id='%s'\n", typestr);
348 for(i=0;i<eb->nsub;i++)
350 t_enxsubblock *sb=&(eb->sub[i]);
351 printf(" Sub block %3d (%5d elems, type=%s) values:\n",
352 i, sb->nr, xdr_datatype_names[sb->type]);
356 case xdr_datatype_float:
357 for(j=0;j<sb->nr;j++)
358 printf("%14d %8.4f\n",j, sb->fval[j]);
360 case xdr_datatype_double:
361 for(j=0;j<sb->nr;j++)
362 printf("%14d %10.6f\n",j, sb->dval[j]);
364 case xdr_datatype_int:
365 for(j=0;j<sb->nr;j++)
366 printf("%14d %10d\n",j, sb->ival[j]);
368 case xdr_datatype_large_int:
369 for(j=0;j<sb->nr;j++)
371 j, gmx_step_str(sb->lval[j],buf));
373 case xdr_datatype_char:
374 for(j=0;j<sb->nr;j++)
375 printf("%14d %1c\n",j, sb->cval[j]);
377 case xdr_datatype_string:
378 for(j=0;j<sb->nr;j++)
379 printf("%14d %80s\n",j, sb->sval[j]);
382 gmx_incons("Unknown subblock type");
396 static void list_mtx(const char *fn)
399 real *full=NULL,value;
400 gmx_sparsematrix_t * sparse=NULL;
402 gmx_mtxio_read(fn,&nrow,&ncol,&full,&sparse);
405 snew(full,nrow*ncol);
406 for(i=0;i<nrow*ncol;i++) {
410 for(i=0;i<sparse->nrow;i++) {
411 for(j=0;j<sparse->ndata[i];j++) {
412 k = sparse->data[i][j].col;
413 value = sparse->data[i][j].value;
414 full[i*ncol+k] = value;
415 full[k*ncol+i] = value;
418 gmx_sparsematrix_destroy(sparse);
421 printf("%d %d\n",nrow,ncol);
422 for(i=0; i<nrow; i++) {
423 for(j=0; j<ncol; j++) {
424 printf(" %g",full[i*ncol+j]);
432 int cmain(int argc,char *argv[])
434 const char *desc[] = {
435 "[TT]gmxdump[tt] reads a run input file ([TT].tpa[tt]/[TT].tpr[tt]/[TT].tpb[tt]),",
436 "a trajectory ([TT].trj[tt]/[TT].trr[tt]/[TT].xtc[tt]), an energy",
437 "file ([TT].ene[tt]/[TT].edr[tt]), or a checkpoint file ([TT].cpt[tt])",
438 "and prints that to standard output in a readable format.",
439 "This program is essential for checking your run input file in case of",
441 "The program can also preprocess a topology to help finding problems.",
442 "Note that currently setting [TT]GMXLIB[tt] is the only way to customize",
443 "directories used for searching include files.",
446 { efTPX, "-s", NULL, ffOPTRD },
447 { efTRX, "-f", NULL, ffOPTRD },
448 { efEDR, "-e", NULL, ffOPTRD },
449 { efCPT, NULL, NULL, ffOPTRD },
450 { efTOP, "-p", NULL, ffOPTRD },
451 { efMTX, "-mtx", "hessian", ffOPTRD },
452 { efMDP, "-om", NULL, ffOPTWR }
454 #define NFILE asize(fnm)
457 /* Command line options */
458 static gmx_bool bXVG=FALSE;
459 static gmx_bool bShowNumbers=TRUE;
460 static gmx_bool bSysTop=FALSE;
462 { "-xvg", FALSE, etBOOL, {&bXVG}, "HIDDENXVG layout for xtc" },
463 { "-nr",FALSE, etBOOL, {&bShowNumbers},"Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" },
464 { "-sys", FALSE, etBOOL, {&bSysTop}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" }
467 CopyRight(stderr,argv[0]);
468 parse_common_args(&argc,argv,0,NFILE,fnm,asize(pa),pa,
469 asize(desc),desc,0,NULL,&oenv);
472 if (ftp2bSet(efTPX,NFILE,fnm))
473 list_tpx(ftp2fn(efTPX,NFILE,fnm),bShowNumbers,
474 ftp2fn_null(efMDP,NFILE,fnm),bSysTop);
475 else if (ftp2bSet(efTRX,NFILE,fnm))
476 list_trx(ftp2fn(efTRX,NFILE,fnm),bXVG);
477 else if (ftp2bSet(efEDR,NFILE,fnm))
478 list_ene(ftp2fn(efEDR,NFILE,fnm));
479 else if (ftp2bSet(efCPT,NFILE,fnm))
480 list_checkpoint(ftp2fn(efCPT,NFILE,fnm),stdout);
481 else if (ftp2bSet(efTOP,NFILE,fnm))
482 list_top(ftp2fn(efTOP,NFILE,fnm));
483 else if (ftp2bSet(efMTX,NFILE,fnm))
484 list_mtx(ftp2fn(efMTX,NFILE,fnm));