Apply re-formatting to C++ in src/ tree.
[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(
370                 top, pbcType, natoms, box, x, isize[0], index[0], &sg, &sk, *nslicex, *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(
405                 "sg_ang_mesh", "S\\sg\\N Angle Order Parameter / Meshpoint", "(nm)", "S\\sg\\N", oenv);
406         fpsk = xvgropen(
407                 "sk_dist_mesh", "S\\sk\\N Distance Order Parameter / Meshpoint", "(nm)", "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,
419                                 "%4f  %4f  %4f  %8f\n",
420                                 (i + 0.5) * box[XX][XX] / (*nslicex),
421                                 (j + 0.5) * box[YY][YY] / (*nslicey),
422                                 (k + 0.5) * box[ZZ][ZZ] / nslicez,
423                                 sg_4d[n][i][j][k]);
424                         fprintf(fpsk,
425                                 "%4f  %4f  %4f  %8f\n",
426                                 (i + 0.5) * box[XX][XX] / (*nslicex),
427                                 (j + 0.5) * box[YY][YY] / (*nslicey),
428                                 (k + 0.5) * box[ZZ][ZZ] / nslicez,
429                                 sk_4d[n][i][j][k]);
430                     }
431                 }
432             }
433         }
434         xvgrclose(fpsg);
435         xvgrclose(fpsk);
436     }
437
438
439     /* Find positions of interface z by scanning orderparam for each frame and for each xy-mesh cylinder along z*/
440
441     /*Simple trial: assume interface is in the middle of -sgang1 and sgang2*/
442     sgintf = 0.5 * (sgang1 + sgang2);
443
444
445     /*Allocate memory for interface arrays; */
446     snew((*intfpos), 2);
447     snew((*intfpos)[0], *nframes);
448     snew((*intfpos)[1], *nframes);
449
450     bins = (*nslicex) * (*nslicey);
451
452
453     snew(perm, nslicez); /*permutation array for sorting along normal coordinate*/
454
455
456     for (n = 0; n < *nframes; n++)
457     {
458         snew((*intfpos)[0][n], bins);
459         snew((*intfpos)[1][n], bins);
460         for (i = 0; i < *nslicex; i++)
461         {
462             for (j = 0; j < *nslicey; j++)
463             {
464                 rangeArray(perm, nslicez); /*reset permutation array to identity*/
465                 /*Binsearch returns 2 bin-numbers where the order param is  <= setpoint sgintf*/
466                 ndx1 = start_binsearch(sg_4d[n][i][j], perm, 0, nslicez / 2 - 1, sgintf, 1);
467                 ndx2 = start_binsearch(sg_4d[n][i][j], perm, nslicez / 2, nslicez - 1, sgintf, -1);
468                 /*Use linear interpolation to smooth out the interface position*/
469
470                 /*left interface (0)*/
471                 /*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){
472                    pos=( (sgintf-sg_4d[n][i][j][perm[ndx1]])*perm[ndx1+1]+(sg_4d[n][i][j][perm[ndx1+1]]-sgintf)*perm[ndx1 ])*/
473                 (*intfpos)[0][n][j + *nslicey * i] = (perm[ndx1] + onehalf) * binw;
474                 /*right interface (1)*/
475                 /*alpha=(sgintf-sg_4d[n][i][j][perm[ndx2]])/(sg_4d[n][i][j][perm[ndx2]+1]-sg_4d[n][i][j][perm[ndx2]]);*/
476                 /*(*intfpos)[1][n][j+*nslicey*i]=((1-alpha)*perm[ndx2]+alpha*(perm[ndx2]+1)+onehalf)*box[ZZ][ZZ]/nslicez;*/
477                 (*intfpos)[1][n][j + *nslicey * i] = (perm[ndx2] + onehalf) * binw;
478             }
479         }
480     }
481
482
483     sfree(sk_4d);
484     sfree(sg_4d);
485 }
486
487 static void writesurftoxpms(real***                          surf,
488                             int                              tblocks,
489                             int                              xbins,
490                             int                              ybins,
491                             real                             bw,
492                             gmx::ArrayRef<const std::string> outfiles,
493                             int                              maplevels)
494 {
495
496     char   numbuf[STRLEN];
497     int    n, i, j;
498     real **profile1, **profile2;
499     real   max1, max2, min1, min2, *xticks, *yticks;
500     t_rgb  lo = { 1, 1, 1 };
501     t_rgb  hi = { 0, 0, 0 };
502     FILE * xpmfile1, *xpmfile2;
503
504     /*Prepare xpm structures for output*/
505
506     /*Allocate memory to tick's and matrices*/
507     snew(xticks, xbins + 1);
508     snew(yticks, ybins + 1);
509
510     profile1 = mk_matrix(xbins, ybins, FALSE);
511     profile2 = mk_matrix(xbins, ybins, FALSE);
512
513     for (i = 0; i < xbins + 1; i++)
514     {
515         xticks[i] += bw;
516     }
517     for (j = 0; j < ybins + 1; j++)
518     {
519         yticks[j] += bw;
520     }
521
522     xpmfile1 = gmx_ffopen(outfiles[0], "w");
523     xpmfile2 = gmx_ffopen(outfiles[1], "w");
524
525     max1 = max2 = 0.0;
526     min1 = min2 = 1000.00;
527
528     for (n = 0; n < tblocks; n++)
529     {
530         sprintf(numbuf, "%5d", n);
531         /*Filling matrices for inclusion in xpm-files*/
532         for (i = 0; i < xbins; i++)
533         {
534             for (j = 0; j < ybins; j++)
535             {
536                 profile1[i][j] = (surf[0][n][j + ybins * i]);
537                 profile2[i][j] = (surf[1][n][j + ybins * i]);
538                 /*Finding max and min values*/
539                 if (profile1[i][j] > max1)
540                 {
541                     max1 = profile1[i][j];
542                 }
543                 if (profile1[i][j] < min1)
544                 {
545                     min1 = profile1[i][j];
546                 }
547                 if (profile2[i][j] > max2)
548                 {
549                     max2 = profile2[i][j];
550                 }
551                 if (profile2[i][j] < min2)
552                 {
553                     min2 = profile2[i][j];
554                 }
555             }
556         }
557
558         write_xpm(xpmfile1, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile1, min1, max1, lo, hi, &maplevels);
559         write_xpm(xpmfile2, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile2, min2, max2, lo, hi, &maplevels);
560     }
561
562     gmx_ffclose(xpmfile1);
563     gmx_ffclose(xpmfile2);
564
565
566     sfree(profile1);
567     sfree(profile2);
568     sfree(xticks);
569     sfree(yticks);
570 }
571
572
573 static void writeraw(real*** surf, int tblocks, int xbins, int ybins, gmx::ArrayRef<const std::string> fnms)
574 {
575     FILE *raw1, *raw2;
576     int   i, j, n;
577
578     raw1 = gmx_ffopen(fnms[0], "w");
579     raw2 = gmx_ffopen(fnms[1], "w");
580     fprintf(raw1, "#Legend\n#TBlock\n#Xbin Ybin Z t\n");
581     fprintf(raw2, "#Legend\n#TBlock\n#Xbin Ybin Z t\n");
582     for (n = 0; n < tblocks; n++)
583     {
584         fprintf(raw1, "%5d\n", n);
585         fprintf(raw2, "%5d\n", n);
586         for (i = 0; i < xbins; i++)
587         {
588             for (j = 0; j < ybins; j++)
589             {
590                 fprintf(raw1, "%i  %i  %8.5f\n", i, j, (surf[0][n][j + ybins * i]));
591                 fprintf(raw2, "%i  %i  %8.5f\n", i, j, (surf[1][n][j + ybins * i]));
592             }
593         }
594     }
595
596     gmx_ffclose(raw1);
597     gmx_ffclose(raw2);
598 }
599
600
601 int gmx_hydorder(int argc, char* argv[])
602 {
603     static const char* desc[] = {
604         "[THISMODULE] computes the tetrahedrality order parameters around a ",
605         "given atom. Both angle an distance order parameters are calculated. See",
606         "P.-L. Chau and A.J. Hardwick, Mol. Phys., 93, (1998), 511-518.",
607         "for more details.[PAR]",
608         "[THISMODULE] calculates the order parameter in a 3d-mesh in the box, and",
609         "with 2 phases in the box gives the user the option to define a 2D interface in time",
610         "separating the faces by specifying parameters [TT]-sgang1[tt] and",
611         "[TT]-sgang2[tt] (it is important to select these judiciously)."
612     };
613
614     int                axis      = 0;
615     static int         nsttblock = 1;
616     static int         nlevels   = 100;
617     static real        binwidth  = 1.0; /* binwidth in mesh           */
618     static real        sg1       = 1;
619     static real        sg2       = 1; /* order parameters for bulk phases */
620     static gmx_bool    bFourier  = FALSE;
621     static gmx_bool    bRawOut   = FALSE;
622     int                frames, xslices, yslices; /* Dimensions of interface arrays*/
623     real***            intfpos; /* Interface arrays (intfnr,t,xy) -potentially large */
624     static const char* normal_axis[] = { nullptr, "z", "x", "y", nullptr };
625
626     t_pargs pa[] = {
627         { "-d", FALSE, etENUM, { normal_axis }, "Direction of the normal on the membrane" },
628         { "-bw", FALSE, etREAL, { &binwidth }, "Binwidth of box mesh" },
629         { "-sgang1", FALSE, etREAL, { &sg1 }, "tetrahedral angle parameter in Phase 1 (bulk)" },
630         { "-sgang2", FALSE, etREAL, { &sg2 }, "tetrahedral angle parameter in Phase 2 (bulk)" },
631         { "-tblock", FALSE, etINT, { &nsttblock }, "Number of frames in one time-block average" },
632         { "-nlevel", FALSE, etINT, { &nlevels }, "Number of Height levels in 2D - XPixMaps" }
633     };
634
635     t_filenm fnm[] = {
636         /* files for g_order    */
637         { efTRX, "-f", nullptr, ffREAD },              /* trajectory file              */
638         { efNDX, "-n", nullptr, ffREAD },              /* index file           */
639         { efTPR, "-s", nullptr, ffREAD },              /* topology file                */
640         { efXPM, "-o", "intf", ffWRMULT },             /* XPM- surface maps     */
641         { efOUT, "-or", "raw", ffOPTWRMULT },          /* xvgr output file           */
642         { efOUT, "-Spect", "intfspect", ffOPTWRMULT }, /* Fourier spectrum interfaces */
643     };
644 #define NFILE asize(fnm)
645
646     /*Filenames*/
647     const char *      ndxfnm, *tpsfnm, *trxfnm;
648     gmx_output_env_t* oenv;
649
650     if (!parse_common_args(
651                 &argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
652     {
653         return 0;
654     }
655     bFourier = opt2bSet("-Spect", NFILE, fnm);
656     bRawOut  = opt2bSet("-or", NFILE, fnm);
657
658     if (binwidth < 0.0)
659     {
660         gmx_fatal(FARGS, "Can not have binwidth < 0");
661     }
662
663     ndxfnm = ftp2fn(efNDX, NFILE, fnm);
664     tpsfnm = ftp2fn(efTPR, NFILE, fnm);
665     trxfnm = ftp2fn(efTRX, NFILE, fnm);
666
667     /* Calculate axis */
668     GMX_RELEASE_ASSERT(normal_axis[0] != nullptr,
669                        "Option setting inconsistency; normal_axis[0] is NULL");
670     if (std::strcmp(normal_axis[0], "x") == 0)
671     {
672         axis = XX;
673     }
674     else if (std::strcmp(normal_axis[0], "y") == 0)
675     {
676         axis = YY;
677     }
678     else if (std::strcmp(normal_axis[0], "z") == 0)
679     {
680         axis = ZZ;
681     }
682     else
683     {
684         gmx_fatal(FARGS, "Invalid axis, use x, y or z");
685     }
686
687     switch (axis)
688     {
689         case 0: fprintf(stderr, "Taking x axis as normal to the membrane\n"); break;
690         case 1: fprintf(stderr, "Taking y axis as normal to the membrane\n"); break;
691         case 2: fprintf(stderr, "Taking z axis as normal to the membrane\n"); break;
692     }
693
694     /* tetraheder order parameter */
695     /* If either of the options is set we compute both */
696     gmx::ArrayRef<const std::string> intfn = opt2fns("-o", NFILE, fnm);
697     if (intfn.size() != 2)
698     {
699         gmx_fatal(FARGS, "No or not correct number (2) of output-files: %td", intfn.ssize());
700     }
701     calc_tetra_order_interface(
702             ndxfnm, tpsfnm, trxfnm, binwidth, nsttblock, &frames, &xslices, &yslices, sg1, sg2, &intfpos, oenv);
703     writesurftoxpms(intfpos, frames, xslices, yslices, binwidth, intfn, nlevels);
704
705     if (bFourier)
706     {
707         gmx::ArrayRef<const std::string> spectra = opt2fns("-Spect", NFILE, fnm);
708         if (spectra.size() != 2)
709         {
710             gmx_fatal(FARGS, "No or not correct number (2) of output-files: %td", spectra.ssize());
711         }
712         powerspectavg(intfpos, frames, xslices, yslices, spectra);
713     }
714
715     if (bRawOut)
716     {
717         gmx::ArrayRef<const std::string> raw = opt2fns("-or", NFILE, fnm);
718         if (raw.size() != 2)
719         {
720             gmx_fatal(FARGS, "No or not correct number (2) of output-files: %td", raw.ssize());
721         }
722         writeraw(intfpos, frames, xslices, yslices, raw);
723     }
724
725     return 0;
726 }