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, 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.
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/fileio/checkpoint.h"
48 #include "gromacs/fileio/enxio.h"
49 #include "gromacs/fileio/gmxfio.h"
50 #include "gromacs/fileio/mtxio.h"
51 #include "gromacs/fileio/tngio.h"
52 #include "gromacs/fileio/tngio_for_tools.h"
53 #include "gromacs/fileio/tpxio.h"
54 #include "gromacs/fileio/trrio.h"
55 #include "gromacs/fileio/xtcio.h"
56 #include "gromacs/gmxpreprocess/gmxcpp.h"
57 #include "gromacs/linearalgebra/sparsematrix.h"
58 #include "gromacs/math/vecdump.h"
59 #include "gromacs/mdrunutility/mdmodules.h"
60 #include "gromacs/mdtypes/forcerec.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/mdtypes/md_enums.h"
63 #include "gromacs/mdtypes/state.h"
64 #include "gromacs/topology/mtop_util.h"
65 #include "gromacs/topology/topology.h"
66 #include "gromacs/trajectory/trajectoryframe.h"
67 #include "gromacs/utility/arraysize.h"
68 #include "gromacs/utility/basedefinitions.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/futil.h"
71 #include "gromacs/utility/smalloc.h"
72 #include "gromacs/utility/txtdump.h"
74 static void list_tpx(const char *fn, gmx_bool bShowNumbers, const char *mdpfn,
78 int indent, i, j, **gcount, atot;
80 t_inputrec *ir = nullptr;
86 read_tpxheader(fn, &tpx, TRUE);
87 gmx::MDModules mdModules;
90 ir = mdModules.inputrec();
95 tpx.bTop ? &mtop : NULL);
99 gp = gmx_fio_fopen(mdpfn, "w");
100 pr_inputrec(gp, 0, NULL, ir, TRUE);
108 top = gmx_mtop_t_to_t_topology(&mtop, false);
111 if (available(stdout, &tpx, 0, fn))
114 pr_title(stdout, indent, fn);
115 pr_inputrec(stdout, 0, "inputrec", ir, FALSE);
117 pr_tpxheader(stdout, indent, "header", &(tpx));
121 pr_mtop(stdout, indent, "topology", &(mtop), bShowNumbers);
125 pr_top(stdout, indent, "topology", &(top), bShowNumbers);
128 pr_rvecs(stdout, indent, "box", tpx.bBox ? state.box : NULL, DIM);
129 pr_rvecs(stdout, indent, "box_rel", tpx.bBox ? state.box_rel : NULL, DIM);
130 pr_rvecs(stdout, indent, "boxv", tpx.bBox ? state.boxv : NULL, DIM);
131 pr_rvecs(stdout, indent, "pres_prev", tpx.bBox ? state.pres_prev : NULL, DIM);
132 pr_rvecs(stdout, indent, "svir_prev", tpx.bBox ? state.svir_prev : NULL, DIM);
133 pr_rvecs(stdout, indent, "fvir_prev", tpx.bBox ? state.fvir_prev : NULL, DIM);
134 /* leave nosehoover_xi in for now to match the tpr version */
135 pr_doubles(stdout, indent, "nosehoover_xi", state.nosehoover_xi.data(), state.ngtc);
136 /*pr_doubles(stdout,indent,"nosehoover_vxi",state.nosehoover_vxi,state.ngtc);*/
137 /*pr_doubles(stdout,indent,"therm_integral",state.therm_integral,state.ngtc);*/
138 pr_rvecs(stdout, indent, "x", tpx.bX ? as_rvec_array(state.x.data()) : NULL, state.natoms);
139 pr_rvecs(stdout, indent, "v", tpx.bV ? as_rvec_array(state.v.data()) : NULL, state.natoms);
142 groups = &mtop.groups;
145 for (i = 0; (i < egcNR); i++)
147 snew(gcount[i], groups->grps[i].nr);
150 for (i = 0; (i < mtop.natoms); i++)
152 for (j = 0; (j < egcNR); j++)
154 gcount[j][ggrpnr(groups, j, i)]++;
157 printf("Group statistics\n");
158 for (i = 0; (i < egcNR); i++)
161 printf("%-12s: ", gtypes[i]);
162 for (j = 0; (j < groups->grps[i].nr); j++)
164 printf(" %5d", gcount[i][j]);
165 atot += gcount[i][j];
167 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_trr(const char *fn)
218 gmx_trr_header_t trrheader;
221 fpread = gmx_trr_open(fn, "r");
224 while (gmx_trr_read_frame_header(fpread, &trrheader, &bOK))
226 snew(x, trrheader.natoms);
227 snew(v, trrheader.natoms);
228 snew(f, trrheader.natoms);
229 if (gmx_trr_read_frame_data(fpread, &trrheader,
230 trrheader.box_size ? box : NULL,
231 trrheader.x_size ? x : NULL,
232 trrheader.v_size ? v : NULL,
233 trrheader.f_size ? f : NULL))
235 sprintf(buf, "%s frame %d", fn, nframe);
237 indent = pr_title(stdout, indent, buf);
238 pr_indent(stdout, indent);
239 fprintf(stdout, "natoms=%10d step=%10" GMX_PRId64 " time=%12.7e lambda=%10g\n",
240 trrheader.natoms, trrheader.step, trrheader.t, trrheader.lambda);
241 if (trrheader.box_size)
243 pr_rvecs(stdout, indent, "box", box, DIM);
245 if (trrheader.x_size)
247 pr_rvecs(stdout, indent, "x", x, trrheader.natoms);
249 if (trrheader.v_size)
251 pr_rvecs(stdout, indent, "v", v, trrheader.natoms);
253 if (trrheader.f_size)
255 pr_rvecs(stdout, indent, "f", f, trrheader.natoms);
260 fprintf(stderr, "\nWARNING: Incomplete frame: nr %d, t=%g\n",
261 nframe, trrheader.t);
271 fprintf(stderr, "\nWARNING: Incomplete frame header: nr %d, t=%g\n",
272 nframe, trrheader.t);
274 gmx_trr_close(fpread);
277 void list_xtc(const char *fn)
289 xd = open_xtc(fn, "r");
290 read_first_xtc(xd, &natoms, &step, &time, box, &x, &prec, &bOK);
295 sprintf(buf, "%s frame %d", fn, nframe);
297 indent = pr_title(stdout, indent, buf);
298 pr_indent(stdout, indent);
299 fprintf(stdout, "natoms=%10d step=%10" GMX_PRId64 " time=%12.7e prec=%10g\n",
300 natoms, step, time, prec);
301 pr_rvecs(stdout, indent, "box", box, DIM);
302 pr_rvecs(stdout, indent, "x", x, natoms);
305 while (read_next_xtc(xd, natoms, &step, &time, box, x, &prec, &bOK));
308 fprintf(stderr, "\nWARNING: Incomplete frame at time %g\n", time);
314 /*! \brief Callback used by list_tng_for_gmx_dump. */
315 static void list_tng_inner(const char *fn,
316 gmx_bool bFirstFrame,
320 gmx_int64_t n_values_per_frame,
331 sprintf(buf, "%s frame %" GMX_PRId64, fn, nframe);
333 indent = pr_title(stdout, indent, buf);
334 pr_indent(stdout, indent);
335 fprintf(stdout, "natoms=%10" GMX_PRId64 " step=%10" GMX_PRId64 " time=%12.7e",
336 n_atoms, step, frame_time);
339 fprintf(stdout, " prec=%10g", prec);
341 fprintf(stdout, "\n");
343 pr_reals_of_dim(stdout, indent, block_name, values, n_atoms, n_values_per_frame);
346 static void list_tng(const char gmx_unused *fn)
349 tng_trajectory_t tng;
350 gmx_int64_t nframe = 0;
351 gmx_int64_t i, *block_ids = NULL, step, ndatablocks;
354 gmx_tng_open(fn, 'r', &tng);
355 gmx_print_tng_molecule_system(tng, stdout);
357 bOK = gmx_get_tng_data_block_types_of_next_frame(tng, -1,
364 for (i = 0; i < ndatablocks; i++)
367 real prec, *values = NULL;
368 gmx_int64_t n_values_per_frame, n_atoms;
369 char block_name[STRLEN];
371 gmx_get_tng_data_next_frame_of_block_type(tng, block_ids[i], &values,
373 &n_values_per_frame, &n_atoms,
375 block_name, STRLEN, &bOK);
378 /* Can't write any output because we don't know what
380 fprintf(stderr, "\nWARNING: Incomplete frame at time %g, will not write output\n", frame_time);
384 list_tng_inner(fn, (0 == i), values, step, frame_time,
385 n_values_per_frame, n_atoms, prec, nframe, block_name);
390 while (gmx_get_tng_data_block_types_of_next_frame(tng, step,
406 void list_trx(const char *fn)
420 fprintf(stderr, "File %s is of an unsupported type. Try using the command\n 'less %s'\n",
425 void list_ene(const char *fn)
429 gmx_enxnm_t *enm = NULL;
434 printf("gmx dump: %s\n", fn);
435 in = open_enx(fn, "r");
436 do_enxnms(in, &nre, &enm);
439 printf("energy components:\n");
440 for (i = 0; (i < nre); i++)
442 printf("%5d %-24s (%s)\n", i, enm[i].name, enm[i].unit);
448 bCont = do_enx(in, fr);
452 printf("\n%24s %12.5e %12s %12s\n", "time:",
453 fr->t, "step:", gmx_step_str(fr->step, buf));
454 printf("%24s %12s %12s %12s\n",
455 "", "", "nsteps:", gmx_step_str(fr->nsteps, buf));
456 printf("%24s %12.5e %12s %12s\n",
457 "delta_t:", fr->dt, "sum steps:", gmx_step_str(fr->nsum, buf));
460 printf("%24s %12s %12s %12s\n",
461 "Component", "Energy", "Av. Energy", "Sum Energy");
464 for (i = 0; (i < nre); i++)
466 printf("%24s %12.5e %12.5e %12.5e\n",
467 enm[i].name, fr->ener[i].e, fr->ener[i].eav,
473 for (i = 0; (i < nre); i++)
475 printf("%24s %12.5e\n",
476 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",
486 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",
497 i, sb->nr, 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++)
523 j, gmx_step_str(sb->lval[j], buf));
526 case xdr_datatype_char:
527 for (j = 0; j < sb->nr; j++)
529 printf("%14d %1c\n", j, sb->cval[j]);
532 case xdr_datatype_string:
533 for (j = 0; j < sb->nr; j++)
535 printf("%14d %80s\n", j, sb->sval[j]);
539 gmx_incons("Unknown subblock type");
554 static void list_mtx(const char *fn)
556 int nrow, ncol, i, j, k;
557 real *full = NULL, value;
558 gmx_sparsematrix_t * sparse = NULL;
560 gmx_mtxio_read(fn, &nrow, &ncol, &full, &sparse);
564 snew(full, nrow*ncol);
565 for (i = 0; i < nrow*ncol; i++)
570 for (i = 0; i < sparse->nrow; i++)
572 for (j = 0; j < sparse->ndata[i]; j++)
574 k = sparse->data[i][j].col;
575 value = sparse->data[i][j].value;
576 full[i*ncol+k] = value;
577 full[k*ncol+i] = value;
580 gmx_sparsematrix_destroy(sparse);
583 printf("%d %d\n", nrow, ncol);
584 for (i = 0; i < nrow; i++)
586 for (j = 0; j < ncol; j++)
588 printf(" %g", full[i*ncol+j]);
596 int gmx_dump(int argc, char *argv[])
598 const char *desc[] = {
599 "[THISMODULE] reads a run input file ([REF].tpr[ref]),",
600 "a trajectory ([REF].trr[ref]/[REF].xtc[ref]/[TT]/tng[tt]), an energy",
601 "file ([REF].edr[ref]) or a checkpoint file ([REF].cpt[ref])",
602 "and prints that to standard output in a readable format.",
603 "This program is essential for checking your run input file in case of",
605 "The program can also preprocess a topology to help finding problems.",
606 "Note that currently setting [TT]GMXLIB[tt] is the only way to customize",
607 "directories used for searching include files.",
609 const char *bugs[] = {
610 "Position restraint output from -sys -s is broken"
613 { efTPR, "-s", NULL, ffOPTRD },
614 { efTRX, "-f", NULL, ffOPTRD },
615 { efEDR, "-e", NULL, ffOPTRD },
616 { efCPT, NULL, NULL, ffOPTRD },
617 { efTOP, "-p", NULL, ffOPTRD },
618 { efMTX, "-mtx", "hessian", ffOPTRD },
619 { efMDP, "-om", NULL, ffOPTWR }
621 #define NFILE asize(fnm)
623 gmx_output_env_t *oenv;
624 /* Command line options */
625 static gmx_bool bShowNumbers = TRUE;
626 static gmx_bool bSysTop = FALSE;
628 { "-nr", FALSE, etBOOL, {&bShowNumbers}, "Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" },
629 { "-sys", FALSE, etBOOL, {&bSysTop}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" }
632 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
633 asize(desc), desc, asize(bugs), bugs, &oenv))
639 if (ftp2bSet(efTPR, NFILE, fnm))
641 list_tpx(ftp2fn(efTPR, NFILE, fnm), bShowNumbers,
642 ftp2fn_null(efMDP, NFILE, fnm), bSysTop);
644 else if (ftp2bSet(efTRX, NFILE, fnm))
646 list_trx(ftp2fn(efTRX, NFILE, fnm));
648 else if (ftp2bSet(efEDR, NFILE, fnm))
650 list_ene(ftp2fn(efEDR, NFILE, fnm));
652 else if (ftp2bSet(efCPT, NFILE, fnm))
654 list_checkpoint(ftp2fn(efCPT, NFILE, fnm), stdout);
656 else if (ftp2bSet(efTOP, NFILE, fnm))
658 list_top(ftp2fn(efTOP, NFILE, fnm));
660 else if (ftp2bSet(efMTX, NFILE, fnm))
662 list_mtx(ftp2fn(efMTX, NFILE, fnm));