Apply re-formatting to C++ in src/ tree.
[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 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include <cmath>
41 #include <cstring>
42
43 #include "gromacs/commandline/filenm.h"
44 #include "gromacs/commandline/pargs.h"
45 #include "gromacs/fileio/trxio.h"
46 #include "gromacs/fileio/xvgr.h"
47 #include "gromacs/gmxana/gmx_ana.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/pbcutil/pbc.h"
50 #include "gromacs/topology/topology.h"
51 #include "gromacs/utility/arraysize.h"
52 #include "gromacs/utility/cstringutil.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/smalloc.h"
56
57 typedef struct
58 {
59     char* label;
60     int   cg;
61     real  q;
62 } t_charge;
63
64 static t_charge* mk_charge(const t_atoms* atoms, int* nncg)
65 {
66     t_charge* cg = nullptr;
67     char      buf[32];
68     int       i, ncg, resnr, anr;
69     real      qq;
70
71     /* Find the charged groups */
72     ncg = 0;
73     for (i = 0; i < atoms->nr; i++)
74     {
75         qq = atoms->atom[i].q;
76         if (std::abs(qq) > 1.0e-5)
77         {
78             srenew(cg, ncg + 1);
79             cg[ncg].q  = qq;
80             cg[ncg].cg = i;
81             anr        = i;
82             resnr      = atoms->atom[anr].resind;
83             sprintf(buf, "%s%d-%d", *(atoms->resinfo[resnr].name), atoms->resinfo[resnr].nr, anr + 1);
84             cg[ncg].label = gmx_strdup(buf);
85             ncg++;
86         }
87     }
88     *nncg = ncg;
89
90     for (i = 0; (i < ncg); i++)
91     {
92         printf("CG: %10s Q: %6g  Atoms:", cg[i].label, cg[i].q);
93         printf(" %4d", cg[i].cg);
94         printf("\n");
95     }
96
97     return cg;
98 }
99
100 int gmx_saltbr(int argc, char* argv[])
101 {
102     const char* desc[] = {
103         "[THISMODULE] plots the distance between all combination of charged groups",
104         "as a function of time. The groups are combined in different ways.",
105         "A minimum distance can be given (i.e. a cut-off), such that groups",
106         "that are never closer than that distance will not be plotted.[PAR]",
107         "Output will be in a number of fixed filenames, [TT]min-min.xvg[tt], [TT]plus-min.xvg[tt]",
108         "and [TT]plus-plus.xvg[tt], or files for every individual ion pair if the [TT]-sep[tt]",
109         "option is selected. In this case, files are named as ",
110         "[TT]sb-(Resname)(Resnr)-(Atomnr)[tt].",
111         "There may be [BB]many[bb] such files."
112     };
113     static gmx_bool bSep     = FALSE;
114     static real     truncate = 1000.0;
115     t_pargs         pa[]     = { { "-t",
116                        FALSE,
117                        etREAL,
118                        { &truncate },
119                        "Groups that are never closer than this distance are not plotted" },
120                      { "-sep",
121                        FALSE,
122                        etBOOL,
123                        { &bSep },
124                        "Use separate files for each interaction (may be MANY)" } };
125     t_filenm        fnm[]    = {
126         { efTRX, "-f", nullptr, ffREAD },
127         { efTPR, nullptr, nullptr, ffREAD },
128     };
129 #define NFILE asize(fnm)
130
131     FILE *             out[3], *fp;
132     static const char* title[3] = { "Distance between positively charged groups",
133                                     "Distance between negatively charged groups",
134                                     "Distance between oppositely charged groups" };
135     static const char* fn[3]    = { "plus-plus.xvg", "min-min.xvg", "plus-min.xvg" };
136     int                nset[3]  = { 0, 0, 0 };
137
138     t_topology*  top;
139     PbcType      pbcType;
140     char*        buf;
141     t_trxstatus* status;
142     int          i, j, k, m, nnn, teller, ncg;
143     real         t, *time, qi, qj;
144     t_charge*    cg;
145     real***      cgdist;
146     int**        nWithin;
147
148     t_pbc             pbc;
149     rvec*             x;
150     matrix            box;
151     gmx_output_env_t* oenv;
152
153     if (!parse_common_args(
154                 &argc, argv, PCA_CAN_TIME, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
155     {
156         return 0;
157     }
158
159     top = read_top(ftp2fn(efTPR, NFILE, fnm), &pbcType);
160     cg  = mk_charge(&top->atoms, &ncg);
161     snew(cgdist, ncg);
162     snew(nWithin, ncg);
163     for (i = 0; (i < ncg); i++)
164     {
165         snew(cgdist[i], ncg);
166         snew(nWithin[i], ncg);
167     }
168
169     read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
170
171     teller = 0;
172     time   = nullptr;
173     do
174     {
175         srenew(time, teller + 1);
176         time[teller] = t;
177
178         set_pbc(&pbc, pbcType, box);
179
180         for (i = 0; (i < ncg); i++)
181         {
182             for (j = i + 1; (j < ncg); j++)
183             {
184                 srenew(cgdist[i][j], teller + 1);
185                 rvec dx;
186                 pbc_dx(&pbc, x[cg[i].cg], x[cg[j].cg], dx);
187                 cgdist[i][j][teller] = norm(dx);
188                 if (cgdist[i][j][teller] < truncate)
189                 {
190                     nWithin[i][j] = 1;
191                 }
192             }
193         }
194
195         teller++;
196     } while (read_next_x(oenv, status, &t, x, box));
197     fprintf(stderr, "\n");
198     close_trx(status);
199
200     if (bSep)
201     {
202         snew(buf, 256);
203         for (i = 0; (i < ncg); i++)
204         {
205             for (j = i + 1; (j < ncg); j++)
206             {
207                 if (nWithin[i][j])
208                 {
209                     sprintf(buf, "sb-%s:%s.xvg", cg[i].label, cg[j].label);
210                     fp = xvgropen(buf, buf, "Time (ps)", "Distance (nm)", oenv);
211                     for (k = 0; (k < teller); k++)
212                     {
213                         fprintf(fp, "%10g  %10g\n", time[k], cgdist[i][j][k]);
214                     }
215                     xvgrclose(fp);
216                 }
217             }
218         }
219         sfree(buf);
220     }
221     else
222     {
223
224         for (m = 0; (m < 3); m++)
225         {
226             out[m] = xvgropen(fn[m], title[m], "Time (ps)", "Distance (nm)", oenv);
227         }
228
229         snew(buf, 256);
230         for (i = 0; (i < ncg); i++)
231         {
232             qi = cg[i].q;
233             for (j = i + 1; (j < ncg); j++)
234             {
235                 qj = cg[j].q;
236                 if (nWithin[i][j])
237                 {
238                     sprintf(buf, "%s:%s", cg[i].label, cg[j].label);
239                     if (qi * qj < 0)
240                     {
241                         nnn = 2;
242                     }
243                     else if (qi + qj > 0)
244                     {
245                         nnn = 0;
246                     }
247                     else
248                     {
249                         nnn = 1;
250                     }
251
252                     if (nset[nnn] == 0)
253                     {
254                         xvgr_legend(out[nnn], 1, &buf, oenv);
255                     }
256                     else
257                     {
258                         if (output_env_get_xvg_format(oenv) == XvgFormat::Xmgr)
259                         {
260                             fprintf(out[nnn], "@ legend string %d \"%s\"\n", nset[nnn], buf);
261                         }
262                         else if (output_env_get_xvg_format(oenv) == XvgFormat::Xmgrace)
263                         {
264                             fprintf(out[nnn], "@ s%d legend \"%s\"\n", nset[nnn], buf);
265                         }
266                     }
267                     nset[nnn]++;
268                     nWithin[i][j] = nnn + 1;
269                 }
270             }
271         }
272         for (k = 0; (k < teller); k++)
273         {
274             for (m = 0; (m < 3); m++)
275             {
276                 fprintf(out[m], "%10g", time[k]);
277             }
278
279             for (i = 0; (i < ncg); i++)
280             {
281                 for (j = i + 1; (j < ncg); j++)
282                 {
283                     nnn = nWithin[i][j];
284                     if (nnn > 0)
285                     {
286                         fprintf(out[nnn - 1], "  %10g", cgdist[i][j][k]);
287                     }
288                 }
289             }
290             for (m = 0; (m < 3); m++)
291             {
292                 fprintf(out[m], "\n");
293             }
294         }
295         for (m = 0; (m < 3); m++)
296         {
297             xvgrclose(out[m]);
298             if (nset[m] == 0)
299             {
300                 remove(fn[m]);
301             }
302         }
303     }
304
305     return 0;
306 }