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