Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_hydorder.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  * $Id: densorder.c,v 0.9
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  *          GROningen MAchine for Chemical Simulations
9  *
10  *                        VERSION 3.0
11  *
12  * Copyright (c) 1991-2001
13  * BIOSON Research Institute, Dept. of Biophysical Chemistry
14  * University of Groningen, The Netherlands
15  *
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  *
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  *
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  *
31  * Do check out http://www.gromacs.org , or mail us at gromacs@gromacs.org .
32  *
33  * And Hey:
34  * Gyas ROwers Mature At Cryogenic Speed
35  */
36
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42 #include <ctype.h>
43
44 #include "sysstuff.h"
45 #include <string.h>
46 #include "typedefs.h"
47 #include "statutil.h"
48 #include "smalloc.h"
49 #include "macros.h"
50 #include "gstat.h"
51 #include "vec.h"
52 #include "xvgr.h"
53 #include "pbc.h"
54 #include "copyrite.h"
55 #include "futil.h"
56 #include "statutil.h"
57 #include "index.h"
58 #include "tpxio.h"
59 #include "matio.h"
60 #include "binsearch.h"
61 #include "powerspect.h"
62
63 /* Print name of first atom in all groups in index file */
64 static void print_types(atom_id index[], atom_id a[], int ngrps,
65                         char *groups[], t_topology *top)
66 {
67     int i;
68
69     fprintf(stderr, "Using following groups: \n");
70     for (i = 0; i < ngrps; i++)
71     {
72         fprintf(stderr, "Groupname: %s First atomname: %s First atomnr %u\n",
73                 groups[i], *(top->atoms.atomname[a[index[i]]]), a[index[i]]);
74     }
75     fprintf(stderr, "\n");
76 }
77
78 static void check_length(real length, int a, int b)
79 {
80     if (length > 0.3)
81     {
82         fprintf(stderr, "WARNING: distance between atoms %d and "
83                 "%d > 0.3 nm (%f). Index file might be corrupt.\n",
84                 a, b, length);
85     }
86 }
87
88 static void find_tetra_order_grid(t_topology top, int ePBC,
89                                   int natoms, matrix box,
90                                   rvec x[], int maxidx, atom_id index[],
91                                   real time, real *sgmean, real *skmean,
92                                   int nslicex, int nslicey, int nslicez,
93                                   real ***sggrid, real ***skgrid)
94 {
95     int         ix, jx, i, j, k, l, n, *nn[4];
96     rvec        dx, rj, rk, urk, urj;
97     real        cost, cost2, *sgmol, *skmol, rmean, rmean2, r2, box2, *r_nn[4];
98     t_pbc       pbc;
99     int         slindex_x, slindex_y, slindex_z;
100     int      ***sl_count;
101     real        onethird = 1.0/3.0;
102     gmx_rmpbc_t gpbc;
103
104     /*  dmat = init_mat(maxidx, FALSE); */
105
106     box2 = box[XX][XX] * box[XX][XX];
107
108     /* Initialize expanded sl_count array */
109     snew(sl_count, nslicex);
110     for (i = 0; i < nslicex; i++)
111     {
112         snew(sl_count[i], nslicey);
113         for (j = 0; j < nslicey; j++)
114         {
115             snew(sl_count[i][j], nslicez);
116         }
117     }
118
119
120     for (i = 0; (i < 4); i++)
121     {
122         snew(r_nn[i], natoms);
123         snew(nn[i], natoms);
124
125         for (j = 0; (j < natoms); j++)
126         {
127             r_nn[i][j] = box2;
128         }
129     }
130
131     snew(sgmol, maxidx);
132     snew(skmol, maxidx);
133
134     /* Must init pbc every step because of pressure coupling */
135     set_pbc(&pbc,ePBC,box);
136     gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms, box);
137     gmx_rmpbc(gpbc, natoms, box, x);
138
139     *sgmean = 0.0;
140     *skmean = 0.0;
141     l       = 0;
142     for (i = 0; (i < maxidx); i++)
143     {   /* loop over index file */
144         ix = index[i];
145         for (j = 0; (j < maxidx); j++)
146         {
147
148             if (i == j)
149             {
150                 continue;
151             }
152
153             jx = index[j];
154
155             pbc_dx(&pbc, x[ix], x[jx], dx);
156             r2 = iprod(dx, dx);
157
158             /* set_mat_entry(dmat,i,j,r2); */
159
160             /* determine the nearest neighbours */
161             if (r2 < r_nn[0][i])
162             {
163                 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
164                 r_nn[2][i] = r_nn[1][i]; nn[2][i] = nn[1][i];
165                 r_nn[1][i] = r_nn[0][i]; nn[1][i] = nn[0][i];
166                 r_nn[0][i] = r2;         nn[0][i] = j;
167             }
168             else if (r2 < r_nn[1][i])
169             {
170                 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
171                 r_nn[2][i] = r_nn[1][i]; nn[2][i] = nn[1][i];
172                 r_nn[1][i] = r2;         nn[1][i] = j;
173             }
174             else if (r2 < r_nn[2][i])
175             {
176                 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
177                 r_nn[2][i] = r2;         nn[2][i] = j;
178             }
179             else if (r2 < r_nn[3][i])
180             {
181                 r_nn[3][i] = r2;         nn[3][i] = j;
182             }
183         }
184
185
186         /* calculate mean distance between nearest neighbours */
187         rmean = 0;
188         for (j = 0; (j < 4); j++)
189         {
190             r_nn[j][i] = sqrt(r_nn[j][i]);
191             rmean     += r_nn[j][i];
192         }
193         rmean /= 4;
194
195         n        = 0;
196         sgmol[i] = 0.0;
197         skmol[i] = 0.0;
198
199         /* Chau1998a eqn 3 */
200         /* angular part tetrahedrality order parameter per atom */
201         for (j = 0; (j < 3); j++)
202         {
203             for (k = j+1; (k < 4); k++)
204             {
205                 pbc_dx(&pbc, x[ix], x[index[nn[k][i]]], rk);
206                 pbc_dx(&pbc, x[ix], x[index[nn[j][i]]], rj);
207
208                 unitv(rk, urk);
209                 unitv(rj, urj);
210
211                 cost  = iprod(urk, urj) + onethird;
212                 cost2 = cost * cost;
213
214                 sgmol[i] += cost2;
215                 l++;
216                 n++;
217             }
218         }
219         /* normalize sgmol between 0.0 and 1.0 */
220         sgmol[i] = 3*sgmol[i]/32;
221         *sgmean += sgmol[i];
222
223         /* distance part tetrahedrality order parameter per atom */
224         rmean2 = 4 * 3 * rmean * rmean;
225         for (j = 0; (j < 4); j++)
226         {
227             skmol[i] += (rmean - r_nn[j][i]) * (rmean - r_nn[j][i]) / rmean2;
228             /*      printf("%d %f (%f %f %f %f) \n",
229                     i, skmol[i], rmean, rmean2, r_nn[j][i], (rmean - r_nn[j][i]) );
230              */
231         }
232
233         *skmean += skmol[i];
234
235         /* Compute sliced stuff in x y z*/
236         slindex_x = gmx_nint((1+x[i][XX]/box[XX][XX])*nslicex) % nslicex;
237         slindex_y = gmx_nint((1+x[i][YY]/box[YY][YY])*nslicey) % nslicey;
238         slindex_z = gmx_nint((1+x[i][ZZ]/box[ZZ][ZZ])*nslicez) % nslicez;
239         sggrid[slindex_x][slindex_y][slindex_z] += sgmol[i];
240         skgrid[slindex_x][slindex_y][slindex_z] += skmol[i];
241         (sl_count[slindex_x][slindex_y][slindex_z])++;
242     } /* loop over entries in index file */
243
244     *sgmean /= maxidx;
245     *skmean /= maxidx;
246
247     for (i = 0; (i < nslicex); i++)
248     {
249         for (j = 0; j < nslicey; j++)
250         {
251             for (k = 0; k < nslicez; k++)
252             {
253                 if (sl_count[i][j][k] > 0)
254                 {
255                     sggrid[i][j][k] /= sl_count[i][j][k];
256                     skgrid[i][j][k] /= sl_count[i][j][k];
257                 }
258             }
259         }
260     }
261
262     sfree(sl_count);
263     sfree(sgmol);
264     sfree(skmol);
265     for (i = 0; (i < 4); i++)
266     {
267         sfree(r_nn[i]);
268         sfree(nn[i]);
269     }
270 }
271
272 /*Determines interface from tetrahedral order parameter in box with specified binwidth.  */
273 /*Outputs interface positions(bins), the number of timeframes, and the number of surface-mesh points in xy*/
274
275 static void calc_tetra_order_interface(const char *fnNDX, const char *fnTPS, const char *fnTRX, real binw, int tblock,
276                                        int *nframes,  int *nslicex, int *nslicey,
277                                        real sgang1, real sgang2, real ****intfpos,
278                                        output_env_t oenv)
279 {
280     FILE         *fpsg   = NULL, *fpsk = NULL;
281     char         *sgslfn = "sg_ang_mesh";  /* Hardcoded filenames for debugging*/
282     char         *skslfn = "sk_dist_mesh";
283     t_topology    top;
284     int           ePBC;
285     char          title[STRLEN], subtitle[STRLEN];
286     t_trxstatus  *status;
287     int           natoms;
288     real          t;
289     rvec         *xtop, *x;
290     matrix        box;
291     real          sg, sk, sgintf, pos;
292     atom_id     **index   = NULL;
293     char        **grpname = NULL;
294     int           i, j, k, n, *isize, ng, nslicez, framenr;
295     real       ***sg_grid = NULL, ***sk_grid = NULL, ***sg_fravg = NULL, ***sk_fravg = NULL, ****sk_4d = NULL, ****sg_4d = NULL;
296     int          *perm;
297     int           ndx1, ndx2;
298     int           bins;
299     const real    onehalf = 1.0/2.0;
300     /* real   ***intfpos[2]; pointers to arrays of two interface positions zcoord(framenr,xbin,ybin): intfpos[interface_index][t][nslicey*x+y]
301      * i.e 1D Row-major order in (t,x,y) */
302
303
304     read_tps_conf(fnTPS, title, &top, &ePBC, &xtop, NULL, box, FALSE);
305
306     *nslicex = (int)(box[XX][XX]/binw + onehalf); /*Calculate slicenr from binwidth*/
307     *nslicey = (int)(box[YY][YY]/binw + onehalf);
308     nslicez  = (int)(box[ZZ][ZZ]/binw +  onehalf);
309
310
311
312     ng = 1;
313     /* get index groups */
314     printf("Select the group that contains the atoms you want to use for the tetrahedrality order parameter calculation:\n");
315     snew(grpname, ng);
316     snew(index, ng);
317     snew(isize, ng);
318     get_index(&top.atoms, fnNDX, ng, isize, index, grpname);
319
320     /* Analyze trajectory */
321     natoms = read_first_x(oenv, &status, fnTRX, &t, &x, box);
322     if (natoms > top.atoms.nr)
323     {
324         gmx_fatal(FARGS, "Topology (%d atoms) does not match trajectory (%d atoms)",
325                   top.atoms.nr, natoms);
326     }
327     check_index(NULL, ng, index[0], NULL, natoms);
328
329
330     /*Prepare structures for temporary storage of frame info*/
331     snew(sg_grid, *nslicex);
332     snew(sk_grid, *nslicex);
333     for (i = 0; i < *nslicex; i++)
334     {
335         snew(sg_grid[i], *nslicey);
336         snew(sk_grid[i], *nslicey);
337         for (j = 0; j < *nslicey; j++)
338         {
339             snew(sg_grid[i][j], nslicez);
340             snew(sk_grid[i][j], nslicez);
341         }
342     }
343
344     sg_4d    = NULL;
345     sk_4d    = NULL;
346     *nframes = 0;
347     framenr  = 0;
348
349 /* Loop over frames*/
350     do
351     {
352         /*Initialize box meshes (temporary storage for each tblock frame -reinitialise every tblock steps */
353         if (framenr%tblock == 0)
354         {
355             srenew(sk_4d, *nframes+1);
356             srenew(sg_4d, *nframes+1);
357             snew(sg_fravg, *nslicex);
358             snew(sk_fravg, *nslicex);
359             for (i = 0; i < *nslicex; i++)
360             {
361                 snew(sg_fravg[i], *nslicey);
362                 snew(sk_fravg[i], *nslicey);
363                 for (j = 0; j < *nslicey; j++)
364                 {
365                     snew(sg_fravg[i][j], nslicez);
366                     snew(sk_fravg[i][j], nslicez);
367                 }
368             }
369         }
370
371         find_tetra_order_grid(top, ePBC, natoms, box, x, isize[0], index[0], t,
372                               &sg, &sk, *nslicex, *nslicey, nslicez, sg_grid, sk_grid);
373         for (i = 0; i < *nslicex; i++)
374         {
375             for (j = 0; j < *nslicey; j++)
376             {
377                 for (k = 0; k < nslicez; k++)
378                 {
379                     sk_fravg[i][j][k] += sk_grid[i][j][k]/tblock;
380                     sg_fravg[i][j][k] += sg_grid[i][j][k]/tblock;
381                 }
382             }
383         }
384
385         framenr++;
386
387         if (framenr%tblock == 0)
388         {
389             sk_4d[*nframes] = sk_fravg;
390             sg_4d[*nframes] = sg_fravg;
391             (*nframes)++;
392         }
393
394     }
395     while (read_next_x(oenv, status, &t, natoms, x, box));
396     close_trj(status);
397
398     sfree(grpname);
399     sfree(index);
400     sfree(isize);
401
402     /*Debugging for printing out the entire order parameter meshes.*/
403     if (debug)
404     {
405         fpsg = xvgropen(sgslfn, "S\\sg\\N Angle Order Parameter / Meshpoint", "(nm)", "S\\sg\\N", oenv);
406         fpsk = xvgropen(skslfn, "S\\sk\\N Distance Order Parameter / Meshpoint", "(nm)", "S\\sk\\N", oenv);
407         for (n = 0; n < (*nframes); n++)
408         {
409             fprintf(fpsg, "%i\n", n);
410             fprintf(fpsk, "%i\n", n);
411             for (i = 0; (i < *nslicex); i++)
412             {
413                 for (j = 0; j < *nslicey; j++)
414                 {
415                     for (k = 0; k < nslicez; k++)
416                     {
417                         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]);
418                         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]);
419                     }
420                 }
421             }
422         }
423         xvgrclose(fpsg);
424         xvgrclose(fpsk);
425     }
426
427
428     /* Find positions of interface z by scanning orderparam for each frame and for each xy-mesh cylinder along z*/
429
430     /*Simple trial: assume interface is in the middle of -sgang1 and sgang2*/
431     sgintf = 0.5*(sgang1+sgang2);
432
433
434     /*Allocate memory for interface arrays; */
435     snew((*intfpos), 2);
436     snew((*intfpos)[0], *nframes);
437     snew((*intfpos)[1], *nframes);
438
439     bins = (*nslicex)*(*nslicey);
440
441
442     snew(perm, nslicez);  /*permutation array for sorting along normal coordinate*/
443
444
445     for (n = 0; n < *nframes; n++)
446     {
447         snew((*intfpos)[0][n], bins);
448         snew((*intfpos)[1][n], bins);
449         for (i = 0; i < *nslicex; i++)
450         {
451             for (j = 0; j < *nslicey; j++)
452             {
453                 rangeArray(perm, nslicez); /*reset permutation array to identity*/
454                 /*Binsearch returns 2 bin-numbers where the order param is  <= setpoint sgintf*/
455                 ndx1 = start_binsearch(sg_4d[n][i][j], perm, 0, nslicez/2-1, sgintf, 1);
456                 ndx2 = start_binsearch(sg_4d[n][i][j], perm, nslicez/2, nslicez-1, sgintf, -1);
457                 /*Use linear interpolation to smooth out the interface position*/
458
459                 /*left interface (0)*/
460                 /*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){
461                    pos=( (sgintf-sg_4d[n][i][j][perm[ndx1]])*perm[ndx1+1]+(sg_4d[n][i][j][perm[ndx1+1]]-sgintf)*perm[ndx1 ])*/
462                 (*intfpos)[0][n][j+*nslicey*i] = (perm[ndx1]+onehalf)*binw;
463                 /*right interface (1)*/
464                 /*alpha=(sgintf-sg_4d[n][i][j][perm[ndx2]])/(sg_4d[n][i][j][perm[ndx2]+1]-sg_4d[n][i][j][perm[ndx2]]);*/
465                 /*(*intfpos)[1][n][j+*nslicey*i]=((1-alpha)*perm[ndx2]+alpha*(perm[ndx2]+1)+onehalf)*box[ZZ][ZZ]/nslicez;*/
466                 (*intfpos)[1][n][j+*nslicey*i] = (perm[ndx2]+onehalf)*binw;
467             }
468         }
469     }
470
471
472     /*sfree(perm);*/
473     sfree(sk_4d);
474     sfree(sg_4d);
475     /*sfree(sg_grid);*/
476     /*sfree(sk_grid);*/
477
478
479 }
480
481 static void writesurftoxpms(real ***surf, int tblocks, int xbins, int ybins, real bw, char **outfiles, int maplevels )
482 {
483
484     char   numbuf[8];
485     int    n, i, j;
486     real **profile1, **profile2;
487     real   max1, max2, min1, min2, *xticks, *yticks;
488     t_rgb  lo = {1, 1, 1};
489     t_rgb  hi = {0, 0, 0};
490     FILE  *xpmfile1, *xpmfile2;
491
492 /*Prepare xpm structures for output*/
493
494 /*Allocate memory to tick's and matrices*/
495     snew (xticks, xbins+1);
496     snew (yticks, ybins+1);
497
498     profile1 = mk_matrix(xbins, ybins, FALSE);
499     profile2 = mk_matrix(xbins, ybins, FALSE);
500
501     for (i = 0; i < xbins+1; i++)
502     {
503         xticks[i] += bw;
504     }
505     for (j = 0; j < ybins+1; j++)
506     {
507         yticks[j] += bw;
508     }
509
510     xpmfile1 = ffopen(outfiles[0], "w");
511     xpmfile2 = ffopen(outfiles[1], "w");
512
513     max1 = max2 = 0.0;
514     min1 = min2 = 1000.00;
515
516     for (n = 0; n < tblocks; n++)
517     {
518         sprintf(numbuf, "%5d", n);
519         /*Filling matrices for inclusion in xpm-files*/
520         for (i = 0; i < xbins; i++)
521         {
522             for (j = 0; j < ybins; j++)
523             {
524                 profile1[i][j] = (surf[0][n][j+ybins*i]);
525                 profile2[i][j] = (surf[1][n][j+ybins*i]);
526                 /*Finding max and min values*/
527                 if (profile1[i][j] > max1)
528                 {
529                     max1 = profile1[i][j];
530                 }
531                 if (profile1[i][j] < min1)
532                 {
533                     min1 = profile1[i][j];
534                 }
535                 if (profile2[i][j] > max2)
536                 {
537                     max2 = profile2[i][j];
538                 }
539                 if (profile2[i][j] < min2)
540                 {
541                     min2 = profile2[i][j];
542                 }
543             }
544         }
545
546         write_xpm(xpmfile1, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile1, min1, max1, lo, hi, &maplevels);
547         write_xpm(xpmfile2, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile2, min2, max2, lo, hi, &maplevels);
548     }
549
550     ffclose(xpmfile1);
551     ffclose(xpmfile2);
552
553
554
555     sfree(profile1);
556     sfree(profile2);
557     sfree(xticks);
558     sfree(yticks);
559 }
560
561
562 static void writeraw(real ***surf, int tblocks, int xbins, int ybins, char **fnms)
563 {
564     FILE *raw1, *raw2;
565     int   i, j, n;
566
567     raw1 = ffopen(fnms[0], "w");
568     raw2 = ffopen(fnms[1], "w");
569     fprintf(raw1, "#Legend\n#TBlock\n#Xbin Ybin Z t\n");
570     fprintf(raw2, "#Legend\n#TBlock\n#Xbin Ybin Z t\n");
571     for (n = 0; n < tblocks; n++)
572     {
573         fprintf(raw1, "%5d\n", n);
574         fprintf(raw2, "%5d\n", n);
575         for (i = 0; i < xbins; i++)
576         {
577             for (j = 0; j < ybins; j++)
578             {
579                 fprintf(raw1, "%i  %i  %8.5f\n", i, j, (surf[0][n][j+ybins*i]));
580                 fprintf(raw2, "%i  %i  %8.5f\n", i, j, (surf[1][n][j+ybins*i]));
581             }
582         }
583     }
584
585     ffclose(raw1);
586     ffclose(raw2);
587 }
588
589
590
591 int gmx_hydorder(int argc, char *argv[])
592 {
593     static const char *desc[] = {
594         "g_hydorder computes the tetrahedrality order parameters around a ",
595         "given atom. Both angle an distance order parameters are calculated. See",
596         "P.-L. Chau and A.J. Hardwick, Mol. Phys., 93, (1998), 511-518.",
597         "for more details.[BR]"
598         "This application calculates the orderparameter in a 3d-mesh in the box, and",
599         "with 2 phases in the box gives the user the option to define a 2D interface in time",
600         "separating the faces by specifying parameters -sgang1 and -sgang2 (It is important",
601         "to select these judiciously)"
602     };
603
604     int                axis      = 0;
605     static int         nsttblock = 1;
606     static int         nlevels   = 100;
607     static real        binwidth  = 1.0; /* binwidth in mesh           */
608     static real        sg1       = 1;
609     static real        sg2       = 1;   /* order parameters for bulk phases */
610     static gmx_bool    bFourier  = FALSE;
611     static gmx_bool    bRawOut   = FALSE;
612     int                frames, xslices, yslices; /* Dimensions of interface arrays*/
613     real            ***intfpos;                  /* Interface arrays (intfnr,t,xy) -potentially large */
614     static char       *normal_axis[] = { NULL, "z", "x", "y", NULL };
615
616     t_pargs            pa[] = {
617         { "-d",   FALSE, etENUM, {normal_axis},
618           "Direction of the normal on the membrane" },
619         { "-bw",  FALSE, etREAL, {&binwidth},
620           "Binwidth of box mesh" },
621         { "-sgang1", FALSE, etREAL, {&sg1},
622           "tetrahedral angle parameter in Phase 1 (bulk)" },
623         { "-sgang2", FALSE, etREAL, {&sg2},
624           "tetrahedral angle parameter in Phase 2 (bulk)" },
625         { "-tblock", FALSE, etINT, {&nsttblock},
626           "Number of frames in one time-block average"},
627         { "-nlevel", FALSE, etINT, {&nlevels},
628           "Number of Height levels in 2D - XPixMaps"}
629     };
630
631     t_filenm           fnm[] = {                      /* files for g_order    */
632         { efTRX, "-f", NULL,  ffREAD },               /* trajectory file              */
633         { efNDX, "-n", NULL,  ffREAD },               /* index file           */
634         { efTPX, "-s", NULL,  ffREAD },               /* topology file                */
635         { efXPM, "-o", "intf",  ffWRMULT},            /* XPM- surface maps      */
636         { efOUT, "-or", "raw", ffOPTWRMULT },         /* xvgr output file           */
637         { efOUT, "-Spect", "intfspect", ffOPTWRMULT}, /* Fourier spectrum interfaces */
638     };
639 #define NFILE asize(fnm)
640
641     /*Filenames*/
642     const char  *ndxfnm, *tpsfnm, *trxfnm;
643     char       **spectra, **intfn, **raw;
644     int          nfspect, nfxpm, nfraw;
645     output_env_t oenv;
646
647     parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
648                       NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
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(efTPX, 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
723
724     thanx(stderr);
725
726     return 0;
727 }