9e45a5ba7d157842597cdb64ab9d28dde809a6f1
[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-2018, The GROMACS development team.
5  * Copyright (c) 2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 #include "gmxpre.h"
38
39 #include <cmath>
40 #include <cstring>
41
42 #include "gromacs/commandline/pargs.h"
43 #include "gromacs/fileio/confio.h"
44 #include "gromacs/fileio/matio.h"
45 #include "gromacs/fileio/trxio.h"
46 #include "gromacs/fileio/xvgr.h"
47 #include "gromacs/gmxana/binsearch.h"
48 #include "gromacs/gmxana/gmx_ana.h"
49 #include "gromacs/gmxana/gstat.h"
50 #include "gromacs/gmxana/powerspect.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/topology/topology.h"
56 #include "gromacs/utility/arraysize.h"
57 #include "gromacs/utility/cstringutil.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/futil.h"
60 #include "gromacs/utility/gmxassert.h"
61 #include "gromacs/utility/smalloc.h"
62
63 static void find_tetra_order_grid(t_topology top,
64                                   PbcType    pbcType,
65                                   int        natoms,
66                                   matrix     box,
67                                   rvec       x[],
68                                   int        maxidx,
69                                   const int  index[],
70                                   real*      sgmean,
71                                   real*      skmean,
72                                   int        nslicex,
73                                   int        nslicey,
74                                   int        nslicez,
75                                   real***    sggrid,
76                                   real***    skgrid)
77 {
78     int         ix, jx, i, j, k, l, n, *nn[4];
79     rvec        dx, rj, rk, urk, urj;
80     real        cost, cost2, *sgmol, *skmol, rmean, rmean2, r2, box2, *r_nn[4];
81     t_pbc       pbc;
82     int         slindex_x, slindex_y, slindex_z;
83     int***      sl_count;
84     real        onethird = 1.0 / 3.0;
85     gmx_rmpbc_t gpbc;
86
87     /*  dmat = init_mat(maxidx, FALSE); */
88
89     box2 = box[XX][XX] * box[XX][XX];
90
91     /* Initialize expanded sl_count array */
92     snew(sl_count, nslicex);
93     for (i = 0; i < nslicex; i++)
94     {
95         snew(sl_count[i], nslicey);
96         for (j = 0; j < nslicey; j++)
97         {
98             snew(sl_count[i][j], nslicez);
99         }
100     }
101
102
103     for (i = 0; (i < 4); i++)
104     {
105         snew(r_nn[i], natoms);
106         snew(nn[i], natoms);
107
108         for (j = 0; (j < natoms); j++)
109         {
110             r_nn[i][j] = box2;
111         }
112     }
113
114     snew(sgmol, maxidx);
115     snew(skmol, maxidx);
116
117     /* Must init pbc every step because of pressure coupling */
118     set_pbc(&pbc, pbcType, box);
119     gpbc = gmx_rmpbc_init(&top.idef, pbcType, natoms);
120     gmx_rmpbc(gpbc, natoms, box, x);
121
122     *sgmean = 0.0;
123     *skmean = 0.0;
124     l       = 0;
125     for (i = 0; (i < maxidx); i++)
126     { /* loop over index file */
127         ix = index[i];
128         for (j = 0; (j < maxidx); j++)
129         {
130
131             if (i == j)
132             {
133                 continue;
134             }
135
136             jx = index[j];
137
138             pbc_dx(&pbc, x[ix], x[jx], dx);
139             r2 = iprod(dx, dx);
140
141             /* set_mat_entry(dmat,i,j,r2); */
142
143             /* determine the nearest neighbours */
144             if (r2 < r_nn[0][i])
145             {
146                 r_nn[3][i] = r_nn[2][i];
147                 nn[3][i]   = nn[2][i];
148                 r_nn[2][i] = r_nn[1][i];
149                 nn[2][i]   = nn[1][i];
150                 r_nn[1][i] = r_nn[0][i];
151                 nn[1][i]   = nn[0][i];
152                 r_nn[0][i] = r2;
153                 nn[0][i]   = j;
154             }
155             else if (r2 < r_nn[1][i])
156             {
157                 r_nn[3][i] = r_nn[2][i];
158                 nn[3][i]   = nn[2][i];
159                 r_nn[2][i] = r_nn[1][i];
160                 nn[2][i]   = nn[1][i];
161                 r_nn[1][i] = r2;
162                 nn[1][i]   = j;
163             }
164             else if (r2 < r_nn[2][i])
165             {
166                 r_nn[3][i] = r_nn[2][i];
167                 nn[3][i]   = nn[2][i];
168                 r_nn[2][i] = r2;
169                 nn[2][i]   = j;
170             }
171             else if (r2 < r_nn[3][i])
172             {
173                 r_nn[3][i] = r2;
174                 nn[3][i]   = j;
175             }
176         }
177
178
179         /* calculate mean distance between nearest neighbours */
180         rmean = 0;
181         for (j = 0; (j < 4); j++)
182         {
183             r_nn[j][i] = std::sqrt(r_nn[j][i]);
184             rmean += r_nn[j][i];
185         }
186         rmean /= 4;
187
188         n        = 0;
189         sgmol[i] = 0.0;
190         skmol[i] = 0.0;
191
192         /* Chau1998a eqn 3 */
193         /* angular part tetrahedrality order parameter per atom */
194         for (j = 0; (j < 3); j++)
195         {
196             for (k = j + 1; (k < 4); k++)
197             {
198                 pbc_dx(&pbc, x[ix], x[index[nn[k][i]]], rk);
199                 pbc_dx(&pbc, x[ix], x[index[nn[j][i]]], rj);
200
201                 unitv(rk, urk);
202                 unitv(rj, urj);
203
204                 cost  = iprod(urk, urj) + onethird;
205                 cost2 = cost * cost;
206
207                 sgmol[i] += cost2;
208                 l++;
209                 n++;
210             }
211         }
212         /* normalize sgmol between 0.0 and 1.0 */
213         sgmol[i] = 3 * sgmol[i] / 32;
214         *sgmean += sgmol[i];
215
216         /* distance part tetrahedrality order parameter per atom */
217         rmean2 = 4 * 3 * rmean * rmean;
218         for (j = 0; (j < 4); j++)
219         {
220             skmol[i] += (rmean - r_nn[j][i]) * (rmean - r_nn[j][i]) / rmean2;
221             /*      printf("%d %f (%f %f %f %f) \n",
222                     i, skmol[i], rmean, rmean2, r_nn[j][i], (rmean - r_nn[j][i]) );
223              */
224         }
225
226         *skmean += skmol[i];
227
228         /* Compute sliced stuff in x y z*/
229         slindex_x = static_cast<int>(std::round((1 + x[i][XX] / box[XX][XX]) * nslicex)) % nslicex;
230         slindex_y = static_cast<int>(std::round((1 + x[i][YY] / box[YY][YY]) * nslicey)) % nslicey;
231         slindex_z = static_cast<int>(std::round((1 + x[i][ZZ] / box[ZZ][ZZ]) * nslicez)) % nslicez;
232         sggrid[slindex_x][slindex_y][slindex_z] += sgmol[i];
233         skgrid[slindex_x][slindex_y][slindex_z] += skmol[i];
234         (sl_count[slindex_x][slindex_y][slindex_z])++;
235     } /* loop over entries in index file */
236
237     *sgmean /= maxidx;
238     *skmean /= maxidx;
239
240     for (i = 0; (i < nslicex); i++)
241     {
242         for (j = 0; j < nslicey; j++)
243         {
244             for (k = 0; k < nslicez; k++)
245             {
246                 if (sl_count[i][j][k] > 0)
247                 {
248                     sggrid[i][j][k] /= sl_count[i][j][k];
249                     skgrid[i][j][k] /= sl_count[i][j][k];
250                 }
251             }
252         }
253     }
254
255     sfree(sl_count);
256     sfree(sgmol);
257     sfree(skmol);
258     for (i = 0; (i < 4); i++)
259     {
260         sfree(r_nn[i]);
261         sfree(nn[i]);
262     }
263 }
264
265 /*Determines interface from tetrahedral order parameter in box with specified binwidth.  */
266 /*Outputs interface positions(bins), the number of timeframes, and the number of surface-mesh points in xy*/
267
268 static void calc_tetra_order_interface(const char*       fnNDX,
269                                        const char*       fnTPS,
270                                        const char*       fnTRX,
271                                        real              binw,
272                                        int               tblock,
273                                        int*              nframes,
274                                        int*              nslicex,
275                                        int*              nslicey,
276                                        real              sgang1,
277                                        real              sgang2,
278                                        real****          intfpos,
279                                        gmx_output_env_t* oenv)
280 {
281     FILE *       fpsg = nullptr, *fpsk = nullptr;
282     t_topology   top;
283     PbcType      pbcType;
284     t_trxstatus* status;
285     int          natoms;
286     real         t;
287     rvec *       xtop, *x;
288     matrix       box;
289     real         sg, sk, sgintf;
290     int**        index   = nullptr;
291     char**       grpname = nullptr;
292     int          i, j, k, n, *isize, ng, nslicez, framenr;
293     real ***sg_grid = nullptr, ***sk_grid = nullptr, ***sg_fravg = nullptr, ***sk_fravg = nullptr,
294          ****sk_4d = nullptr, ****sg_4d = nullptr;
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, &top, &pbcType, &xtop, nullptr, box, FALSE);
304
305     *nslicex = static_cast<int>(box[XX][XX] / binw + onehalf); /*Calculate slicenr from binwidth*/
306     *nslicey = static_cast<int>(box[YY][YY] / binw + onehalf);
307     nslicez  = static_cast<int>(box[ZZ][ZZ] / binw + onehalf);
308
309
310     ng = 1;
311     /* get index groups */
312     printf("Select the group that contains the atoms you want to use for the tetrahedrality order "
313            "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)", top.atoms.nr, natoms);
324     }
325     check_index(nullptr, ng, index[0], nullptr, natoms);
326
327
328     /*Prepare structures for temporary storage of frame info*/
329     snew(sg_grid, *nslicex);
330     snew(sk_grid, *nslicex);
331     for (i = 0; i < *nslicex; i++)
332     {
333         snew(sg_grid[i], *nslicey);
334         snew(sk_grid[i], *nslicey);
335         for (j = 0; j < *nslicey; j++)
336         {
337             snew(sg_grid[i][j], nslicez);
338             snew(sk_grid[i][j], nslicez);
339         }
340     }
341
342     sg_4d    = nullptr;
343     sk_4d    = nullptr;
344     *nframes = 0;
345     framenr  = 0;
346
347     /* Loop over frames*/
348     do
349     {
350         /*Initialize box meshes (temporary storage for each tblock frame -reinitialise every tblock steps */
351         if (framenr % tblock == 0)
352         {
353             srenew(sk_4d, *nframes + 1);
354             srenew(sg_4d, *nframes + 1);
355             snew(sg_fravg, *nslicex);
356             snew(sk_fravg, *nslicex);
357             for (i = 0; i < *nslicex; i++)
358             {
359                 snew(sg_fravg[i], *nslicey);
360                 snew(sk_fravg[i], *nslicey);
361                 for (j = 0; j < *nslicey; j++)
362                 {
363                     snew(sg_fravg[i][j], nslicez);
364                     snew(sk_fravg[i][j], nslicez);
365                 }
366             }
367         }
368
369         find_tetra_order_grid(top, pbcType, natoms, box, x, isize[0], index[0], &sg, &sk, *nslicex,
370                               *nslicey, nslicez, sg_grid, sk_grid);
371         GMX_RELEASE_ASSERT(sk_fravg != nullptr, "Trying to dereference NULL sk_fravg pointer");
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             GMX_RELEASE_ASSERT(sk_4d != nullptr, "Trying to dereference NULL sk_4d pointer");
389             sk_4d[*nframes] = sk_fravg;
390             sg_4d[*nframes] = sg_fravg;
391             (*nframes)++;
392         }
393
394     } while (read_next_x(oenv, status, &t, x, box));
395     close_trx(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)",
405                         "S\\sg\\N", oenv);
406         fpsk = xvgropen("sk_dist_mesh", "S\\sk\\N Distance Order Parameter / Meshpoint", "(nm)",
407                         "S\\sk\\N", oenv);
408         for (n = 0; n < (*nframes); n++)
409         {
410             fprintf(fpsg, "%i\n", n);
411             fprintf(fpsk, "%i\n", n);
412             for (i = 0; (i < *nslicex); i++)
413             {
414                 for (j = 0; j < *nslicey; j++)
415                 {
416                     for (k = 0; k < nslicez; k++)
417                     {
418                         fprintf(fpsg, "%4f  %4f  %4f  %8f\n", (i + 0.5) * box[XX][XX] / (*nslicex),
419                                 (j + 0.5) * box[YY][YY] / (*nslicey),
420                                 (k + 0.5) * box[ZZ][ZZ] / nslicez, sg_4d[n][i][j][k]);
421                         fprintf(fpsk, "%4f  %4f  %4f  %8f\n", (i + 0.5) * box[XX][XX] / (*nslicex),
422                                 (j + 0.5) * box[YY][YY] / (*nslicey),
423                                 (k + 0.5) * box[ZZ][ZZ] / nslicez, sk_4d[n][i][j][k]);
424                     }
425                 }
426             }
427         }
428         xvgrclose(fpsg);
429         xvgrclose(fpsk);
430     }
431
432
433     /* Find positions of interface z by scanning orderparam for each frame and for each xy-mesh cylinder along z*/
434
435     /*Simple trial: assume interface is in the middle of -sgang1 and sgang2*/
436     sgintf = 0.5 * (sgang1 + sgang2);
437
438
439     /*Allocate memory for interface arrays; */
440     snew((*intfpos), 2);
441     snew((*intfpos)[0], *nframes);
442     snew((*intfpos)[1], *nframes);
443
444     bins = (*nslicex) * (*nslicey);
445
446
447     snew(perm, nslicez); /*permutation array for sorting along normal coordinate*/
448
449
450     for (n = 0; n < *nframes; n++)
451     {
452         snew((*intfpos)[0][n], bins);
453         snew((*intfpos)[1][n], bins);
454         for (i = 0; i < *nslicex; i++)
455         {
456             for (j = 0; j < *nslicey; j++)
457             {
458                 rangeArray(perm, nslicez); /*reset permutation array to identity*/
459                 /*Binsearch returns 2 bin-numbers where the order param is  <= setpoint sgintf*/
460                 ndx1 = start_binsearch(sg_4d[n][i][j], perm, 0, nslicez / 2 - 1, sgintf, 1);
461                 ndx2 = start_binsearch(sg_4d[n][i][j], perm, nslicez / 2, nslicez - 1, sgintf, -1);
462                 /*Use linear interpolation to smooth out the interface position*/
463
464                 /*left interface (0)*/
465                 /*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){
466                    pos=( (sgintf-sg_4d[n][i][j][perm[ndx1]])*perm[ndx1+1]+(sg_4d[n][i][j][perm[ndx1+1]]-sgintf)*perm[ndx1 ])*/
467                 (*intfpos)[0][n][j + *nslicey * i] = (perm[ndx1] + onehalf) * binw;
468                 /*right interface (1)*/
469                 /*alpha=(sgintf-sg_4d[n][i][j][perm[ndx2]])/(sg_4d[n][i][j][perm[ndx2]+1]-sg_4d[n][i][j][perm[ndx2]]);*/
470                 /*(*intfpos)[1][n][j+*nslicey*i]=((1-alpha)*perm[ndx2]+alpha*(perm[ndx2]+1)+onehalf)*box[ZZ][ZZ]/nslicez;*/
471                 (*intfpos)[1][n][j + *nslicey * i] = (perm[ndx2] + onehalf) * binw;
472             }
473         }
474     }
475
476
477     sfree(sk_4d);
478     sfree(sg_4d);
479 }
480
481 static void writesurftoxpms(real***                          surf,
482                             int                              tblocks,
483                             int                              xbins,
484                             int                              ybins,
485                             real                             bw,
486                             gmx::ArrayRef<const std::string> outfiles,
487                             int                              maplevels)
488 {
489
490     char   numbuf[STRLEN];
491     int    n, i, j;
492     real **profile1, **profile2;
493     real   max1, max2, min1, min2, *xticks, *yticks;
494     t_rgb  lo = { 1, 1, 1 };
495     t_rgb  hi = { 0, 0, 0 };
496     FILE * xpmfile1, *xpmfile2;
497
498     /*Prepare xpm structures for output*/
499
500     /*Allocate memory to tick's and matrices*/
501     snew(xticks, xbins + 1);
502     snew(yticks, ybins + 1);
503
504     profile1 = mk_matrix(xbins, ybins, FALSE);
505     profile2 = mk_matrix(xbins, ybins, FALSE);
506
507     for (i = 0; i < xbins + 1; i++)
508     {
509         xticks[i] += bw;
510     }
511     for (j = 0; j < ybins + 1; j++)
512     {
513         yticks[j] += bw;
514     }
515
516     xpmfile1 = gmx_ffopen(outfiles[0], "w");
517     xpmfile2 = gmx_ffopen(outfiles[1], "w");
518
519     max1 = max2 = 0.0;
520     min1 = min2 = 1000.00;
521
522     for (n = 0; n < tblocks; n++)
523     {
524         sprintf(numbuf, "%5d", n);
525         /*Filling matrices for inclusion in xpm-files*/
526         for (i = 0; i < xbins; i++)
527         {
528             for (j = 0; j < ybins; j++)
529             {
530                 profile1[i][j] = (surf[0][n][j + ybins * i]);
531                 profile2[i][j] = (surf[1][n][j + ybins * i]);
532                 /*Finding max and min values*/
533                 if (profile1[i][j] > max1)
534                 {
535                     max1 = profile1[i][j];
536                 }
537                 if (profile1[i][j] < min1)
538                 {
539                     min1 = profile1[i][j];
540                 }
541                 if (profile2[i][j] > max2)
542                 {
543                     max2 = profile2[i][j];
544                 }
545                 if (profile2[i][j] < min2)
546                 {
547                     min2 = profile2[i][j];
548                 }
549             }
550         }
551
552         write_xpm(xpmfile1, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks,
553                   profile1, min1, max1, lo, hi, &maplevels);
554         write_xpm(xpmfile2, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks,
555                   profile2, min2, max2, lo, hi, &maplevels);
556     }
557
558     gmx_ffclose(xpmfile1);
559     gmx_ffclose(xpmfile2);
560
561
562     sfree(profile1);
563     sfree(profile2);
564     sfree(xticks);
565     sfree(yticks);
566 }
567
568
569 static void writeraw(real*** surf, int tblocks, int xbins, int ybins, gmx::ArrayRef<const std::string> fnms)
570 {
571     FILE *raw1, *raw2;
572     int   i, j, n;
573
574     raw1 = gmx_ffopen(fnms[0], "w");
575     raw2 = gmx_ffopen(fnms[1], "w");
576     fprintf(raw1, "#Legend\n#TBlock\n#Xbin Ybin Z t\n");
577     fprintf(raw2, "#Legend\n#TBlock\n#Xbin Ybin Z t\n");
578     for (n = 0; n < tblocks; n++)
579     {
580         fprintf(raw1, "%5d\n", n);
581         fprintf(raw2, "%5d\n", n);
582         for (i = 0; i < xbins; i++)
583         {
584             for (j = 0; j < ybins; j++)
585             {
586                 fprintf(raw1, "%i  %i  %8.5f\n", i, j, (surf[0][n][j + ybins * i]));
587                 fprintf(raw2, "%i  %i  %8.5f\n", i, j, (surf[1][n][j + ybins * i]));
588             }
589         }
590     }
591
592     gmx_ffclose(raw1);
593     gmx_ffclose(raw2);
594 }
595
596
597 int gmx_hydorder(int argc, char* argv[])
598 {
599     static const char* desc[] = {
600         "[THISMODULE] computes the tetrahedrality order parameters around a ",
601         "given atom. Both angle an distance order parameters are calculated. See",
602         "P.-L. Chau and A.J. Hardwick, Mol. Phys., 93, (1998), 511-518.",
603         "for more details.[PAR]",
604         "[THISMODULE] calculates the order parameter in a 3d-mesh in the box, and",
605         "with 2 phases in the box gives the user the option to define a 2D interface in time",
606         "separating the faces by specifying parameters [TT]-sgang1[tt] and",
607         "[TT]-sgang2[tt] (it is important to select these judiciously)."
608     };
609
610     int                axis      = 0;
611     static int         nsttblock = 1;
612     static int         nlevels   = 100;
613     static real        binwidth  = 1.0; /* binwidth in mesh           */
614     static real        sg1       = 1;
615     static real        sg2       = 1; /* order parameters for bulk phases */
616     static gmx_bool    bFourier  = FALSE;
617     static gmx_bool    bRawOut   = FALSE;
618     int                frames, xslices, yslices; /* Dimensions of interface arrays*/
619     real***            intfpos; /* Interface arrays (intfnr,t,xy) -potentially large */
620     static const char* normal_axis[] = { nullptr, "z", "x", "y", nullptr };
621
622     t_pargs pa[] = {
623         { "-d", FALSE, etENUM, { normal_axis }, "Direction of the normal on the membrane" },
624         { "-bw", FALSE, etREAL, { &binwidth }, "Binwidth of box mesh" },
625         { "-sgang1", FALSE, etREAL, { &sg1 }, "tetrahedral angle parameter in Phase 1 (bulk)" },
626         { "-sgang2", FALSE, etREAL, { &sg2 }, "tetrahedral angle parameter in Phase 2 (bulk)" },
627         { "-tblock", FALSE, etINT, { &nsttblock }, "Number of frames in one time-block average" },
628         { "-nlevel", FALSE, etINT, { &nlevels }, "Number of Height levels in 2D - XPixMaps" }
629     };
630
631     t_filenm fnm[] = {
632         /* files for g_order    */
633         { efTRX, "-f", nullptr, ffREAD },              /* trajectory file              */
634         { efNDX, "-n", nullptr, ffREAD },              /* index file           */
635         { efTPR, "-s", nullptr, ffREAD },              /* topology file                */
636         { efXPM, "-o", "intf", ffWRMULT },             /* XPM- surface maps     */
637         { efOUT, "-or", "raw", ffOPTWRMULT },          /* xvgr output file           */
638         { efOUT, "-Spect", "intfspect", ffOPTWRMULT }, /* Fourier spectrum interfaces */
639     };
640 #define NFILE asize(fnm)
641
642     /*Filenames*/
643     const char *      ndxfnm, *tpsfnm, *trxfnm;
644     gmx_output_env_t* oenv;
645
646     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME, NFILE, fnm, asize(pa), pa,
647                            asize(desc), desc, 0, nullptr, &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(efTPR, NFILE, fnm);
661     trxfnm = ftp2fn(efTRX, NFILE, fnm);
662
663     /* Calculate axis */
664     GMX_RELEASE_ASSERT(normal_axis[0] != nullptr,
665                        "Option setting inconsistency; normal_axis[0] is NULL");
666     if (std::strcmp(normal_axis[0], "x") == 0)
667     {
668         axis = XX;
669     }
670     else if (std::strcmp(normal_axis[0], "y") == 0)
671     {
672         axis = YY;
673     }
674     else if (std::strcmp(normal_axis[0], "z") == 0)
675     {
676         axis = ZZ;
677     }
678     else
679     {
680         gmx_fatal(FARGS, "Invalid axis, use x, y or z");
681     }
682
683     switch (axis)
684     {
685         case 0: fprintf(stderr, "Taking x axis as normal to the membrane\n"); break;
686         case 1: fprintf(stderr, "Taking y axis as normal to the membrane\n"); break;
687         case 2: fprintf(stderr, "Taking z axis as normal to the membrane\n"); break;
688     }
689
690     /* tetraheder order parameter */
691     /* If either of the options is set we compute both */
692     gmx::ArrayRef<const std::string> intfn = opt2fns("-o", NFILE, fnm);
693     if (intfn.size() != 2)
694     {
695         gmx_fatal(FARGS, "No or not correct number (2) of output-files: %td", intfn.ssize());
696     }
697     calc_tetra_order_interface(ndxfnm, tpsfnm, trxfnm, binwidth, nsttblock, &frames, &xslices,
698                                &yslices, sg1, sg2, &intfpos, oenv);
699     writesurftoxpms(intfpos, frames, xslices, yslices, binwidth, intfn, nlevels);
700
701     if (bFourier)
702     {
703         gmx::ArrayRef<const std::string> spectra = opt2fns("-Spect", NFILE, fnm);
704         if (spectra.size() != 2)
705         {
706             gmx_fatal(FARGS, "No or not correct number (2) of output-files: %td", spectra.ssize());
707         }
708         powerspectavg(intfpos, frames, xslices, yslices, spectra);
709     }
710
711     if (bRawOut)
712     {
713         gmx::ArrayRef<const std::string> raw = opt2fns("-or", NFILE, fnm);
714         if (raw.size() != 2)
715         {
716             gmx_fatal(FARGS, "No or not correct number (2) of output-files: %td", raw.ssize());
717         }
718         writeraw(intfpos, frames, xslices, yslices, raw);
719     }
720
721     return 0;
722 }