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