Split txtdump.*, part 1
[gromacs.git] / src / gromacs / tools / dump.cpp
blob30309902c2332f9d4b16726e3454043092f1cf2b
1 /*
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.
37 #include "gmxpre.h"
39 #include "config.h"
41 #include <cassert>
42 #include <cmath>
43 #include <cstdio>
44 #include <cstring>
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/txtdump.h"
57 #include "gromacs/fileio/xtcio.h"
58 #include "gromacs/gmxpreprocess/gmxcpp.h"
59 #include "gromacs/linearalgebra/sparsematrix.h"
60 #include "gromacs/mdtypes/forcerec.h"
61 #include "gromacs/mdtypes/md_enums.h"
62 #include "gromacs/mdtypes/state.h"
63 #include "gromacs/topology/mtop_util.h"
64 #include "gromacs/topology/topology.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"
71 static void list_tpx(const char *fn, gmx_bool bShowNumbers, const char *mdpfn,
72 gmx_bool bSysTop)
74 FILE *gp;
75 int indent, i, j, **gcount, atot;
76 t_state state;
77 t_inputrec ir;
78 t_tpxheader tpx;
79 gmx_mtop_t mtop;
80 gmx_groups_t *groups;
81 t_topology top;
83 read_tpxheader(fn, &tpx, TRUE, NULL, NULL);
85 read_tpx_state(fn,
86 tpx.bIr ? &ir : NULL,
87 &state,
88 tpx.bTop ? &mtop : NULL);
90 if (mdpfn && tpx.bIr)
92 gp = gmx_fio_fopen(mdpfn, "w");
93 pr_inputrec(gp, 0, NULL, &(ir), TRUE);
94 gmx_fio_fclose(gp);
97 if (!mdpfn)
99 if (bSysTop)
101 top = gmx_mtop_t_to_t_topology(&mtop);
104 if (available(stdout, &tpx, 0, fn))
106 indent = 0;
107 pr_title(stdout, indent, fn);
108 pr_inputrec(stdout, 0, "inputrec", tpx.bIr ? &(ir) : NULL, FALSE);
110 pr_tpxheader(stdout, indent, "header", &(tpx));
112 if (!bSysTop)
114 pr_mtop(stdout, indent, "topology", &(mtop), bShowNumbers);
116 else
118 pr_top(stdout, indent, "topology", &(top), bShowNumbers);
121 pr_rvecs(stdout, indent, "box", tpx.bBox ? state.box : NULL, DIM);
122 pr_rvecs(stdout, indent, "box_rel", tpx.bBox ? state.box_rel : NULL, DIM);
123 pr_rvecs(stdout, indent, "boxv", tpx.bBox ? state.boxv : NULL, DIM);
124 pr_rvecs(stdout, indent, "pres_prev", tpx.bBox ? state.pres_prev : NULL, DIM);
125 pr_rvecs(stdout, indent, "svir_prev", tpx.bBox ? state.svir_prev : NULL, DIM);
126 pr_rvecs(stdout, indent, "fvir_prev", tpx.bBox ? state.fvir_prev : NULL, DIM);
127 /* leave nosehoover_xi in for now to match the tpr version */
128 pr_doubles(stdout, indent, "nosehoover_xi", state.nosehoover_xi, state.ngtc);
129 /*pr_doubles(stdout,indent,"nosehoover_vxi",state.nosehoover_vxi,state.ngtc);*/
130 /*pr_doubles(stdout,indent,"therm_integral",state.therm_integral,state.ngtc);*/
131 pr_rvecs(stdout, indent, "x", tpx.bX ? state.x : NULL, state.natoms);
132 pr_rvecs(stdout, indent, "v", tpx.bV ? state.v : NULL, state.natoms);
135 groups = &mtop.groups;
137 snew(gcount, egcNR);
138 for (i = 0; (i < egcNR); i++)
140 snew(gcount[i], groups->grps[i].nr);
143 for (i = 0; (i < mtop.natoms); i++)
145 for (j = 0; (j < egcNR); j++)
147 gcount[j][ggrpnr(groups, j, i)]++;
150 printf("Group statistics\n");
151 for (i = 0; (i < egcNR); i++)
153 atot = 0;
154 printf("%-12s: ", gtypes[i]);
155 for (j = 0; (j < groups->grps[i].nr); j++)
157 printf(" %5d", gcount[i][j]);
158 atot += gcount[i][j];
160 printf(" (total %d atoms)\n", atot);
161 sfree(gcount[i]);
163 sfree(gcount);
165 done_state(&state);
168 static void list_top(const char *fn)
170 int status, done;
171 #define BUFLEN 256
172 char buf[BUFLEN];
173 gmx_cpp_t handle;
174 char *cppopts[] = { NULL };
176 status = cpp_open_file(fn, &handle, cppopts);
177 if (status != 0)
179 gmx_fatal(FARGS, cpp_error(&handle, status));
183 status = cpp_read_line(&handle, BUFLEN, buf);
184 done = (status == eCPP_EOF);
185 if (!done)
187 if (status != eCPP_OK)
189 gmx_fatal(FARGS, cpp_error(&handle, status));
191 else
193 printf("%s\n", buf);
197 while (!done);
198 status = cpp_close_file(&handle);
199 if (status != eCPP_OK)
201 gmx_fatal(FARGS, cpp_error(&handle, status));
205 static void list_trr(const char *fn)
207 t_fileio *fpread;
208 int nframe, indent;
209 char buf[256];
210 rvec *x, *v, *f;
211 matrix box;
212 gmx_trr_header_t trrheader;
213 gmx_bool bOK;
215 fpread = gmx_trr_open(fn, "r");
217 nframe = 0;
218 while (gmx_trr_read_frame_header(fpread, &trrheader, &bOK))
220 snew(x, trrheader.natoms);
221 snew(v, trrheader.natoms);
222 snew(f, trrheader.natoms);
223 if (gmx_trr_read_frame_data(fpread, &trrheader,
224 trrheader.box_size ? box : NULL,
225 trrheader.x_size ? x : NULL,
226 trrheader.v_size ? v : NULL,
227 trrheader.f_size ? f : NULL))
229 sprintf(buf, "%s frame %d", fn, nframe);
230 indent = 0;
231 indent = pr_title(stdout, indent, buf);
232 pr_indent(stdout, indent);
233 fprintf(stdout, "natoms=%10d step=%10d time=%12.7e lambda=%10g\n",
234 trrheader.natoms, trrheader.step, trrheader.t, trrheader.lambda);
235 if (trrheader.box_size)
237 pr_rvecs(stdout, indent, "box", box, DIM);
239 if (trrheader.x_size)
241 pr_rvecs(stdout, indent, "x", x, trrheader.natoms);
243 if (trrheader.v_size)
245 pr_rvecs(stdout, indent, "v", v, trrheader.natoms);
247 if (trrheader.f_size)
249 pr_rvecs(stdout, indent, "f", f, trrheader.natoms);
252 else
254 fprintf(stderr, "\nWARNING: Incomplete frame: nr %d, t=%g\n",
255 nframe, trrheader.t);
258 sfree(x);
259 sfree(v);
260 sfree(f);
261 nframe++;
263 if (!bOK)
265 fprintf(stderr, "\nWARNING: Incomplete frame header: nr %d, t=%g\n",
266 nframe, trrheader.t);
268 gmx_trr_close(fpread);
271 void list_xtc(const char *fn)
273 t_fileio *xd;
274 int indent;
275 char buf[256];
276 rvec *x;
277 matrix box;
278 int nframe, natoms, step;
279 real prec, time;
280 gmx_bool bOK;
282 xd = open_xtc(fn, "r");
283 read_first_xtc(xd, &natoms, &step, &time, box, &x, &prec, &bOK);
285 nframe = 0;
288 sprintf(buf, "%s frame %d", fn, nframe);
289 indent = 0;
290 indent = pr_title(stdout, indent, buf);
291 pr_indent(stdout, indent);
292 fprintf(stdout, "natoms=%10d step=%10d time=%12.7e prec=%10g\n",
293 natoms, step, time, prec);
294 pr_rvecs(stdout, indent, "box", box, DIM);
295 pr_rvecs(stdout, indent, "x", x, natoms);
296 nframe++;
298 while (read_next_xtc(xd, natoms, &step, &time, box, x, &prec, &bOK));
299 if (!bOK)
301 fprintf(stderr, "\nWARNING: Incomplete frame at time %g\n", time);
303 sfree(x);
304 close_xtc(xd);
307 /*! \brief Callback used by list_tng_for_gmx_dump. */
308 static void list_tng_inner(const char *fn,
309 gmx_bool bFirstFrame,
310 real *values,
311 gmx_int64_t step,
312 double frame_time,
313 gmx_int64_t n_values_per_frame,
314 gmx_int64_t n_atoms,
315 real prec,
316 gmx_int64_t nframe,
317 char *block_name)
319 char buf[256];
320 int indent = 0;
322 if (bFirstFrame)
324 sprintf(buf, "%s frame %" GMX_PRId64, fn, nframe);
325 indent = 0;
326 indent = pr_title(stdout, indent, buf);
327 pr_indent(stdout, indent);
328 fprintf(stdout, "natoms=%10" GMX_PRId64 " step=%10" GMX_PRId64 " time=%12.7e",
329 n_atoms, step, frame_time);
330 if (prec > 0)
332 fprintf(stdout, " prec=%10g", prec);
334 fprintf(stdout, "\n");
336 pr_reals_of_dim(stdout, indent, block_name, values, n_atoms, n_values_per_frame);
339 static void list_tng(const char gmx_unused *fn)
341 #ifdef GMX_USE_TNG
342 tng_trajectory_t tng;
343 gmx_int64_t nframe = 0;
344 gmx_int64_t i, *block_ids = NULL, step, ndatablocks;
345 gmx_bool bOK;
347 gmx_tng_open(fn, 'r', &tng);
348 gmx_print_tng_molecule_system(tng, stdout);
350 bOK = gmx_get_tng_data_block_types_of_next_frame(tng, -1,
352 NULL,
353 &step, &ndatablocks,
354 &block_ids);
357 for (i = 0; i < ndatablocks; i++)
359 double frame_time;
360 real prec, *values = NULL;
361 gmx_int64_t n_values_per_frame, n_atoms;
362 char block_name[STRLEN];
364 gmx_get_tng_data_next_frame_of_block_type(tng, block_ids[i], &values,
365 &step, &frame_time,
366 &n_values_per_frame, &n_atoms,
367 &prec,
368 block_name, STRLEN, &bOK);
369 if (!bOK)
371 /* Can't write any output because we don't know what
372 arrays are valid. */
373 fprintf(stderr, "\nWARNING: Incomplete frame at time %g, will not write output\n", frame_time);
375 else
377 list_tng_inner(fn, (0 == i), values, step, frame_time,
378 n_values_per_frame, n_atoms, prec, nframe, block_name);
381 nframe++;
383 while (gmx_get_tng_data_block_types_of_next_frame(tng, step,
385 NULL,
386 &step,
387 &ndatablocks,
388 &block_ids));
390 if (block_ids)
392 sfree(block_ids);
395 gmx_tng_close(&tng);
396 #endif
399 void list_trx(const char *fn)
401 switch (fn2ftp(fn))
403 case efXTC:
404 list_xtc(fn);
405 break;
406 case efTRR:
407 list_trr(fn);
408 break;
409 case efTNG:
410 list_tng(fn);
411 break;
412 default:
413 fprintf(stderr, "File %s is of an unsupported type. Try using the command\n 'less %s'\n",
414 fn, fn);
418 void list_ene(const char *fn)
420 ener_file_t in;
421 gmx_bool bCont;
422 gmx_enxnm_t *enm = NULL;
423 t_enxframe *fr;
424 int i, j, nre, b;
425 char buf[22];
427 printf("gmx dump: %s\n", fn);
428 in = open_enx(fn, "r");
429 do_enxnms(in, &nre, &enm);
430 assert(enm);
432 printf("energy components:\n");
433 for (i = 0; (i < nre); i++)
435 printf("%5d %-24s (%s)\n", i, enm[i].name, enm[i].unit);
438 snew(fr, 1);
441 bCont = do_enx(in, fr);
443 if (bCont)
445 printf("\n%24s %12.5e %12s %12s\n", "time:",
446 fr->t, "step:", gmx_step_str(fr->step, buf));
447 printf("%24s %12s %12s %12s\n",
448 "", "", "nsteps:", gmx_step_str(fr->nsteps, buf));
449 printf("%24s %12.5e %12s %12s\n",
450 "delta_t:", fr->dt, "sum steps:", gmx_step_str(fr->nsum, buf));
451 if (fr->nre == nre)
453 printf("%24s %12s %12s %12s\n",
454 "Component", "Energy", "Av. Energy", "Sum Energy");
455 if (fr->nsum > 0)
457 for (i = 0; (i < nre); i++)
459 printf("%24s %12.5e %12.5e %12.5e\n",
460 enm[i].name, fr->ener[i].e, fr->ener[i].eav,
461 fr->ener[i].esum);
464 else
466 for (i = 0; (i < nre); i++)
468 printf("%24s %12.5e\n",
469 enm[i].name, fr->ener[i].e);
473 for (b = 0; b < fr->nblock; b++)
475 const char *typestr = "";
477 t_enxblock *eb = &(fr->block[b]);
478 printf("Block data %2d (%3d subblocks, id=%d)\n",
479 b, eb->nsub, eb->id);
481 if (eb->id < enxNR)
483 typestr = enx_block_id_name[eb->id];
485 printf(" id='%s'\n", typestr);
486 for (i = 0; i < eb->nsub; i++)
488 t_enxsubblock *sb = &(eb->sub[i]);
489 printf(" Sub block %3d (%5d elems, type=%s) values:\n",
490 i, sb->nr, xdr_datatype_names[sb->type]);
492 switch (sb->type)
494 case xdr_datatype_float:
495 for (j = 0; j < sb->nr; j++)
497 printf("%14d %8.4f\n", j, sb->fval[j]);
499 break;
500 case xdr_datatype_double:
501 for (j = 0; j < sb->nr; j++)
503 printf("%14d %10.6f\n", j, sb->dval[j]);
505 break;
506 case xdr_datatype_int:
507 for (j = 0; j < sb->nr; j++)
509 printf("%14d %10d\n", j, sb->ival[j]);
511 break;
512 case xdr_datatype_int64:
513 for (j = 0; j < sb->nr; j++)
515 printf("%14d %s\n",
516 j, gmx_step_str(sb->lval[j], buf));
518 break;
519 case xdr_datatype_char:
520 for (j = 0; j < sb->nr; j++)
522 printf("%14d %1c\n", j, sb->cval[j]);
524 break;
525 case xdr_datatype_string:
526 for (j = 0; j < sb->nr; j++)
528 printf("%14d %80s\n", j, sb->sval[j]);
530 break;
531 default:
532 gmx_incons("Unknown subblock type");
538 while (bCont);
540 close_enx(in);
542 free_enxframe(fr);
543 sfree(fr);
544 sfree(enm);
547 static void list_mtx(const char *fn)
549 int nrow, ncol, i, j, k;
550 real *full = NULL, value;
551 gmx_sparsematrix_t * sparse = NULL;
553 gmx_mtxio_read(fn, &nrow, &ncol, &full, &sparse);
555 if (full == NULL)
557 snew(full, nrow*ncol);
558 for (i = 0; i < nrow*ncol; i++)
560 full[i] = 0;
563 for (i = 0; i < sparse->nrow; i++)
565 for (j = 0; j < sparse->ndata[i]; j++)
567 k = sparse->data[i][j].col;
568 value = sparse->data[i][j].value;
569 full[i*ncol+k] = value;
570 full[k*ncol+i] = value;
573 gmx_sparsematrix_destroy(sparse);
576 printf("%d %d\n", nrow, ncol);
577 for (i = 0; i < nrow; i++)
579 for (j = 0; j < ncol; j++)
581 printf(" %g", full[i*ncol+j]);
583 printf("\n");
586 sfree(full);
589 int gmx_dump(int argc, char *argv[])
591 const char *desc[] = {
592 "[THISMODULE] reads a run input file ([REF].tpr[ref]),",
593 "a trajectory ([REF].trr[ref]/[REF].xtc[ref]/[TT]/tng[tt]), an energy",
594 "file ([REF].edr[ref]) or a checkpoint file ([REF].cpt[ref])",
595 "and prints that to standard output in a readable format.",
596 "This program is essential for checking your run input file in case of",
597 "problems.[PAR]",
598 "The program can also preprocess a topology to help finding problems.",
599 "Note that currently setting [TT]GMXLIB[tt] is the only way to customize",
600 "directories used for searching include files.",
602 const char *bugs[] = {
603 "Position restraint output from -sys -s is broken"
605 t_filenm fnm[] = {
606 { efTPR, "-s", NULL, ffOPTRD },
607 { efTRX, "-f", NULL, ffOPTRD },
608 { efEDR, "-e", NULL, ffOPTRD },
609 { efCPT, NULL, NULL, ffOPTRD },
610 { efTOP, "-p", NULL, ffOPTRD },
611 { efMTX, "-mtx", "hessian", ffOPTRD },
612 { efMDP, "-om", NULL, ffOPTWR }
614 #define NFILE asize(fnm)
616 gmx_output_env_t *oenv;
617 /* Command line options */
618 static gmx_bool bShowNumbers = TRUE;
619 static gmx_bool bSysTop = FALSE;
620 t_pargs pa[] = {
621 { "-nr", FALSE, etBOOL, {&bShowNumbers}, "Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" },
622 { "-sys", FALSE, etBOOL, {&bSysTop}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" }
625 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
626 asize(desc), desc, asize(bugs), bugs, &oenv))
628 return 0;
632 if (ftp2bSet(efTPR, NFILE, fnm))
634 list_tpx(ftp2fn(efTPR, NFILE, fnm), bShowNumbers,
635 ftp2fn_null(efMDP, NFILE, fnm), bSysTop);
637 else if (ftp2bSet(efTRX, NFILE, fnm))
639 list_trx(ftp2fn(efTRX, NFILE, fnm));
641 else if (ftp2bSet(efEDR, NFILE, fnm))
643 list_ene(ftp2fn(efEDR, NFILE, fnm));
645 else if (ftp2bSet(efCPT, NFILE, fnm))
647 list_checkpoint(ftp2fn(efCPT, NFILE, fnm), stdout);
649 else if (ftp2bSet(efTOP, NFILE, fnm))
651 list_top(ftp2fn(efTOP, NFILE, fnm));
653 else if (ftp2bSet(efMTX, NFILE, fnm))
655 list_mtx(ftp2fn(efMTX, NFILE, fnm));
658 return 0;