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