Merge branch 'release-4-6'
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_sorient.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include "macros.h"
40 #include "statutil.h"
41 #include "smalloc.h"
42 #include "gstat.h"
43 #include "vec.h"
44 #include "xvgr.h"
45 #include "pbc.h"
46 #include "index.h"
47 #include "tpxio.h"
48 #include "gmx_ana.h"
49
50
51 static void calc_com_pbc(int nrefat, t_topology *top, rvec x[], t_pbc *pbc,
52                          atom_id index[], rvec xref, gmx_bool bPBC)
53 {
54     const real tol = 1e-4;
55     gmx_bool   bChanged;
56     int        m, j, ai, iter;
57     real       mass, mtot;
58     rvec       dx, xtest;
59
60     /* First simple calculation */
61     clear_rvec(xref);
62     mtot = 0;
63     for (m = 0; (m < nrefat); m++)
64     {
65         ai   = index[m];
66         mass = top->atoms.atom[ai].m;
67         for (j = 0; (j < DIM); j++)
68         {
69             xref[j] += mass*x[ai][j];
70         }
71         mtot += mass;
72     }
73     svmul(1/mtot, xref, xref);
74     /* Now check if any atom is more than half the box from the COM */
75     if (bPBC)
76     {
77         iter = 0;
78         do
79         {
80             bChanged = FALSE;
81             for (m = 0; (m < nrefat); m++)
82             {
83                 ai   = index[m];
84                 mass = top->atoms.atom[ai].m/mtot;
85                 pbc_dx(pbc, x[ai], xref, dx);
86                 rvec_add(xref, dx, xtest);
87                 for (j = 0; (j < DIM); j++)
88                 {
89                     if (fabs(xtest[j]-x[ai][j]) > tol)
90                     {
91                         /* Here we have used the wrong image for contributing to the COM */
92                         xref[j] += mass*(xtest[j]-x[ai][j]);
93                         x[ai][j] = xtest[j];
94                         bChanged = TRUE;
95                     }
96                 }
97             }
98             if (bChanged)
99             {
100                 printf("COM: %8.3f  %8.3f  %8.3f  iter = %d\n", xref[XX], xref[YY], xref[ZZ], iter);
101             }
102             iter++;
103         }
104         while (bChanged);
105     }
106 }
107
108 int gmx_sorient(int argc, char *argv[])
109 {
110     t_topology      top;
111     int             ePBC = -1;
112     char            title[STRLEN];
113     t_trxstatus    *status;
114     int             natoms;
115     real            t;
116     rvec           *xtop, *x;
117     matrix          box;
118
119     FILE           *fp;
120     int             i, j, p, sa0, sa1, sa2, n, ntot, nf, m, *hist1, *hist2, *histn, nbin1, nbin2, nrbin;
121     real           *histi1, *histi2, invbw, invrbw;
122     double          sum1, sum2;
123     int            *isize, nrefgrp, nrefat;
124     atom_id       **index;
125     char          **grpname;
126     real            inp, outp, two_pi, nav, normfac, rmin2, rmax2, rcut, rcut2, r2, r, mass, mtot;
127     real            c1, c2;
128     char            str[STRLEN];
129     gmx_bool        bTPS;
130     rvec            xref, dx, dxh1, dxh2, outer;
131     gmx_rmpbc_t     gpbc = NULL;
132     t_pbc           pbc;
133     const char     *legr[] = {
134         "<cos(\\8q\\4\\s1\\N)>",
135         "<3cos\\S2\\N(\\8q\\4\\s2\\N)-1>"
136     };
137     const char     *legc[] = {
138         "cos(\\8q\\4\\s1\\N)",
139         "3cos\\S2\\N(\\8q\\4\\s2\\N)-1"
140     };
141
142     const char     *desc[] = {
143         "[TT]g_sorient[tt] analyzes solvent orientation around solutes.",
144         "It calculates two angles between the vector from one or more",
145         "reference positions to the first atom of each solvent molecule:[PAR]",
146         "[GRK]theta[grk][SUB]1[sub]: the angle with the vector from the first atom of the solvent",
147         "molecule to the midpoint between atoms 2 and 3.[BR]",
148         "[GRK]theta[grk][SUB]2[sub]: the angle with the normal of the solvent plane, defined by the",
149         "same three atoms, or, when the option [TT]-v23[tt] is set, ",
150         "the angle with the vector between atoms 2 and 3.[PAR]",
151         "The reference can be a set of atoms or",
152         "the center of mass of a set of atoms. The group of solvent atoms should",
153         "consist of 3 atoms per solvent molecule.",
154         "Only solvent molecules between [TT]-rmin[tt] and [TT]-rmax[tt] are",
155         "considered for [TT]-o[tt] and [TT]-no[tt] each frame.[PAR]",
156         "[TT]-o[tt]: distribtion of [MATH][COS][GRK]theta[grk][SUB]1[sub][cos][math] for rmin<=r<=rmax.[PAR]",
157         "[TT]-no[tt]: distribution of [MATH][COS][GRK]theta[grk][SUB]2[sub][cos][math] for rmin<=r<=rmax.[PAR]",
158         "[TT]-ro[tt]: [MATH][CHEVRON][COS][GRK]theta[grk][SUB]1[sub][cos][chevron][math] and [MATH][CHEVRON]3[COS]^2[GRK]theta[grk][SUB]2[sub][cos]-1[chevron][math] as a function of the",
159         "distance.[PAR]",
160         "[TT]-co[tt]: the sum over all solvent molecules within distance r",
161         "of [MATH][COS][GRK]theta[grk][SUB]1[sub][cos][math] and [MATH]3[COS]^2([GRK]theta[grk][SUB]2[sub])-1[cos][math] as a function of r.[PAR]",
162         "[TT]-rc[tt]: the distribution of the solvent molecules as a function of r"
163     };
164
165     output_env_t    oenv;
166     static gmx_bool bCom = FALSE, bVec23 = FALSE, bPBC = FALSE;
167     static real     rmin = 0.0, rmax = 0.5, binwidth = 0.02, rbinw = 0.02;
168     t_pargs         pa[] = {
169         { "-com",  FALSE, etBOOL,  {&bCom},
170           "Use the center of mass as the reference postion" },
171         { "-v23",  FALSE, etBOOL,  {&bVec23},
172           "Use the vector between atoms 2 and 3" },
173         { "-rmin",  FALSE, etREAL, {&rmin}, "Minimum distance (nm)" },
174         { "-rmax",  FALSE, etREAL, {&rmax}, "Maximum distance (nm)" },
175         { "-cbin",  FALSE, etREAL, {&binwidth}, "Binwidth for the cosine" },
176         { "-rbin",  FALSE, etREAL, {&rbinw}, "Binwidth for r (nm)" },
177         { "-pbc",   FALSE, etBOOL, {&bPBC}, "Check PBC for the center of mass calculation. Only necessary when your reference group consists of several molecules." }
178     };
179
180     t_filenm        fnm[] = {
181         { efTRX, NULL,  NULL,  ffREAD },
182         { efTPS, NULL,  NULL,  ffREAD },
183         { efNDX, NULL,  NULL,  ffOPTRD },
184         { efXVG, NULL,  "sori.xvg",  ffWRITE },
185         { efXVG, "-no", "snor.xvg",  ffWRITE },
186         { efXVG, "-ro", "sord.xvg",  ffWRITE },
187         { efXVG, "-co", "scum.xvg",  ffWRITE },
188         { efXVG, "-rc", "scount.xvg",  ffWRITE }
189     };
190 #define NFILE asize(fnm)
191
192     if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
193                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
194     {
195         return 0;
196     }
197
198     two_pi = 2/M_PI;
199
200     bTPS = (opt2bSet("-s", NFILE, fnm) || !opt2bSet("-n", NFILE, fnm) || bCom);
201     if (bTPS)
202     {
203         read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, box,
204                       bCom);
205     }
206
207     /* get index groups */
208     printf("Select a group of reference particles and a solvent group:\n");
209     snew(grpname, 2);
210     snew(index, 2);
211     snew(isize, 2);
212     if (bTPS)
213     {
214         get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 2, isize, index, grpname);
215     }
216     else
217     {
218         get_index(NULL, ftp2fn(efNDX, NFILE, fnm), 2, isize, index, grpname);
219     }
220
221     if (bCom)
222     {
223         nrefgrp = 1;
224         nrefat  = isize[0];
225     }
226     else
227     {
228         nrefgrp = isize[0];
229         nrefat  = 1;
230     }
231
232     if (isize[1] % 3)
233     {
234         gmx_fatal(FARGS, "The number of solvent atoms (%d) is not a multiple of 3",
235                   isize[1]);
236     }
237
238     /* initialize reading trajectory:                         */
239     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
240
241     rmin2 = sqr(rmin);
242     rmax2 = sqr(rmax);
243     rcut  = 0.99*sqrt(max_cutoff2(guess_ePBC(box), box));
244     if (rcut == 0)
245     {
246         rcut = 10*rmax;
247     }
248     rcut2 = sqr(rcut);
249
250     invbw = 1/binwidth;
251     nbin1 = 1+(int)(2*invbw + 0.5);
252     nbin2 = 1+(int)(invbw + 0.5);
253
254     invrbw = 1/rbinw;
255
256     snew(hist1, nbin1);
257     snew(hist2, nbin2);
258     nrbin = 1+(int)(rcut/rbinw);
259     if (nrbin == 0)
260     {
261         nrbin = 1;
262     }
263     snew(histi1, nrbin);
264     snew(histi2, nrbin);
265     snew(histn, nrbin);
266
267     ntot = 0;
268     nf   = 0;
269     sum1 = 0;
270     sum2 = 0;
271
272     if (bTPS)
273     {
274         /* make molecules whole again */
275         gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
276     }
277     /* start analysis of trajectory */
278     do
279     {
280         if (bTPS)
281         {
282             /* make molecules whole again */
283             gmx_rmpbc(gpbc, natoms, box, x);
284         }
285
286         set_pbc(&pbc, ePBC, box);
287         n    = 0;
288         inp  = 0;
289         outp = 0;
290         for (p = 0; (p < nrefgrp); p++)
291         {
292             if (bCom)
293             {
294                 calc_com_pbc(nrefat, &top, x, &pbc, index[0], xref, bPBC);
295             }
296             else
297             {
298                 copy_rvec(x[index[0][p]], xref);
299             }
300
301             for (m = 0; m < isize[1]; m += 3)
302             {
303                 sa0 = index[1][m];
304                 sa1 = index[1][m+1];
305                 sa2 = index[1][m+2];
306                 range_check(sa0, 0, natoms);
307                 range_check(sa1, 0, natoms);
308                 range_check(sa2, 0, natoms);
309                 pbc_dx(&pbc, x[sa0], xref, dx);
310                 r2  = norm2(dx);
311                 if (r2 < rcut2)
312                 {
313                     r = sqrt(r2);
314                     if (!bVec23)
315                     {
316                         /* Determine the normal to the plain */
317                         rvec_sub(x[sa1], x[sa0], dxh1);
318                         rvec_sub(x[sa2], x[sa0], dxh2);
319                         rvec_inc(dxh1, dxh2);
320                         svmul(1/r, dx, dx);
321                         unitv(dxh1, dxh1);
322                         inp = iprod(dx, dxh1);
323                         cprod(dxh1, dxh2, outer);
324                         unitv(outer, outer);
325                         outp = iprod(dx, outer);
326                     }
327                     else
328                     {
329                         /* Use the vector between the 2nd and 3rd atom */
330                         rvec_sub(x[sa2], x[sa1], dxh2);
331                         unitv(dxh2, dxh2);
332                         outp = iprod(dx, dxh2)/r;
333                     }
334                     {
335                         int ii = (int)(invrbw*r);
336                         range_check(ii, 0, nrbin);
337                         histi1[ii] += inp;
338                         histi2[ii] += 3*sqr(outp) - 1;
339                         histn[ii]++;
340                     }
341                     if ((r2 >= rmin2) && (r2 < rmax2))
342                     {
343                         int ii1 = (int)(invbw*(inp + 1));
344                         int ii2 = (int)(invbw*fabs(outp));
345
346                         range_check(ii1, 0, nbin1);
347                         range_check(ii2, 0, nbin2);
348                         hist1[ii1]++;
349                         hist2[ii2]++;
350                         sum1 += inp;
351                         sum2 += outp;
352                         n++;
353                     }
354                 }
355             }
356         }
357         ntot += n;
358         nf++;
359
360     }
361     while (read_next_x(oenv, status, &t, x, box));
362
363     /* clean up */
364     sfree(x);
365     close_trj(status);
366     gmx_rmpbc_done(gpbc);
367
368     /* Add the bin for the exact maximum to the previous bin */
369     hist1[nbin1-1] += hist1[nbin1];
370     hist2[nbin2-1] += hist2[nbin2];
371
372     nav     = (real)ntot/(nrefgrp*nf);
373     normfac = invbw/ntot;
374
375     fprintf(stderr,  "Average nr of molecules between %g and %g nm: %.1f\n",
376             rmin, rmax, nav);
377     if (ntot > 0)
378     {
379         sum1 /= ntot;
380         sum2 /= ntot;
381         fprintf(stderr, "Average cos(theta1)     between %g and %g nm: %6.3f\n",
382                 rmin, rmax, sum1);
383         fprintf(stderr, "Average 3cos2(theta2)-1 between %g and %g nm: %6.3f\n",
384                 rmin, rmax, sum2);
385     }
386
387     sprintf(str, "Solvent orientation between %g and %g nm", rmin, rmax);
388     fp = xvgropen(opt2fn("-o", NFILE, fnm), str, "cos(\\8q\\4\\s1\\N)", "", oenv);
389     if (output_env_get_print_xvgr_codes(oenv))
390     {
391         fprintf(fp, "@ subtitle \"average shell size %.1f molecules\"\n", nav);
392     }
393     for (i = 0; i < nbin1; i++)
394     {
395         fprintf(fp, "%g %g\n", (i+0.5)*binwidth-1, 2*normfac*hist1[i]);
396     }
397     ffclose(fp);
398
399     sprintf(str, "Solvent normal orientation between %g and %g nm", rmin, rmax);
400     fp = xvgropen(opt2fn("-no", NFILE, fnm), str, "cos(\\8q\\4\\s2\\N)", "", oenv);
401     if (output_env_get_print_xvgr_codes(oenv))
402     {
403         fprintf(fp, "@ subtitle \"average shell size %.1f molecules\"\n", nav);
404     }
405     for (i = 0; i < nbin2; i++)
406     {
407         fprintf(fp, "%g %g\n", (i+0.5)*binwidth, normfac*hist2[i]);
408     }
409     ffclose(fp);
410
411
412     sprintf(str, "Solvent orientation");
413     fp = xvgropen(opt2fn("-ro", NFILE, fnm), str, "r (nm)", "", oenv);
414     if (output_env_get_print_xvgr_codes(oenv))
415     {
416         fprintf(fp, "@ subtitle \"as a function of distance\"\n");
417     }
418     xvgr_legend(fp, 2, legr, oenv);
419     for (i = 0; i < nrbin; i++)
420     {
421         fprintf(fp, "%g %g %g\n", (i+0.5)*rbinw,
422                 histn[i] ? histi1[i]/histn[i] : 0,
423                 histn[i] ? histi2[i]/histn[i] : 0);
424     }
425     ffclose(fp);
426
427     sprintf(str, "Cumulative solvent orientation");
428     fp = xvgropen(opt2fn("-co", NFILE, fnm), str, "r (nm)", "", oenv);
429     if (output_env_get_print_xvgr_codes(oenv))
430     {
431         fprintf(fp, "@ subtitle \"as a function of distance\"\n");
432     }
433     xvgr_legend(fp, 2, legc, oenv);
434     normfac = 1.0/(nrefgrp*nf);
435     c1      = 0;
436     c2      = 0;
437     fprintf(fp, "%g %g %g\n", 0.0, c1, c2);
438     for (i = 0; i < nrbin; i++)
439     {
440         c1 += histi1[i]*normfac;
441         c2 += histi2[i]*normfac;
442         fprintf(fp, "%g %g %g\n", (i+1)*rbinw, c1, c2);
443     }
444     ffclose(fp);
445
446     sprintf(str, "Solvent distribution");
447     fp = xvgropen(opt2fn("-rc", NFILE, fnm), str, "r (nm)", "molecules/nm", oenv);
448     if (output_env_get_print_xvgr_codes(oenv))
449     {
450         fprintf(fp, "@ subtitle \"as a function of distance\"\n");
451     }
452     normfac = 1.0/(rbinw*nf);
453     for (i = 0; i < nrbin; i++)
454     {
455         fprintf(fp, "%g %g\n", (i+0.5)*rbinw, histn[i]*normfac);
456     }
457     ffclose(fp);
458
459     do_view(oenv, opt2fn("-o", NFILE, fnm), NULL);
460     do_view(oenv, opt2fn("-no", NFILE, fnm), NULL);
461     do_view(oenv, opt2fn("-ro", NFILE, fnm), "-nxy");
462     do_view(oenv, opt2fn("-co", NFILE, fnm), "-nxy");
463
464     return 0;
465 }