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
45 #include "gromacs/fileio/filenm.h"
46 #include "gromacs/fileio/trxio.h"
48 #include "gromacs/fileio/futil.h"
49 #include "gmx_fatal.h"
62 static t_charge *mk_charge(t_atoms *atoms, t_block *cgs, int *nncg)
66 int i, j, ncg, resnr, anr;
69 /* Find the charged groups */
71 for (i = 0; (i < cgs->nr); i++)
74 for (j = cgs->index[i]; (j < cgs->index[i+1]); j++)
76 qq += atoms->atom[j].q;
78 if (fabs(qq) > 1.0e-5)
84 resnr = atoms->atom[anr].resind;
85 sprintf(buf, "%s%d-%d",
86 *(atoms->resinfo[resnr].name),
87 atoms->resinfo[resnr].nr,
89 cg[ncg].label = strdup(buf);
95 for (i = 0; (i < ncg); i++)
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++)
109 static real calc_dist(t_pbc *pbc, rvec x[], t_block *cgs, int icg, int jcg)
113 real d2, mindist2 = 1000;
115 for (i = cgs->index[icg]; (i < cgs->index[icg+1]); i++)
117 for (j = cgs->index[jcg]; (j < cgs->index[jcg+1]); j++)
119 pbc_dx(pbc, x[i], x[j], dx);
127 return sqrt(mindist2);
130 int gmx_saltbr(int argc, char *argv[])
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."
142 static gmx_bool bSep = FALSE;
143 static real truncate = 1000.0;
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)" }
151 { efTRX, "-f", NULL, ffREAD },
152 { efTPX, NULL, NULL, ffREAD },
154 #define NFILE asize(fnm)
157 static const char *title[3] = {
158 "Distance between positively charged groups",
159 "Distance between negatively charged groups",
160 "Distance between oppositely charged groups"
162 static const char *fn[3] = {
167 int nset[3] = {0, 0, 0};
173 int i, j, k, m, nnn, teller, ncg, n1, n2, n3, natoms;
174 real t, *time, qi, qj;
186 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_BE_NICE,
187 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
192 top = read_top(ftp2fn(efTPX, NFILE, fnm), &ePBC);
193 cg = mk_charge(&top->atoms, &(top->cgs), &ncg);
196 for (i = 0; (i < ncg); i++)
198 snew(cgdist[i], ncg);
199 snew(nWithin[i], ncg);
202 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
208 srenew(time, teller+1);
211 set_pbc(&pbc, ePBC, box);
213 for (i = 0; (i < ncg); i++)
215 for (j = i+1; (j < ncg); j++)
217 srenew(cgdist[i][j], teller+1);
218 cgdist[i][j][teller] =
219 calc_dist(&pbc, x, &(top->cgs), cg[i].cg, cg[j].cg);
220 if (cgdist[i][j][teller] < truncate)
229 while (read_next_x(oenv, status, &t, x, box));
230 fprintf(stderr, "\n");
236 for (i = 0; (i < ncg); i++)
238 for (j = i+1; (j < ncg); j++)
242 sprintf(buf, "sb-%s:%s.xvg", cg[i].label, cg[j].label);
243 fp = xvgropen(buf, buf, "Time (ps)", "Distance (nm)", oenv);
244 for (k = 0; (k < teller); k++)
246 fprintf(fp, "%10g %10g\n", time[k], cgdist[i][j][k]);
257 for (m = 0; (m < 3); m++)
259 out[m] = xvgropen(fn[m], title[m], "Time (ps)", "Distance (nm)", oenv);
263 for (i = 0; (i < ncg); i++)
266 for (j = i+1; (j < ncg); j++)
271 sprintf(buf, "%s:%s", cg[i].label, cg[j].label);
287 xvgr_legend(out[nnn], 1, (const char**)&buf, oenv);
291 if (output_env_get_xvg_format(oenv) == exvgXMGR)
293 fprintf(out[nnn], "@ legend string %d \"%s\"\n", nset[nnn], buf);
295 else if (output_env_get_xvg_format(oenv) == exvgXMGRACE)
297 fprintf(out[nnn], "@ s%d legend \"%s\"\n", nset[nnn], buf);
301 nWithin[i][j] = nnn+1;
305 for (k = 0; (k < teller); k++)
307 for (m = 0; (m < 3); m++)
309 fprintf(out[m], "%10g", time[k]);
312 for (i = 0; (i < ncg); i++)
314 for (j = i+1; (j < ncg); j++)
319 fprintf(out[nnn-1], " %10g", cgdist[i][j][k]);
323 for (m = 0; (m < 3); m++)
325 fprintf(out[m], "\n");
328 for (m = 0; (m < 3); m++)