6232cd1a7fe4592d76a4f2cbed0b5fd6341ca0d1
[alexxy/gromacs.git] / src / programs / gmx / gmxdump.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2013, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <stdio.h>
40 #include <string.h>
41 #include <math.h>
42 #include <assert.h>
43 #include "main.h"
44 #include "macros.h"
45 #include "gromacs/fileio/futil.h"
46 #include "statutil.h"
47 #include "sysstuff.h"
48 #include "txtdump.h"
49 #include "gmx_fatal.h"
50 #include "gromacs/fileio/xtcio.h"
51 #include "gromacs/fileio/enxio.h"
52 #include "smalloc.h"
53 #include "names.h"
54 #include "gromacs/fileio/gmxfio.h"
55 #include "gromacs/fileio/tpxio.h"
56 #include "gromacs/fileio/trnio.h"
57 #include "txtdump.h"
58 #include "gmxcpp.h"
59 #include "checkpoint.h"
60 #include "mtop_util.h"
61
62 #include "gromacs/linearalgebra/mtxio.h"
63 #include "gromacs/linearalgebra/sparsematrix.h"
64
65
66 static void list_tpx(const char *fn, gmx_bool bShowNumbers, const char *mdpfn,
67                      gmx_bool bSysTop)
68 {
69     FILE         *gp;
70     int           fp, indent, i, j, **gcount, atot;
71     t_state       state;
72     rvec         *f = NULL;
73     t_inputrec    ir;
74     t_tpxheader   tpx;
75     gmx_mtop_t    mtop;
76     gmx_groups_t *groups;
77     t_topology    top;
78
79     read_tpxheader(fn, &tpx, TRUE, NULL, NULL);
80
81     read_tpx_state(fn,
82                    tpx.bIr  ? &ir : NULL,
83                    &state, tpx.bF ? f : NULL,
84                    tpx.bTop ? &mtop : NULL);
85
86     if (mdpfn && tpx.bIr)
87     {
88         gp = gmx_fio_fopen(mdpfn, "w");
89         pr_inputrec(gp, 0, NULL, &(ir), TRUE);
90         gmx_fio_fclose(gp);
91     }
92
93     if (!mdpfn)
94     {
95         if (bSysTop)
96         {
97             top = gmx_mtop_t_to_t_topology(&mtop);
98         }
99
100         if (available(stdout, &tpx, 0, fn))
101         {
102             indent = 0;
103             indent = pr_title(stdout, indent, fn);
104             pr_inputrec(stdout, 0, "inputrec", tpx.bIr ? &(ir) : NULL, FALSE);
105
106             indent = 0;
107             pr_header(stdout, indent, "header", &(tpx));
108
109             if (!bSysTop)
110             {
111                 pr_mtop(stdout, indent, "topology", &(mtop), bShowNumbers);
112             }
113             else
114             {
115                 pr_top(stdout, indent, "topology", &(top), bShowNumbers);
116             }
117
118             pr_rvecs(stdout, indent, "box", tpx.bBox ? state.box : NULL, DIM);
119             pr_rvecs(stdout, indent, "box_rel", tpx.bBox ? state.box_rel : NULL, DIM);
120             pr_rvecs(stdout, indent, "boxv", tpx.bBox ? state.boxv : NULL, DIM);
121             pr_rvecs(stdout, indent, "pres_prev", tpx.bBox ? state.pres_prev : NULL, DIM);
122             pr_rvecs(stdout, indent, "svir_prev", tpx.bBox ? state.svir_prev : NULL, DIM);
123             pr_rvecs(stdout, indent, "fvir_prev", tpx.bBox ? state.fvir_prev : NULL, DIM);
124             /* leave nosehoover_xi in for now to match the tpr version */
125             pr_doubles(stdout, indent, "nosehoover_xi", state.nosehoover_xi, state.ngtc);
126             /*pr_doubles(stdout,indent,"nosehoover_vxi",state.nosehoover_vxi,state.ngtc);*/
127             /*pr_doubles(stdout,indent,"therm_integral",state.therm_integral,state.ngtc);*/
128             pr_rvecs(stdout, indent, "x", tpx.bX ? state.x : NULL, state.natoms);
129             pr_rvecs(stdout, indent, "v", tpx.bV ? state.v : NULL, state.natoms);
130             if (tpx.bF)
131             {
132                 pr_rvecs(stdout, indent, "f", f, state.natoms);
133             }
134         }
135
136         groups = &mtop.groups;
137
138         snew(gcount, egcNR);
139         for (i = 0; (i < egcNR); i++)
140         {
141             snew(gcount[i], groups->grps[i].nr);
142         }
143
144         for (i = 0; (i < mtop.natoms); i++)
145         {
146             for (j = 0; (j < egcNR); j++)
147             {
148                 gcount[j][ggrpnr(groups, j, i)]++;
149             }
150         }
151         printf("Group statistics\n");
152         for (i = 0; (i < egcNR); i++)
153         {
154             atot = 0;
155             printf("%-12s: ", gtypes[i]);
156             for (j = 0; (j < groups->grps[i].nr); j++)
157             {
158                 printf("  %5d", gcount[i][j]);
159                 atot += gcount[i][j];
160             }
161             printf("  (total %d atoms)\n", atot);
162             sfree(gcount[i]);
163         }
164         sfree(gcount);
165     }
166     done_state(&state);
167     sfree(f);
168 }
169
170 static void list_top(const char *fn)
171 {
172     int       status, done;
173 #define BUFLEN 256
174     char      buf[BUFLEN];
175     gmx_cpp_t handle;
176     char     *cppopts[] = { NULL };
177
178     status = cpp_open_file(fn, &handle, cppopts);
179     if (status != 0)
180     {
181         gmx_fatal(FARGS, cpp_error(&handle, status));
182     }
183     do
184     {
185         status = cpp_read_line(&handle, BUFLEN, buf);
186         done   = (status == eCPP_EOF);
187         if (!done)
188         {
189             if (status != eCPP_OK)
190             {
191                 gmx_fatal(FARGS, cpp_error(&handle, status));
192             }
193             else
194             {
195                 printf("%s\n", buf);
196             }
197         }
198     }
199     while (!done);
200     status = cpp_close_file(&handle);
201     if (status != eCPP_OK)
202     {
203         gmx_fatal(FARGS, cpp_error(&handle, status));
204     }
205 }
206
207 static void list_trn(const char *fn)
208 {
209     t_fileio       *fpread, *fpwrite;
210     int             nframe, indent;
211     char            buf[256];
212     rvec           *x, *v, *f;
213     matrix          box;
214     t_trnheader     trn;
215     gmx_bool        bOK;
216
217     fpread  = open_trn(fn, "r");
218     fpwrite = open_tpx(NULL, "w");
219     gmx_fio_setdebug(fpwrite, TRUE);
220
221     nframe = 0;
222     while (fread_trnheader(fpread, &trn, &bOK))
223     {
224         snew(x, trn.natoms);
225         snew(v, trn.natoms);
226         snew(f, trn.natoms);
227         if (fread_htrn(fpread, &trn,
228                        trn.box_size ? box : NULL,
229                        trn.x_size   ? x : NULL,
230                        trn.v_size   ? v : NULL,
231                        trn.f_size   ? f : NULL))
232         {
233             sprintf(buf, "%s frame %d", fn, nframe);
234             indent = 0;
235             indent = pr_title(stdout, indent, buf);
236             pr_indent(stdout, indent);
237             fprintf(stdout, "natoms=%10d  step=%10d  time=%12.7e  lambda=%10g\n",
238                     trn.natoms, trn.step, trn.t, trn.lambda);
239             if (trn.box_size)
240             {
241                 pr_rvecs(stdout, indent, "box", box, DIM);
242             }
243             if (trn.x_size)
244             {
245                 pr_rvecs(stdout, indent, "x", x, trn.natoms);
246             }
247             if (trn.v_size)
248             {
249                 pr_rvecs(stdout, indent, "v", v, trn.natoms);
250             }
251             if (trn.f_size)
252             {
253                 pr_rvecs(stdout, indent, "f", f, trn.natoms);
254             }
255         }
256         else
257         {
258             fprintf(stderr, "\nWARNING: Incomplete frame: nr %d, t=%g\n",
259                     nframe, trn.t);
260         }
261
262         sfree(x);
263         sfree(v);
264         sfree(f);
265         nframe++;
266     }
267     if (!bOK)
268     {
269         fprintf(stderr, "\nWARNING: Incomplete frame header: nr %d, t=%g\n",
270                 nframe, trn.t);
271     }
272     close_tpx(fpwrite);
273     close_trn(fpread);
274 }
275
276 void list_xtc(const char *fn, gmx_bool bXVG)
277 {
278     t_fileio  *xd;
279     int        indent;
280     char       buf[256];
281     rvec      *x;
282     matrix     box;
283     int        nframe, natoms, step;
284     real       prec, time;
285     gmx_bool   bOK;
286
287     xd = open_xtc(fn, "r");
288     read_first_xtc(xd, &natoms, &step, &time, box, &x, &prec, &bOK);
289
290     nframe = 0;
291     do
292     {
293         if (bXVG)
294         {
295             int i, d;
296
297             fprintf(stdout, "%g", time);
298             for (i = 0; i < natoms; i++)
299             {
300                 for (d = 0; d < DIM; d++)
301                 {
302                     fprintf(stdout, " %g", x[i][d]);
303                 }
304             }
305             fprintf(stdout, "\n");
306         }
307         else
308         {
309             sprintf(buf, "%s frame %d", fn, nframe);
310             indent = 0;
311             indent = pr_title(stdout, indent, buf);
312             pr_indent(stdout, indent);
313             fprintf(stdout, "natoms=%10d  step=%10d  time=%12.7e  prec=%10g\n",
314                     natoms, step, time, prec);
315             pr_rvecs(stdout, indent, "box", box, DIM);
316             pr_rvecs(stdout, indent, "x", x, natoms);
317         }
318         nframe++;
319     }
320     while (read_next_xtc(xd, natoms, &step, &time, box, x, &prec, &bOK));
321     if (!bOK)
322     {
323         fprintf(stderr, "\nWARNING: Incomplete frame at time %g\n", time);
324     }
325     close_xtc(xd);
326 }
327
328 void list_trx(const char *fn, gmx_bool bXVG)
329 {
330     int ftp;
331
332     ftp = fn2ftp(fn);
333     if (ftp == efXTC)
334     {
335         list_xtc(fn, bXVG);
336     }
337     else if ((ftp == efTRR) || (ftp == efTRJ))
338     {
339         list_trn(fn);
340     }
341     else
342     {
343         fprintf(stderr, "File %s is of an unsupported type. Try using the command\n 'less %s'\n",
344                 fn, fn);
345     }
346 }
347
348 void list_ene(const char *fn)
349 {
350     int            ndr;
351     ener_file_t    in;
352     gmx_bool       bCont;
353     gmx_enxnm_t   *enm = NULL;
354     t_enxframe    *fr;
355     int            i, j, nre, b;
356     real           rav, minthird;
357     char           buf[22];
358
359     printf("gmx dump: %s\n", fn);
360     in = open_enx(fn, "r");
361     do_enxnms(in, &nre, &enm);
362     assert(enm);
363
364     printf("energy components:\n");
365     for (i = 0; (i < nre); i++)
366     {
367         printf("%5d  %-24s (%s)\n", i, enm[i].name, enm[i].unit);
368     }
369
370     minthird = -1.0/3.0;
371     snew(fr, 1);
372     do
373     {
374         bCont = do_enx(in, fr);
375
376         if (bCont)
377         {
378             printf("\n%24s  %12.5e  %12s  %12s\n", "time:",
379                    fr->t, "step:", gmx_step_str(fr->step, buf));
380             printf("%24s  %12s  %12s  %12s\n",
381                    "", "", "nsteps:", gmx_step_str(fr->nsteps, buf));
382             printf("%24s  %12.5e  %12s  %12s\n",
383                    "delta_t:", fr->dt, "sum steps:", gmx_step_str(fr->nsum, buf));
384             if (fr->nre == nre)
385             {
386                 printf("%24s  %12s  %12s  %12s\n",
387                        "Component", "Energy", "Av. Energy", "Sum Energy");
388                 if (fr->nsum > 0)
389                 {
390                     for (i = 0; (i < nre); i++)
391                     {
392                         printf("%24s  %12.5e  %12.5e  %12.5e\n",
393                                enm[i].name, fr->ener[i].e, fr->ener[i].eav,
394                                fr->ener[i].esum);
395                     }
396                 }
397                 else
398                 {
399                     for (i = 0; (i < nre); i++)
400                     {
401                         printf("%24s  %12.5e\n",
402                                enm[i].name, fr->ener[i].e);
403                     }
404                 }
405             }
406             for (b = 0; b < fr->nblock; b++)
407             {
408                 const char *typestr = "";
409
410                 t_enxblock *eb = &(fr->block[b]);
411                 printf("Block data %2d (%3d subblocks, id=%d)\n",
412                        b, eb->nsub, eb->id);
413
414                 if (eb->id < enxNR)
415                 {
416                     typestr = enx_block_id_name[eb->id];
417                 }
418                 printf("  id='%s'\n", typestr);
419                 for (i = 0; i < eb->nsub; i++)
420                 {
421                     t_enxsubblock *sb = &(eb->sub[i]);
422                     printf("  Sub block %3d (%5d elems, type=%s) values:\n",
423                            i, sb->nr, xdr_datatype_names[sb->type]);
424
425                     switch (sb->type)
426                     {
427                         case xdr_datatype_float:
428                             for (j = 0; j < sb->nr; j++)
429                             {
430                                 printf("%14d   %8.4f\n", j, sb->fval[j]);
431                             }
432                             break;
433                         case xdr_datatype_double:
434                             for (j = 0; j < sb->nr; j++)
435                             {
436                                 printf("%14d   %10.6f\n", j, sb->dval[j]);
437                             }
438                             break;
439                         case xdr_datatype_int:
440                             for (j = 0; j < sb->nr; j++)
441                             {
442                                 printf("%14d %10d\n", j, sb->ival[j]);
443                             }
444                             break;
445                         case xdr_datatype_large_int:
446                             for (j = 0; j < sb->nr; j++)
447                             {
448                                 printf("%14d %s\n",
449                                        j, gmx_step_str(sb->lval[j], buf));
450                             }
451                             break;
452                         case xdr_datatype_char:
453                             for (j = 0; j < sb->nr; j++)
454                             {
455                                 printf("%14d %1c\n", j, sb->cval[j]);
456                             }
457                             break;
458                         case xdr_datatype_string:
459                             for (j = 0; j < sb->nr; j++)
460                             {
461                                 printf("%14d %80s\n", j, sb->sval[j]);
462                             }
463                             break;
464                         default:
465                             gmx_incons("Unknown subblock type");
466                     }
467                 }
468             }
469         }
470     }
471     while (bCont);
472
473     close_enx(in);
474
475     free_enxframe(fr);
476     sfree(fr);
477     sfree(enm);
478 }
479
480 static void list_mtx(const char *fn)
481 {
482     int                  nrow, ncol, i, j, k;
483     real                *full   = NULL, value;
484     gmx_sparsematrix_t * sparse = NULL;
485
486     gmx_mtxio_read(fn, &nrow, &ncol, &full, &sparse);
487
488     if (full == NULL)
489     {
490         snew(full, nrow*ncol);
491         for (i = 0; i < nrow*ncol; i++)
492         {
493             full[i] = 0;
494         }
495
496         for (i = 0; i < sparse->nrow; i++)
497         {
498             for (j = 0; j < sparse->ndata[i]; j++)
499             {
500                 k              = sparse->data[i][j].col;
501                 value          = sparse->data[i][j].value;
502                 full[i*ncol+k] = value;
503                 full[k*ncol+i] = value;
504             }
505         }
506         gmx_sparsematrix_destroy(sparse);
507     }
508
509     printf("%d %d\n", nrow, ncol);
510     for (i = 0; i < nrow; i++)
511     {
512         for (j = 0; j < ncol; j++)
513         {
514             printf(" %g", full[i*ncol+j]);
515         }
516         printf("\n");
517     }
518
519     sfree(full);
520 }
521
522 int gmx_gmxdump(int argc, char *argv[])
523 {
524     const char *desc[] = {
525         "[TT]gmx dump[tt] reads a run input file ([TT].tpa[tt]/[TT].tpr[tt]/[TT].tpb[tt]),",
526         "a trajectory ([TT].trj[tt]/[TT].trr[tt]/[TT].xtc[tt]), an energy",
527         "file ([TT].ene[tt]/[TT].edr[tt]), or a checkpoint file ([TT].cpt[tt])",
528         "and prints that to standard output in a readable format.",
529         "This program is essential for checking your run input file in case of",
530         "problems.[PAR]",
531         "The program can also preprocess a topology to help finding problems.",
532         "Note that currently setting [TT]GMXLIB[tt] is the only way to customize",
533         "directories used for searching include files.",
534     };
535     const char *bugs[] = {
536         "Position restraint output from -sys -s is broken"
537     };
538     t_filenm    fnm[] = {
539         { efTPX, "-s", NULL, ffOPTRD },
540         { efTRX, "-f", NULL, ffOPTRD },
541         { efEDR, "-e", NULL, ffOPTRD },
542         { efCPT, NULL, NULL, ffOPTRD },
543         { efTOP, "-p", NULL, ffOPTRD },
544         { efMTX, "-mtx", "hessian", ffOPTRD },
545         { efMDP, "-om", NULL, ffOPTWR }
546     };
547 #define NFILE asize(fnm)
548
549     output_env_t    oenv;
550     /* Command line options */
551     static gmx_bool bXVG         = FALSE;
552     static gmx_bool bShowNumbers = TRUE;
553     static gmx_bool bSysTop      = FALSE;
554     t_pargs         pa[]         = {
555         { "-xvg", FALSE, etBOOL, {&bXVG}, "HIDDENXVG layout for xtc" },
556         { "-nr", FALSE, etBOOL, {&bShowNumbers}, "Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" },
557         { "-sys", FALSE, etBOOL, {&bSysTop}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" }
558     };
559
560     if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
561                            asize(desc), desc, asize(bugs), bugs, &oenv))
562     {
563         return 0;
564     }
565
566
567     if (ftp2bSet(efTPX, NFILE, fnm))
568     {
569         list_tpx(ftp2fn(efTPX, NFILE, fnm), bShowNumbers,
570                  ftp2fn_null(efMDP, NFILE, fnm), bSysTop);
571     }
572     else if (ftp2bSet(efTRX, NFILE, fnm))
573     {
574         list_trx(ftp2fn(efTRX, NFILE, fnm), bXVG);
575     }
576     else if (ftp2bSet(efEDR, NFILE, fnm))
577     {
578         list_ene(ftp2fn(efEDR, NFILE, fnm));
579     }
580     else if (ftp2bSet(efCPT, NFILE, fnm))
581     {
582         list_checkpoint(ftp2fn(efCPT, NFILE, fnm), stdout);
583     }
584     else if (ftp2bSet(efTOP, NFILE, fnm))
585     {
586         list_top(ftp2fn(efTOP, NFILE, fnm));
587     }
588     else if (ftp2bSet(efMTX, NFILE, fnm))
589     {
590         list_mtx(ftp2fn(efMTX, NFILE, fnm));
591     }
592
593     return 0;
594 }