2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2013, The GROMACS development team.
6 * Copyright (c) 2013, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
47 #include "gromacs/fileio/futil.h"
51 #include "gmx_fatal.h"
52 #include "gromacs/fileio/xtcio.h"
53 #include "gromacs/fileio/enxio.h"
56 #include "gromacs/fileio/gmxfio.h"
57 #include "gromacs/fileio/tpxio.h"
58 #include "gromacs/fileio/trnio.h"
61 #include "checkpoint.h"
62 #include "mtop_util.h"
64 #include "gromacs/linearalgebra/mtxio.h"
65 #include "gromacs/linearalgebra/sparsematrix.h"
68 static void list_tpx(const char *fn, gmx_bool bShowNumbers, const char *mdpfn,
72 int fp, indent, i, j, **gcount, atot;
81 read_tpxheader(fn, &tpx, TRUE, NULL, NULL);
85 &state, tpx.bF ? f : NULL,
86 tpx.bTop ? &mtop : NULL);
90 gp = gmx_fio_fopen(mdpfn, "w");
91 pr_inputrec(gp, 0, NULL, &(ir), TRUE);
99 top = gmx_mtop_t_to_t_topology(&mtop);
102 if (available(stdout, &tpx, 0, fn))
105 indent = pr_title(stdout, indent, fn);
106 pr_inputrec(stdout, 0, "inputrec", tpx.bIr ? &(ir) : NULL, FALSE);
109 pr_header(stdout, indent, "header", &(tpx));
113 pr_mtop(stdout, indent, "topology", &(mtop), bShowNumbers);
117 pr_top(stdout, indent, "topology", &(top), bShowNumbers);
120 pr_rvecs(stdout, indent, "box", tpx.bBox ? state.box : NULL, DIM);
121 pr_rvecs(stdout, indent, "box_rel", tpx.bBox ? state.box_rel : NULL, DIM);
122 pr_rvecs(stdout, indent, "boxv", tpx.bBox ? state.boxv : NULL, DIM);
123 pr_rvecs(stdout, indent, "pres_prev", tpx.bBox ? state.pres_prev : NULL, DIM);
124 pr_rvecs(stdout, indent, "svir_prev", tpx.bBox ? state.svir_prev : NULL, DIM);
125 pr_rvecs(stdout, indent, "fvir_prev", tpx.bBox ? state.fvir_prev : NULL, DIM);
126 /* leave nosehoover_xi in for now to match the tpr version */
127 pr_doubles(stdout, indent, "nosehoover_xi", state.nosehoover_xi, state.ngtc);
128 /*pr_doubles(stdout,indent,"nosehoover_vxi",state.nosehoover_vxi,state.ngtc);*/
129 /*pr_doubles(stdout,indent,"therm_integral",state.therm_integral,state.ngtc);*/
130 pr_rvecs(stdout, indent, "x", tpx.bX ? state.x : NULL, state.natoms);
131 pr_rvecs(stdout, indent, "v", tpx.bV ? state.v : NULL, state.natoms);
134 pr_rvecs(stdout, indent, "f", f, state.natoms);
138 groups = &mtop.groups;
141 for (i = 0; (i < egcNR); i++)
143 snew(gcount[i], groups->grps[i].nr);
146 for (i = 0; (i < mtop.natoms); i++)
148 for (j = 0; (j < egcNR); j++)
150 gcount[j][ggrpnr(groups, j, i)]++;
153 printf("Group statistics\n");
154 for (i = 0; (i < egcNR); i++)
157 printf("%-12s: ", gtypes[i]);
158 for (j = 0; (j < groups->grps[i].nr); j++)
160 printf(" %5d", gcount[i][j]);
161 atot += gcount[i][j];
163 printf(" (total %d atoms)\n", atot);
172 static void list_top(const char *fn)
178 char *cppopts[] = { NULL };
180 status = cpp_open_file(fn, &handle, cppopts);
183 gmx_fatal(FARGS, cpp_error(&handle, status));
187 status = cpp_read_line(&handle, BUFLEN, buf);
188 done = (status == eCPP_EOF);
191 if (status != eCPP_OK)
193 gmx_fatal(FARGS, cpp_error(&handle, status));
202 status = cpp_close_file(&handle);
203 if (status != eCPP_OK)
205 gmx_fatal(FARGS, cpp_error(&handle, status));
209 static void list_trn(const char *fn)
211 t_fileio *fpread, *fpwrite;
219 fpread = open_trn(fn, "r");
220 fpwrite = open_tpx(NULL, "w");
221 gmx_fio_setdebug(fpwrite, TRUE);
224 while (fread_trnheader(fpread, &trn, &bOK))
229 if (fread_htrn(fpread, &trn,
230 trn.box_size ? box : NULL,
231 trn.x_size ? x : NULL,
232 trn.v_size ? v : NULL,
233 trn.f_size ? f : NULL))
235 sprintf(buf, "%s frame %d", fn, nframe);
237 indent = pr_title(stdout, indent, buf);
238 pr_indent(stdout, indent);
239 fprintf(stdout, "natoms=%10d step=%10d time=%12.7e lambda=%10g\n",
240 trn.natoms, trn.step, trn.t, trn.lambda);
243 pr_rvecs(stdout, indent, "box", box, DIM);
247 pr_rvecs(stdout, indent, "x", x, trn.natoms);
251 pr_rvecs(stdout, indent, "v", v, trn.natoms);
255 pr_rvecs(stdout, indent, "f", f, trn.natoms);
260 fprintf(stderr, "\nWARNING: Incomplete frame: nr %d, t=%g\n",
271 fprintf(stderr, "\nWARNING: Incomplete frame header: nr %d, t=%g\n",
278 void list_xtc(const char *fn, gmx_bool bXVG)
285 int nframe, natoms, step;
289 xd = open_xtc(fn, "r");
290 read_first_xtc(xd, &natoms, &step, &time, box, &x, &prec, &bOK);
299 fprintf(stdout, "%g", time);
300 for (i = 0; i < natoms; i++)
302 for (d = 0; d < DIM; d++)
304 fprintf(stdout, " %g", x[i][d]);
307 fprintf(stdout, "\n");
311 sprintf(buf, "%s frame %d", fn, nframe);
313 indent = pr_title(stdout, indent, buf);
314 pr_indent(stdout, indent);
315 fprintf(stdout, "natoms=%10d step=%10d time=%12.7e prec=%10g\n",
316 natoms, step, time, prec);
317 pr_rvecs(stdout, indent, "box", box, DIM);
318 pr_rvecs(stdout, indent, "x", x, natoms);
322 while (read_next_xtc(xd, natoms, &step, &time, box, x, &prec, &bOK));
325 fprintf(stderr, "\nWARNING: Incomplete frame at time %g\n", time);
331 void list_trx(const char *fn, gmx_bool bXVG)
340 else if ((ftp == efTRR) || (ftp == efTRJ))
346 fprintf(stderr, "File %s is of an unsupported type. Try using the command\n 'less %s'\n",
351 void list_ene(const char *fn)
356 gmx_enxnm_t *enm = NULL;
362 printf("gmx dump: %s\n", fn);
363 in = open_enx(fn, "r");
364 do_enxnms(in, &nre, &enm);
367 printf("energy components:\n");
368 for (i = 0; (i < nre); i++)
370 printf("%5d %-24s (%s)\n", i, enm[i].name, enm[i].unit);
377 bCont = do_enx(in, fr);
381 printf("\n%24s %12.5e %12s %12s\n", "time:",
382 fr->t, "step:", gmx_step_str(fr->step, buf));
383 printf("%24s %12s %12s %12s\n",
384 "", "", "nsteps:", gmx_step_str(fr->nsteps, buf));
385 printf("%24s %12.5e %12s %12s\n",
386 "delta_t:", fr->dt, "sum steps:", gmx_step_str(fr->nsum, buf));
389 printf("%24s %12s %12s %12s\n",
390 "Component", "Energy", "Av. Energy", "Sum Energy");
393 for (i = 0; (i < nre); i++)
395 printf("%24s %12.5e %12.5e %12.5e\n",
396 enm[i].name, fr->ener[i].e, fr->ener[i].eav,
402 for (i = 0; (i < nre); i++)
404 printf("%24s %12.5e\n",
405 enm[i].name, fr->ener[i].e);
409 for (b = 0; b < fr->nblock; b++)
411 const char *typestr = "";
413 t_enxblock *eb = &(fr->block[b]);
414 printf("Block data %2d (%3d subblocks, id=%d)\n",
415 b, eb->nsub, eb->id);
419 typestr = enx_block_id_name[eb->id];
421 printf(" id='%s'\n", typestr);
422 for (i = 0; i < eb->nsub; i++)
424 t_enxsubblock *sb = &(eb->sub[i]);
425 printf(" Sub block %3d (%5d elems, type=%s) values:\n",
426 i, sb->nr, xdr_datatype_names[sb->type]);
430 case xdr_datatype_float:
431 for (j = 0; j < sb->nr; j++)
433 printf("%14d %8.4f\n", j, sb->fval[j]);
436 case xdr_datatype_double:
437 for (j = 0; j < sb->nr; j++)
439 printf("%14d %10.6f\n", j, sb->dval[j]);
442 case xdr_datatype_int:
443 for (j = 0; j < sb->nr; j++)
445 printf("%14d %10d\n", j, sb->ival[j]);
448 case xdr_datatype_large_int:
449 for (j = 0; j < sb->nr; j++)
452 j, gmx_step_str(sb->lval[j], buf));
455 case xdr_datatype_char:
456 for (j = 0; j < sb->nr; j++)
458 printf("%14d %1c\n", j, sb->cval[j]);
461 case xdr_datatype_string:
462 for (j = 0; j < sb->nr; j++)
464 printf("%14d %80s\n", j, sb->sval[j]);
468 gmx_incons("Unknown subblock type");
483 static void list_mtx(const char *fn)
485 int nrow, ncol, i, j, k;
486 real *full = NULL, value;
487 gmx_sparsematrix_t * sparse = NULL;
489 gmx_mtxio_read(fn, &nrow, &ncol, &full, &sparse);
493 snew(full, nrow*ncol);
494 for (i = 0; i < nrow*ncol; i++)
499 for (i = 0; i < sparse->nrow; i++)
501 for (j = 0; j < sparse->ndata[i]; j++)
503 k = sparse->data[i][j].col;
504 value = sparse->data[i][j].value;
505 full[i*ncol+k] = value;
506 full[k*ncol+i] = value;
509 gmx_sparsematrix_destroy(sparse);
512 printf("%d %d\n", nrow, ncol);
513 for (i = 0; i < nrow; i++)
515 for (j = 0; j < ncol; j++)
517 printf(" %g", full[i*ncol+j]);
525 int gmx_gmxdump(int argc, char *argv[])
527 const char *desc[] = {
528 "[THISMODULE] reads a run input file ([TT].tpa[tt]/[TT].tpr[tt]/[TT].tpb[tt]),",
529 "a trajectory ([TT].trj[tt]/[TT].trr[tt]/[TT].xtc[tt]), an energy",
530 "file ([TT].ene[tt]/[TT].edr[tt]), or a checkpoint file ([TT].cpt[tt])",
531 "and prints that to standard output in a readable format.",
532 "This program is essential for checking your run input file in case of",
534 "The program can also preprocess a topology to help finding problems.",
535 "Note that currently setting [TT]GMXLIB[tt] is the only way to customize",
536 "directories used for searching include files.",
538 const char *bugs[] = {
539 "Position restraint output from -sys -s is broken"
542 { efTPX, "-s", NULL, ffOPTRD },
543 { efTRX, "-f", NULL, ffOPTRD },
544 { efEDR, "-e", NULL, ffOPTRD },
545 { efCPT, NULL, NULL, ffOPTRD },
546 { efTOP, "-p", NULL, ffOPTRD },
547 { efMTX, "-mtx", "hessian", ffOPTRD },
548 { efMDP, "-om", NULL, ffOPTWR }
550 #define NFILE asize(fnm)
553 /* Command line options */
554 static gmx_bool bXVG = FALSE;
555 static gmx_bool bShowNumbers = TRUE;
556 static gmx_bool bSysTop = FALSE;
558 { "-xvg", FALSE, etBOOL, {&bXVG}, "HIDDENXVG layout for xtc" },
559 { "-nr", FALSE, etBOOL, {&bShowNumbers}, "Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" },
560 { "-sys", FALSE, etBOOL, {&bSysTop}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" }
563 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
564 asize(desc), desc, asize(bugs), bugs, &oenv))
570 if (ftp2bSet(efTPX, NFILE, fnm))
572 list_tpx(ftp2fn(efTPX, NFILE, fnm), bShowNumbers,
573 ftp2fn_null(efMDP, NFILE, fnm), bSysTop);
575 else if (ftp2bSet(efTRX, NFILE, fnm))
577 list_trx(ftp2fn(efTRX, NFILE, fnm), bXVG);
579 else if (ftp2bSet(efEDR, NFILE, fnm))
581 list_ene(ftp2fn(efEDR, NFILE, fnm));
583 else if (ftp2bSet(efCPT, NFILE, fnm))
585 list_checkpoint(ftp2fn(efCPT, NFILE, fnm), stdout);
587 else if (ftp2bSet(efTOP, NFILE, fnm))
589 list_top(ftp2fn(efTOP, NFILE, fnm));
591 else if (ftp2bSet(efMTX, NFILE, fnm))
593 list_mtx(ftp2fn(efMTX, NFILE, fnm));