Sort all includes in src/gromacs
[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 #include "gmxpre.h"
36
37 #include <ctype.h>
38 #include <math.h>
39 #include <string.h>
40
41 #include "gromacs/commandline/pargs.h"
42 #include "gromacs/fileio/matio.h"
43 #include "gromacs/fileio/tpxio.h"
44 #include "gromacs/fileio/trxio.h"
45 #include "gromacs/fileio/xvgr.h"
46 #include "gromacs/gmxana/binsearch.h"
47 #include "gromacs/gmxana/dens_filter.h"
48 #include "gromacs/gmxana/gmx_ana.h"
49 #include "gromacs/gmxana/gstat.h"
50 #include "gromacs/gmxana/powerspect.h"
51 #include "gromacs/legacyheaders/copyrite.h"
52 #include "gromacs/legacyheaders/macros.h"
53 #include "gromacs/legacyheaders/typedefs.h"
54 #include "gromacs/math/units.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/pbcutil/pbc.h"
57 #include "gromacs/pbcutil/rmpbc.h"
58 #include "gromacs/topology/index.h"
59 #include "gromacs/utility/cstringutil.h"
60 #include "gromacs/utility/exceptions.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/futil.h"
63 #include "gromacs/utility/smalloc.h"
64
65 #ifdef GMX_DOUBLE
66 #define FLOOR(x) ((int) floor(x))
67 #else
68 #define FLOOR(x) ((int) floorf(x))
69 #endif
70
71 enum {
72     methSEL, methBISECT, methFUNCFIT, methNR
73 };
74
75 static void center_coords(t_atoms *atoms, matrix box, rvec x0[], int axis)
76 {
77     int  i, m;
78     real tmass, mm;
79     rvec com, shift, box_center;
80
81     tmass = 0;
82     clear_rvec(com);
83     for (i = 0; (i < atoms->nr); i++)
84     {
85         mm     = atoms->atom[i].m;
86         tmass += mm;
87         for (m = 0; (m < DIM); m++)
88         {
89             com[m] += mm*x0[i][m];
90         }
91     }
92     for (m = 0; (m < DIM); m++)
93     {
94         com[m] /= tmass;
95     }
96     calc_box_center(ecenterDEF, box, box_center);
97     rvec_sub(box_center, com, shift);
98     shift[axis] -= box_center[axis];
99
100     for (i = 0; (i < atoms->nr); i++)
101     {
102         rvec_dec(x0[i], shift);
103     }
104 }
105
106
107 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
109 {
110 /*
111  * *****Densdevel pointer to array of density values in slices and frame-blocks Densdevel[*nsttblock][*xslices][*yslices][*zslices]
112  * Densslice[x][y][z]
113  * nsttblock - nr of frames in each time-block
114  * bw  widths of normal slices
115  *
116  * axis  - axis direction (normal to slices)
117  * nndx - number ot atoms in **index
118  * grpn  - group number in index
119  */
120     t_trxstatus *status;
121     gmx_rmpbc_t  gpbc = NULL;
122     matrix       box;                    /* Box - 3x3 -each step*/
123     rvec        *x0;                     /* List of Coord without PBC*/
124     int          i, j,                   /* loop indices, checks etc*/
125                  ax1     = 0, ax2 = 0,   /* tangent directions */
126                  framenr = 0,            /* frame number in trajectory*/
127                  slicex, slicey, slicez; /*slice # of x y z position */
128     real ***Densslice = NULL;            /* Density-slice in one frame*/
129     real    dscale;                      /*physical scaling factor*/
130     real    t, x, y, z;                  /* time and coordinates*/
131     rvec    bbww;
132
133     *tblock = 0; /* blocknr in block average - initialise to 0*/
134     /* Axis: X=0, Y=1,Z=2 */
135     switch (axis)
136     {
137         case 0:
138             ax1 = YY; ax2 = ZZ; /*Surface: YZ*/
139             break;
140         case 1:
141             ax1 = ZZ; ax2 = XX; /* Surface : XZ*/
142             break;
143         case 2:
144             ax1 = XX; ax2 = YY; /* Surface XY*/
145             break;
146         default:
147             gmx_fatal(FARGS, "Invalid axes. Terminating\n");
148     }
149
150     if (read_first_x(oenv, &status, fn, &t, &x0, box) == 0)
151     {
152         gmx_fatal(FARGS, "Could not read coordinates from file"); /* Open trajectory for read*/
153
154
155     }
156     *zslices = 1+FLOOR(box[axis][axis]/bwz);
157     *yslices = 1+FLOOR(box[ax2][ax2]/bw);
158     *xslices = 1+FLOOR(box[ax1][ax1]/bw);
159     if (bps1d)
160     {
161         if (*xslices < *yslices)
162         {
163             *xslices = 1;
164         }
165         else
166         {
167             *yslices = 1;
168         }
169     }
170     fprintf(stderr,
171             "\nDividing the box in %5d x %5d x %5d slices with binw %f along axis %d\n", *xslices, *yslices, *zslices, bw, axis );
172
173
174     /****Start trajectory processing***/
175
176     /*Initialize Densdevel and PBC-remove*/
177     gpbc = gmx_rmpbc_init(&top->idef, ePBC, top->atoms.nr);
178
179     *Densdevel = NULL;
180
181     do
182     {
183         bbww[XX] = box[ax1][ax1]/ *xslices;
184         bbww[YY] = box[ax2][ax2]/ *yslices;
185         bbww[ZZ] = box[axis][axis]/ *zslices;
186         gmx_rmpbc(gpbc, top->atoms.nr, box, x0);
187         /*Reset Densslice every nsttblock steps*/
188         /* The first conditional is for clang to understand that this branch is
189          * always taken the first time. */
190         if (Densslice == NULL || framenr % nsttblock == 0)
191         {
192             snew(Densslice, *xslices);
193             for (i = 0; i < *xslices; i++)
194             {
195                 snew(Densslice[i], *yslices);
196                 for (j = 0; j < *yslices; j++)
197                 {
198                     snew(Densslice[i][j], *zslices);
199                 }
200             }
201
202             /* Allocate Memory to  extra frame in Densdevel -  rather stupid approach:
203              * A single frame each time, although only every nsttblock steps.
204              */
205             srenew(*Densdevel, *tblock+1);
206             (*Densdevel)[*tblock] = Densslice;
207         }
208
209         dscale = (*xslices)*(*yslices)*(*zslices)*AMU/ (box[ax1][ax1]*box[ax2][ax2]*box[axis][axis]*nsttblock*(NANO*NANO*NANO));
210
211         if (bCenter)
212         {
213             center_coords(&top->atoms, box, x0, axis);
214         }
215
216
217         for (j = 0; j < gnx[0]; j++)
218         {   /*Loop over all atoms in selected index*/
219             x = x0[index[0][j]][ax1];
220             y = x0[index[0][j]][ax2];
221             z = x0[index[0][j]][axis];
222             while (x < 0)
223             {
224                 x += box[ax1][ax1];
225             }
226             while (x > box[ax1][ax1])
227             {
228                 x -= box[ax1][ax1];
229             }
230
231             while (y < 0)
232             {
233                 y += box[ax2][ax2];
234             }
235             while (y > box[ax2][ax2])
236             {
237                 y -= box[ax2][ax2];
238             }
239
240             while (z < 0)
241             {
242                 z += box[axis][axis];
243             }
244             while (z > box[axis][axis])
245             {
246                 z -= box[axis][axis];
247             }
248
249             slicex = ((int) (x/bbww[XX])) % *xslices;
250             slicey = ((int) (y/bbww[YY])) % *yslices;
251             slicez = ((int) (z/bbww[ZZ])) % *zslices;
252             Densslice[slicex][slicey][slicez] += (top->atoms.atom[index[0][j]].m*dscale);
253         }
254
255         framenr++;
256
257         if (framenr % nsttblock == 0)
258         {
259             /*Implicit incrementation of Densdevel via renewal of Densslice*/
260             /*only every nsttblock steps*/
261             (*tblock)++;
262         }
263
264     }
265     while (read_next_x(oenv, status, &t, x0, box));
266
267
268     /*Free memory we no longer need and exit.*/
269     gmx_rmpbc_done(gpbc);
270     close_trj(status);
271
272     if (0)
273     {
274         FILE *fp;
275         fp = fopen("koko.xvg", "w");
276         for (j = 0; (j < *zslices); j++)
277         {
278             fprintf(fp, "%5d", j);
279             for (i = 0; (i < *tblock); i++)
280             {
281                 fprintf(fp, "  %10g", (*Densdevel)[i][9][1][j]);
282             }
283             fprintf(fp, "\n");
284         }
285         fclose(fp);
286     }
287
288 }
289
290 static void outputfield(const char *fldfn, real ****Densmap,
291                         int xslices, int yslices, int zslices, int tdim)
292 {
293 /*Debug-filename and filehandle*/
294     FILE *fldH;
295     int   n, i, j, k;
296     int   dim[4];
297     real  totdens = 0;
298
299     dim[0] = tdim;
300     dim[1] = xslices;
301     dim[2] = yslices;
302     dim[3] = zslices;
303
304     fldH = gmx_ffopen(fldfn, "w");
305     fwrite(dim, sizeof(int), 4, fldH);
306     for (n = 0; n < tdim; n++)
307     {
308         for (i = 0; i < xslices; i++)
309         {
310             for (j = 0; j < yslices; j++)
311             {
312                 for (k = 0; k < zslices; k++)
313                 {
314                     fwrite(&(Densmap[n][i][j][k]), sizeof(real), 1, fldH);
315                     totdens += (Densmap[n][i][j][k]);
316                 }
317             }
318         }
319     }
320     totdens /= (xslices*yslices*zslices*tdim);
321     fprintf(stderr, "Total density [kg/m^3]  %8f", totdens);
322     gmx_ffclose(fldH);
323 }
324
325 static void filterdensmap(real ****Densmap, int xslices, int yslices, int zslices, int tblocks, int ftsize)
326 {
327     real *kernel;
328     real  std, var;
329     int   i, j, n, order;
330     order = ftsize/2;
331     std   = ((real)order/2.0);
332     var   = std*std;
333     snew(kernel, ftsize);
334     gausskernel(kernel, ftsize, var);
335     for (n = 0; n < tblocks; n++)
336     {
337         for (i = 0; i < xslices; i++)
338         {
339             for (j = 0; j < yslices; j++)
340             {
341                 periodic_convolution(zslices, Densmap[n][i][j], ftsize, kernel);
342             }
343         }
344     }
345 }
346
347
348
349
350 static void interfaces_txy (real ****Densmap, int xslices, int yslices, int zslices,
351                             int tblocks, real binwidth, int method,
352                             real dens1, real dens2, t_interf ****intf1,
353                             t_interf ****intf2, const output_env_t oenv)
354 {
355     /*Returns two pointers to 3D arrays of t_interf structs containing (position,thickness) of the interface(s)*/
356     FILE       *xvg;
357     real       *zDensavg; /* zDensavg[z]*/
358     int         i, j, k, n;
359     int         xysize;
360     int         ndx1, ndx2, *zperm;
361     real        densmid;
362     real        splitpoint, startpoint, endpoint;
363     real       *sigma1, *sigma2;
364     real        beginfit1[4];
365     real        beginfit2[4];
366     real       *fit1 = NULL, *fit2 = NULL;
367     const real *avgfit1;
368     const real *avgfit2;
369     const real  onehalf = 1.00/2.00;
370     t_interf ***int1    = NULL, ***int2 = NULL; /*Interface matrices [t][x,y] - last index in row-major order*/
371     /*Create int1(t,xy) and int2(t,xy) arrays with correct number of interf_t elements*/
372     xysize = xslices*yslices;
373     snew(int1, tblocks);
374     snew(int2, tblocks);
375     for (i = 0; i < tblocks; i++)
376     {
377         snew(int1[i], xysize);
378         snew(int2[i], xysize);
379         for (j = 0; j < xysize; j++)
380         {
381             snew(int1[i][j], 1);
382             snew(int2[i][j], 1);
383             init_interf(int1[i][j]);
384             init_interf(int2[i][j]);
385         }
386     }
387
388     if (method == methBISECT)
389     {
390         densmid = onehalf*(dens1+dens2);
391         snew(zperm, zslices);
392         for (n = 0; n < tblocks; n++)
393         {
394             for (i = 0; i < xslices; i++)
395             {
396                 for (j = 0; j < yslices; j++)
397                 {
398                     rangeArray(zperm, zslices); /*reset permutation array to identity*/
399                     /*Binsearch returns slice-nr where the order param is  <= setpoint sgmid*/
400                     ndx1 = start_binsearch(Densmap[n][i][j], zperm, 0, zslices/2-1, densmid, 1);
401                     ndx2 = start_binsearch(Densmap[n][i][j], zperm, zslices/2, zslices-1, densmid, -1);
402
403                     /* Linear interpolation (for use later if time allows)
404                      * rho_1s= Densmap[n][i][j][zperm[ndx1]]
405                      * rho_1e =Densmap[n][i][j][zperm[ndx1+1]] - in worst case might be far off
406                      * rho_2s =Densmap[n][i][j][zperm[ndx2+1]]
407                      * rho_2e =Densmap[n][i][j][zperm[ndx2]]
408                      * For 1st interface we have:
409                        densl= Densmap[n][i][j][zperm[ndx1]];
410                        densr= Densmap[n][i][j][zperm[ndx1+1]];
411                        alpha=(densmid-densl)/(densr-densl);
412                        deltandx=zperm[ndx1+1]-zperm[ndx1];
413
414                        if(debug){
415                        printf("Alpha, Deltandx  %f %i\n", alpha,deltandx);
416                        }
417                        if(abs(alpha)>1.0 || abs(deltandx)>3){
418                        pos=zperm[ndx1];
419                        spread=-1;
420                        }
421                        else {
422                        pos=zperm[ndx1]+alpha*deltandx;
423                        spread=binwidth*deltandx;
424                        }
425                      * For the 2nd interface  can use the same formulation, since alpha should become negative ie:
426                      * alpha=(densmid-Densmap[n][i][j][zperm[ndx2]])/(Densmap[n][i][j][zperm[nxd2+1]]-Densmap[n][i][j][zperm[ndx2]]);
427                      * deltandx=zperm[ndx2+1]-zperm[ndx2];
428                      * pos=zperm[ndx2]+alpha*deltandx;   */
429
430                     /*After filtering we use the direct approach        */
431                     int1[n][j+(i*yslices)]->Z = (zperm[ndx1]+onehalf)*binwidth;
432                     int1[n][j+(i*yslices)]->t = binwidth;
433                     int2[n][j+(i*yslices)]->Z = (zperm[ndx2]+onehalf)*binwidth;
434                     int2[n][j+(i*yslices)]->t = binwidth;
435                 }
436             }
437         }
438     }
439
440     if (method == methFUNCFIT)
441     {
442         /*Assume a box divided in 2 along midpoint of z for starters*/
443         startpoint = 0.0;
444         endpoint   = binwidth*zslices;
445         splitpoint = (startpoint+endpoint)/2.0;
446         /*Initial fit proposals*/
447         beginfit1[0] = dens1;
448         beginfit1[1] = dens2;
449         beginfit1[2] = (splitpoint/2);
450         beginfit1[3] = 0.5;
451
452         beginfit2[0] = dens2;
453         beginfit2[1] = dens1;
454         beginfit2[2] = (3*splitpoint/2);
455         beginfit2[3] = 0.5;
456
457         snew(zDensavg, zslices);
458         snew(sigma1, zslices);
459         snew(sigma2, zslices);
460
461         for (k = 0; k < zslices; k++)
462         {
463             sigma1[k] = sigma2[k] = 1;
464         }
465         /*Calculate average density along z - avoid smoothing by using coarse-grained-mesh*/
466         for (k = 0; k < zslices; k++)
467         {
468             for (n = 0; n < tblocks; n++)
469             {
470                 for (i = 0; i < xslices; i++)
471                 {
472                     for (j = 0; j < yslices; j++)
473                     {
474                         zDensavg[k] += (Densmap[n][i][j][k]/(xslices*yslices*tblocks));
475                     }
476                 }
477             }
478         }
479
480         if (debug)
481         {
482             xvg = xvgropen("DensprofileonZ.xvg", "Averaged Densityprofile on Z", "z[nm]", "Density[kg/m^3]", oenv);
483             for (k = 0; k < zslices; k++)
484             {
485                 fprintf(xvg, "%4f.3   %8f.4\n", k*binwidth, zDensavg[k]);
486             }
487             xvgrclose(xvg);
488         }
489
490         /*Fit average density in z over whole trajectory to obtain tentative fit-parameters in fit1 and fit2*/
491
492         /*Fit 1st half of box*/
493         do_lmfit(zslices, zDensavg, sigma1, binwidth, NULL, startpoint, splitpoint, oenv, FALSE, effnERF, beginfit1, 3);
494         /*Fit 2nd half of box*/
495         do_lmfit(zslices, zDensavg, sigma2, binwidth, NULL, splitpoint, endpoint, oenv, FALSE, effnERF, beginfit2, 3);
496
497         /*Initialise the const arrays for storing the average fit parameters*/
498         avgfit1 = beginfit1;
499         avgfit2 = beginfit2;
500
501
502
503         /*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*/
504         for (n = 0; n < tblocks; n++)
505         {
506             for (i = 0; i < xslices; i++)
507             {
508                 for (j = 0; j < yslices; j++)
509                 {
510                     /*Reinitialise fit for each mesh-point*/
511                     srenew(fit1, 4);
512                     srenew(fit2, 4);
513                     for (k = 0; k < 4; k++)
514                     {
515                         fit1[k] = avgfit1[k];
516                         fit2[k] = avgfit2[k];
517                     }
518                     /*Now fit and store in structures in row-major order int[n][i][j]*/
519                     do_lmfit(zslices, Densmap[n][i][j], sigma1, binwidth, NULL, startpoint, splitpoint, oenv, FALSE, effnERF, fit1, 1);
520                     int1[n][j+(yslices*i)]->Z = fit1[2];
521                     int1[n][j+(yslices*i)]->t = fit1[3];
522                     do_lmfit(zslices, Densmap[n][i][j], sigma2, binwidth, NULL, splitpoint, endpoint, oenv, FALSE, effnERF, fit2, 2);
523                     int2[n][j+(yslices*i)]->Z = fit2[2];
524                     int2[n][j+(yslices*i)]->t = fit2[3];
525                 }
526             }
527         }
528     }
529
530
531     *intf1 = int1;
532     *intf2 = int2;
533
534 }
535
536 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 )
537 {
538     char   numbuf[13];
539     int    n, i, j;
540     real **profile1, **profile2;
541     real   max1, max2, min1, min2, *xticks, *yticks;
542     t_rgb  lo = {0, 0, 0};
543     t_rgb  hi = {1, 1, 1};
544     FILE  *xpmfile1, *xpmfile2;
545
546 /*Prepare xpm structures for output*/
547
548 /*Allocate memory to tick's and matrices*/
549     snew (xticks, xbins+1);
550     snew (yticks, ybins+1);
551
552     profile1 = mk_matrix(xbins, ybins, FALSE);
553     profile2 = mk_matrix(xbins, ybins, FALSE);
554
555     for (i = 0; i < xbins+1; i++)
556     {
557         xticks[i] += bw;
558     }
559     for (j = 0; j < ybins+1; j++)
560     {
561         yticks[j] += bw;
562     }
563
564     xpmfile1 = gmx_ffopen(outfiles[0], "w");
565     xpmfile2 = gmx_ffopen(outfiles[1], "w");
566
567     max1 = max2 = 0.0;
568     min1 = min2 = zbins*bwz;
569
570     for (n = 0; n < tblocks; n++)
571     {
572         sprintf(numbuf, "tblock: %4i", n);
573 /*Filling matrices for inclusion in xpm-files*/
574         for (i = 0; i < xbins; i++)
575         {
576             for (j = 0; j < ybins; j++)
577             {
578                 profile1[i][j] = (surf1[n][j+ybins*i])->Z;
579                 profile2[i][j] = (surf2[n][j+ybins*i])->Z;
580                 /*Finding max and min values*/
581                 if (profile1[i][j] > max1)
582                 {
583                     max1 = profile1[i][j];
584                 }
585                 if (profile1[i][j] < min1)
586                 {
587                     min1 = profile1[i][j];
588                 }
589                 if (profile2[i][j] > max2)
590                 {
591                     max2 = profile2[i][j];
592                 }
593                 if (profile2[i][j] < min2)
594                 {
595                     min2 = profile2[i][j];
596                 }
597             }
598         }
599
600         write_xpm(xpmfile1, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile1, min1, max1, lo, hi, &maplevels);
601         write_xpm(xpmfile2, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile2, min2, max2, lo, hi, &maplevels);
602     }
603
604     gmx_ffclose(xpmfile1);
605     gmx_ffclose(xpmfile2);
606
607
608     sfree(profile1);
609     sfree(profile2);
610     sfree(xticks);
611     sfree(yticks);
612 }
613
614 static void writeraw(t_interf ***int1, t_interf ***int2, int tblocks,
615                      int xbins, int ybins, char **fnms,
616                      const output_env_t oenv)
617 {
618     FILE *raw1, *raw2;
619     int   i, j, n;
620
621     raw1 = gmx_ffopen(fnms[0], "w");
622     raw2 = gmx_ffopen(fnms[1], "w");
623     try
624     {
625         gmx::BinaryInformationSettings settings;
626         settings.generatedByHeader(true);
627         settings.linePrefix("# ");
628         gmx::printBinaryInformation(raw1, output_env_get_program_context(oenv),
629                                     settings);
630         gmx::printBinaryInformation(raw2, output_env_get_program_context(oenv),
631                                     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         { efTPR, "-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(efTPR, 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, oenv);
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 }