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