1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2 * $Id: densorder.c,v 0.9
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
12 * Copyright (c) 1991-2001
13 * BIOSON Research Institute, Dept. of Biophysical Chemistry
14 * University of Groningen, The Netherlands
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 * Do check out http://www.gromacs.org , or mail us at gromacs@gromacs.org .
34 * Gyas ROwers Mature At Cryogenic Speed
60 #include "binsearch.h"
61 #include "powerspect.h"
64 /* Print name of first atom in all groups in index file */
65 static void print_types(atom_id index[], atom_id a[], int ngrps,
66 char *groups[], t_topology *top)
70 fprintf(stderr,"Using following groups: \n");
71 for(i = 0; i < ngrps; i++)
72 fprintf(stderr,"Groupname: %s First atomname: %s First atomnr %u\n",
73 groups[i], *(top->atoms.atomname[a[index[i]]]), a[index[i]]);
77 static void check_length(real length, int a, int b)
80 fprintf(stderr,"WARNING: distance between atoms %d and "
81 "%d > 0.3 nm (%f). Index file might be corrupt.\n",
85 static void find_tetra_order_grid(t_topology top, int ePBC,
86 int natoms, matrix box,
87 rvec x[],int maxidx,atom_id index[],
88 real time, real *sgmean, real *skmean,
89 int nslicex, int nslicey, int nslicez,
90 real ***sggrid,real ***skgrid)
92 int ix,jx,i,j,k,l,n,*nn[4];
93 rvec dx,rj,rk,urk,urj;
94 real cost,cost2,*sgmol,*skmol,rmean,rmean2,r2,box2,*r_nn[4];
96 int slindex_x,slindex_y,slindex_z;
98 real onethird=1.0/3.0;
101 /* dmat = init_mat(maxidx, FALSE); */
103 box2 = box[XX][XX] * box[XX][XX];
105 /* Initialize expanded sl_count array */
106 snew(sl_count,nslicex);
107 for(i=0;i<nslicex;i++)
109 snew(sl_count[i],nslicey);
110 for(j=0;j<nslicey;j++)
112 snew(sl_count[i][j],nslicez);
117 for (i=0; (i<4); i++)
119 snew(r_nn[i],natoms);
122 for (j=0; (j<natoms); j++)
131 /* Must init pbc every step because of pressure coupling */
132 gpbc = gmx_rmpbc_init(&top.idef,ePBC,natoms,box);
133 gmx_rmpbc(gpbc,natoms,box,x);
138 for (i=0; (i<maxidx); i++)
139 { /* loop over index file */
141 for (j=0; (j<maxidx); j++)
144 if (i == j) continue;
148 pbc_dx(&pbc,x[ix],x[jx],dx);
151 /* set_mat_entry(dmat,i,j,r2); */
153 /* determine the nearest neighbours */
156 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
157 r_nn[2][i] = r_nn[1][i]; nn[2][i] = nn[1][i];
158 r_nn[1][i] = r_nn[0][i]; nn[1][i] = nn[0][i];
159 r_nn[0][i] = r2; nn[0][i] = j;
160 } else if (r2 < r_nn[1][i])
162 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
163 r_nn[2][i] = r_nn[1][i]; nn[2][i] = nn[1][i];
164 r_nn[1][i] = r2; nn[1][i] = j;
165 } else if (r2 < r_nn[2][i])
167 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
168 r_nn[2][i] = r2; nn[2][i] = j;
169 } else if (r2 < r_nn[3][i])
171 r_nn[3][i] = r2; nn[3][i] = j;
176 /* calculate mean distance between nearest neighbours */
178 for (j=0; (j<4); j++)
180 r_nn[j][i] = sqrt(r_nn[j][i]);
189 /* Chau1998a eqn 3 */
190 /* angular part tetrahedrality order parameter per atom */
191 for (j=0; (j<3); j++)
193 for (k=j+1; (k<4); k++)
195 pbc_dx(&pbc,x[ix],x[index[nn[k][i]]],rk);
196 pbc_dx(&pbc,x[ix],x[index[nn[j][i]]],rj);
201 cost = iprod(urk,urj) + onethird;
209 /* normalize sgmol between 0.0 and 1.0 */
210 sgmol[i] = 3*sgmol[i]/32;
213 /* distance part tetrahedrality order parameter per atom */
214 rmean2 = 4 * 3 * rmean * rmean;
215 for (j=0; (j<4); j++)
217 skmol[i] += (rmean - r_nn[j][i]) * (rmean - r_nn[j][i]) / rmean2;
218 /* printf("%d %f (%f %f %f %f) \n",
219 i, skmol[i], rmean, rmean2, r_nn[j][i], (rmean - r_nn[j][i]) );
225 /* Compute sliced stuff in x y z*/
226 slindex_x = gmx_nint((1+x[i][XX]/box[XX][XX])*nslicex) % nslicex;
227 slindex_y = gmx_nint((1+x[i][YY]/box[YY][YY])*nslicey) % nslicey;
228 slindex_z = gmx_nint((1+x[i][ZZ]/box[ZZ][ZZ])*nslicez) % nslicez;
229 sggrid[slindex_x][slindex_y][slindex_z] += sgmol[i];
230 skgrid[slindex_x][slindex_y][slindex_z] += skmol[i];
231 (sl_count[slindex_x][slindex_y][slindex_z])++;
232 } /* loop over entries in index file */
237 for(i=0; (i<nslicex); i++)
239 for(j=0; j<nslicey; j++)
241 for(k=0;k<nslicez;k++)
243 if (sl_count[i][j][k] > 0)
245 sggrid[i][j][k] /= sl_count[i][j][k];
246 skgrid[i][j][k] /= sl_count[i][j][k];
255 for (i=0; (i<4); i++)
262 /*Determines interface from tetrahedral order parameter in box with specified binwidth. */
263 /*Outputs interface positions(bins), the number of timeframes, and the number of surface-mesh points in xy*/
265 static void calc_tetra_order_interface(const char *fnNDX,const char *fnTPS,const char *fnTRX, real binw, int tblock,
266 int *nframes, int *nslicex, int *nslicey,
267 real sgang1,real sgang2,real ****intfpos,
270 FILE *fpsg=NULL,*fpsk=NULL;
271 char *sgslfn="sg_ang_mesh"; /* Hardcoded filenames for debugging*/
272 char *skslfn="sk_dist_mesh";
275 char title[STRLEN],subtitle[STRLEN];
281 real sg,sk, sgintf, pos;
282 atom_id **index=NULL;
284 int i,j,k,n,*isize,ng, nslicez, framenr;
285 real ***sg_grid=NULL,***sk_grid=NULL, ***sg_fravg=NULL, ***sk_fravg=NULL, ****sk_4d=NULL, ****sg_4d=NULL;
289 const real onehalf=1.0/2.0;
290 /* real ***intfpos[2]; pointers to arrays of two interface positions zcoord(framenr,xbin,ybin): intfpos[interface_index][t][nslicey*x+y]
291 * i.e 1D Row-major order in (t,x,y) */
294 read_tps_conf(fnTPS,title,&top,&ePBC,&xtop,NULL,box,FALSE);
296 *nslicex= (int)(box[XX][XX]/binw + onehalf); /*Calculate slicenr from binwidth*/
297 *nslicey= (int)(box[YY][YY]/binw + onehalf);
298 nslicez= (int)(box[ZZ][ZZ]/binw + onehalf);
303 /* get index groups */
304 printf("Select the group that contains the atoms you want to use for the tetrahedrality order parameter calculation:\n");
308 get_index(&top.atoms,fnNDX,ng,isize,index,grpname);
310 /* Analyze trajectory */
311 natoms=read_first_x(oenv,&status,fnTRX,&t,&x,box);
312 if ( natoms > top.atoms.nr )
313 gmx_fatal(FARGS,"Topology (%d atoms) does not match trajectory (%d atoms)",
314 top.atoms.nr,natoms);
315 check_index(NULL,ng,index[0],NULL,natoms);
318 /*Prepare structures for temporary storage of frame info*/
319 snew(sg_grid,*nslicex);
320 snew(sk_grid,*nslicex);
321 for(i=0;i<*nslicex;i++)
323 snew(sg_grid[i],*nslicey);
324 snew(sk_grid[i],*nslicey);
325 for(j=0;j<*nslicey;j++)
327 snew(sg_grid[i][j],nslicez);
328 snew(sk_grid[i][j],nslicez);
337 /* Loop over frames*/
340 /*Initialize box meshes (temporary storage for each tblock frame -reinitialise every tblock steps */
341 if(framenr%tblock ==0)
343 srenew(sk_4d,*nframes+1);
344 srenew(sg_4d,*nframes+1);
345 snew(sg_fravg,*nslicex);
346 snew(sk_fravg,*nslicex);
347 for(i=0;i<*nslicex;i++)
349 snew(sg_fravg[i],*nslicey);
350 snew(sk_fravg[i],*nslicey);
351 for(j=0;j<*nslicey;j++)
353 snew(sg_fravg[i][j],nslicez);
354 snew(sk_fravg[i][j],nslicez);
359 find_tetra_order_grid(top,ePBC,natoms,box,x,isize[0],index[0],t,
360 &sg,&sk,*nslicex,*nslicey,nslicez,sg_grid,sk_grid);
361 for(i=0;i<*nslicex;i++)
363 for(j=0;j<*nslicey;j++)
365 for(k=0;k<nslicez;k++)
367 sk_fravg[i][j][k]+=sk_grid[i][j][k]/tblock;
368 sg_fravg[i][j][k]+=sg_grid[i][j][k]/tblock;
375 if(framenr%tblock== 0)
377 sk_4d[*nframes] = sk_fravg;
378 sg_4d[*nframes] = sg_fravg;
382 } while (read_next_x(oenv,status,&t,natoms,x,box));
389 /*Debugging for printing out the entire order parameter meshes.*/
392 fpsg = xvgropen(sgslfn,"S\\sg\\N Angle Order Parameter / Meshpoint","(nm)","S\\sg\\N",oenv);
393 fpsk = xvgropen(skslfn,"S\\sk\\N Distance Order Parameter / Meshpoint","(nm)","S\\sk\\N",oenv);
394 for(n=0;n<(*nframes);n++)
396 fprintf(fpsg,"%i\n",n);
397 fprintf(fpsk,"%i\n",n);
398 for(i=0; (i<*nslicex); i++)
400 for(j=0;j<*nslicey;j++)
402 for(k=0;k<nslicez;k++)
404 fprintf(fpsg,"%4f %4f %4f %8f\n",(i+0.5)*box[XX][XX]/(*nslicex),(j+0.5)*box[YY][YY]/(*nslicey),(k+0.5)*box[ZZ][ZZ]/nslicez,sg_4d[n][i][j][k]);
405 fprintf(fpsk,"%4f %4f %4f %8f\n",(i+0.5)*box[XX][XX]/(*nslicex),(j+0.5)*box[YY][YY]/(*nslicey),(k+0.5)*box[ZZ][ZZ]/nslicez,sk_4d[n][i][j][k]);
415 /* Find positions of interface z by scanning orderparam for each frame and for each xy-mesh cylinder along z*/
417 /*Simple trial: assume interface is in the middle of -sgang1 and sgang2*/
418 sgintf=0.5*(sgang1+sgang2);
421 /*Allocate memory for interface arrays; */
423 snew((*intfpos)[0],*nframes);
424 snew((*intfpos)[1],*nframes);
426 bins=(*nslicex)*(*nslicey);
429 snew(perm,nslicez); /*permutation array for sorting along normal coordinate*/
432 for (n=0;n<*nframes;n++)
434 snew((*intfpos)[0][n],bins);
435 snew((*intfpos)[1][n],bins);
436 for(i=0;i<*nslicex;i++)
438 for(j=0;j<*nslicey;j++)
440 rangeArray(perm,nslicez); /*reset permutation array to identity*/
441 /*Binsearch returns 2 bin-numbers where the order param is <= setpoint sgintf*/
442 ndx1=start_binsearch(sg_4d[n][i][j],perm,0,nslicez/2-1,sgintf,1);
443 ndx2=start_binsearch(sg_4d[n][i][j],perm,nslicez/2,nslicez-1,sgintf,-1);
444 /*Use linear interpolation to smooth out the interface position*/
446 /*left interface (0)*/
447 /*if((sg_4d[n][i][j][perm[ndx1+1]]-sg_4d[n][i][j][perm[ndx1]])/sg_4d[n][i][j][perm[ndx1]] > 0.01){
448 pos=( (sgintf-sg_4d[n][i][j][perm[ndx1]])*perm[ndx1+1]+(sg_4d[n][i][j][perm[ndx1+1]]-sgintf)*perm[ndx1 ])*/
449 (*intfpos)[0][n][j+*nslicey*i]=(perm[ndx1]+onehalf)*binw;
450 /*right interface (1)*/
451 /*alpha=(sgintf-sg_4d[n][i][j][perm[ndx2]])/(sg_4d[n][i][j][perm[ndx2]+1]-sg_4d[n][i][j][perm[ndx2]]);*/
452 /*(*intfpos)[1][n][j+*nslicey*i]=((1-alpha)*perm[ndx2]+alpha*(perm[ndx2]+1)+onehalf)*box[ZZ][ZZ]/nslicez;*/
453 (*intfpos)[1][n][j+*nslicey*i]=(perm[ndx2]+onehalf)*binw;
468 static void writesurftoxpms(real ***surf, int tblocks,int xbins, int ybins, real bw,char **outfiles,int maplevels )
473 real **profile1, **profile2;
474 real max1, max2, min1, min2, *xticks, *yticks;
477 FILE *xpmfile1, *xpmfile2;
479 /*Prepare xpm structures for output*/
481 /*Allocate memory to tick's and matrices*/
482 snew (xticks,xbins+1);
483 snew (yticks,ybins+1);
485 profile1=mk_matrix(xbins,ybins,FALSE);
486 profile2=mk_matrix(xbins,ybins,FALSE);
488 for (i=0;i<xbins+1;i++) xticks[i]+=bw;
489 for (j=0;j<ybins+1;j++) yticks[j]+=bw;
491 xpmfile1 = ffopen(outfiles[0],"w");
492 xpmfile2 = ffopen(outfiles[1],"w");
497 for(n=0;n<tblocks;n++)
499 sprintf(numbuf,"%5d",n);
500 /*Filling matrices for inclusion in xpm-files*/
505 profile1[i][j]=(surf[0][n][j+ybins*i]);
506 profile2[i][j]=(surf[1][n][j+ybins*i]);
507 /*Finding max and min values*/
508 if(profile1[i][j]>max1) max1=profile1[i][j];
509 if(profile1[i][j]<min1) min1=profile1[i][j];
510 if(profile2[i][j]>max2) max2=profile2[i][j];
511 if(profile2[i][j]<min2) min2=profile2[i][j];
515 write_xpm(xpmfile1,3,numbuf,"Height","x[nm]","y[nm]",xbins,ybins,xticks,yticks,profile1,min1,max1,lo,hi,&maplevels);
516 write_xpm(xpmfile2,3,numbuf,"Height","x[nm]","y[nm]",xbins,ybins,xticks,yticks,profile2,min2,max2,lo,hi,&maplevels);
531 static void writeraw(real ***surf, int tblocks,int xbins, int ybins,char **fnms){
535 raw1=ffopen(fnms[0],"w");
536 raw2=ffopen(fnms[1],"w");
537 fprintf(raw1,"#Legend\n#TBlock\n#Xbin Ybin Z t\n");
538 fprintf(raw2,"#Legend\n#TBlock\n#Xbin Ybin Z t\n");
539 for (n=0;n<tblocks;n++)
541 fprintf(raw1,"%5d\n",n);
542 fprintf(raw2,"%5d\n",n);
547 fprintf(raw1,"%i %i %8.5f\n",i,j,(surf[0][n][j+ybins*i]));
548 fprintf(raw2,"%i %i %8.5f\n",i,j,(surf[1][n][j+ybins*i]));
559 int gmx_hydorder(int argc,char *argv[])
561 static const char *desc[] = {
562 "g_hydorder computes the tetrahedrality order parameters around a ",
563 "given atom. Both angle an distance order parameters are calculated. See",
564 "P.-L. Chau and A.J. Hardwick, Mol. Phys., 93, (1998), 511-518.",
565 "for more details.[BR]"
566 "This application calculates the orderparameter in a 3d-mesh in the box, and",
567 "with 2 phases in the box gives the user the option to define a 2D interface in time",
568 "separating the faces by specifying parameters -sgang1 and -sgang2 (It is important",
569 "to select these judiciously)"};
572 static int nsttblock=1;
573 static int nlevels=100;
574 static real binwidth = 1.0; /* binwidth in mesh */
576 static real sg2 = 1; /* order parameters for bulk phases */
577 static gmx_bool bFourier = FALSE;
578 static gmx_bool bRawOut = FALSE;
579 int frames,xslices,yslices; /* Dimensions of interface arrays*/
580 real ***intfpos; /* Interface arrays (intfnr,t,xy) -potentially large */
581 static char *normal_axis[] = { NULL, "z", "x", "y", NULL };
584 { "-d", FALSE, etENUM, {normal_axis},
585 "Direction of the normal on the membrane" },
586 { "-bw", FALSE, etREAL, {&binwidth},
587 "Binwidth of box mesh" },
588 { "-sgang1",FALSE,etREAL, {&sg1},
589 "tetrahedral angle parameter in Phase 1 (bulk)" },
590 { "-sgang2",FALSE,etREAL, {&sg2},
591 "tetrahedral angle parameter in Phase 2 (bulk)" },
592 { "-tblock",FALSE,etINT,{&nsttblock},
593 "Number of frames in one time-block average"},
594 { "-nlevel", FALSE,etINT, {&nlevels},
595 "Number of Height levels in 2D - XPixMaps"}
598 t_filenm fnm[] = { /* files for g_order */
599 { efTRX, "-f", NULL, ffREAD }, /* trajectory file */
600 { efNDX, "-n", NULL, ffREAD }, /* index file */
601 { efTPX, "-s", NULL, ffREAD }, /* topology file */
602 { efXPM, "-o", "intf", ffWRMULT}, /* XPM- surface maps */
603 { efOUT,"-or","raw", ffOPTWRMULT }, /* xvgr output file */
604 { efOUT,"-Spect","intfspect",ffOPTWRMULT}, /* Fourier spectrum interfaces */
606 #define NFILE asize(fnm)
609 const char *ndxfnm,*tpsfnm,*trxfnm;
610 char **spectra,**intfn, **raw;
611 int nfspect,nfxpm, nfraw;
614 CopyRight(stderr,argv[0]);
616 parse_common_args(&argc,argv,PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
617 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
618 bFourier= opt2bSet("-Spect",NFILE,fnm);
619 bRawOut = opt2bSet("-or",NFILE,fnm);
622 gmx_fatal(FARGS,"Can not have binwidth < 0");
624 ndxfnm = ftp2fn(efNDX,NFILE,fnm);
625 tpsfnm = ftp2fn(efTPX,NFILE,fnm);
626 trxfnm = ftp2fn(efTRX,NFILE,fnm);
629 if (strcmp(normal_axis[0],"x") == 0) axis = XX;
630 else if (strcmp(normal_axis[0],"y") == 0) axis = YY;
631 else if (strcmp(normal_axis[0],"z") == 0) axis = ZZ;
632 else gmx_fatal(FARGS,"Invalid axis, use x, y or z");
636 fprintf(stderr,"Taking x axis as normal to the membrane\n");
639 fprintf(stderr,"Taking y axis as normal to the membrane\n");
642 fprintf(stderr,"Taking z axis as normal to the membrane\n");
646 /* tetraheder order parameter */
647 /* If either of the options is set we compute both */
648 nfxpm=opt2fns(&intfn,"-o",NFILE,fnm);
651 gmx_fatal(FARGS,"No or not correct number (2) of output-files: %d",nfxpm);
653 calc_tetra_order_interface(ndxfnm,tpsfnm,trxfnm,binwidth,nsttblock,&frames,&xslices,&yslices,sg1,sg2,&intfpos,oenv);
654 writesurftoxpms(intfpos,frames,xslices,yslices,binwidth,intfn,nlevels);
658 nfspect=opt2fns(&spectra,"-Spect",NFILE,fnm);
661 gmx_fatal(FARGS,"No or not correct number (2) of output-files: %d",nfspect);
663 powerspectavg(intfpos, frames,xslices,yslices,spectra);
668 nfraw=opt2fns(&raw,"-or",NFILE,fnm);
671 gmx_fatal(FARGS,"No or not correct number (2) of output-files: %d",nfraw);
673 writeraw(intfpos,frames,xslices,yslices,raw);