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,2015,2016,2017,2018,2019, 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.
48 #include "gromacs/commandline/pargs.h"
49 #include "gromacs/fileio/checkpoint.h"
50 #include "gromacs/fileio/enxio.h"
51 #include "gromacs/fileio/gmxfio.h"
52 #include "gromacs/fileio/mtxio.h"
53 #include "gromacs/fileio/tngio.h"
54 #include "gromacs/fileio/tpxio.h"
55 #include "gromacs/fileio/trrio.h"
56 #include "gromacs/fileio/xtcio.h"
57 #include "gromacs/gmxpreprocess/gmxcpp.h"
58 #include "gromacs/linearalgebra/sparsematrix.h"
59 #include "gromacs/math/vecdump.h"
60 #include "gromacs/mdrunutility/mdmodules.h"
61 #include "gromacs/mdtypes/forcerec.h"
62 #include "gromacs/mdtypes/inputrec.h"
63 #include "gromacs/mdtypes/md_enums.h"
64 #include "gromacs/mdtypes/state.h"
65 #include "gromacs/topology/mtop_util.h"
66 #include "gromacs/topology/topology.h"
67 #include "gromacs/trajectory/energyframe.h"
68 #include "gromacs/trajectory/trajectoryframe.h"
69 #include "gromacs/utility/arraysize.h"
70 #include "gromacs/utility/basedefinitions.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/futil.h"
73 #include "gromacs/utility/smalloc.h"
74 #include "gromacs/utility/txtdump.h"
83 void list_tpx(const char *fn,
84 gmx_bool bShowNumbers,
85 gmx_bool bShowParameters,
88 gmx_bool bOriginalInputrec)
91 int indent, i, j, **gcount, atot;
97 read_tpxheader(fn, &tpx, TRUE);
100 tpx.bIr ? &ir : nullptr,
102 tpx.bTop ? &mtop : nullptr);
103 if (tpx.bIr && !bOriginalInputrec)
105 MDModules().adjustInputrecBasedOnModules(&ir);
108 if (mdpfn && tpx.bIr)
110 gp = gmx_fio_fopen(mdpfn, "w");
111 pr_inputrec(gp, 0, nullptr, &ir, TRUE);
119 top = gmx_mtop_t_to_t_topology(&mtop, false);
122 if (available(stdout, &tpx, 0, fn))
125 pr_title(stdout, indent, fn);
126 pr_inputrec(stdout, 0, "inputrec", tpx.bIr ? &ir : nullptr, FALSE);
128 pr_tpxheader(stdout, indent, "header", &(tpx));
132 pr_mtop(stdout, indent, "topology", &(mtop), bShowNumbers, bShowParameters);
136 pr_top(stdout, indent, "topology", &(top), bShowNumbers, bShowParameters);
139 pr_rvecs(stdout, indent, "box", tpx.bBox ? state.box : nullptr, DIM);
140 pr_rvecs(stdout, indent, "box_rel", tpx.bBox ? state.box_rel : nullptr, DIM);
141 pr_rvecs(stdout, indent, "boxv", tpx.bBox ? state.boxv : nullptr, DIM);
142 pr_rvecs(stdout, indent, "pres_prev", tpx.bBox ? state.pres_prev : nullptr, DIM);
143 pr_rvecs(stdout, indent, "svir_prev", tpx.bBox ? state.svir_prev : nullptr, DIM);
144 pr_rvecs(stdout, indent, "fvir_prev", tpx.bBox ? state.fvir_prev : nullptr, DIM);
145 /* leave nosehoover_xi in for now to match the tpr version */
146 pr_doubles(stdout, indent, "nosehoover_xi", state.nosehoover_xi.data(), state.ngtc);
147 /*pr_doubles(stdout,indent,"nosehoover_vxi",state.nosehoover_vxi,state.ngtc);*/
148 /*pr_doubles(stdout,indent,"therm_integral",state.therm_integral,state.ngtc);*/
149 pr_rvecs(stdout, indent, "x", tpx.bX ? state.x.rvec_array() : nullptr, state.natoms);
150 pr_rvecs(stdout, indent, "v", tpx.bV ? state.v.rvec_array() : nullptr, state.natoms);
153 const gmx_groups_t &groups = mtop.groups;
156 for (i = 0; (i < egcNR); i++)
158 snew(gcount[i], groups.grps[i].nr);
161 for (i = 0; (i < mtop.natoms); i++)
163 for (j = 0; (j < egcNR); j++)
165 gcount[j][getGroupType(groups, j, i)]++;
168 printf("Group statistics\n");
169 for (i = 0; (i < egcNR); i++)
172 printf("%-12s: ", gtypes[i]);
173 for (j = 0; (j < groups.grps[i].nr); j++)
175 printf(" %5d", gcount[i][j]);
176 atot += gcount[i][j];
178 printf(" (total %d atoms)\n", atot);
185 //! Dump a topology file
186 void list_top(const char *fn)
192 char *cppopts[] = { nullptr };
194 status = cpp_open_file(fn, &handle, cppopts);
197 gmx_fatal(FARGS, "%s", cpp_error(&handle, status));
201 status = cpp_read_line(&handle, BUFLEN, buf);
202 done = static_cast<int>(status == eCPP_EOF);
205 if (status != eCPP_OK)
207 gmx_fatal(FARGS, "%s", cpp_error(&handle, status));
216 status = cpp_close_file(&handle);
217 if (status != eCPP_OK)
219 gmx_fatal(FARGS, "%s", cpp_error(&handle, status));
224 void list_trr(const char *fn)
231 gmx_trr_header_t trrheader;
234 fpread = gmx_trr_open(fn, "r");
237 while (gmx_trr_read_frame_header(fpread, &trrheader, &bOK))
239 snew(x, trrheader.natoms);
240 snew(v, trrheader.natoms);
241 snew(f, trrheader.natoms);
242 if (gmx_trr_read_frame_data(fpread, &trrheader,
243 trrheader.box_size ? box : nullptr,
244 trrheader.x_size ? x : nullptr,
245 trrheader.v_size ? v : nullptr,
246 trrheader.f_size ? f : nullptr))
248 sprintf(buf, "%s frame %d", fn, nframe);
250 indent = pr_title(stdout, indent, buf);
251 pr_indent(stdout, indent);
252 fprintf(stdout, "natoms=%10d step=%10" PRId64 " time=%12.7e lambda=%10g\n",
253 trrheader.natoms, trrheader.step, trrheader.t, trrheader.lambda);
254 if (trrheader.box_size)
256 pr_rvecs(stdout, indent, "box", box, DIM);
258 if (trrheader.x_size)
260 pr_rvecs(stdout, indent, "x", x, trrheader.natoms);
262 if (trrheader.v_size)
264 pr_rvecs(stdout, indent, "v", v, trrheader.natoms);
266 if (trrheader.f_size)
268 pr_rvecs(stdout, indent, "f", f, trrheader.natoms);
273 fprintf(stderr, "\nWARNING: Incomplete frame: nr %d, t=%g\n",
274 nframe, trrheader.t);
284 fprintf(stderr, "\nWARNING: Incomplete frame header: nr %d, t=%g\n",
285 nframe, trrheader.t);
287 gmx_trr_close(fpread);
291 void list_xtc(const char *fn)
303 xd = open_xtc(fn, "r");
304 read_first_xtc(xd, &natoms, &step, &time, box, &x, &prec, &bOK);
309 sprintf(buf, "%s frame %d", fn, nframe);
311 indent = pr_title(stdout, indent, buf);
312 pr_indent(stdout, indent);
313 fprintf(stdout, "natoms=%10d step=%10" PRId64 " time=%12.7e prec=%10g\n",
314 natoms, step, time, prec);
315 pr_rvecs(stdout, indent, "box", box, DIM);
316 pr_rvecs(stdout, indent, "x", x, natoms);
319 while (read_next_xtc(xd, natoms, &step, &time, box, x, &prec, &bOK) != 0);
322 fprintf(stderr, "\nWARNING: Incomplete frame at time %g\n", time);
330 /*! \brief Callback used by list_tng_for_gmx_dump. */
331 void list_tng_inner(const char *fn,
332 gmx_bool bFirstFrame,
336 int64_t n_values_per_frame,
347 sprintf(buf, "%s frame %" PRId64, fn, nframe);
349 indent = pr_title(stdout, indent, buf);
350 pr_indent(stdout, indent);
351 fprintf(stdout, "natoms=%10" PRId64 " step=%10" PRId64 " time=%12.7e",
352 n_atoms, step, frame_time);
355 fprintf(stdout, " prec=%10g", prec);
357 fprintf(stdout, "\n");
359 pr_reals_of_dim(stdout, indent, block_name, values, n_atoms, n_values_per_frame);
365 void list_tng(const char *fn)
368 gmx_tng_trajectory_t tng;
370 int64_t i, *block_ids = nullptr, step, ndatablocks;
372 real *values = nullptr;
374 gmx_tng_open(fn, 'r', &tng);
375 gmx_print_tng_molecule_system(tng, stdout);
377 bOK = gmx_get_tng_data_block_types_of_next_frame(tng, -1,
384 for (i = 0; i < ndatablocks; i++)
388 int64_t n_values_per_frame, n_atoms;
389 char block_name[STRLEN];
391 gmx_get_tng_data_next_frame_of_block_type(tng, block_ids[i], &values,
393 &n_values_per_frame, &n_atoms,
395 block_name, STRLEN, &bOK);
398 /* Can't write any output because we don't know what
400 fprintf(stderr, "\nWARNING: Incomplete frame at time %g, will not write output\n", frame_time);
404 list_tng_inner(fn, (0 == i), values, step, frame_time,
405 n_values_per_frame, n_atoms, prec, nframe, block_name);
410 while (gmx_get_tng_data_block_types_of_next_frame(tng, step,
424 GMX_UNUSED_VALUE(fn);
428 //! Dump a trajectory file
429 void list_trx(const char *fn)
443 fprintf(stderr, "File %s is of an unsupported type. Try using the command\n 'less %s'\n",
448 //! Dump an energy file
449 void list_ene(const char *fn)
453 gmx_enxnm_t *enm = nullptr;
458 printf("gmx dump: %s\n", fn);
459 in = open_enx(fn, "r");
460 do_enxnms(in, &nre, &enm);
463 printf("energy components:\n");
464 for (i = 0; (i < nre); i++)
466 printf("%5d %-24s (%s)\n", i, enm[i].name, enm[i].unit);
472 bCont = do_enx(in, fr);
476 printf("\n%24s %12.5e %12s %12s\n", "time:",
477 fr->t, "step:", gmx_step_str(fr->step, buf));
478 printf("%24s %12s %12s %12s\n",
479 "", "", "nsteps:", gmx_step_str(fr->nsteps, buf));
480 printf("%24s %12.5e %12s %12s\n",
481 "delta_t:", fr->dt, "sum steps:", gmx_step_str(fr->nsum, buf));
484 printf("%24s %12s %12s %12s\n",
485 "Component", "Energy", "Av. Energy", "Sum Energy");
488 for (i = 0; (i < nre); i++)
490 printf("%24s %12.5e %12.5e %12.5e\n",
491 enm[i].name, fr->ener[i].e, fr->ener[i].eav,
497 for (i = 0; (i < nre); i++)
499 printf("%24s %12.5e\n",
500 enm[i].name, fr->ener[i].e);
504 for (b = 0; b < fr->nblock; b++)
506 const char *typestr = "";
508 t_enxblock *eb = &(fr->block[b]);
509 printf("Block data %2d (%3d subblocks, id=%d)\n",
510 b, eb->nsub, eb->id);
514 typestr = enx_block_id_name[eb->id];
516 printf(" id='%s'\n", typestr);
517 for (i = 0; i < eb->nsub; i++)
519 t_enxsubblock *sb = &(eb->sub[i]);
520 printf(" Sub block %3d (%5d elems, type=%s) values:\n",
521 i, sb->nr, xdr_datatype_names[sb->type]);
525 case xdr_datatype_float:
526 for (j = 0; j < sb->nr; j++)
528 printf("%14d %8.4f\n", j, sb->fval[j]);
531 case xdr_datatype_double:
532 for (j = 0; j < sb->nr; j++)
534 printf("%14d %10.6f\n", j, sb->dval[j]);
537 case xdr_datatype_int:
538 for (j = 0; j < sb->nr; j++)
540 printf("%14d %10d\n", j, sb->ival[j]);
543 case xdr_datatype_int64:
544 for (j = 0; j < sb->nr; j++)
547 j, gmx_step_str(sb->lval[j], buf));
550 case xdr_datatype_char:
551 for (j = 0; j < sb->nr; j++)
553 printf("%14d %1c\n", j, sb->cval[j]);
556 case xdr_datatype_string:
557 for (j = 0; j < sb->nr; j++)
559 printf("%14d %80s\n", j, sb->sval[j]);
563 gmx_incons("Unknown subblock type");
578 //! Dump a (Hessian) matrix file
579 void list_mtx(const char *fn)
581 int nrow, ncol, i, j, k;
582 real *full = nullptr, value;
583 gmx_sparsematrix_t * sparse = nullptr;
585 gmx_mtxio_read(fn, &nrow, &ncol, &full, &sparse);
589 snew(full, nrow*ncol);
590 for (i = 0; i < nrow*ncol; i++)
595 for (i = 0; i < sparse->nrow; i++)
597 for (j = 0; j < sparse->ndata[i]; j++)
599 k = sparse->data[i][j].col;
600 value = sparse->data[i][j].value;
601 full[i*ncol+k] = value;
602 full[k*ncol+i] = value;
605 gmx_sparsematrix_destroy(sparse);
608 printf("%d %d\n", nrow, ncol);
609 for (i = 0; i < nrow; i++)
611 for (j = 0; j < ncol; j++)
613 printf(" %g", full[i*ncol+j]);
623 int gmx_dump(int argc, char *argv[])
625 const char *desc[] = {
626 "[THISMODULE] reads a run input file ([REF].tpr[ref]),",
627 "a trajectory ([REF].trr[ref]/[REF].xtc[ref]/[TT]/tng[tt]), an energy",
628 "file ([REF].edr[ref]) or a checkpoint file ([REF].cpt[ref])",
629 "and prints that to standard output in a readable format.",
630 "This program is essential for checking your run input file in case of",
632 "The program can also preprocess a topology to help finding problems.",
633 "Note that currently setting [TT]GMXLIB[tt] is the only way to customize",
634 "directories used for searching include files.",
636 const char *bugs[] = {
637 "Position restraint output from -sys -s is broken"
640 { efTPR, "-s", nullptr, ffOPTRD },
641 { efTRX, "-f", nullptr, ffOPTRD },
642 { efEDR, "-e", nullptr, ffOPTRD },
643 { efCPT, nullptr, nullptr, ffOPTRD },
644 { efTOP, "-p", nullptr, ffOPTRD },
645 { efMTX, "-mtx", "hessian", ffOPTRD },
646 { efMDP, "-om", nullptr, ffOPTWR }
648 #define NFILE asize(fnm)
650 gmx_output_env_t *oenv;
651 /* Command line options */
652 gmx_bool bShowNumbers = TRUE;
653 gmx_bool bShowParams = FALSE;
654 gmx_bool bSysTop = FALSE;
655 gmx_bool bOriginalInputrec = FALSE;
657 { "-nr", FALSE, etBOOL, {&bShowNumbers}, "Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" },
658 { "-param", FALSE, etBOOL, {&bShowParams}, "Show parameters for each bonded interaction (for comparing dumps, it is useful to combine this with -nonr)" },
659 { "-sys", FALSE, etBOOL, {&bSysTop}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" },
660 { "-orgir", FALSE, etBOOL, {&bOriginalInputrec}, "Show input parameters from tpr as they were written by the version that produced the file, instead of how the current version reads them" }
663 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
664 asize(desc), desc, asize(bugs), bugs, &oenv))
670 if (ftp2bSet(efTPR, NFILE, fnm))
672 list_tpx(ftp2fn(efTPR, NFILE, fnm), bShowNumbers, bShowParams,
673 ftp2fn_null(efMDP, NFILE, fnm), bSysTop, bOriginalInputrec);
675 else if (ftp2bSet(efTRX, NFILE, fnm))
677 list_trx(ftp2fn(efTRX, NFILE, fnm));
679 else if (ftp2bSet(efEDR, NFILE, fnm))
681 list_ene(ftp2fn(efEDR, NFILE, fnm));
683 else if (ftp2bSet(efCPT, NFILE, fnm))
685 list_checkpoint(ftp2fn(efCPT, NFILE, fnm), stdout);
687 else if (ftp2bSet(efTOP, NFILE, fnm))
689 list_top(ftp2fn(efTOP, NFILE, fnm));
691 else if (ftp2bSet(efMTX, NFILE, fnm))
693 list_mtx(ftp2fn(efMTX, NFILE, fnm));