Merge branch release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_saltbr.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38 #include <math.h>
39 #include <string.h>
40
41 #include "macros.h"
42 #include "vec.h"
43 #include "sysstuff.h"
44 #include "typedefs.h"
45 #include "filenm.h"
46 #include "statutil.h"
47 #include "futil.h"
48 #include "gmx_fatal.h"
49 #include "smalloc.h"
50 #include "pbc.h"
51 #include "xvgr.h"
52 #include "gmx_ana.h"
53
54
55 typedef struct {
56     char *label;
57     int   cg;
58     real  q;
59 } t_charge;
60
61 static t_charge *mk_charge(t_atoms *atoms, t_block *cgs, int *nncg)
62 {
63     t_charge *cg = NULL;
64     char      buf[32];
65     int       i, j, ncg, resnr, anr;
66     real      qq;
67
68     /* Find the charged groups */
69     ncg = 0;
70     for (i = 0; (i < cgs->nr); i++)
71     {
72         qq = 0.0;
73         for (j = cgs->index[i]; (j < cgs->index[i+1]); j++)
74         {
75             qq += atoms->atom[j].q;
76         }
77         if (fabs(qq) > 1.0e-5)
78         {
79             srenew(cg, ncg+1);
80             cg[ncg].q  = qq;
81             cg[ncg].cg = i;
82             anr        = cgs->index[i];
83             resnr      = atoms->atom[anr].resind;
84             sprintf(buf, "%s%d-%d",
85                     *(atoms->resinfo[resnr].name),
86                     atoms->resinfo[resnr].nr,
87                     anr+1);
88             cg[ncg].label = strdup(buf);
89             ncg++;
90         }
91     }
92     *nncg = ncg;
93
94     for (i = 0; (i < ncg); i++)
95     {
96         printf("CG: %10s Q: %6g  Atoms:",
97                cg[i].label, cg[i].q);
98         for (j = cgs->index[cg[i].cg]; (j < cgs->index[cg[i].cg+1]); j++)
99         {
100             printf(" %4d", j);
101         }
102         printf("\n");
103     }
104
105     return cg;
106 }
107
108 static real calc_dist(t_pbc *pbc, rvec x[], t_block *cgs, int icg, int jcg)
109 {
110     int  i, j;
111     rvec dx;
112     real d2, mindist2 = 1000;
113
114     for (i = cgs->index[icg]; (i < cgs->index[icg+1]); i++)
115     {
116         for (j = cgs->index[jcg]; (j < cgs->index[jcg+1]); j++)
117         {
118             pbc_dx(pbc, x[i], x[j], dx);
119             d2 = norm2(dx);
120             if (d2 < mindist2)
121             {
122                 mindist2 = d2;
123             }
124         }
125     }
126     return sqrt(mindist2);
127 }
128
129 int gmx_saltbr(int argc, char *argv[])
130 {
131     const char     *desc[] = {
132         "[TT]g_saltbr[tt] plots the distance between all combination of charged groups",
133         "as a function of time. The groups are combined in different ways.",
134         "A minimum distance can be given (i.e. a cut-off), such that groups",
135         "that are never closer than that distance will not be plotted.[PAR]",
136         "Output will be in a number of fixed filenames, [TT]min-min.xvg[tt], [TT]plus-min.xvg[tt]",
137         "and [TT]plus-plus.xvg[tt], or files for every individual ion pair if the [TT]-sep[tt]",
138         "option is selected. In this case, files are named as [TT]sb-(Resname)(Resnr)-(Atomnr)[tt].",
139         "There may be [BB]many[bb] such files."
140     };
141     static gmx_bool bSep     = FALSE;
142     static real     truncate = 1000.0;
143     t_pargs         pa[]     = {
144         { "-t",   FALSE, etREAL, {&truncate},
145           "Groups that are never closer than this distance are not plotted" },
146         { "-sep", FALSE, etBOOL, {&bSep},
147           "Use separate files for each interaction (may be MANY)" }
148     };
149     t_filenm        fnm[] = {
150         { efTRX, "-f",  NULL, ffREAD },
151         { efTPX, NULL,  NULL, ffREAD },
152     };
153 #define NFILE asize(fnm)
154
155     FILE              *out[3], *fp;
156     static const char *title[3] = {
157         "Distance between positively charged groups",
158         "Distance between negatively charged groups",
159         "Distance between oppositely charged groups"
160     };
161     static const char *fn[3] = {
162         "plus-plus.xvg",
163         "min-min.xvg",
164         "plus-min.xvg"
165     };
166     int                nset[3] = {0, 0, 0};
167
168     t_topology        *top;
169     int                ePBC;
170     char              *buf;
171     t_trxstatus       *status;
172     int                i, j, k, m, nnn, teller, ncg, n1, n2, n3, natoms;
173     real               t, *time, qi, qj;
174     t_charge          *cg;
175     real            ***cgdist;
176     int              **nWithin;
177
178     double             t0, dt;
179     char               label[234];
180     t_pbc              pbc;
181     rvec              *x;
182     matrix             box;
183     output_env_t       oenv;
184
185     parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_BE_NICE,
186                       NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
187
188     top = read_top(ftp2fn(efTPX, NFILE, fnm), &ePBC);
189     cg  = mk_charge(&top->atoms, &(top->cgs), &ncg);
190     snew(cgdist, ncg);
191     snew(nWithin, ncg);
192     for (i = 0; (i < ncg); i++)
193     {
194         snew(cgdist[i], ncg);
195         snew(nWithin[i], ncg);
196     }
197
198     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
199
200     teller = 0;
201     time   = NULL;
202     do
203     {
204         srenew(time, teller+1);
205         time[teller] = t;
206
207         set_pbc(&pbc, ePBC, box);
208
209         for (i = 0; (i < ncg); i++)
210         {
211             for (j = i+1; (j < ncg); j++)
212             {
213                 srenew(cgdist[i][j], teller+1);
214                 cgdist[i][j][teller] =
215                     calc_dist(&pbc, x, &(top->cgs), cg[i].cg, cg[j].cg);
216                 if (cgdist[i][j][teller] < truncate)
217                 {
218                     nWithin[i][j] = 1;
219                 }
220             }
221         }
222
223         teller++;
224     }
225     while (read_next_x(oenv, status, &t, x, box));
226     fprintf(stderr, "\n");
227     close_trj(status);
228
229     if (bSep)
230     {
231         snew(buf, 256);
232         for (i = 0; (i < ncg); i++)
233         {
234             for (j = i+1; (j < ncg); j++)
235             {
236                 if (nWithin[i][j])
237                 {
238                     sprintf(buf, "sb-%s:%s.xvg", cg[i].label, cg[j].label);
239                     fp = xvgropen(buf, buf, "Time (ps)", "Distance (nm)", oenv);
240                     for (k = 0; (k < teller); k++)
241                     {
242                         fprintf(fp, "%10g  %10g\n", time[k], cgdist[i][j][k]);
243                     }
244                     ffclose(fp);
245                 }
246             }
247         }
248         sfree(buf);
249     }
250     else
251     {
252
253         for (m = 0; (m < 3); m++)
254         {
255             out[m] = xvgropen(fn[m], title[m], "Time (ps)", "Distance (nm)", oenv);
256         }
257
258         snew(buf, 256);
259         for (i = 0; (i < ncg); i++)
260         {
261             qi = cg[i].q;
262             for (j = i+1; (j < ncg); j++)
263             {
264                 qj = cg[j].q;
265                 if (nWithin[i][j])
266                 {
267                     sprintf(buf, "%s:%s", cg[i].label, cg[j].label);
268                     if (qi*qj < 0)
269                     {
270                         nnn = 2;
271                     }
272                     else if (qi+qj > 0)
273                     {
274                         nnn = 0;
275                     }
276                     else
277                     {
278                         nnn = 1;
279                     }
280
281                     if (nset[nnn] == 0)
282                     {
283                         xvgr_legend(out[nnn], 1, (const char**)&buf, oenv);
284                     }
285                     else
286                     {
287                         if (output_env_get_xvg_format(oenv) == exvgXMGR)
288                         {
289                             fprintf(out[nnn], "@ legend string %d \"%s\"\n", nset[nnn], buf);
290                         }
291                         else if (output_env_get_xvg_format(oenv) == exvgXMGRACE)
292                         {
293                             fprintf(out[nnn], "@ s%d legend \"%s\"\n", nset[nnn], buf);
294                         }
295                     }
296                     nset[nnn]++;
297                     nWithin[i][j] = nnn+1;
298                 }
299             }
300         }
301         for (k = 0; (k < teller); k++)
302         {
303             for (m = 0; (m < 3); m++)
304             {
305                 fprintf(out[m], "%10g", time[k]);
306             }
307
308             for (i = 0; (i < ncg); i++)
309             {
310                 for (j = i+1; (j < ncg); j++)
311                 {
312                     nnn = nWithin[i][j];
313                     if (nnn > 0)
314                     {
315                         fprintf(out[nnn-1], "  %10g", cgdist[i][j][k]);
316                     }
317                 }
318             }
319             for (m = 0; (m < 3); m++)
320             {
321                 fprintf(out[m], "\n");
322             }
323         }
324         for (m = 0; (m < 3); m++)
325         {
326             ffclose(out[m]);
327             if (nset[m] == 0)
328             {
329                 remove(fn[m]);
330             }
331         }
332     }
333
334     return 0;
335 }