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