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
57 #include "dens_filter.h"
58 #include "binsearch.h"
59 #include "powerspect.h"
63 #define FLOOR(x) ((int) floor(x))
65 #define FLOOR(x) ((int) floorf(x))
69 methSEL, methBISECT, methFUNCFIT, methNR
72 static void center_coords(t_atoms *atoms, matrix box, rvec x0[], int axis)
76 rvec com, shift, box_center;
80 for (i = 0; (i < atoms->nr); i++)
82 mm = atoms->atom[i].m;
84 for (m = 0; (m < DIM); m++)
86 com[m] += mm*x0[i][m];
89 for (m = 0; (m < DIM); m++)
93 calc_box_center(ecenterDEF, box, box_center);
94 rvec_sub(box_center, com, shift);
95 shift[axis] -= box_center[axis];
97 for (i = 0; (i < atoms->nr); i++)
99 rvec_dec(x0[i], shift);
104 static void density_in_time (const char *fn, atom_id **index, int gnx[], real bw, real bwz, int nsttblock, real *****Densdevel, int *xslices, int *yslices, int *zslices, int *tblock, t_topology *top, int ePBC, int axis, gmx_bool bCenter, gmx_bool bps1d, const output_env_t oenv)
108 * *****Densdevel pointer to array of density values in slices and frame-blocks Densdevel[*nsttblock][*xslices][*yslices][*zslices]
110 * nsttblock - nr of frames in each time-block
111 * bw widths of normal slices
113 * axis - axis direction (normal to slices)
114 * nndx - number ot atoms in **index
115 * grpn - group number in index
118 gmx_rmpbc_t gpbc = NULL;
119 matrix box; /* Box - 3x3 -each step*/
120 rvec *x0; /* List of Coord without PBC*/
121 int natoms, i, j, k, n, /* loop indices, checks etc*/
122 ax1 = 0, ax2 = 0, /* tangent directions */
123 framenr = 0, /* frame number in trajectory*/
124 slicex, slicey, slicez; /*slice # of x y z position */
125 real ***Densslice = NULL; /* Density-slice in one frame*/
126 real dscale; /*physical scaling factor*/
127 real t, x, y, z; /* time and coordinates*/
130 *tblock = 0; /* blocknr in block average - initialise to 0*/
131 /* Axis: X=0, Y=1,Z=2 */
135 ax1 = YY; ax2 = ZZ; /*Surface: YZ*/
138 ax1 = ZZ; ax2 = XX; /* Surface : XZ*/
141 ax1 = XX; ax2 = YY; /* Surface XY*/
144 gmx_fatal(FARGS, "Invalid axes. Terminating\n");
147 if ( (natoms = read_first_x(oenv, &status, fn, &t, &x0, box)) == 0)
149 gmx_fatal(FARGS, "Could not read coordinates from file"); /* Open trajectory for read*/
153 *zslices = 1+FLOOR(box[axis][axis]/bwz);
154 *yslices = 1+FLOOR(box[ax2][ax2]/bw);
155 *xslices = 1+FLOOR(box[ax1][ax1]/bw);
158 if (*xslices < *yslices)
168 "\nDividing the box in %5d x %5d x %5d slices with binw %f along axis %d\n", *xslices, *yslices, *zslices, bw, axis );
171 /****Start trajectory processing***/
173 /*Initialize Densdevel and PBC-remove*/
174 gpbc = gmx_rmpbc_init(&top->idef, ePBC, top->atoms.nr, box);
180 bbww[XX] = box[ax1][ax1]/ *xslices;
181 bbww[YY] = box[ax2][ax2]/ *yslices;
182 bbww[ZZ] = box[axis][axis]/ *zslices;
183 gmx_rmpbc(gpbc, top->atoms.nr, box, x0);
184 /*Reset Densslice every nsttblock steps*/
185 if (framenr % nsttblock == 0)
187 snew(Densslice, *xslices);
188 for (i = 0; i < *xslices; i++)
190 snew(Densslice[i], *yslices);
191 for (j = 0; j < *yslices; j++)
193 snew(Densslice[i][j], *zslices);
197 /*Allocate Memory to extra frame in Densdevel - rather stupid approach: *A single frame each time, although only every nsttblock steps.*/
198 srenew(*Densdevel, *tblock+1);
203 dscale = (*xslices)*(*yslices)*(*zslices)*AMU/ (box[ax1][ax1]*box[ax2][ax2]*box[axis][axis]*nsttblock*(NANO*NANO*NANO));
207 center_coords(&top->atoms, box, x0, axis);
211 for (j = 0; j < gnx[0]; j++)
212 { /*Loop over all atoms in selected index*/
213 x = x0[index[0][j]][ax1];
214 y = x0[index[0][j]][ax2];
215 z = x0[index[0][j]][axis];
220 while (x > box[ax1][ax1])
229 while (y > box[ax2][ax2])
236 z += box[axis][axis];
238 while (z > box[axis][axis])
240 z -= box[axis][axis];
243 slicex = ((int) (x/bbww[XX])) % *xslices;
244 slicey = ((int) (y/bbww[YY])) % *yslices;
245 slicez = ((int) (z/bbww[ZZ])) % *zslices;
246 Densslice[slicex][slicey][slicez] += (top->atoms.atom[index[0][j]].m*dscale);
253 if (framenr % nsttblock == 0)
255 /*Implicit incrementation of Densdevel via renewal of Densslice*/
256 /*only every nsttblock steps*/
257 (*Densdevel)[*tblock] = Densslice;
262 while (read_next_x(oenv, status, &t, natoms, x0, box));
265 /*Free memory we no longer need and exit.*/
266 gmx_rmpbc_done(gpbc);
272 fp = fopen("koko.xvg", "w");
273 for (j = 0; (j < *zslices); j++)
275 fprintf(fp, "%5d", j);
276 for (i = 0; (i < *tblock); i++)
278 fprintf(fp, " %10g", (*Densdevel)[i][9][1][j]);
287 static void outputfield(const char *fldfn, real ****Densmap,
288 int xslices, int yslices, int zslices, int tdim)
290 /*Debug-filename and filehandle*/
301 fldH = ffopen(fldfn, "w");
302 fwrite(dim, sizeof(int), 4, fldH);
303 for (n = 0; n < tdim; n++)
305 for (i = 0; i < xslices; i++)
307 for (j = 0; j < yslices; j++)
309 for (k = 0; k < zslices; k++)
311 fwrite(&(Densmap[n][i][j][k]), sizeof(real), 1, fldH);
312 totdens += (Densmap[n][i][j][k]);
317 totdens /= (xslices*yslices*zslices*tdim);
318 fprintf(stderr, "Total density [kg/m^3] %8f", totdens);
322 static void filterdensmap(real ****Densmap, int xslices, int yslices, int zslices, int tblocks, int ftsize)
327 int i, j, k, n, order;
329 std = ((real)order/2.0);
331 snew(kernel, ftsize);
332 gausskernel(kernel, ftsize, var);
333 for (n = 0; n < tblocks; n++)
335 for (i = 0; i < xslices; i++)
337 for (j = 0; j < yslices; j++)
339 periodic_convolution(zslices, Densmap[n][i][j], ftsize, kernel);
348 static void interfaces_txy (real ****Densmap, int xslices, int yslices, int zslices,
349 int tblocks, real binwidth, int method,
350 real dens1, real dens2, t_interf ****intf1,
351 t_interf ****intf2, const output_env_t oenv)
353 /*Returns two pointers to 3D arrays of t_interf structs containing (position,thickness) of the interface(s)*/
355 real *zDensavg; /* zDensavg[z]*/
358 int ndx1, ndx2, deltandx, *zperm;
359 real densmid, densl, densr, alpha, pos, spread;
360 real splitpoint, startpoint, endpoint;
361 real *sigma1, *sigma2;
364 real *fit1 = NULL, *fit2 = NULL;
367 const real onehalf = 1.00/2.00;
368 t_interf ***int1 = NULL, ***int2 = NULL; /*Interface matrices [t][x,y] - last index in row-major order*/
369 /*Create int1(t,xy) and int2(t,xy) arrays with correct number of interf_t elements*/
370 xysize = xslices*yslices;
373 for (i = 0; i < tblocks; i++)
375 snew(int1[i], xysize);
376 snew(int2[i], xysize);
377 for (j = 0; j < xysize; j++)
381 init_interf(int1[i][j]);
382 init_interf(int2[i][j]);
386 if (method == methBISECT)
388 densmid = onehalf*(dens1+dens2);
389 snew(zperm, zslices);
390 for (n = 0; n < tblocks; n++)
392 for (i = 0; i < xslices; i++)
394 for (j = 0; j < yslices; j++)
396 rangeArray(zperm, zslices); /*reset permutation array to identity*/
397 /*Binsearch returns slice-nr where the order param is <= setpoint sgmid*/
398 ndx1 = start_binsearch(Densmap[n][i][j], zperm, 0, zslices/2-1, densmid, 1);
399 ndx2 = start_binsearch(Densmap[n][i][j], zperm, zslices/2, zslices-1, densmid, -1);
401 /* Linear interpolation (for use later if time allows)
402 * rho_1s= Densmap[n][i][j][zperm[ndx1]]
403 * rho_1e =Densmap[n][i][j][zperm[ndx1+1]] - in worst case might be far off
404 * rho_2s =Densmap[n][i][j][zperm[ndx2+1]]
405 * rho_2e =Densmap[n][i][j][zperm[ndx2]]
406 * For 1st interface we have:
407 densl= Densmap[n][i][j][zperm[ndx1]];
408 densr= Densmap[n][i][j][zperm[ndx1+1]];
409 alpha=(densmid-densl)/(densr-densl);
410 deltandx=zperm[ndx1+1]-zperm[ndx1];
413 printf("Alpha, Deltandx %f %i\n", alpha,deltandx);
415 if(abs(alpha)>1.0 || abs(deltandx)>3){
420 pos=zperm[ndx1]+alpha*deltandx;
421 spread=binwidth*deltandx;
423 * For the 2nd interface can use the same formulation, since alpha should become negative ie:
424 * alpha=(densmid-Densmap[n][i][j][zperm[ndx2]])/(Densmap[n][i][j][zperm[nxd2+1]]-Densmap[n][i][j][zperm[ndx2]]);
425 * deltandx=zperm[ndx2+1]-zperm[ndx2];
426 * pos=zperm[ndx2]+alpha*deltandx; */
428 /*After filtering we use the direct approach */
429 int1[n][j+(i*yslices)]->Z = (zperm[ndx1]+onehalf)*binwidth;
430 int1[n][j+(i*yslices)]->t = binwidth;
431 int2[n][j+(i*yslices)]->Z = (zperm[ndx2]+onehalf)*binwidth;
432 int2[n][j+(i*yslices)]->t = binwidth;
438 if (method == methFUNCFIT)
440 /*Assume a box divided in 2 along midpoint of z for starters*/
442 endpoint = binwidth*zslices;
443 splitpoint = (startpoint+endpoint)/2.0;
444 /*Initial fit proposals*/
445 beginfit1[0] = dens1;
446 beginfit1[1] = dens2;
447 beginfit1[2] = (splitpoint/2);
450 beginfit2[0] = dens2;
451 beginfit2[1] = dens1;
452 beginfit2[2] = (3*splitpoint/2);
455 snew(zDensavg, zslices);
456 snew(sigma1, zslices);
457 snew(sigma2, zslices);
459 for (k = 0; k < zslices; k++)
461 sigma1[k] = sigma2[k] = 1;
463 /*Calculate average density along z - avoid smoothing by using coarse-grained-mesh*/
464 for (k = 0; k < zslices; k++)
466 for (n = 0; n < tblocks; n++)
468 for (i = 0; i < xslices; i++)
470 for (j = 0; j < yslices; j++)
472 zDensavg[k] += (Densmap[n][i][j][k]/(xslices*yslices*tblocks));
480 xvg = xvgropen("DensprofileonZ.xvg", "Averaged Densityprofile on Z", "z[nm]", "Density[kg/m^3]", oenv);
481 for (k = 0; k < zslices; k++)
483 fprintf(xvg, "%4f.3 %8f.4\n", k*binwidth, zDensavg[k]);
488 /*Fit average density in z over whole trajectory to obtain tentative fit-parameters in fit1 and fit2*/
490 /*Fit 1st half of box*/
491 do_lmfit(zslices, zDensavg, sigma1, binwidth, NULL, startpoint, splitpoint, oenv, FALSE, effnERF, beginfit1, 3);
492 /*Fit 2nd half of box*/
493 do_lmfit(zslices, zDensavg, sigma2, binwidth, NULL, splitpoint, endpoint, oenv, FALSE, effnERF, beginfit2, 3);
495 /*Initialise the const arrays for storing the average fit parameters*/
501 /*Now do fit over each x y and t slice to get Zint(x,y,t) - loop is very large, we potentially should average over time directly*/
502 for (n = 0; n < tblocks; n++)
504 for (i = 0; i < xslices; i++)
506 for (j = 0; j < yslices; j++)
508 /*Reinitialise fit for each mesh-point*/
511 for (k = 0; k < 4; k++)
513 fit1[k] = avgfit1[k];
514 fit2[k] = avgfit2[k];
516 /*Now fit and store in structures in row-major order int[n][i][j]*/
517 do_lmfit(zslices, Densmap[n][i][j], sigma1, binwidth, NULL, startpoint, splitpoint, oenv, FALSE, effnERF, fit1, 1);
518 int1[n][j+(yslices*i)]->Z = fit1[2];
519 int1[n][j+(yslices*i)]->t = fit1[3];
520 do_lmfit(zslices, Densmap[n][i][j], sigma2, binwidth, NULL, splitpoint, endpoint, oenv, FALSE, effnERF, fit2, 2);
521 int2[n][j+(yslices*i)]->Z = fit2[2];
522 int2[n][j+(yslices*i)]->t = fit2[3];
534 static void writesurftoxpms(t_interf ***surf1, t_interf ***surf2, int tblocks, int xbins, int ybins, int zbins, real bw, real bwz, char **outfiles, int maplevels )
538 real **profile1, **profile2;
539 real max1, max2, min1, min2, *xticks, *yticks;
540 t_rgb lo = {0, 0, 0};
541 t_rgb hi = {1, 1, 1};
542 FILE *xpmfile1, *xpmfile2;
544 /*Prepare xpm structures for output*/
546 /*Allocate memory to tick's and matrices*/
547 snew (xticks, xbins+1);
548 snew (yticks, ybins+1);
550 profile1 = mk_matrix(xbins, ybins, FALSE);
551 profile2 = mk_matrix(xbins, ybins, FALSE);
553 for (i = 0; i < xbins+1; i++)
557 for (j = 0; j < ybins+1; j++)
562 xpmfile1 = ffopen(outfiles[0], "w");
563 xpmfile2 = ffopen(outfiles[1], "w");
566 min1 = min2 = zbins*bwz;
568 for (n = 0; n < tblocks; n++)
570 sprintf(numbuf, "tblock: %4i", n);
571 /*Filling matrices for inclusion in xpm-files*/
572 for (i = 0; i < xbins; i++)
574 for (j = 0; j < ybins; j++)
576 profile1[i][j] = (surf1[n][j+ybins*i])->Z;
577 profile2[i][j] = (surf2[n][j+ybins*i])->Z;
578 /*Finding max and min values*/
579 if (profile1[i][j] > max1)
581 max1 = profile1[i][j];
583 if (profile1[i][j] < min1)
585 min1 = profile1[i][j];
587 if (profile2[i][j] > max2)
589 max2 = profile2[i][j];
591 if (profile2[i][j] < min2)
593 min2 = profile2[i][j];
598 write_xpm(xpmfile1, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile1, min1, max1, lo, hi, &maplevels);
599 write_xpm(xpmfile2, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile2, min2, max2, lo, hi, &maplevels);
612 static void writeraw(t_interf ***int1, t_interf ***int2, int tblocks, int xbins, int ybins, char **fnms)
617 raw1 = ffopen(fnms[0], "w");
618 raw2 = ffopen(fnms[1], "w");
619 fprintf(raw1, "#Produced by: %s #\n", command_line());
620 fprintf(raw2, "#Produced by: %s #\n", command_line());
621 fprintf(raw1, "#Legend: nt nx ny\n#Xbin Ybin Z t\n");
622 fprintf(raw2, "#Legend: nt nx ny\n#Xbin Ybin Z t\n");
623 fprintf(raw1, "%i %i %i\n", tblocks, xbins, ybins);
624 fprintf(raw2, "%i %i %i\n", tblocks, xbins, ybins);
625 for (n = 0; n < tblocks; n++)
627 for (i = 0; i < xbins; i++)
629 for (j = 0; j < ybins; j++)
631 fprintf(raw1, "%i %i %8.5f %6.4f\n", i, j, (int1[n][j+ybins*i])->Z, (int1[n][j+ybins*i])->t);
632 fprintf(raw2, "%i %i %8.5f %6.4f\n", i, j, (int2[n][j+ybins*i])->Z, (int2[n][j+ybins*i])->t);
643 int gmx_densorder(int argc, char *argv[])
645 static const char *desc[] = {
646 "A small program to reduce a two-phase density distribution",
647 "along an axis, computed over a MD trajectory",
648 "to 2D surfaces fluctuating in time, by a fit to",
649 "a functional profile for interfacial densities",
650 "A time-averaged spatial representation of the",
651 "interfaces can be output with the option -tavg"
654 /* Extra arguments - but note how you always get the begin/end
655 * options when running the program, without mentioning them here!
660 char title[STRLEN], **grpname;
662 static real binw = 0.2;
663 static real binwz = 0.05;
664 static real dens1 = 0.00;
665 static real dens2 = 1000.00;
666 static int ftorder = 0;
667 static int nsttblock = 100;
669 static char *axtitle = "Z";
670 atom_id **index; /* Index list for single group*/
671 int xslices, yslices, zslices, tblock;
672 static gmx_bool bGraph = FALSE;
673 static gmx_bool bCenter = FALSE;
674 static gmx_bool bFourier = FALSE;
675 static gmx_bool bRawOut = FALSE;
676 static gmx_bool bOut = FALSE;
677 static gmx_bool b1d = FALSE;
678 static int nlevels = 100;
679 /*Densitymap - Densmap[t][x][y][z]*/
680 real ****Densmap = NULL;
681 /* Surfaces surf[t][surf_x,surf_y]*/
682 t_interf ***surf1, ***surf2;
684 static const char *meth[] = {NULL, "bisect", "functional", NULL};
687 char **graphfiles, **rawfiles, **spectra; /* Filenames for xpm-surface maps, rawdata and powerspectra */
688 int nfxpm = -1, nfraw, nfspect; /* # files for interface maps and spectra = # interfaces */
691 { "-1d", FALSE, etBOOL, {&b1d},
692 "Pseudo-1d interface geometry"},
693 { "-bw", FALSE, etREAL, {&binw},
694 "Binwidth of density distribution tangential to interface"},
695 { "-bwn", FALSE, etREAL, {&binwz},
696 "Binwidth of density distribution normal to interface"},
697 { "-order", FALSE, etINT, {&ftorder},
698 "Order of Gaussian filter, order 0 equates to NO filtering"},
699 {"-axis", FALSE, etSTR, {&axtitle},
700 "Axis Direction - X, Y or Z"},
701 {"-method", FALSE, etENUM, {meth},
702 "Interface location method"},
703 {"-d1", FALSE, etREAL, {&dens1},
704 "Bulk density phase 1 (at small z)"},
705 {"-d2", FALSE, etREAL, {&dens2},
706 "Bulk density phase 2 (at large z)"},
707 { "-tblock", FALSE, etINT, {&nsttblock},
708 "Number of frames in one time-block average"},
709 { "-nlevel", FALSE, etINT, {&nlevels},
710 "Number of Height levels in 2D - XPixMaps"}
715 { efTPX, "-s", NULL, ffREAD }, /* this is for the topology */
716 { efTRX, "-f", NULL, ffREAD }, /* and this for the trajectory */
717 { efNDX, "-n", NULL, ffREAD}, /* this is to select groups */
718 { efDAT, "-o", "Density4D", ffOPTWR}, /* This is for outputting the entire 4D densityfield in binary format */
719 { efOUT, "-or", NULL, ffOPTWRMULT}, /* This is for writing out the entire information in the t_interf arrays */
720 { efXPM, "-og", "interface", ffOPTWRMULT}, /* This is for writing out the interface meshes - one xpm-file per tblock*/
721 { efOUT, "-Spect", "intfspect", ffOPTWRMULT}, /* This is for the trajectory averaged Fourier-spectra*/
724 #define NFILE asize(fnm)
726 /* This is the routine responsible for adding default options,
727 * calling the X/motif interface, etc. */
728 parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
729 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
733 bFourier = opt2bSet("-Spect", NFILE, fnm);
734 bRawOut = opt2bSet("-or", NFILE, fnm);
735 bGraph = opt2bSet("-og", NFILE, fnm);
736 bOut = opt2bSet("-o", NFILE, fnm);
737 top = read_top(ftp2fn(efTPX, NFILE, fnm), &ePBC);
743 axis = toupper(axtitle[0]) - 'X';
745 get_index(&top->atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, ngx, index, grpname);
747 density_in_time(ftp2fn(efTRX, NFILE, fnm), index, ngx, binw, binwz, nsttblock, &Densmap, &xslices, &yslices, &zslices, &tblock, top, ePBC, axis, bCenter, b1d, oenv);
751 filterdensmap(Densmap, xslices, yslices, zslices, tblock, 2*ftorder+1);
756 outputfield(opt2fn("-o", NFILE, fnm), Densmap, xslices, yslices, zslices, tblock);
759 interfaces_txy(Densmap, xslices, yslices, zslices, tblock, binwz, eMeth, dens1, dens2, &surf1, &surf2, oenv);
764 /*Output surface-xpms*/
765 nfxpm = opt2fns(&graphfiles, "-og", NFILE, fnm);
768 gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfxpm);
770 writesurftoxpms(surf1, surf2, tblock, xslices, yslices, zslices, binw, binwz, graphfiles, zslices);
780 nfraw = opt2fns(&rawfiles, "-or", NFILE, fnm);
783 gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfxpm);
785 writeraw(surf1, surf2, tblock, xslices, yslices, rawfiles);
792 nfspect = opt2fns(&spectra, "-Spect", NFILE, fnm);
795 gmx_fatal(FARGS, "No or not correct number (2) of output-file-series: %d",
798 powerspectavg_intf(surf1, surf2, tblock, xslices, yslices, spectra);
802 if (bGraph || bFourier || bRawOut)