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