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,2014, 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.
50 #include "gromacs/commandline/pargs.h"
51 #include "gromacs/fileio/enxio.h"
52 #include "gromacs/fileio/gmxfio.h"
53 #include "gromacs/fileio/mtxio.h"
54 #include "gromacs/fileio/tngio.h"
55 #include "gromacs/fileio/tngio_for_tools.h"
56 #include "gromacs/fileio/tpxio.h"
57 #include "gromacs/fileio/trnio.h"
58 #include "gromacs/fileio/xtcio.h"
59 #include "gromacs/gmxpreprocess/gmxcpp.h"
60 #include "gromacs/legacyheaders/checkpoint.h"
61 #include "gromacs/legacyheaders/macros.h"
62 #include "gromacs/legacyheaders/names.h"
63 #include "gromacs/legacyheaders/txtdump.h"
64 #include "gromacs/linearalgebra/sparsematrix.h"
65 #include "gromacs/topology/mtop_util.h"
66 #include "gromacs/utility/fatalerror.h"
67 #include "gromacs/utility/futil.h"
68 #include "gromacs/utility/smalloc.h"
70 static void list_tpx(const char *fn, gmx_bool bShowNumbers, const char *mdpfn,
74 int fp, indent, i, j, **gcount, atot;
83 read_tpxheader(fn, &tpx, TRUE, NULL, NULL);
87 &state, tpx.bF ? f : NULL,
88 tpx.bTop ? &mtop : NULL);
92 gp = gmx_fio_fopen(mdpfn, "w");
93 pr_inputrec(gp, 0, NULL, &(ir), TRUE);
101 top = gmx_mtop_t_to_t_topology(&mtop);
104 if (available(stdout, &tpx, 0, fn))
107 indent = pr_title(stdout, indent, fn);
108 pr_inputrec(stdout, 0, "inputrec", tpx.bIr ? &(ir) : NULL, FALSE);
111 pr_header(stdout, indent, "header", &(tpx));
115 pr_mtop(stdout, indent, "topology", &(mtop), bShowNumbers);
119 pr_top(stdout, indent, "topology", &(top), bShowNumbers);
122 pr_rvecs(stdout, indent, "box", tpx.bBox ? state.box : NULL, DIM);
123 pr_rvecs(stdout, indent, "box_rel", tpx.bBox ? state.box_rel : NULL, DIM);
124 pr_rvecs(stdout, indent, "boxv", tpx.bBox ? state.boxv : NULL, DIM);
125 pr_rvecs(stdout, indent, "pres_prev", tpx.bBox ? state.pres_prev : NULL, DIM);
126 pr_rvecs(stdout, indent, "svir_prev", tpx.bBox ? state.svir_prev : NULL, DIM);
127 pr_rvecs(stdout, indent, "fvir_prev", tpx.bBox ? state.fvir_prev : NULL, DIM);
128 /* leave nosehoover_xi in for now to match the tpr version */
129 pr_doubles(stdout, indent, "nosehoover_xi", state.nosehoover_xi, state.ngtc);
130 /*pr_doubles(stdout,indent,"nosehoover_vxi",state.nosehoover_vxi,state.ngtc);*/
131 /*pr_doubles(stdout,indent,"therm_integral",state.therm_integral,state.ngtc);*/
132 pr_rvecs(stdout, indent, "x", tpx.bX ? state.x : NULL, state.natoms);
133 pr_rvecs(stdout, indent, "v", tpx.bV ? state.v : NULL, state.natoms);
136 pr_rvecs(stdout, indent, "f", f, state.natoms);
140 groups = &mtop.groups;
143 for (i = 0; (i < egcNR); i++)
145 snew(gcount[i], groups->grps[i].nr);
148 for (i = 0; (i < mtop.natoms); i++)
150 for (j = 0; (j < egcNR); j++)
152 gcount[j][ggrpnr(groups, j, i)]++;
155 printf("Group statistics\n");
156 for (i = 0; (i < egcNR); i++)
159 printf("%-12s: ", gtypes[i]);
160 for (j = 0; (j < groups->grps[i].nr); j++)
162 printf(" %5d", gcount[i][j]);
163 atot += gcount[i][j];
165 printf(" (total %d atoms)\n", atot);
174 static void list_top(const char *fn)
180 char *cppopts[] = { NULL };
182 status = cpp_open_file(fn, &handle, cppopts);
185 gmx_fatal(FARGS, cpp_error(&handle, status));
189 status = cpp_read_line(&handle, BUFLEN, buf);
190 done = (status == eCPP_EOF);
193 if (status != eCPP_OK)
195 gmx_fatal(FARGS, cpp_error(&handle, status));
204 status = cpp_close_file(&handle);
205 if (status != eCPP_OK)
207 gmx_fatal(FARGS, cpp_error(&handle, status));
211 static void list_trn(const char *fn)
213 t_fileio *fpread, *fpwrite;
221 fpread = open_trn(fn, "r");
222 fpwrite = open_tpx(NULL, "w");
223 gmx_fio_setdebug(fpwrite, TRUE);
226 while (fread_trnheader(fpread, &trn, &bOK))
231 if (fread_htrn(fpread, &trn,
232 trn.box_size ? box : NULL,
233 trn.x_size ? x : NULL,
234 trn.v_size ? v : NULL,
235 trn.f_size ? f : NULL))
237 sprintf(buf, "%s frame %d", fn, nframe);
239 indent = pr_title(stdout, indent, buf);
240 pr_indent(stdout, indent);
241 fprintf(stdout, "natoms=%10d step=%10d time=%12.7e lambda=%10g\n",
242 trn.natoms, trn.step, trn.t, trn.lambda);
245 pr_rvecs(stdout, indent, "box", box, DIM);
249 pr_rvecs(stdout, indent, "x", x, trn.natoms);
253 pr_rvecs(stdout, indent, "v", v, trn.natoms);
257 pr_rvecs(stdout, indent, "f", f, trn.natoms);
262 fprintf(stderr, "\nWARNING: Incomplete frame: nr %d, t=%g\n",
273 fprintf(stderr, "\nWARNING: Incomplete frame header: nr %d, t=%g\n",
280 void list_xtc(const char *fn)
287 int nframe, natoms, step;
291 xd = open_xtc(fn, "r");
292 read_first_xtc(xd, &natoms, &step, &time, box, &x, &prec, &bOK);
297 sprintf(buf, "%s frame %d", fn, nframe);
299 indent = pr_title(stdout, indent, buf);
300 pr_indent(stdout, indent);
301 fprintf(stdout, "natoms=%10d step=%10d time=%12.7e prec=%10g\n",
302 natoms, step, time, prec);
303 pr_rvecs(stdout, indent, "box", box, DIM);
304 pr_rvecs(stdout, indent, "x", x, natoms);
307 while (read_next_xtc(xd, natoms, &step, &time, box, x, &prec, &bOK));
310 fprintf(stderr, "\nWARNING: Incomplete frame at time %g\n", time);
316 /*! \brief Callback used by list_tng_for_gmx_dump. */
317 static void list_tng_inner(const char *fn,
318 gmx_bool bFirstFrame,
322 gmx_int64_t n_values_per_frame,
333 sprintf(buf, "%s frame %" GMX_PRId64, fn, nframe);
335 indent = pr_title(stdout, indent, buf);
336 pr_indent(stdout, indent);
337 fprintf(stdout, "natoms=%10" GMX_PRId64 " step=%10" GMX_PRId64 " time=%12.7e",
338 n_atoms, step, frame_time);
341 fprintf(stdout, " prec=%10g", prec);
343 fprintf(stdout, "\n");
345 pr_reals_of_dim(stdout, indent, block_name, values, n_atoms, n_values_per_frame);
348 static void list_tng(const char gmx_unused *fn)
351 tng_trajectory_t tng;
352 gmx_int64_t nframe = 0;
353 gmx_int64_t i, *block_ids = NULL, step, ndatablocks;
356 gmx_tng_open(fn, 'r', &tng);
357 gmx_print_tng_molecule_system(tng, stdout);
359 bOK = gmx_get_tng_data_block_types_of_next_frame(tng, -1,
366 for (i = 0; i < ndatablocks; i++)
369 real prec, *values = NULL;
370 gmx_int64_t n_values_per_frame, n_atoms;
371 char block_name[STRLEN];
373 gmx_get_tng_data_next_frame_of_block_type(tng, block_ids[i], &values,
375 &n_values_per_frame, &n_atoms,
377 block_name, STRLEN, &bOK);
380 /* Can't write any output because we don't know what
382 fprintf(stderr, "\nWARNING: Incomplete frame at time %g, will not write output\n", frame_time);
386 list_tng_inner(fn, (0 == i), values, step, frame_time,
387 n_values_per_frame, n_atoms, prec, nframe, block_name);
392 while (gmx_get_tng_data_block_types_of_next_frame(tng, step,
408 void list_trx(const char *fn)
422 fprintf(stderr, "File %s is of an unsupported type. Try using the command\n 'less %s'\n",
427 void list_ene(const char *fn)
432 gmx_enxnm_t *enm = NULL;
438 printf("gmx dump: %s\n", fn);
439 in = open_enx(fn, "r");
440 do_enxnms(in, &nre, &enm);
443 printf("energy components:\n");
444 for (i = 0; (i < nre); i++)
446 printf("%5d %-24s (%s)\n", i, enm[i].name, enm[i].unit);
453 bCont = do_enx(in, fr);
457 printf("\n%24s %12.5e %12s %12s\n", "time:",
458 fr->t, "step:", gmx_step_str(fr->step, buf));
459 printf("%24s %12s %12s %12s\n",
460 "", "", "nsteps:", gmx_step_str(fr->nsteps, buf));
461 printf("%24s %12.5e %12s %12s\n",
462 "delta_t:", fr->dt, "sum steps:", gmx_step_str(fr->nsum, buf));
465 printf("%24s %12s %12s %12s\n",
466 "Component", "Energy", "Av. Energy", "Sum Energy");
469 for (i = 0; (i < nre); i++)
471 printf("%24s %12.5e %12.5e %12.5e\n",
472 enm[i].name, fr->ener[i].e, fr->ener[i].eav,
478 for (i = 0; (i < nre); i++)
480 printf("%24s %12.5e\n",
481 enm[i].name, fr->ener[i].e);
485 for (b = 0; b < fr->nblock; b++)
487 const char *typestr = "";
489 t_enxblock *eb = &(fr->block[b]);
490 printf("Block data %2d (%3d subblocks, id=%d)\n",
491 b, eb->nsub, eb->id);
495 typestr = enx_block_id_name[eb->id];
497 printf(" id='%s'\n", typestr);
498 for (i = 0; i < eb->nsub; i++)
500 t_enxsubblock *sb = &(eb->sub[i]);
501 printf(" Sub block %3d (%5d elems, type=%s) values:\n",
502 i, sb->nr, xdr_datatype_names[sb->type]);
506 case xdr_datatype_float:
507 for (j = 0; j < sb->nr; j++)
509 printf("%14d %8.4f\n", j, sb->fval[j]);
512 case xdr_datatype_double:
513 for (j = 0; j < sb->nr; j++)
515 printf("%14d %10.6f\n", j, sb->dval[j]);
518 case xdr_datatype_int:
519 for (j = 0; j < sb->nr; j++)
521 printf("%14d %10d\n", j, sb->ival[j]);
524 case xdr_datatype_int64:
525 for (j = 0; j < sb->nr; j++)
528 j, gmx_step_str(sb->lval[j], buf));
531 case xdr_datatype_char:
532 for (j = 0; j < sb->nr; j++)
534 printf("%14d %1c\n", j, sb->cval[j]);
537 case xdr_datatype_string:
538 for (j = 0; j < sb->nr; j++)
540 printf("%14d %80s\n", j, sb->sval[j]);
544 gmx_incons("Unknown subblock type");
559 static void list_mtx(const char *fn)
561 int nrow, ncol, i, j, k;
562 real *full = NULL, value;
563 gmx_sparsematrix_t * sparse = NULL;
565 gmx_mtxio_read(fn, &nrow, &ncol, &full, &sparse);
569 snew(full, nrow*ncol);
570 for (i = 0; i < nrow*ncol; i++)
575 for (i = 0; i < sparse->nrow; i++)
577 for (j = 0; j < sparse->ndata[i]; j++)
579 k = sparse->data[i][j].col;
580 value = sparse->data[i][j].value;
581 full[i*ncol+k] = value;
582 full[k*ncol+i] = value;
585 gmx_sparsematrix_destroy(sparse);
588 printf("%d %d\n", nrow, ncol);
589 for (i = 0; i < nrow; i++)
591 for (j = 0; j < ncol; j++)
593 printf(" %g", full[i*ncol+j]);
601 int gmx_dump(int argc, char *argv[])
603 const char *desc[] = {
604 "[THISMODULE] reads a run input file ([TT].tpr[tt]),",
605 "a trajectory ([TT].trr[tt]/[TT].xtc[tt]/[TT]/tng[tt]), an energy",
606 "file ([TT].edr[tt]) or a checkpoint file ([TT].cpt[tt])",
607 "and prints that to standard output in a readable format.",
608 "This program is essential for checking your run input file in case of",
610 "The program can also preprocess a topology to help finding problems.",
611 "Note that currently setting [TT]GMXLIB[tt] is the only way to customize",
612 "directories used for searching include files.",
614 const char *bugs[] = {
615 "Position restraint output from -sys -s is broken"
618 { efTPR, "-s", NULL, ffOPTRD },
619 { efTRX, "-f", NULL, ffOPTRD },
620 { efEDR, "-e", NULL, ffOPTRD },
621 { efCPT, NULL, NULL, ffOPTRD },
622 { efTOP, "-p", NULL, ffOPTRD },
623 { efMTX, "-mtx", "hessian", ffOPTRD },
624 { efMDP, "-om", NULL, ffOPTWR }
626 #define NFILE asize(fnm)
629 /* Command line options */
630 static gmx_bool bShowNumbers = TRUE;
631 static gmx_bool bSysTop = FALSE;
633 { "-nr", FALSE, etBOOL, {&bShowNumbers}, "Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" },
634 { "-sys", FALSE, etBOOL, {&bSysTop}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" }
637 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
638 asize(desc), desc, asize(bugs), bugs, &oenv))
644 if (ftp2bSet(efTPR, NFILE, fnm))
646 list_tpx(ftp2fn(efTPR, NFILE, fnm), bShowNumbers,
647 ftp2fn_null(efMDP, NFILE, fnm), bSysTop);
649 else if (ftp2bSet(efTRX, NFILE, fnm))
651 list_trx(ftp2fn(efTRX, NFILE, fnm));
653 else if (ftp2bSet(efEDR, NFILE, fnm))
655 list_ene(ftp2fn(efEDR, NFILE, fnm));
657 else if (ftp2bSet(efCPT, NFILE, fnm))
659 list_checkpoint(ftp2fn(efCPT, NFILE, fnm), stdout);
661 else if (ftp2bSet(efTOP, NFILE, fnm))
663 list_top(ftp2fn(efTOP, NFILE, fnm));
665 else if (ftp2bSet(efMTX, NFILE, fnm))
667 list_mtx(ftp2fn(efMTX, NFILE, fnm));