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