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