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