Removing scripts, introducing wrapper binaries instead.
[gromacs.git] / src / tools / g_msd.c
blob02c1016ae95908f8953b4bcb148063e17f1710d5
1 /*
2 * $Id$
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * Green Red Orange Magenta Azure Cyan Skyblue
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #include <string.h>
41 #include <ctype.h>
42 #include <math.h>
44 #include "sysstuff.h"
45 #include "smalloc.h"
46 #include "macros.h"
47 #include "statutil.h"
48 #include "maths.h"
49 #include "futil.h"
50 #include "index.h"
51 #include "copyrite.h"
52 #include "typedefs.h"
53 #include "xvgr.h"
54 #include "gstat.h"
55 #include "tpxio.h"
56 #include "pbc.h"
57 #include "vec.h"
59 #define FACTOR 1000.0 /* Convert nm^2/ps to 10e-5 cm^2/s */
60 /* NORMAL = total diffusion coefficient (default). X,Y,Z is diffusion
61 coefficient in X,Y,Z direction. LATERAL is diffusion coefficient in
62 plane perpendicular to axis
64 enum { NOT_USED, NORMAL, X, Y, Z, LATERAL };
66 typedef struct {
67 real t0,delta_t,dim_factor;
68 real **data,*time,*data_x,*data_y,*data_z,*data_xy,*mass;
69 rvec **x0;
70 t_block *mols;
71 t_lsq **lsq;
72 int type,axis,natoms,nrestart,nnx,nframes,nlast,ngrp;
73 int *n_offs;
74 int **ndata;
75 } t_corr;
77 typedef real t_calc_func(t_corr *,int,atom_id[],int,rvec[]);
78 typedef void t_prep_data_func(t_corr *this,int gnx,atom_id index[],
79 rvec xcur[],rvec xprev[],matrix box);
81 static real thistime(t_corr *this)
83 return this->time[this->nframes];
86 static bool in_data(t_corr *this,int nx00)
88 return this->nframes-this->n_offs[nx00];
91 t_corr *init_corr(int nrgrp,int type,int axis,real dim_factor,
92 bool bMass,bool bMol,real dt,t_topology *top)
94 t_corr *this;
95 t_atoms *atoms;
96 int i;
98 snew(this,1);
99 this->type = type;
100 this->axis = axis;
101 this->ngrp = nrgrp;
102 this->nrestart = 0;
103 this->delta_t = dt;
104 this->x0 = NULL;
105 this->n_offs = NULL;
106 this->nframes = 0;
107 this->nlast = 0;
108 this->dim_factor = dim_factor;
110 snew(this->ndata,nrgrp);
111 snew(this->data,nrgrp);
112 for(i=0; (i<nrgrp); i++) {
113 this->ndata[i] = NULL;
114 this->data[i] = NULL;
116 this->time = NULL;
117 this->lsq = NULL;
118 if (bMass) {
119 atoms=&top->atoms;
120 snew(this->mass,atoms->nr);
121 for(i=0; (i<atoms->nr); i++) {
122 this->mass[i]=atoms->atom[i].m;
125 if (bMol)
126 this->mols = &(top->blocks[ebMOLS]);
128 return this;
131 static void done_corr(t_corr *this)
133 int i;
135 sfree(this->n_offs);
136 for(i=0; (i<this->nrestart); i++)
137 sfree(this->x0[i]);
138 sfree(this->x0);
141 static void init_restart(t_corr *this)
143 int i;
145 this->nrestart++;
147 srenew(this->x0,this->nrestart);
148 snew(this->x0[this->nrestart-1],this->natoms);
149 srenew(this->n_offs,this->nrestart);
152 static void corr_print(t_corr *this,char *fn,char *title,char *yaxis,
153 real beginfit,real endfit,
154 real *DD,real *SigmaD,char *grpname[])
156 FILE *out;
157 int i,j;
159 out=xvgropen(fn,title,xvgr_tlabel(),yaxis);
160 if (DD) {
161 fprintf(out,"# Diffusion constants fitted from time %g to %g (%s)\n",
162 beginfit,endfit,time_unit());
163 for(i=0; i<this->ngrp; i++)
164 fprintf(out,"# D[%10s] = %.3f (+/- %.3f) (1e-5 cm^2/s)\n",
165 grpname[i],DD[i],SigmaD[i]);
167 for(i=0; i<this->nframes; i++) {
168 fprintf(out,"%10g",convert_time(this->time[i]));
169 for(j=0; j<this->ngrp; j++)
170 fprintf(out," %10g",this->data[j][i]);
171 fprintf(out,"\n");
173 fclose(out);
176 /* called from corr_loop, to do the main calculations */
177 static void calc_corr(t_corr *this,int nr,int nx,atom_id index[],rvec xc[],
178 t_calc_func *calc1)
180 int nx0;
181 real g;
183 /* Check for new starting point */
184 if (this->nlast < this->nrestart) {
185 if ((thistime(this) >= (this->nlast*this->delta_t)) && (nr==0)) {
186 memcpy(this->x0[this->nlast],xc,this->natoms*sizeof(xc[0]));
187 this->n_offs[this->nlast]=this->nframes;
188 this->nlast++;
192 /* nx0 appears to be the number of new starting points,
193 * so for all starting points, call calc1.
195 for(nx0=0; (nx0<this->nlast); nx0++) {
196 g = calc1(this,nx,index,nx0,xc);
197 #ifdef DEBUG2
198 printf("g[%d]=%g\n",nx0,g);
199 #endif
200 this->data [nr][in_data(this,nx0)]+=g;
201 this->ndata[nr][in_data(this,nx0)]++;
205 static real calc1_norm(t_corr *this,int nx,atom_id index[],int nx0,rvec xc[])
207 int i,ix,m;
208 real g,r,r2;
210 g=0.0;
212 for(i=0; (i<nx); i++) {
213 ix=index[i];
214 r2=0.0;
215 switch (this->type) {
216 case NORMAL:
217 for(m=0; (m<DIM); m++) {
218 r = this->x0[nx0][ix][m]-xc[ix][m];
219 r2 += r*r;
221 break;
222 case X:
223 case Y:
224 case Z:
225 r = this->x0[nx0][ix][this->type-X]-xc[ix][this->type-X];
226 break;
227 case LATERAL:
228 for(m=0; (m<DIM); m++) {
229 if (m != this->axis) {
230 r = this->x0[nx0][ix][m]-xc[ix][m];
231 r2 += r*r;
234 break;
235 default:
236 fatal_error(0,"Error: did not expect option value %d",this->type);
238 g+=r2;
240 g/=nx;
242 return g;
245 /* this is the mean loop for the correlation type functions
246 * fx and nx are file pointers to things like read_first_x and
247 * read_next_x
249 void corr_loop(t_corr *this,char *fn,int gnx[],atom_id *index[],
250 t_calc_func *calc1,t_prep_data_func *prep1,real dt)
252 rvec *x[2];
253 real t;
254 int i,j,status,cur=0,maxframes=0;
255 #define prev (1-cur)
256 matrix box;
258 this->natoms=read_first_x(&status,fn,&this->t0,&(x[prev]),box);
259 #ifdef DEBUG
260 fprintf(stderr,"Read %d atoms for first frame\n",this->natoms);
261 #endif
263 snew(x[cur],this->natoms);
265 /* init_restart(this,nrestart,dt); */
266 memcpy(x[cur],x[prev],this->natoms*sizeof(x[prev][0]));
267 t=this->t0;
268 do {
269 if (bRmod(t,this->t0,dt))
270 init_restart(this);
271 if (this->nframes >= maxframes-1) {
272 if (maxframes==0) {
273 for(i=0; (i<this->ngrp); i++) {
274 this->ndata[i] = NULL;
275 this->data[i] = NULL;
277 this->time = NULL;
279 maxframes+=10;
280 for(i=0; (i<this->ngrp); i++) {
281 srenew(this->ndata[i],maxframes);
282 srenew(this->data[i],maxframes);
283 for(j=maxframes-10; j<maxframes; j++) {
284 this->ndata[i][j]=0;
285 this->data[i][j]=0;
288 srenew(this->time,maxframes);
291 this->time[this->nframes] = t - this->t0;
293 /* loop over all groups in index file */
294 for(i=0; (i<this->ngrp); i++) {
295 /* nice for putting things in boxes and such */
296 prep1(this,gnx[i],index[i],x[cur],x[prev],box);
297 /* calculate something useful, like mean square displacements */
298 calc_corr(this,i,gnx[i],index[i],x[cur],calc1);
300 cur=prev;
302 this->nframes++;
303 } while (read_next_x(status,&t,this->natoms,x[cur],box));
304 fprintf(stderr,"\nUsed %d restart points spaced %g %s over %g %s\n\n",
305 this->nrestart,
306 convert_time(dt), time_unit(),
307 convert_time(this->time[this->nframes-1]), time_unit() );
309 close_trj(status);
312 static void prep_data_mol(t_corr *this,int gnx,atom_id index[],
313 rvec xcur[],rvec xprev[],matrix box)
315 int i,j,k,m,d,ind;
316 rvec hbox;
318 /* Remove periodicity */
319 for(m=0; (m<DIM); m++)
320 hbox[m]=0.5*box[m][m];
321 for(i=0; i<gnx; i++) {
322 ind=index[i];
323 for(j=this->mols->index[ind]; j<this->mols->index[ind+1]; j++) {
324 k=this->mols->a[j];
325 for(m=DIM-1; m>=0; m--) {
326 while (xcur[k][m]-xprev[k][m] <= -hbox[m])
327 for(d=0; d<=m; d++)
328 xcur[k][d] += box[m][d];
329 while (xcur[k][m]-xprev[k][m] > hbox[m])
330 for(d=0; d<=m; d++)
331 xcur[k][d] -= box[m][d];
337 static real calc_one_mw(t_corr *this,int ix,int nx0,rvec xc[],real *tm)
339 real r2,r,mm;
340 int m;
342 mm=this->mass[ix];
343 if (mm < 1)
344 return 0;
345 (*tm)+=this->mass[ix];
346 r2=0.0;
347 switch (this->type) {
348 case NORMAL:
349 for(m=0; (m<DIM); m++) {
350 r = this->x0[nx0][ix][m]-xc[ix][m];
351 r2 += this->mass[ix]*r*r;
353 break;
354 case X:
355 case Y:
356 case Z:
357 r = this->x0[nx0][ix][this->type-X]-xc[ix][this->type-X];
358 r2 = this->mass[ix]*r*r;
359 break;
360 case LATERAL:
361 for(m=0; (m<DIM); m++) {
362 if (m != this->axis) {
363 r = this->x0[nx0][ix][m]-xc[ix][m];
364 r2 += this->mass[ix]*r*r;
367 break;
368 default:
369 fatal_error(0,"Options got screwed. Did not expect value %d\n",this->type);
370 } /* end switch */
371 return r2;
374 static real calc1_mw(t_corr *this,int nx,atom_id index[],int nx0,rvec xc[])
376 int i;
377 real g,tm;
379 g=tm=0.0;
380 for(i=0; (i<nx); i++)
381 g+=calc_one_mw(this,index[i],nx0,xc,&tm);
383 g/=tm;
385 return g;
388 static void prep_data_norm(t_corr *this,int gnx,atom_id index[],
389 rvec xcur[],rvec xprev[],matrix box)
391 int i,m,ind;
392 rvec hbox;
394 /* Remove periodicity */
395 for(m=0; (m<DIM); m++)
396 hbox[m]=0.5*box[m][m];
397 for(i=0; (i<gnx); i++) {
398 ind=index[i];
399 for(m=DIM-1; m>=0; m--) {
400 while(xcur[ind][m]-xprev[ind][m] <= hbox[m])
401 rvec_inc(xcur[ind],box[m]);
402 while(xcur[ind][m]-xprev[ind][m] > hbox[m])
403 rvec_dec(xcur[ind],box[m]);
408 static real calc1_mol(t_corr *this,int nx,atom_id index[],int nx0,rvec xc[])
410 int i,ii,j;
411 real g,mm,gtot,tt;
413 if (this->lsq == NULL) {
414 this->nnx = nx;
415 snew(this->lsq,this->nrestart);
416 for(i=0; (i<this->nrestart); i++)
417 snew(this->lsq[i],nx);
419 tt=this->time[in_data(this,nx0)];
420 gtot=0;
421 for(i=0; (i<nx); i++) {
422 ii=index[i];
423 g=mm=0;
424 for(j=this->mols->index[ii]; (j<this->mols->index[ii+1]); j++)
425 g += calc_one_mw(this,this->mols->a[j],nx0,xc,&mm);
427 g/=mm;
428 gtot+=g;
429 add_lsq(&(this->lsq[nx0][i]),tt,g);
431 return gtot/nx;
434 void printmol(t_corr *this,char *fn)
436 #define NDIST 100
437 FILE *out;
438 t_lsq lsq1;
439 int i,j;
440 real a,b,*D,Dav,D2av,VarD;
442 out=xvgropen(fn,"Diffusion Coefficients / Molecule","Molecule","D");
444 snew(D,this->nnx);
445 Dav = D2av = 0;
446 for(i=0; (i<this->nnx); i++) {
447 init_lsq(&lsq1);
448 for(j=0; (j<this->nrestart); j++) {
449 lsq1.sx+=this->lsq[j][i].sx;
450 lsq1.sy+=this->lsq[j][i].sy;
451 lsq1.xx+=this->lsq[j][i].xx;
452 lsq1.yx+=this->lsq[j][i].yx;
453 lsq1.np+=this->lsq[j][i].np;
455 get_lsq_ab(&lsq1,&a,&b);
456 D[i] = a*FACTOR/this->dim_factor;
457 Dav += D[i];
458 D2av += sqr(D[i]);
459 fprintf(out,"%10d %10g\n",i,D[i]);
461 fclose(out);
462 do_view(fn,"-graphtype bar");
464 /* Compute variance, stddev and error */
465 Dav /= this->nnx;
466 D2av /= this->nnx;
467 VarD = D2av - sqr(Dav);
468 printf("<D> = %.3f Std. Dev. = %.3f Error = %.3f\n",
469 Dav,sqrt(VarD),sqrt(VarD/this->nnx));
471 sfree(D);
474 void do_corr(char *trx_file, char *ndx_file, char *msd_file, char *mol_file,
475 int nrgrp, t_topology *top,bool bMW,
476 int type,real dim_factor,int axis,
477 real dt,real beginfit,real endfit)
479 t_corr *msd;
480 int *gnx;
481 atom_id **index;
482 char **grpname;
483 int i,i0,i1,j,N;
484 real *DD,*SigmaD,a,a2,b;
486 snew(gnx,nrgrp);
487 snew(index,nrgrp);
488 snew(grpname,nrgrp);
490 if (mol_file && !ndx_file) {
491 gnx[0] = top->blocks[ebMOLS].nr;
492 snew(index[0],gnx[0]);
493 grpname[0] = "Molecules";
494 for(i=0; (i<gnx[0]); i++)
495 index[0][i] = i;
497 else
498 get_index(&top->atoms,ndx_file,nrgrp,gnx,index,grpname);
500 msd = init_corr(nrgrp,type,axis,dim_factor,bMW,(mol_file!=NULL),dt,top);
502 corr_loop(msd,trx_file,gnx,index,
503 (mol_file!=NULL) ? calc1_mol : (bMW ? calc1_mw : calc1_norm),
504 (mol_file!=NULL) ? prep_data_mol : prep_data_norm,dt);
506 /* Correct for the number of points */
507 for(j=0; (j<msd->ngrp); j++)
508 for(i=0; (i<msd->nframes); i++)
509 msd->data[j][i] /= msd->ndata[j][i];
511 if (mol_file)
512 printmol(msd,mol_file);
514 DD = NULL;
515 SigmaD = NULL;
516 for(i0=0; i0<msd->nframes && msd->time[i0]<beginfit; i0++)
518 if (endfit == -1) {
519 i1 = msd->nframes;
520 endfit = msd->time[i1-1];
521 } else
522 for(i1=i0; i1<msd->nframes && msd->time[i1]<=endfit; i1++)
524 N = i1-i0;
525 if (N <= 2) {
526 fprintf(stdout,"Not enough points for fitting (%d).\n"
527 "Can not determine the diffusion constant.\n",N);
528 } else {
529 snew(DD,msd->ngrp);
530 snew(SigmaD,msd->ngrp);
531 for(j=0; j<msd->ngrp; j++) {
532 if (N >= 4) {
533 lsq_y_ax_b(N/2,&(msd->time[i0]),&(msd->data[j][i0]),&a,&b);
534 lsq_y_ax_b(N/2,&(msd->time[i0+N/2]),&(msd->data[j][i0+N/2]),&a2,&b);
535 SigmaD[j] = fabs(a-a2);
536 } else
537 SigmaD[j] = 0;
538 lsq_y_ax_b(N,&(msd->time[i0]),&(msd->data[j][i0]),&(DD[j]),&b);
539 DD[j] *= FACTOR/msd->dim_factor;
540 SigmaD[j] *= FACTOR/msd->dim_factor;
541 fprintf(stdout,"D[%10s] %.3f (+/- %.3f) 1e-5 cm^2/s\n",
542 grpname[j],DD[j],SigmaD[j]);
545 /* Print mean square displacement */
546 corr_print(msd,msd_file,
547 "Mean Square Displacement",
548 "MSD (nm\\S2\\N)",
549 beginfit,endfit,DD,SigmaD,grpname);
552 int gmx_msd(int argc,char *argv[])
554 static char *desc[] = {
555 "g_msd computes the mean square displacement (MSD) of atoms from",
556 "their initial positions. This provides an easy way to compute",
557 "the diffusion constant using the Einstein relation.",
558 "The diffusion constant is calculated by least squares fitting a",
559 "straight line through the MSD from [TT]-beginfit[tt] to",
560 "[TT]-endfit[tt]. An error estimate given, which is the difference",
561 "of the diffusion coefficients obtained from fits over the two halfs",
562 "of the fit interval.[PAR]",
563 "Option [TT]-mol[tt] plots the MSD for molecules, this implies",
564 "[TT]-mw[tt], i.e. for each inidividual molecule an diffusion constant",
565 "is computed. When using an index file, it should contain molecule",
566 "numbers instead of atom numbers.",
567 "Using this option one also gets an accurate error estimate",
568 "based on the statistics between individual molecules. Since one usually",
569 "is interested in self-diffusion at infinite dilution this is probably",
570 "the most useful number.[PAR]",
572 static char *normtype[]= { NULL,"no","x","y","z",NULL };
573 static char *axtitle[] = { NULL,"no","x","y","z",NULL };
574 static int ngroup = 1;
575 static real dt = 0;
576 static real beginfit = 0;
577 static real endfit = -1;
578 static bool bMW = TRUE;
579 t_pargs pa[] = {
580 { "-type", FALSE, etENUM, {normtype},
581 "Compute diffusion coefficient in one direction" },
582 { "-lateral", FALSE, etENUM, {axtitle},
583 "Calculate the lateral diffusion in a plane perpendicular to" },
584 { "-ngroup", FALSE, etINT, {&ngroup},
585 "Number of groups to calculate MSD for" },
586 { "-mw", FALSE, etBOOL, {&bMW},
587 "Mass weighted MSD" },
588 { "-trestart",FALSE, etTIME, {&dt},
589 "Time between restarting points in trajectory (%t)" },
590 { "-beginfit",FALSE, etTIME, {&beginfit},
591 "Start time for fitting the MSD (%t)" },
592 { "-endfit",FALSE, etTIME, {&endfit},
593 "End time for fitting the MSD (%t), -1 is till end" }
596 t_filenm fnm[] = {
597 { efTRX, NULL, NULL, ffREAD },
598 { efTPS, NULL, NULL, ffREAD },
599 { efNDX, NULL, NULL, ffOPTRD },
600 { efXVG, NULL, "msd", ffWRITE },
601 { efXVG, "-mol", "diff_mol",ffOPTWR }
603 #define NFILE asize(fnm)
605 t_topology top;
606 matrix box;
607 char title[256];
608 char *trx_file, *tps_file, *ndx_file, *msd_file, *mol_file;
609 rvec *xdum;
610 bool bTop;
611 int axis,type;
612 real dim_factor;
614 CopyRight(stderr,argv[0]);
616 parse_common_args(&argc,argv,PCA_CAN_VIEW | PCA_CAN_TIME | PCA_TIME_UNIT | PCA_BE_NICE,
617 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL);
618 trx_file = ftp2fn_null(efTRX,NFILE,fnm);
619 tps_file = ftp2fn_null(efTPS,NFILE,fnm);
620 ndx_file = ftp2fn_null(efNDX,NFILE,fnm);
621 msd_file = ftp2fn_null(efXVG,NFILE,fnm);
622 mol_file = opt2fn_null("-mol",NFILE,fnm);
624 if (ngroup < 1)
625 fatal_error(0,"Must have at least 1 group (now %d)",ngroup);
627 if (mol_file) {
628 bMW = TRUE;
629 fprintf(stderr,"Calculating diffusion coefficients for molecules.\n");
632 if (normtype[0][0]!='n') {
633 type = normtype[0][0] - 'x' + X; /* See defines above */
634 dim_factor = 2.0;
636 else {
637 type = NORMAL;
638 dim_factor = 6.0;
640 if ((type==NORMAL) && (axtitle[0][0]!='n')) {
641 type=LATERAL;
642 dim_factor = 4.0;
643 axis = (axtitle[0][0] - 'x'); /* See defines above */
645 else
646 axis = 0;
647 bTop=read_tps_conf(tps_file,title,&top,&xdum,NULL,box,bMW);
648 if (mol_file && !bTop)
649 fatal_error(0,"Could not read a topology from %s. Try a tpr file instead.",
650 tps_file);
652 do_corr(trx_file,ndx_file,msd_file,mol_file,ngroup,
653 &top,bMW,type,dim_factor,axis,dt,beginfit,endfit);
655 view_all(NFILE, fnm);
657 thanx(stderr);
659 return 0;