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