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.
38 * \brief Implements gmx dump utility.
40 * \ingroup module_tools
53 #include "gromacs/commandline/cmdlineoptionsmodule.h"
54 #include "gromacs/fileio/checkpoint.h"
55 #include "gromacs/fileio/enxio.h"
56 #include "gromacs/fileio/filetypes.h"
57 #include "gromacs/fileio/gmxfio.h"
58 #include "gromacs/fileio/mtxio.h"
59 #include "gromacs/fileio/tngio.h"
60 #include "gromacs/fileio/tpxio.h"
61 #include "gromacs/fileio/trrio.h"
62 #include "gromacs/fileio/xtcio.h"
63 #include "gromacs/gmxpreprocess/gmxcpp.h"
64 #include "gromacs/math/vecdump.h"
65 #include "gromacs/mdrun/mdmodules.h"
66 #include "gromacs/mdtypes/forcerec.h"
67 #include "gromacs/mdtypes/inputrec.h"
68 #include "gromacs/mdtypes/md_enums.h"
69 #include "gromacs/mdtypes/state.h"
70 #include "gromacs/options/basicoptions.h"
71 #include "gromacs/options/filenameoption.h"
72 #include "gromacs/options/ioptionscontainer.h"
73 #include "gromacs/topology/mtop_util.h"
74 #include "gromacs/topology/topology.h"
75 #include "gromacs/trajectory/energyframe.h"
76 #include "gromacs/trajectory/trajectoryframe.h"
77 #include "gromacs/utility/arraysize.h"
78 #include "gromacs/utility/basedefinitions.h"
79 #include "gromacs/utility/fatalerror.h"
80 #include "gromacs/utility/futil.h"
81 #include "gromacs/utility/smalloc.h"
82 #include "gromacs/utility/txtdump.h"
91 void list_tpr(const char* fn,
92 gmx_bool bShowNumbers,
93 gmx_bool bShowParameters,
96 gmx_bool bOriginalInputrec)
104 TpxFileHeader tpx = readTpxHeader(fn, true);
107 read_tpx_state(fn, tpx.bIr ? &ir : nullptr, &state, tpx.bTop ? &mtop : nullptr);
108 if (tpx.bIr && !bOriginalInputrec)
110 MDModules().adjustInputrecBasedOnModules(&ir);
113 if (mdpfn && tpx.bIr)
115 gp = gmx_fio_fopen(mdpfn, "w");
116 pr_inputrec(gp, 0, nullptr, &ir, TRUE);
124 top = gmx_mtop_t_to_t_topology(&mtop, false);
127 if (available(stdout, &tpx, 0, fn))
130 pr_title(stdout, indent, fn);
131 pr_inputrec(stdout, 0, "inputrec", tpx.bIr ? &ir : nullptr, FALSE);
133 pr_tpxheader(stdout, indent, "header", &(tpx));
137 pr_mtop(stdout, indent, "topology", &(mtop), bShowNumbers, bShowParameters);
141 pr_top(stdout, indent, "topology", &(top), bShowNumbers, bShowParameters);
144 pr_rvecs(stdout, indent, "box", tpx.bBox ? state.box : nullptr, DIM);
145 pr_rvecs(stdout, indent, "box_rel", tpx.bBox ? state.box_rel : nullptr, DIM);
146 pr_rvecs(stdout, indent, "boxv", tpx.bBox ? state.boxv : nullptr, DIM);
147 pr_rvecs(stdout, indent, "pres_prev", tpx.bBox ? state.pres_prev : nullptr, DIM);
148 pr_rvecs(stdout, indent, "svir_prev", tpx.bBox ? state.svir_prev : nullptr, DIM);
149 pr_rvecs(stdout, indent, "fvir_prev", tpx.bBox ? state.fvir_prev : nullptr, DIM);
150 /* leave nosehoover_xi in for now to match the tpr version */
151 pr_doubles(stdout, indent, "nosehoover_xi", state.nosehoover_xi.data(), state.ngtc);
152 /*pr_doubles(stdout,indent,"nosehoover_vxi",state.nosehoover_vxi,state.ngtc);*/
153 /*pr_doubles(stdout,indent,"therm_integral",state.therm_integral,state.ngtc);*/
154 pr_rvecs(stdout, indent, "x", tpx.bX ? state.x.rvec_array() : nullptr, state.natoms);
155 pr_rvecs(stdout, indent, "v", tpx.bV ? state.v.rvec_array() : nullptr, state.natoms);
158 const SimulationGroups& groups = mtop.groups;
160 gmx::EnumerationArray<SimulationAtomGroupType, std::vector<int>> gcount;
161 for (auto group : keysOf(gcount))
163 gcount[group].resize(groups.groups[group].size());
166 for (int i = 0; (i < mtop.natoms); i++)
168 for (auto group : keysOf(gcount))
170 gcount[group][getGroupType(groups, group, i)]++;
173 printf("Group statistics\n");
174 for (auto group : keysOf(gcount))
177 printf("%-12s: ", shortName(group));
178 for (const auto& entry : gcount[group])
180 printf(" %5d", entry);
183 printf(" (total %d atoms)\n", atot);
188 //! Dump a topology file
189 void list_top(const char* fn)
192 // Legacy string length macro
195 char* cppopts[] = { nullptr };
197 status = cpp_open_file(fn, &handle, cppopts);
200 gmx_fatal(FARGS, "%s", cpp_error(&handle, status));
204 status = cpp_read_line(&handle, STRLEN, buf);
205 done = static_cast<int>(status == eCPP_EOF);
208 if (status != eCPP_OK)
210 gmx_fatal(FARGS, "%s", cpp_error(&handle, status));
218 status = cpp_close_file(&handle);
219 if (status != eCPP_OK)
221 gmx_fatal(FARGS, "%s", cpp_error(&handle, status));
226 void list_trr(const char* fn)
233 gmx_trr_header_t trrheader;
236 fpread = gmx_trr_open(fn, "r");
239 while (gmx_trr_read_frame_header(fpread, &trrheader, &bOK))
241 snew(x, trrheader.natoms);
242 snew(v, trrheader.natoms);
243 snew(f, trrheader.natoms);
244 if (gmx_trr_read_frame_data(fpread, &trrheader, trrheader.box_size ? box : nullptr,
245 trrheader.x_size ? x : nullptr, 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", nframe, trrheader.t);
283 fprintf(stderr, "\nWARNING: Incomplete frame header: nr %d, t=%g\n", nframe, trrheader.t);
285 gmx_trr_close(fpread);
289 void list_xtc(const char* fn)
301 xd = open_xtc(fn, "r");
302 read_first_xtc(xd, &natoms, &step, &time, box, &x, &prec, &bOK);
307 sprintf(buf, "%s frame %d", fn, nframe);
309 indent = pr_title(stdout, indent, buf);
310 pr_indent(stdout, indent);
311 fprintf(stdout, "natoms=%10d step=%10" PRId64 " time=%12.7e prec=%10g\n", natoms, step,
313 pr_rvecs(stdout, indent, "box", box, DIM);
314 pr_rvecs(stdout, indent, "x", x, natoms);
316 } while (read_next_xtc(xd, natoms, &step, &time, box, x, &prec, &bOK) != 0);
319 fprintf(stderr, "\nWARNING: Incomplete frame at time %g\n", time);
327 /*! \brief Callback used by list_tng_for_gmx_dump. */
328 void list_tng_inner(const char* fn,
329 gmx_bool bFirstFrame,
333 int64_t n_values_per_frame,
344 sprintf(buf, "%s frame %" PRId64, fn, nframe);
346 indent = pr_title(stdout, indent, buf);
347 pr_indent(stdout, indent);
348 fprintf(stdout, "natoms=%10" PRId64 " step=%10" PRId64 " time=%12.7e", n_atoms, step, frame_time);
351 fprintf(stdout, " prec=%10g", prec);
353 fprintf(stdout, "\n");
355 pr_reals_of_dim(stdout, indent, block_name, values, n_atoms, n_values_per_frame);
361 void list_tng(const char* fn)
364 gmx_tng_trajectory_t tng;
366 int64_t i, *block_ids = nullptr, step, ndatablocks;
368 real* values = nullptr;
370 gmx_tng_open(fn, 'r', &tng);
371 gmx_print_tng_molecule_system(tng, stdout);
373 bOK = gmx_get_tng_data_block_types_of_next_frame(tng, -1, 0, nullptr, &step, &ndatablocks, &block_ids);
376 for (i = 0; i < ndatablocks; i++)
380 int64_t n_values_per_frame, n_atoms;
381 char block_name[STRLEN];
383 gmx_get_tng_data_next_frame_of_block_type(tng, block_ids[i], &values, &step,
384 &frame_time, &n_values_per_frame, &n_atoms,
385 &prec, block_name, STRLEN, &bOK);
388 /* Can't write any output because we don't know what
390 fprintf(stderr, "\nWARNING: Incomplete frame at time %g, will not write output\n",
395 list_tng_inner(fn, (0 == i), values, step, frame_time, n_values_per_frame, n_atoms,
396 prec, nframe, block_name);
400 } while (gmx_get_tng_data_block_types_of_next_frame(tng, step, 0, nullptr, &step, &ndatablocks,
410 GMX_UNUSED_VALUE(fn);
414 //! Dump a trajectory file
415 void list_trx(const char* fn)
419 case efXTC: list_xtc(fn); break;
420 case efTRR: list_trr(fn); break;
421 case efTNG: list_tng(fn); break;
423 fprintf(stderr, "File %s is of an unsupported type. Try using the command\n 'less %s'\n",
428 //! Dump an energy file
429 void list_ene(const char* fn)
433 gmx_enxnm_t* enm = nullptr;
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);
452 bCont = do_enx(in, fr);
456 printf("\n%24s %12.5e %12s %12s\n", "time:", fr->t, "step:", gmx_step_str(fr->step, buf));
457 printf("%24s %12s %12s %12s\n", "", "", "nsteps:", gmx_step_str(fr->nsteps, buf));
458 printf("%24s %12.5e %12s %12s\n", "delta_t:", fr->dt,
459 "sum steps:", gmx_step_str(fr->nsum, buf));
462 printf("%24s %12s %12s %12s\n", "Component", "Energy", "Av. Energy",
466 for (i = 0; (i < nre); i++)
468 printf("%24s %12.5e %12.5e %12.5e\n", enm[i].name, fr->ener[i].e,
469 fr->ener[i].eav, fr->ener[i].esum);
474 for (i = 0; (i < nre); i++)
476 printf("%24s %12.5e\n", enm[i].name, fr->ener[i].e);
480 for (b = 0; b < fr->nblock; b++)
482 const char* typestr = "";
484 t_enxblock* eb = &(fr->block[b]);
485 printf("Block data %2d (%3d subblocks, id=%d)\n", b, eb->nsub, eb->id);
489 typestr = enx_block_id_name[eb->id];
491 printf(" id='%s'\n", typestr);
492 for (i = 0; i < eb->nsub; i++)
494 t_enxsubblock* sb = &(eb->sub[i]);
495 printf(" Sub block %3d (%5d elems, type=%s) values:\n", i, sb->nr,
496 xdr_datatype_names[sb->type]);
500 case xdr_datatype_float:
501 for (j = 0; j < sb->nr; j++)
503 printf("%14d %8.4f\n", j, sb->fval[j]);
506 case xdr_datatype_double:
507 for (j = 0; j < sb->nr; j++)
509 printf("%14d %10.6f\n", j, sb->dval[j]);
512 case xdr_datatype_int:
513 for (j = 0; j < sb->nr; j++)
515 printf("%14d %10d\n", j, sb->ival[j]);
518 case xdr_datatype_int64:
519 for (j = 0; j < sb->nr; j++)
521 printf("%14d %s\n", j, gmx_step_str(sb->lval[j], buf));
524 case xdr_datatype_char:
525 for (j = 0; j < sb->nr; j++)
527 printf("%14d %1c\n", j, sb->cval[j]);
530 case xdr_datatype_string:
531 for (j = 0; j < sb->nr; j++)
533 printf("%14d %80s\n", j, sb->sval[j]);
536 default: gmx_incons("Unknown subblock type");
550 //! Dump a (Hessian) matrix file
551 void list_mtx(const char* fn)
553 int nrow, ncol, i, j, k;
554 real * full = nullptr, value;
555 gmx_sparsematrix_t* sparse = nullptr;
557 gmx_mtxio_read(fn, &nrow, &ncol, &full, &sparse);
561 snew(full, nrow * ncol);
562 for (i = 0; i < nrow * ncol; i++)
567 for (i = 0; i < sparse->nrow; i++)
569 for (j = 0; j < sparse->ndata[i]; j++)
571 k = sparse->data[i][j].col;
572 value = sparse->data[i][j].value;
573 full[i * ncol + k] = value;
574 full[k * ncol + i] = value;
577 gmx_sparsematrix_destroy(sparse);
580 printf("%d %d\n", nrow, ncol);
581 for (i = 0; i < nrow; i++)
583 for (j = 0; j < ncol; j++)
585 printf(" %g", full[i * ncol + j]);
593 class Dump : public ICommandLineOptionsModule
598 // From ICommandLineOptionsModule
599 void init(CommandLineModuleSettings* /*settings*/) override {}
601 void initOptions(IOptionsContainer* options, ICommandLineOptionsModuleSettings* settings) override;
603 void optionsFinished() override;
608 //! Commandline options
610 bool bShowNumbers_ = true;
611 bool bShowParams_ = false;
612 bool bSysTop_ = false;
613 bool bOriginalInputrec_ = false;
615 //! Commandline file options
617 std::string inputTprFilename_;
618 std::string inputTrajectoryFilename_;
619 std::string inputEnergyFilename_;
620 std::string inputCheckpointFilename_;
621 std::string inputTopologyFilename_;
622 std::string inputMatrixFilename_;
623 std::string outputMdpFilename_;
627 void Dump::initOptions(IOptionsContainer* options, ICommandLineOptionsModuleSettings* settings)
629 const char* desc[] = { "[THISMODULE] reads a run input file ([REF].tpr[ref]),",
630 "a trajectory ([REF].trr[ref]/[REF].xtc[ref]/[TT]tng[tt]), an energy",
631 "file ([REF].edr[ref]), a checkpoint file ([REF].cpt[ref])",
632 "or topology file ([REF].top[ref])",
633 "and prints that to standard output in a readable format.",
634 "This program is essential for checking your run input file in case of",
636 settings->setHelpText(desc);
638 const char* bugs[] = {
639 "The [REF].mdp[ref] file produced by [TT]-om[tt] can not be read by grompp."
641 settings->setBugText(bugs);
642 // TODO If this ancient note acknowledging a bug is still true,
643 // fix it or block that run path:
644 // Position restraint output from -sys -s is broken
647 FileNameOption("s").filetype(eftRunInput).inputFile().store(&inputTprFilename_).description("Run input file to dump"));
648 options->addOption(FileNameOption("f")
649 .filetype(eftTrajectory)
651 .store(&inputTrajectoryFilename_)
652 .description("Trajectory file to dump"));
654 FileNameOption("e").filetype(eftEnergy).inputFile().store(&inputEnergyFilename_).description("Energy file to dump"));
656 FileNameOption("cp").legacyType(efCPT).inputFile().store(&inputCheckpointFilename_).description("Checkpoint file to dump"));
658 FileNameOption("p").legacyType(efTOP).inputFile().store(&inputTopologyFilename_).description("Topology file to dump"));
660 FileNameOption("mtx").legacyType(efMTX).inputFile().store(&inputMatrixFilename_).description("Hessian matrix to dump"));
661 options->addOption(FileNameOption("om")
664 .store(&outputMdpFilename_)
665 .description("grompp input file from run input file"));
668 BooleanOption("nr").store(&bShowNumbers_).defaultValue(true).description("Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)"));
670 BooleanOption("param").store(&bShowParams_).defaultValue(false).description("Show parameters for each bonded interaction (for comparing dumps, it is useful to combine this with -nonr)"));
672 BooleanOption("sys").store(&bShowParams_).defaultValue(false).description("List the atoms and bonded interactions for the whole system instead of for each molecule type"));
674 BooleanOption("orgir").store(&bShowParams_).defaultValue(false).description("Show input parameters from tpr as they were written by the version that produced the file, instead of how the current version reads them"));
677 void Dump::optionsFinished()
679 // TODO Currently gmx dump ignores user input that seeks to dump
680 // multiple files. Here, we could enforce that the user only asks
686 if (!inputTprFilename_.empty())
688 list_tpr(inputTprFilename_.c_str(), bShowNumbers_, bShowParams_,
689 outputMdpFilename_.empty() ? nullptr : outputMdpFilename_.c_str(), bSysTop_,
692 else if (!inputTrajectoryFilename_.empty())
694 list_trx(inputTrajectoryFilename_.c_str());
696 else if (!inputEnergyFilename_.empty())
698 list_ene(inputEnergyFilename_.c_str());
700 else if (!inputCheckpointFilename_.empty())
702 list_checkpoint(inputCheckpointFilename_.c_str(), stdout);
704 else if (!inputTopologyFilename_.empty())
706 list_top(inputTopologyFilename_.c_str());
708 else if (!inputMatrixFilename_.empty())
710 list_mtx(inputMatrixFilename_.c_str());
718 const char DumpInfo::name[] = "dump";
719 const char DumpInfo::shortDescription[] = "Make binary files human readable";
720 ICommandLineOptionsModulePointer DumpInfo::create()
722 return std::make_unique<Dump>();