Parallel vs. sequentiual code: I get binary identical result (apart from comment...
[gromacs/qmmm-gamess-us.git] / src / tools / gmx_bond.c
blob0b0ec2eb1b0b74dd326ba7c44ab7328a1e677ca9
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 * Green Red Orange Magenta Azure Cyan Skyblue
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38 #include <math.h>
39 #include <string.h>
40 #include "sysstuff.h"
41 #include "typedefs.h"
42 #include "smalloc.h"
43 #include "macros.h"
44 #include "vec.h"
45 #include "pbc.h"
46 #include "xvgr.h"
47 #include "copyrite.h"
48 #include "gmx_fatal.h"
49 #include "futil.h"
50 #include "statutil.h"
51 #include "index.h"
52 #include "gmx_statistics.h"
53 #include "tpxio.h"
54 #include "gmx_ana.h"
57 static void make_dist_leg(FILE *fp,int gnx,atom_id index[],t_atoms *atoms,
58 const output_env_t oenv)
60 char **leg;
61 int i;
63 snew(leg,gnx/2);
64 for(i=0; i<gnx; i+=2) {
65 snew(leg[i/2],256);
66 sprintf(leg[i/2],"%s %s%d - %s %s%d",
67 *(atoms->atomname[index[i]]),
68 *(atoms->resinfo[atoms->atom[index[i]].resind].name),
69 atoms->resinfo[atoms->atom[index[i]].resind].nr,
70 *(atoms->atomname[index[i+1]]),
71 *(atoms->resinfo[atoms->atom[index[i+1]].resind].name),
72 atoms->resinfo[atoms->atom[index[i+1]].resind].nr);
74 xvgr_legend(fp,gnx/2,leg,oenv);
75 for(i=0; i<gnx/2; i++)
76 sfree(leg[i]);
77 sfree(leg);
80 static void do_bonds(FILE *log,const char *fn,const char *fbond,
81 const char *fdist, int gnx,atom_id index[],
82 real blen,real tol,bool bAver,
83 t_topology *top,int ePBC,bool bAverDist,
84 const output_env_t oenv)
86 #define MAXTAB 1000
87 FILE *out,*outd=NULL;
88 int *btab=NULL;
89 real b0=0,b1,db=0;
90 real bond,bav;
91 gmx_stats_t b_one=NULL,*b_all=NULL;
92 /*real mean, mean2, sqrdev2, sigma2;
93 int counter;*/
94 rvec *x;
95 rvec dx;
96 int status,natoms;
97 matrix box;
98 real t,fac;
99 int bind,i,nframes,i0,i1;
100 t_pbc pbc;
101 int N;
102 real aver,sigma,error;
104 if (!bAver) {
105 snew(b_all,gnx/2);
106 for(i=0; (i<gnx/2); i++)
107 b_all[i] = gmx_stats_init();
109 else {
110 b_one = gmx_stats_init();
111 snew(btab,MAXTAB+1);
114 natoms=read_first_x(oenv,&status,fn,&t,&x,box);
115 if (natoms == 0)
116 gmx_fatal(FARGS,"No atoms in trajectory!");
118 if (fdist) {
119 outd = xvgropen(fdist,bAverDist ? "Average distance" : "Distances",
120 "Time (ps)","Distance (nm)",oenv);
121 if (!bAverDist)
122 make_dist_leg(outd,gnx,index,&(top->atoms),oenv);
125 nframes=0;
126 do {
127 set_pbc(&pbc,ePBC,box);
128 if (fdist)
129 fprintf(outd," %8.4f",t);
130 nframes++; /* count frames */
131 bav = 0.0;
132 for(i=0; (i<gnx); i+=2) {
133 pbc_dx(&pbc,x[index[i]],x[index[i+1]],dx);
134 bond = norm(dx);
135 if (bAverDist)
136 bav += bond;
137 else if (fdist)
138 fprintf(outd," %.3f",bond);
139 if (bAver) {
140 gmx_stats_add_point(b_one,t,bond,0,0);
141 if (db == 0) {
142 if (blen == -1) {
143 b0 = 0;
144 b1 = 0.2;
145 db = (b1-b0)/MAXTAB;
147 else {
148 b0 = (1.0-tol)*blen;
149 b1 = (1.0+tol)*blen;
150 db = (2.0*(b1-b0))/MAXTAB;
153 bind = (int)((bond-b0)/db+0.5);
154 if ((bind >= 0) && (bind <= MAXTAB))
155 btab[bind]++;
156 else {
158 printf("bond: %4d-%4d bond=%10.5e, dx=(%10.5e,%10.5e,%10.5e)\n",
159 index[i],index[i+1],bond,dx[XX],dx[YY],dx[ZZ]);
163 else {
164 gmx_stats_add_point(b_all[i/2],t,bond,0,0);
167 if (bAverDist)
168 fprintf(outd," %.5f",bav*2.0/gnx);
169 if (fdist)
170 fprintf(outd,"\n");
171 } while (read_next_x(oenv,status,&t,natoms,x,box));
172 close_trj(status);
174 if (fdist)
175 ffclose(outd);
178 mean = mean / counter;
179 mean2 = mean2 / counter;
180 sqrdev2 = (mean2 - mean*mean);
181 sigma2 = sqrdev2*counter / (counter - 1);
183 /* For definitions see "Weet wat je meet" */
184 if (bAver) {
185 printf("\n");
186 gmx_stats_get_npoints(b_one,&N);
187 printf("Total number of samples : %d\n",N);
188 gmx_stats_get_ase(b_one,&aver,&sigma,&error);
189 printf("Mean : %g\n",aver);
190 printf("Standard deviation of the distribution: %g\n",sigma);
191 printf("Standard deviation of the mean : %g\n",error);
192 gmx_stats_done(b_one);
193 sfree(b_one);
195 out=xvgropen(fbond,"Bond Stretching Distribution",
196 "Bond Length (nm)","",oenv);
198 for(i0=0; ((i0 < MAXTAB) && (btab[i0]==0)); i0++)
200 i0=max(0,i0-1);
201 for(i1=MAXTAB; ((i1 > 0) && (btab[i1]==0)); i1--)
203 i1=min(MAXTAB,i1+1);
205 if (i0 >= i1)
206 gmx_fatal(FARGS,"No distribution... (i0 = %d, i1 = %d)? ? ! ! ? !",i0,i1);
208 fac=2.0/(nframes*gnx*db);
209 for(i=i0; (i<=i1); i++)
210 fprintf(out,"%8.5f %8.5f\n",b0+i*db,btab[i]*fac);
211 ffclose(out);
213 else {
214 fprintf(log,"%5s %5s %8s %8s\n","i","j","b_aver","sigma");
215 for(i=0; (i<gnx/2); i++) {
216 gmx_stats_get_ase(b_all[i],&aver,&sigma,NULL);
217 fprintf(log,"%5u %5u %8.5f %8.5f\n",1+index[2*i],1+index[2*i+1],
218 aver,sigma);
219 gmx_stats_done(b_all[i]);
220 sfree(b_all[i]);
222 sfree(b_all);
226 int gmx_bond(int argc,char *argv[])
228 const char *desc[] = {
229 "g_bond makes a distribution of bond lengths. If all is well a",
230 "gaussian distribution should be made when using a harmonic potential.",
231 "bonds are read from a single group in the index file in order i1-j1",
232 "i2-j2 thru in-jn.[PAR]",
233 "[TT]-tol[tt] gives the half-width of the distribution as a fraction",
234 "of the bondlength ([TT]-blen[tt]). That means, for a bond of 0.2",
235 "a tol of 0.1 gives a distribution from 0.18 to 0.22.[PAR]",
236 "Option [TT]-d[tt] plots all the distances as a function of time.",
237 "This requires a structure file for the atom and residue names in",
238 "the output. If however the option [TT]-averdist[tt] is given (as well",
239 "or separately) the average bond length is plotted instead."
241 const char *bugs[] = {
242 "It should be possible to get bond information from the topology."
244 static real blen=-1.0,tol=0.1;
245 static bool bAver=TRUE,bAverDist=TRUE;
246 t_pargs pa[] = {
247 { "-blen", FALSE, etREAL, {&blen},
248 "Bond length. By default length of first bond" },
249 { "-tol", FALSE, etREAL, {&tol},
250 "Half width of distribution as fraction of blen" },
251 { "-aver", FALSE, etBOOL, {&bAver},
252 "Average bond length distributions" },
253 { "-averdist", FALSE, etBOOL, {&bAverDist},
254 "Average distances (turns on -d)" }
256 FILE *fp;
257 char *grpname;
258 const char *fdist;
259 int gnx;
260 atom_id *index;
261 char title[STRLEN];
262 t_topology top;
263 int ePBC=-1;
264 rvec *x;
265 matrix box;
266 output_env_t oenv;
268 t_filenm fnm[] = {
269 { efTRX, "-f", NULL, ffREAD },
270 { efNDX, NULL, NULL, ffREAD },
271 { efTPS, NULL, NULL, ffOPTRD },
272 { efXVG, "-o", "bonds", ffWRITE },
273 { efLOG, NULL, "bonds", ffOPTWR },
274 { efXVG, "-d", "distance", ffOPTWR }
276 #define NFILE asize(fnm)
278 CopyRight(stderr,argv[0]);
279 parse_common_args(&argc,argv,PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE ,
280 NFILE,fnm,asize(pa),pa,asize(desc),desc,asize(bugs),bugs,
281 &oenv);
283 if (bAverDist)
284 fdist = opt2fn("-d",NFILE,fnm);
285 else {
286 fdist = opt2fn_null("-d",NFILE,fnm);
287 if (fdist)
288 read_tps_conf(ftp2fn(efTPS,NFILE,fnm),title,&top,&ePBC,&x,NULL,box,
289 FALSE);
292 rd_index(ftp2fn(efNDX,NFILE,fnm),1,&gnx,&index,&grpname);
293 if ( !even(gnx) )
294 fprintf(stderr,"WARNING: odd number of atoms (%d) in group!\n",gnx);
295 fprintf(stderr,"Will gather information on %d bonds\n",gnx/2);
297 if (!bAver)
298 fp = ftp2FILE(efLOG,NFILE,fnm,"w");
299 else
300 fp = NULL;
302 do_bonds(fp,ftp2fn(efTRX,NFILE,fnm),opt2fn("-o",NFILE,fnm),fdist,gnx,index,
303 blen,tol,bAver,&top,ePBC,bAverDist,oenv);
305 do_view(oenv,opt2fn("-o",NFILE,fnm),"-nxy");
306 do_view(oenv,opt2fn_null("-d",NFILE,fnm),"-nxy");
308 thanx(stderr);
310 return 0;