Parallel vs. sequentiual code: I get binary identical result (apart from comment...
[gromacs/qmmm-gamess-us.git] / src / tools / gmx_msd.c
blob6e28ebbcaf211a4da417d8ea32a2dce58c0a8257
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
39 #include <string.h>
40 #include <ctype.h>
41 #include <math.h>
43 #include "sysstuff.h"
44 #include "smalloc.h"
45 #include "macros.h"
46 #include "statutil.h"
47 #include "maths.h"
48 #include "futil.h"
49 #include "index.h"
50 #include "copyrite.h"
51 #include "typedefs.h"
52 #include "xvgr.h"
53 #include "gstat.h"
54 #include "gmx_statistics.h"
55 #include "tpxio.h"
56 #include "pbc.h"
57 #include "vec.h"
58 #include "confio.h"
59 #include "gmx_ana.h"
62 #define FACTOR 1000.0 /* Convert nm^2/ps to 10e-5 cm^2/s */
63 /* NORMAL = total diffusion coefficient (default). X,Y,Z is diffusion
64 coefficient in X,Y,Z direction. LATERAL is diffusion coefficient in
65 plane perpendicular to axis
67 typedef enum { NOT_USED, NORMAL, X, Y, Z, LATERAL } msd_type;
69 typedef struct {
70 real t0; /* start time and time increment between */
71 real delta_t; /* time between restart points */
72 real beginfit, /* the begin/end time for fits as reals between */
73 endfit; /* 0 and 1 */
74 real dim_factor; /* the dimensionality factor for the diffusion
75 constant */
76 real **data; /* the displacement data. First index is the group
77 number, second is frame number */
78 real *time; /* frame time */
79 real *mass; /* masses for mass-weighted msd */
80 matrix **datam;
81 rvec **x0; /* original positions */
82 rvec *com; /* center of mass correction for each frame */
83 gmx_stats_t **lsq; /* fitting stats for individual molecule msds */
84 msd_type type; /* the type of msd to calculate (lateral, etc.)*/
85 int axis; /* the axis along which to calculate */
86 int ncoords;
87 int nrestart; /* number of restart points */
88 int nmol; /* number of molecules (for bMol) */
89 int nframes; /* number of frames */
90 int nlast;
91 int ngrp; /* number of groups to use for msd calculation */
92 int *n_offs;
93 int **ndata; /* the number of msds (particles/mols) per data
94 point. */
95 } t_corr;
97 typedef real t_calc_func(t_corr *,int,atom_id[],int,rvec[],rvec,bool,matrix,
98 const output_env_t oenv);
100 static real thistime(t_corr *curr)
102 return curr->time[curr->nframes];
105 static bool in_data(t_corr *curr,int nx00)
107 return curr->nframes-curr->n_offs[nx00];
110 t_corr *init_corr(int nrgrp,int type,int axis,real dim_factor,
111 int nmol,bool bTen,bool bMass,real dt,t_topology *top,
112 real beginfit,real endfit)
114 t_corr *curr;
115 t_atoms *atoms;
116 int i;
118 snew(curr,1);
119 curr->type = (msd_type)type;
120 curr->axis = axis;
121 curr->ngrp = nrgrp;
122 curr->nrestart = 0;
123 curr->delta_t = dt;
124 curr->beginfit = (1 - 2*GMX_REAL_EPS)*beginfit;
125 curr->endfit = (1 + 2*GMX_REAL_EPS)*endfit;
126 curr->x0 = NULL;
127 curr->n_offs = NULL;
128 curr->nframes = 0;
129 curr->nlast = 0;
130 curr->dim_factor = dim_factor;
132 snew(curr->ndata,nrgrp);
133 snew(curr->data,nrgrp);
134 if (bTen)
135 snew(curr->datam,nrgrp);
136 for(i=0; (i<nrgrp); i++) {
137 curr->ndata[i] = NULL;
138 curr->data[i] = NULL;
139 if (bTen)
140 curr->datam[i] = NULL;
142 curr->time = NULL;
143 curr->lsq = NULL;
144 curr->nmol = nmol;
145 if (curr->nmol > 0) {
146 snew(curr->mass,curr->nmol);
147 for(i=0; i<curr->nmol; i++)
148 curr->mass[i] = 1;
150 else {
151 if (bMass) {
152 atoms = &top->atoms;
153 snew(curr->mass,atoms->nr);
154 for(i=0; (i<atoms->nr); i++) {
155 curr->mass[i] = atoms->atom[i].m;
160 return curr;
163 static void done_corr(t_corr *curr)
165 int i;
167 sfree(curr->n_offs);
168 for(i=0; (i<curr->nrestart); i++)
169 sfree(curr->x0[i]);
170 sfree(curr->x0);
173 static void corr_print(t_corr *curr,bool bTen,const char *fn,const char *title,
174 const char *yaxis,
175 real msdtime,real beginfit,real endfit,
176 real *DD,real *SigmaD,char *grpname[],
177 const output_env_t oenv)
179 FILE *out;
180 int i,j;
182 out=xvgropen(fn,title,output_env_get_xvgr_tlabel(oenv),yaxis,oenv);
183 if (DD) {
184 fprintf(out,"# MSD gathered over %g %s with %d restarts\n",
185 msdtime,output_env_get_time_unit(oenv),curr->nrestart);
186 fprintf(out,"# Diffusion constants fitted from time %g to %g %s\n",
187 beginfit,endfit,output_env_get_time_unit(oenv));
188 for(i=0; i<curr->ngrp; i++)
189 fprintf(out,"# D[%10s] = %.4f (+/- %.4f) (1e-5 cm^2/s)\n",
190 grpname[i],DD[i],SigmaD[i]);
192 for(i=0; i<curr->nframes; i++) {
193 fprintf(out,"%10g",output_env_conv_time(oenv,curr->time[i]));
194 for(j=0; j<curr->ngrp; j++) {
195 fprintf(out," %10g",curr->data[j][i]);
196 if (bTen) {
197 fprintf(out," %10g %10g %10g %10g %10g %10g",
198 curr->datam[j][i][XX][XX],
199 curr->datam[j][i][YY][YY],
200 curr->datam[j][i][ZZ][ZZ],
201 curr->datam[j][i][YY][XX],
202 curr->datam[j][i][ZZ][XX],
203 curr->datam[j][i][ZZ][YY]);
206 fprintf(out,"\n");
208 ffclose(out);
211 /* called from corr_loop, to do the main calculations */
212 static void calc_corr(t_corr *curr,int nr,int nx,atom_id index[],rvec xc[],
213 bool bRmCOMM,rvec com,t_calc_func *calc1,bool bTen,
214 const output_env_t oenv)
216 int nx0;
217 real g;
218 matrix mat;
219 rvec dcom;
221 /* Check for new starting point */
222 if (curr->nlast < curr->nrestart) {
223 if ((thistime(curr) >= (curr->nlast*curr->delta_t)) && (nr==0)) {
224 memcpy(curr->x0[curr->nlast],xc,curr->ncoords*sizeof(xc[0]));
225 curr->n_offs[curr->nlast]=curr->nframes;
226 copy_rvec(com,curr->com[curr->nlast]);
227 curr->nlast++;
231 /* nx0 appears to be the number of new starting points,
232 * so for all starting points, call calc1.
234 for(nx0=0; (nx0<curr->nlast); nx0++) {
235 if (bRmCOMM) {
236 rvec_sub(com,curr->com[nx0],dcom);
237 } else {
238 clear_rvec(dcom);
240 g = calc1(curr,nx,index,nx0,xc,dcom,bTen,mat,oenv);
241 #ifdef DEBUG2
242 printf("g[%d]=%g\n",nx0,g);
243 #endif
244 curr->data[nr][in_data(curr,nx0)] += g;
245 if (bTen) {
246 m_add(curr->datam[nr][in_data(curr,nx0)],mat,
247 curr->datam[nr][in_data(curr,nx0)]);
249 curr->ndata[nr][in_data(curr,nx0)]++;
253 /* the non-mass-weighted mean-squared displacement calcuation */
254 static real calc1_norm(t_corr *curr,int nx,atom_id index[],int nx0,rvec xc[],
255 rvec dcom,bool bTen,matrix mat, const output_env_t oenv)
257 int i,ix,m,m2;
258 real g,r,r2;
259 rvec rv;
261 g=0.0;
262 clear_mat(mat);
264 for(i=0; (i<nx); i++) {
265 ix=index[i];
266 r2=0.0;
267 switch (curr->type) {
268 case NORMAL:
269 for(m=0; (m<DIM); m++) {
270 rv[m] = xc[ix][m] - curr->x0[nx0][ix][m] - dcom[m];
271 r2 += rv[m]*rv[m];
272 if (bTen) {
273 for(m2=0; m2<=m; m2++)
274 mat[m][m2] += rv[m]*rv[m2];
277 break;
278 case X:
279 case Y:
280 case Z:
281 r = xc[ix][curr->type-X] - curr->x0[nx0][ix][curr->type-X] -
282 dcom[curr->type-X];
283 r2 += r*r;
284 break;
285 case LATERAL:
286 for(m=0; (m<DIM); m++) {
287 if (m != curr->axis) {
288 r = xc[ix][m] - curr->x0[nx0][ix][m] - dcom[m];
289 r2 += r*r;
292 break;
293 default:
294 gmx_fatal(FARGS,"Error: did not expect option value %d",curr->type);
296 g+=r2;
298 g/=nx;
299 msmul(mat,1.0/nx,mat);
301 return g;
304 /* calculate the com of molecules in x and put it into xa */
305 static void calc_mol_com(int nmol,int *molindex,t_block *mols,t_atoms *atoms,
306 rvec *x,rvec *xa)
308 int m,mol,i,d;
309 rvec xm;
310 real mass,mtot;
312 for(m=0; m<nmol; m++) {
313 mol = molindex[m];
314 clear_rvec(xm);
315 mtot = 0;
316 for(i=mols->index[mol]; i<mols->index[mol+1]; i++) {
317 mass = atoms->atom[i].m;
318 for(d=0; d<DIM; d++)
319 xm[d] += mass*x[i][d];
320 mtot += mass;
322 svmul(1/mtot,xm,xa[m]);
326 static real calc_one_mw(t_corr *curr,int ix,int nx0,rvec xc[],real *tm,
327 rvec dcom,bool bTen,matrix mat)
329 real r2,r,mm;
330 rvec rv;
331 int m,m2;
333 mm=curr->mass[ix];
334 if (mm == 0)
335 return 0;
336 (*tm) += mm;
337 r2 = 0.0;
338 switch (curr->type) {
339 case NORMAL:
340 for(m=0; (m<DIM); m++) {
341 rv[m] = xc[ix][m] - curr->x0[nx0][ix][m] - dcom[m];
342 r2 += mm*rv[m]*rv[m];
343 if (bTen) {
344 for(m2=0; m2<=m; m2++)
345 mat[m][m2] += mm*rv[m]*rv[m2];
348 break;
349 case X:
350 case Y:
351 case Z:
352 r = xc[ix][curr->type-X] - curr->x0[nx0][ix][curr->type-X] -
353 dcom[curr->type-X];
354 r2 = mm*r*r;
355 break;
356 case LATERAL:
357 for(m=0; (m<DIM); m++) {
358 if (m != curr->axis) {
359 r = xc[ix][m] - curr->x0[nx0][ix][m] - dcom[m];
360 r2 += mm*r*r;
363 break;
364 default:
365 gmx_fatal(FARGS,"Options got screwed. Did not expect value %d\n",curr->type);
366 } /* end switch */
367 return r2;
370 /* the normal, mass-weighted mean-squared displacement calcuation */
371 static real calc1_mw(t_corr *curr,int nx,atom_id index[],int nx0,rvec xc[],
372 rvec dcom,bool bTen,matrix mat,const output_env_t oenv)
374 int i;
375 real g,tm;
377 g=tm=0.0;
378 clear_mat(mat);
379 for(i=0; (i<nx); i++)
380 g += calc_one_mw(curr,index[i],nx0,xc,&tm,dcom,bTen,mat);
382 g/=tm;
383 if (bTen)
384 msmul(mat,1/tm,mat);
386 return g;
389 /* prepare the coordinates by removing periodic boundary crossings.
390 gnx = the number of atoms/molecules
391 index = the indices
392 xcur = the current coordinates
393 xprev = the previous coordinates
394 box = the box matrix */
395 static void prep_data(bool bMol,int gnx,atom_id index[],
396 rvec xcur[],rvec xprev[],matrix box)
398 int i,m,ind;
399 rvec hbox;
401 /* Remove periodicity */
402 for(m=0; (m<DIM); m++)
403 hbox[m]=0.5*box[m][m];
405 for(i=0; (i<gnx); i++) {
406 if (bMol)
407 ind = i;
408 else
409 ind = index[i];
411 for(m=DIM-1; m>=0; m--)
413 while(xcur[ind][m]-xprev[ind][m] <= -hbox[m])
414 rvec_inc(xcur[ind],box[m]);
415 while(xcur[ind][m]-xprev[ind][m] > hbox[m])
416 rvec_dec(xcur[ind],box[m]);
421 /* calculate the center of mass for a group
422 gnx = the number of atoms/molecules
423 index = the indices
424 xcur = the current coordinates
425 xprev = the previous coordinates
426 box = the box matrix
427 atoms = atom data (for mass)
428 com(output) = center of mass */
429 static void calc_com(bool bMol, int gnx, atom_id index[],
430 rvec xcur[],rvec xprev[],matrix box, t_atoms *atoms,
431 rvec com)
433 int i,m,ind;
434 real mass;
435 double tmass;
436 dvec sx;
438 clear_dvec(sx);
439 tmass = 0;
440 mass = 1;
442 prep_data(bMol, gnx, index, xcur, xprev, box);
443 for(i=0; (i<gnx); i++)
445 if (bMol)
446 ind = i;
447 else
448 ind = index[i];
451 mass = atoms->atom[ind].m;
452 for(m=0; m<DIM; m++)
453 sx[m] += mass*xcur[ind][m];
454 tmass += mass;
456 for(m=0; m<DIM; m++)
457 com[m] = sx[m]/tmass;
461 static real calc1_mol(t_corr *curr,int nx,atom_id index[],int nx0,rvec xc[],
462 rvec dcom,bool bTen,matrix mat, const output_env_t oenv)
464 int i;
465 real g,tm,gtot,tt;
467 tt = curr->time[in_data(curr,nx0)];
468 gtot = 0;
469 tm = 0;
470 clear_mat(mat);
471 for(i=0; (i<nx); i++) {
472 g = calc_one_mw(curr,i,nx0,xc,&tm,dcom,bTen,mat);
473 /* We don't need to normalize as the mass was set to 1 */
474 gtot += g;
475 if (tt >= curr->beginfit && (curr->endfit < 0 || tt <= curr->endfit))
476 gmx_stats_add_point(curr->lsq[nx0][i],tt,g,0,0);
478 msmul(mat,1.0/nx,mat);
480 return gtot/nx;
483 void printmol(t_corr *curr,const char *fn,
484 const char *fn_pdb,int *molindex,t_topology *top,
485 rvec *x,int ePBC,matrix box, const output_env_t oenv)
487 #define NDIST 100
488 FILE *out;
489 gmx_stats_t lsq1;
490 int i,j;
491 real a,b,D,Dav,D2av,VarD,sqrtD,sqrtD_max,scale;
492 t_pdbinfo *pdbinfo=NULL;
493 int *mol2a=NULL;
495 out=xvgropen(fn,"Diffusion Coefficients / Molecule","Molecule","D",oenv);
497 if (fn_pdb) {
498 if (top->atoms.pdbinfo == NULL)
499 snew(top->atoms.pdbinfo,top->atoms.nr);
500 pdbinfo = top->atoms.pdbinfo;
501 mol2a = top->mols.index;
504 Dav = D2av = 0;
505 sqrtD_max = 0;
506 for(i=0; (i<curr->nmol); i++) {
507 lsq1 = gmx_stats_init();
508 for(j=0; (j<curr->nrestart); j++) {
509 real xx,yy,dx,dy;
511 while(gmx_stats_get_point(curr->lsq[j][i],&xx,&yy,&dx,&dy) == estatsOK)
512 gmx_stats_add_point(lsq1,xx,yy,dx,dy);
514 gmx_stats_get_ab(lsq1,elsqWEIGHT_NONE,&a,&b,NULL,NULL,NULL,NULL);
515 gmx_stats_done(lsq1);
516 sfree(lsq1);
517 D = a*FACTOR/curr->dim_factor;
518 if (D < 0)
519 D = 0;
520 Dav += D;
521 D2av += sqr(D);
522 fprintf(out,"%10d %10g\n",i,D);
523 if (pdbinfo) {
524 sqrtD = sqrt(D);
525 if (sqrtD > sqrtD_max)
526 sqrtD_max = sqrtD;
527 for(j=mol2a[molindex[i]]; j<mol2a[molindex[i]+1]; j++)
528 pdbinfo[j].bfac = sqrtD;
531 ffclose(out);
532 do_view(oenv,fn,"-graphtype bar");
534 /* Compute variance, stddev and error */
535 Dav /= curr->nmol;
536 D2av /= curr->nmol;
537 VarD = D2av - sqr(Dav);
538 printf("<D> = %.4f Std. Dev. = %.4f Error = %.4f\n",
539 Dav,sqrt(VarD),sqrt(VarD/curr->nmol));
541 if (fn_pdb && x) {
542 scale = 1;
543 while (scale*sqrtD_max > 10)
544 scale *= 0.1;
545 while (scale*sqrtD_max < 0.1)
546 scale *= 10;
547 for(i=0; i<top->atoms.nr; i++)
548 pdbinfo[i].bfac *= scale;
549 write_sto_conf(fn_pdb,"molecular MSD",&top->atoms,x,NULL,ePBC,box);
553 /* this is the main loop for the correlation type functions
554 * fx and nx are file pointers to things like read_first_x and
555 * read_next_x
557 int corr_loop(t_corr *curr,const char *fn,t_topology *top,int ePBC,
558 bool bMol,int gnx[],atom_id *index[],
559 t_calc_func *calc1,bool bTen, int *gnx_com, atom_id *index_com[],
560 real dt, real t_pdb,rvec **x_pdb,matrix box_pdb,
561 const output_env_t oenv)
563 rvec *x[2]; /* the coordinates to read */
564 rvec *xa[2]; /* the coordinates to calculate displacements for */
565 rvec com={0};
566 real t,t_prev=0;
567 int natoms,i,j,status,cur=0,maxframes=0;
568 #define prev (1-cur)
569 matrix box;
570 bool bFirst;
572 natoms = read_first_x(oenv,&status,fn,&curr->t0,&(x[cur]),box);
573 #ifdef DEBUG
574 fprintf(stderr,"Read %d atoms for first frame\n",natoms);
575 #endif
576 if ((gnx_com!=NULL) && natoms < top->atoms.nr) {
577 fprintf(stderr,"WARNING: The trajectory only contains part of the system (%d of %d atoms) and therefore the COM motion of only this part of the system will be removed\n",natoms,top->atoms.nr);
580 snew(x[prev],natoms);
582 if (bMol) {
583 curr->ncoords = curr->nmol;
584 snew(xa[0],curr->ncoords);
585 snew(xa[1],curr->ncoords);
586 } else {
587 curr->ncoords = natoms;
588 xa[0] = x[0];
589 xa[1] = x[1];
592 bFirst = TRUE;
593 t=curr->t0;
594 if (x_pdb)
595 *x_pdb = NULL;
597 /* the loop over all frames */
600 if (x_pdb && ((bFirst && t_pdb < t) ||
601 (!bFirst &&
602 t_pdb > t - 0.5*(t - t_prev) &&
603 t_pdb < t + 0.5*(t - t_prev))))
605 if (*x_pdb == NULL)
606 snew(*x_pdb,natoms);
607 for(i=0; i<natoms; i++)
608 copy_rvec(x[cur][i],(*x_pdb)[i]);
609 copy_mat(box,box_pdb);
613 /* check whether we've reached a restart point */
614 if (bRmod(t,curr->t0,dt)) {
615 curr->nrestart++;
617 srenew(curr->x0,curr->nrestart);
618 snew(curr->x0[curr->nrestart-1],curr->ncoords);
619 srenew(curr->com,curr->nrestart);
620 srenew(curr->n_offs,curr->nrestart);
621 srenew(curr->lsq,curr->nrestart);
622 snew(curr->lsq[curr->nrestart-1],curr->nmol);
623 for(i=0;i<curr->nmol;i++)
624 curr->lsq[curr->nrestart-1][i] = gmx_stats_init();
626 if (debug)
627 fprintf(debug,"Extended data structures because of new restart %d\n",
628 curr->nrestart);
630 /* create or extend the frame-based arrays */
631 if (curr->nframes >= maxframes-1) {
632 if (maxframes==0) {
633 for(i=0; (i<curr->ngrp); i++) {
634 curr->ndata[i] = NULL;
635 curr->data[i] = NULL;
636 if (bTen)
637 curr->datam[i] = NULL;
639 curr->time = NULL;
641 maxframes+=10;
642 for(i=0; (i<curr->ngrp); i++) {
643 srenew(curr->ndata[i],maxframes);
644 srenew(curr->data[i],maxframes);
645 if (bTen)
646 srenew(curr->datam[i],maxframes);
647 for(j=maxframes-10; j<maxframes; j++) {
648 curr->ndata[i][j]=0;
649 curr->data[i][j]=0;
650 if (bTen)
651 clear_mat(curr->datam[i][j]);
654 srenew(curr->time,maxframes);
657 /* set the time */
658 curr->time[curr->nframes] = t - curr->t0;
660 /* for the first frame, the previous frame is a copy of the first frame */
661 if (bFirst) {
662 memcpy(xa[prev],xa[cur],curr->ncoords*sizeof(xa[prev][0]));
663 bFirst = FALSE;
666 /* make the molecules whole */
667 if (bMol)
668 rm_pbc(&top->idef,ePBC,natoms,box,x[cur],x[cur]);
670 /* first remove the periodic boundary condition crossings */
671 for(i=0;i<curr->ngrp;i++)
673 prep_data(bMol, gnx[i], index[i], xa[cur], xa[prev], box);
676 /* calculate the center of mass */
677 if (gnx_com)
679 prep_data(bMol, gnx_com[0], index_com[0], xa[cur], xa[prev], box);
680 calc_com(bMol, gnx_com[0], index_com[0], xa[cur], xa[prev], box,
681 &top->atoms, com);
684 /* calculate the molecules' centers of masses and put them into xa */
685 if (bMol)
686 calc_mol_com(gnx[0],index[0],&top->mols,&top->atoms, x[cur],xa[cur]);
688 /* loop over all groups in index file */
689 for(i=0; (i<curr->ngrp); i++)
691 /* calculate something useful, like mean square displacements */
692 calc_corr(curr,i,gnx[i],index[i],xa[cur], (gnx_com!=NULL),com,
693 calc1,bTen,oenv);
695 cur=prev;
696 t_prev = t;
698 curr->nframes++;
699 } while (read_next_x(oenv,status,&t,natoms,x[cur],box));
700 fprintf(stderr,"\nUsed %d restart points spaced %g %s over %g %s\n\n",
701 curr->nrestart,
702 output_env_conv_time(oenv,dt), output_env_get_time_unit(oenv),
703 output_env_conv_time(oenv,curr->time[curr->nframes-1]), output_env_get_time_unit(oenv) );
705 close_trj(status);
707 return natoms;
710 static void index_atom2mol(int *n,int *index,t_block *mols)
712 int nat,i,nmol,mol,j;
714 nat = *n;
715 i = 0;
716 nmol = 0;
717 mol = 0;
718 while (i < nat) {
719 while (index[i] > mols->index[mol]) {
720 mol++;
721 if (mol >= mols->nr)
722 gmx_fatal(FARGS,"Atom index out of range: %d",index[i]+1);
724 for(j=mols->index[mol]; j<mols->index[mol+1]; j++) {
725 if (i >= nat || index[i] != j)
726 gmx_fatal(FARGS,"The index group does not consist of whole molecules");
727 i++;
729 index[nmol++] = mol;
732 fprintf(stderr,"Split group of %d atoms into %d molecules\n",nat,nmol);
734 *n = nmol;
737 void do_corr(const char *trx_file, const char *ndx_file, const char *msd_file,
738 const char *mol_file, const char *pdb_file,real t_pdb,
739 int nrgrp, t_topology *top,int ePBC,
740 bool bTen,bool bMW,bool bRmCOMM,
741 int type,real dim_factor,int axis,
742 real dt,real beginfit,real endfit,const output_env_t oenv)
744 t_corr *msd;
745 int *gnx; /* the selected groups' sizes */
746 atom_id **index; /* selected groups' indices */
747 char **grpname;
748 int i,i0,i1,j,N,nat_trx;
749 real *DD,*SigmaD,a,a2,b,r,chi2;
750 rvec *x;
751 matrix box;
752 int *gnx_com=NULL; /* the COM removal group size */
753 atom_id **index_com=NULL; /* the COM removal group atom indices */
754 char **grpname_com=NULL; /* the COM removal group name */
756 snew(gnx,nrgrp);
757 snew(index,nrgrp);
758 snew(grpname,nrgrp);
760 fprintf(stderr, "\nSelect a group to calculate mean squared displacement for:\n");
761 get_index(&top->atoms,ndx_file,nrgrp,gnx,index,grpname);
763 if (bRmCOMM)
765 snew(gnx_com,1);
766 snew(index_com,1);
767 snew(grpname_com,1);
769 fprintf(stderr, "\nNow select a group for center of mass removal:\n");
770 get_index(&top->atoms, ndx_file, 1, gnx_com, index_com, grpname_com);
773 if (mol_file)
774 index_atom2mol(&gnx[0],index[0],&top->mols);
776 msd = init_corr(nrgrp,type,axis,dim_factor,
777 mol_file==NULL ? 0 : gnx[0],bTen,bMW,dt,top,
778 beginfit,endfit);
780 nat_trx =
781 corr_loop(msd,trx_file,top,ePBC,mol_file ? gnx[0] : 0,gnx,index,
782 (mol_file!=NULL) ? calc1_mol : (bMW ? calc1_mw : calc1_norm),
783 bTen, gnx_com, index_com, dt,t_pdb,
784 pdb_file ? &x : NULL,box,oenv);
786 /* Correct for the number of points */
787 for(j=0; (j<msd->ngrp); j++) {
788 for(i=0; (i<msd->nframes); i++) {
789 msd->data[j][i] /= msd->ndata[j][i];
790 if (bTen)
791 msmul(msd->datam[j][i],1.0/msd->ndata[j][i],msd->datam[j][i]);
795 if (mol_file) {
796 if (pdb_file && x == NULL) {
797 fprintf(stderr,"\nNo frame found need time tpdb = %g ps\n"
798 "Can not write %s\n\n",t_pdb,pdb_file);
800 i = top->atoms.nr;
801 top->atoms.nr = nat_trx;
802 printmol(msd,mol_file,pdb_file,index[0],top,x,ePBC,box,oenv);
803 top->atoms.nr = i;
806 DD = NULL;
807 SigmaD = NULL;
809 if (beginfit == -1) {
810 i0 = (int)(0.1*(msd->nframes - 1) + 0.5);
811 beginfit = msd->time[i0];
812 } else
813 for(i0=0; i0<msd->nframes && msd->time[i0]<beginfit; i0++)
816 if (endfit == -1) {
817 i1 = (int)(0.9*(msd->nframes - 1) + 0.5) + 1;
818 endfit = msd->time[i1-1];
819 } else
820 for(i1=i0; i1<msd->nframes && msd->time[i1]<=endfit; i1++)
822 fprintf(stdout,"Fitting from %g to %g %s\n\n",beginfit,endfit,
823 output_env_get_time_unit(oenv));
825 N = i1-i0;
826 if (N <= 2) {
827 fprintf(stdout,"Not enough points for fitting (%d).\n"
828 "Can not determine the diffusion constant.\n",N);
829 } else {
830 snew(DD,msd->ngrp);
831 snew(SigmaD,msd->ngrp);
832 for(j=0; j<msd->ngrp; j++) {
833 if (N >= 4) {
834 lsq_y_ax_b(N/2,&(msd->time[i0]),&(msd->data[j][i0]),&a,&b,&r,&chi2);
835 lsq_y_ax_b(N/2,&(msd->time[i0+N/2]),&(msd->data[j][i0+N/2]),&a2,&b,&r,&chi2);
836 SigmaD[j] = fabs(a-a2);
837 } else
838 SigmaD[j] = 0;
839 lsq_y_ax_b(N,&(msd->time[i0]),&(msd->data[j][i0]),&(DD[j]),&b,&r,&chi2);
840 DD[j] *= FACTOR/msd->dim_factor;
841 SigmaD[j] *= FACTOR/msd->dim_factor;
842 if (DD[j] > 0.01 && DD[j] < 1e4)
843 fprintf(stdout,"D[%10s] %.4f (+/- %.4f) 1e-5 cm^2/s\n",
844 grpname[j],DD[j],SigmaD[j]);
845 else
846 fprintf(stdout,"D[%10s] %.4g (+/- %.4g) 1e-5 cm^2/s\n",
847 grpname[j],DD[j], SigmaD[j]);
850 /* Print mean square displacement */
851 corr_print(msd,bTen,msd_file,
852 "Mean Square Displacement",
853 "MSD (nm\\S2\\N)",
854 msd->time[msd->nframes-1],beginfit,endfit,DD,SigmaD,grpname,oenv);
857 int gmx_msd(int argc,char *argv[])
859 const char *desc[] = {
860 "g_msd computes the mean square displacement (MSD) of atoms from",
861 "a set of initial positions. This provides an easy way to compute",
862 "the diffusion constant using the Einstein relation.",
863 "The time between the reference points for the MSD calculation",
864 "is set with [TT]-trestart[tt].",
865 "The diffusion constant is calculated by least squares fitting a",
866 "straight line (D*t + c) through the MSD(t) from [TT]-beginfit[tt] to",
867 "[TT]-endfit[tt] (note that t is time from the reference positions,",
868 "not simulation time). An error estimate given, which is the difference",
869 "of the diffusion coefficients obtained from fits over the two halves",
870 "of the fit interval.[PAR]",
871 "There are three, mutually exclusive, options to determine different",
872 "types of mean square displacement: [TT]-type[tt], [TT]-lateral[tt]",
873 "and [TT]-ten[tt]. Option [TT]-ten[tt] writes the full MSD tensor for",
874 "each group, the order in the output is: trace xx yy zz yx zx zy.[PAR]",
875 "If [TT]-mol[tt] is set, g_msd plots the MSD for individual molecules: ",
876 "for each individual molecule a diffusion constant is computed for ",
877 "its center of mass. The chosen index group will be split into ",
878 "molecules.[PAR]",
879 "The default way to calculate a MSD is by using mass-weighted averages.",
880 "This can be turned off with [TT]-nomw[tt].[PAR]",
881 "With the option [TT]-rmcomm[tt], the center of mass motion of a ",
882 "specific group can be removed. For trajectories produced with ",
883 "GROMACS this is usually not necessary, ",
884 "as mdrun usually already removes the center of mass motion.",
885 "When you use this option be sure that the whole system is stored",
886 "in the trajectory file.[PAR]",
887 "The diffusion coefficient is determined by linear regression of the MSD,",
888 "where, unlike for the normal output of D, the times are weighted",
889 "according to the number of reference points, i.e. short times have",
890 "a higher weight. Also when [TT]-beginfit[tt]=-1,fitting starts at 10%",
891 "and when [TT]-endfit[tt]=-1, fitting goes to 90%.",
892 "Using this option one also gets an accurate error estimate",
893 "based on the statistics between individual molecules.",
894 "Note that this diffusion coefficient and error estimate are only",
895 "accurate when the MSD is completely linear between",
896 "[TT]-beginfit[tt] and [TT]-endfit[tt].[PAR]",
897 "Option [TT]-pdb[tt] writes a pdb file with the coordinates of the frame",
898 "at time [TT]-tpdb[tt] with in the B-factor field the square root of",
899 "the diffusion coefficient of the molecule.",
900 "This option implies option [TT]-mol[tt]."
902 static const char *normtype[]= { NULL,"no","x","y","z",NULL };
903 static const char *axtitle[] = { NULL,"no","x","y","z",NULL };
904 static int ngroup = 1;
905 static real dt = 10;
906 static real t_pdb = 0;
907 static real beginfit = -1;
908 static real endfit = -1;
909 static bool bTen = FALSE;
910 static bool bMW = TRUE;
911 static bool bRmCOMM = FALSE;
912 t_pargs pa[] = {
913 { "-type", FALSE, etENUM, {normtype},
914 "Compute diffusion coefficient in one direction" },
915 { "-lateral", FALSE, etENUM, {axtitle},
916 "Calculate the lateral diffusion in a plane perpendicular to" },
917 { "-ten", FALSE, etBOOL, {&bTen},
918 "Calculate the full tensor" },
919 { "-ngroup", FALSE, etINT, {&ngroup},
920 "Number of groups to calculate MSD for" },
921 { "-mw", FALSE, etBOOL, {&bMW},
922 "Mass weighted MSD" },
923 { "-rmcomm", FALSE, etBOOL, {&bRmCOMM},
924 "Remove center of mass motion" },
925 { "-tpdb",FALSE, etTIME, {&t_pdb},
926 "The frame to use for option -pdb (%t)" },
927 { "-trestart",FALSE, etTIME, {&dt},
928 "Time between restarting points in trajectory (%t)" },
929 { "-beginfit",FALSE, etTIME, {&beginfit},
930 "Start time for fitting the MSD (%t), -1 is 10%" },
931 { "-endfit",FALSE, etTIME, {&endfit},
932 "End time for fitting the MSD (%t), -1 is 90%" }
935 t_filenm fnm[] = {
936 { efTRX, NULL, NULL, ffREAD },
937 { efTPS, NULL, NULL, ffREAD },
938 { efNDX, NULL, NULL, ffOPTRD },
939 { efXVG, NULL, "msd", ffWRITE },
940 { efXVG, "-mol", "diff_mol",ffOPTWR },
941 { efPDB, "-pdb", "diff_mol", ffOPTWR }
943 #define NFILE asize(fnm)
945 t_topology top;
946 int ePBC;
947 matrix box;
948 char title[256];
949 const char *trx_file, *tps_file, *ndx_file, *msd_file, *mol_file, *pdb_file;
950 rvec *xdum;
951 bool bTop;
952 int axis,type;
953 real dim_factor;
954 output_env_t oenv;
956 CopyRight(stderr,argv[0]);
958 parse_common_args(&argc,argv,
959 PCA_CAN_VIEW | PCA_CAN_TIME | PCA_TIME_UNIT | PCA_BE_NICE,
960 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
961 trx_file = ftp2fn_null(efTRX,NFILE,fnm);
962 tps_file = ftp2fn_null(efTPS,NFILE,fnm);
963 ndx_file = ftp2fn_null(efNDX,NFILE,fnm);
964 msd_file = ftp2fn_null(efXVG,NFILE,fnm);
965 pdb_file = opt2fn_null("-pdb",NFILE,fnm);
966 if (pdb_file)
967 mol_file = opt2fn("-mol",NFILE,fnm);
968 else
969 mol_file = opt2fn_null("-mol",NFILE,fnm);
971 if (ngroup < 1)
972 gmx_fatal(FARGS,"Must have at least 1 group (now %d)",ngroup);
973 if (mol_file && ngroup > 1)
974 gmx_fatal(FARGS,"With molecular msd can only have 1 group (now %d)",
975 ngroup);
978 if (mol_file) {
979 bMW = TRUE;
980 fprintf(stderr,"Calculating diffusion coefficients for molecules.\n");
983 if (normtype[0][0]!='n') {
984 type = normtype[0][0] - 'x' + X; /* See defines above */
985 dim_factor = 2.0;
987 else {
988 type = NORMAL;
989 dim_factor = 6.0;
991 if ((type==NORMAL) && (axtitle[0][0]!='n')) {
992 type=LATERAL;
993 dim_factor = 4.0;
994 axis = (axtitle[0][0] - 'x'); /* See defines above */
996 else
997 axis = 0;
999 if (bTen && type != NORMAL)
1000 gmx_fatal(FARGS,"Can only calculate the full tensor for 3D msd");
1002 bTop = read_tps_conf(tps_file,title,&top,&ePBC,&xdum,NULL,box,bMW||bRmCOMM);
1003 if (mol_file && !bTop)
1004 gmx_fatal(FARGS,
1005 "Could not read a topology from %s. Try a tpr file instead.",
1006 tps_file);
1008 do_corr(trx_file,ndx_file,msd_file,mol_file,pdb_file,t_pdb,ngroup,
1009 &top,ePBC,bTen,bMW,bRmCOMM,type,dim_factor,axis,dt,beginfit,endfit,
1010 oenv);
1012 view_all(oenv,NFILE, fnm);
1014 thanx(stderr);
1016 return 0;