Merge release-4-6 into release-5-0
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_densorder.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2011,2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <ctype.h>
40 #include <math.h>
41 #include <string.h>
42
43 #include "sysstuff.h"
44 #include "gromacs/utility/cstringutil.h"
45 #include "typedefs.h"
46 #include "gromacs/utility/smalloc.h"
47 #include "macros.h"
48 #include "gstat.h"
49 #include "vec.h"
50 #include "xvgr.h"
51 #include "pbc.h"
52 #include "gromacs/fileio/futil.h"
53 #include "gromacs/commandline/pargs.h"
54 #include "index.h"
55 #include "gromacs/fileio/tpxio.h"
56 #include "gromacs/fileio/trxio.h"
57 #include "physics.h"
58 #include "gromacs/fileio/matio.h"
59 #include "dens_filter.h"
60 #include "binsearch.h"
61 #include "powerspect.h"
62 #include "gmx_ana.h"
63 #include "copyrite.h"
64
65 #include "gromacs/utility/exceptions.h"
66 #include "gromacs/legacyheaders/gmx_fatal.h"
67 #include "gromacs/utility/programcontext.h"
68
69 #ifdef GMX_DOUBLE
70 #define FLOOR(x) ((int) floor(x))
71 #else
72 #define FLOOR(x) ((int) floorf(x))
73 #endif
74
75 enum {
76     methSEL, methBISECT, methFUNCFIT, methNR
77 };
78
79 static void center_coords(t_atoms *atoms, matrix box, rvec x0[], int axis)
80 {
81     int  i, m;
82     real tmass, mm;
83     rvec com, shift, box_center;
84
85     tmass = 0;
86     clear_rvec(com);
87     for (i = 0; (i < atoms->nr); i++)
88     {
89         mm     = atoms->atom[i].m;
90         tmass += mm;
91         for (m = 0; (m < DIM); m++)
92         {
93             com[m] += mm*x0[i][m];
94         }
95     }
96     for (m = 0; (m < DIM); m++)
97     {
98         com[m] /= tmass;
99     }
100     calc_box_center(ecenterDEF, box, box_center);
101     rvec_sub(box_center, com, shift);
102     shift[axis] -= box_center[axis];
103
104     for (i = 0; (i < atoms->nr); i++)
105     {
106         rvec_dec(x0[i], shift);
107     }
108 }
109
110
111 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)
112
113 {
114 /*
115  * *****Densdevel pointer to array of density values in slices and frame-blocks Densdevel[*nsttblock][*xslices][*yslices][*zslices]
116  * Densslice[x][y][z]
117  * nsttblock - nr of frames in each time-block
118  * bw  widths of normal slices
119  *
120  * axis  - axis direction (normal to slices)
121  * nndx - number ot atoms in **index
122  * grpn  - group number in index
123  */
124     t_trxstatus *status;
125     gmx_rmpbc_t  gpbc = NULL;
126     matrix       box;                    /* Box - 3x3 -each step*/
127     rvec        *x0;                     /* List of Coord without PBC*/
128     int          i, j,                   /* loop indices, checks etc*/
129                  ax1     = 0, ax2 = 0,   /* tangent directions */
130                  framenr = 0,            /* frame number in trajectory*/
131                  slicex, slicey, slicez; /*slice # of x y z position */
132     real ***Densslice = NULL;            /* Density-slice in one frame*/
133     real    dscale;                      /*physical scaling factor*/
134     real    t, x, y, z;                  /* time and coordinates*/
135     rvec    bbww;
136
137     *tblock = 0; /* blocknr in block average - initialise to 0*/
138     /* Axis: X=0, Y=1,Z=2 */
139     switch (axis)
140     {
141         case 0:
142             ax1 = YY; ax2 = ZZ; /*Surface: YZ*/
143             break;
144         case 1:
145             ax1 = ZZ; ax2 = XX; /* Surface : XZ*/
146             break;
147         case 2:
148             ax1 = XX; ax2 = YY; /* Surface XY*/
149             break;
150         default:
151             gmx_fatal(FARGS, "Invalid axes. Terminating\n");
152     }
153
154     if (read_first_x(oenv, &status, fn, &t, &x0, box) == 0)
155     {
156         gmx_fatal(FARGS, "Could not read coordinates from file"); /* Open trajectory for read*/
157
158
159     }
160     *zslices = 1+FLOOR(box[axis][axis]/bwz);
161     *yslices = 1+FLOOR(box[ax2][ax2]/bw);
162     *xslices = 1+FLOOR(box[ax1][ax1]/bw);
163     if (bps1d)
164     {
165         if (*xslices < *yslices)
166         {
167             *xslices = 1;
168         }
169         else
170         {
171             *yslices = 1;
172         }
173     }
174     fprintf(stderr,
175             "\nDividing the box in %5d x %5d x %5d slices with binw %f along axis %d\n", *xslices, *yslices, *zslices, bw, axis );
176
177
178     /****Start trajectory processing***/
179
180     /*Initialize Densdevel and PBC-remove*/
181     gpbc = gmx_rmpbc_init(&top->idef, ePBC, top->atoms.nr);
182
183     *Densdevel = NULL;
184
185     do
186     {
187         bbww[XX] = box[ax1][ax1]/ *xslices;
188         bbww[YY] = box[ax2][ax2]/ *yslices;
189         bbww[ZZ] = box[axis][axis]/ *zslices;
190         gmx_rmpbc(gpbc, top->atoms.nr, box, x0);
191         /*Reset Densslice every nsttblock steps*/
192         /* The first conditional is for clang to understand that this branch is
193          * always taken the first time. */
194         if (Densslice == NULL || framenr % nsttblock == 0)
195         {
196             snew(Densslice, *xslices);
197             for (i = 0; i < *xslices; i++)
198             {
199                 snew(Densslice[i], *yslices);
200                 for (j = 0; j < *yslices; j++)
201                 {
202                     snew(Densslice[i][j], *zslices);
203                 }
204             }
205
206             /* Allocate Memory to  extra frame in Densdevel -  rather stupid approach:
207              * A single frame each time, although only every nsttblock steps.
208              */
209             srenew(*Densdevel, *tblock+1);
210             (*Densdevel)[*tblock] = Densslice;
211         }
212
213         dscale = (*xslices)*(*yslices)*(*zslices)*AMU/ (box[ax1][ax1]*box[ax2][ax2]*box[axis][axis]*nsttblock*(NANO*NANO*NANO));
214
215         if (bCenter)
216         {
217             center_coords(&top->atoms, box, x0, axis);
218         }
219
220
221         for (j = 0; j < gnx[0]; j++)
222         {   /*Loop over all atoms in selected index*/
223             x = x0[index[0][j]][ax1];
224             y = x0[index[0][j]][ax2];
225             z = x0[index[0][j]][axis];
226             while (x < 0)
227             {
228                 x += box[ax1][ax1];
229             }
230             while (x > box[ax1][ax1])
231             {
232                 x -= box[ax1][ax1];
233             }
234
235             while (y < 0)
236             {
237                 y += box[ax2][ax2];
238             }
239             while (y > box[ax2][ax2])
240             {
241                 y -= box[ax2][ax2];
242             }
243
244             while (z < 0)
245             {
246                 z += box[axis][axis];
247             }
248             while (z > box[axis][axis])
249             {
250                 z -= box[axis][axis];
251             }
252
253             slicex = ((int) (x/bbww[XX])) % *xslices;
254             slicey = ((int) (y/bbww[YY])) % *yslices;
255             slicez = ((int) (z/bbww[ZZ])) % *zslices;
256             Densslice[slicex][slicey][slicez] += (top->atoms.atom[index[0][j]].m*dscale);
257         }
258
259         framenr++;
260
261         if (framenr % nsttblock == 0)
262         {
263             /*Implicit incrementation of Densdevel via renewal of Densslice*/
264             /*only every nsttblock steps*/
265             (*tblock)++;
266         }
267
268     }
269     while (read_next_x(oenv, status, &t, x0, box));
270
271
272     /*Free memory we no longer need and exit.*/
273     gmx_rmpbc_done(gpbc);
274     close_trj(status);
275
276     if (0)
277     {
278         FILE *fp;
279         fp = fopen("koko.xvg", "w");
280         for (j = 0; (j < *zslices); j++)
281         {
282             fprintf(fp, "%5d", j);
283             for (i = 0; (i < *tblock); i++)
284             {
285                 fprintf(fp, "  %10g", (*Densdevel)[i][9][1][j]);
286             }
287             fprintf(fp, "\n");
288         }
289         fclose(fp);
290     }
291
292 }
293
294 static void outputfield(const char *fldfn, real ****Densmap,
295                         int xslices, int yslices, int zslices, int tdim)
296 {
297 /*Debug-filename and filehandle*/
298     FILE *fldH;
299     int   n, i, j, k;
300     int   dim[4];
301     real  totdens = 0;
302
303     dim[0] = tdim;
304     dim[1] = xslices;
305     dim[2] = yslices;
306     dim[3] = zslices;
307
308     fldH = gmx_ffopen(fldfn, "w");
309     fwrite(dim, sizeof(int), 4, fldH);
310     for (n = 0; n < tdim; n++)
311     {
312         for (i = 0; i < xslices; i++)
313         {
314             for (j = 0; j < yslices; j++)
315             {
316                 for (k = 0; k < zslices; k++)
317                 {
318                     fwrite(&(Densmap[n][i][j][k]), sizeof(real), 1, fldH);
319                     totdens += (Densmap[n][i][j][k]);
320                 }
321             }
322         }
323     }
324     totdens /= (xslices*yslices*zslices*tdim);
325     fprintf(stderr, "Total density [kg/m^3]  %8f", totdens);
326     gmx_ffclose(fldH);
327 }
328
329 static void filterdensmap(real ****Densmap, int xslices, int yslices, int zslices, int tblocks, int ftsize)
330 {
331     real *kernel;
332     real  std, var;
333     int   i, j, n, order;
334     order = ftsize/2;
335     std   = ((real)order/2.0);
336     var   = std*std;
337     snew(kernel, ftsize);
338     gausskernel(kernel, ftsize, var);
339     for (n = 0; n < tblocks; n++)
340     {
341         for (i = 0; i < xslices; i++)
342         {
343             for (j = 0; j < yslices; j++)
344             {
345                 periodic_convolution(zslices, Densmap[n][i][j], ftsize, kernel);
346             }
347         }
348     }
349 }
350
351
352
353
354 static void interfaces_txy (real ****Densmap, int xslices, int yslices, int zslices,
355                             int tblocks, real binwidth, int method,
356                             real dens1, real dens2, t_interf ****intf1,
357                             t_interf ****intf2, const output_env_t oenv)
358 {
359     /*Returns two pointers to 3D arrays of t_interf structs containing (position,thickness) of the interface(s)*/
360     FILE       *xvg;
361     real       *zDensavg; /* zDensavg[z]*/
362     int         i, j, k, n;
363     int         xysize;
364     int         ndx1, ndx2, *zperm;
365     real        densmid;
366     real        splitpoint, startpoint, endpoint;
367     real       *sigma1, *sigma2;
368     real        beginfit1[4];
369     real        beginfit2[4];
370     real       *fit1 = NULL, *fit2 = NULL;
371     const real *avgfit1;
372     const real *avgfit2;
373     const real  onehalf = 1.00/2.00;
374     t_interf ***int1    = NULL, ***int2 = NULL; /*Interface matrices [t][x,y] - last index in row-major order*/
375     /*Create int1(t,xy) and int2(t,xy) arrays with correct number of interf_t elements*/
376     xysize = xslices*yslices;
377     snew(int1, tblocks);
378     snew(int2, tblocks);
379     for (i = 0; i < tblocks; i++)
380     {
381         snew(int1[i], xysize);
382         snew(int2[i], xysize);
383         for (j = 0; j < xysize; j++)
384         {
385             snew(int1[i][j], 1);
386             snew(int2[i][j], 1);
387             init_interf(int1[i][j]);
388             init_interf(int2[i][j]);
389         }
390     }
391
392     if (method == methBISECT)
393     {
394         densmid = onehalf*(dens1+dens2);
395         snew(zperm, zslices);
396         for (n = 0; n < tblocks; n++)
397         {
398             for (i = 0; i < xslices; i++)
399             {
400                 for (j = 0; j < yslices; j++)
401                 {
402                     rangeArray(zperm, zslices); /*reset permutation array to identity*/
403                     /*Binsearch returns slice-nr where the order param is  <= setpoint sgmid*/
404                     ndx1 = start_binsearch(Densmap[n][i][j], zperm, 0, zslices/2-1, densmid, 1);
405                     ndx2 = start_binsearch(Densmap[n][i][j], zperm, zslices/2, zslices-1, densmid, -1);
406
407                     /* Linear interpolation (for use later if time allows)
408                      * rho_1s= Densmap[n][i][j][zperm[ndx1]]
409                      * rho_1e =Densmap[n][i][j][zperm[ndx1+1]] - in worst case might be far off
410                      * rho_2s =Densmap[n][i][j][zperm[ndx2+1]]
411                      * rho_2e =Densmap[n][i][j][zperm[ndx2]]
412                      * For 1st interface we have:
413                        densl= Densmap[n][i][j][zperm[ndx1]];
414                        densr= Densmap[n][i][j][zperm[ndx1+1]];
415                        alpha=(densmid-densl)/(densr-densl);
416                        deltandx=zperm[ndx1+1]-zperm[ndx1];
417
418                        if(debug){
419                        printf("Alpha, Deltandx  %f %i\n", alpha,deltandx);
420                        }
421                        if(abs(alpha)>1.0 || abs(deltandx)>3){
422                        pos=zperm[ndx1];
423                        spread=-1;
424                        }
425                        else {
426                        pos=zperm[ndx1]+alpha*deltandx;
427                        spread=binwidth*deltandx;
428                        }
429                      * For the 2nd interface  can use the same formulation, since alpha should become negative ie:
430                      * alpha=(densmid-Densmap[n][i][j][zperm[ndx2]])/(Densmap[n][i][j][zperm[nxd2+1]]-Densmap[n][i][j][zperm[ndx2]]);
431                      * deltandx=zperm[ndx2+1]-zperm[ndx2];
432                      * pos=zperm[ndx2]+alpha*deltandx;   */
433
434                     /*After filtering we use the direct approach        */
435                     int1[n][j+(i*yslices)]->Z = (zperm[ndx1]+onehalf)*binwidth;
436                     int1[n][j+(i*yslices)]->t = binwidth;
437                     int2[n][j+(i*yslices)]->Z = (zperm[ndx2]+onehalf)*binwidth;
438                     int2[n][j+(i*yslices)]->t = binwidth;
439                 }
440             }
441         }
442     }
443
444     if (method == methFUNCFIT)
445     {
446         /*Assume a box divided in 2 along midpoint of z for starters*/
447         startpoint = 0.0;
448         endpoint   = binwidth*zslices;
449         splitpoint = (startpoint+endpoint)/2.0;
450         /*Initial fit proposals*/
451         beginfit1[0] = dens1;
452         beginfit1[1] = dens2;
453         beginfit1[2] = (splitpoint/2);
454         beginfit1[3] = 0.5;
455
456         beginfit2[0] = dens2;
457         beginfit2[1] = dens1;
458         beginfit2[2] = (3*splitpoint/2);
459         beginfit2[3] = 0.5;
460
461         snew(zDensavg, zslices);
462         snew(sigma1, zslices);
463         snew(sigma2, zslices);
464
465         for (k = 0; k < zslices; k++)
466         {
467             sigma1[k] = sigma2[k] = 1;
468         }
469         /*Calculate average density along z - avoid smoothing by using coarse-grained-mesh*/
470         for (k = 0; k < zslices; k++)
471         {
472             for (n = 0; n < tblocks; n++)
473             {
474                 for (i = 0; i < xslices; i++)
475                 {
476                     for (j = 0; j < yslices; j++)
477                     {
478                         zDensavg[k] += (Densmap[n][i][j][k]/(xslices*yslices*tblocks));
479                     }
480                 }
481             }
482         }
483
484         if (debug)
485         {
486             xvg = xvgropen("DensprofileonZ.xvg", "Averaged Densityprofile on Z", "z[nm]", "Density[kg/m^3]", oenv);
487             for (k = 0; k < zslices; k++)
488             {
489                 fprintf(xvg, "%4f.3   %8f.4\n", k*binwidth, zDensavg[k]);
490             }
491             xvgrclose(xvg);
492         }
493
494         /*Fit average density in z over whole trajectory to obtain tentative fit-parameters in fit1 and fit2*/
495
496         /*Fit 1st half of box*/
497         do_lmfit(zslices, zDensavg, sigma1, binwidth, NULL, startpoint, splitpoint, oenv, FALSE, effnERF, beginfit1, 3);
498         /*Fit 2nd half of box*/
499         do_lmfit(zslices, zDensavg, sigma2, binwidth, NULL, splitpoint, endpoint, oenv, FALSE, effnERF, beginfit2, 3);
500
501         /*Initialise the const arrays for storing the average fit parameters*/
502         avgfit1 = beginfit1;
503         avgfit2 = beginfit2;
504
505
506
507         /*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*/
508         for (n = 0; n < tblocks; n++)
509         {
510             for (i = 0; i < xslices; i++)
511             {
512                 for (j = 0; j < yslices; j++)
513                 {
514                     /*Reinitialise fit for each mesh-point*/
515                     srenew(fit1, 4);
516                     srenew(fit2, 4);
517                     for (k = 0; k < 4; k++)
518                     {
519                         fit1[k] = avgfit1[k];
520                         fit2[k] = avgfit2[k];
521                     }
522                     /*Now fit and store in structures in row-major order int[n][i][j]*/
523                     do_lmfit(zslices, Densmap[n][i][j], sigma1, binwidth, NULL, startpoint, splitpoint, oenv, FALSE, effnERF, fit1, 1);
524                     int1[n][j+(yslices*i)]->Z = fit1[2];
525                     int1[n][j+(yslices*i)]->t = fit1[3];
526                     do_lmfit(zslices, Densmap[n][i][j], sigma2, binwidth, NULL, splitpoint, endpoint, oenv, FALSE, effnERF, fit2, 2);
527                     int2[n][j+(yslices*i)]->Z = fit2[2];
528                     int2[n][j+(yslices*i)]->t = fit2[3];
529                 }
530             }
531         }
532     }
533
534
535     *intf1 = int1;
536     *intf2 = int2;
537
538 }
539
540 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 )
541 {
542     char   numbuf[13];
543     int    n, i, j;
544     real **profile1, **profile2;
545     real   max1, max2, min1, min2, *xticks, *yticks;
546     t_rgb  lo = {0, 0, 0};
547     t_rgb  hi = {1, 1, 1};
548     FILE  *xpmfile1, *xpmfile2;
549
550 /*Prepare xpm structures for output*/
551
552 /*Allocate memory to tick's and matrices*/
553     snew (xticks, xbins+1);
554     snew (yticks, ybins+1);
555
556     profile1 = mk_matrix(xbins, ybins, FALSE);
557     profile2 = mk_matrix(xbins, ybins, FALSE);
558
559     for (i = 0; i < xbins+1; i++)
560     {
561         xticks[i] += bw;
562     }
563     for (j = 0; j < ybins+1; j++)
564     {
565         yticks[j] += bw;
566     }
567
568     xpmfile1 = gmx_ffopen(outfiles[0], "w");
569     xpmfile2 = gmx_ffopen(outfiles[1], "w");
570
571     max1 = max2 = 0.0;
572     min1 = min2 = zbins*bwz;
573
574     for (n = 0; n < tblocks; n++)
575     {
576         sprintf(numbuf, "tblock: %4i", n);
577 /*Filling matrices for inclusion in xpm-files*/
578         for (i = 0; i < xbins; i++)
579         {
580             for (j = 0; j < ybins; j++)
581             {
582                 profile1[i][j] = (surf1[n][j+ybins*i])->Z;
583                 profile2[i][j] = (surf2[n][j+ybins*i])->Z;
584                 /*Finding max and min values*/
585                 if (profile1[i][j] > max1)
586                 {
587                     max1 = profile1[i][j];
588                 }
589                 if (profile1[i][j] < min1)
590                 {
591                     min1 = profile1[i][j];
592                 }
593                 if (profile2[i][j] > max2)
594                 {
595                     max2 = profile2[i][j];
596                 }
597                 if (profile2[i][j] < min2)
598                 {
599                     min2 = profile2[i][j];
600                 }
601             }
602         }
603
604         write_xpm(xpmfile1, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile1, min1, max1, lo, hi, &maplevels);
605         write_xpm(xpmfile2, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile2, min2, max2, lo, hi, &maplevels);
606     }
607
608     gmx_ffclose(xpmfile1);
609     gmx_ffclose(xpmfile2);
610
611
612     sfree(profile1);
613     sfree(profile2);
614     sfree(xticks);
615     sfree(yticks);
616 }
617
618 static void writeraw(t_interf ***int1, t_interf ***int2, int tblocks, int xbins, int ybins, char **fnms)
619 {
620     FILE *raw1, *raw2;
621     int   i, j, n;
622
623     raw1 = gmx_ffopen(fnms[0], "w");
624     raw2 = gmx_ffopen(fnms[1], "w");
625     try
626     {
627         gmx::BinaryInformationSettings settings;
628         settings.generatedByHeader(true);
629         settings.linePrefix("# ");
630         gmx::printBinaryInformation(raw1, gmx::getProgramContext(), settings);
631         gmx::printBinaryInformation(raw2, gmx::getProgramContext(), settings);
632     }
633     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
634     fprintf(raw1, "# Legend: nt nx ny\n# Xbin Ybin Z t\n");
635     fprintf(raw2, "# Legend: nt nx ny\n# Xbin Ybin Z t\n");
636     fprintf(raw1, "%i %i %i\n", tblocks, xbins, ybins);
637     fprintf(raw2, "%i %i %i\n", tblocks, xbins, ybins);
638     for (n = 0; n < tblocks; n++)
639     {
640         for (i = 0; i < xbins; i++)
641         {
642             for (j = 0; j < ybins; j++)
643             {
644                 fprintf(raw1, "%i  %i  %8.5f  %6.4f\n", i, j, (int1[n][j+ybins*i])->Z, (int1[n][j+ybins*i])->t);
645                 fprintf(raw2, "%i  %i  %8.5f  %6.4f\n", i, j, (int2[n][j+ybins*i])->Z, (int2[n][j+ybins*i])->t);
646             }
647         }
648     }
649
650     gmx_ffclose(raw1);
651     gmx_ffclose(raw2);
652 }
653
654
655
656 int gmx_densorder(int argc, char *argv[])
657 {
658     static const char *desc[] = {
659         "[THISMODULE] reduces a two-phase density distribution",
660         "along an axis, computed over a MD trajectory,",
661         "to 2D surfaces fluctuating in time, by a fit to",
662         "a functional profile for interfacial densities.",
663         "A time-averaged spatial representation of the",
664         "interfaces can be output with the option [TT]-tavg[tt]."
665     };
666
667     /* Extra arguments - but note how you always get the begin/end
668      * options when running the program, without mentioning them here!
669      */
670
671     output_env_t       oenv;
672     t_topology        *top;
673     char             **grpname;
674     int                ePBC, *ngx;
675     static real        binw      = 0.2;
676     static real        binwz     = 0.05;
677     static real        dens1     = 0.00;
678     static real        dens2     = 1000.00;
679     static int         ftorder   = 0;
680     static int         nsttblock = 100;
681     static int         axis      = 2;
682     static const char *axtitle   = "Z";
683     atom_id          **index; /* Index list for single group*/
684     int                xslices, yslices, zslices, tblock;
685     static gmx_bool    bGraph   = FALSE;
686     static gmx_bool    bCenter  = FALSE;
687     static gmx_bool    bFourier = FALSE;
688     static gmx_bool    bRawOut  = FALSE;
689     static gmx_bool    bOut     = FALSE;
690     static gmx_bool    b1d      = FALSE;
691     static int         nlevels  = 100;
692     /*Densitymap - Densmap[t][x][y][z]*/
693     real           ****Densmap = NULL;
694     /* Surfaces surf[t][surf_x,surf_y]*/
695     t_interf        ***surf1, ***surf2;
696
697     static const char *meth[] = {NULL, "bisect", "functional", NULL};
698     int                eMeth;
699
700     char             **graphfiles, **rawfiles, **spectra; /* Filenames for xpm-surface maps, rawdata and powerspectra */
701     int                nfxpm = -1, nfraw, nfspect;        /* # files for interface maps and spectra = # interfaces */
702
703     t_pargs            pa[] = {
704         { "-1d", FALSE, etBOOL, {&b1d},
705           "Pseudo-1d interface geometry"},
706         { "-bw", FALSE, etREAL, {&binw},
707           "Binwidth of density distribution tangential to interface"},
708         { "-bwn", FALSE, etREAL, {&binwz},
709           "Binwidth of density distribution normal to interface"},
710         { "-order", FALSE, etINT, {&ftorder},
711           "Order of Gaussian filter, order 0 equates to NO filtering"},
712         {"-axis", FALSE, etSTR, {&axtitle},
713          "Axis Direction - X, Y or Z"},
714         {"-method", FALSE, etENUM, {meth},
715          "Interface location method"},
716         {"-d1", FALSE, etREAL, {&dens1},
717          "Bulk density phase 1 (at small z)"},
718         {"-d2", FALSE, etREAL, {&dens2},
719          "Bulk density phase 2 (at large z)"},
720         { "-tblock", FALSE, etINT, {&nsttblock},
721           "Number of frames in one time-block average"},
722         { "-nlevel", FALSE, etINT, {&nlevels},
723           "Number of Height levels in 2D - XPixMaps"}
724     };
725
726
727     t_filenm fnm[] = {
728         { efTPX, "-s",  NULL, ffREAD },               /* this is for the topology */
729         { efTRX, "-f", NULL, ffREAD },                /* and this for the trajectory */
730         { efNDX, "-n", NULL, ffREAD},                 /* this is to select groups */
731         { efDAT, "-o", "Density4D", ffOPTWR},         /* This is for outputting the entire 4D densityfield in binary format */
732         { efOUT, "-or", NULL, ffOPTWRMULT},           /* This is for writing out the entire information in the t_interf arrays */
733         { efXPM, "-og", "interface", ffOPTWRMULT},    /* This is for writing out the interface meshes - one xpm-file per tblock*/
734         { efOUT, "-Spect", "intfspect", ffOPTWRMULT}, /* This is for the trajectory averaged Fourier-spectra*/
735     };
736
737 #define NFILE asize(fnm)
738
739     /* This is the routine responsible for adding default options,
740      * calling the X/motif interface, etc. */
741     if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
742                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
743     {
744         return 0;
745     }
746
747
748     eMeth    = nenum(meth);
749     bFourier = opt2bSet("-Spect", NFILE, fnm);
750     bRawOut  = opt2bSet("-or", NFILE, fnm);
751     bGraph   = opt2bSet("-og", NFILE, fnm);
752     bOut     = opt2bSet("-o", NFILE, fnm);
753     top      = read_top(ftp2fn(efTPX, NFILE, fnm), &ePBC);
754     snew(grpname, 1);
755     snew(index, 1);
756     snew(ngx, 1);
757
758 /* Calculate axis */
759     axis = toupper(axtitle[0]) - 'X';
760
761     get_index(&top->atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, ngx, index, grpname);
762
763     density_in_time(ftp2fn(efTRX, NFILE, fnm), index, ngx, binw, binwz, nsttblock, &Densmap, &xslices, &yslices, &zslices, &tblock, top, ePBC, axis, bCenter, b1d, oenv);
764
765     if (ftorder > 0)
766     {
767         filterdensmap(Densmap, xslices, yslices, zslices, tblock, 2*ftorder+1);
768     }
769
770     if (bOut)
771     {
772         outputfield(opt2fn("-o", NFILE, fnm), Densmap, xslices, yslices, zslices, tblock);
773     }
774
775     interfaces_txy(Densmap, xslices, yslices, zslices, tblock, binwz, eMeth, dens1, dens2, &surf1, &surf2, oenv);
776
777     if (bGraph)
778     {
779
780         /*Output surface-xpms*/
781         nfxpm = opt2fns(&graphfiles, "-og", NFILE, fnm);
782         if (nfxpm != 2)
783         {
784             gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfxpm);
785         }
786         writesurftoxpms(surf1, surf2, tblock, xslices, yslices, zslices, binw, binwz, graphfiles, zslices);
787     }
788
789
790
791
792
793 /*Output raw-data*/
794     if (bRawOut)
795     {
796         nfraw = opt2fns(&rawfiles, "-or", NFILE, fnm);
797         if (nfraw != 2)
798         {
799             gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfxpm);
800         }
801         writeraw(surf1, surf2, tblock, xslices, yslices, rawfiles);
802     }
803
804
805
806     if (bFourier)
807     {
808         nfspect = opt2fns(&spectra, "-Spect", NFILE, fnm);
809         if (nfspect != 2)
810         {
811             gmx_fatal(FARGS, "No or not correct number (2) of output-file-series: %d",
812                       nfspect);
813         }
814         powerspectavg_intf(surf1, surf2, tblock, xslices, yslices, spectra);
815     }
816
817     sfree(Densmap);
818     if (bGraph || bFourier || bRawOut)
819     {
820         sfree(surf1);
821         sfree(surf2);
822     }
823
824     return 0;
825 }