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/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/mdtypes/forcerec.h"
60 #include "gromacs/mdtypes/md_enums.h"
61 #include "gromacs/mdtypes/state.h"
62 #include "gromacs/topology/mtop_util.h"
63 #include "gromacs/topology/topology.h"
64 #include "gromacs/trajectory/trajectoryframe.h"
65 #include "gromacs/utility/arraysize.h"
66 #include "gromacs/utility/basedefinitions.h"
67 #include "gromacs/utility/fatalerror.h"
68 #include "gromacs/utility/futil.h"
69 #include "gromacs/utility/smalloc.h"
70 #include "gromacs/utility/txtdump.h"
72 static void list_tpx(const char *fn
, gmx_bool bShowNumbers
, const char *mdpfn
,
76 int indent
, i
, j
, **gcount
, atot
;
84 read_tpxheader(fn
, &tpx
, TRUE
, NULL
, NULL
);
89 tpx
.bTop
? &mtop
: NULL
);
93 gp
= gmx_fio_fopen(mdpfn
, "w");
94 pr_inputrec(gp
, 0, NULL
, &(ir
), TRUE
);
102 top
= gmx_mtop_t_to_t_topology(&mtop
);
105 if (available(stdout
, &tpx
, 0, fn
))
108 pr_title(stdout
, indent
, fn
);
109 pr_inputrec(stdout
, 0, "inputrec", tpx
.bIr
? &(ir
) : NULL
, FALSE
);
111 pr_tpxheader(stdout
, indent
, "header", &(tpx
));
115 pr_mtop(stdout
, indent
, "topology", &(mtop
), bShowNumbers
);
119 pr_top(stdout
, indent
, "topology", &(top
), bShowNumbers
);
122 pr_rvecs(stdout
, indent
, "box", tpx
.bBox
? state
.box
: NULL
, DIM
);
123 pr_rvecs(stdout
, indent
, "box_rel", tpx
.bBox
? state
.box_rel
: NULL
, DIM
);
124 pr_rvecs(stdout
, indent
, "boxv", tpx
.bBox
? state
.boxv
: NULL
, DIM
);
125 pr_rvecs(stdout
, indent
, "pres_prev", tpx
.bBox
? state
.pres_prev
: NULL
, DIM
);
126 pr_rvecs(stdout
, indent
, "svir_prev", tpx
.bBox
? state
.svir_prev
: NULL
, DIM
);
127 pr_rvecs(stdout
, indent
, "fvir_prev", tpx
.bBox
? state
.fvir_prev
: NULL
, DIM
);
128 /* leave nosehoover_xi in for now to match the tpr version */
129 pr_doubles(stdout
, indent
, "nosehoover_xi", state
.nosehoover_xi
, state
.ngtc
);
130 /*pr_doubles(stdout,indent,"nosehoover_vxi",state.nosehoover_vxi,state.ngtc);*/
131 /*pr_doubles(stdout,indent,"therm_integral",state.therm_integral,state.ngtc);*/
132 pr_rvecs(stdout
, indent
, "x", tpx
.bX
? state
.x
: NULL
, state
.natoms
);
133 pr_rvecs(stdout
, indent
, "v", tpx
.bV
? state
.v
: NULL
, 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
);
169 static void list_top(const char *fn
)
175 char *cppopts
[] = { NULL
};
177 status
= cpp_open_file(fn
, &handle
, cppopts
);
180 gmx_fatal(FARGS
, cpp_error(&handle
, status
));
184 status
= cpp_read_line(&handle
, BUFLEN
, buf
);
185 done
= (status
== eCPP_EOF
);
188 if (status
!= eCPP_OK
)
190 gmx_fatal(FARGS
, cpp_error(&handle
, status
));
199 status
= cpp_close_file(&handle
);
200 if (status
!= eCPP_OK
)
202 gmx_fatal(FARGS
, cpp_error(&handle
, status
));
206 static void list_trr(const char *fn
)
213 gmx_trr_header_t trrheader
;
216 fpread
= gmx_trr_open(fn
, "r");
219 while (gmx_trr_read_frame_header(fpread
, &trrheader
, &bOK
))
221 snew(x
, trrheader
.natoms
);
222 snew(v
, trrheader
.natoms
);
223 snew(f
, trrheader
.natoms
);
224 if (gmx_trr_read_frame_data(fpread
, &trrheader
,
225 trrheader
.box_size
? box
: NULL
,
226 trrheader
.x_size
? x
: NULL
,
227 trrheader
.v_size
? v
: NULL
,
228 trrheader
.f_size
? f
: NULL
))
230 sprintf(buf
, "%s frame %d", fn
, nframe
);
232 indent
= pr_title(stdout
, indent
, buf
);
233 pr_indent(stdout
, indent
);
234 fprintf(stdout
, "natoms=%10d step=%10d time=%12.7e lambda=%10g\n",
235 trrheader
.natoms
, trrheader
.step
, trrheader
.t
, trrheader
.lambda
);
236 if (trrheader
.box_size
)
238 pr_rvecs(stdout
, indent
, "box", box
, DIM
);
240 if (trrheader
.x_size
)
242 pr_rvecs(stdout
, indent
, "x", x
, trrheader
.natoms
);
244 if (trrheader
.v_size
)
246 pr_rvecs(stdout
, indent
, "v", v
, trrheader
.natoms
);
248 if (trrheader
.f_size
)
250 pr_rvecs(stdout
, indent
, "f", f
, trrheader
.natoms
);
255 fprintf(stderr
, "\nWARNING: Incomplete frame: nr %d, t=%g\n",
256 nframe
, trrheader
.t
);
266 fprintf(stderr
, "\nWARNING: Incomplete frame header: nr %d, t=%g\n",
267 nframe
, trrheader
.t
);
269 gmx_trr_close(fpread
);
272 void list_xtc(const char *fn
)
279 int nframe
, natoms
, step
;
283 xd
= open_xtc(fn
, "r");
284 read_first_xtc(xd
, &natoms
, &step
, &time
, box
, &x
, &prec
, &bOK
);
289 sprintf(buf
, "%s frame %d", fn
, nframe
);
291 indent
= pr_title(stdout
, indent
, buf
);
292 pr_indent(stdout
, indent
);
293 fprintf(stdout
, "natoms=%10d step=%10d time=%12.7e prec=%10g\n",
294 natoms
, step
, time
, prec
);
295 pr_rvecs(stdout
, indent
, "box", box
, DIM
);
296 pr_rvecs(stdout
, indent
, "x", x
, natoms
);
299 while (read_next_xtc(xd
, natoms
, &step
, &time
, box
, x
, &prec
, &bOK
));
302 fprintf(stderr
, "\nWARNING: Incomplete frame at time %g\n", time
);
308 /*! \brief Callback used by list_tng_for_gmx_dump. */
309 static void list_tng_inner(const char *fn
,
310 gmx_bool bFirstFrame
,
314 gmx_int64_t n_values_per_frame
,
325 sprintf(buf
, "%s frame %" GMX_PRId64
, fn
, nframe
);
327 indent
= pr_title(stdout
, indent
, buf
);
328 pr_indent(stdout
, indent
);
329 fprintf(stdout
, "natoms=%10" GMX_PRId64
" step=%10" GMX_PRId64
" time=%12.7e",
330 n_atoms
, step
, frame_time
);
333 fprintf(stdout
, " prec=%10g", prec
);
335 fprintf(stdout
, "\n");
337 pr_reals_of_dim(stdout
, indent
, block_name
, values
, n_atoms
, n_values_per_frame
);
340 static void list_tng(const char gmx_unused
*fn
)
343 tng_trajectory_t tng
;
344 gmx_int64_t nframe
= 0;
345 gmx_int64_t i
, *block_ids
= NULL
, step
, ndatablocks
;
348 gmx_tng_open(fn
, 'r', &tng
);
349 gmx_print_tng_molecule_system(tng
, stdout
);
351 bOK
= gmx_get_tng_data_block_types_of_next_frame(tng
, -1,
358 for (i
= 0; i
< ndatablocks
; i
++)
361 real prec
, *values
= NULL
;
362 gmx_int64_t n_values_per_frame
, n_atoms
;
363 char block_name
[STRLEN
];
365 gmx_get_tng_data_next_frame_of_block_type(tng
, block_ids
[i
], &values
,
367 &n_values_per_frame
, &n_atoms
,
369 block_name
, STRLEN
, &bOK
);
372 /* Can't write any output because we don't know what
374 fprintf(stderr
, "\nWARNING: Incomplete frame at time %g, will not write output\n", frame_time
);
378 list_tng_inner(fn
, (0 == i
), values
, step
, frame_time
,
379 n_values_per_frame
, n_atoms
, prec
, nframe
, block_name
);
384 while (gmx_get_tng_data_block_types_of_next_frame(tng
, step
,
400 void list_trx(const char *fn
)
414 fprintf(stderr
, "File %s is of an unsupported type. Try using the command\n 'less %s'\n",
419 void list_ene(const char *fn
)
423 gmx_enxnm_t
*enm
= NULL
;
428 printf("gmx dump: %s\n", fn
);
429 in
= open_enx(fn
, "r");
430 do_enxnms(in
, &nre
, &enm
);
433 printf("energy components:\n");
434 for (i
= 0; (i
< nre
); i
++)
436 printf("%5d %-24s (%s)\n", i
, enm
[i
].name
, enm
[i
].unit
);
442 bCont
= do_enx(in
, fr
);
446 printf("\n%24s %12.5e %12s %12s\n", "time:",
447 fr
->t
, "step:", gmx_step_str(fr
->step
, buf
));
448 printf("%24s %12s %12s %12s\n",
449 "", "", "nsteps:", gmx_step_str(fr
->nsteps
, buf
));
450 printf("%24s %12.5e %12s %12s\n",
451 "delta_t:", fr
->dt
, "sum steps:", gmx_step_str(fr
->nsum
, buf
));
454 printf("%24s %12s %12s %12s\n",
455 "Component", "Energy", "Av. Energy", "Sum Energy");
458 for (i
= 0; (i
< nre
); i
++)
460 printf("%24s %12.5e %12.5e %12.5e\n",
461 enm
[i
].name
, fr
->ener
[i
].e
, fr
->ener
[i
].eav
,
467 for (i
= 0; (i
< nre
); i
++)
469 printf("%24s %12.5e\n",
470 enm
[i
].name
, fr
->ener
[i
].e
);
474 for (b
= 0; b
< fr
->nblock
; b
++)
476 const char *typestr
= "";
478 t_enxblock
*eb
= &(fr
->block
[b
]);
479 printf("Block data %2d (%3d subblocks, id=%d)\n",
480 b
, eb
->nsub
, eb
->id
);
484 typestr
= enx_block_id_name
[eb
->id
];
486 printf(" id='%s'\n", typestr
);
487 for (i
= 0; i
< eb
->nsub
; i
++)
489 t_enxsubblock
*sb
= &(eb
->sub
[i
]);
490 printf(" Sub block %3d (%5d elems, type=%s) values:\n",
491 i
, sb
->nr
, xdr_datatype_names
[sb
->type
]);
495 case xdr_datatype_float
:
496 for (j
= 0; j
< sb
->nr
; j
++)
498 printf("%14d %8.4f\n", j
, sb
->fval
[j
]);
501 case xdr_datatype_double
:
502 for (j
= 0; j
< sb
->nr
; j
++)
504 printf("%14d %10.6f\n", j
, sb
->dval
[j
]);
507 case xdr_datatype_int
:
508 for (j
= 0; j
< sb
->nr
; j
++)
510 printf("%14d %10d\n", j
, sb
->ival
[j
]);
513 case xdr_datatype_int64
:
514 for (j
= 0; j
< sb
->nr
; j
++)
517 j
, gmx_step_str(sb
->lval
[j
], buf
));
520 case xdr_datatype_char
:
521 for (j
= 0; j
< sb
->nr
; j
++)
523 printf("%14d %1c\n", j
, sb
->cval
[j
]);
526 case xdr_datatype_string
:
527 for (j
= 0; j
< sb
->nr
; j
++)
529 printf("%14d %80s\n", j
, sb
->sval
[j
]);
533 gmx_incons("Unknown subblock type");
548 static void list_mtx(const char *fn
)
550 int nrow
, ncol
, i
, j
, k
;
551 real
*full
= NULL
, value
;
552 gmx_sparsematrix_t
* sparse
= NULL
;
554 gmx_mtxio_read(fn
, &nrow
, &ncol
, &full
, &sparse
);
558 snew(full
, nrow
*ncol
);
559 for (i
= 0; i
< nrow
*ncol
; i
++)
564 for (i
= 0; i
< sparse
->nrow
; i
++)
566 for (j
= 0; j
< sparse
->ndata
[i
]; j
++)
568 k
= sparse
->data
[i
][j
].col
;
569 value
= sparse
->data
[i
][j
].value
;
570 full
[i
*ncol
+k
] = value
;
571 full
[k
*ncol
+i
] = value
;
574 gmx_sparsematrix_destroy(sparse
);
577 printf("%d %d\n", nrow
, ncol
);
578 for (i
= 0; i
< nrow
; i
++)
580 for (j
= 0; j
< ncol
; j
++)
582 printf(" %g", full
[i
*ncol
+j
]);
590 int gmx_dump(int argc
, char *argv
[])
592 const char *desc
[] = {
593 "[THISMODULE] reads a run input file ([REF].tpr[ref]),",
594 "a trajectory ([REF].trr[ref]/[REF].xtc[ref]/[TT]/tng[tt]), an energy",
595 "file ([REF].edr[ref]) or a checkpoint file ([REF].cpt[ref])",
596 "and prints that to standard output in a readable format.",
597 "This program is essential for checking your run input file in case of",
599 "The program can also preprocess a topology to help finding problems.",
600 "Note that currently setting [TT]GMXLIB[tt] is the only way to customize",
601 "directories used for searching include files.",
603 const char *bugs
[] = {
604 "Position restraint output from -sys -s is broken"
607 { efTPR
, "-s", NULL
, ffOPTRD
},
608 { efTRX
, "-f", NULL
, ffOPTRD
},
609 { efEDR
, "-e", NULL
, ffOPTRD
},
610 { efCPT
, NULL
, NULL
, ffOPTRD
},
611 { efTOP
, "-p", NULL
, ffOPTRD
},
612 { efMTX
, "-mtx", "hessian", ffOPTRD
},
613 { efMDP
, "-om", NULL
, ffOPTWR
}
615 #define NFILE asize(fnm)
617 gmx_output_env_t
*oenv
;
618 /* Command line options */
619 static gmx_bool bShowNumbers
= TRUE
;
620 static gmx_bool bSysTop
= FALSE
;
622 { "-nr", FALSE
, etBOOL
, {&bShowNumbers
}, "Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" },
623 { "-sys", FALSE
, etBOOL
, {&bSysTop
}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" }
626 if (!parse_common_args(&argc
, argv
, 0, NFILE
, fnm
, asize(pa
), pa
,
627 asize(desc
), desc
, asize(bugs
), bugs
, &oenv
))
633 if (ftp2bSet(efTPR
, NFILE
, fnm
))
635 list_tpx(ftp2fn(efTPR
, NFILE
, fnm
), bShowNumbers
,
636 ftp2fn_null(efMDP
, NFILE
, fnm
), bSysTop
);
638 else if (ftp2bSet(efTRX
, NFILE
, fnm
))
640 list_trx(ftp2fn(efTRX
, NFILE
, fnm
));
642 else if (ftp2bSet(efEDR
, NFILE
, fnm
))
644 list_ene(ftp2fn(efEDR
, NFILE
, fnm
));
646 else if (ftp2bSet(efCPT
, NFILE
, fnm
))
648 list_checkpoint(ftp2fn(efCPT
, NFILE
, fnm
), stdout
);
650 else if (ftp2bSet(efTOP
, NFILE
, fnm
))
652 list_top(ftp2fn(efTOP
, NFILE
, fnm
));
654 else if (ftp2bSet(efMTX
, NFILE
, fnm
))
656 list_mtx(ftp2fn(efMTX
, NFILE
, fnm
));