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.
45 #include "gromacs/utility/cstringutil.h"
46 #include "gromacs/utility/smalloc.h"
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/commandline/pargs.h"
49 #include "gromacs/pbcutil/pbc.h"
50 #include "gromacs/legacyheaders/force.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/futil.h"
53 #include "gromacs/math/utilities.h"
54 #include "gromacs/legacyheaders/macros.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/fileio/tpxio.h"
57 #include "gromacs/legacyheaders/mdrun.h"
58 #include "gromacs/random/random.h"
59 #include "gromacs/topology/index.h"
62 static void insert_ion(int nsa, int *nwater,
63 gmx_bool bSet[], int repl[], atom_id index[],
65 int sign, int q, const char *ionname,
67 real rmin, gmx_rng_t rng)
81 ei = nw*gmx_rng_uniform_real(rng);
84 while (bSet[ei] && (maxrand > 0));
87 gmx_fatal(FARGS, "No more replaceable solvent!");
90 fprintf(stderr, "Replacing solvent molecule %d (atom %d) with %s\n",
91 ei, index[nsa*ei], ionname);
93 /* Replace solvent molecule charges with ion charge */
97 atoms->atom[index[nsa*ei]].q = q;
98 for (i = 1; i < nsa; i++)
100 atoms->atom[index[nsa*ei+i]].q = 0;
103 /* Mark all solvent molecules within rmin as unavailable for substitution */
107 for (i = 0; (i < nw); i++)
111 pbc_dx(pbc, x[index[nsa*ei]], x[index[nsa*i]], dx);
112 if (iprod(dx, dx) < rmin2)
122 static char *aname(const char *mname)
127 str = gmx_strdup(mname);
129 while (i > 1 && (isdigit(str[i]) || (str[i] == '+') || (str[i] == '-')))
138 void sort_ions(int nsa, int nw, int repl[], atom_id index[],
139 t_atoms *atoms, rvec x[],
140 const char *p_name, const char *n_name)
142 int i, j, k, r, np, nn, starta, startr, npi, nni;
144 char **pptr = NULL, **nptr = NULL, **paptr = NULL, **naptr = NULL;
148 /* Put all the solvent in front and count the added ions */
152 for (i = 0; i < nw; i++)
157 for (k = 0; k < nsa; k++)
159 copy_rvec(x[index[nsa*i+k]], xt[j++]);
174 /* Put the positive and negative ions at the end */
175 starta = index[nsa*(nw - np - nn)];
176 startr = atoms->atom[starta].resind;
181 pptr[0] = gmx_strdup(p_name);
183 paptr[0] = aname(p_name);
188 nptr[0] = gmx_strdup(n_name);
190 naptr[0] = aname(n_name);
194 for (i = 0; i < nw; i++)
201 copy_rvec(x[index[nsa*i]], xt[j]);
202 atoms->atomname[j] = paptr;
203 atoms->atom[j].resind = k;
204 atoms->resinfo[k].name = pptr;
211 copy_rvec(x[index[nsa*i]], xt[j]);
212 atoms->atomname[j] = naptr;
213 atoms->atom[j].resind = k;
214 atoms->resinfo[k].name = nptr;
218 for (i = index[nsa*nw-1]+1; i < atoms->nr; i++)
220 j = i-(nsa-1)*(np+nn);
221 atoms->atomname[j] = atoms->atomname[i];
222 atoms->atom[j] = atoms->atom[i];
223 copy_rvec(x[i], xt[j]);
225 atoms->nr -= (nsa-1)*(np+nn);
227 /* Copy the new positions back */
228 for (i = index[0]; i < atoms->nr; i++)
230 copy_rvec(xt[i], x[i]);
236 static void update_topol(const char *topinout, int p_num, int n_num,
237 const char *p_name, const char *n_name, char *grpname)
239 #define TEMP_FILENM "temp.top"
241 char buf[STRLEN], buf2[STRLEN], *temp, **mol_line = NULL;
242 int line, i, nsol, nmol_line, sol_line, nsol_last;
245 printf("\nProcessing topology\n");
246 fpin = gmx_ffopen(topinout, "r");
247 fpout = gmx_ffopen(TEMP_FILENM, "w");
254 while (fgets(buf, STRLEN, fpin))
258 if ((temp = strchr(buf2, '\n')) != NULL)
266 if ((temp = strchr(buf2, '\n')) != NULL)
271 if (buf2[strlen(buf2)-1] == ']')
273 buf2[strlen(buf2)-1] = '\0';
276 bMolecules = (gmx_strcasecmp(buf2, "molecules") == 0);
278 fprintf(fpout, "%s", buf);
280 else if (!bMolecules)
282 fprintf(fpout, "%s", buf);
286 /* Check if this is a line with solvent molecules */
287 sscanf(buf, "%s", buf2);
288 if (gmx_strcasecmp(buf2, grpname) == 0)
290 sol_line = nmol_line;
291 sscanf(buf, "%*s %d", &nsol_last);
293 /* Store this molecules section line */
294 srenew(mol_line, nmol_line+1);
295 mol_line[nmol_line] = gmx_strdup(buf);
304 gmx_fatal(FARGS, "No line with moleculetype '%s' found the [ molecules ] section of file '%s'", grpname, topinout);
306 if (nsol_last < p_num+n_num)
309 gmx_fatal(FARGS, "The last entry for moleculetype '%s' in the [ molecules ] section of file '%s' has less solvent molecules (%d) than were replaced (%d)", grpname, topinout, nsol_last, p_num+n_num);
312 /* Print all the molecule entries */
313 for (i = 0; i < nmol_line; i++)
317 fprintf(fpout, "%s", mol_line[i]);
321 printf("Replacing %d solute molecules in topology file (%s) "
322 " by %d %s and %d %s ions.\n",
323 p_num+n_num, topinout, p_num, p_name, n_num, n_name);
324 nsol_last -= p_num + n_num;
327 fprintf(fpout, "%-10s %d\n", grpname, nsol_last);
331 fprintf(fpout, "%-15s %d\n", p_name, p_num);
335 fprintf(fpout, "%-15s %d\n", n_name, n_num);
340 /* use gmx_ffopen to generate backup of topinout */
341 fpout = gmx_ffopen(topinout, "w");
343 rename(TEMP_FILENM, topinout);
347 int gmx_genion(int argc, char *argv[])
349 const char *desc[] = {
350 "[THISMODULE] randomly replaces solvent molecules with monoatomic ions.",
351 "The group of solvent molecules should be continuous and all molecules",
352 "should have the same number of atoms.",
353 "The user should add the ion molecules to the topology file or use",
354 "the [TT]-p[tt] option to automatically modify the topology.[PAR]",
355 "The ion molecule type, residue and atom names in all force fields",
356 "are the capitalized element names without sign. This molecule name",
357 "should be given with [TT]-pname[tt] or [TT]-nname[tt], and the",
358 "[TT][molecules][tt] section of your topology updated accordingly,",
359 "either by hand or with [TT]-p[tt]. Do not use an atom name instead!",
360 "[PAR]Ions which can have multiple charge states get the multiplicity",
361 "added, without sign, for the uncommon states only.[PAR]",
362 "For larger ions, e.g. sulfate we recommended using [gmx-insert-molecules]."
364 const char *bugs[] = {
365 "If you specify a salt concentration existing ions are not taken into "
366 "account. In effect you therefore specify the amount of salt to be added.",
368 static int p_num = 0, n_num = 0, p_q = 1, n_q = -1;
369 static const char *p_name = "NA", *n_name = "CL";
370 static real rmin = 0.6, conc = 0;
371 static int seed = 1993;
372 static gmx_bool bNeutral = FALSE;
373 static t_pargs pa[] = {
374 { "-np", FALSE, etINT, {&p_num}, "Number of positive ions" },
375 { "-pname", FALSE, etSTR, {&p_name}, "Name of the positive ion" },
376 { "-pq", FALSE, etINT, {&p_q}, "Charge of the positive ion" },
377 { "-nn", FALSE, etINT, {&n_num}, "Number of negative ions" },
378 { "-nname", FALSE, etSTR, {&n_name}, "Name of the negative ion" },
379 { "-nq", FALSE, etINT, {&n_q}, "Charge of the negative ion" },
380 { "-rmin", FALSE, etREAL, {&rmin}, "Minimum distance between ions" },
381 { "-seed", FALSE, etINT, {&seed}, "Seed for random number generator" },
382 { "-conc", FALSE, etREAL, {&conc},
383 "Specify salt concentration (mol/liter). This will add sufficient ions to reach up to the specified concentration as computed from the volume of the cell in the input [TT].tpr[tt] file. Overrides the [TT]-np[tt] and [TT]-nn[tt] options." },
384 { "-neutral", FALSE, etBOOL, {&bNeutral}, "This option will add enough ions to neutralize the system. These ions are added on top of those specified with [TT]-np[tt]/[TT]-nn[tt] or [TT]-conc[tt]. "}
394 char *grpname, title[STRLEN];
396 int i, nw, nwa, nsa, nsalt, iqtot;
400 { efTPX, NULL, NULL, ffREAD },
401 { efNDX, NULL, NULL, ffOPTRD },
402 { efSTO, "-o", NULL, ffWRITE },
403 { efTOP, "-p", "topol", ffOPTRW }
405 #define NFILE asize(fnm)
407 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
408 asize(desc), desc, asize(bugs), bugs, &oenv))
413 /* Check input for something sensible */
414 if ((p_num < 0) || (n_num < 0))
416 gmx_fatal(FARGS, "Negative number of ions to add?");
419 if (conc > 0 && (p_num > 0 || n_num > 0))
421 fprintf(stderr, "WARNING: -conc specified, overriding -nn and -np.\n");
424 /* Read atom positions and charges */
425 read_tps_conf(ftp2fn(efTPX, NFILE, fnm), title, &top, &ePBC, &x, &v, box, FALSE);
428 /* Compute total charge */
430 for (i = 0; (i < atoms.nr); i++)
432 qtot += atoms.atom[i].q;
434 iqtot = gmx_nint(qtot);
439 /* Compute number of ions to be added */
441 nsalt = gmx_nint(conc*vol*AVOGADRO/1e24);
442 p_num = abs(nsalt*n_q);
443 n_num = abs(nsalt*p_q);
447 int qdelta = p_num*p_q + n_num*n_q + iqtot;
449 /* Check if the system is neutralizable
450 * is (qdelta == p_q*p_num + n_q*n_num) solvable for p_num and n_num? */
451 int gcd = gmx_greatest_common_divisor(n_q, p_q);
452 if ((qdelta % gcd) != 0)
454 gmx_fatal(FARGS, "Can't neutralize this system using -nq %d and"
455 " -pq %d.\n", n_q, p_q);
473 if ((p_num == 0) && (n_num == 0))
475 fprintf(stderr, "No ions to add, will just copy input configuration.\n");
479 printf("Will try to add %d %s ions and %d %s ions.\n",
480 p_num, p_name, n_num, n_name);
481 printf("Select a continuous group of solvent molecules\n");
482 get_index(&atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nwa, &index, &grpname);
483 for (i = 1; i < nwa; i++)
485 if (index[i] != index[i-1]+1)
487 gmx_fatal(FARGS, "The solvent group %s is not continuous: "
488 "index[%d]=%d, index[%d]=%d",
489 grpname, i, index[i-1]+1, i+1, index[i]+1);
493 while ((nsa < nwa) &&
494 (atoms.atom[index[nsa]].resind ==
495 atoms.atom[index[nsa-1]].resind))
501 gmx_fatal(FARGS, "Your solvent group size (%d) is not a multiple of %d",
505 fprintf(stderr, "Number of (%d-atomic) solvent molecules: %d\n", nsa, nw);
506 if (p_num+n_num > nw)
508 gmx_fatal(FARGS, "Not enough solvent for adding ions");
511 if (opt2bSet("-p", NFILE, fnm))
513 update_topol(opt2fn("-p", NFILE, fnm), p_num, n_num, p_name, n_name, grpname);
520 snew(atoms.pdbinfo, atoms.nr);
522 set_pbc(&pbc, ePBC, box);
526 rng = gmx_rng_init(gmx_rng_make_seed());
530 rng = gmx_rng_init(seed);
532 /* Now loop over the ions that have to be placed */
535 insert_ion(nsa, &nw, bSet, repl, index, x, &pbc,
536 1, p_q, p_name, &atoms, rmin, rng);
540 insert_ion(nsa, &nw, bSet, repl, index, x, &pbc,
541 -1, n_q, n_name, &atoms, rmin, rng);
543 gmx_rng_destroy(rng);
544 fprintf(stderr, "\n");
548 sort_ions(nsa, nw, repl, index, &atoms, x, p_name, n_name);
551 sfree(atoms.pdbinfo);
552 atoms.pdbinfo = NULL;
554 write_sto_conf(ftp2fn(efSTO, NFILE, fnm), *top.name, &atoms, x, NULL, ePBC, box);