2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
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"
56 #include "gromacs/utility/fatalerror.h"
64 static t_charge *mk_charge(t_atoms *atoms, t_block *cgs, int *nncg)
68 int i, j, ncg, resnr, anr;
71 /* Find the charged groups */
73 for (i = 0; (i < cgs->nr); i++)
76 for (j = cgs->index[i]; (j < cgs->index[i+1]); j++)
78 qq += atoms->atom[j].q;
80 if (fabs(qq) > 1.0e-5)
86 resnr = atoms->atom[anr].resind;
87 sprintf(buf, "%s%d-%d",
88 *(atoms->resinfo[resnr].name),
89 atoms->resinfo[resnr].nr,
91 cg[ncg].label = gmx_strdup(buf);
97 for (i = 0; (i < ncg); i++)
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++)
111 static real calc_dist(t_pbc *pbc, rvec x[], t_block *cgs, int icg, int jcg)
115 real d2, mindist2 = 1000;
117 for (i = cgs->index[icg]; (i < cgs->index[icg+1]); i++)
119 for (j = cgs->index[jcg]; (j < cgs->index[jcg+1]); j++)
121 pbc_dx(pbc, x[i], x[j], dx);
129 return sqrt(mindist2);
132 int gmx_saltbr(int argc, char *argv[])
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."
144 static gmx_bool bSep = FALSE;
145 static real truncate = 1000.0;
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)" }
153 { efTRX, "-f", NULL, ffREAD },
154 { efTPX, NULL, NULL, ffREAD },
156 #define NFILE asize(fnm)
159 static const char *title[3] = {
160 "Distance between positively charged groups",
161 "Distance between negatively charged groups",
162 "Distance between oppositely charged groups"
164 static const char *fn[3] = {
169 int nset[3] = {0, 0, 0};
175 int i, j, k, m, nnn, teller, ncg, n1, n2, n3, natoms;
176 real t, *time, qi, qj;
188 if (!parse_common_args(&argc, argv, PCA_CAN_TIME,
189 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
194 top = read_top(ftp2fn(efTPX, NFILE, fnm), &ePBC);
195 cg = mk_charge(&top->atoms, &(top->cgs), &ncg);
198 for (i = 0; (i < ncg); i++)
200 snew(cgdist[i], ncg);
201 snew(nWithin[i], ncg);
204 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
210 srenew(time, teller+1);
213 set_pbc(&pbc, ePBC, box);
215 for (i = 0; (i < ncg); i++)
217 for (j = i+1; (j < ncg); j++)
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)
231 while (read_next_x(oenv, status, &t, x, box));
232 fprintf(stderr, "\n");
238 for (i = 0; (i < ncg); i++)
240 for (j = i+1; (j < ncg); j++)
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++)
248 fprintf(fp, "%10g %10g\n", time[k], cgdist[i][j][k]);
259 for (m = 0; (m < 3); m++)
261 out[m] = xvgropen(fn[m], title[m], "Time (ps)", "Distance (nm)", oenv);
265 for (i = 0; (i < ncg); i++)
268 for (j = i+1; (j < ncg); j++)
273 sprintf(buf, "%s:%s", cg[i].label, cg[j].label);
289 xvgr_legend(out[nnn], 1, (const char**)&buf, oenv);
293 if (output_env_get_xvg_format(oenv) == exvgXMGR)
295 fprintf(out[nnn], "@ legend string %d \"%s\"\n", nset[nnn], buf);
297 else if (output_env_get_xvg_format(oenv) == exvgXMGRACE)
299 fprintf(out[nnn], "@ s%d legend \"%s\"\n", nset[nnn], buf);
303 nWithin[i][j] = nnn+1;
307 for (k = 0; (k < teller); k++)
309 for (m = 0; (m < 3); m++)
311 fprintf(out[m], "%10g", time[k]);
314 for (i = 0; (i < ncg); i++)
316 for (j = i+1; (j < ncg); j++)
321 fprintf(out[nnn-1], " %10g", cgdist[i][j][k]);
325 for (m = 0; (m < 3); m++)
327 fprintf(out[m], "\n");
330 for (m = 0; (m < 3); m++)