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