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