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, 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/enxio.h"
48 #include "gromacs/fileio/gmxfio.h"
49 #include "gromacs/fileio/mtxio.h"
50 #include "gromacs/fileio/tngio.h"
51 #include "gromacs/fileio/tngio_for_tools.h"
52 #include "gromacs/fileio/tpxio.h"
53 #include "gromacs/fileio/trnio.h"
54 #include "gromacs/fileio/xtcio.h"
55 #include "gromacs/gmxpreprocess/gmxcpp.h"
56 #include "gromacs/legacyheaders/checkpoint.h"
57 #include "gromacs/legacyheaders/macros.h"
58 #include "gromacs/legacyheaders/names.h"
59 #include "gromacs/legacyheaders/txtdump.h"
60 #include "gromacs/linearalgebra/sparsematrix.h"
61 #include "gromacs/topology/mtop_util.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/futil.h"
64 #include "gromacs/utility/smalloc.h"
66 static void list_tpx(const char *fn, gmx_bool bShowNumbers, const char *mdpfn,
70 int fp, indent, i, j, **gcount, atot;
79 read_tpxheader(fn, &tpx, TRUE, NULL, NULL);
83 &state, tpx.bF ? f : NULL,
84 tpx.bTop ? &mtop : NULL);
88 gp = gmx_fio_fopen(mdpfn, "w");
89 pr_inputrec(gp, 0, NULL, &(ir), TRUE);
97 top = gmx_mtop_t_to_t_topology(&mtop);
100 if (available(stdout, &tpx, 0, fn))
103 indent = pr_title(stdout, indent, fn);
104 pr_inputrec(stdout, 0, "inputrec", tpx.bIr ? &(ir) : NULL, FALSE);
107 pr_header(stdout, indent, "header", &(tpx));
111 pr_mtop(stdout, indent, "topology", &(mtop), bShowNumbers);
115 pr_top(stdout, indent, "topology", &(top), bShowNumbers);
118 pr_rvecs(stdout, indent, "box", tpx.bBox ? state.box : NULL, DIM);
119 pr_rvecs(stdout, indent, "box_rel", tpx.bBox ? state.box_rel : NULL, DIM);
120 pr_rvecs(stdout, indent, "boxv", tpx.bBox ? state.boxv : NULL, DIM);
121 pr_rvecs(stdout, indent, "pres_prev", tpx.bBox ? state.pres_prev : NULL, DIM);
122 pr_rvecs(stdout, indent, "svir_prev", tpx.bBox ? state.svir_prev : NULL, DIM);
123 pr_rvecs(stdout, indent, "fvir_prev", tpx.bBox ? state.fvir_prev : NULL, DIM);
124 /* leave nosehoover_xi in for now to match the tpr version */
125 pr_doubles(stdout, indent, "nosehoover_xi", state.nosehoover_xi, state.ngtc);
126 /*pr_doubles(stdout,indent,"nosehoover_vxi",state.nosehoover_vxi,state.ngtc);*/
127 /*pr_doubles(stdout,indent,"therm_integral",state.therm_integral,state.ngtc);*/
128 pr_rvecs(stdout, indent, "x", tpx.bX ? state.x : NULL, state.natoms);
129 pr_rvecs(stdout, indent, "v", tpx.bV ? state.v : NULL, state.natoms);
132 pr_rvecs(stdout, indent, "f", f, state.natoms);
136 groups = &mtop.groups;
139 for (i = 0; (i < egcNR); i++)
141 snew(gcount[i], groups->grps[i].nr);
144 for (i = 0; (i < mtop.natoms); i++)
146 for (j = 0; (j < egcNR); j++)
148 gcount[j][ggrpnr(groups, j, i)]++;
151 printf("Group statistics\n");
152 for (i = 0; (i < egcNR); i++)
155 printf("%-12s: ", gtypes[i]);
156 for (j = 0; (j < groups->grps[i].nr); j++)
158 printf(" %5d", gcount[i][j]);
159 atot += gcount[i][j];
161 printf(" (total %d atoms)\n", atot);
170 static void list_top(const char *fn)
176 char *cppopts[] = { NULL };
178 status = cpp_open_file(fn, &handle, cppopts);
181 gmx_fatal(FARGS, cpp_error(&handle, status));
185 status = cpp_read_line(&handle, BUFLEN, buf);
186 done = (status == eCPP_EOF);
189 if (status != eCPP_OK)
191 gmx_fatal(FARGS, cpp_error(&handle, status));
200 status = cpp_close_file(&handle);
201 if (status != eCPP_OK)
203 gmx_fatal(FARGS, cpp_error(&handle, status));
207 static void list_trn(const char *fn)
217 fpread = open_trn(fn, "r");
220 while (fread_trnheader(fpread, &trn, &bOK))
225 if (fread_htrn(fpread, &trn,
226 trn.box_size ? box : NULL,
227 trn.x_size ? x : NULL,
228 trn.v_size ? v : NULL,
229 trn.f_size ? f : NULL))
231 sprintf(buf, "%s frame %d", fn, nframe);
233 indent = pr_title(stdout, indent, buf);
234 pr_indent(stdout, indent);
235 fprintf(stdout, "natoms=%10d step=%10d time=%12.7e lambda=%10g\n",
236 trn.natoms, trn.step, trn.t, trn.lambda);
239 pr_rvecs(stdout, indent, "box", box, DIM);
243 pr_rvecs(stdout, indent, "x", x, trn.natoms);
247 pr_rvecs(stdout, indent, "v", v, trn.natoms);
251 pr_rvecs(stdout, indent, "f", f, trn.natoms);
256 fprintf(stderr, "\nWARNING: Incomplete frame: nr %d, t=%g\n",
267 fprintf(stderr, "\nWARNING: Incomplete frame header: nr %d, t=%g\n",
273 void list_xtc(const char *fn)
280 int nframe, natoms, step;
284 xd = open_xtc(fn, "r");
285 read_first_xtc(xd, &natoms, &step, &time, box, &x, &prec, &bOK);
290 sprintf(buf, "%s frame %d", fn, nframe);
292 indent = pr_title(stdout, indent, buf);
293 pr_indent(stdout, indent);
294 fprintf(stdout, "natoms=%10d step=%10d time=%12.7e prec=%10g\n",
295 natoms, step, time, prec);
296 pr_rvecs(stdout, indent, "box", box, DIM);
297 pr_rvecs(stdout, indent, "x", x, natoms);
300 while (read_next_xtc(xd, natoms, &step, &time, box, x, &prec, &bOK));
303 fprintf(stderr, "\nWARNING: Incomplete frame at time %g\n", time);
309 /*! \brief Callback used by list_tng_for_gmx_dump. */
310 static void list_tng_inner(const char *fn,
311 gmx_bool bFirstFrame,
315 gmx_int64_t n_values_per_frame,
326 sprintf(buf, "%s frame %" GMX_PRId64, fn, nframe);
328 indent = pr_title(stdout, indent, buf);
329 pr_indent(stdout, indent);
330 fprintf(stdout, "natoms=%10" GMX_PRId64 " step=%10" GMX_PRId64 " time=%12.7e",
331 n_atoms, step, frame_time);
334 fprintf(stdout, " prec=%10g", prec);
336 fprintf(stdout, "\n");
338 pr_reals_of_dim(stdout, indent, block_name, values, n_atoms, n_values_per_frame);
341 static void list_tng(const char gmx_unused *fn)
344 tng_trajectory_t tng;
345 gmx_int64_t nframe = 0;
346 gmx_int64_t i, *block_ids = NULL, step, ndatablocks;
349 gmx_tng_open(fn, 'r', &tng);
350 gmx_print_tng_molecule_system(tng, stdout);
352 bOK = gmx_get_tng_data_block_types_of_next_frame(tng, -1,
359 for (i = 0; i < ndatablocks; i++)
362 real prec, *values = NULL;
363 gmx_int64_t n_values_per_frame, n_atoms;
364 char block_name[STRLEN];
366 gmx_get_tng_data_next_frame_of_block_type(tng, block_ids[i], &values,
368 &n_values_per_frame, &n_atoms,
370 block_name, STRLEN, &bOK);
373 /* Can't write any output because we don't know what
375 fprintf(stderr, "\nWARNING: Incomplete frame at time %g, will not write output\n", frame_time);
379 list_tng_inner(fn, (0 == i), values, step, frame_time,
380 n_values_per_frame, n_atoms, prec, nframe, block_name);
385 while (gmx_get_tng_data_block_types_of_next_frame(tng, step,
401 void list_trx(const char *fn)
415 fprintf(stderr, "File %s is of an unsupported type. Try using the command\n 'less %s'\n",
420 void list_ene(const char *fn)
425 gmx_enxnm_t *enm = NULL;
431 printf("gmx dump: %s\n", fn);
432 in = open_enx(fn, "r");
433 do_enxnms(in, &nre, &enm);
436 printf("energy components:\n");
437 for (i = 0; (i < nre); i++)
439 printf("%5d %-24s (%s)\n", i, enm[i].name, enm[i].unit);
446 bCont = do_enx(in, fr);
450 printf("\n%24s %12.5e %12s %12s\n", "time:",
451 fr->t, "step:", gmx_step_str(fr->step, buf));
452 printf("%24s %12s %12s %12s\n",
453 "", "", "nsteps:", gmx_step_str(fr->nsteps, buf));
454 printf("%24s %12.5e %12s %12s\n",
455 "delta_t:", fr->dt, "sum steps:", gmx_step_str(fr->nsum, buf));
458 printf("%24s %12s %12s %12s\n",
459 "Component", "Energy", "Av. Energy", "Sum Energy");
462 for (i = 0; (i < nre); i++)
464 printf("%24s %12.5e %12.5e %12.5e\n",
465 enm[i].name, fr->ener[i].e, fr->ener[i].eav,
471 for (i = 0; (i < nre); i++)
473 printf("%24s %12.5e\n",
474 enm[i].name, fr->ener[i].e);
478 for (b = 0; b < fr->nblock; b++)
480 const char *typestr = "";
482 t_enxblock *eb = &(fr->block[b]);
483 printf("Block data %2d (%3d subblocks, id=%d)\n",
484 b, eb->nsub, eb->id);
488 typestr = enx_block_id_name[eb->id];
490 printf(" id='%s'\n", typestr);
491 for (i = 0; i < eb->nsub; i++)
493 t_enxsubblock *sb = &(eb->sub[i]);
494 printf(" Sub block %3d (%5d elems, type=%s) values:\n",
495 i, sb->nr, xdr_datatype_names[sb->type]);
499 case xdr_datatype_float:
500 for (j = 0; j < sb->nr; j++)
502 printf("%14d %8.4f\n", j, sb->fval[j]);
505 case xdr_datatype_double:
506 for (j = 0; j < sb->nr; j++)
508 printf("%14d %10.6f\n", j, sb->dval[j]);
511 case xdr_datatype_int:
512 for (j = 0; j < sb->nr; j++)
514 printf("%14d %10d\n", j, sb->ival[j]);
517 case xdr_datatype_int64:
518 for (j = 0; j < sb->nr; j++)
521 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]);
537 gmx_incons("Unknown subblock type");
552 static void list_mtx(const char *fn)
554 int nrow, ncol, i, j, k;
555 real *full = NULL, value;
556 gmx_sparsematrix_t * sparse = NULL;
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 int gmx_dump(int argc, char *argv[])
596 const char *desc[] = {
597 "[THISMODULE] reads a run input file ([REF].tpr[ref]),",
598 "a trajectory ([REF].trr[ref]/[REF].xtc[ref]/[TT]/tng[tt]), an energy",
599 "file ([REF].edr[ref]) or a checkpoint file ([REF].cpt[ref])",
600 "and prints that to standard output in a readable format.",
601 "This program is essential for checking your run input file in case of",
603 "The program can also preprocess a topology to help finding problems.",
604 "Note that currently setting [TT]GMXLIB[tt] is the only way to customize",
605 "directories used for searching include files.",
607 const char *bugs[] = {
608 "Position restraint output from -sys -s is broken"
611 { efTPR, "-s", NULL, ffOPTRD },
612 { efTRX, "-f", NULL, ffOPTRD },
613 { efEDR, "-e", NULL, ffOPTRD },
614 { efCPT, NULL, NULL, ffOPTRD },
615 { efTOP, "-p", NULL, ffOPTRD },
616 { efMTX, "-mtx", "hessian", ffOPTRD },
617 { efMDP, "-om", NULL, ffOPTWR }
619 #define NFILE asize(fnm)
622 /* Command line options */
623 static gmx_bool bShowNumbers = TRUE;
624 static gmx_bool bSysTop = FALSE;
626 { "-nr", FALSE, etBOOL, {&bShowNumbers}, "Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" },
627 { "-sys", FALSE, etBOOL, {&bSysTop}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" }
630 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
631 asize(desc), desc, asize(bugs), bugs, &oenv))
637 if (ftp2bSet(efTPR, NFILE, fnm))
639 list_tpx(ftp2fn(efTPR, NFILE, fnm), bShowNumbers,
640 ftp2fn_null(efMDP, NFILE, fnm), bSysTop);
642 else if (ftp2bSet(efTRX, NFILE, fnm))
644 list_trx(ftp2fn(efTRX, NFILE, fnm));
646 else if (ftp2bSet(efEDR, NFILE, fnm))
648 list_ene(ftp2fn(efEDR, NFILE, fnm));
650 else if (ftp2bSet(efCPT, NFILE, fnm))
652 list_checkpoint(ftp2fn(efCPT, NFILE, fnm), stdout);
654 else if (ftp2bSet(efTOP, NFILE, fnm))
656 list_top(ftp2fn(efTOP, NFILE, fnm));
658 else if (ftp2bSet(efMTX, NFILE, fnm))
660 list_mtx(ftp2fn(efMTX, NFILE, fnm));