Moved first batch of analysis tool source to C++
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_saltbr.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,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 <math.h>
40 #include <string.h>
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 (fabs(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 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, n1, n2, n3, natoms;
174     real               t, *time, qi, qj;
175     t_charge          *cg;
176     real            ***cgdist;
177     int              **nWithin;
178
179     double             t0, dt;
180     char               label[234];
181     t_pbc              pbc;
182     rvec              *x;
183     matrix             box;
184     output_env_t       oenv;
185
186     if (!parse_common_args(&argc, argv, PCA_CAN_TIME,
187                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
188     {
189         return 0;
190     }
191
192     top = read_top(ftp2fn(efTPR, NFILE, fnm), &ePBC);
193     cg  = mk_charge(&top->atoms, &(top->cgs), &ncg);
194     snew(cgdist, ncg);
195     snew(nWithin, ncg);
196     for (i = 0; (i < ncg); i++)
197     {
198         snew(cgdist[i], ncg);
199         snew(nWithin[i], ncg);
200     }
201
202     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
203
204     teller = 0;
205     time   = NULL;
206     do
207     {
208         srenew(time, teller+1);
209         time[teller] = t;
210
211         set_pbc(&pbc, ePBC, box);
212
213         for (i = 0; (i < ncg); i++)
214         {
215             for (j = i+1; (j < ncg); j++)
216             {
217                 srenew(cgdist[i][j], teller+1);
218                 cgdist[i][j][teller] =
219                     calc_dist(&pbc, x, &(top->cgs), cg[i].cg, cg[j].cg);
220                 if (cgdist[i][j][teller] < truncate)
221                 {
222                     nWithin[i][j] = 1;
223                 }
224             }
225         }
226
227         teller++;
228     }
229     while (read_next_x(oenv, status, &t, x, box));
230     fprintf(stderr, "\n");
231     close_trj(status);
232
233     if (bSep)
234     {
235         snew(buf, 256);
236         for (i = 0; (i < ncg); i++)
237         {
238             for (j = i+1; (j < ncg); j++)
239             {
240                 if (nWithin[i][j])
241                 {
242                     sprintf(buf, "sb-%s:%s.xvg", cg[i].label, cg[j].label);
243                     fp = xvgropen(buf, buf, "Time (ps)", "Distance (nm)", oenv);
244                     for (k = 0; (k < teller); k++)
245                     {
246                         fprintf(fp, "%10g  %10g\n", time[k], cgdist[i][j][k]);
247                     }
248                     xvgrclose(fp);
249                 }
250             }
251         }
252         sfree(buf);
253     }
254     else
255     {
256
257         for (m = 0; (m < 3); m++)
258         {
259             out[m] = xvgropen(fn[m], title[m], "Time (ps)", "Distance (nm)", oenv);
260         }
261
262         snew(buf, 256);
263         for (i = 0; (i < ncg); i++)
264         {
265             qi = cg[i].q;
266             for (j = i+1; (j < ncg); j++)
267             {
268                 qj = cg[j].q;
269                 if (nWithin[i][j])
270                 {
271                     sprintf(buf, "%s:%s", cg[i].label, cg[j].label);
272                     if (qi*qj < 0)
273                     {
274                         nnn = 2;
275                     }
276                     else if (qi+qj > 0)
277                     {
278                         nnn = 0;
279                     }
280                     else
281                     {
282                         nnn = 1;
283                     }
284
285                     if (nset[nnn] == 0)
286                     {
287                         xvgr_legend(out[nnn], 1, (const char**)&buf, oenv);
288                     }
289                     else
290                     {
291                         if (output_env_get_xvg_format(oenv) == exvgXMGR)
292                         {
293                             fprintf(out[nnn], "@ legend string %d \"%s\"\n", nset[nnn], buf);
294                         }
295                         else if (output_env_get_xvg_format(oenv) == exvgXMGRACE)
296                         {
297                             fprintf(out[nnn], "@ s%d legend \"%s\"\n", nset[nnn], buf);
298                         }
299                     }
300                     nset[nnn]++;
301                     nWithin[i][j] = nnn+1;
302                 }
303             }
304         }
305         for (k = 0; (k < teller); k++)
306         {
307             for (m = 0; (m < 3); m++)
308             {
309                 fprintf(out[m], "%10g", time[k]);
310             }
311
312             for (i = 0; (i < ncg); i++)
313             {
314                 for (j = i+1; (j < ncg); j++)
315                 {
316                     nnn = nWithin[i][j];
317                     if (nnn > 0)
318                     {
319                         fprintf(out[nnn-1], "  %10g", cgdist[i][j][k]);
320                     }
321                 }
322             }
323             for (m = 0; (m < 3); m++)
324             {
325                 fprintf(out[m], "\n");
326             }
327         }
328         for (m = 0; (m < 3); m++)
329         {
330             xvgrclose(out[m]);
331             if (nset[m] == 0)
332             {
333                 remove(fn[m]);
334             }
335         }
336     }
337
338     return 0;
339 }