Moved copyrite.* to fileio from gmxlib and legacyheaders.
[gromacs.git] / src / gromacs / gmxana / gmx_hydorder.cpp
blob0038e70d0eff7823bd61c771eb45c59c27f33eb2
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010,2011,2012,2013,2014,2015, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
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.
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.
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.
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.
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.
36 #include "gmxpre.h"
38 #include <cmath>
39 #include <cstring>
41 #include "gromacs/commandline/pargs.h"
42 #include "gromacs/fileio/confio.h"
43 #include "gromacs/fileio/matio.h"
44 #include "gromacs/fileio/trxio.h"
45 #include "gromacs/fileio/xvgr.h"
46 #include "gromacs/gmxana/binsearch.h"
47 #include "gromacs/gmxana/gmx_ana.h"
48 #include "gromacs/gmxana/gstat.h"
49 #include "gromacs/gmxana/powerspect.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/pbcutil/pbc.h"
52 #include "gromacs/pbcutil/rmpbc.h"
53 #include "gromacs/topology/index.h"
54 #include "gromacs/topology/topology.h"
55 #include "gromacs/utility/arraysize.h"
56 #include "gromacs/utility/cstringutil.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/futil.h"
59 #include "gromacs/utility/gmxassert.h"
60 #include "gromacs/utility/smalloc.h"
62 /* Print name of first atom in all groups in index file */
63 static void print_types(atom_id index[], atom_id a[], int ngrps,
64 char *groups[], t_topology *top)
66 int i;
68 fprintf(stderr, "Using following groups: \n");
69 for (i = 0; i < ngrps; i++)
71 fprintf(stderr, "Groupname: %s First atomname: %s First atomnr %d\n",
72 groups[i], *(top->atoms.atomname[a[index[i]]]), a[index[i]]);
74 fprintf(stderr, "\n");
77 static void check_length(real length, int a, int b)
79 if (length > 0.3)
81 fprintf(stderr, "WARNING: distance between atoms %d and "
82 "%d > 0.3 nm (%f). Index file might be corrupt.\n",
83 a, b, length);
87 static void find_tetra_order_grid(t_topology top, int ePBC,
88 int natoms, matrix box,
89 rvec x[], int maxidx, atom_id index[],
90 real *sgmean, real *skmean,
91 int nslicex, int nslicey, int nslicez,
92 real ***sggrid, real ***skgrid)
94 int ix, jx, i, j, k, l, n, *nn[4];
95 rvec dx, rj, rk, urk, urj;
96 real cost, cost2, *sgmol, *skmol, rmean, rmean2, r2, box2, *r_nn[4];
97 t_pbc pbc;
98 int slindex_x, slindex_y, slindex_z;
99 int ***sl_count;
100 real onethird = 1.0/3.0;
101 gmx_rmpbc_t gpbc;
103 /* dmat = init_mat(maxidx, FALSE); */
105 box2 = box[XX][XX] * box[XX][XX];
107 /* Initialize expanded sl_count array */
108 snew(sl_count, nslicex);
109 for (i = 0; i < nslicex; i++)
111 snew(sl_count[i], nslicey);
112 for (j = 0; j < nslicey; j++)
114 snew(sl_count[i][j], nslicez);
119 for (i = 0; (i < 4); i++)
121 snew(r_nn[i], natoms);
122 snew(nn[i], natoms);
124 for (j = 0; (j < natoms); j++)
126 r_nn[i][j] = box2;
130 snew(sgmol, maxidx);
131 snew(skmol, maxidx);
133 /* Must init pbc every step because of pressure coupling */
134 set_pbc(&pbc, ePBC, box);
135 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
136 gmx_rmpbc(gpbc, natoms, box, x);
138 *sgmean = 0.0;
139 *skmean = 0.0;
140 l = 0;
141 for (i = 0; (i < maxidx); i++)
142 { /* loop over index file */
143 ix = index[i];
144 for (j = 0; (j < maxidx); j++)
147 if (i == j)
149 continue;
152 jx = index[j];
154 pbc_dx(&pbc, x[ix], x[jx], dx);
155 r2 = iprod(dx, dx);
157 /* set_mat_entry(dmat,i,j,r2); */
159 /* determine the nearest neighbours */
160 if (r2 < r_nn[0][i])
162 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
163 r_nn[2][i] = r_nn[1][i]; nn[2][i] = nn[1][i];
164 r_nn[1][i] = r_nn[0][i]; nn[1][i] = nn[0][i];
165 r_nn[0][i] = r2; nn[0][i] = j;
167 else if (r2 < r_nn[1][i])
169 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
170 r_nn[2][i] = r_nn[1][i]; nn[2][i] = nn[1][i];
171 r_nn[1][i] = r2; nn[1][i] = j;
173 else if (r2 < r_nn[2][i])
175 r_nn[3][i] = r_nn[2][i]; nn[3][i] = nn[2][i];
176 r_nn[2][i] = r2; nn[2][i] = j;
178 else if (r2 < r_nn[3][i])
180 r_nn[3][i] = r2; nn[3][i] = j;
185 /* calculate mean distance between nearest neighbours */
186 rmean = 0;
187 for (j = 0; (j < 4); j++)
189 r_nn[j][i] = std::sqrt(r_nn[j][i]);
190 rmean += r_nn[j][i];
192 rmean /= 4;
194 n = 0;
195 sgmol[i] = 0.0;
196 skmol[i] = 0.0;
198 /* Chau1998a eqn 3 */
199 /* angular part tetrahedrality order parameter per atom */
200 for (j = 0; (j < 3); j++)
202 for (k = j+1; (k < 4); k++)
204 pbc_dx(&pbc, x[ix], x[index[nn[k][i]]], rk);
205 pbc_dx(&pbc, x[ix], x[index[nn[j][i]]], rj);
207 unitv(rk, urk);
208 unitv(rj, urj);
210 cost = iprod(urk, urj) + onethird;
211 cost2 = cost * cost;
213 sgmol[i] += cost2;
214 l++;
215 n++;
218 /* normalize sgmol between 0.0 and 1.0 */
219 sgmol[i] = 3*sgmol[i]/32;
220 *sgmean += sgmol[i];
222 /* distance part tetrahedrality order parameter per atom */
223 rmean2 = 4 * 3 * rmean * rmean;
224 for (j = 0; (j < 4); j++)
226 skmol[i] += (rmean - r_nn[j][i]) * (rmean - r_nn[j][i]) / rmean2;
227 /* printf("%d %f (%f %f %f %f) \n",
228 i, skmol[i], rmean, rmean2, r_nn[j][i], (rmean - r_nn[j][i]) );
232 *skmean += skmol[i];
234 /* Compute sliced stuff in x y z*/
235 slindex_x = static_cast<int>(std::round((1+x[i][XX]/box[XX][XX])*nslicex)) % nslicex;
236 slindex_y = static_cast<int>(std::round((1+x[i][YY]/box[YY][YY])*nslicey)) % nslicey;
237 slindex_z = static_cast<int>(std::round((1+x[i][ZZ]/box[ZZ][ZZ])*nslicez)) % nslicez;
238 sggrid[slindex_x][slindex_y][slindex_z] += sgmol[i];
239 skgrid[slindex_x][slindex_y][slindex_z] += skmol[i];
240 (sl_count[slindex_x][slindex_y][slindex_z])++;
241 } /* loop over entries in index file */
243 *sgmean /= maxidx;
244 *skmean /= maxidx;
246 for (i = 0; (i < nslicex); i++)
248 for (j = 0; j < nslicey; j++)
250 for (k = 0; k < nslicez; k++)
252 if (sl_count[i][j][k] > 0)
254 sggrid[i][j][k] /= sl_count[i][j][k];
255 skgrid[i][j][k] /= sl_count[i][j][k];
261 sfree(sl_count);
262 sfree(sgmol);
263 sfree(skmol);
264 for (i = 0; (i < 4); i++)
266 sfree(r_nn[i]);
267 sfree(nn[i]);
271 /*Determines interface from tetrahedral order parameter in box with specified binwidth. */
272 /*Outputs interface positions(bins), the number of timeframes, and the number of surface-mesh points in xy*/
274 static void calc_tetra_order_interface(const char *fnNDX, const char *fnTPS, const char *fnTRX, real binw, int tblock,
275 int *nframes, int *nslicex, int *nslicey,
276 real sgang1, real sgang2, real ****intfpos,
277 gmx_output_env_t *oenv)
279 FILE *fpsg = NULL, *fpsk = NULL;
280 t_topology top;
281 int ePBC;
282 t_trxstatus *status;
283 int natoms;
284 real t;
285 rvec *xtop, *x;
286 matrix box;
287 real sg, sk, sgintf;
288 atom_id **index = NULL;
289 char **grpname = NULL;
290 int i, j, k, n, *isize, ng, nslicez, framenr;
291 real ***sg_grid = NULL, ***sk_grid = NULL, ***sg_fravg = NULL, ***sk_fravg = NULL, ****sk_4d = NULL, ****sg_4d = NULL;
292 int *perm;
293 int ndx1, ndx2;
294 int bins;
295 const real onehalf = 1.0/2.0;
296 /* real ***intfpos[2]; pointers to arrays of two interface positions zcoord(framenr,xbin,ybin): intfpos[interface_index][t][nslicey*x+y]
297 * i.e 1D Row-major order in (t,x,y) */
300 read_tps_conf(fnTPS, &top, &ePBC, &xtop, NULL, box, FALSE);
302 *nslicex = static_cast<int>(box[XX][XX]/binw + onehalf); /*Calculate slicenr from binwidth*/
303 *nslicey = static_cast<int>(box[YY][YY]/binw + onehalf);
304 nslicez = static_cast<int>(box[ZZ][ZZ]/binw + onehalf);
308 ng = 1;
309 /* get index groups */
310 printf("Select the group that contains the atoms you want to use for the tetrahedrality order parameter calculation:\n");
311 snew(grpname, ng);
312 snew(index, ng);
313 snew(isize, ng);
314 get_index(&top.atoms, fnNDX, ng, isize, index, grpname);
316 /* Analyze trajectory */
317 natoms = read_first_x(oenv, &status, fnTRX, &t, &x, box);
318 if (natoms > top.atoms.nr)
320 gmx_fatal(FARGS, "Topology (%d atoms) does not match trajectory (%d atoms)",
321 top.atoms.nr, natoms);
323 check_index(NULL, ng, index[0], NULL, natoms);
326 /*Prepare structures for temporary storage of frame info*/
327 snew(sg_grid, *nslicex);
328 snew(sk_grid, *nslicex);
329 for (i = 0; i < *nslicex; i++)
331 snew(sg_grid[i], *nslicey);
332 snew(sk_grid[i], *nslicey);
333 for (j = 0; j < *nslicey; j++)
335 snew(sg_grid[i][j], nslicez);
336 snew(sk_grid[i][j], nslicez);
340 sg_4d = NULL;
341 sk_4d = NULL;
342 *nframes = 0;
343 framenr = 0;
345 /* Loop over frames*/
348 /*Initialize box meshes (temporary storage for each tblock frame -reinitialise every tblock steps */
349 if (framenr%tblock == 0)
351 srenew(sk_4d, *nframes+1);
352 srenew(sg_4d, *nframes+1);
353 snew(sg_fravg, *nslicex);
354 snew(sk_fravg, *nslicex);
355 for (i = 0; i < *nslicex; i++)
357 snew(sg_fravg[i], *nslicey);
358 snew(sk_fravg[i], *nslicey);
359 for (j = 0; j < *nslicey; j++)
361 snew(sg_fravg[i][j], nslicez);
362 snew(sk_fravg[i][j], nslicez);
367 find_tetra_order_grid(top, ePBC, natoms, box, x, isize[0], index[0],
368 &sg, &sk, *nslicex, *nslicey, nslicez, sg_grid, sk_grid);
369 GMX_RELEASE_ASSERT(sk_fravg != NULL, "Trying to dereference NULL sk_fravg pointer");
370 for (i = 0; i < *nslicex; i++)
372 for (j = 0; j < *nslicey; j++)
374 for (k = 0; k < nslicez; k++)
376 sk_fravg[i][j][k] += sk_grid[i][j][k]/tblock;
377 sg_fravg[i][j][k] += sg_grid[i][j][k]/tblock;
382 framenr++;
384 if (framenr%tblock == 0)
386 GMX_RELEASE_ASSERT(sk_4d != NULL, "Trying to dereference NULL sk_4d pointer");
387 sk_4d[*nframes] = sk_fravg;
388 sg_4d[*nframes] = sg_fravg;
389 (*nframes)++;
393 while (read_next_x(oenv, status, &t, x, box));
394 close_trj(status);
396 sfree(grpname);
397 sfree(index);
398 sfree(isize);
400 /*Debugging for printing out the entire order parameter meshes.*/
401 if (debug)
403 fpsg = xvgropen("sg_ang_mesh", "S\\sg\\N Angle Order Parameter / Meshpoint", "(nm)", "S\\sg\\N", oenv);
404 fpsk = xvgropen("sk_dist_mesh", "S\\sk\\N Distance Order Parameter / Meshpoint", "(nm)", "S\\sk\\N", oenv);
405 for (n = 0; n < (*nframes); n++)
407 fprintf(fpsg, "%i\n", n);
408 fprintf(fpsk, "%i\n", n);
409 for (i = 0; (i < *nslicex); i++)
411 for (j = 0; j < *nslicey; j++)
413 for (k = 0; k < nslicez; k++)
415 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]);
416 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]);
421 xvgrclose(fpsg);
422 xvgrclose(fpsk);
426 /* Find positions of interface z by scanning orderparam for each frame and for each xy-mesh cylinder along z*/
428 /*Simple trial: assume interface is in the middle of -sgang1 and sgang2*/
429 sgintf = 0.5*(sgang1+sgang2);
432 /*Allocate memory for interface arrays; */
433 snew((*intfpos), 2);
434 snew((*intfpos)[0], *nframes);
435 snew((*intfpos)[1], *nframes);
437 bins = (*nslicex)*(*nslicey);
440 snew(perm, nslicez); /*permutation array for sorting along normal coordinate*/
443 for (n = 0; n < *nframes; n++)
445 snew((*intfpos)[0][n], bins);
446 snew((*intfpos)[1][n], bins);
447 for (i = 0; i < *nslicex; i++)
449 for (j = 0; j < *nslicey; j++)
451 rangeArray(perm, nslicez); /*reset permutation array to identity*/
452 /*Binsearch returns 2 bin-numbers where the order param is <= setpoint sgintf*/
453 ndx1 = start_binsearch(sg_4d[n][i][j], perm, 0, nslicez/2-1, sgintf, 1);
454 ndx2 = start_binsearch(sg_4d[n][i][j], perm, nslicez/2, nslicez-1, sgintf, -1);
455 /*Use linear interpolation to smooth out the interface position*/
457 /*left interface (0)*/
458 /*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){
459 pos=( (sgintf-sg_4d[n][i][j][perm[ndx1]])*perm[ndx1+1]+(sg_4d[n][i][j][perm[ndx1+1]]-sgintf)*perm[ndx1 ])*/
460 (*intfpos)[0][n][j+*nslicey*i] = (perm[ndx1]+onehalf)*binw;
461 /*right interface (1)*/
462 /*alpha=(sgintf-sg_4d[n][i][j][perm[ndx2]])/(sg_4d[n][i][j][perm[ndx2]+1]-sg_4d[n][i][j][perm[ndx2]]);*/
463 /*(*intfpos)[1][n][j+*nslicey*i]=((1-alpha)*perm[ndx2]+alpha*(perm[ndx2]+1)+onehalf)*box[ZZ][ZZ]/nslicez;*/
464 (*intfpos)[1][n][j+*nslicey*i] = (perm[ndx2]+onehalf)*binw;
470 sfree(sk_4d);
471 sfree(sg_4d);
474 static void writesurftoxpms(real ***surf, int tblocks, int xbins, int ybins, real bw, char **outfiles, int maplevels )
477 char numbuf[8];
478 int n, i, j;
479 real **profile1, **profile2;
480 real max1, max2, min1, min2, *xticks, *yticks;
481 t_rgb lo = {1, 1, 1};
482 t_rgb hi = {0, 0, 0};
483 FILE *xpmfile1, *xpmfile2;
485 /*Prepare xpm structures for output*/
487 /*Allocate memory to tick's and matrices*/
488 snew (xticks, xbins+1);
489 snew (yticks, ybins+1);
491 profile1 = mk_matrix(xbins, ybins, FALSE);
492 profile2 = mk_matrix(xbins, ybins, FALSE);
494 for (i = 0; i < xbins+1; i++)
496 xticks[i] += bw;
498 for (j = 0; j < ybins+1; j++)
500 yticks[j] += bw;
503 xpmfile1 = gmx_ffopen(outfiles[0], "w");
504 xpmfile2 = gmx_ffopen(outfiles[1], "w");
506 max1 = max2 = 0.0;
507 min1 = min2 = 1000.00;
509 for (n = 0; n < tblocks; n++)
511 sprintf(numbuf, "%5d", n);
512 /*Filling matrices for inclusion in xpm-files*/
513 for (i = 0; i < xbins; i++)
515 for (j = 0; j < ybins; j++)
517 profile1[i][j] = (surf[0][n][j+ybins*i]);
518 profile2[i][j] = (surf[1][n][j+ybins*i]);
519 /*Finding max and min values*/
520 if (profile1[i][j] > max1)
522 max1 = profile1[i][j];
524 if (profile1[i][j] < min1)
526 min1 = profile1[i][j];
528 if (profile2[i][j] > max2)
530 max2 = profile2[i][j];
532 if (profile2[i][j] < min2)
534 min2 = profile2[i][j];
539 write_xpm(xpmfile1, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile1, min1, max1, lo, hi, &maplevels);
540 write_xpm(xpmfile2, 3, numbuf, "Height", "x[nm]", "y[nm]", xbins, ybins, xticks, yticks, profile2, min2, max2, lo, hi, &maplevels);
543 gmx_ffclose(xpmfile1);
544 gmx_ffclose(xpmfile2);
548 sfree(profile1);
549 sfree(profile2);
550 sfree(xticks);
551 sfree(yticks);
555 static void writeraw(real ***surf, int tblocks, int xbins, int ybins, char **fnms)
557 FILE *raw1, *raw2;
558 int i, j, n;
560 raw1 = gmx_ffopen(fnms[0], "w");
561 raw2 = gmx_ffopen(fnms[1], "w");
562 fprintf(raw1, "#Legend\n#TBlock\n#Xbin Ybin Z t\n");
563 fprintf(raw2, "#Legend\n#TBlock\n#Xbin Ybin Z t\n");
564 for (n = 0; n < tblocks; n++)
566 fprintf(raw1, "%5d\n", n);
567 fprintf(raw2, "%5d\n", n);
568 for (i = 0; i < xbins; i++)
570 for (j = 0; j < ybins; j++)
572 fprintf(raw1, "%i %i %8.5f\n", i, j, (surf[0][n][j+ybins*i]));
573 fprintf(raw2, "%i %i %8.5f\n", i, j, (surf[1][n][j+ybins*i]));
578 gmx_ffclose(raw1);
579 gmx_ffclose(raw2);
584 int gmx_hydorder(int argc, char *argv[])
586 static const char *desc[] = {
587 "[THISMODULE] computes the tetrahedrality order parameters around a ",
588 "given atom. Both angle an distance order parameters are calculated. See",
589 "P.-L. Chau and A.J. Hardwick, Mol. Phys., 93, (1998), 511-518.",
590 "for more details.[PAR]"
591 "[THISMODULE] calculates the order parameter in a 3d-mesh in the box, and",
592 "with 2 phases in the box gives the user the option to define a 2D interface in time",
593 "separating the faces by specifying parameters [TT]-sgang1[tt] and",
594 "[TT]-sgang2[tt] (it is important to select these judiciously)."
597 int axis = 0;
598 static int nsttblock = 1;
599 static int nlevels = 100;
600 static real binwidth = 1.0; /* binwidth in mesh */
601 static real sg1 = 1;
602 static real sg2 = 1; /* order parameters for bulk phases */
603 static gmx_bool bFourier = FALSE;
604 static gmx_bool bRawOut = FALSE;
605 int frames, xslices, yslices; /* Dimensions of interface arrays*/
606 real ***intfpos; /* Interface arrays (intfnr,t,xy) -potentially large */
607 static const char *normal_axis[] = { NULL, "z", "x", "y", NULL };
609 t_pargs pa[] = {
610 { "-d", FALSE, etENUM, {normal_axis},
611 "Direction of the normal on the membrane" },
612 { "-bw", FALSE, etREAL, {&binwidth},
613 "Binwidth of box mesh" },
614 { "-sgang1", FALSE, etREAL, {&sg1},
615 "tetrahedral angle parameter in Phase 1 (bulk)" },
616 { "-sgang2", FALSE, etREAL, {&sg2},
617 "tetrahedral angle parameter in Phase 2 (bulk)" },
618 { "-tblock", FALSE, etINT, {&nsttblock},
619 "Number of frames in one time-block average"},
620 { "-nlevel", FALSE, etINT, {&nlevels},
621 "Number of Height levels in 2D - XPixMaps"}
624 t_filenm fnm[] = { /* files for g_order */
625 { efTRX, "-f", NULL, ffREAD }, /* trajectory file */
626 { efNDX, "-n", NULL, ffREAD }, /* index file */
627 { efTPR, "-s", NULL, ffREAD }, /* topology file */
628 { efXPM, "-o", "intf", ffWRMULT}, /* XPM- surface maps */
629 { efOUT, "-or", "raw", ffOPTWRMULT }, /* xvgr output file */
630 { efOUT, "-Spect", "intfspect", ffOPTWRMULT}, /* Fourier spectrum interfaces */
632 #define NFILE asize(fnm)
634 /*Filenames*/
635 const char *ndxfnm, *tpsfnm, *trxfnm;
636 char **spectra, **intfn, **raw;
637 int nfspect, nfxpm, nfraw;
638 gmx_output_env_t *oenv;
640 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME,
641 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
643 return 0;
645 bFourier = opt2bSet("-Spect", NFILE, fnm);
646 bRawOut = opt2bSet("-or", NFILE, fnm);
648 if (binwidth < 0.0)
650 gmx_fatal(FARGS, "Can not have binwidth < 0");
653 ndxfnm = ftp2fn(efNDX, NFILE, fnm);
654 tpsfnm = ftp2fn(efTPR, NFILE, fnm);
655 trxfnm = ftp2fn(efTRX, NFILE, fnm);
657 /* Calculate axis */
658 GMX_RELEASE_ASSERT(normal_axis[0] != NULL, "Option setting inconsistency; normal_axis[0] is NULL");
659 if (std::strcmp(normal_axis[0], "x") == 0)
661 axis = XX;
663 else if (std::strcmp(normal_axis[0], "y") == 0)
665 axis = YY;
667 else if (std::strcmp(normal_axis[0], "z") == 0)
669 axis = ZZ;
671 else
673 gmx_fatal(FARGS, "Invalid axis, use x, y or z");
676 switch (axis)
678 case 0:
679 fprintf(stderr, "Taking x axis as normal to the membrane\n");
680 break;
681 case 1:
682 fprintf(stderr, "Taking y axis as normal to the membrane\n");
683 break;
684 case 2:
685 fprintf(stderr, "Taking z axis as normal to the membrane\n");
686 break;
689 /* tetraheder order parameter */
690 /* If either of the options is set we compute both */
691 nfxpm = opt2fns(&intfn, "-o", NFILE, fnm);
692 if (nfxpm != 2)
694 gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfxpm);
696 calc_tetra_order_interface(ndxfnm, tpsfnm, trxfnm, binwidth, nsttblock, &frames, &xslices, &yslices, sg1, sg2, &intfpos, oenv);
697 writesurftoxpms(intfpos, frames, xslices, yslices, binwidth, intfn, nlevels);
699 if (bFourier)
701 nfspect = opt2fns(&spectra, "-Spect", NFILE, fnm);
702 if (nfspect != 2)
704 gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfspect);
706 powerspectavg(intfpos, frames, xslices, yslices, spectra);
709 if (bRawOut)
711 nfraw = opt2fns(&raw, "-or", NFILE, fnm);
712 if (nfraw != 2)
714 gmx_fatal(FARGS, "No or not correct number (2) of output-files: %d", nfraw);
716 writeraw(intfpos, frames, xslices, yslices, raw);
719 return 0;