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/trx.h"
56 #include "gromacs/fileio/xtcio.h"
57 #include "gromacs/gmxpreprocess/gmxcpp.h"
58 #include "gromacs/legacyheaders/names.h"
59 #include "gromacs/legacyheaders/txtdump.h"
60 #include "gromacs/linearalgebra/sparsematrix.h"
61 #include "gromacs/mdtypes/state.h"
62 #include "gromacs/topology/mtop_util.h"
63 #include "gromacs/utility/arraysize.h"
64 #include "gromacs/utility/fatalerror.h"
65 #include "gromacs/utility/futil.h"
66 #include "gromacs/utility/smalloc.h"
68 static void list_tpx(const char *fn
, gmx_bool bShowNumbers
, const char *mdpfn
,
72 int indent
, i
, j
, **gcount
, atot
;
80 read_tpxheader(fn
, &tpx
, TRUE
, NULL
, NULL
);
85 tpx
.bTop
? &mtop
: NULL
);
89 gp
= gmx_fio_fopen(mdpfn
, "w");
90 pr_inputrec(gp
, 0, NULL
, &(ir
), TRUE
);
98 top
= gmx_mtop_t_to_t_topology(&mtop
);
101 if (available(stdout
, &tpx
, 0, fn
))
104 pr_title(stdout
, indent
, fn
);
105 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 groups
= &mtop
.groups
;
135 for (i
= 0; (i
< egcNR
); i
++)
137 snew(gcount
[i
], groups
->grps
[i
].nr
);
140 for (i
= 0; (i
< mtop
.natoms
); i
++)
142 for (j
= 0; (j
< egcNR
); j
++)
144 gcount
[j
][ggrpnr(groups
, j
, i
)]++;
147 printf("Group statistics\n");
148 for (i
= 0; (i
< egcNR
); i
++)
151 printf("%-12s: ", gtypes
[i
]);
152 for (j
= 0; (j
< groups
->grps
[i
].nr
); j
++)
154 printf(" %5d", gcount
[i
][j
]);
155 atot
+= gcount
[i
][j
];
157 printf(" (total %d atoms)\n", atot
);
165 static void list_top(const char *fn
)
171 char *cppopts
[] = { NULL
};
173 status
= cpp_open_file(fn
, &handle
, cppopts
);
176 gmx_fatal(FARGS
, cpp_error(&handle
, status
));
180 status
= cpp_read_line(&handle
, BUFLEN
, buf
);
181 done
= (status
== eCPP_EOF
);
184 if (status
!= eCPP_OK
)
186 gmx_fatal(FARGS
, cpp_error(&handle
, status
));
195 status
= cpp_close_file(&handle
);
196 if (status
!= eCPP_OK
)
198 gmx_fatal(FARGS
, cpp_error(&handle
, status
));
202 static void list_trr(const char *fn
)
209 gmx_trr_header_t trrheader
;
212 fpread
= gmx_trr_open(fn
, "r");
215 while (gmx_trr_read_frame_header(fpread
, &trrheader
, &bOK
))
217 snew(x
, trrheader
.natoms
);
218 snew(v
, trrheader
.natoms
);
219 snew(f
, trrheader
.natoms
);
220 if (gmx_trr_read_frame_data(fpread
, &trrheader
,
221 trrheader
.box_size
? box
: NULL
,
222 trrheader
.x_size
? x
: NULL
,
223 trrheader
.v_size
? v
: NULL
,
224 trrheader
.f_size
? f
: NULL
))
226 sprintf(buf
, "%s frame %d", fn
, nframe
);
228 indent
= pr_title(stdout
, indent
, buf
);
229 pr_indent(stdout
, indent
);
230 fprintf(stdout
, "natoms=%10d step=%10d time=%12.7e lambda=%10g\n",
231 trrheader
.natoms
, trrheader
.step
, trrheader
.t
, trrheader
.lambda
);
232 if (trrheader
.box_size
)
234 pr_rvecs(stdout
, indent
, "box", box
, DIM
);
236 if (trrheader
.x_size
)
238 pr_rvecs(stdout
, indent
, "x", x
, trrheader
.natoms
);
240 if (trrheader
.v_size
)
242 pr_rvecs(stdout
, indent
, "v", v
, trrheader
.natoms
);
244 if (trrheader
.f_size
)
246 pr_rvecs(stdout
, indent
, "f", f
, trrheader
.natoms
);
251 fprintf(stderr
, "\nWARNING: Incomplete frame: nr %d, t=%g\n",
252 nframe
, trrheader
.t
);
262 fprintf(stderr
, "\nWARNING: Incomplete frame header: nr %d, t=%g\n",
263 nframe
, trrheader
.t
);
265 gmx_trr_close(fpread
);
268 void list_xtc(const char *fn
)
275 int nframe
, natoms
, step
;
279 xd
= open_xtc(fn
, "r");
280 read_first_xtc(xd
, &natoms
, &step
, &time
, box
, &x
, &prec
, &bOK
);
285 sprintf(buf
, "%s frame %d", fn
, nframe
);
287 indent
= pr_title(stdout
, indent
, buf
);
288 pr_indent(stdout
, indent
);
289 fprintf(stdout
, "natoms=%10d step=%10d time=%12.7e prec=%10g\n",
290 natoms
, step
, time
, prec
);
291 pr_rvecs(stdout
, indent
, "box", box
, DIM
);
292 pr_rvecs(stdout
, indent
, "x", x
, natoms
);
295 while (read_next_xtc(xd
, natoms
, &step
, &time
, box
, x
, &prec
, &bOK
));
298 fprintf(stderr
, "\nWARNING: Incomplete frame at time %g\n", time
);
304 /*! \brief Callback used by list_tng_for_gmx_dump. */
305 static void list_tng_inner(const char *fn
,
306 gmx_bool bFirstFrame
,
310 gmx_int64_t n_values_per_frame
,
321 sprintf(buf
, "%s frame %" GMX_PRId64
, fn
, nframe
);
323 indent
= pr_title(stdout
, indent
, buf
);
324 pr_indent(stdout
, indent
);
325 fprintf(stdout
, "natoms=%10" GMX_PRId64
" step=%10" GMX_PRId64
" time=%12.7e",
326 n_atoms
, step
, frame_time
);
329 fprintf(stdout
, " prec=%10g", prec
);
331 fprintf(stdout
, "\n");
333 pr_reals_of_dim(stdout
, indent
, block_name
, values
, n_atoms
, n_values_per_frame
);
336 static void list_tng(const char gmx_unused
*fn
)
339 tng_trajectory_t tng
;
340 gmx_int64_t nframe
= 0;
341 gmx_int64_t i
, *block_ids
= NULL
, step
, ndatablocks
;
344 gmx_tng_open(fn
, 'r', &tng
);
345 gmx_print_tng_molecule_system(tng
, stdout
);
347 bOK
= gmx_get_tng_data_block_types_of_next_frame(tng
, -1,
354 for (i
= 0; i
< ndatablocks
; i
++)
357 real prec
, *values
= NULL
;
358 gmx_int64_t n_values_per_frame
, n_atoms
;
359 char block_name
[STRLEN
];
361 gmx_get_tng_data_next_frame_of_block_type(tng
, block_ids
[i
], &values
,
363 &n_values_per_frame
, &n_atoms
,
365 block_name
, STRLEN
, &bOK
);
368 /* Can't write any output because we don't know what
370 fprintf(stderr
, "\nWARNING: Incomplete frame at time %g, will not write output\n", frame_time
);
374 list_tng_inner(fn
, (0 == i
), values
, step
, frame_time
,
375 n_values_per_frame
, n_atoms
, prec
, nframe
, block_name
);
380 while (gmx_get_tng_data_block_types_of_next_frame(tng
, step
,
396 void list_trx(const char *fn
)
410 fprintf(stderr
, "File %s is of an unsupported type. Try using the command\n 'less %s'\n",
415 void list_ene(const char *fn
)
419 gmx_enxnm_t
*enm
= NULL
;
424 printf("gmx dump: %s\n", fn
);
425 in
= open_enx(fn
, "r");
426 do_enxnms(in
, &nre
, &enm
);
429 printf("energy components:\n");
430 for (i
= 0; (i
< nre
); i
++)
432 printf("%5d %-24s (%s)\n", i
, enm
[i
].name
, enm
[i
].unit
);
438 bCont
= do_enx(in
, fr
);
442 printf("\n%24s %12.5e %12s %12s\n", "time:",
443 fr
->t
, "step:", gmx_step_str(fr
->step
, buf
));
444 printf("%24s %12s %12s %12s\n",
445 "", "", "nsteps:", gmx_step_str(fr
->nsteps
, buf
));
446 printf("%24s %12.5e %12s %12s\n",
447 "delta_t:", fr
->dt
, "sum steps:", gmx_step_str(fr
->nsum
, buf
));
450 printf("%24s %12s %12s %12s\n",
451 "Component", "Energy", "Av. Energy", "Sum Energy");
454 for (i
= 0; (i
< nre
); i
++)
456 printf("%24s %12.5e %12.5e %12.5e\n",
457 enm
[i
].name
, fr
->ener
[i
].e
, fr
->ener
[i
].eav
,
463 for (i
= 0; (i
< nre
); i
++)
465 printf("%24s %12.5e\n",
466 enm
[i
].name
, fr
->ener
[i
].e
);
470 for (b
= 0; b
< fr
->nblock
; b
++)
472 const char *typestr
= "";
474 t_enxblock
*eb
= &(fr
->block
[b
]);
475 printf("Block data %2d (%3d subblocks, id=%d)\n",
476 b
, eb
->nsub
, eb
->id
);
480 typestr
= enx_block_id_name
[eb
->id
];
482 printf(" id='%s'\n", typestr
);
483 for (i
= 0; i
< eb
->nsub
; i
++)
485 t_enxsubblock
*sb
= &(eb
->sub
[i
]);
486 printf(" Sub block %3d (%5d elems, type=%s) values:\n",
487 i
, sb
->nr
, xdr_datatype_names
[sb
->type
]);
491 case xdr_datatype_float
:
492 for (j
= 0; j
< sb
->nr
; j
++)
494 printf("%14d %8.4f\n", j
, sb
->fval
[j
]);
497 case xdr_datatype_double
:
498 for (j
= 0; j
< sb
->nr
; j
++)
500 printf("%14d %10.6f\n", j
, sb
->dval
[j
]);
503 case xdr_datatype_int
:
504 for (j
= 0; j
< sb
->nr
; j
++)
506 printf("%14d %10d\n", j
, sb
->ival
[j
]);
509 case xdr_datatype_int64
:
510 for (j
= 0; j
< sb
->nr
; j
++)
513 j
, gmx_step_str(sb
->lval
[j
], buf
));
516 case xdr_datatype_char
:
517 for (j
= 0; j
< sb
->nr
; j
++)
519 printf("%14d %1c\n", j
, sb
->cval
[j
]);
522 case xdr_datatype_string
:
523 for (j
= 0; j
< sb
->nr
; j
++)
525 printf("%14d %80s\n", j
, sb
->sval
[j
]);
529 gmx_incons("Unknown subblock type");
544 static void list_mtx(const char *fn
)
546 int nrow
, ncol
, i
, j
, k
;
547 real
*full
= NULL
, value
;
548 gmx_sparsematrix_t
* sparse
= NULL
;
550 gmx_mtxio_read(fn
, &nrow
, &ncol
, &full
, &sparse
);
554 snew(full
, nrow
*ncol
);
555 for (i
= 0; i
< nrow
*ncol
; i
++)
560 for (i
= 0; i
< sparse
->nrow
; i
++)
562 for (j
= 0; j
< sparse
->ndata
[i
]; j
++)
564 k
= sparse
->data
[i
][j
].col
;
565 value
= sparse
->data
[i
][j
].value
;
566 full
[i
*ncol
+k
] = value
;
567 full
[k
*ncol
+i
] = value
;
570 gmx_sparsematrix_destroy(sparse
);
573 printf("%d %d\n", nrow
, ncol
);
574 for (i
= 0; i
< nrow
; i
++)
576 for (j
= 0; j
< ncol
; j
++)
578 printf(" %g", full
[i
*ncol
+j
]);
586 int gmx_dump(int argc
, char *argv
[])
588 const char *desc
[] = {
589 "[THISMODULE] reads a run input file ([REF].tpr[ref]),",
590 "a trajectory ([REF].trr[ref]/[REF].xtc[ref]/[TT]/tng[tt]), an energy",
591 "file ([REF].edr[ref]) or a checkpoint file ([REF].cpt[ref])",
592 "and prints that to standard output in a readable format.",
593 "This program is essential for checking your run input file in case of",
595 "The program can also preprocess a topology to help finding problems.",
596 "Note that currently setting [TT]GMXLIB[tt] is the only way to customize",
597 "directories used for searching include files.",
599 const char *bugs
[] = {
600 "Position restraint output from -sys -s is broken"
603 { efTPR
, "-s", NULL
, ffOPTRD
},
604 { efTRX
, "-f", NULL
, ffOPTRD
},
605 { efEDR
, "-e", NULL
, ffOPTRD
},
606 { efCPT
, NULL
, NULL
, ffOPTRD
},
607 { efTOP
, "-p", NULL
, ffOPTRD
},
608 { efMTX
, "-mtx", "hessian", ffOPTRD
},
609 { efMDP
, "-om", NULL
, ffOPTWR
}
611 #define NFILE asize(fnm)
613 gmx_output_env_t
*oenv
;
614 /* Command line options */
615 static gmx_bool bShowNumbers
= TRUE
;
616 static gmx_bool bSysTop
= FALSE
;
618 { "-nr", FALSE
, etBOOL
, {&bShowNumbers
}, "Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" },
619 { "-sys", FALSE
, etBOOL
, {&bSysTop
}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" }
622 if (!parse_common_args(&argc
, argv
, 0, NFILE
, fnm
, asize(pa
), pa
,
623 asize(desc
), desc
, asize(bugs
), bugs
, &oenv
))
629 if (ftp2bSet(efTPR
, NFILE
, fnm
))
631 list_tpx(ftp2fn(efTPR
, NFILE
, fnm
), bShowNumbers
,
632 ftp2fn_null(efMDP
, NFILE
, fnm
), bSysTop
);
634 else if (ftp2bSet(efTRX
, NFILE
, fnm
))
636 list_trx(ftp2fn(efTRX
, NFILE
, fnm
));
638 else if (ftp2bSet(efEDR
, NFILE
, fnm
))
640 list_ene(ftp2fn(efEDR
, NFILE
, fnm
));
642 else if (ftp2bSet(efCPT
, NFILE
, fnm
))
644 list_checkpoint(ftp2fn(efCPT
, NFILE
, fnm
), stdout
);
646 else if (ftp2bSet(efTOP
, NFILE
, fnm
))
648 list_top(ftp2fn(efTOP
, NFILE
, fnm
));
650 else if (ftp2bSet(efMTX
, NFILE
, fnm
))
652 list_mtx(ftp2fn(efMTX
, NFILE
, fnm
));