85e05de0e8ada5a44a612635ea529a5155b1d401
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_saltbr.cpp
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,2015, 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 <cmath>
40 #include <cstring>
41
42 #include "gromacs/commandline/pargs.h"
43 #include "gromacs/fileio/filenm.h"
44 #include "gromacs/fileio/trxio.h"
45 #include "gromacs/fileio/xvgr.h"
46 #include "gromacs/gmxana/gmx_ana.h"
47 #include "gromacs/legacyheaders/macros.h"
48 #include "gromacs/legacyheaders/typedefs.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/pbcutil/pbc.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/futil.h"
54 #include "gromacs/utility/smalloc.h"
55
56 typedef struct {
57     char *label;
58     int   cg;
59     real  q;
60 } t_charge;
61
62 static t_charge *mk_charge(t_atoms *atoms, t_block *cgs, int *nncg)
63 {
64     t_charge *cg = NULL;
65     char      buf[32];
66     int       i, j, ncg, resnr, anr;
67     real      qq;
68
69     /* Find the charged groups */
70     ncg = 0;
71     for (i = 0; (i < cgs->nr); i++)
72     {
73         qq = 0.0;
74         for (j = cgs->index[i]; (j < cgs->index[i+1]); j++)
75         {
76             qq += atoms->atom[j].q;
77         }
78         if (std::abs(qq) > 1.0e-5)
79         {
80             srenew(cg, ncg+1);
81             cg[ncg].q  = qq;
82             cg[ncg].cg = i;
83             anr        = cgs->index[i];
84             resnr      = atoms->atom[anr].resind;
85             sprintf(buf, "%s%d-%d",
86                     *(atoms->resinfo[resnr].name),
87                     atoms->resinfo[resnr].nr,
88                     anr+1);
89             cg[ncg].label = gmx_strdup(buf);
90             ncg++;
91         }
92     }
93     *nncg = ncg;
94
95     for (i = 0; (i < ncg); i++)
96     {
97         printf("CG: %10s Q: %6g  Atoms:",
98                cg[i].label, cg[i].q);
99         for (j = cgs->index[cg[i].cg]; (j < cgs->index[cg[i].cg+1]); j++)
100         {
101             printf(" %4d", j);
102         }
103         printf("\n");
104     }
105
106     return cg;
107 }
108
109 static real calc_dist(t_pbc *pbc, rvec x[], t_block *cgs, int icg, int jcg)
110 {
111     int  i, j;
112     rvec dx;
113     real d2, mindist2 = 1000;
114
115     for (i = cgs->index[icg]; (i < cgs->index[icg+1]); i++)
116     {
117         for (j = cgs->index[jcg]; (j < cgs->index[jcg+1]); j++)
118         {
119             pbc_dx(pbc, x[i], x[j], dx);
120             d2 = norm2(dx);
121             if (d2 < mindist2)
122             {
123                 mindist2 = d2;
124             }
125         }
126     }
127     return std::sqrt(mindist2);
128 }
129
130 int gmx_saltbr(int argc, char *argv[])
131 {
132     const char     *desc[] = {
133         "[THISMODULE] plots the distance between all combination of charged groups",
134         "as a function of time. The groups are combined in different ways.",
135         "A minimum distance can be given (i.e. a cut-off), such that groups",
136         "that are never closer than that distance will not be plotted.[PAR]",
137         "Output will be in a number of fixed filenames, [TT]min-min.xvg[tt], [TT]plus-min.xvg[tt]",
138         "and [TT]plus-plus.xvg[tt], or files for every individual ion pair if the [TT]-sep[tt]",
139         "option is selected. In this case, files are named as [TT]sb-(Resname)(Resnr)-(Atomnr)[tt].",
140         "There may be [BB]many[bb] such files."
141     };
142     static gmx_bool bSep     = FALSE;
143     static real     truncate = 1000.0;
144     t_pargs         pa[]     = {
145         { "-t",   FALSE, etREAL, {&truncate},
146           "Groups that are never closer than this distance are not plotted" },
147         { "-sep", FALSE, etBOOL, {&bSep},
148           "Use separate files for each interaction (may be MANY)" }
149     };
150     t_filenm        fnm[] = {
151         { efTRX, "-f",  NULL, ffREAD },
152         { efTPR, NULL,  NULL, ffREAD },
153     };
154 #define NFILE asize(fnm)
155
156     FILE              *out[3], *fp;
157     static const char *title[3] = {
158         "Distance between positively charged groups",
159         "Distance between negatively charged groups",
160         "Distance between oppositely charged groups"
161     };
162     static const char *fn[3] = {
163         "plus-plus.xvg",
164         "min-min.xvg",
165         "plus-min.xvg"
166     };
167     int                nset[3] = {0, 0, 0};
168
169     t_topology        *top;
170     int                ePBC;
171     char              *buf;
172     t_trxstatus       *status;
173     int                i, j, k, m, nnn, teller, ncg;
174     real               t, *time, qi, qj;
175     t_charge          *cg;
176     real            ***cgdist;
177     int              **nWithin;
178
179     t_pbc              pbc;
180     rvec              *x;
181     matrix             box;
182     output_env_t       oenv;
183
184     if (!parse_common_args(&argc, argv, PCA_CAN_TIME,
185                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
186     {
187         return 0;
188     }
189
190     top = read_top(ftp2fn(efTPR, NFILE, fnm), &ePBC);
191     cg  = mk_charge(&top->atoms, &(top->cgs), &ncg);
192     snew(cgdist, ncg);
193     snew(nWithin, ncg);
194     for (i = 0; (i < ncg); i++)
195     {
196         snew(cgdist[i], ncg);
197         snew(nWithin[i], ncg);
198     }
199
200     read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
201
202     teller = 0;
203     time   = NULL;
204     do
205     {
206         srenew(time, teller+1);
207         time[teller] = t;
208
209         set_pbc(&pbc, ePBC, box);
210
211         for (i = 0; (i < ncg); i++)
212         {
213             for (j = i+1; (j < ncg); j++)
214             {
215                 srenew(cgdist[i][j], teller+1);
216                 cgdist[i][j][teller] =
217                     calc_dist(&pbc, x, &(top->cgs), cg[i].cg, cg[j].cg);
218                 if (cgdist[i][j][teller] < truncate)
219                 {
220                     nWithin[i][j] = 1;
221                 }
222             }
223         }
224
225         teller++;
226     }
227     while (read_next_x(oenv, status, &t, x, box));
228     fprintf(stderr, "\n");
229     close_trj(status);
230
231     if (bSep)
232     {
233         snew(buf, 256);
234         for (i = 0; (i < ncg); i++)
235         {
236             for (j = i+1; (j < ncg); j++)
237             {
238                 if (nWithin[i][j])
239                 {
240                     sprintf(buf, "sb-%s:%s.xvg", cg[i].label, cg[j].label);
241                     fp = xvgropen(buf, buf, "Time (ps)", "Distance (nm)", oenv);
242                     for (k = 0; (k < teller); k++)
243                     {
244                         fprintf(fp, "%10g  %10g\n", time[k], cgdist[i][j][k]);
245                     }
246                     xvgrclose(fp);
247                 }
248             }
249         }
250         sfree(buf);
251     }
252     else
253     {
254
255         for (m = 0; (m < 3); m++)
256         {
257             out[m] = xvgropen(fn[m], title[m], "Time (ps)", "Distance (nm)", oenv);
258         }
259
260         snew(buf, 256);
261         for (i = 0; (i < ncg); i++)
262         {
263             qi = cg[i].q;
264             for (j = i+1; (j < ncg); j++)
265             {
266                 qj = cg[j].q;
267                 if (nWithin[i][j])
268                 {
269                     sprintf(buf, "%s:%s", cg[i].label, cg[j].label);
270                     if (qi*qj < 0)
271                     {
272                         nnn = 2;
273                     }
274                     else if (qi+qj > 0)
275                     {
276                         nnn = 0;
277                     }
278                     else
279                     {
280                         nnn = 1;
281                     }
282
283                     if (nset[nnn] == 0)
284                     {
285                         xvgr_legend(out[nnn], 1, (const char**)&buf, oenv);
286                     }
287                     else
288                     {
289                         if (output_env_get_xvg_format(oenv) == exvgXMGR)
290                         {
291                             fprintf(out[nnn], "@ legend string %d \"%s\"\n", nset[nnn], buf);
292                         }
293                         else if (output_env_get_xvg_format(oenv) == exvgXMGRACE)
294                         {
295                             fprintf(out[nnn], "@ s%d legend \"%s\"\n", nset[nnn], buf);
296                         }
297                     }
298                     nset[nnn]++;
299                     nWithin[i][j] = nnn+1;
300                 }
301             }
302         }
303         for (k = 0; (k < teller); k++)
304         {
305             for (m = 0; (m < 3); m++)
306             {
307                 fprintf(out[m], "%10g", time[k]);
308             }
309
310             for (i = 0; (i < ncg); i++)
311             {
312                 for (j = i+1; (j < ncg); j++)
313                 {
314                     nnn = nWithin[i][j];
315                     if (nnn > 0)
316                     {
317                         fprintf(out[nnn-1], "  %10g", cgdist[i][j][k]);
318                     }
319                 }
320             }
321             for (m = 0; (m < 3); m++)
322             {
323                 fprintf(out[m], "\n");
324             }
325         }
326         for (m = 0; (m < 3); m++)
327         {
328             xvgrclose(out[m]);
329             if (nset[m] == 0)
330             {
331                 remove(fn[m]);
332             }
333         }
334     }
335
336     return 0;
337 }