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 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms, box);
139 gmx_rmpbc(gpbc, natoms, box, x);
144 for (i = 0; (i < maxidx); i++)
145 { /* loop over index file */
147 for (j = 0; (j < maxidx); j++)
157 pbc_dx(&pbc, x[ix], x[jx], dx);
160 /* set_mat_entry(dmat,i,j,r2); */
162 /* determine the nearest neighbours */
165 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
166 r_nn[2][i] = r_nn[1][i]; nn[2][i] = nn[1][i];
167 r_nn[1][i] = r_nn[0][i]; nn[1][i] = nn[0][i];
168 r_nn[0][i] = r2; nn[0][i] = j;
170 else if (r2 < r_nn[1][i])
172 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
173 r_nn[2][i] = r_nn[1][i]; nn[2][i] = nn[1][i];
174 r_nn[1][i] = r2; nn[1][i] = j;
176 else if (r2 < r_nn[2][i])
178 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
179 r_nn[2][i] = r2; nn[2][i] = j;
181 else if (r2 < r_nn[3][i])
183 r_nn[3][i] = r2; nn[3][i] = j;
188 /* calculate mean distance between nearest neighbours */
190 for (j = 0; (j < 4); j++)
192 r_nn[j][i] = sqrt(r_nn[j][i]);
201 /* Chau1998a eqn 3 */
202 /* angular part tetrahedrality order parameter per atom */
203 for (j = 0; (j < 3); j++)
205 for (k = j+1; (k < 4); k++)
207 pbc_dx(&pbc, x[ix], x[index[nn[k][i]]], rk);
208 pbc_dx(&pbc, x[ix], x[index[nn[j][i]]], rj);
213 cost = iprod(urk, urj) + onethird;
221 /* normalize sgmol between 0.0 and 1.0 */
222 sgmol[i] = 3*sgmol[i]/32;
225 /* distance part tetrahedrality order parameter per atom */
226 rmean2 = 4 * 3 * rmean * rmean;
227 for (j = 0; (j < 4); j++)
229 skmol[i] += (rmean - r_nn[j][i]) * (rmean - r_nn[j][i]) / rmean2;
230 /* printf("%d %f (%f %f %f %f) \n",
231 i, skmol[i], rmean, rmean2, r_nn[j][i], (rmean - r_nn[j][i]) );
237 /* Compute sliced stuff in x y z*/
238 slindex_x = gmx_nint((1+x[i][XX]/box[XX][XX])*nslicex) % nslicex;
239 slindex_y = gmx_nint((1+x[i][YY]/box[YY][YY])*nslicey) % nslicey;
240 slindex_z = gmx_nint((1+x[i][ZZ]/box[ZZ][ZZ])*nslicez) % nslicez;
241 sggrid[slindex_x][slindex_y][slindex_z] += sgmol[i];
242 skgrid[slindex_x][slindex_y][slindex_z] += skmol[i];
243 (sl_count[slindex_x][slindex_y][slindex_z])++;
244 } /* loop over entries in index file */
249 for (i = 0; (i < nslicex); i++)
251 for (j = 0; j < nslicey; j++)
253 for (k = 0; k < nslicez; k++)
255 if (sl_count[i][j][k] > 0)
257 sggrid[i][j][k] /= sl_count[i][j][k];
258 skgrid[i][j][k] /= sl_count[i][j][k];
267 for (i = 0; (i < 4); i++)
274 /*Determines interface from tetrahedral order parameter in box with specified binwidth. */
275 /*Outputs interface positions(bins), the number of timeframes, and the number of surface-mesh points in xy*/
277 static void calc_tetra_order_interface(const char *fnNDX, const char *fnTPS, const char *fnTRX, real binw, int tblock,
278 int *nframes, int *nslicex, int *nslicey,
279 real sgang1, real sgang2, real ****intfpos,
282 FILE *fpsg = NULL, *fpsk = NULL;
283 char *sgslfn = "sg_ang_mesh"; /* Hardcoded filenames for debugging*/
284 char *skslfn = "sk_dist_mesh";
287 char title[STRLEN], subtitle[STRLEN];
293 real sg, sk, sgintf, pos;
294 atom_id **index = NULL;
295 char **grpname = NULL;
296 int i, j, k, n, *isize, ng, nslicez, framenr;
297 real ***sg_grid = NULL, ***sk_grid = NULL, ***sg_fravg = NULL, ***sk_fravg = NULL, ****sk_4d = NULL, ****sg_4d = NULL;
301 const real onehalf = 1.0/2.0;
302 /* real ***intfpos[2]; pointers to arrays of two interface positions zcoord(framenr,xbin,ybin): intfpos[interface_index][t][nslicey*x+y]
303 * i.e 1D Row-major order in (t,x,y) */
306 read_tps_conf(fnTPS, title, &top, &ePBC, &xtop, NULL, box, FALSE);
308 *nslicex = (int)(box[XX][XX]/binw + onehalf); /*Calculate slicenr from binwidth*/
309 *nslicey = (int)(box[YY][YY]/binw + onehalf);
310 nslicez = (int)(box[ZZ][ZZ]/binw + onehalf);
315 /* get index groups */
316 printf("Select the group that contains the atoms you want to use for the tetrahedrality order parameter calculation:\n");
320 get_index(&top.atoms, fnNDX, ng, isize, index, grpname);
322 /* Analyze trajectory */
323 natoms = read_first_x(oenv, &status, fnTRX, &t, &x, box);
324 if (natoms > top.atoms.nr)
326 gmx_fatal(FARGS, "Topology (%d atoms) does not match trajectory (%d atoms)",
327 top.atoms.nr, natoms);
329 check_index(NULL, ng, index[0], NULL, natoms);
332 /*Prepare structures for temporary storage of frame info*/
333 snew(sg_grid, *nslicex);
334 snew(sk_grid, *nslicex);
335 for (i = 0; i < *nslicex; i++)
337 snew(sg_grid[i], *nslicey);
338 snew(sk_grid[i], *nslicey);
339 for (j = 0; j < *nslicey; j++)
341 snew(sg_grid[i][j], nslicez);
342 snew(sk_grid[i][j], nslicez);
351 /* Loop over frames*/
354 /*Initialize box meshes (temporary storage for each tblock frame -reinitialise every tblock steps */
355 if (framenr%tblock == 0)
357 srenew(sk_4d, *nframes+1);
358 srenew(sg_4d, *nframes+1);
359 snew(sg_fravg, *nslicex);
360 snew(sk_fravg, *nslicex);
361 for (i = 0; i < *nslicex; i++)
363 snew(sg_fravg[i], *nslicey);
364 snew(sk_fravg[i], *nslicey);
365 for (j = 0; j < *nslicey; j++)
367 snew(sg_fravg[i][j], nslicez);
368 snew(sk_fravg[i][j], nslicez);
373 find_tetra_order_grid(top, ePBC, natoms, box, x, isize[0], index[0], t,
374 &sg, &sk, *nslicex, *nslicey, nslicez, sg_grid, sk_grid);
375 for (i = 0; i < *nslicex; i++)
377 for (j = 0; j < *nslicey; j++)
379 for (k = 0; k < nslicez; k++)
381 sk_fravg[i][j][k] += sk_grid[i][j][k]/tblock;
382 sg_fravg[i][j][k] += sg_grid[i][j][k]/tblock;
389 if (framenr%tblock == 0)
391 sk_4d[*nframes] = sk_fravg;
392 sg_4d[*nframes] = sg_fravg;
397 while (read_next_x(oenv, status, &t, natoms, x, box));
404 /*Debugging for printing out the entire order parameter meshes.*/
407 fpsg = xvgropen(sgslfn, "S\\sg\\N Angle Order Parameter / Meshpoint", "(nm)", "S\\sg\\N", oenv);
408 fpsk = xvgropen(skslfn, "S\\sk\\N Distance Order Parameter / Meshpoint", "(nm)", "S\\sk\\N", oenv);
409 for (n = 0; n < (*nframes); n++)
411 fprintf(fpsg, "%i\n", n);
412 fprintf(fpsk, "%i\n", n);
413 for (i = 0; (i < *nslicex); i++)
415 for (j = 0; j < *nslicey; j++)
417 for (k = 0; k < nslicez; k++)
419 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]);
420 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]);
430 /* Find positions of interface z by scanning orderparam for each frame and for each xy-mesh cylinder along z*/
432 /*Simple trial: assume interface is in the middle of -sgang1 and sgang2*/
433 sgintf = 0.5*(sgang1+sgang2);
436 /*Allocate memory for interface arrays; */
438 snew((*intfpos)[0], *nframes);
439 snew((*intfpos)[1], *nframes);
441 bins = (*nslicex)*(*nslicey);
444 snew(perm, nslicez); /*permutation array for sorting along normal coordinate*/
447 for (n = 0; n < *nframes; n++)
449 snew((*intfpos)[0][n], bins);
450 snew((*intfpos)[1][n], bins);
451 for (i = 0; i < *nslicex; i++)
453 for (j = 0; j < *nslicey; j++)
455 rangeArray(perm, nslicez); /*reset permutation array to identity*/
456 /*Binsearch returns 2 bin-numbers where the order param is <= setpoint sgintf*/
457 ndx1 = start_binsearch(sg_4d[n][i][j], perm, 0, nslicez/2-1, sgintf, 1);
458 ndx2 = start_binsearch(sg_4d[n][i][j], perm, nslicez/2, nslicez-1, sgintf, -1);
459 /*Use linear interpolation to smooth out the interface position*/
461 /*left interface (0)*/
462 /*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){
463 pos=( (sgintf-sg_4d[n][i][j][perm[ndx1]])*perm[ndx1+1]+(sg_4d[n][i][j][perm[ndx1+1]]-sgintf)*perm[ndx1 ])*/
464 (*intfpos)[0][n][j+*nslicey*i] = (perm[ndx1]+onehalf)*binw;
465 /*right interface (1)*/
466 /*alpha=(sgintf-sg_4d[n][i][j][perm[ndx2]])/(sg_4d[n][i][j][perm[ndx2]+1]-sg_4d[n][i][j][perm[ndx2]]);*/
467 /*(*intfpos)[1][n][j+*nslicey*i]=((1-alpha)*perm[ndx2]+alpha*(perm[ndx2]+1)+onehalf)*box[ZZ][ZZ]/nslicez;*/
468 (*intfpos)[1][n][j+*nslicey*i] = (perm[ndx2]+onehalf)*binw;
483 static void writesurftoxpms(real ***surf, int tblocks, int xbins, int ybins, real bw, char **outfiles, int maplevels )
488 real **profile1, **profile2;
489 real max1, max2, min1, min2, *xticks, *yticks;
490 t_rgb lo = {1, 1, 1};
491 t_rgb hi = {0, 0, 0};
492 FILE *xpmfile1, *xpmfile2;
494 /*Prepare xpm structures for output*/
496 /*Allocate memory to tick's and matrices*/
497 snew (xticks, xbins+1);
498 snew (yticks, ybins+1);
500 profile1 = mk_matrix(xbins, ybins, FALSE);
501 profile2 = mk_matrix(xbins, ybins, FALSE);
503 for (i = 0; i < xbins+1; i++)
507 for (j = 0; j < ybins+1; j++)
512 xpmfile1 = ffopen(outfiles[0], "w");
513 xpmfile2 = ffopen(outfiles[1], "w");
516 min1 = min2 = 1000.00;
518 for (n = 0; n < tblocks; n++)
520 sprintf(numbuf, "%5d", n);
521 /*Filling matrices for inclusion in xpm-files*/
522 for (i = 0; i < xbins; i++)
524 for (j = 0; j < ybins; j++)
526 profile1[i][j] = (surf[0][n][j+ybins*i]);
527 profile2[i][j] = (surf[1][n][j+ybins*i]);
528 /*Finding max and min values*/
529 if (profile1[i][j] > max1)
531 max1 = profile1[i][j];
533 if (profile1[i][j] < min1)
535 min1 = profile1[i][j];
537 if (profile2[i][j] > max2)
539 max2 = profile2[i][j];
541 if (profile2[i][j] < min2)
543 min2 = profile2[i][j];
548 write_xpm(xpmfile1, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile1, min1, max1, lo, hi, &maplevels);
549 write_xpm(xpmfile2, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile2, min2, max2, lo, hi, &maplevels);
564 static void writeraw(real ***surf, int tblocks, int xbins, int ybins, char **fnms)
569 raw1 = ffopen(fnms[0], "w");
570 raw2 = ffopen(fnms[1], "w");
571 fprintf(raw1, "#Legend\n#TBlock\n#Xbin Ybin Z t\n");
572 fprintf(raw2, "#Legend\n#TBlock\n#Xbin Ybin Z t\n");
573 for (n = 0; n < tblocks; n++)
575 fprintf(raw1, "%5d\n", n);
576 fprintf(raw2, "%5d\n", n);
577 for (i = 0; i < xbins; i++)
579 for (j = 0; j < ybins; j++)
581 fprintf(raw1, "%i %i %8.5f\n", i, j, (surf[0][n][j+ybins*i]));
582 fprintf(raw2, "%i %i %8.5f\n", i, j, (surf[1][n][j+ybins*i]));
593 int gmx_hydorder(int argc, char *argv[])
595 static const char *desc[] = {
596 "g_hydorder computes the tetrahedrality order parameters around a ",
597 "given atom. Both angle an distance order parameters are calculated. See",
598 "P.-L. Chau and A.J. Hardwick, Mol. Phys., 93, (1998), 511-518.",
599 "for more details.[BR]"
600 "This application calculates the orderparameter in a 3d-mesh in the box, and",
601 "with 2 phases in the box gives the user the option to define a 2D interface in time",
602 "separating the faces by specifying parameters -sgang1 and -sgang2 (It is important",
603 "to select these judiciously)"
607 static int nsttblock = 1;
608 static int nlevels = 100;
609 static real binwidth = 1.0; /* binwidth in mesh */
611 static real sg2 = 1; /* order parameters for bulk phases */
612 static gmx_bool bFourier = FALSE;
613 static gmx_bool bRawOut = FALSE;
614 int frames, xslices, yslices; /* Dimensions of interface arrays*/
615 real ***intfpos; /* Interface arrays (intfnr,t,xy) -potentially large */
616 static char *normal_axis[] = { NULL, "z", "x", "y", NULL };
619 { "-d", FALSE, etENUM, {normal_axis},
620 "Direction of the normal on the membrane" },
621 { "-bw", FALSE, etREAL, {&binwidth},
622 "Binwidth of box mesh" },
623 { "-sgang1", FALSE, etREAL, {&sg1},
624 "tetrahedral angle parameter in Phase 1 (bulk)" },
625 { "-sgang2", FALSE, etREAL, {&sg2},
626 "tetrahedral angle parameter in Phase 2 (bulk)" },
627 { "-tblock", FALSE, etINT, {&nsttblock},
628 "Number of frames in one time-block average"},
629 { "-nlevel", FALSE, etINT, {&nlevels},
630 "Number of Height levels in 2D - XPixMaps"}
633 t_filenm fnm[] = { /* files for g_order */
634 { efTRX, "-f", NULL, ffREAD }, /* trajectory file */
635 { efNDX, "-n", NULL, ffREAD }, /* index file */
636 { efTPX, "-s", NULL, ffREAD }, /* topology file */
637 { efXPM, "-o", "intf", ffWRMULT}, /* XPM- surface maps */
638 { efOUT, "-or", "raw", ffOPTWRMULT }, /* xvgr output file */
639 { efOUT, "-Spect", "intfspect", ffOPTWRMULT}, /* Fourier spectrum interfaces */
641 #define NFILE asize(fnm)
644 const char *ndxfnm, *tpsfnm, *trxfnm;
645 char **spectra, **intfn, **raw;
646 int nfspect, nfxpm, nfraw;
649 CopyRight(stderr, argv[0]);
651 parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
652 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
653 bFourier = opt2bSet("-Spect", NFILE, fnm);
654 bRawOut = opt2bSet("-or", NFILE, fnm);
658 gmx_fatal(FARGS, "Can not have binwidth < 0");
661 ndxfnm = ftp2fn(efNDX, NFILE, fnm);
662 tpsfnm = ftp2fn(efTPX, NFILE, fnm);
663 trxfnm = ftp2fn(efTRX, NFILE, fnm);
666 if (strcmp(normal_axis[0], "x") == 0)
670 else if (strcmp(normal_axis[0], "y") == 0)
674 else if (strcmp(normal_axis[0], "z") == 0)
680 gmx_fatal(FARGS, "Invalid axis, use x, y or z");
686 fprintf(stderr, "Taking x axis as normal to the membrane\n");
689 fprintf(stderr, "Taking y axis as normal to the membrane\n");
692 fprintf(stderr, "Taking z axis as normal to the membrane\n");
696 /* tetraheder order parameter */
697 /* If either of the options is set we compute both */
698 nfxpm = opt2fns(&intfn, "-o", NFILE, fnm);
701 gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfxpm);
703 calc_tetra_order_interface(ndxfnm, tpsfnm, trxfnm, binwidth, nsttblock, &frames, &xslices, &yslices, sg1, sg2, &intfpos, oenv);
704 writesurftoxpms(intfpos, frames, xslices, yslices, binwidth, intfn, nlevels);
708 nfspect = opt2fns(&spectra, "-Spect", NFILE, fnm);
711 gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfspect);
713 powerspectavg(intfpos, frames, xslices, yslices, spectra);
718 nfraw = opt2fns(&raw, "-or", NFILE, fnm);
721 gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfraw);
723 writeraw(intfpos, frames, xslices, yslices, raw);