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