Converted 2xnn kernel to C++
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_hydorder.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2011,2012,2013,2014,2015, 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 <cmath>
39 #include <cstring>
40
41 #include "gromacs/commandline/pargs.h"
42 #include "gromacs/fileio/confio.h"
43 #include "gromacs/fileio/matio.h"
44 #include "gromacs/fileio/trxio.h"
45 #include "gromacs/fileio/xvgr.h"
46 #include "gromacs/gmxana/binsearch.h"
47 #include "gromacs/gmxana/gmx_ana.h"
48 #include "gromacs/gmxana/gstat.h"
49 #include "gromacs/gmxana/powerspect.h"
50 #include "gromacs/legacyheaders/macros.h"
51 #include "gromacs/legacyheaders/typedefs.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/pbcutil/pbc.h"
54 #include "gromacs/pbcutil/rmpbc.h"
55 #include "gromacs/topology/index.h"
56 #include "gromacs/utility/cstringutil.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/futil.h"
59 #include "gromacs/utility/gmxassert.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] = std::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     t_topology    top;
281     int           ePBC;
282     char          title[STRLEN];
283     t_trxstatus  *status;
284     int           natoms;
285     real          t;
286     rvec         *xtop, *x;
287     matrix        box;
288     real          sg, sk, sgintf;
289     atom_id     **index   = NULL;
290     char        **grpname = NULL;
291     int           i, j, k, n, *isize, ng, nslicez, framenr;
292     real       ***sg_grid = NULL, ***sk_grid = NULL, ***sg_fravg = NULL, ***sk_fravg = NULL, ****sk_4d = NULL, ****sg_4d = NULL;
293     int          *perm;
294     int           ndx1, ndx2;
295     int           bins;
296     const real    onehalf = 1.0/2.0;
297     /* real   ***intfpos[2]; pointers to arrays of two interface positions zcoord(framenr,xbin,ybin): intfpos[interface_index][t][nslicey*x+y]
298      * i.e 1D Row-major order in (t,x,y) */
299
300
301     read_tps_conf(fnTPS, title, &top, &ePBC, &xtop, NULL, box, FALSE);
302
303     *nslicex = static_cast<int>(box[XX][XX]/binw + onehalf); /*Calculate slicenr from binwidth*/
304     *nslicey = static_cast<int>(box[YY][YY]/binw + onehalf);
305     nslicez  = static_cast<int>(box[ZZ][ZZ]/binw +  onehalf);
306
307
308
309     ng = 1;
310     /* get index groups */
311     printf("Select the group that contains the atoms you want to use for the tetrahedrality order parameter calculation:\n");
312     snew(grpname, ng);
313     snew(index, ng);
314     snew(isize, ng);
315     get_index(&top.atoms, fnNDX, ng, isize, index, grpname);
316
317     /* Analyze trajectory */
318     natoms = read_first_x(oenv, &status, fnTRX, &t, &x, box);
319     if (natoms > top.atoms.nr)
320     {
321         gmx_fatal(FARGS, "Topology (%d atoms) does not match trajectory (%d atoms)",
322                   top.atoms.nr, natoms);
323     }
324     check_index(NULL, ng, index[0], NULL, natoms);
325
326
327     /*Prepare structures for temporary storage of frame info*/
328     snew(sg_grid, *nslicex);
329     snew(sk_grid, *nslicex);
330     for (i = 0; i < *nslicex; i++)
331     {
332         snew(sg_grid[i], *nslicey);
333         snew(sk_grid[i], *nslicey);
334         for (j = 0; j < *nslicey; j++)
335         {
336             snew(sg_grid[i][j], nslicez);
337             snew(sk_grid[i][j], nslicez);
338         }
339     }
340
341     sg_4d    = NULL;
342     sk_4d    = NULL;
343     *nframes = 0;
344     framenr  = 0;
345
346 /* Loop over frames*/
347     do
348     {
349         /*Initialize box meshes (temporary storage for each tblock frame -reinitialise every tblock steps */
350         if (framenr%tblock == 0)
351         {
352             srenew(sk_4d, *nframes+1);
353             srenew(sg_4d, *nframes+1);
354             snew(sg_fravg, *nslicex);
355             snew(sk_fravg, *nslicex);
356             for (i = 0; i < *nslicex; i++)
357             {
358                 snew(sg_fravg[i], *nslicey);
359                 snew(sk_fravg[i], *nslicey);
360                 for (j = 0; j < *nslicey; j++)
361                 {
362                     snew(sg_fravg[i][j], nslicez);
363                     snew(sk_fravg[i][j], nslicez);
364                 }
365             }
366         }
367
368         find_tetra_order_grid(top, ePBC, natoms, box, x, isize[0], index[0],
369                               &sg, &sk, *nslicex, *nslicey, nslicez, sg_grid, sk_grid);
370         GMX_RELEASE_ASSERT(sk_fravg != NULL, "Trying to dereference NULL sk_fravg pointer");
371         for (i = 0; i < *nslicex; i++)
372         {
373             for (j = 0; j < *nslicey; j++)
374             {
375                 for (k = 0; k < nslicez; k++)
376                 {
377                     sk_fravg[i][j][k] += sk_grid[i][j][k]/tblock;
378                     sg_fravg[i][j][k] += sg_grid[i][j][k]/tblock;
379                 }
380             }
381         }
382
383         framenr++;
384
385         if (framenr%tblock == 0)
386         {
387             GMX_RELEASE_ASSERT(sk_4d != NULL, "Trying to dereference NULL sk_4d pointer");
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("sg_ang_mesh", "S\\sg\\N Angle Order Parameter / Meshpoint", "(nm)", "S\\sg\\N", oenv);
405         fpsk = xvgropen("sk_dist_mesh", "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(sk_4d);
472     sfree(sg_4d);
473 }
474
475 static void writesurftoxpms(real ***surf, int tblocks, int xbins, int ybins, real bw, char **outfiles, int maplevels )
476 {
477
478     char   numbuf[8];
479     int    n, i, j;
480     real **profile1, **profile2;
481     real   max1, max2, min1, min2, *xticks, *yticks;
482     t_rgb  lo = {1, 1, 1};
483     t_rgb  hi = {0, 0, 0};
484     FILE  *xpmfile1, *xpmfile2;
485
486 /*Prepare xpm structures for output*/
487
488 /*Allocate memory to tick's and matrices*/
489     snew (xticks, xbins+1);
490     snew (yticks, ybins+1);
491
492     profile1 = mk_matrix(xbins, ybins, FALSE);
493     profile2 = mk_matrix(xbins, ybins, FALSE);
494
495     for (i = 0; i < xbins+1; i++)
496     {
497         xticks[i] += bw;
498     }
499     for (j = 0; j < ybins+1; j++)
500     {
501         yticks[j] += bw;
502     }
503
504     xpmfile1 = gmx_ffopen(outfiles[0], "w");
505     xpmfile2 = gmx_ffopen(outfiles[1], "w");
506
507     max1 = max2 = 0.0;
508     min1 = min2 = 1000.00;
509
510     for (n = 0; n < tblocks; n++)
511     {
512         sprintf(numbuf, "%5d", n);
513         /*Filling matrices for inclusion in xpm-files*/
514         for (i = 0; i < xbins; i++)
515         {
516             for (j = 0; j < ybins; j++)
517             {
518                 profile1[i][j] = (surf[0][n][j+ybins*i]);
519                 profile2[i][j] = (surf[1][n][j+ybins*i]);
520                 /*Finding max and min values*/
521                 if (profile1[i][j] > max1)
522                 {
523                     max1 = profile1[i][j];
524                 }
525                 if (profile1[i][j] < min1)
526                 {
527                     min1 = profile1[i][j];
528                 }
529                 if (profile2[i][j] > max2)
530                 {
531                     max2 = profile2[i][j];
532                 }
533                 if (profile2[i][j] < min2)
534                 {
535                     min2 = profile2[i][j];
536                 }
537             }
538         }
539
540         write_xpm(xpmfile1, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile1, min1, max1, lo, hi, &maplevels);
541         write_xpm(xpmfile2, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile2, min2, max2, lo, hi, &maplevels);
542     }
543
544     gmx_ffclose(xpmfile1);
545     gmx_ffclose(xpmfile2);
546
547
548
549     sfree(profile1);
550     sfree(profile2);
551     sfree(xticks);
552     sfree(yticks);
553 }
554
555
556 static void writeraw(real ***surf, int tblocks, int xbins, int ybins, char **fnms)
557 {
558     FILE *raw1, *raw2;
559     int   i, j, n;
560
561     raw1 = gmx_ffopen(fnms[0], "w");
562     raw2 = gmx_ffopen(fnms[1], "w");
563     fprintf(raw1, "#Legend\n#TBlock\n#Xbin Ybin Z t\n");
564     fprintf(raw2, "#Legend\n#TBlock\n#Xbin Ybin Z t\n");
565     for (n = 0; n < tblocks; n++)
566     {
567         fprintf(raw1, "%5d\n", n);
568         fprintf(raw2, "%5d\n", n);
569         for (i = 0; i < xbins; i++)
570         {
571             for (j = 0; j < ybins; j++)
572             {
573                 fprintf(raw1, "%i  %i  %8.5f\n", i, j, (surf[0][n][j+ybins*i]));
574                 fprintf(raw2, "%i  %i  %8.5f\n", i, j, (surf[1][n][j+ybins*i]));
575             }
576         }
577     }
578
579     gmx_ffclose(raw1);
580     gmx_ffclose(raw2);
581 }
582
583
584
585 int gmx_hydorder(int argc, char *argv[])
586 {
587     static const char *desc[] = {
588         "[THISMODULE] computes the tetrahedrality order parameters around a ",
589         "given atom. Both angle an distance order parameters are calculated. See",
590         "P.-L. Chau and A.J. Hardwick, Mol. Phys., 93, (1998), 511-518.",
591         "for more details.[PAR]"
592         "[THISMODULE] calculates the order parameter in a 3d-mesh in the box, and",
593         "with 2 phases in the box gives the user the option to define a 2D interface in time",
594         "separating the faces by specifying parameters [TT]-sgang1[tt] and",
595         "[TT]-sgang2[tt] (it is important to select these judiciously)."
596     };
597
598     int                axis      = 0;
599     static int         nsttblock = 1;
600     static int         nlevels   = 100;
601     static real        binwidth  = 1.0; /* binwidth in mesh           */
602     static real        sg1       = 1;
603     static real        sg2       = 1;   /* order parameters for bulk phases */
604     static gmx_bool    bFourier  = FALSE;
605     static gmx_bool    bRawOut   = FALSE;
606     int                frames, xslices, yslices; /* Dimensions of interface arrays*/
607     real            ***intfpos;                  /* Interface arrays (intfnr,t,xy) -potentially large */
608     static const char *normal_axis[] = { NULL, "z", "x", "y", NULL };
609
610     t_pargs            pa[] = {
611         { "-d",   FALSE, etENUM, {normal_axis},
612           "Direction of the normal on the membrane" },
613         { "-bw",  FALSE, etREAL, {&binwidth},
614           "Binwidth of box mesh" },
615         { "-sgang1", FALSE, etREAL, {&sg1},
616           "tetrahedral angle parameter in Phase 1 (bulk)" },
617         { "-sgang2", FALSE, etREAL, {&sg2},
618           "tetrahedral angle parameter in Phase 2 (bulk)" },
619         { "-tblock", FALSE, etINT, {&nsttblock},
620           "Number of frames in one time-block average"},
621         { "-nlevel", FALSE, etINT, {&nlevels},
622           "Number of Height levels in 2D - XPixMaps"}
623     };
624
625     t_filenm           fnm[] = {                      /* files for g_order    */
626         { efTRX, "-f", NULL,  ffREAD },               /* trajectory file              */
627         { efNDX, "-n", NULL,  ffREAD },               /* index file           */
628         { efTPR, "-s", NULL,  ffREAD },               /* topology file                */
629         { efXPM, "-o", "intf",  ffWRMULT},            /* XPM- surface maps      */
630         { efOUT, "-or", "raw", ffOPTWRMULT },         /* xvgr output file           */
631         { efOUT, "-Spect", "intfspect", ffOPTWRMULT}, /* Fourier spectrum interfaces */
632     };
633 #define NFILE asize(fnm)
634
635     /*Filenames*/
636     const char  *ndxfnm, *tpsfnm, *trxfnm;
637     char       **spectra, **intfn, **raw;
638     int          nfspect, nfxpm, nfraw;
639     output_env_t oenv;
640
641     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME,
642                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
643     {
644         return 0;
645     }
646     bFourier = opt2bSet("-Spect", NFILE, fnm);
647     bRawOut  = opt2bSet("-or", NFILE, fnm);
648
649     if (binwidth < 0.0)
650     {
651         gmx_fatal(FARGS, "Can not have binwidth < 0");
652     }
653
654     ndxfnm = ftp2fn(efNDX, NFILE, fnm);
655     tpsfnm = ftp2fn(efTPR, NFILE, fnm);
656     trxfnm = ftp2fn(efTRX, NFILE, fnm);
657
658     /* Calculate axis */
659     GMX_RELEASE_ASSERT(normal_axis[0] != NULL, "Option setting inconsistency; normal_axis[0] is NULL");
660     if (std::strcmp(normal_axis[0], "x") == 0)
661     {
662         axis = XX;
663     }
664     else if (std::strcmp(normal_axis[0], "y") == 0)
665     {
666         axis = YY;
667     }
668     else if (std::strcmp(normal_axis[0], "z") == 0)
669     {
670         axis = ZZ;
671     }
672     else
673     {
674         gmx_fatal(FARGS, "Invalid axis, use x, y or z");
675     }
676
677     switch (axis)
678     {
679         case 0:
680             fprintf(stderr, "Taking x axis as normal to the membrane\n");
681             break;
682         case 1:
683             fprintf(stderr, "Taking y axis as normal to the membrane\n");
684             break;
685         case 2:
686             fprintf(stderr, "Taking z axis as normal to the membrane\n");
687             break;
688     }
689
690     /* tetraheder order parameter */
691     /* If either of the options is set we compute both */
692     nfxpm = opt2fns(&intfn, "-o", NFILE, fnm);
693     if (nfxpm != 2)
694     {
695         gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfxpm);
696     }
697     calc_tetra_order_interface(ndxfnm, tpsfnm, trxfnm, binwidth, nsttblock, &frames, &xslices, &yslices, sg1, sg2, &intfpos, oenv);
698     writesurftoxpms(intfpos, frames, xslices, yslices, binwidth, intfn, nlevels);
699
700     if (bFourier)
701     {
702         nfspect = opt2fns(&spectra, "-Spect", NFILE, fnm);
703         if (nfspect != 2)
704         {
705             gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfspect);
706         }
707         powerspectavg(intfpos, frames, xslices, yslices, spectra);
708     }
709
710     if (bRawOut)
711     {
712         nfraw = opt2fns(&raw, "-or", NFILE, fnm);
713         if (nfraw != 2)
714         {
715             gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfraw);
716         }
717         writeraw(intfpos, frames, xslices, yslices, raw);
718     }
719
720     return 0;
721 }