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