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
49 #include "gmx_fatal.h"
59 #include "checkpoint.h"
60 #include "mtop_util.h"
61 #include "sparsematrix.h"
64 static void dump_top(FILE *fp,t_topology *top,char *tpr)
68 fprintf(fp,"; Topology generated from %s by program %s\n",tpr,Program());
69 fprintf(fp,"[ defaults ]\n 1 1 no 1.0 1.0\n\n");
70 fprintf(fp,"[ atomtypes ]\n");
71 fprintf(fp,";name at.num mass charge ptype c6 c12\n");
72 snew(types,top->atomtypes.nr);
73 for(i=0; (i<top->atomtypes.nr); i++) {
74 for(j=0; (j<top->atoms.nr) && (top->atoms.atom[j].type != i); j++)
76 if (j<top->atoms.nr) {
78 fprintf(fp,"%5s %4d %8.4f %8.4f %2s %8.3f %8.3f\n",
79 *(top->atoms.atomtype[j]),top->atomtypes.atomnumber[i],
83 fprintf(fp,"[ nonbonded_params ]\n");
84 for(i=k=0; (i<top->idef.ntypes); i++) {
85 for(j=0; (j<top->idef.ntypes); j++,k++) {
86 fprintf(fp,"%12s %12s 1 %12.5e %12.5e\n",
87 *(top->atoms.atomtype[types[i]]),
88 *(top->atoms.atomtype[types[j]]),
89 top->idef.iparams[k].lj.c12,top->idef.iparams[k].lj.c6);
95 static void list_tpx(const char *fn, gmx_bool bShowNumbers,const char *mdpfn,
99 int fp,indent,i,j,**gcount,atot;
105 gmx_groups_t *groups;
108 read_tpxheader(fn,&tpx,TRUE,NULL,NULL);
111 tpx.bIr ? &ir : NULL,
112 &state,tpx.bF ? f : NULL,
113 tpx.bTop ? &mtop: NULL);
115 if (mdpfn && tpx.bIr) {
116 gp = gmx_fio_fopen(mdpfn,"w");
117 pr_inputrec(gp,0,NULL,&(ir),TRUE);
123 top = gmx_mtop_t_to_t_topology(&mtop);
125 if (available(stdout,&tpx,0,fn)) {
127 indent=pr_title(stdout,indent,fn);
128 pr_inputrec(stdout,0,"inputrec",tpx.bIr ? &(ir) : NULL,FALSE);
131 pr_header(stdout,indent,"header",&(tpx));
134 pr_mtop(stdout,indent,"topology",&(mtop),bShowNumbers);
136 pr_top(stdout,indent,"topology",&(top),bShowNumbers);
138 pr_rvecs(stdout,indent,"box",tpx.bBox ? state.box : NULL,DIM);
139 pr_rvecs(stdout,indent,"box_rel",tpx.bBox ? state.box_rel : NULL,DIM);
140 pr_rvecs(stdout,indent,"boxv",tpx.bBox ? state.boxv : NULL,DIM);
141 pr_rvecs(stdout,indent,"pres_prev",tpx.bBox ? state.pres_prev : NULL,DIM);
142 pr_rvecs(stdout,indent,"svir_prev",tpx.bBox ? state.svir_prev : NULL,DIM);
143 pr_rvecs(stdout,indent,"fvir_prev",tpx.bBox ? state.fvir_prev : NULL,DIM);
144 /* leave nosehoover_xi in for now to match the tpr version */
145 pr_doubles(stdout,indent,"nosehoover_xi",state.nosehoover_xi,state.ngtc);
146 /*pr_doubles(stdout,indent,"nosehoover_vxi",state.nosehoover_vxi,state.ngtc);*/
147 /*pr_doubles(stdout,indent,"therm_integral",state.therm_integral,state.ngtc);*/
148 pr_rvecs(stdout,indent,"x",tpx.bX ? state.x : NULL,state.natoms);
149 pr_rvecs(stdout,indent,"v",tpx.bV ? state.v : NULL,state.natoms);
151 pr_rvecs(stdout,indent,"f",f,state.natoms);
155 groups = &mtop.groups;
158 for(i=0; (i<egcNR); i++)
159 snew(gcount[i],groups->grps[i].nr);
161 for(i=0; (i<mtop.natoms); i++) {
162 for(j=0; (j<egcNR); j++)
163 gcount[j][ggrpnr(groups,j,i)]++;
165 printf("Group statistics\n");
166 for(i=0; (i<egcNR); i++) {
168 printf("%-12s: ",gtypes[i]);
169 for(j=0; (j<groups->grps[i].nr); j++) {
170 printf(" %5d",gcount[i][j]);
173 printf(" (total %d atoms)\n",atot);
182 static void list_top(const char *fn)
188 char *cppopts[] = { NULL };
190 status = cpp_open_file(fn,&handle,cppopts);
192 gmx_fatal(FARGS,cpp_error(&handle,status));
194 status = cpp_read_line(&handle,BUFLEN,buf);
195 done = (status == eCPP_EOF);
197 if (status != eCPP_OK)
198 gmx_fatal(FARGS,cpp_error(&handle,status));
203 status = cpp_close_file(&handle);
204 if (status != eCPP_OK)
205 gmx_fatal(FARGS,cpp_error(&handle,status));
208 static void list_trn(const char *fn)
210 t_fileio *fpread, *fpwrite;
218 fpread = open_trn(fn,"r");
219 fpwrite = open_tpx(NULL,"w");
220 gmx_fio_setdebug(fpwrite,TRUE);
223 while (fread_trnheader(fpread,&trn,&bOK)) {
227 if (fread_htrn(fpread,&trn,
228 trn.box_size ? box : NULL,
229 trn.x_size ? x : NULL,
230 trn.v_size ? v : NULL,
231 trn.f_size ? f : NULL)) {
232 sprintf(buf,"%s frame %d",fn,nframe);
234 indent=pr_title(stdout,indent,buf);
235 pr_indent(stdout,indent);
236 fprintf(stdout,"natoms=%10d step=%10d time=%12.7e lambda=%10g\n",
237 trn.natoms,trn.step,trn.t,trn.lambda);
239 pr_rvecs(stdout,indent,"box",box,DIM);
241 pr_rvecs(stdout,indent,"x",x,trn.natoms);
243 pr_rvecs(stdout,indent,"v",v,trn.natoms);
245 pr_rvecs(stdout,indent,"f",f,trn.natoms);
248 fprintf(stderr,"\nWARNING: Incomplete frame: nr %d, t=%g\n",
257 fprintf(stderr,"\nWARNING: Incomplete frame header: nr %d, t=%g\n",
263 void list_xtc(const char *fn, gmx_bool bXVG)
270 int nframe,natoms,step;
274 xd = open_xtc(fn,"r");
275 read_first_xtc(xd,&natoms,&step,&time,box,&x,&prec,&bOK);
282 fprintf(stdout,"%g",time);
283 for(i=0; i<natoms; i++)
285 fprintf(stdout," %g",x[i][d]);
286 fprintf(stdout,"\n");
288 sprintf(buf,"%s frame %d",fn,nframe);
290 indent=pr_title(stdout,indent,buf);
291 pr_indent(stdout,indent);
292 fprintf(stdout,"natoms=%10d step=%10d time=%12.7e prec=%10g\n",
293 natoms,step,time,prec);
294 pr_rvecs(stdout,indent,"box",box,DIM);
295 pr_rvecs(stdout,indent,"x",x,natoms);
298 } while (read_next_xtc(xd,natoms,&step,&time,box,x,&prec,&bOK));
300 fprintf(stderr,"\nWARNING: Incomplete frame at time %g\n",time);
304 void list_trx(const char *fn,gmx_bool bXVG)
311 else if ((ftp == efTRR) || (ftp == efTRJ))
314 fprintf(stderr,"File %s is of an unsupported type. Try using the command\n 'less %s'\n",
318 void list_ene(const char *fn)
323 gmx_enxnm_t *enm=NULL;
329 printf("gmxdump: %s\n",fn);
330 in = open_enx(fn,"r");
331 do_enxnms(in,&nre,&enm);
333 printf("energy components:\n");
334 for(i=0; (i<nre); i++)
335 printf("%5d %-24s (%s)\n",i,enm[i].name,enm[i].unit);
344 printf("\n%24s %12.5e %12s %12s\n","time:",
345 fr->t,"step:",gmx_step_str(fr->step,buf));
346 printf("%24s %12s %12s %12s\n",
347 "","","nsteps:",gmx_step_str(fr->nsteps,buf));
348 printf("%24s %12.5e %12s %12s\n",
349 "delta_t:",fr->dt,"sum steps:",gmx_step_str(fr->nsum,buf));
350 if (fr->nre == nre) {
351 printf("%24s %12s %12s %12s\n",
352 "Component","Energy","Av. Energy","Sum Energy");
354 for(i=0; (i<nre); i++)
355 printf("%24s %12.5e %12.5e %12.5e\n",
356 enm[i].name,fr->ener[i].e,fr->ener[i].eav,
359 for(i=0; (i<nre); i++)
360 printf("%24s %12.5e\n",
361 enm[i].name,fr->ener[i].e);
364 for(b=0; b<fr->nblock; b++)
366 const char *typestr="";
368 t_enxblock *eb=&(fr->block[b]);
369 printf("Block data %2d (%3d subblocks, id=%d)\n",
370 b, eb->nsub, eb->id);
373 typestr=enx_block_id_name[eb->id];
374 printf(" id='%s'\n", typestr);
375 for(i=0;i<eb->nsub;i++)
377 t_enxsubblock *sb=&(eb->sub[i]);
378 printf(" Sub block %3d (%5d elems, type=%s) values:\n",
379 i, sb->nr, xdr_datatype_names[sb->type]);
383 case xdr_datatype_float:
384 for(j=0;j<sb->nr;j++)
385 printf("%14d %8.4f\n",j, sb->fval[j]);
387 case xdr_datatype_double:
388 for(j=0;j<sb->nr;j++)
389 printf("%14d %10.6f\n",j, sb->dval[j]);
391 case xdr_datatype_int:
392 for(j=0;j<sb->nr;j++)
393 printf("%14d %10d\n",j, sb->ival[j]);
395 case xdr_datatype_large_int:
396 for(j=0;j<sb->nr;j++)
398 j, gmx_step_str(sb->lval[j],buf));
400 case xdr_datatype_char:
401 for(j=0;j<sb->nr;j++)
402 printf("%14d %1c\n",j, sb->cval[j]);
404 case xdr_datatype_string:
405 for(j=0;j<sb->nr;j++)
406 printf("%14d %80s\n",j, sb->sval[j]);
409 gmx_incons("Unknown subblock type");
423 static void list_mtx(const char *fn)
426 real *full=NULL,value;
427 gmx_sparsematrix_t * sparse=NULL;
429 gmx_mtxio_read(fn,&nrow,&ncol,&full,&sparse);
432 snew(full,nrow*ncol);
433 for(i=0;i<nrow*ncol;i++) {
437 for(i=0;i<sparse->nrow;i++) {
438 for(j=0;j<sparse->ndata[i];j++) {
439 k = sparse->data[i][j].col;
440 value = sparse->data[i][j].value;
441 full[i*ncol+k] = value;
442 full[k*ncol+i] = value;
445 gmx_sparsematrix_destroy(sparse);
448 printf("%d %d\n",nrow,ncol);
449 for(i=0; i<nrow; i++) {
450 for(j=0; j<ncol; j++) {
451 printf(" %g",full[i*ncol+j]);
459 int main(int argc,char *argv[])
461 const char *desc[] = {
462 "[TT]gmxdump[tt] reads a run input file ([TT].tpa[tt]/[TT].tpr[tt]/[TT].tpb[tt]),",
463 "a trajectory ([TT].trj[tt]/[TT].trr[tt]/[TT].xtc[tt]), an energy",
464 "file ([TT].ene[tt]/[TT].edr[tt]), or a checkpoint file ([TT].cpt[tt])",
465 "and prints that to standard output in a readable format.",
466 "This program is essential for checking your run input file in case of",
468 "The program can also preprocess a topology to help finding problems.",
469 "Note that currently setting [TT]GMXLIB[tt] is the only way to customize",
470 "directories used for searching include files.",
473 { efTPX, "-s", NULL, ffOPTRD },
474 { efTRX, "-f", NULL, ffOPTRD },
475 { efEDR, "-e", NULL, ffOPTRD },
476 { efCPT, NULL, NULL, ffOPTRD },
477 { efTOP, "-p", NULL, ffOPTRD },
478 { efMTX, "-mtx", "hessian", ffOPTRD },
479 { efMDP, "-om", NULL, ffOPTWR }
481 #define NFILE asize(fnm)
484 /* Command line options */
485 static gmx_bool bXVG=FALSE;
486 static gmx_bool bShowNumbers=TRUE;
487 static gmx_bool bSysTop=FALSE;
489 { "-xvg", FALSE, etBOOL, {&bXVG}, "HIDDENXVG layout for xtc" },
490 { "-nr",FALSE, etBOOL, {&bShowNumbers},"Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" },
491 { "-sys", FALSE, etBOOL, {&bSysTop}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" }
494 CopyRight(stderr,argv[0]);
495 parse_common_args(&argc,argv,0,NFILE,fnm,asize(pa),pa,
496 asize(desc),desc,0,NULL,&oenv);
499 if (ftp2bSet(efTPX,NFILE,fnm))
500 list_tpx(ftp2fn(efTPX,NFILE,fnm),bShowNumbers,
501 ftp2fn_null(efMDP,NFILE,fnm),bSysTop);
502 else if (ftp2bSet(efTRX,NFILE,fnm))
503 list_trx(ftp2fn(efTRX,NFILE,fnm),bXVG);
504 else if (ftp2bSet(efEDR,NFILE,fnm))
505 list_ene(ftp2fn(efEDR,NFILE,fnm));
506 else if (ftp2bSet(efCPT,NFILE,fnm))
507 list_checkpoint(ftp2fn(efCPT,NFILE,fnm),stdout);
508 else if (ftp2bSet(efTOP,NFILE,fnm))
509 list_top(ftp2fn(efTOP,NFILE,fnm));
510 else if (ftp2bSet(efMTX,NFILE,fnm))
511 list_mtx(ftp2fn(efMTX,NFILE,fnm));