3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
48 #include "gmx_fatal.h"
61 static t_charge *mk_charge(t_atoms *atoms, t_block *cgs, int *nncg)
65 int i, j, ncg, resnr, anr;
68 /* Find the charged groups */
70 for (i = 0; (i < cgs->nr); i++)
73 for (j = cgs->index[i]; (j < cgs->index[i+1]); j++)
75 qq += atoms->atom[j].q;
77 if (fabs(qq) > 1.0e-5)
83 resnr = atoms->atom[anr].resind;
84 sprintf(buf, "%s%d-%d",
85 *(atoms->resinfo[resnr].name),
86 atoms->resinfo[resnr].nr,
88 cg[ncg].label = strdup(buf);
94 for (i = 0; (i < ncg); i++)
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++)
108 static real calc_dist(t_pbc *pbc, rvec x[], t_block *cgs, int icg, int jcg)
112 real d2, mindist2 = 1000;
114 for (i = cgs->index[icg]; (i < cgs->index[icg+1]); i++)
116 for (j = cgs->index[jcg]; (j < cgs->index[jcg+1]); j++)
118 pbc_dx(pbc, x[i], x[j], dx);
126 return sqrt(mindist2);
129 int gmx_saltbr(int argc, char *argv[])
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."
141 static gmx_bool bSep = FALSE;
142 static real truncate = 1000.0;
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)" }
150 { efTRX, "-f", NULL, ffREAD },
151 { efTPX, NULL, NULL, ffREAD },
153 #define NFILE asize(fnm)
156 static const char *title[3] = {
157 "Distance between positively charged groups",
158 "Distance between negatively charged groups",
159 "Distance between oppositely charged groups"
161 static const char *fn[3] = {
166 int nset[3] = {0, 0, 0};
172 int i, j, k, m, nnn, teller, ncg, n1, n2, n3, natoms;
173 real t, *time, qi, qj;
185 parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_BE_NICE,
186 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
188 top = read_top(ftp2fn(efTPX, NFILE, fnm), &ePBC);
189 cg = mk_charge(&top->atoms, &(top->cgs), &ncg);
192 for (i = 0; (i < ncg); i++)
194 snew(cgdist[i], ncg);
195 snew(nWithin[i], ncg);
198 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
204 srenew(time, teller+1);
207 set_pbc(&pbc, ePBC, box);
209 for (i = 0; (i < ncg); i++)
211 for (j = i+1; (j < ncg); j++)
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)
225 while (read_next_x(oenv, status, &t, x, box));
226 fprintf(stderr, "\n");
232 for (i = 0; (i < ncg); i++)
234 for (j = i+1; (j < ncg); j++)
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++)
242 fprintf(fp, "%10g %10g\n", time[k], cgdist[i][j][k]);
253 for (m = 0; (m < 3); m++)
255 out[m] = xvgropen(fn[m], title[m], "Time (ps)", "Distance (nm)", oenv);
259 for (i = 0; (i < ncg); i++)
262 for (j = i+1; (j < ncg); j++)
267 sprintf(buf, "%s:%s", cg[i].label, cg[j].label);
283 xvgr_legend(out[nnn], 1, (const char**)&buf, oenv);
287 if (output_env_get_xvg_format(oenv) == exvgXMGR)
289 fprintf(out[nnn], "@ legend string %d \"%s\"\n", nset[nnn], buf);
291 else if (output_env_get_xvg_format(oenv) == exvgXMGRACE)
293 fprintf(out[nnn], "@ s%d legend \"%s\"\n", nset[nnn], buf);
297 nWithin[i][j] = nnn+1;
301 for (k = 0; (k < teller); k++)
303 for (m = 0; (m < 3); m++)
305 fprintf(out[m], "%10g", time[k]);
308 for (i = 0; (i < ncg); i++)
310 for (j = i+1; (j < ncg); j++)
315 fprintf(out[nnn-1], " %10g", cgdist[i][j][k]);
319 for (m = 0; (m < 3); m++)
321 fprintf(out[m], "\n");
324 for (m = 0; (m < 3); m++)