Merge branch release-5-1 into release-2016
[alexxy/gromacs.git] / src / gromacs / tools / dump.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "config.h"
40
41 #include <cassert>
42 #include <cmath>
43 #include <cstdio>
44 #include <cstring>
45
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/inputrec.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/trajectory/trajectoryframe.h"
66 #include "gromacs/utility/arraysize.h"
67 #include "gromacs/utility/basedefinitions.h"
68 #include "gromacs/utility/fatalerror.h"
69 #include "gromacs/utility/futil.h"
70 #include "gromacs/utility/smalloc.h"
71 #include "gromacs/utility/txtdump.h"
72
73 static void list_tpx(const char *fn, gmx_bool bShowNumbers, const char *mdpfn,
74                      gmx_bool bSysTop)
75 {
76     FILE         *gp;
77     int           indent, i, j, **gcount, atot;
78     t_state       state;
79     t_inputrec    ir;
80     t_tpxheader   tpx;
81     gmx_mtop_t    mtop;
82     gmx_groups_t *groups;
83     t_topology    top;
84
85     read_tpxheader(fn, &tpx, TRUE);
86
87     read_tpx_state(fn,
88                    tpx.bIr  ? &ir : NULL,
89                    &state,
90                    tpx.bTop ? &mtop : NULL);
91
92     if (mdpfn && tpx.bIr)
93     {
94         gp = gmx_fio_fopen(mdpfn, "w");
95         pr_inputrec(gp, 0, NULL, &(ir), TRUE);
96         gmx_fio_fclose(gp);
97     }
98
99     if (!mdpfn)
100     {
101         if (bSysTop)
102         {
103             top = gmx_mtop_t_to_t_topology(&mtop);
104         }
105
106         if (available(stdout, &tpx, 0, fn))
107         {
108             indent = 0;
109             pr_title(stdout, indent, fn);
110             pr_inputrec(stdout, 0, "inputrec", tpx.bIr ? &(ir) : NULL, FALSE);
111
112             pr_tpxheader(stdout, indent, "header", &(tpx));
113
114             if (!bSysTop)
115             {
116                 pr_mtop(stdout, indent, "topology", &(mtop), bShowNumbers);
117             }
118             else
119             {
120                 pr_top(stdout, indent, "topology", &(top), bShowNumbers);
121             }
122
123             pr_rvecs(stdout, indent, "box", tpx.bBox ? state.box : NULL, DIM);
124             pr_rvecs(stdout, indent, "box_rel", tpx.bBox ? state.box_rel : NULL, DIM);
125             pr_rvecs(stdout, indent, "boxv", tpx.bBox ? state.boxv : NULL, DIM);
126             pr_rvecs(stdout, indent, "pres_prev", tpx.bBox ? state.pres_prev : NULL, DIM);
127             pr_rvecs(stdout, indent, "svir_prev", tpx.bBox ? state.svir_prev : NULL, DIM);
128             pr_rvecs(stdout, indent, "fvir_prev", tpx.bBox ? state.fvir_prev : NULL, DIM);
129             /* leave nosehoover_xi in for now to match the tpr version */
130             pr_doubles(stdout, indent, "nosehoover_xi", state.nosehoover_xi, state.ngtc);
131             /*pr_doubles(stdout,indent,"nosehoover_vxi",state.nosehoover_vxi,state.ngtc);*/
132             /*pr_doubles(stdout,indent,"therm_integral",state.therm_integral,state.ngtc);*/
133             pr_rvecs(stdout, indent, "x", tpx.bX ? state.x : NULL, state.natoms);
134             pr_rvecs(stdout, indent, "v", tpx.bV ? state.v : NULL, state.natoms);
135         }
136
137         groups = &mtop.groups;
138
139         snew(gcount, egcNR);
140         for (i = 0; (i < egcNR); i++)
141         {
142             snew(gcount[i], groups->grps[i].nr);
143         }
144
145         for (i = 0; (i < mtop.natoms); i++)
146         {
147             for (j = 0; (j < egcNR); j++)
148             {
149                 gcount[j][ggrpnr(groups, j, i)]++;
150             }
151         }
152         printf("Group statistics\n");
153         for (i = 0; (i < egcNR); i++)
154         {
155             atot = 0;
156             printf("%-12s: ", gtypes[i]);
157             for (j = 0; (j < groups->grps[i].nr); j++)
158             {
159                 printf("  %5d", gcount[i][j]);
160                 atot += gcount[i][j];
161             }
162             printf("  (total %d atoms)\n", atot);
163             sfree(gcount[i]);
164         }
165         sfree(gcount);
166     }
167     done_state(&state);
168 }
169
170 static void list_top(const char *fn)
171 {
172     int       status, done;
173 #define BUFLEN 256
174     char      buf[BUFLEN];
175     gmx_cpp_t handle;
176     char     *cppopts[] = { NULL };
177
178     status = cpp_open_file(fn, &handle, cppopts);
179     if (status != 0)
180     {
181         gmx_fatal(FARGS, cpp_error(&handle, status));
182     }
183     do
184     {
185         status = cpp_read_line(&handle, BUFLEN, buf);
186         done   = (status == eCPP_EOF);
187         if (!done)
188         {
189             if (status != eCPP_OK)
190             {
191                 gmx_fatal(FARGS, cpp_error(&handle, status));
192             }
193             else
194             {
195                 printf("%s\n", buf);
196             }
197         }
198     }
199     while (!done);
200     status = cpp_close_file(&handle);
201     if (status != eCPP_OK)
202     {
203         gmx_fatal(FARGS, cpp_error(&handle, status));
204     }
205 }
206
207 static void list_trr(const char *fn)
208 {
209     t_fileio         *fpread;
210     int               nframe, indent;
211     char              buf[256];
212     rvec             *x, *v, *f;
213     matrix            box;
214     gmx_trr_header_t  trrheader;
215     gmx_bool          bOK;
216
217     fpread  = gmx_trr_open(fn, "r");
218
219     nframe = 0;
220     while (gmx_trr_read_frame_header(fpread, &trrheader, &bOK))
221     {
222         snew(x, trrheader.natoms);
223         snew(v, trrheader.natoms);
224         snew(f, trrheader.natoms);
225         if (gmx_trr_read_frame_data(fpread, &trrheader,
226                                     trrheader.box_size ? box : NULL,
227                                     trrheader.x_size   ? x : NULL,
228                                     trrheader.v_size   ? v : NULL,
229                                     trrheader.f_size   ? f : NULL))
230         {
231             sprintf(buf, "%s frame %d", fn, nframe);
232             indent = 0;
233             indent = pr_title(stdout, indent, buf);
234             pr_indent(stdout, indent);
235             fprintf(stdout, "natoms=%10d  step=%10" GMX_PRId64 "  time=%12.7e  lambda=%10g\n",
236                     trrheader.natoms, trrheader.step, trrheader.t, trrheader.lambda);
237             if (trrheader.box_size)
238             {
239                 pr_rvecs(stdout, indent, "box", box, DIM);
240             }
241             if (trrheader.x_size)
242             {
243                 pr_rvecs(stdout, indent, "x", x, trrheader.natoms);
244             }
245             if (trrheader.v_size)
246             {
247                 pr_rvecs(stdout, indent, "v", v, trrheader.natoms);
248             }
249             if (trrheader.f_size)
250             {
251                 pr_rvecs(stdout, indent, "f", f, trrheader.natoms);
252             }
253         }
254         else
255         {
256             fprintf(stderr, "\nWARNING: Incomplete frame: nr %d, t=%g\n",
257                     nframe, trrheader.t);
258         }
259
260         sfree(x);
261         sfree(v);
262         sfree(f);
263         nframe++;
264     }
265     if (!bOK)
266     {
267         fprintf(stderr, "\nWARNING: Incomplete frame header: nr %d, t=%g\n",
268                 nframe, trrheader.t);
269     }
270     gmx_trr_close(fpread);
271 }
272
273 void list_xtc(const char *fn)
274 {
275     t_fileio   *xd;
276     int         indent;
277     char        buf[256];
278     rvec       *x;
279     matrix      box;
280     int         nframe, natoms;
281     gmx_int64_t step;
282     real        prec, time;
283     gmx_bool    bOK;
284
285     xd = open_xtc(fn, "r");
286     read_first_xtc(xd, &natoms, &step, &time, box, &x, &prec, &bOK);
287
288     nframe = 0;
289     do
290     {
291         sprintf(buf, "%s frame %d", fn, nframe);
292         indent = 0;
293         indent = pr_title(stdout, indent, buf);
294         pr_indent(stdout, indent);
295         fprintf(stdout, "natoms=%10d  step=%10" GMX_PRId64 "  time=%12.7e  prec=%10g\n",
296                 natoms, step, time, prec);
297         pr_rvecs(stdout, indent, "box", box, DIM);
298         pr_rvecs(stdout, indent, "x", x, natoms);
299         nframe++;
300     }
301     while (read_next_xtc(xd, natoms, &step, &time, box, x, &prec, &bOK));
302     if (!bOK)
303     {
304         fprintf(stderr, "\nWARNING: Incomplete frame at time %g\n", time);
305     }
306     sfree(x);
307     close_xtc(xd);
308 }
309
310 /*! \brief Callback used by list_tng_for_gmx_dump. */
311 static void list_tng_inner(const char *fn,
312                            gmx_bool    bFirstFrame,
313                            real       *values,
314                            gmx_int64_t step,
315                            double      frame_time,
316                            gmx_int64_t n_values_per_frame,
317                            gmx_int64_t n_atoms,
318                            real        prec,
319                            gmx_int64_t nframe,
320                            char       *block_name)
321 {
322     char                 buf[256];
323     int                  indent = 0;
324
325     if (bFirstFrame)
326     {
327         sprintf(buf, "%s frame %" GMX_PRId64, fn, nframe);
328         indent = 0;
329         indent = pr_title(stdout, indent, buf);
330         pr_indent(stdout, indent);
331         fprintf(stdout, "natoms=%10" GMX_PRId64 "  step=%10" GMX_PRId64 "  time=%12.7e",
332                 n_atoms, step, frame_time);
333         if (prec > 0)
334         {
335             fprintf(stdout, "  prec=%10g", prec);
336         }
337         fprintf(stdout, "\n");
338     }
339     pr_reals_of_dim(stdout, indent, block_name, values, n_atoms, n_values_per_frame);
340 }
341
342 static void list_tng(const char gmx_unused *fn)
343 {
344 #ifdef GMX_USE_TNG
345     tng_trajectory_t     tng;
346     gmx_int64_t          nframe = 0;
347     gmx_int64_t          i, *block_ids = NULL, step, ndatablocks;
348     gmx_bool             bOK;
349     real                *values = NULL;
350
351     gmx_tng_open(fn, 'r', &tng);
352     gmx_print_tng_molecule_system(tng, stdout);
353
354     bOK    = gmx_get_tng_data_block_types_of_next_frame(tng, -1,
355                                                         0,
356                                                         NULL,
357                                                         &step, &ndatablocks,
358                                                         &block_ids);
359     do
360     {
361         for (i = 0; i < ndatablocks; i++)
362         {
363             double               frame_time;
364             real                 prec;
365             gmx_int64_t          n_values_per_frame, n_atoms;
366             char                 block_name[STRLEN];
367
368             gmx_get_tng_data_next_frame_of_block_type(tng, block_ids[i], &values,
369                                                       &step, &frame_time,
370                                                       &n_values_per_frame, &n_atoms,
371                                                       &prec,
372                                                       block_name, STRLEN, &bOK);
373             if (!bOK)
374             {
375                 /* Can't write any output because we don't know what
376                    arrays are valid. */
377                 fprintf(stderr, "\nWARNING: Incomplete frame at time %g, will not write output\n", frame_time);
378             }
379             else
380             {
381                 list_tng_inner(fn, (0 == i), values, step, frame_time,
382                                n_values_per_frame, n_atoms, prec, nframe, block_name);
383             }
384         }
385         nframe++;
386     }
387     while (gmx_get_tng_data_block_types_of_next_frame(tng, step,
388                                                       0,
389                                                       NULL,
390                                                       &step,
391                                                       &ndatablocks,
392                                                       &block_ids));
393
394     if (block_ids)
395     {
396         sfree(block_ids);
397     }
398     sfree(values);
399     gmx_tng_close(&tng);
400 #endif
401 }
402
403 void list_trx(const char *fn)
404 {
405     switch (fn2ftp(fn))
406     {
407         case efXTC:
408             list_xtc(fn);
409             break;
410         case efTRR:
411             list_trr(fn);
412             break;
413         case efTNG:
414             list_tng(fn);
415             break;
416         default:
417             fprintf(stderr, "File %s is of an unsupported type. Try using the command\n 'less %s'\n",
418                     fn, fn);
419     }
420 }
421
422 void list_ene(const char *fn)
423 {
424     ener_file_t    in;
425     gmx_bool       bCont;
426     gmx_enxnm_t   *enm = NULL;
427     t_enxframe    *fr;
428     int            i, j, nre, b;
429     char           buf[22];
430
431     printf("gmx dump: %s\n", fn);
432     in = open_enx(fn, "r");
433     do_enxnms(in, &nre, &enm);
434     assert(enm);
435
436     printf("energy components:\n");
437     for (i = 0; (i < nre); i++)
438     {
439         printf("%5d  %-24s (%s)\n", i, enm[i].name, enm[i].unit);
440     }
441
442     snew(fr, 1);
443     do
444     {
445         bCont = do_enx(in, fr);
446
447         if (bCont)
448         {
449             printf("\n%24s  %12.5e  %12s  %12s\n", "time:",
450                    fr->t, "step:", gmx_step_str(fr->step, buf));
451             printf("%24s  %12s  %12s  %12s\n",
452                    "", "", "nsteps:", gmx_step_str(fr->nsteps, buf));
453             printf("%24s  %12.5e  %12s  %12s\n",
454                    "delta_t:", fr->dt, "sum steps:", gmx_step_str(fr->nsum, buf));
455             if (fr->nre == nre)
456             {
457                 printf("%24s  %12s  %12s  %12s\n",
458                        "Component", "Energy", "Av. Energy", "Sum Energy");
459                 if (fr->nsum > 0)
460                 {
461                     for (i = 0; (i < nre); i++)
462                     {
463                         printf("%24s  %12.5e  %12.5e  %12.5e\n",
464                                enm[i].name, fr->ener[i].e, fr->ener[i].eav,
465                                fr->ener[i].esum);
466                     }
467                 }
468                 else
469                 {
470                     for (i = 0; (i < nre); i++)
471                     {
472                         printf("%24s  %12.5e\n",
473                                enm[i].name, fr->ener[i].e);
474                     }
475                 }
476             }
477             for (b = 0; b < fr->nblock; b++)
478             {
479                 const char *typestr = "";
480
481                 t_enxblock *eb = &(fr->block[b]);
482                 printf("Block data %2d (%3d subblocks, id=%d)\n",
483                        b, eb->nsub, eb->id);
484
485                 if (eb->id < enxNR)
486                 {
487                     typestr = enx_block_id_name[eb->id];
488                 }
489                 printf("  id='%s'\n", typestr);
490                 for (i = 0; i < eb->nsub; i++)
491                 {
492                     t_enxsubblock *sb = &(eb->sub[i]);
493                     printf("  Sub block %3d (%5d elems, type=%s) values:\n",
494                            i, sb->nr, xdr_datatype_names[sb->type]);
495
496                     switch (sb->type)
497                     {
498                         case xdr_datatype_float:
499                             for (j = 0; j < sb->nr; j++)
500                             {
501                                 printf("%14d   %8.4f\n", j, sb->fval[j]);
502                             }
503                             break;
504                         case xdr_datatype_double:
505                             for (j = 0; j < sb->nr; j++)
506                             {
507                                 printf("%14d   %10.6f\n", j, sb->dval[j]);
508                             }
509                             break;
510                         case xdr_datatype_int:
511                             for (j = 0; j < sb->nr; j++)
512                             {
513                                 printf("%14d %10d\n", j, sb->ival[j]);
514                             }
515                             break;
516                         case xdr_datatype_int64:
517                             for (j = 0; j < sb->nr; j++)
518                             {
519                                 printf("%14d %s\n",
520                                        j, gmx_step_str(sb->lval[j], buf));
521                             }
522                             break;
523                         case xdr_datatype_char:
524                             for (j = 0; j < sb->nr; j++)
525                             {
526                                 printf("%14d %1c\n", j, sb->cval[j]);
527                             }
528                             break;
529                         case xdr_datatype_string:
530                             for (j = 0; j < sb->nr; j++)
531                             {
532                                 printf("%14d %80s\n", j, sb->sval[j]);
533                             }
534                             break;
535                         default:
536                             gmx_incons("Unknown subblock type");
537                     }
538                 }
539             }
540         }
541     }
542     while (bCont);
543
544     close_enx(in);
545
546     free_enxframe(fr);
547     sfree(fr);
548     sfree(enm);
549 }
550
551 static void list_mtx(const char *fn)
552 {
553     int                  nrow, ncol, i, j, k;
554     real                *full   = NULL, value;
555     gmx_sparsematrix_t * sparse = NULL;
556
557     gmx_mtxio_read(fn, &nrow, &ncol, &full, &sparse);
558
559     if (full == NULL)
560     {
561         snew(full, nrow*ncol);
562         for (i = 0; i < nrow*ncol; i++)
563         {
564             full[i] = 0;
565         }
566
567         for (i = 0; i < sparse->nrow; i++)
568         {
569             for (j = 0; j < sparse->ndata[i]; j++)
570             {
571                 k              = sparse->data[i][j].col;
572                 value          = sparse->data[i][j].value;
573                 full[i*ncol+k] = value;
574                 full[k*ncol+i] = value;
575             }
576         }
577         gmx_sparsematrix_destroy(sparse);
578     }
579
580     printf("%d %d\n", nrow, ncol);
581     for (i = 0; i < nrow; i++)
582     {
583         for (j = 0; j < ncol; j++)
584         {
585             printf(" %g", full[i*ncol+j]);
586         }
587         printf("\n");
588     }
589
590     sfree(full);
591 }
592
593 int gmx_dump(int argc, char *argv[])
594 {
595     const char *desc[] = {
596         "[THISMODULE] reads a run input file ([REF].tpr[ref]),",
597         "a trajectory ([REF].trr[ref]/[REF].xtc[ref]/[TT]/tng[tt]), an energy",
598         "file ([REF].edr[ref]) or a checkpoint file ([REF].cpt[ref])",
599         "and prints that to standard output in a readable format.",
600         "This program is essential for checking your run input file in case of",
601         "problems.[PAR]",
602         "The program can also preprocess a topology to help finding problems.",
603         "Note that currently setting [TT]GMXLIB[tt] is the only way to customize",
604         "directories used for searching include files.",
605     };
606     const char *bugs[] = {
607         "Position restraint output from -sys -s is broken"
608     };
609     t_filenm    fnm[] = {
610         { efTPR, "-s", NULL, ffOPTRD },
611         { efTRX, "-f", NULL, ffOPTRD },
612         { efEDR, "-e", NULL, ffOPTRD },
613         { efCPT, NULL, NULL, ffOPTRD },
614         { efTOP, "-p", NULL, ffOPTRD },
615         { efMTX, "-mtx", "hessian", ffOPTRD },
616         { efMDP, "-om", NULL, ffOPTWR }
617     };
618 #define NFILE asize(fnm)
619
620     gmx_output_env_t *oenv;
621     /* Command line options */
622     static gmx_bool   bShowNumbers = TRUE;
623     static gmx_bool   bSysTop      = FALSE;
624     t_pargs           pa[]         = {
625         { "-nr", FALSE, etBOOL, {&bShowNumbers}, "Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" },
626         { "-sys", FALSE, etBOOL, {&bSysTop}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" }
627     };
628
629     if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
630                            asize(desc), desc, asize(bugs), bugs, &oenv))
631     {
632         return 0;
633     }
634
635
636     if (ftp2bSet(efTPR, NFILE, fnm))
637     {
638         list_tpx(ftp2fn(efTPR, NFILE, fnm), bShowNumbers,
639                  ftp2fn_null(efMDP, NFILE, fnm), bSysTop);
640     }
641     else if (ftp2bSet(efTRX, NFILE, fnm))
642     {
643         list_trx(ftp2fn(efTRX, NFILE, fnm));
644     }
645     else if (ftp2bSet(efEDR, NFILE, fnm))
646     {
647         list_ene(ftp2fn(efEDR, NFILE, fnm));
648     }
649     else if (ftp2bSet(efCPT, NFILE, fnm))
650     {
651         list_checkpoint(ftp2fn(efCPT, NFILE, fnm), stdout);
652     }
653     else if (ftp2bSet(efTOP, NFILE, fnm))
654     {
655         list_top(ftp2fn(efTOP, NFILE, fnm));
656     }
657     else if (ftp2bSet(efMTX, NFILE, fnm))
658     {
659         list_mtx(ftp2fn(efMTX, NFILE, fnm));
660     }
661
662     return 0;
663 }