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