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