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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief Implements gmx dump utility.
41 * \ingroup module_tools
54 #include "gromacs/commandline/cmdlineoptionsmodule.h"
55 #include "gromacs/fileio/checkpoint.h"
56 #include "gromacs/fileio/enxio.h"
57 #include "gromacs/fileio/filetypes.h"
58 #include "gromacs/fileio/gmxfio.h"
59 #include "gromacs/fileio/mtxio.h"
60 #include "gromacs/fileio/tngio.h"
61 #include "gromacs/fileio/tpxio.h"
62 #include "gromacs/fileio/trrio.h"
63 #include "gromacs/fileio/xtcio.h"
64 #include "gromacs/gmxpreprocess/gmxcpp.h"
65 #include "gromacs/math/vecdump.h"
66 #include "gromacs/mdrun/mdmodules.h"
67 #include "gromacs/mdtypes/forcerec.h"
68 #include "gromacs/mdtypes/inputrec.h"
69 #include "gromacs/mdtypes/md_enums.h"
70 #include "gromacs/mdtypes/state.h"
71 #include "gromacs/options/basicoptions.h"
72 #include "gromacs/options/filenameoption.h"
73 #include "gromacs/options/ioptionscontainer.h"
74 #include "gromacs/topology/mtop_util.h"
75 #include "gromacs/topology/topology.h"
76 #include "gromacs/trajectory/energyframe.h"
77 #include "gromacs/trajectory/trajectoryframe.h"
78 #include "gromacs/utility/arraysize.h"
79 #include "gromacs/utility/basedefinitions.h"
80 #include "gromacs/utility/fatalerror.h"
81 #include "gromacs/utility/futil.h"
82 #include "gromacs/utility/smalloc.h"
83 #include "gromacs/utility/txtdump.h"
92 void list_tpr(const char* fn
,
93 gmx_bool bShowNumbers
,
94 gmx_bool bShowParameters
,
97 gmx_bool bOriginalInputrec
)
105 TpxFileHeader tpx
= readTpxHeader(fn
, true);
108 read_tpx_state(fn
, tpx
.bIr
? &ir
: nullptr, &state
, tpx
.bTop
? &mtop
: nullptr);
109 if (tpx
.bIr
&& !bOriginalInputrec
)
111 MDModules().adjustInputrecBasedOnModules(&ir
);
114 if (mdpfn
&& tpx
.bIr
)
116 gp
= gmx_fio_fopen(mdpfn
, "w");
117 pr_inputrec(gp
, 0, nullptr, &ir
, TRUE
);
125 top
= gmx_mtop_t_to_t_topology(&mtop
, false);
128 if (available(stdout
, &tpx
, 0, fn
))
131 pr_title(stdout
, indent
, fn
);
132 pr_inputrec(stdout
, 0, "inputrec", tpx
.bIr
? &ir
: nullptr, FALSE
);
134 pr_tpxheader(stdout
, indent
, "header", &(tpx
));
138 pr_mtop(stdout
, indent
, "topology", &(mtop
), bShowNumbers
, bShowParameters
);
142 pr_top(stdout
, indent
, "topology", &(top
), bShowNumbers
, bShowParameters
);
145 pr_rvecs(stdout
, indent
, "box", tpx
.bBox
? state
.box
: nullptr, DIM
);
146 pr_rvecs(stdout
, indent
, "box_rel", tpx
.bBox
? state
.box_rel
: nullptr, DIM
);
147 pr_rvecs(stdout
, indent
, "boxv", tpx
.bBox
? state
.boxv
: nullptr, DIM
);
148 pr_rvecs(stdout
, indent
, "pres_prev", tpx
.bBox
? state
.pres_prev
: nullptr, DIM
);
149 pr_rvecs(stdout
, indent
, "svir_prev", tpx
.bBox
? state
.svir_prev
: nullptr, DIM
);
150 pr_rvecs(stdout
, indent
, "fvir_prev", tpx
.bBox
? state
.fvir_prev
: nullptr, DIM
);
151 /* leave nosehoover_xi in for now to match the tpr version */
152 pr_doubles(stdout
, indent
, "nosehoover_xi", state
.nosehoover_xi
.data(), state
.ngtc
);
153 /*pr_doubles(stdout,indent,"nosehoover_vxi",state.nosehoover_vxi,state.ngtc);*/
154 /*pr_doubles(stdout,indent,"therm_integral",state.therm_integral,state.ngtc);*/
155 pr_rvecs(stdout
, indent
, "x", tpx
.bX
? state
.x
.rvec_array() : nullptr, state
.natoms
);
156 pr_rvecs(stdout
, indent
, "v", tpx
.bV
? state
.v
.rvec_array() : nullptr, state
.natoms
);
159 const SimulationGroups
& groups
= mtop
.groups
;
161 gmx::EnumerationArray
<SimulationAtomGroupType
, std::vector
<int>> gcount
;
162 for (auto group
: keysOf(gcount
))
164 gcount
[group
].resize(groups
.groups
[group
].size());
167 for (int i
= 0; (i
< mtop
.natoms
); i
++)
169 for (auto group
: keysOf(gcount
))
171 gcount
[group
][getGroupType(groups
, group
, i
)]++;
174 printf("Group statistics\n");
175 for (auto group
: keysOf(gcount
))
178 printf("%-12s: ", shortName(group
));
179 for (const auto& entry
: gcount
[group
])
181 printf(" %5d", entry
);
184 printf(" (total %d atoms)\n", atot
);
189 //! Dump a topology file
190 void list_top(const char* fn
)
193 // Legacy string length macro
196 char* cppopts
[] = { nullptr };
198 status
= cpp_open_file(fn
, &handle
, cppopts
);
201 gmx_fatal(FARGS
, "%s", cpp_error(&handle
, status
));
205 status
= cpp_read_line(&handle
, STRLEN
, buf
);
206 done
= static_cast<int>(status
== eCPP_EOF
);
209 if (status
!= eCPP_OK
)
211 gmx_fatal(FARGS
, "%s", cpp_error(&handle
, status
));
219 status
= cpp_close_file(&handle
);
220 if (status
!= eCPP_OK
)
222 gmx_fatal(FARGS
, "%s", cpp_error(&handle
, status
));
227 void list_trr(const char* fn
)
234 gmx_trr_header_t trrheader
;
237 fpread
= gmx_trr_open(fn
, "r");
240 while (gmx_trr_read_frame_header(fpread
, &trrheader
, &bOK
))
242 snew(x
, trrheader
.natoms
);
243 snew(v
, trrheader
.natoms
);
244 snew(f
, trrheader
.natoms
);
245 if (gmx_trr_read_frame_data(fpread
, &trrheader
, trrheader
.box_size
? box
: nullptr,
246 trrheader
.x_size
? x
: nullptr, trrheader
.v_size
? v
: nullptr,
247 trrheader
.f_size
? f
: nullptr))
249 sprintf(buf
, "%s frame %d", fn
, nframe
);
251 indent
= pr_title(stdout
, indent
, buf
);
252 pr_indent(stdout
, indent
);
253 fprintf(stdout
, "natoms=%10d step=%10" PRId64
" time=%12.7e lambda=%10g\n",
254 trrheader
.natoms
, trrheader
.step
, trrheader
.t
, trrheader
.lambda
);
255 if (trrheader
.box_size
)
257 pr_rvecs(stdout
, indent
, "box", box
, DIM
);
259 if (trrheader
.x_size
)
261 pr_rvecs(stdout
, indent
, "x", x
, trrheader
.natoms
);
263 if (trrheader
.v_size
)
265 pr_rvecs(stdout
, indent
, "v", v
, trrheader
.natoms
);
267 if (trrheader
.f_size
)
269 pr_rvecs(stdout
, indent
, "f", f
, trrheader
.natoms
);
274 fprintf(stderr
, "\nWARNING: Incomplete frame: nr %d, t=%g\n", nframe
, trrheader
.t
);
284 fprintf(stderr
, "\nWARNING: Incomplete frame header: nr %d, t=%g\n", nframe
, trrheader
.t
);
286 gmx_trr_close(fpread
);
290 void list_xtc(const char* fn
)
302 xd
= open_xtc(fn
, "r");
303 read_first_xtc(xd
, &natoms
, &step
, &time
, box
, &x
, &prec
, &bOK
);
308 sprintf(buf
, "%s frame %d", fn
, nframe
);
310 indent
= pr_title(stdout
, indent
, buf
);
311 pr_indent(stdout
, indent
);
312 fprintf(stdout
, "natoms=%10d step=%10" PRId64
" time=%12.7e prec=%10g\n", natoms
, step
,
314 pr_rvecs(stdout
, indent
, "box", box
, DIM
);
315 pr_rvecs(stdout
, indent
, "x", x
, natoms
);
317 } while (read_next_xtc(xd
, natoms
, &step
, &time
, box
, x
, &prec
, &bOK
) != 0);
320 fprintf(stderr
, "\nWARNING: Incomplete frame at time %g\n", time
);
328 /*! \brief Callback used by list_tng_for_gmx_dump. */
329 void list_tng_inner(const char* fn
,
330 gmx_bool bFirstFrame
,
334 int64_t n_values_per_frame
,
345 sprintf(buf
, "%s frame %" PRId64
, fn
, nframe
);
347 indent
= pr_title(stdout
, indent
, buf
);
348 pr_indent(stdout
, indent
);
349 fprintf(stdout
, "natoms=%10" PRId64
" step=%10" PRId64
" time=%12.7e", n_atoms
, step
, frame_time
);
352 fprintf(stdout
, " prec=%10g", prec
);
354 fprintf(stdout
, "\n");
356 pr_reals_of_dim(stdout
, indent
, block_name
, values
, n_atoms
, n_values_per_frame
);
362 void list_tng(const char* fn
)
365 gmx_tng_trajectory_t tng
;
367 int64_t i
, *block_ids
= nullptr, step
, ndatablocks
;
369 real
* values
= nullptr;
371 gmx_tng_open(fn
, 'r', &tng
);
372 gmx_print_tng_molecule_system(tng
, stdout
);
374 bOK
= gmx_get_tng_data_block_types_of_next_frame(tng
, -1, 0, nullptr, &step
, &ndatablocks
, &block_ids
);
377 for (i
= 0; i
< ndatablocks
; i
++)
381 int64_t n_values_per_frame
, n_atoms
;
382 char block_name
[STRLEN
];
384 gmx_get_tng_data_next_frame_of_block_type(tng
, block_ids
[i
], &values
, &step
,
385 &frame_time
, &n_values_per_frame
, &n_atoms
,
386 &prec
, block_name
, STRLEN
, &bOK
);
389 /* Can't write any output because we don't know what
391 fprintf(stderr
, "\nWARNING: Incomplete frame at time %g, will not write output\n",
396 list_tng_inner(fn
, (0 == i
), values
, step
, frame_time
, n_values_per_frame
, n_atoms
,
397 prec
, nframe
, block_name
);
401 } while (gmx_get_tng_data_block_types_of_next_frame(tng
, step
, 0, nullptr, &step
, &ndatablocks
,
411 GMX_UNUSED_VALUE(fn
);
415 //! Dump a trajectory file
416 void list_trx(const char* fn
)
420 case efXTC
: list_xtc(fn
); break;
421 case efTRR
: list_trr(fn
); break;
422 case efTNG
: list_tng(fn
); break;
424 fprintf(stderr
, "File %s is of an unsupported type. Try using the command\n 'less %s'\n",
429 //! Dump an energy file
430 void list_ene(const char* fn
)
434 gmx_enxnm_t
* enm
= nullptr;
439 printf("gmx dump: %s\n", fn
);
440 in
= open_enx(fn
, "r");
441 do_enxnms(in
, &nre
, &enm
);
444 printf("energy components:\n");
445 for (i
= 0; (i
< nre
); i
++)
447 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:", fr
->t
, "step:", gmx_step_str(fr
->step
, buf
));
458 printf("%24s %12s %12s %12s\n", "", "", "nsteps:", gmx_step_str(fr
->nsteps
, buf
));
459 printf("%24s %12.5e %12s %12s\n", "delta_t:", fr
->dt
,
460 "sum steps:", gmx_step_str(fr
->nsum
, buf
));
463 printf("%24s %12s %12s %12s\n", "Component", "Energy", "Av. Energy",
467 for (i
= 0; (i
< nre
); i
++)
469 printf("%24s %12.5e %12.5e %12.5e\n", enm
[i
].name
, fr
->ener
[i
].e
,
470 fr
->ener
[i
].eav
, fr
->ener
[i
].esum
);
475 for (i
= 0; (i
< nre
); i
++)
477 printf("%24s %12.5e\n", enm
[i
].name
, fr
->ener
[i
].e
);
481 for (b
= 0; b
< fr
->nblock
; b
++)
483 const char* typestr
= "";
485 t_enxblock
* eb
= &(fr
->block
[b
]);
486 printf("Block data %2d (%3d subblocks, id=%d)\n", b
, eb
->nsub
, eb
->id
);
490 typestr
= enx_block_id_name
[eb
->id
];
492 printf(" id='%s'\n", typestr
);
493 for (i
= 0; i
< eb
->nsub
; i
++)
495 t_enxsubblock
* sb
= &(eb
->sub
[i
]);
496 printf(" Sub block %3d (%5d elems, type=%s) values:\n", i
, sb
->nr
,
497 xdr_datatype_names
[sb
->type
]);
501 case xdr_datatype_float
:
502 for (j
= 0; j
< sb
->nr
; j
++)
504 printf("%14d %8.4f\n", j
, sb
->fval
[j
]);
507 case xdr_datatype_double
:
508 for (j
= 0; j
< sb
->nr
; j
++)
510 printf("%14d %10.6f\n", j
, sb
->dval
[j
]);
513 case xdr_datatype_int
:
514 for (j
= 0; j
< sb
->nr
; j
++)
516 printf("%14d %10d\n", j
, sb
->ival
[j
]);
519 case xdr_datatype_int64
:
520 for (j
= 0; j
< sb
->nr
; j
++)
522 printf("%14d %s\n", j
, gmx_step_str(sb
->lval
[j
], buf
));
525 case xdr_datatype_char
:
526 for (j
= 0; j
< sb
->nr
; j
++)
528 printf("%14d %1c\n", j
, sb
->cval
[j
]);
531 case xdr_datatype_string
:
532 for (j
= 0; j
< sb
->nr
; j
++)
534 printf("%14d %80s\n", j
, sb
->sval
[j
]);
537 default: gmx_incons("Unknown subblock type");
551 //! Dump a (Hessian) matrix file
552 void list_mtx(const char* fn
)
554 int nrow
, ncol
, i
, j
, k
;
555 real
* full
= nullptr, value
;
556 gmx_sparsematrix_t
* sparse
= nullptr;
558 gmx_mtxio_read(fn
, &nrow
, &ncol
, &full
, &sparse
);
562 snew(full
, nrow
* ncol
);
563 for (i
= 0; i
< nrow
* ncol
; i
++)
568 for (i
= 0; i
< sparse
->nrow
; i
++)
570 for (j
= 0; j
< sparse
->ndata
[i
]; j
++)
572 k
= sparse
->data
[i
][j
].col
;
573 value
= sparse
->data
[i
][j
].value
;
574 full
[i
* ncol
+ k
] = value
;
575 full
[k
* ncol
+ i
] = value
;
578 gmx_sparsematrix_destroy(sparse
);
581 printf("%d %d\n", nrow
, ncol
);
582 for (i
= 0; i
< nrow
; i
++)
584 for (j
= 0; j
< ncol
; j
++)
586 printf(" %g", full
[i
* ncol
+ j
]);
594 class Dump
: public ICommandLineOptionsModule
599 // From ICommandLineOptionsModule
600 void init(CommandLineModuleSettings
* /*settings*/) override
{}
602 void initOptions(IOptionsContainer
* options
, ICommandLineOptionsModuleSettings
* settings
) override
;
604 void optionsFinished() override
;
609 //! Commandline options
611 bool bShowNumbers_
= true;
612 bool bShowParams_
= false;
613 bool bSysTop_
= false;
614 bool bOriginalInputrec_
= false;
616 //! Commandline file options
618 std::string inputTprFilename_
;
619 std::string inputTrajectoryFilename_
;
620 std::string inputEnergyFilename_
;
621 std::string inputCheckpointFilename_
;
622 std::string inputTopologyFilename_
;
623 std::string inputMatrixFilename_
;
624 std::string outputMdpFilename_
;
628 void Dump::initOptions(IOptionsContainer
* options
, ICommandLineOptionsModuleSettings
* settings
)
630 const char* desc
[] = { "[THISMODULE] reads a run input file ([REF].tpr[ref]),",
631 "a trajectory ([REF].trr[ref]/[REF].xtc[ref]/[TT]tng[tt]), an energy",
632 "file ([REF].edr[ref]), a checkpoint file ([REF].cpt[ref])",
633 "or topology file ([REF].top[ref])",
634 "and prints that to standard output in a readable format.",
635 "This program is essential for checking your run input file in case of",
637 settings
->setHelpText(desc
);
639 const char* bugs
[] = {
640 "The [REF].mdp[ref] file produced by [TT]-om[tt] can not be read by grompp."
642 settings
->setBugText(bugs
);
643 // TODO If this ancient note acknowledging a bug is still true,
644 // fix it or block that run path:
645 // Position restraint output from -sys -s is broken
648 FileNameOption("s").filetype(eftRunInput
).inputFile().store(&inputTprFilename_
).description("Run input file to dump"));
649 options
->addOption(FileNameOption("f")
650 .filetype(eftTrajectory
)
652 .store(&inputTrajectoryFilename_
)
653 .description("Trajectory file to dump"));
655 FileNameOption("e").filetype(eftEnergy
).inputFile().store(&inputEnergyFilename_
).description("Energy file to dump"));
657 FileNameOption("cp").legacyType(efCPT
).inputFile().store(&inputCheckpointFilename_
).description("Checkpoint file to dump"));
659 FileNameOption("p").legacyType(efTOP
).inputFile().store(&inputTopologyFilename_
).description("Topology file to dump"));
661 FileNameOption("mtx").legacyType(efMTX
).inputFile().store(&inputMatrixFilename_
).description("Hessian matrix to dump"));
662 options
->addOption(FileNameOption("om")
665 .store(&outputMdpFilename_
)
666 .description("grompp input file from run input file"));
669 BooleanOption("nr").store(&bShowNumbers_
).defaultValue(true).description("Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)"));
671 BooleanOption("param").store(&bShowParams_
).defaultValue(false).description("Show parameters for each bonded interaction (for comparing dumps, it is useful to combine this with -nonr)"));
673 BooleanOption("sys").store(&bShowParams_
).defaultValue(false).description("List the atoms and bonded interactions for the whole system instead of for each molecule type"));
675 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"));
678 void Dump::optionsFinished()
680 // TODO Currently gmx dump ignores user input that seeks to dump
681 // multiple files. Here, we could enforce that the user only asks
687 if (!inputTprFilename_
.empty())
689 list_tpr(inputTprFilename_
.c_str(), bShowNumbers_
, bShowParams_
,
690 outputMdpFilename_
.empty() ? nullptr : outputMdpFilename_
.c_str(), bSysTop_
,
693 else if (!inputTrajectoryFilename_
.empty())
695 list_trx(inputTrajectoryFilename_
.c_str());
697 else if (!inputEnergyFilename_
.empty())
699 list_ene(inputEnergyFilename_
.c_str());
701 else if (!inputCheckpointFilename_
.empty())
703 list_checkpoint(inputCheckpointFilename_
.c_str(), stdout
);
705 else if (!inputTopologyFilename_
.empty())
707 list_top(inputTopologyFilename_
.c_str());
709 else if (!inputMatrixFilename_
.empty())
711 list_mtx(inputMatrixFilename_
.c_str());
719 const char DumpInfo::name
[] = "dump";
720 const char DumpInfo::shortDescription
[] = "Make binary files human readable";
721 ICommandLineOptionsModulePointer
DumpInfo::create()
723 return std::make_unique
<Dump
>();