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,2015,2017,2018, 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/math/functions.h"
44 #include "gromacs/math/units.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/mdtypes/forcerec.h"
47 #include "gromacs/mdtypes/inputrec.h"
48 #include "gromacs/mdtypes/md_enums.h"
49 #include "gromacs/mdtypes/mdatom.h"
50 #include "gromacs/pbcutil/ishift.h"
51 #include "gromacs/pbcutil/mshift.h"
52 #include "gromacs/pbcutil/pbc.h"
53 #include "gromacs/topology/topology.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/pleasecite.h"
57 real RF_excl_correction(const t_forcerec *fr,
59 const t_mdatoms *mdatoms,
61 bool usingDomainDecomposition,
69 /* Calculate the reaction-field energy correction for this node:
70 * epsfac q_i q_j (k_rf r_ij^2 - c_rf)
71 * and force correction for all excluded pairs, including self pairs.
73 int i, j, j1, j2, k, ki;
74 double q2sumA, q2sumB, ener;
75 const real *chargeA, *chargeB;
76 real ek, ec, L1, qiA, qiB, qqA, qqB, qqL, v;
81 int end = mdatoms->homenr;
83 gmx_bool bMolPBC = fr->bMolPBC;
87 /* For test particle insertion we only correct for the test molecule */
88 start = mdatoms->nr - fr->n_tpi;
91 ek = fr->ic->epsfac*fr->ic->k_rf;
92 ec = fr->ic->epsfac*fr->ic->c_rf;
93 chargeA = mdatoms->chargeA;
94 chargeB = mdatoms->chargeB;
98 if (usingDomainDecomposition)
110 if (mdatoms->nChargePerturbed == 0)
112 for (i = start; i < niat; i++)
119 /* Do the exclusions */
121 j2 = excl->index[i+1];
122 for (j = j1; j < j2; j++)
127 qqA = qiA*chargeA[k];
132 rvec_sub(x[i], x[k], dx);
133 ivec_sub(SHIFT_IVEC(g, i), SHIFT_IVEC(g, k), dt);
138 ki = pbc_dx_aiuc(pbc, x[i], x[k], dx);
142 rvec_sub(x[i], x[k], dx);
144 ener += qqA*(ek*norm2(dx) - ec);
145 svmul(-2*qqA*ek, dx, df);
148 rvec_inc(fshift[ki], df);
149 rvec_dec(fshift[CENTRAL], df);
154 ener += -0.5*ec*q2sumA;
159 for (i = start; i < niat; i++)
168 /* Do the exclusions */
170 j2 = excl->index[i+1];
171 for (j = j1; j < j2; j++)
176 qqA = qiA*chargeA[k];
177 qqB = qiB*chargeB[k];
178 if (qqA != 0 || qqB != 0)
180 qqL = L1*qqA + lambda*qqB;
183 rvec_sub(x[i], x[k], dx);
184 ivec_sub(SHIFT_IVEC(g, i), SHIFT_IVEC(g, k), dt);
189 ki = pbc_dx_aiuc(pbc, x[i], x[k], dx);
193 rvec_sub(x[i], x[k], dx);
195 v = ek*norm2(dx) - ec;
197 svmul(-2*qqL*ek, dx, df);
200 rvec_inc(fshift[ki], df);
201 rvec_dec(fshift[CENTRAL], df);
202 *dvdlambda += (qqB - qqA)*v;
207 ener += -0.5*ec*(L1*q2sumA + lambda*q2sumB);
208 *dvdlambda += -0.5*ec*(q2sumB - q2sumA);
213 fprintf(debug, "RF exclusion energy: %g\n", ener);
219 void calc_rffac(FILE *fplog, int eel, real eps_r, real eps_rf, real Rc, real Temp,
220 real zsq, matrix box,
221 real *krf, real *crf)
223 /* Compute constants for Generalized reaction field */
224 real kappa, k1, k2, I, rmin;
231 /* Consistency check */
234 gmx_fatal(FARGS, "Temperature is %f while using"
235 " Generalized Reaction Field\n", Temp);
237 /* Ionic strength (only needed for eelGRF */
240 kappa = std::sqrt(2*I/(EPSILON0*eps_rf*BOLTZ*Temp));
248 /* eps == 0 signals infinite dielectric */
251 *krf = 1/(2*Rc*Rc*Rc);
256 k2 = eps_rf*gmx::square(static_cast<real>(kappa*Rc));
258 *krf = ((eps_rf - eps_r)*k1 + 0.5*k2)/((2*eps_rf + eps_r)*k1 + k2)/(Rc*Rc*Rc);
260 *crf = 1/Rc + *krf*Rc*Rc;
261 // Make sure we don't lose resolution in pow() by casting real arg to double
262 rmin = gmx::invcbrt(static_cast<double>(*krf*2.0));
268 please_cite(fplog, "Tironi95a");
269 fprintf(fplog, "%s:\n"
270 "epsRF = %10g, I = %10g, volume = %10g, kappa = %10g\n"
271 "rc = %10g, krf = %10g, crf = %10g, epsfac = %10g\n",
272 eel_names[eel], eps_rf, I, vol, kappa, Rc, *krf, *crf,
277 fprintf(fplog, "%s:\n"
278 "epsRF = %g, rc = %g, krf = %g, crf = %g, epsfac = %g\n",
279 eel_names[eel], eps_rf, Rc, *krf, *crf, ONE_4PI_EPS0/eps_r);
282 "The electrostatics potential has its minimum at r = %g\n",
288 void init_generalized_rf(FILE *fplog,
289 const gmx_mtop_t *mtop, const t_inputrec *ir,
293 real q, zsq, nrdf, T;
294 const gmx_moltype_t *molt;
297 if (ir->efep != efepNO && fplog)
299 fprintf(fplog, "\nWARNING: the generalized reaction field constants are determined from topology A only\n\n");
302 for (const gmx_molblock_t &molb : mtop->molblock)
304 molt = &mtop->moltype[molb.type];
306 for (i = 0; (i < cgs->nr); i++)
309 for (j = cgs->index[i]; (j < cgs->index[i+1]); j++)
311 q += molt->atoms.atom[j].q;
313 zsq += molb.nmol*q*q;
320 for (i = 0; (i < ir->opts.ngtc); i++)
322 nrdf += ir->opts.nrdf[i];
323 T += (ir->opts.nrdf[i] * ir->opts.ref_t[i]);
327 gmx_fatal(FARGS, "No degrees of freedom!");