2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2001
5 * BIOSON Research Institute, Dept. of Biophysical Chemistry
6 * University of Groningen, The Netherlands
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
62 #include "binsearch.h"
63 #include "powerspect.h"
66 /* Print name of first atom in all groups in index file */
67 static void print_types(atom_id index[], atom_id a[], int ngrps,
68 char *groups[], t_topology *top)
72 fprintf(stderr, "Using following groups: \n");
73 for (i = 0; i < ngrps; i++)
75 fprintf(stderr, "Groupname: %s First atomname: %s First atomnr %u\n",
76 groups[i], *(top->atoms.atomname[a[index[i]]]), a[index[i]]);
78 fprintf(stderr, "\n");
81 static void check_length(real length, int a, int b)
85 fprintf(stderr, "WARNING: distance between atoms %d and "
86 "%d > 0.3 nm (%f). Index file might be corrupt.\n",
91 static void find_tetra_order_grid(t_topology top, int ePBC,
92 int natoms, matrix box,
93 rvec x[], int maxidx, atom_id index[],
94 real time, real *sgmean, real *skmean,
95 int nslicex, int nslicey, int nslicez,
96 real ***sggrid, real ***skgrid)
98 int ix, jx, i, j, k, l, n, *nn[4];
99 rvec dx, rj, rk, urk, urj;
100 real cost, cost2, *sgmol, *skmol, rmean, rmean2, r2, box2, *r_nn[4];
102 int slindex_x, slindex_y, slindex_z;
104 real onethird = 1.0/3.0;
107 /* dmat = init_mat(maxidx, FALSE); */
109 box2 = box[XX][XX] * box[XX][XX];
111 /* Initialize expanded sl_count array */
112 snew(sl_count, nslicex);
113 for (i = 0; i < nslicex; i++)
115 snew(sl_count[i], nslicey);
116 for (j = 0; j < nslicey; j++)
118 snew(sl_count[i][j], nslicez);
123 for (i = 0; (i < 4); i++)
125 snew(r_nn[i], natoms);
128 for (j = 0; (j < natoms); j++)
137 /* Must init pbc every step because of pressure coupling */
138 set_pbc(&pbc,ePBC,box);
139 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms, box);
140 gmx_rmpbc(gpbc, natoms, box, x);
145 for (i = 0; (i < maxidx); i++)
146 { /* loop over index file */
148 for (j = 0; (j < maxidx); j++)
158 pbc_dx(&pbc, x[ix], x[jx], dx);
161 /* set_mat_entry(dmat,i,j,r2); */
163 /* determine the nearest neighbours */
166 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
167 r_nn[2][i] = r_nn[1][i]; nn[2][i] = nn[1][i];
168 r_nn[1][i] = r_nn[0][i]; nn[1][i] = nn[0][i];
169 r_nn[0][i] = r2; nn[0][i] = j;
171 else if (r2 < r_nn[1][i])
173 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
174 r_nn[2][i] = r_nn[1][i]; nn[2][i] = nn[1][i];
175 r_nn[1][i] = r2; nn[1][i] = j;
177 else if (r2 < r_nn[2][i])
179 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
180 r_nn[2][i] = r2; nn[2][i] = j;
182 else if (r2 < r_nn[3][i])
184 r_nn[3][i] = r2; nn[3][i] = j;
189 /* calculate mean distance between nearest neighbours */
191 for (j = 0; (j < 4); j++)
193 r_nn[j][i] = sqrt(r_nn[j][i]);
202 /* Chau1998a eqn 3 */
203 /* angular part tetrahedrality order parameter per atom */
204 for (j = 0; (j < 3); j++)
206 for (k = j+1; (k < 4); k++)
208 pbc_dx(&pbc, x[ix], x[index[nn[k][i]]], rk);
209 pbc_dx(&pbc, x[ix], x[index[nn[j][i]]], rj);
214 cost = iprod(urk, urj) + onethird;
222 /* normalize sgmol between 0.0 and 1.0 */
223 sgmol[i] = 3*sgmol[i]/32;
226 /* distance part tetrahedrality order parameter per atom */
227 rmean2 = 4 * 3 * rmean * rmean;
228 for (j = 0; (j < 4); j++)
230 skmol[i] += (rmean - r_nn[j][i]) * (rmean - r_nn[j][i]) / rmean2;
231 /* printf("%d %f (%f %f %f %f) \n",
232 i, skmol[i], rmean, rmean2, r_nn[j][i], (rmean - r_nn[j][i]) );
238 /* Compute sliced stuff in x y z*/
239 slindex_x = gmx_nint((1+x[i][XX]/box[XX][XX])*nslicex) % nslicex;
240 slindex_y = gmx_nint((1+x[i][YY]/box[YY][YY])*nslicey) % nslicey;
241 slindex_z = gmx_nint((1+x[i][ZZ]/box[ZZ][ZZ])*nslicez) % nslicez;
242 sggrid[slindex_x][slindex_y][slindex_z] += sgmol[i];
243 skgrid[slindex_x][slindex_y][slindex_z] += skmol[i];
244 (sl_count[slindex_x][slindex_y][slindex_z])++;
245 } /* loop over entries in index file */
250 for (i = 0; (i < nslicex); i++)
252 for (j = 0; j < nslicey; j++)
254 for (k = 0; k < nslicez; k++)
256 if (sl_count[i][j][k] > 0)
258 sggrid[i][j][k] /= sl_count[i][j][k];
259 skgrid[i][j][k] /= sl_count[i][j][k];
268 for (i = 0; (i < 4); i++)
275 /*Determines interface from tetrahedral order parameter in box with specified binwidth. */
276 /*Outputs interface positions(bins), the number of timeframes, and the number of surface-mesh points in xy*/
278 static void calc_tetra_order_interface(const char *fnNDX, const char *fnTPS, const char *fnTRX, real binw, int tblock,
279 int *nframes, int *nslicex, int *nslicey,
280 real sgang1, real sgang2, real ****intfpos,
283 FILE *fpsg = NULL, *fpsk = NULL;
284 char *sgslfn = "sg_ang_mesh"; /* Hardcoded filenames for debugging*/
285 char *skslfn = "sk_dist_mesh";
288 char title[STRLEN], subtitle[STRLEN];
294 real sg, sk, sgintf, pos;
295 atom_id **index = NULL;
296 char **grpname = NULL;
297 int i, j, k, n, *isize, ng, nslicez, framenr;
298 real ***sg_grid = NULL, ***sk_grid = NULL, ***sg_fravg = NULL, ***sk_fravg = NULL, ****sk_4d = NULL, ****sg_4d = NULL;
302 const real onehalf = 1.0/2.0;
303 /* real ***intfpos[2]; pointers to arrays of two interface positions zcoord(framenr,xbin,ybin): intfpos[interface_index][t][nslicey*x+y]
304 * i.e 1D Row-major order in (t,x,y) */
307 read_tps_conf(fnTPS, title, &top, &ePBC, &xtop, NULL, box, FALSE);
309 *nslicex = (int)(box[XX][XX]/binw + onehalf); /*Calculate slicenr from binwidth*/
310 *nslicey = (int)(box[YY][YY]/binw + onehalf);
311 nslicez = (int)(box[ZZ][ZZ]/binw + onehalf);
316 /* get index groups */
317 printf("Select the group that contains the atoms you want to use for the tetrahedrality order parameter calculation:\n");
321 get_index(&top.atoms, fnNDX, ng, isize, index, grpname);
323 /* Analyze trajectory */
324 natoms = read_first_x(oenv, &status, fnTRX, &t, &x, box);
325 if (natoms > top.atoms.nr)
327 gmx_fatal(FARGS, "Topology (%d atoms) does not match trajectory (%d atoms)",
328 top.atoms.nr, natoms);
330 check_index(NULL, ng, index[0], NULL, natoms);
333 /*Prepare structures for temporary storage of frame info*/
334 snew(sg_grid, *nslicex);
335 snew(sk_grid, *nslicex);
336 for (i = 0; i < *nslicex; i++)
338 snew(sg_grid[i], *nslicey);
339 snew(sk_grid[i], *nslicey);
340 for (j = 0; j < *nslicey; j++)
342 snew(sg_grid[i][j], nslicez);
343 snew(sk_grid[i][j], nslicez);
352 /* Loop over frames*/
355 /*Initialize box meshes (temporary storage for each tblock frame -reinitialise every tblock steps */
356 if (framenr%tblock == 0)
358 srenew(sk_4d, *nframes+1);
359 srenew(sg_4d, *nframes+1);
360 snew(sg_fravg, *nslicex);
361 snew(sk_fravg, *nslicex);
362 for (i = 0; i < *nslicex; i++)
364 snew(sg_fravg[i], *nslicey);
365 snew(sk_fravg[i], *nslicey);
366 for (j = 0; j < *nslicey; j++)
368 snew(sg_fravg[i][j], nslicez);
369 snew(sk_fravg[i][j], nslicez);
374 find_tetra_order_grid(top, ePBC, natoms, box, x, isize[0], index[0], t,
375 &sg, &sk, *nslicex, *nslicey, nslicez, sg_grid, sk_grid);
376 for (i = 0; i < *nslicex; i++)
378 for (j = 0; j < *nslicey; j++)
380 for (k = 0; k < nslicez; k++)
382 sk_fravg[i][j][k] += sk_grid[i][j][k]/tblock;
383 sg_fravg[i][j][k] += sg_grid[i][j][k]/tblock;
390 if (framenr%tblock == 0)
392 sk_4d[*nframes] = sk_fravg;
393 sg_4d[*nframes] = sg_fravg;
398 while (read_next_x(oenv, status, &t, natoms, x, box));
405 /*Debugging for printing out the entire order parameter meshes.*/
408 fpsg = xvgropen(sgslfn, "S\\sg\\N Angle Order Parameter / Meshpoint", "(nm)", "S\\sg\\N", oenv);
409 fpsk = xvgropen(skslfn, "S\\sk\\N Distance Order Parameter / Meshpoint", "(nm)", "S\\sk\\N", oenv);
410 for (n = 0; n < (*nframes); n++)
412 fprintf(fpsg, "%i\n", n);
413 fprintf(fpsk, "%i\n", n);
414 for (i = 0; (i < *nslicex); i++)
416 for (j = 0; j < *nslicey; j++)
418 for (k = 0; k < nslicez; k++)
420 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]);
421 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]);
431 /* Find positions of interface z by scanning orderparam for each frame and for each xy-mesh cylinder along z*/
433 /*Simple trial: assume interface is in the middle of -sgang1 and sgang2*/
434 sgintf = 0.5*(sgang1+sgang2);
437 /*Allocate memory for interface arrays; */
439 snew((*intfpos)[0], *nframes);
440 snew((*intfpos)[1], *nframes);
442 bins = (*nslicex)*(*nslicey);
445 snew(perm, nslicez); /*permutation array for sorting along normal coordinate*/
448 for (n = 0; n < *nframes; n++)
450 snew((*intfpos)[0][n], bins);
451 snew((*intfpos)[1][n], bins);
452 for (i = 0; i < *nslicex; i++)
454 for (j = 0; j < *nslicey; j++)
456 rangeArray(perm, nslicez); /*reset permutation array to identity*/
457 /*Binsearch returns 2 bin-numbers where the order param is <= setpoint sgintf*/
458 ndx1 = start_binsearch(sg_4d[n][i][j], perm, 0, nslicez/2-1, sgintf, 1);
459 ndx2 = start_binsearch(sg_4d[n][i][j], perm, nslicez/2, nslicez-1, sgintf, -1);
460 /*Use linear interpolation to smooth out the interface position*/
462 /*left interface (0)*/
463 /*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){
464 pos=( (sgintf-sg_4d[n][i][j][perm[ndx1]])*perm[ndx1+1]+(sg_4d[n][i][j][perm[ndx1+1]]-sgintf)*perm[ndx1 ])*/
465 (*intfpos)[0][n][j+*nslicey*i] = (perm[ndx1]+onehalf)*binw;
466 /*right interface (1)*/
467 /*alpha=(sgintf-sg_4d[n][i][j][perm[ndx2]])/(sg_4d[n][i][j][perm[ndx2]+1]-sg_4d[n][i][j][perm[ndx2]]);*/
468 /*(*intfpos)[1][n][j+*nslicey*i]=((1-alpha)*perm[ndx2]+alpha*(perm[ndx2]+1)+onehalf)*box[ZZ][ZZ]/nslicez;*/
469 (*intfpos)[1][n][j+*nslicey*i] = (perm[ndx2]+onehalf)*binw;
484 static void writesurftoxpms(real ***surf, int tblocks, int xbins, int ybins, real bw, char **outfiles, int maplevels )
489 real **profile1, **profile2;
490 real max1, max2, min1, min2, *xticks, *yticks;
491 t_rgb lo = {1, 1, 1};
492 t_rgb hi = {0, 0, 0};
493 FILE *xpmfile1, *xpmfile2;
495 /*Prepare xpm structures for output*/
497 /*Allocate memory to tick's and matrices*/
498 snew (xticks, xbins+1);
499 snew (yticks, ybins+1);
501 profile1 = mk_matrix(xbins, ybins, FALSE);
502 profile2 = mk_matrix(xbins, ybins, FALSE);
504 for (i = 0; i < xbins+1; i++)
508 for (j = 0; j < ybins+1; j++)
513 xpmfile1 = ffopen(outfiles[0], "w");
514 xpmfile2 = ffopen(outfiles[1], "w");
517 min1 = min2 = 1000.00;
519 for (n = 0; n < tblocks; n++)
521 sprintf(numbuf, "%5d", n);
522 /*Filling matrices for inclusion in xpm-files*/
523 for (i = 0; i < xbins; i++)
525 for (j = 0; j < ybins; j++)
527 profile1[i][j] = (surf[0][n][j+ybins*i]);
528 profile2[i][j] = (surf[1][n][j+ybins*i]);
529 /*Finding max and min values*/
530 if (profile1[i][j] > max1)
532 max1 = profile1[i][j];
534 if (profile1[i][j] < min1)
536 min1 = profile1[i][j];
538 if (profile2[i][j] > max2)
540 max2 = profile2[i][j];
542 if (profile2[i][j] < min2)
544 min2 = profile2[i][j];
549 write_xpm(xpmfile1, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile1, min1, max1, lo, hi, &maplevels);
550 write_xpm(xpmfile2, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile2, min2, max2, lo, hi, &maplevels);
565 static void writeraw(real ***surf, int tblocks, int xbins, int ybins, char **fnms)
570 raw1 = ffopen(fnms[0], "w");
571 raw2 = ffopen(fnms[1], "w");
572 fprintf(raw1, "#Legend\n#TBlock\n#Xbin Ybin Z t\n");
573 fprintf(raw2, "#Legend\n#TBlock\n#Xbin Ybin Z t\n");
574 for (n = 0; n < tblocks; n++)
576 fprintf(raw1, "%5d\n", n);
577 fprintf(raw2, "%5d\n", n);
578 for (i = 0; i < xbins; i++)
580 for (j = 0; j < ybins; j++)
582 fprintf(raw1, "%i %i %8.5f\n", i, j, (surf[0][n][j+ybins*i]));
583 fprintf(raw2, "%i %i %8.5f\n", i, j, (surf[1][n][j+ybins*i]));
594 int gmx_hydorder(int argc, char *argv[])
596 static const char *desc[] = {
597 "g_hydorder computes the tetrahedrality order parameters around a ",
598 "given atom. Both angle an distance order parameters are calculated. See",
599 "P.-L. Chau and A.J. Hardwick, Mol. Phys., 93, (1998), 511-518.",
600 "for more details.[BR]"
601 "This application calculates the orderparameter in a 3d-mesh in the box, and",
602 "with 2 phases in the box gives the user the option to define a 2D interface in time",
603 "separating the faces by specifying parameters -sgang1 and -sgang2 (It is important",
604 "to select these judiciously)"
608 static int nsttblock = 1;
609 static int nlevels = 100;
610 static real binwidth = 1.0; /* binwidth in mesh */
612 static real sg2 = 1; /* order parameters for bulk phases */
613 static gmx_bool bFourier = FALSE;
614 static gmx_bool bRawOut = FALSE;
615 int frames, xslices, yslices; /* Dimensions of interface arrays*/
616 real ***intfpos; /* Interface arrays (intfnr,t,xy) -potentially large */
617 static char *normal_axis[] = { NULL, "z", "x", "y", NULL };
620 { "-d", FALSE, etENUM, {normal_axis},
621 "Direction of the normal on the membrane" },
622 { "-bw", FALSE, etREAL, {&binwidth},
623 "Binwidth of box mesh" },
624 { "-sgang1", FALSE, etREAL, {&sg1},
625 "tetrahedral angle parameter in Phase 1 (bulk)" },
626 { "-sgang2", FALSE, etREAL, {&sg2},
627 "tetrahedral angle parameter in Phase 2 (bulk)" },
628 { "-tblock", FALSE, etINT, {&nsttblock},
629 "Number of frames in one time-block average"},
630 { "-nlevel", FALSE, etINT, {&nlevels},
631 "Number of Height levels in 2D - XPixMaps"}
634 t_filenm fnm[] = { /* files for g_order */
635 { efTRX, "-f", NULL, ffREAD }, /* trajectory file */
636 { efNDX, "-n", NULL, ffREAD }, /* index file */
637 { efTPX, "-s", NULL, ffREAD }, /* topology file */
638 { efXPM, "-o", "intf", ffWRMULT}, /* XPM- surface maps */
639 { efOUT, "-or", "raw", ffOPTWRMULT }, /* xvgr output file */
640 { efOUT, "-Spect", "intfspect", ffOPTWRMULT}, /* Fourier spectrum interfaces */
642 #define NFILE asize(fnm)
645 const char *ndxfnm, *tpsfnm, *trxfnm;
646 char **spectra, **intfn, **raw;
647 int nfspect, nfxpm, nfraw;
650 CopyRight(stderr, argv[0]);
652 parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
653 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
654 bFourier = opt2bSet("-Spect", NFILE, fnm);
655 bRawOut = opt2bSet("-or", NFILE, fnm);
659 gmx_fatal(FARGS, "Can not have binwidth < 0");
662 ndxfnm = ftp2fn(efNDX, NFILE, fnm);
663 tpsfnm = ftp2fn(efTPX, NFILE, fnm);
664 trxfnm = ftp2fn(efTRX, NFILE, fnm);
667 if (strcmp(normal_axis[0], "x") == 0)
671 else if (strcmp(normal_axis[0], "y") == 0)
675 else if (strcmp(normal_axis[0], "z") == 0)
681 gmx_fatal(FARGS, "Invalid axis, use x, y or z");
687 fprintf(stderr, "Taking x axis as normal to the membrane\n");
690 fprintf(stderr, "Taking y axis as normal to the membrane\n");
693 fprintf(stderr, "Taking z axis as normal to the membrane\n");
697 /* tetraheder order parameter */
698 /* If either of the options is set we compute both */
699 nfxpm = opt2fns(&intfn, "-o", NFILE, fnm);
702 gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfxpm);
704 calc_tetra_order_interface(ndxfnm, tpsfnm, trxfnm, binwidth, nsttblock, &frames, &xslices, &yslices, sg1, sg2, &intfpos, oenv);
705 writesurftoxpms(intfpos, frames, xslices, yslices, binwidth, intfn, nlevels);
709 nfspect = opt2fns(&spectra, "-Spect", NFILE, fnm);
712 gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfspect);
714 powerspectavg(intfpos, frames, xslices, yslices, spectra);
719 nfraw = opt2fns(&raw, "-or", NFILE, fnm);
722 gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfraw);
724 writeraw(intfpos, frames, xslices, yslices, raw);