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