OpenMM: added combination rule check, disabled restraint check for now as it's buggy
[gromacs/qmmm-gamess-us.git] / src / kernel / gmxdump.c
blobc6f023b2154eb515230a0ab659546ae6a652cd75
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
39 #include <stdio.h>
40 #include <string.h>
41 #include <math.h>
42 #include "main.h"
43 #include "macros.h"
44 #include "futil.h"
45 #include "statutil.h"
46 #include "copyrite.h"
47 #include "sysstuff.h"
48 #include "txtdump.h"
49 #include "gmx_fatal.h"
50 #include "xtcio.h"
51 #include "enxio.h"
52 #include "smalloc.h"
53 #include "names.h"
54 #include "gmxfio.h"
55 #include "tpxio.h"
56 #include "trnio.h"
57 #include "txtdump.h"
58 #include "gmxcpp.h"
59 #include "checkpoint.h"
60 #include "mtop_util.h"
61 #include "sparsematrix.h"
62 #include "mtxio.h"
64 static void dump_top(FILE *fp,t_topology *top,char *tpr)
66 int i,j,k,*types;
68 fprintf(fp,"; Topology generated from %s by program %s\n",tpr,Program());
69 fprintf(fp,"[ defaults ]\n 1 1 no 1.0 1.0\n\n");
70 fprintf(fp,"[ atomtypes ]\n");
71 fprintf(fp,";name at.num mass charge ptype c6 c12\n");
72 snew(types,top->atomtypes.nr);
73 for(i=0; (i<top->atomtypes.nr); i++) {
74 for(j=0; (j<top->atoms.nr) && (top->atoms.atom[j].type != i); j++)
76 if (j<top->atoms.nr) {
77 types[i] = j;
78 fprintf(fp,"%5s %4d %8.4f %8.4f %2s %8.3f %8.3f\n",
79 *(top->atoms.atomtype[j]),top->atomtypes.atomnumber[i],
80 0.0,0.0,"A",0.0,0.0);
83 fprintf(fp,"[ nonbonded_params ]\n");
84 for(i=k=0; (i<top->idef.ntypes); i++) {
85 for(j=0; (j<top->idef.ntypes); j++,k++) {
86 fprintf(fp,"%12s %12s 1 %12.5e %12.5e\n",
87 *(top->atoms.atomtype[types[i]]),
88 *(top->atoms.atomtype[types[j]]),
89 top->idef.iparams[k].lj.c12,top->idef.iparams[k].lj.c6);
92 sfree(types);
95 static void list_tpx(const char *fn, bool bShowNumbers,const char *mdpfn,
96 bool bSysTop)
98 FILE *gp;
99 int fp,indent,i,j,**gcount,atot;
100 t_state state;
101 rvec *f=NULL;
102 t_inputrec ir;
103 t_tpxheader tpx;
104 gmx_mtop_t mtop;
105 gmx_groups_t *groups;
106 t_topology top;
108 read_tpxheader(fn,&tpx,TRUE,NULL,NULL);
110 read_tpx_state(fn,
111 tpx.bIr ? &ir : NULL,
112 &state,tpx.bF ? f : NULL,
113 tpx.bTop ? &mtop: NULL);
115 if (mdpfn && tpx.bIr) {
116 gp = gmx_fio_fopen(mdpfn,"w");
117 pr_inputrec(gp,0,NULL,&(ir),TRUE);
118 gmx_fio_fclose(gp);
121 if (!mdpfn) {
122 if (bSysTop)
123 top = gmx_mtop_t_to_t_topology(&mtop);
125 if (available(stdout,&tpx,0,fn)) {
126 indent=0;
127 indent=pr_title(stdout,indent,fn);
128 pr_inputrec(stdout,0,"inputrec",tpx.bIr ? &(ir) : NULL,FALSE);
130 indent = 0;
131 pr_header(stdout,indent,"header",&(tpx));
133 if (!bSysTop)
134 pr_mtop(stdout,indent,"topology",&(mtop),bShowNumbers);
135 else
136 pr_top(stdout,indent,"topology",&(top),bShowNumbers);
138 pr_rvecs(stdout,indent,"box",tpx.bBox ? state.box : NULL,DIM);
139 pr_rvecs(stdout,indent,"box_rel",tpx.bBox ? state.box_rel : NULL,DIM);
140 pr_rvecs(stdout,indent,"boxv",tpx.bBox ? state.boxv : NULL,DIM);
141 pr_rvecs(stdout,indent,"pres_prev",tpx.bBox ? state.pres_prev : NULL,DIM);
142 pr_rvecs(stdout,indent,"svir_prev",tpx.bBox ? state.svir_prev : NULL,DIM);
143 pr_rvecs(stdout,indent,"fvir_prev",tpx.bBox ? state.fvir_prev : NULL,DIM);
144 /* leave nosehoover_xi in for now to match the tpr version */
145 pr_doubles(stdout,indent,"nosehoover_xi",state.nosehoover_xi,state.ngtc);
146 /*pr_doubles(stdout,indent,"nosehoover_vxi",state.nosehoover_vxi,state.ngtc);*/
147 /*pr_doubles(stdout,indent,"therm_integral",state.therm_integral,state.ngtc);*/
148 pr_rvecs(stdout,indent,"x",tpx.bX ? state.x : NULL,state.natoms);
149 pr_rvecs(stdout,indent,"v",tpx.bV ? state.v : NULL,state.natoms);
150 if (state,tpx.bF) {
151 pr_rvecs(stdout,indent,"f",f,state.natoms);
155 groups = &mtop.groups;
157 snew(gcount,egcNR);
158 for(i=0; (i<egcNR); i++)
159 snew(gcount[i],groups->grps[i].nr);
161 for(i=0; (i<mtop.natoms); i++) {
162 for(j=0; (j<egcNR); j++)
163 gcount[j][ggrpnr(groups,j,i)]++;
165 printf("Group statistics\n");
166 for(i=0; (i<egcNR); i++) {
167 atot=0;
168 printf("%-12s: ",gtypes[i]);
169 for(j=0; (j<groups->grps[i].nr); j++) {
170 printf(" %5d",gcount[i][j]);
171 atot+=gcount[i][j];
173 printf(" (total %d atoms)\n",atot);
174 sfree(gcount[i]);
176 sfree(gcount);
178 done_state(&state);
179 sfree(f);
182 static void list_top(const char *fn)
184 int status,done;
185 #define BUFLEN 256
186 char buf[BUFLEN];
187 gmx_cpp_t handle;
188 char *cppopts[] = { NULL };
190 status = cpp_open_file(fn,&handle,cppopts);
191 if (status != 0)
192 gmx_fatal(FARGS,cpp_error(&handle,status));
193 do {
194 status = cpp_read_line(&handle,BUFLEN,buf);
195 done = (status == eCPP_EOF);
196 if (!done) {
197 if (status != eCPP_OK)
198 gmx_fatal(FARGS,cpp_error(&handle,status));
199 else
200 printf("%s\n",buf);
202 } while (!done);
203 status = cpp_close_file(&handle);
204 if (status != eCPP_OK)
205 gmx_fatal(FARGS,cpp_error(&handle,status));
208 static void list_trn(const char *fn)
210 int fpread,fpwrite,nframe,indent;
211 char buf[256];
212 rvec *x,*v,*f;
213 matrix box;
214 t_trnheader trn;
215 bool bOK;
217 fpread = open_trn(fn,"r");
218 fpwrite = open_tpx(NULL,"w");
219 gmx_fio_setdebug(fpwrite,TRUE);
221 nframe = 0;
222 while (fread_trnheader(fpread,&trn,&bOK)) {
223 snew(x,trn.natoms);
224 snew(v,trn.natoms);
225 snew(f,trn.natoms);
226 if (fread_htrn(fpread,&trn,
227 trn.box_size ? box : NULL,
228 trn.x_size ? x : NULL,
229 trn.v_size ? v : NULL,
230 trn.f_size ? f : NULL)) {
231 sprintf(buf,"%s frame %d",fn,nframe);
232 indent=0;
233 indent=pr_title(stdout,indent,buf);
234 pr_indent(stdout,indent);
235 fprintf(stdout,"natoms=%10d step=%10d time=%12.7e lambda=%10g\n",
236 trn.natoms,trn.step,trn.t,trn.lambda);
237 if (trn.box_size)
238 pr_rvecs(stdout,indent,"box",box,DIM);
239 if (trn.x_size)
240 pr_rvecs(stdout,indent,"x",x,trn.natoms);
241 if (trn.v_size)
242 pr_rvecs(stdout,indent,"v",v,trn.natoms);
243 if (trn.f_size)
244 pr_rvecs(stdout,indent,"f",f,trn.natoms);
246 else
247 fprintf(stderr,"\nWARNING: Incomplete frame: nr %d, t=%g\n",
248 nframe,trn.t);
250 sfree(x);
251 sfree(v);
252 sfree(f);
253 nframe++;
255 if (!bOK)
256 fprintf(stderr,"\nWARNING: Incomplete frame header: nr %d, t=%g\n",
257 nframe,trn.t);
258 close_tpx(fpwrite);
259 close_trn(fpread);
262 void list_xtc(const char *fn, bool bXVG)
264 int xd,indent;
265 char buf[256];
266 rvec *x;
267 matrix box;
268 int nframe,natoms,step;
269 real prec,time;
270 bool bOK;
272 xd = open_xtc(fn,"r");
273 read_first_xtc(xd,&natoms,&step,&time,box,&x,&prec,&bOK);
275 nframe=0;
276 do {
277 if (bXVG) {
278 int i,d;
280 fprintf(stdout,"%g",time);
281 for(i=0; i<natoms; i++)
282 for(d=0; d<DIM; d++)
283 fprintf(stdout," %g",x[i][d]);
284 fprintf(stdout,"\n");
285 } else {
286 sprintf(buf,"%s frame %d",fn,nframe);
287 indent=0;
288 indent=pr_title(stdout,indent,buf);
289 pr_indent(stdout,indent);
290 fprintf(stdout,"natoms=%10d step=%10d time=%12.7e prec=%10g\n",
291 natoms,step,time,prec);
292 pr_rvecs(stdout,indent,"box",box,DIM);
293 pr_rvecs(stdout,indent,"x",x,natoms);
295 nframe++;
296 } while (read_next_xtc(xd,natoms,&step,&time,box,x,&prec,&bOK));
297 if (!bOK)
298 fprintf(stderr,"\nWARNING: Incomplete frame at time %g\n",time);
299 close_xtc(xd);
302 void list_trx(const char *fn,bool bXVG)
304 int ftp;
306 ftp = fn2ftp(fn);
307 if (ftp == efXTC)
308 list_xtc(fn,bXVG);
309 else if ((ftp == efTRR) || (ftp == efTRJ))
310 list_trn(fn);
311 else
312 fprintf(stderr,"File %s not supported. Try using more %s\n",
313 fn,fn);
316 void list_ene(const char *fn)
318 int ndr;
319 ener_file_t in;
320 bool bCont;
321 gmx_enxnm_t *enm=NULL;
322 t_enxframe *fr;
323 int i,nre,b;
324 real rav,minthird;
325 char buf[22];
327 printf("gmxdump: %s\n",fn);
328 in = open_enx(fn,"r");
329 do_enxnms(in,&nre,&enm);
331 printf("energy components:\n");
332 for(i=0; (i<nre); i++)
333 printf("%5d %-24s (%s)\n",i,enm[i].name,enm[i].unit);
335 minthird=-1.0/3.0;
336 snew(fr,1);
337 do {
338 bCont=do_enx(in,fr);
340 if (bCont) {
341 printf("\n%24s %12.5e %12s %12s\n","time:",
342 fr->t,"step:",gmx_step_str(fr->step,buf));
343 printf("%24s %12s %12s %12s\n",
344 "","","nsteps:",gmx_step_str(fr->nsteps,buf));
345 printf("%24s %12s %12s %12s\n",
346 "","","sum steps:",gmx_step_str(fr->nsum,buf));
347 if (fr->nre == nre) {
348 printf("%24s %12s %12s %12s\n",
349 "Component","Energy","Av. Energy","Sum Energy");
350 if (fr->nsum > 0) {
351 for(i=0; (i<nre); i++)
352 printf("%24s %12.5e %12.5e %12.5e\n",
353 enm[i].name,fr->ener[i].e,fr->ener[i].eav,fr->ener[i].esum);
354 } else {
355 for(i=0; (i<nre); i++)
356 printf("%24s %12.5e\n",
357 enm[i].name,fr->ener[i].e);
360 if (fr->ndisre > 0) {
361 printf("Distance restraint %8s %8s\n","r(t)","<r^-3>^-3");
362 for(i=0; i<fr->ndisre; i++) {
363 rav=pow(fr->disre_rm3tav[i],minthird);
364 printf("%17d %8.4f %8.4f\n",i,fr->disre_rt[i],rav);
367 for(b=0; b<fr->nblock; b++)
368 if (fr->nr[b] > 0) {
369 printf("Block data %2d (%4d elm.) %8s\n",b,fr->nr[b],"value");
370 for(i=0; i<fr->nr[b]; i++)
371 printf("%24d %8.4f\n",i,fr->block[b][i]);
374 } while (bCont);
376 close_enx(in);
378 free_enxframe(fr);
379 sfree(fr);
380 sfree(enm);
383 static void list_mtx(const char *fn)
385 int nrow,ncol,i,j,k;
386 real *full=NULL,value;
387 gmx_sparsematrix_t * sparse=NULL;
389 gmx_mtxio_read(fn,&nrow,&ncol,&full,&sparse);
391 if (full == NULL) {
392 snew(full,nrow*ncol);
393 for(i=0;i<nrow*ncol;i++) {
394 full[i] = 0;
397 for(i=0;i<sparse->nrow;i++) {
398 for(j=0;j<sparse->ndata[i];j++) {
399 k = sparse->data[i][j].col;
400 value = sparse->data[i][j].value;
401 full[i*ncol+k] = value;
402 full[k*ncol+i] = value;
405 gmx_sparsematrix_destroy(sparse);
408 printf("%d %d\n",nrow,ncol);
409 for(i=0; i<nrow; i++) {
410 for(j=0; j<ncol; j++) {
411 printf(" %g",full[i*ncol+j]);
413 printf("\n");
416 sfree(full);
419 int main(int argc,char *argv[])
421 const char *desc[] = {
422 "gmxdump reads a run input file ([TT].tpa[tt]/[TT].tpr[tt]/[TT].tpb[tt]),",
423 "a trajectory ([TT].trj[tt]/[TT].trr[tt]/[TT].xtc[tt]), an energy",
424 "file ([TT].ene[tt]/[TT].edr[tt]), or a checkpoint file ([TT].cpt[tt])",
425 "and prints that to standard output in a readable format.",
426 "This program is essential for checking your run input file in case of",
427 "problems.[PAR]",
428 "The program can also preprocess a topology to help finding problems.",
429 "Note that currently setting GMXLIB is the only way to customize",
430 "directories used for searching include files.",
432 t_filenm fnm[] = {
433 { efTPX, "-s", NULL, ffOPTRD },
434 { efTRX, "-f", NULL, ffOPTRD },
435 { efEDR, "-e", NULL, ffOPTRD },
436 { efCPT, NULL, NULL, ffOPTRD },
437 { efTOP, "-p", NULL, ffOPTRD },
438 { efMTX, "-mtx", "hessian", ffOPTRD },
439 { efMDP, "-om", NULL, ffOPTWR }
441 #define NFILE asize(fnm)
443 output_env_t oenv;
444 /* Command line options */
445 static bool bXVG=FALSE;
446 static bool bShowNumbers=TRUE;
447 static bool bSysTop=FALSE;
448 t_pargs pa[] = {
449 { "-xvg", FALSE, etBOOL, {&bXVG}, "HIDDENXVG layout for xtc" },
450 { "-nr",FALSE, etBOOL, {&bShowNumbers},"Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" },
451 { "-sys", FALSE, etBOOL, {&bSysTop}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" }
454 CopyRight(stderr,argv[0]);
455 parse_common_args(&argc,argv,0,NFILE,fnm,asize(pa),pa,
456 asize(desc),desc,0,NULL,&oenv);
459 if (ftp2bSet(efTPX,NFILE,fnm))
460 list_tpx(ftp2fn(efTPX,NFILE,fnm),bShowNumbers,
461 ftp2fn_null(efMDP,NFILE,fnm),bSysTop);
462 else if (ftp2bSet(efTRX,NFILE,fnm))
463 list_trx(ftp2fn(efTRX,NFILE,fnm),bXVG);
464 else if (ftp2bSet(efEDR,NFILE,fnm))
465 list_ene(ftp2fn(efEDR,NFILE,fnm));
466 else if (ftp2bSet(efCPT,NFILE,fnm))
467 list_checkpoint(ftp2fn(efCPT,NFILE,fnm),stdout);
468 else if (ftp2bSet(efTOP,NFILE,fnm))
469 list_top(ftp2fn(efTOP,NFILE,fnm));
470 else if (ftp2bSet(efMTX,NFILE,fnm))
471 list_mtx(ftp2fn(efMTX,NFILE,fnm));
473 thanx(stderr);
475 return 0;