284a99bc3158e4b9f1ca5625d57ff97fe1b20aa6
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_hydorder.c
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
36 #include "gmxpre.h"
37
38 #include "config.h"
39
40 #include <math.h>
41 #include <string.h>
42
43 #include "gromacs/legacyheaders/typedefs.h"
44 #include "gromacs/legacyheaders/macros.h"
45 #include "gstat.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/pbcutil/pbc.h"
48 #include "gromacs/utility/futil.h"
49 #include "gromacs/topology/index.h"
50 #include "gromacs/fileio/tpxio.h"
51 #include "gromacs/fileio/trxio.h"
52 #include "binsearch.h"
53 #include "powerspect.h"
54
55 #include "gromacs/commandline/pargs.h"
56 #include "gromacs/fileio/matio.h"
57 #include "gromacs/fileio/xvgr.h"
58 #include "gromacs/pbcutil/rmpbc.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/smalloc.h"
61
62 /* Print name of first atom in all groups in index file */
63 static void print_types(atom_id index[], atom_id a[], int ngrps,
64                         char *groups[], t_topology *top)
65 {
66     int i;
67
68     fprintf(stderr, "Using following groups: \n");
69     for (i = 0; i < ngrps; i++)
70     {
71         fprintf(stderr, "Groupname: %s First atomname: %s First atomnr %d\n",
72                 groups[i], *(top->atoms.atomname[a[index[i]]]), a[index[i]]);
73     }
74     fprintf(stderr, "\n");
75 }
76
77 static void check_length(real length, int a, int b)
78 {
79     if (length > 0.3)
80     {
81         fprintf(stderr, "WARNING: distance between atoms %d and "
82                 "%d > 0.3 nm (%f). Index file might be corrupt.\n",
83                 a, b, length);
84     }
85 }
86
87 static void find_tetra_order_grid(t_topology top, int ePBC,
88                                   int natoms, matrix box,
89                                   rvec x[], int maxidx, atom_id index[],
90                                   real *sgmean, real *skmean,
91                                   int nslicex, int nslicey, int nslicez,
92                                   real ***sggrid, real ***skgrid)
93 {
94     int         ix, jx, i, j, k, l, n, *nn[4];
95     rvec        dx, rj, rk, urk, urj;
96     real        cost, cost2, *sgmol, *skmol, rmean, rmean2, r2, box2, *r_nn[4];
97     t_pbc       pbc;
98     int         slindex_x, slindex_y, slindex_z;
99     int      ***sl_count;
100     real        onethird = 1.0/3.0;
101     gmx_rmpbc_t gpbc;
102
103     /*  dmat = init_mat(maxidx, FALSE); */
104
105     box2 = box[XX][XX] * box[XX][XX];
106
107     /* Initialize expanded sl_count array */
108     snew(sl_count, nslicex);
109     for (i = 0; i < nslicex; i++)
110     {
111         snew(sl_count[i], nslicey);
112         for (j = 0; j < nslicey; j++)
113         {
114             snew(sl_count[i][j], nslicez);
115         }
116     }
117
118
119     for (i = 0; (i < 4); i++)
120     {
121         snew(r_nn[i], natoms);
122         snew(nn[i], natoms);
123
124         for (j = 0; (j < natoms); j++)
125         {
126             r_nn[i][j] = box2;
127         }
128     }
129
130     snew(sgmol, maxidx);
131     snew(skmol, maxidx);
132
133     /* Must init pbc every step because of pressure coupling */
134     set_pbc(&pbc, ePBC, box);
135     gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
136     gmx_rmpbc(gpbc, natoms, box, x);
137
138     *sgmean = 0.0;
139     *skmean = 0.0;
140     l       = 0;
141     for (i = 0; (i < maxidx); i++)
142     {   /* loop over index file */
143         ix = index[i];
144         for (j = 0; (j < maxidx); j++)
145         {
146
147             if (i == j)
148             {
149                 continue;
150             }
151
152             jx = index[j];
153
154             pbc_dx(&pbc, x[ix], x[jx], dx);
155             r2 = iprod(dx, dx);
156
157             /* set_mat_entry(dmat,i,j,r2); */
158
159             /* determine the nearest neighbours */
160             if (r2 < r_nn[0][i])
161             {
162                 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
163                 r_nn[2][i] = r_nn[1][i]; nn[2][i] = nn[1][i];
164                 r_nn[1][i] = r_nn[0][i]; nn[1][i] = nn[0][i];
165                 r_nn[0][i] = r2;         nn[0][i] = j;
166             }
167             else if (r2 < r_nn[1][i])
168             {
169                 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
170                 r_nn[2][i] = r_nn[1][i]; nn[2][i] = nn[1][i];
171                 r_nn[1][i] = r2;         nn[1][i] = j;
172             }
173             else if (r2 < r_nn[2][i])
174             {
175                 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
176                 r_nn[2][i] = r2;         nn[2][i] = j;
177             }
178             else if (r2 < r_nn[3][i])
179             {
180                 r_nn[3][i] = r2;         nn[3][i] = j;
181             }
182         }
183
184
185         /* calculate mean distance between nearest neighbours */
186         rmean = 0;
187         for (j = 0; (j < 4); j++)
188         {
189             r_nn[j][i] = sqrt(r_nn[j][i]);
190             rmean     += r_nn[j][i];
191         }
192         rmean /= 4;
193
194         n        = 0;
195         sgmol[i] = 0.0;
196         skmol[i] = 0.0;
197
198         /* Chau1998a eqn 3 */
199         /* angular part tetrahedrality order parameter per atom */
200         for (j = 0; (j < 3); j++)
201         {
202             for (k = j+1; (k < 4); k++)
203             {
204                 pbc_dx(&pbc, x[ix], x[index[nn[k][i]]], rk);
205                 pbc_dx(&pbc, x[ix], x[index[nn[j][i]]], rj);
206
207                 unitv(rk, urk);
208                 unitv(rj, urj);
209
210                 cost  = iprod(urk, urj) + onethird;
211                 cost2 = cost * cost;
212
213                 sgmol[i] += cost2;
214                 l++;
215                 n++;
216             }
217         }
218         /* normalize sgmol between 0.0 and 1.0 */
219         sgmol[i] = 3*sgmol[i]/32;
220         *sgmean += sgmol[i];
221
222         /* distance part tetrahedrality order parameter per atom */
223         rmean2 = 4 * 3 * rmean * rmean;
224         for (j = 0; (j < 4); j++)
225         {
226             skmol[i] += (rmean - r_nn[j][i]) * (rmean - r_nn[j][i]) / rmean2;
227             /*      printf("%d %f (%f %f %f %f) \n",
228                     i, skmol[i], rmean, rmean2, r_nn[j][i], (rmean - r_nn[j][i]) );
229              */
230         }
231
232         *skmean += skmol[i];
233
234         /* Compute sliced stuff in x y z*/
235         slindex_x = gmx_nint((1+x[i][XX]/box[XX][XX])*nslicex) % nslicex;
236         slindex_y = gmx_nint((1+x[i][YY]/box[YY][YY])*nslicey) % nslicey;
237         slindex_z = gmx_nint((1+x[i][ZZ]/box[ZZ][ZZ])*nslicez) % nslicez;
238         sggrid[slindex_x][slindex_y][slindex_z] += sgmol[i];
239         skgrid[slindex_x][slindex_y][slindex_z] += skmol[i];
240         (sl_count[slindex_x][slindex_y][slindex_z])++;
241     } /* loop over entries in index file */
242
243     *sgmean /= maxidx;
244     *skmean /= maxidx;
245
246     for (i = 0; (i < nslicex); i++)
247     {
248         for (j = 0; j < nslicey; j++)
249         {
250             for (k = 0; k < nslicez; k++)
251             {
252                 if (sl_count[i][j][k] > 0)
253                 {
254                     sggrid[i][j][k] /= sl_count[i][j][k];
255                     skgrid[i][j][k] /= sl_count[i][j][k];
256                 }
257             }
258         }
259     }
260
261     sfree(sl_count);
262     sfree(sgmol);
263     sfree(skmol);
264     for (i = 0; (i < 4); i++)
265     {
266         sfree(r_nn[i]);
267         sfree(nn[i]);
268     }
269 }
270
271 /*Determines interface from tetrahedral order parameter in box with specified binwidth.  */
272 /*Outputs interface positions(bins), the number of timeframes, and the number of surface-mesh points in xy*/
273
274 static void calc_tetra_order_interface(const char *fnNDX, const char *fnTPS, const char *fnTRX, real binw, int tblock,
275                                        int *nframes,  int *nslicex, int *nslicey,
276                                        real sgang1, real sgang2, real ****intfpos,
277                                        output_env_t oenv)
278 {
279     FILE         *fpsg   = NULL, *fpsk = NULL;
280     char         *sgslfn = "sg_ang_mesh";  /* Hardcoded filenames for debugging*/
281     char         *skslfn = "sk_dist_mesh";
282     t_topology    top;
283     int           ePBC;
284     char          title[STRLEN], subtitle[STRLEN];
285     t_trxstatus  *status;
286     int           natoms;
287     real          t;
288     rvec         *xtop, *x;
289     matrix        box;
290     real          sg, sk, sgintf, pos;
291     atom_id     **index   = NULL;
292     char        **grpname = NULL;
293     int           i, j, k, n, *isize, ng, nslicez, framenr;
294     real       ***sg_grid = NULL, ***sk_grid = NULL, ***sg_fravg = NULL, ***sk_fravg = NULL, ****sk_4d = NULL, ****sg_4d = NULL;
295     int          *perm;
296     int           ndx1, ndx2;
297     int           bins;
298     const real    onehalf = 1.0/2.0;
299     /* real   ***intfpos[2]; pointers to arrays of two interface positions zcoord(framenr,xbin,ybin): intfpos[interface_index][t][nslicey*x+y]
300      * i.e 1D Row-major order in (t,x,y) */
301
302
303     read_tps_conf(fnTPS, title, &top, &ePBC, &xtop, NULL, box, FALSE);
304
305     *nslicex = (int)(box[XX][XX]/binw + onehalf); /*Calculate slicenr from binwidth*/
306     *nslicey = (int)(box[YY][YY]/binw + onehalf);
307     nslicez  = (int)(box[ZZ][ZZ]/binw +  onehalf);
308
309
310
311     ng = 1;
312     /* get index groups */
313     printf("Select the group that contains the atoms you want to use for the tetrahedrality order parameter calculation:\n");
314     snew(grpname, ng);
315     snew(index, ng);
316     snew(isize, ng);
317     get_index(&top.atoms, fnNDX, ng, isize, index, grpname);
318
319     /* Analyze trajectory */
320     natoms = read_first_x(oenv, &status, fnTRX, &t, &x, box);
321     if (natoms > top.atoms.nr)
322     {
323         gmx_fatal(FARGS, "Topology (%d atoms) does not match trajectory (%d atoms)",
324                   top.atoms.nr, natoms);
325     }
326     check_index(NULL, ng, index[0], NULL, natoms);
327
328
329     /*Prepare structures for temporary storage of frame info*/
330     snew(sg_grid, *nslicex);
331     snew(sk_grid, *nslicex);
332     for (i = 0; i < *nslicex; i++)
333     {
334         snew(sg_grid[i], *nslicey);
335         snew(sk_grid[i], *nslicey);
336         for (j = 0; j < *nslicey; j++)
337         {
338             snew(sg_grid[i][j], nslicez);
339             snew(sk_grid[i][j], nslicez);
340         }
341     }
342
343     sg_4d    = NULL;
344     sk_4d    = NULL;
345     *nframes = 0;
346     framenr  = 0;
347
348 /* Loop over frames*/
349     do
350     {
351         /*Initialize box meshes (temporary storage for each tblock frame -reinitialise every tblock steps */
352         if (framenr%tblock == 0)
353         {
354             srenew(sk_4d, *nframes+1);
355             srenew(sg_4d, *nframes+1);
356             snew(sg_fravg, *nslicex);
357             snew(sk_fravg, *nslicex);
358             for (i = 0; i < *nslicex; i++)
359             {
360                 snew(sg_fravg[i], *nslicey);
361                 snew(sk_fravg[i], *nslicey);
362                 for (j = 0; j < *nslicey; j++)
363                 {
364                     snew(sg_fravg[i][j], nslicez);
365                     snew(sk_fravg[i][j], nslicez);
366                 }
367             }
368         }
369
370         find_tetra_order_grid(top, ePBC, natoms, box, x, isize[0], index[0],
371                               &sg, &sk, *nslicex, *nslicey, nslicez, sg_grid, sk_grid);
372         for (i = 0; i < *nslicex; i++)
373         {
374             for (j = 0; j < *nslicey; j++)
375             {
376                 for (k = 0; k < nslicez; k++)
377                 {
378                     sk_fravg[i][j][k] += sk_grid[i][j][k]/tblock;
379                     sg_fravg[i][j][k] += sg_grid[i][j][k]/tblock;
380                 }
381             }
382         }
383
384         framenr++;
385
386         if (framenr%tblock == 0)
387         {
388             sk_4d[*nframes] = sk_fravg;
389             sg_4d[*nframes] = sg_fravg;
390             (*nframes)++;
391         }
392
393     }
394     while (read_next_x(oenv, status, &t, x, box));
395     close_trj(status);
396
397     sfree(grpname);
398     sfree(index);
399     sfree(isize);
400
401     /*Debugging for printing out the entire order parameter meshes.*/
402     if (debug)
403     {
404         fpsg = xvgropen(sgslfn, "S\\sg\\N Angle Order Parameter / Meshpoint", "(nm)", "S\\sg\\N", oenv);
405         fpsk = xvgropen(skslfn, "S\\sk\\N Distance Order Parameter / Meshpoint", "(nm)", "S\\sk\\N", oenv);
406         for (n = 0; n < (*nframes); n++)
407         {
408             fprintf(fpsg, "%i\n", n);
409             fprintf(fpsk, "%i\n", n);
410             for (i = 0; (i < *nslicex); i++)
411             {
412                 for (j = 0; j < *nslicey; j++)
413                 {
414                     for (k = 0; k < nslicez; k++)
415                     {
416                         fprintf(fpsg, "%4f  %4f  %4f  %8f\n", (i+0.5)*box[XX][XX]/(*nslicex), (j+0.5)*box[YY][YY]/(*nslicey), (k+0.5)*box[ZZ][ZZ]/nslicez, sg_4d[n][i][j][k]);
417                         fprintf(fpsk, "%4f  %4f  %4f  %8f\n", (i+0.5)*box[XX][XX]/(*nslicex), (j+0.5)*box[YY][YY]/(*nslicey), (k+0.5)*box[ZZ][ZZ]/nslicez, sk_4d[n][i][j][k]);
418                     }
419                 }
420             }
421         }
422         xvgrclose(fpsg);
423         xvgrclose(fpsk);
424     }
425
426
427     /* Find positions of interface z by scanning orderparam for each frame and for each xy-mesh cylinder along z*/
428
429     /*Simple trial: assume interface is in the middle of -sgang1 and sgang2*/
430     sgintf = 0.5*(sgang1+sgang2);
431
432
433     /*Allocate memory for interface arrays; */
434     snew((*intfpos), 2);
435     snew((*intfpos)[0], *nframes);
436     snew((*intfpos)[1], *nframes);
437
438     bins = (*nslicex)*(*nslicey);
439
440
441     snew(perm, nslicez);  /*permutation array for sorting along normal coordinate*/
442
443
444     for (n = 0; n < *nframes; n++)
445     {
446         snew((*intfpos)[0][n], bins);
447         snew((*intfpos)[1][n], bins);
448         for (i = 0; i < *nslicex; i++)
449         {
450             for (j = 0; j < *nslicey; j++)
451             {
452                 rangeArray(perm, nslicez); /*reset permutation array to identity*/
453                 /*Binsearch returns 2 bin-numbers where the order param is  <= setpoint sgintf*/
454                 ndx1 = start_binsearch(sg_4d[n][i][j], perm, 0, nslicez/2-1, sgintf, 1);
455                 ndx2 = start_binsearch(sg_4d[n][i][j], perm, nslicez/2, nslicez-1, sgintf, -1);
456                 /*Use linear interpolation to smooth out the interface position*/
457
458                 /*left interface (0)*/
459                 /*if((sg_4d[n][i][j][perm[ndx1+1]]-sg_4d[n][i][j][perm[ndx1]])/sg_4d[n][i][j][perm[ndx1]] > 0.01){
460                    pos=( (sgintf-sg_4d[n][i][j][perm[ndx1]])*perm[ndx1+1]+(sg_4d[n][i][j][perm[ndx1+1]]-sgintf)*perm[ndx1 ])*/
461                 (*intfpos)[0][n][j+*nslicey*i] = (perm[ndx1]+onehalf)*binw;
462                 /*right interface (1)*/
463                 /*alpha=(sgintf-sg_4d[n][i][j][perm[ndx2]])/(sg_4d[n][i][j][perm[ndx2]+1]-sg_4d[n][i][j][perm[ndx2]]);*/
464                 /*(*intfpos)[1][n][j+*nslicey*i]=((1-alpha)*perm[ndx2]+alpha*(perm[ndx2]+1)+onehalf)*box[ZZ][ZZ]/nslicez;*/
465                 (*intfpos)[1][n][j+*nslicey*i] = (perm[ndx2]+onehalf)*binw;
466             }
467         }
468     }
469
470
471     /*sfree(perm);*/
472     sfree(sk_4d);
473     sfree(sg_4d);
474     /*sfree(sg_grid);*/
475     /*sfree(sk_grid);*/
476
477
478 }
479
480 static void writesurftoxpms(real ***surf, int tblocks, int xbins, int ybins, real bw, char **outfiles, int maplevels )
481 {
482
483     char   numbuf[8];
484     int    n, i, j;
485     real **profile1, **profile2;
486     real   max1, max2, min1, min2, *xticks, *yticks;
487     t_rgb  lo = {1, 1, 1};
488     t_rgb  hi = {0, 0, 0};
489     FILE  *xpmfile1, *xpmfile2;
490
491 /*Prepare xpm structures for output*/
492
493 /*Allocate memory to tick's and matrices*/
494     snew (xticks, xbins+1);
495     snew (yticks, ybins+1);
496
497     profile1 = mk_matrix(xbins, ybins, FALSE);
498     profile2 = mk_matrix(xbins, ybins, FALSE);
499
500     for (i = 0; i < xbins+1; i++)
501     {
502         xticks[i] += bw;
503     }
504     for (j = 0; j < ybins+1; j++)
505     {
506         yticks[j] += bw;
507     }
508
509     xpmfile1 = gmx_ffopen(outfiles[0], "w");
510     xpmfile2 = gmx_ffopen(outfiles[1], "w");
511
512     max1 = max2 = 0.0;
513     min1 = min2 = 1000.00;
514
515     for (n = 0; n < tblocks; n++)
516     {
517         sprintf(numbuf, "%5d", n);
518         /*Filling matrices for inclusion in xpm-files*/
519         for (i = 0; i < xbins; i++)
520         {
521             for (j = 0; j < ybins; j++)
522             {
523                 profile1[i][j] = (surf[0][n][j+ybins*i]);
524                 profile2[i][j] = (surf[1][n][j+ybins*i]);
525                 /*Finding max and min values*/
526                 if (profile1[i][j] > max1)
527                 {
528                     max1 = profile1[i][j];
529                 }
530                 if (profile1[i][j] < min1)
531                 {
532                     min1 = profile1[i][j];
533                 }
534                 if (profile2[i][j] > max2)
535                 {
536                     max2 = profile2[i][j];
537                 }
538                 if (profile2[i][j] < min2)
539                 {
540                     min2 = profile2[i][j];
541                 }
542             }
543         }
544
545         write_xpm(xpmfile1, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile1, min1, max1, lo, hi, &maplevels);
546         write_xpm(xpmfile2, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile2, min2, max2, lo, hi, &maplevels);
547     }
548
549     gmx_ffclose(xpmfile1);
550     gmx_ffclose(xpmfile2);
551
552
553
554     sfree(profile1);
555     sfree(profile2);
556     sfree(xticks);
557     sfree(yticks);
558 }
559
560
561 static void writeraw(real ***surf, int tblocks, int xbins, int ybins, char **fnms)
562 {
563     FILE *raw1, *raw2;
564     int   i, j, n;
565
566     raw1 = gmx_ffopen(fnms[0], "w");
567     raw2 = gmx_ffopen(fnms[1], "w");
568     fprintf(raw1, "#Legend\n#TBlock\n#Xbin Ybin Z t\n");
569     fprintf(raw2, "#Legend\n#TBlock\n#Xbin Ybin Z t\n");
570     for (n = 0; n < tblocks; n++)
571     {
572         fprintf(raw1, "%5d\n", n);
573         fprintf(raw2, "%5d\n", n);
574         for (i = 0; i < xbins; i++)
575         {
576             for (j = 0; j < ybins; j++)
577             {
578                 fprintf(raw1, "%i  %i  %8.5f\n", i, j, (surf[0][n][j+ybins*i]));
579                 fprintf(raw2, "%i  %i  %8.5f\n", i, j, (surf[1][n][j+ybins*i]));
580             }
581         }
582     }
583
584     gmx_ffclose(raw1);
585     gmx_ffclose(raw2);
586 }
587
588
589
590 int gmx_hydorder(int argc, char *argv[])
591 {
592     static const char *desc[] = {
593         "[THISMODULE] computes the tetrahedrality order parameters around a ",
594         "given atom. Both angle an distance order parameters are calculated. See",
595         "P.-L. Chau and A.J. Hardwick, Mol. Phys., 93, (1998), 511-518.",
596         "for more details.[PAR]"
597         "[THISMODULE] calculates the order parameter in a 3d-mesh in the box, and",
598         "with 2 phases in the box gives the user the option to define a 2D interface in time",
599         "separating the faces by specifying parameters [TT]-sgang1[tt] and",
600         "[TT]-sgang2[tt] (it is important to select these judiciously)."
601     };
602
603     int                axis      = 0;
604     static int         nsttblock = 1;
605     static int         nlevels   = 100;
606     static real        binwidth  = 1.0; /* binwidth in mesh           */
607     static real        sg1       = 1;
608     static real        sg2       = 1;   /* order parameters for bulk phases */
609     static gmx_bool    bFourier  = FALSE;
610     static gmx_bool    bRawOut   = FALSE;
611     int                frames, xslices, yslices; /* Dimensions of interface arrays*/
612     real            ***intfpos;                  /* Interface arrays (intfnr,t,xy) -potentially large */
613     static char       *normal_axis[] = { NULL, "z", "x", "y", NULL };
614
615     t_pargs            pa[] = {
616         { "-d",   FALSE, etENUM, {normal_axis},
617           "Direction of the normal on the membrane" },
618         { "-bw",  FALSE, etREAL, {&binwidth},
619           "Binwidth of box mesh" },
620         { "-sgang1", FALSE, etREAL, {&sg1},
621           "tetrahedral angle parameter in Phase 1 (bulk)" },
622         { "-sgang2", FALSE, etREAL, {&sg2},
623           "tetrahedral angle parameter in Phase 2 (bulk)" },
624         { "-tblock", FALSE, etINT, {&nsttblock},
625           "Number of frames in one time-block average"},
626         { "-nlevel", FALSE, etINT, {&nlevels},
627           "Number of Height levels in 2D - XPixMaps"}
628     };
629
630     t_filenm           fnm[] = {                      /* files for g_order    */
631         { efTRX, "-f", NULL,  ffREAD },               /* trajectory file              */
632         { efNDX, "-n", NULL,  ffREAD },               /* index file           */
633         { efTPX, "-s", NULL,  ffREAD },               /* topology file                */
634         { efXPM, "-o", "intf",  ffWRMULT},            /* XPM- surface maps      */
635         { efOUT, "-or", "raw", ffOPTWRMULT },         /* xvgr output file           */
636         { efOUT, "-Spect", "intfspect", ffOPTWRMULT}, /* Fourier spectrum interfaces */
637     };
638 #define NFILE asize(fnm)
639
640     /*Filenames*/
641     const char  *ndxfnm, *tpsfnm, *trxfnm;
642     char       **spectra, **intfn, **raw;
643     int          nfspect, nfxpm, nfraw;
644     output_env_t oenv;
645
646     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME,
647                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
648     {
649         return 0;
650     }
651     bFourier = opt2bSet("-Spect", NFILE, fnm);
652     bRawOut  = opt2bSet("-or", NFILE, fnm);
653
654     if (binwidth < 0.0)
655     {
656         gmx_fatal(FARGS, "Can not have binwidth < 0");
657     }
658
659     ndxfnm = ftp2fn(efNDX, NFILE, fnm);
660     tpsfnm = ftp2fn(efTPX, NFILE, fnm);
661     trxfnm = ftp2fn(efTRX, NFILE, fnm);
662
663     /* Calculate axis */
664     if (strcmp(normal_axis[0], "x") == 0)
665     {
666         axis = XX;
667     }
668     else if (strcmp(normal_axis[0], "y") == 0)
669     {
670         axis = YY;
671     }
672     else if (strcmp(normal_axis[0], "z") == 0)
673     {
674         axis = ZZ;
675     }
676     else
677     {
678         gmx_fatal(FARGS, "Invalid axis, use x, y or z");
679     }
680
681     switch (axis)
682     {
683         case 0:
684             fprintf(stderr, "Taking x axis as normal to the membrane\n");
685             break;
686         case 1:
687             fprintf(stderr, "Taking y axis as normal to the membrane\n");
688             break;
689         case 2:
690             fprintf(stderr, "Taking z axis as normal to the membrane\n");
691             break;
692     }
693
694     /* tetraheder order parameter */
695     /* If either of the options is set we compute both */
696     nfxpm = opt2fns(&intfn, "-o", NFILE, fnm);
697     if (nfxpm != 2)
698     {
699         gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfxpm);
700     }
701     calc_tetra_order_interface(ndxfnm, tpsfnm, trxfnm, binwidth, nsttblock, &frames, &xslices, &yslices, sg1, sg2, &intfpos, oenv);
702     writesurftoxpms(intfpos, frames, xslices, yslices, binwidth, intfn, nlevels);
703
704     if (bFourier)
705     {
706         nfspect = opt2fns(&spectra, "-Spect", NFILE, fnm);
707         if (nfspect != 2)
708         {
709             gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfspect);
710         }
711         powerspectavg(intfpos, frames, xslices, yslices, spectra);
712     }
713
714     if (bRawOut)
715     {
716         nfraw = opt2fns(&raw, "-or", NFILE, fnm);
717         if (nfraw != 2)
718         {
719             gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfraw);
720         }
721         writeraw(intfpos, frames, xslices, yslices, raw);
722     }
723
724     return 0;
725 }