541d6da1bfdc2107cbfdbefee978d3fab0d15fe5
[alexxy/gromacs.git] / src / gromacs / tools / dump.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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,2021, 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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 /*! \internal \file
39  * \brief Implements gmx dump utility.
40  *
41  * \ingroup module_tools
42  */
43 #include "gmxpre.h"
44
45 #include "dump.h"
46
47 #include "config.h"
48
49 #include <cassert>
50 #include <cmath>
51 #include <cstdio>
52 #include <cstring>
53
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"
84
85 namespace gmx
86 {
87
88 namespace
89 {
90
91 //! Dump a TPR file
92 void list_tpr(const char* fn,
93               gmx_bool    bShowNumbers,
94               gmx_bool    bShowParameters,
95               const char* mdpfn,
96               gmx_bool    bSysTop,
97               gmx_bool    bOriginalInputrec)
98 {
99     FILE*      gp;
100     int        indent, atot;
101     t_state    state;
102     gmx_mtop_t mtop;
103     t_topology top;
104
105     TpxFileHeader tpx = readTpxHeader(fn, true);
106     t_inputrec    ir;
107
108     read_tpx_state(fn, tpx.bIr ? &ir : nullptr, &state, tpx.bTop ? &mtop : nullptr);
109     if (tpx.bIr && !bOriginalInputrec)
110     {
111         MDModules().adjustInputrecBasedOnModules(&ir);
112     }
113
114     if (mdpfn && tpx.bIr)
115     {
116         gp = gmx_fio_fopen(mdpfn, "w");
117         pr_inputrec(gp, 0, nullptr, &ir, TRUE);
118         gmx_fio_fclose(gp);
119     }
120
121     if (!mdpfn)
122     {
123         if (bSysTop)
124         {
125             top = gmx_mtop_t_to_t_topology(&mtop, false);
126         }
127
128         if (available(stdout, &tpx, 0, fn))
129         {
130             indent = 0;
131             pr_title(stdout, indent, fn);
132             pr_inputrec(stdout, 0, "inputrec", tpx.bIr ? &ir : nullptr, FALSE);
133
134             pr_tpxheader(stdout, indent, "header", &(tpx));
135
136             if (!bSysTop)
137             {
138                 pr_mtop(stdout, indent, "topology", &(mtop), bShowNumbers, bShowParameters);
139             }
140             else
141             {
142                 pr_top(stdout, indent, "topology", &(top), bShowNumbers, bShowParameters);
143             }
144
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);
157         }
158
159         const SimulationGroups& groups = mtop.groups;
160
161         gmx::EnumerationArray<SimulationAtomGroupType, std::vector<int>> gcount;
162         for (auto group : keysOf(gcount))
163         {
164             gcount[group].resize(groups.groups[group].size());
165         }
166
167         for (int i = 0; (i < mtop.natoms); i++)
168         {
169             for (auto group : keysOf(gcount))
170             {
171                 if (group == SimulationAtomGroupType::AccelerationUnused)
172                 {
173                     continue;
174                 }
175                 gcount[group][getGroupType(groups, group, i)]++;
176             }
177         }
178         printf("Group statistics\n");
179         for (auto group : keysOf(gcount))
180         {
181             atot = 0;
182             printf("%-12s: ", shortName(group));
183             for (const auto& entry : gcount[group])
184             {
185                 printf("  %5d", entry);
186                 atot += entry;
187             }
188             printf("  (total %d atoms)\n", atot);
189         }
190     }
191 }
192
193 //! Dump a topology file
194 void list_top(const char* fn)
195 {
196     int status, done;
197     // Legacy string length macro
198     char      buf[STRLEN];
199     gmx_cpp_t handle;
200     char*     cppopts[] = { nullptr };
201
202     status = cpp_open_file(fn, &handle, cppopts);
203     if (status != 0)
204     {
205         gmx_fatal(FARGS, "%s", cpp_error(&handle, status));
206     }
207     do
208     {
209         status = cpp_read_line(&handle, STRLEN, buf);
210         done   = static_cast<int>(status == eCPP_EOF);
211         if (!done)
212         {
213             if (status != eCPP_OK)
214             {
215                 gmx_fatal(FARGS, "%s", cpp_error(&handle, status));
216             }
217             else
218             {
219                 printf("%s\n", buf);
220             }
221         }
222     } while (done == 0);
223     status = cpp_close_file(&handle);
224     if (status != eCPP_OK)
225     {
226         gmx_fatal(FARGS, "%s", cpp_error(&handle, status));
227     }
228 }
229
230 //! Dump a TRR file
231 void list_trr(const char* fn)
232 {
233     t_fileio*        fpread;
234     int              nframe, indent;
235     char             buf[256];
236     rvec *           x, *v, *f;
237     matrix           box;
238     gmx_trr_header_t trrheader;
239     gmx_bool         bOK;
240
241     fpread = gmx_trr_open(fn, "r");
242
243     nframe = 0;
244     while (gmx_trr_read_frame_header(fpread, &trrheader, &bOK))
245     {
246         snew(x, trrheader.natoms);
247         snew(v, trrheader.natoms);
248         snew(f, trrheader.natoms);
249         if (gmx_trr_read_frame_data(fpread,
250                                     &trrheader,
251                                     trrheader.box_size ? box : nullptr,
252                                     trrheader.x_size ? x : nullptr,
253                                     trrheader.v_size ? v : nullptr,
254                                     trrheader.f_size ? f : nullptr))
255         {
256             sprintf(buf, "%s frame %d", fn, nframe);
257             indent = 0;
258             indent = pr_title(stdout, indent, buf);
259             pr_indent(stdout, indent);
260             fprintf(stdout,
261                     "natoms=%10d  step=%10" PRId64 "  time=%12.7e  lambda=%10g\n",
262                     trrheader.natoms,
263                     trrheader.step,
264                     trrheader.t,
265                     trrheader.lambda);
266             if (trrheader.box_size)
267             {
268                 pr_rvecs(stdout, indent, "box", box, DIM);
269             }
270             if (trrheader.x_size)
271             {
272                 pr_rvecs(stdout, indent, "x", x, trrheader.natoms);
273             }
274             if (trrheader.v_size)
275             {
276                 pr_rvecs(stdout, indent, "v", v, trrheader.natoms);
277             }
278             if (trrheader.f_size)
279             {
280                 pr_rvecs(stdout, indent, "f", f, trrheader.natoms);
281             }
282         }
283         else
284         {
285             fprintf(stderr, "\nWARNING: Incomplete frame: nr %d, t=%g\n", nframe, trrheader.t);
286         }
287
288         sfree(x);
289         sfree(v);
290         sfree(f);
291         nframe++;
292     }
293     if (!bOK)
294     {
295         fprintf(stderr, "\nWARNING: Incomplete frame header: nr %d, t=%g\n", nframe, trrheader.t);
296     }
297     gmx_trr_close(fpread);
298 }
299
300 //! Dump an xtc file
301 void list_xtc(const char* fn)
302 {
303     t_fileio* xd;
304     int       indent;
305     char      buf[256];
306     rvec*     x;
307     matrix    box;
308     int       nframe, natoms;
309     int64_t   step;
310     real      prec, time;
311     gmx_bool  bOK;
312
313     xd = open_xtc(fn, "r");
314     read_first_xtc(xd, &natoms, &step, &time, box, &x, &prec, &bOK);
315
316     nframe = 0;
317     do
318     {
319         sprintf(buf, "%s frame %d", fn, nframe);
320         indent = 0;
321         indent = pr_title(stdout, indent, buf);
322         pr_indent(stdout, indent);
323         fprintf(stdout, "natoms=%10d  step=%10" PRId64 "  time=%12.7e  prec=%10g\n", natoms, step, time, prec);
324         pr_rvecs(stdout, indent, "box", box, DIM);
325         pr_rvecs(stdout, indent, "x", x, natoms);
326         nframe++;
327     } while (read_next_xtc(xd, natoms, &step, &time, box, x, &prec, &bOK) != 0);
328     if (!bOK)
329     {
330         fprintf(stderr, "\nWARNING: Incomplete frame at time %g\n", time);
331     }
332     sfree(x);
333     close_xtc(xd);
334 }
335
336 #if GMX_USE_TNG
337
338 /*! \brief Callback used by list_tng_for_gmx_dump. */
339 void list_tng_inner(const char* fn,
340                     gmx_bool    bFirstFrame,
341                     real*       values,
342                     int64_t     step,
343                     double      frame_time,
344                     int64_t     n_values_per_frame,
345                     int64_t     n_atoms,
346                     real        prec,
347                     int64_t     nframe,
348                     char*       block_name)
349 {
350     char buf[256];
351     int  indent = 0;
352
353     if (bFirstFrame)
354     {
355         sprintf(buf, "%s frame %" PRId64, fn, nframe);
356         indent = 0;
357         indent = pr_title(stdout, indent, buf);
358         pr_indent(stdout, indent);
359         fprintf(stdout, "natoms=%10" PRId64 "  step=%10" PRId64 "  time=%12.7e", n_atoms, step, frame_time);
360         if (prec > 0)
361         {
362             fprintf(stdout, "  prec=%10g", prec);
363         }
364         fprintf(stdout, "\n");
365     }
366     pr_reals_of_dim(stdout, indent, block_name, values, n_atoms, n_values_per_frame);
367 }
368
369 #endif
370
371 //! Dump a TNG file
372 void list_tng(const char* fn)
373 {
374 #if GMX_USE_TNG
375     gmx_tng_trajectory_t tng;
376     int64_t              nframe = 0;
377     int64_t              i, *block_ids = nullptr, step, ndatablocks;
378     gmx_bool             bOK;
379     real*                values = nullptr;
380
381     gmx_tng_open(fn, 'r', &tng);
382     gmx_print_tng_molecule_system(tng, stdout);
383
384     bOK = gmx_get_tng_data_block_types_of_next_frame(tng, -1, 0, nullptr, &step, &ndatablocks, &block_ids);
385     do
386     {
387         for (i = 0; i < ndatablocks; i++)
388         {
389             double  frame_time;
390             real    prec;
391             int64_t n_values_per_frame, n_atoms;
392             char    block_name[STRLEN];
393
394             gmx_get_tng_data_next_frame_of_block_type(tng,
395                                                       block_ids[i],
396                                                       &values,
397                                                       &step,
398                                                       &frame_time,
399                                                       &n_values_per_frame,
400                                                       &n_atoms,
401                                                       &prec,
402                                                       block_name,
403                                                       STRLEN,
404                                                       &bOK);
405             if (!bOK)
406             {
407                 /* Can't write any output because we don't know what
408                    arrays are valid. */
409                 fprintf(stderr, "\nWARNING: Incomplete frame at time %g, will not write output\n", frame_time);
410             }
411             else
412             {
413                 list_tng_inner(
414                         fn, (0 == i), values, step, frame_time, n_values_per_frame, n_atoms, prec, nframe, block_name);
415             }
416         }
417         nframe++;
418     } while (gmx_get_tng_data_block_types_of_next_frame(
419             tng, step, 0, nullptr, &step, &ndatablocks, &block_ids));
420
421     if (block_ids)
422     {
423         sfree(block_ids);
424     }
425     sfree(values);
426     gmx_tng_close(&tng);
427 #else
428     GMX_UNUSED_VALUE(fn);
429 #endif
430 }
431
432 //! Dump a trajectory file
433 void list_trx(const char* fn)
434 {
435     switch (fn2ftp(fn))
436     {
437         case efXTC: list_xtc(fn); break;
438         case efTRR: list_trr(fn); break;
439         case efTNG: list_tng(fn); break;
440         default:
441             fprintf(stderr, "File %s is of an unsupported type. Try using the command\n 'less %s'\n", fn, fn);
442     }
443 }
444
445 //! Dump an energy file
446 void list_ene(const char* fn)
447 {
448     ener_file_t  in;
449     gmx_bool     bCont;
450     gmx_enxnm_t* enm = nullptr;
451     t_enxframe*  fr;
452     int          i, j, nre, b;
453     char         buf[22];
454
455     printf("gmx dump: %s\n", fn);
456     in = open_enx(fn, "r");
457     do_enxnms(in, &nre, &enm);
458     assert(enm);
459
460     printf("energy components:\n");
461     for (i = 0; (i < nre); i++)
462     {
463         printf("%5d  %-24s (%s)\n", i, enm[i].name, enm[i].unit);
464     }
465
466     snew(fr, 1);
467     do
468     {
469         bCont = do_enx(in, fr);
470
471         if (bCont)
472         {
473             printf("\n%24s  %12.5e  %12s  %12s\n", "time:", fr->t, "step:", gmx_step_str(fr->step, buf));
474             printf("%24s  %12s  %12s  %12s\n", "", "", "nsteps:", gmx_step_str(fr->nsteps, buf));
475             printf("%24s  %12.5e  %12s  %12s\n", "delta_t:", fr->dt, "sum steps:", gmx_step_str(fr->nsum, buf));
476             if (fr->nre == nre)
477             {
478                 printf("%24s  %12s  %12s  %12s\n",
479                        "Component",
480                        "Energy",
481                        "Av. Energy",
482                        "Sum Energy");
483                 if (fr->nsum > 0)
484                 {
485                     for (i = 0; (i < nre); i++)
486                     {
487                         printf("%24s  %12.5e  %12.5e  %12.5e\n",
488                                enm[i].name,
489                                fr->ener[i].e,
490                                fr->ener[i].eav,
491                                fr->ener[i].esum);
492                     }
493                 }
494                 else
495                 {
496                     for (i = 0; (i < nre); i++)
497                     {
498                         printf("%24s  %12.5e\n", enm[i].name, fr->ener[i].e);
499                     }
500                 }
501             }
502             for (b = 0; b < fr->nblock; b++)
503             {
504                 const char* typestr = "";
505
506                 t_enxblock* eb = &(fr->block[b]);
507                 printf("Block data %2d (%3d subblocks, id=%d)\n", b, eb->nsub, eb->id);
508
509                 if (eb->id < enxNR)
510                 {
511                     typestr = enx_block_id_name[eb->id];
512                 }
513                 printf("  id='%s'\n", typestr);
514                 for (i = 0; i < eb->nsub; i++)
515                 {
516                     t_enxsubblock* sb = &(eb->sub[i]);
517                     printf("  Sub block %3d (%5d elems, type=%s) values:\n",
518                            i,
519                            sb->nr,
520                            enumValueToString(sb->type));
521
522                     switch (sb->type)
523                     {
524                         case XdrDataType::Float:
525                             for (j = 0; j < sb->nr; j++)
526                             {
527                                 printf("%14d   %8.4f\n", j, sb->fval[j]);
528                             }
529                             break;
530                         case XdrDataType::Double:
531                             for (j = 0; j < sb->nr; j++)
532                             {
533                                 printf("%14d   %10.6f\n", j, sb->dval[j]);
534                             }
535                             break;
536                         case XdrDataType::Int:
537                             for (j = 0; j < sb->nr; j++)
538                             {
539                                 printf("%14d %10d\n", j, sb->ival[j]);
540                             }
541                             break;
542                         case XdrDataType::Int64:
543                             for (j = 0; j < sb->nr; j++)
544                             {
545                                 printf("%14d %s\n", j, gmx_step_str(sb->lval[j], buf));
546                             }
547                             break;
548                         case XdrDataType::Char:
549                             for (j = 0; j < sb->nr; j++)
550                             {
551                                 printf("%14d %1c\n", j, sb->cval[j]);
552                             }
553                             break;
554                         case XdrDataType::String:
555                             for (j = 0; j < sb->nr; j++)
556                             {
557                                 printf("%14d %80s\n", j, sb->sval[j]);
558                             }
559                             break;
560                         default: gmx_incons("Unknown subblock type");
561                     }
562                 }
563             }
564         }
565     } while (bCont);
566
567     close_enx(in);
568
569     free_enxframe(fr);
570     sfree(fr);
571     sfree(enm);
572 }
573
574 //! Dump a (Hessian) matrix file
575 void list_mtx(const char* fn)
576 {
577     int                 nrow, ncol, i, j, k;
578     real *              full   = nullptr, value;
579     gmx_sparsematrix_t* sparse = nullptr;
580
581     gmx_mtxio_read(fn, &nrow, &ncol, &full, &sparse);
582
583     if (full == nullptr)
584     {
585         snew(full, nrow * ncol);
586         for (i = 0; i < nrow * ncol; i++)
587         {
588             full[i] = 0;
589         }
590
591         for (i = 0; i < sparse->nrow; i++)
592         {
593             for (j = 0; j < sparse->ndata[i]; j++)
594             {
595                 k                  = sparse->data[i][j].col;
596                 value              = sparse->data[i][j].value;
597                 full[i * ncol + k] = value;
598                 full[k * ncol + i] = value;
599             }
600         }
601         gmx_sparsematrix_destroy(sparse);
602     }
603
604     printf("%d %d\n", nrow, ncol);
605     for (i = 0; i < nrow; i++)
606     {
607         for (j = 0; j < ncol; j++)
608         {
609             printf(" %g", full[i * ncol + j]);
610         }
611         printf("\n");
612     }
613
614     sfree(full);
615 }
616
617 class Dump : public ICommandLineOptionsModule
618 {
619 public:
620     Dump() {}
621
622     // From ICommandLineOptionsModule
623     void init(CommandLineModuleSettings* /*settings*/) override {}
624
625     void initOptions(IOptionsContainer* options, ICommandLineOptionsModuleSettings* settings) override;
626
627     void optionsFinished() override;
628
629     int run() override;
630
631 private:
632     //! Commandline options
633     //! \{
634     bool bShowNumbers_      = true;
635     bool bShowParams_       = false;
636     bool bSysTop_           = false;
637     bool bOriginalInputrec_ = false;
638     //! \}
639     //! Commandline file options
640     //! \{
641     std::string inputTprFilename_;
642     std::string inputTrajectoryFilename_;
643     std::string inputEnergyFilename_;
644     std::string inputCheckpointFilename_;
645     std::string inputTopologyFilename_;
646     std::string inputMatrixFilename_;
647     std::string outputMdpFilename_;
648     //! \}
649 };
650
651 void Dump::initOptions(IOptionsContainer* options, ICommandLineOptionsModuleSettings* settings)
652 {
653     const char* desc[] = { "[THISMODULE] reads a run input file ([REF].tpr[ref]),",
654                            "a trajectory ([REF].trr[ref]/[REF].xtc[ref]/[TT]tng[tt]), an energy",
655                            "file ([REF].edr[ref]), a checkpoint file ([REF].cpt[ref])",
656                            "or topology file ([REF].top[ref])",
657                            "and prints that to standard output in a readable format.",
658                            "This program is essential for checking your run input file in case of",
659                            "problems." };
660     settings->setHelpText(desc);
661
662     const char* bugs[] = {
663         "The [REF].mdp[ref] file produced by [TT]-om[tt] can not be read by grompp."
664     };
665     settings->setBugText(bugs);
666     // TODO If this ancient note acknowledging a bug is still true,
667     // fix it or block that run path:
668     //   Position restraint output from -sys -s is broken
669
670     options->addOption(FileNameOption("s")
671                                .filetype(OptionFileType::RunInput)
672                                .inputFile()
673                                .store(&inputTprFilename_)
674                                .description("Run input file to dump"));
675     options->addOption(FileNameOption("f")
676                                .filetype(OptionFileType::Trajectory)
677                                .inputFile()
678                                .store(&inputTrajectoryFilename_)
679                                .description("Trajectory file to dump"));
680     options->addOption(FileNameOption("e")
681                                .filetype(OptionFileType::Energy)
682                                .inputFile()
683                                .store(&inputEnergyFilename_)
684                                .description("Energy file to dump"));
685     options->addOption(
686             FileNameOption("cp").legacyType(efCPT).inputFile().store(&inputCheckpointFilename_).description("Checkpoint file to dump"));
687     options->addOption(
688             FileNameOption("p").legacyType(efTOP).inputFile().store(&inputTopologyFilename_).description("Topology file to dump"));
689     options->addOption(
690             FileNameOption("mtx").legacyType(efMTX).inputFile().store(&inputMatrixFilename_).description("Hessian matrix to dump"));
691     options->addOption(FileNameOption("om")
692                                .legacyType(efMDP)
693                                .outputFile()
694                                .store(&outputMdpFilename_)
695                                .description("grompp input file from run input file"));
696
697     options->addOption(
698             BooleanOption("nr").store(&bShowNumbers_).defaultValue(true).description("Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)"));
699     options->addOption(
700             BooleanOption("param").store(&bShowParams_).defaultValue(false).description("Show parameters for each bonded interaction (for comparing dumps, it is useful to combine this with -nonr)"));
701     options->addOption(
702             BooleanOption("sys").store(&bShowParams_).defaultValue(false).description("List the atoms and bonded interactions for the whole system instead of for each molecule type"));
703     options->addOption(
704             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"));
705 }
706
707 void Dump::optionsFinished()
708 {
709     // TODO Currently gmx dump ignores user input that seeks to dump
710     // multiple files. Here, we could enforce that the user only asks
711     // to dump one file.
712 }
713
714 int Dump::run()
715 {
716     if (!inputTprFilename_.empty())
717     {
718         list_tpr(inputTprFilename_.c_str(),
719                  bShowNumbers_,
720                  bShowParams_,
721                  outputMdpFilename_.empty() ? nullptr : outputMdpFilename_.c_str(),
722                  bSysTop_,
723                  bOriginalInputrec_);
724     }
725     else if (!inputTrajectoryFilename_.empty())
726     {
727         list_trx(inputTrajectoryFilename_.c_str());
728     }
729     else if (!inputEnergyFilename_.empty())
730     {
731         list_ene(inputEnergyFilename_.c_str());
732     }
733     else if (!inputCheckpointFilename_.empty())
734     {
735         list_checkpoint(inputCheckpointFilename_.c_str(), stdout);
736     }
737     else if (!inputTopologyFilename_.empty())
738     {
739         list_top(inputTopologyFilename_.c_str());
740     }
741     else if (!inputMatrixFilename_.empty())
742     {
743         list_mtx(inputMatrixFilename_.c_str());
744     }
745
746     return 0;
747 }
748
749 } // namespace
750
751 const char                       DumpInfo::name[]             = "dump";
752 const char                       DumpInfo::shortDescription[] = "Make binary files human readable";
753 ICommandLineOptionsModulePointer DumpInfo::create()
754 {
755     return std::make_unique<Dump>();
756 }
757
758 } // namespace gmx