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