More gmxfio simplification
[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;
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
219     nframe = 0;
220     while (fread_trnheader(fpread, &trn, &bOK))
221     {
222         snew(x, trn.natoms);
223         snew(v, trn.natoms);
224         snew(f, trn.natoms);
225         if (fread_htrn(fpread, &trn,
226                        trn.box_size ? box : NULL,
227                        trn.x_size   ? x : NULL,
228                        trn.v_size   ? v : NULL,
229                        trn.f_size   ? f : NULL))
230         {
231             sprintf(buf, "%s frame %d", fn, nframe);
232             indent = 0;
233             indent = pr_title(stdout, indent, buf);
234             pr_indent(stdout, indent);
235             fprintf(stdout, "natoms=%10d  step=%10d  time=%12.7e  lambda=%10g\n",
236                     trn.natoms, trn.step, trn.t, trn.lambda);
237             if (trn.box_size)
238             {
239                 pr_rvecs(stdout, indent, "box", box, DIM);
240             }
241             if (trn.x_size)
242             {
243                 pr_rvecs(stdout, indent, "x", x, trn.natoms);
244             }
245             if (trn.v_size)
246             {
247                 pr_rvecs(stdout, indent, "v", v, trn.natoms);
248             }
249             if (trn.f_size)
250             {
251                 pr_rvecs(stdout, indent, "f", f, trn.natoms);
252             }
253         }
254         else
255         {
256             fprintf(stderr, "\nWARNING: Incomplete frame: nr %d, t=%g\n",
257                     nframe, trn.t);
258         }
259
260         sfree(x);
261         sfree(v);
262         sfree(f);
263         nframe++;
264     }
265     if (!bOK)
266     {
267         fprintf(stderr, "\nWARNING: Incomplete frame header: nr %d, t=%g\n",
268                 nframe, trn.t);
269     }
270     close_trn(fpread);
271 }
272
273 void list_xtc(const char *fn)
274 {
275     t_fileio  *xd;
276     int        indent;
277     char       buf[256];
278     rvec      *x;
279     matrix     box;
280     int        nframe, natoms, step;
281     real       prec, time;
282     gmx_bool   bOK;
283
284     xd = open_xtc(fn, "r");
285     read_first_xtc(xd, &natoms, &step, &time, box, &x, &prec, &bOK);
286
287     nframe = 0;
288     do
289     {
290         sprintf(buf, "%s frame %d", fn, nframe);
291         indent = 0;
292         indent = pr_title(stdout, indent, buf);
293         pr_indent(stdout, indent);
294         fprintf(stdout, "natoms=%10d  step=%10d  time=%12.7e  prec=%10g\n",
295                 natoms, step, time, prec);
296         pr_rvecs(stdout, indent, "box", box, DIM);
297         pr_rvecs(stdout, indent, "x", x, natoms);
298         nframe++;
299     }
300     while (read_next_xtc(xd, natoms, &step, &time, box, x, &prec, &bOK));
301     if (!bOK)
302     {
303         fprintf(stderr, "\nWARNING: Incomplete frame at time %g\n", time);
304     }
305     sfree(x);
306     close_xtc(xd);
307 }
308
309 /*! \brief Callback used by list_tng_for_gmx_dump. */
310 static void list_tng_inner(const char *fn,
311                            gmx_bool    bFirstFrame,
312                            real       *values,
313                            gmx_int64_t step,
314                            double      frame_time,
315                            gmx_int64_t n_values_per_frame,
316                            gmx_int64_t n_atoms,
317                            real        prec,
318                            gmx_int64_t nframe,
319                            char       *block_name)
320 {
321     char                 buf[256];
322     int                  indent = 0;
323
324     if (bFirstFrame)
325     {
326         sprintf(buf, "%s frame %" GMX_PRId64, fn, nframe);
327         indent = 0;
328         indent = pr_title(stdout, indent, buf);
329         pr_indent(stdout, indent);
330         fprintf(stdout, "natoms=%10" GMX_PRId64 "  step=%10" GMX_PRId64 "  time=%12.7e",
331                 n_atoms, step, frame_time);
332         if (prec > 0)
333         {
334             fprintf(stdout, "  prec=%10g", prec);
335         }
336         fprintf(stdout, "\n");
337     }
338     pr_reals_of_dim(stdout, indent, block_name, values, n_atoms, n_values_per_frame);
339 }
340
341 static void list_tng(const char gmx_unused *fn)
342 {
343 #ifdef GMX_USE_TNG
344     tng_trajectory_t     tng;
345     gmx_int64_t          nframe = 0;
346     gmx_int64_t          i, *block_ids = NULL, step, ndatablocks;
347     gmx_bool             bOK;
348
349     gmx_tng_open(fn, 'r', &tng);
350     gmx_print_tng_molecule_system(tng, stdout);
351
352     bOK    = gmx_get_tng_data_block_types_of_next_frame(tng, -1,
353                                                         0,
354                                                         NULL,
355                                                         &step, &ndatablocks,
356                                                         &block_ids);
357     do
358     {
359         for (i = 0; i < ndatablocks; i++)
360         {
361             double               frame_time;
362             real                 prec, *values = NULL;
363             gmx_int64_t          n_values_per_frame, n_atoms;
364             char                 block_name[STRLEN];
365
366             gmx_get_tng_data_next_frame_of_block_type(tng, block_ids[i], &values,
367                                                       &step, &frame_time,
368                                                       &n_values_per_frame, &n_atoms,
369                                                       &prec,
370                                                       block_name, STRLEN, &bOK);
371             if (!bOK)
372             {
373                 /* Can't write any output because we don't know what
374                    arrays are valid. */
375                 fprintf(stderr, "\nWARNING: Incomplete frame at time %g, will not write output\n", frame_time);
376             }
377             else
378             {
379                 list_tng_inner(fn, (0 == i), values, step, frame_time,
380                                n_values_per_frame, n_atoms, prec, nframe, block_name);
381             }
382         }
383         nframe++;
384     }
385     while (gmx_get_tng_data_block_types_of_next_frame(tng, step,
386                                                       0,
387                                                       NULL,
388                                                       &step,
389                                                       &ndatablocks,
390                                                       &block_ids));
391
392     if (block_ids)
393     {
394         sfree(block_ids);
395     }
396
397     gmx_tng_close(&tng);
398 #endif
399 }
400
401 void list_trx(const char *fn)
402 {
403     switch (fn2ftp(fn))
404     {
405         case efXTC:
406             list_xtc(fn);
407             break;
408         case efTRR:
409             list_trn(fn);
410             break;
411         case efTNG:
412             list_tng(fn);
413             break;
414         default:
415             fprintf(stderr, "File %s is of an unsupported type. Try using the command\n 'less %s'\n",
416                     fn, fn);
417     }
418 }
419
420 void list_ene(const char *fn)
421 {
422     int            ndr;
423     ener_file_t    in;
424     gmx_bool       bCont;
425     gmx_enxnm_t   *enm = NULL;
426     t_enxframe    *fr;
427     int            i, j, nre, b;
428     real           rav, minthird;
429     char           buf[22];
430
431     printf("gmx dump: %s\n", fn);
432     in = open_enx(fn, "r");
433     do_enxnms(in, &nre, &enm);
434     assert(enm);
435
436     printf("energy components:\n");
437     for (i = 0; (i < nre); i++)
438     {
439         printf("%5d  %-24s (%s)\n", i, enm[i].name, enm[i].unit);
440     }
441
442     minthird = -1.0/3.0;
443     snew(fr, 1);
444     do
445     {
446         bCont = do_enx(in, fr);
447
448         if (bCont)
449         {
450             printf("\n%24s  %12.5e  %12s  %12s\n", "time:",
451                    fr->t, "step:", gmx_step_str(fr->step, buf));
452             printf("%24s  %12s  %12s  %12s\n",
453                    "", "", "nsteps:", gmx_step_str(fr->nsteps, buf));
454             printf("%24s  %12.5e  %12s  %12s\n",
455                    "delta_t:", fr->dt, "sum steps:", gmx_step_str(fr->nsum, buf));
456             if (fr->nre == nre)
457             {
458                 printf("%24s  %12s  %12s  %12s\n",
459                        "Component", "Energy", "Av. Energy", "Sum Energy");
460                 if (fr->nsum > 0)
461                 {
462                     for (i = 0; (i < nre); i++)
463                     {
464                         printf("%24s  %12.5e  %12.5e  %12.5e\n",
465                                enm[i].name, fr->ener[i].e, fr->ener[i].eav,
466                                fr->ener[i].esum);
467                     }
468                 }
469                 else
470                 {
471                     for (i = 0; (i < nre); i++)
472                     {
473                         printf("%24s  %12.5e\n",
474                                enm[i].name, fr->ener[i].e);
475                     }
476                 }
477             }
478             for (b = 0; b < fr->nblock; b++)
479             {
480                 const char *typestr = "";
481
482                 t_enxblock *eb = &(fr->block[b]);
483                 printf("Block data %2d (%3d subblocks, id=%d)\n",
484                        b, eb->nsub, eb->id);
485
486                 if (eb->id < enxNR)
487                 {
488                     typestr = enx_block_id_name[eb->id];
489                 }
490                 printf("  id='%s'\n", typestr);
491                 for (i = 0; i < eb->nsub; i++)
492                 {
493                     t_enxsubblock *sb = &(eb->sub[i]);
494                     printf("  Sub block %3d (%5d elems, type=%s) values:\n",
495                            i, sb->nr, xdr_datatype_names[sb->type]);
496
497                     switch (sb->type)
498                     {
499                         case xdr_datatype_float:
500                             for (j = 0; j < sb->nr; j++)
501                             {
502                                 printf("%14d   %8.4f\n", j, sb->fval[j]);
503                             }
504                             break;
505                         case xdr_datatype_double:
506                             for (j = 0; j < sb->nr; j++)
507                             {
508                                 printf("%14d   %10.6f\n", j, sb->dval[j]);
509                             }
510                             break;
511                         case xdr_datatype_int:
512                             for (j = 0; j < sb->nr; j++)
513                             {
514                                 printf("%14d %10d\n", j, sb->ival[j]);
515                             }
516                             break;
517                         case xdr_datatype_int64:
518                             for (j = 0; j < sb->nr; j++)
519                             {
520                                 printf("%14d %s\n",
521                                        j, gmx_step_str(sb->lval[j], buf));
522                             }
523                             break;
524                         case xdr_datatype_char:
525                             for (j = 0; j < sb->nr; j++)
526                             {
527                                 printf("%14d %1c\n", j, sb->cval[j]);
528                             }
529                             break;
530                         case xdr_datatype_string:
531                             for (j = 0; j < sb->nr; j++)
532                             {
533                                 printf("%14d %80s\n", j, sb->sval[j]);
534                             }
535                             break;
536                         default:
537                             gmx_incons("Unknown subblock type");
538                     }
539                 }
540             }
541         }
542     }
543     while (bCont);
544
545     close_enx(in);
546
547     free_enxframe(fr);
548     sfree(fr);
549     sfree(enm);
550 }
551
552 static void list_mtx(const char *fn)
553 {
554     int                  nrow, ncol, i, j, k;
555     real                *full   = NULL, value;
556     gmx_sparsematrix_t * sparse = NULL;
557
558     gmx_mtxio_read(fn, &nrow, &ncol, &full, &sparse);
559
560     if (full == NULL)
561     {
562         snew(full, nrow*ncol);
563         for (i = 0; i < nrow*ncol; i++)
564         {
565             full[i] = 0;
566         }
567
568         for (i = 0; i < sparse->nrow; i++)
569         {
570             for (j = 0; j < sparse->ndata[i]; j++)
571             {
572                 k              = sparse->data[i][j].col;
573                 value          = sparse->data[i][j].value;
574                 full[i*ncol+k] = value;
575                 full[k*ncol+i] = value;
576             }
577         }
578         gmx_sparsematrix_destroy(sparse);
579     }
580
581     printf("%d %d\n", nrow, ncol);
582     for (i = 0; i < nrow; i++)
583     {
584         for (j = 0; j < ncol; j++)
585         {
586             printf(" %g", full[i*ncol+j]);
587         }
588         printf("\n");
589     }
590
591     sfree(full);
592 }
593
594 int gmx_dump(int argc, char *argv[])
595 {
596     const char *desc[] = {
597         "[THISMODULE] reads a run input file ([REF].tpr[ref]),",
598         "a trajectory ([REF].trr[ref]/[REF].xtc[ref]/[TT]/tng[tt]), an energy",
599         "file ([REF].edr[ref]) or a checkpoint file ([REF].cpt[ref])",
600         "and prints that to standard output in a readable format.",
601         "This program is essential for checking your run input file in case of",
602         "problems.[PAR]",
603         "The program can also preprocess a topology to help finding problems.",
604         "Note that currently setting [TT]GMXLIB[tt] is the only way to customize",
605         "directories used for searching include files.",
606     };
607     const char *bugs[] = {
608         "Position restraint output from -sys -s is broken"
609     };
610     t_filenm    fnm[] = {
611         { efTPR, "-s", NULL, ffOPTRD },
612         { efTRX, "-f", NULL, ffOPTRD },
613         { efEDR, "-e", NULL, ffOPTRD },
614         { efCPT, NULL, NULL, ffOPTRD },
615         { efTOP, "-p", NULL, ffOPTRD },
616         { efMTX, "-mtx", "hessian", ffOPTRD },
617         { efMDP, "-om", NULL, ffOPTWR }
618     };
619 #define NFILE asize(fnm)
620
621     output_env_t    oenv;
622     /* Command line options */
623     static gmx_bool bShowNumbers = TRUE;
624     static gmx_bool bSysTop      = FALSE;
625     t_pargs         pa[]         = {
626         { "-nr", FALSE, etBOOL, {&bShowNumbers}, "Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" },
627         { "-sys", FALSE, etBOOL, {&bSysTop}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" }
628     };
629
630     if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
631                            asize(desc), desc, asize(bugs), bugs, &oenv))
632     {
633         return 0;
634     }
635
636
637     if (ftp2bSet(efTPR, NFILE, fnm))
638     {
639         list_tpx(ftp2fn(efTPR, NFILE, fnm), bShowNumbers,
640                  ftp2fn_null(efMDP, NFILE, fnm), bSysTop);
641     }
642     else if (ftp2bSet(efTRX, NFILE, fnm))
643     {
644         list_trx(ftp2fn(efTRX, NFILE, fnm));
645     }
646     else if (ftp2bSet(efEDR, NFILE, fnm))
647     {
648         list_ene(ftp2fn(efEDR, NFILE, fnm));
649     }
650     else if (ftp2bSet(efCPT, NFILE, fnm))
651     {
652         list_checkpoint(ftp2fn(efCPT, NFILE, fnm), stdout);
653     }
654     else if (ftp2bSet(efTOP, NFILE, fnm))
655     {
656         list_top(ftp2fn(efTOP, NFILE, fnm));
657     }
658     else if (ftp2bSet(efMTX, NFILE, fnm))
659     {
660         list_mtx(ftp2fn(efMTX, NFILE, fnm));
661     }
662
663     return 0;
664 }