Merge branch 'release-4-6'
[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 "copyrite.h"
48 #include "futil.h"
49 #include "gmx_fatal.h"
50 #include "smalloc.h"
51 #include "pbc.h"
52 #include "xvgr.h"
53 #include "gmx_ana.h"
54
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 = 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(" %4u", 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         "[TT]g_saltbr[tt] 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         { efTPX, 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     parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_BE_NICE,
187                       NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
188
189     top = read_top(ftp2fn(efTPX, NFILE, fnm), &ePBC);
190     cg  = mk_charge(&top->atoms, &(top->cgs), &ncg);
191     snew(cgdist, ncg);
192     snew(nWithin, ncg);
193     for (i = 0; (i < ncg); i++)
194     {
195         snew(cgdist[i], ncg);
196         snew(nWithin[i], ncg);
197     }
198
199     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
200
201     teller = 0;
202     time   = NULL;
203     do
204     {
205         srenew(time, teller+1);
206         time[teller] = t;
207
208         set_pbc(&pbc, ePBC, box);
209
210         for (i = 0; (i < ncg); i++)
211         {
212             for (j = i+1; (j < ncg); j++)
213             {
214                 srenew(cgdist[i][j], teller+1);
215                 cgdist[i][j][teller] =
216                     calc_dist(&pbc, x, &(top->cgs), cg[i].cg, cg[j].cg);
217                 if (cgdist[i][j][teller] < truncate)
218                 {
219                     nWithin[i][j] = 1;
220                 }
221             }
222         }
223
224         teller++;
225     }
226     while (read_next_x(oenv, status, &t, natoms, x, box));
227     fprintf(stderr, "\n");
228     close_trj(status);
229
230     if (bSep)
231     {
232         snew(buf, 256);
233         for (i = 0; (i < ncg); i++)
234         {
235             for (j = i+1; (j < ncg); j++)
236             {
237                 if (nWithin[i][j])
238                 {
239                     sprintf(buf, "sb-%s:%s.xvg", cg[i].label, cg[j].label);
240                     fp = xvgropen(buf, buf, "Time (ps)", "Distance (nm)", oenv);
241                     for (k = 0; (k < teller); k++)
242                     {
243                         fprintf(fp, "%10g  %10g\n", time[k], cgdist[i][j][k]);
244                     }
245                     ffclose(fp);
246                 }
247             }
248         }
249         sfree(buf);
250     }
251     else
252     {
253
254         for (m = 0; (m < 3); m++)
255         {
256             out[m] = xvgropen(fn[m], title[m], "Time (ps)", "Distance (nm)", oenv);
257         }
258
259         snew(buf, 256);
260         for (i = 0; (i < ncg); i++)
261         {
262             qi = cg[i].q;
263             for (j = i+1; (j < ncg); j++)
264             {
265                 qj = cg[j].q;
266                 if (nWithin[i][j])
267                 {
268                     sprintf(buf, "%s:%s", cg[i].label, cg[j].label);
269                     if (qi*qj < 0)
270                     {
271                         nnn = 2;
272                     }
273                     else if (qi+qj > 0)
274                     {
275                         nnn = 0;
276                     }
277                     else
278                     {
279                         nnn = 1;
280                     }
281
282                     if (nset[nnn] == 0)
283                     {
284                         xvgr_legend(out[nnn], 1, (const char**)&buf, oenv);
285                     }
286                     else
287                     {
288                         if (output_env_get_xvg_format(oenv) == exvgXMGR)
289                         {
290                             fprintf(out[nnn], "@ legend string %d \"%s\"\n", nset[nnn], buf);
291                         }
292                         else if (output_env_get_xvg_format(oenv) == exvgXMGRACE)
293                         {
294                             fprintf(out[nnn], "@ s%d legend \"%s\"\n", nset[nnn], buf);
295                         }
296                     }
297                     nset[nnn]++;
298                     nWithin[i][j] = nnn+1;
299                 }
300             }
301         }
302         for (k = 0; (k < teller); k++)
303         {
304             for (m = 0; (m < 3); m++)
305             {
306                 fprintf(out[m], "%10g", time[k]);
307             }
308
309             for (i = 0; (i < ncg); i++)
310             {
311                 for (j = i+1; (j < ncg); j++)
312                 {
313                     nnn = nWithin[i][j];
314                     if (nnn > 0)
315                     {
316                         fprintf(out[nnn-1], "  %10g", cgdist[i][j][k]);
317                     }
318                 }
319             }
320             for (m = 0; (m < 3); m++)
321             {
322                 fprintf(out[m], "\n");
323             }
324         }
325         for (m = 0; (m < 3); m++)
326         {
327             ffclose(out[m]);
328             if (nset[m] == 0)
329             {
330                 remove(fn[m]);
331             }
332         }
333     }
334     thanx(stderr);
335
336     return 0;
337 }