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 * GROwing Monsters And Cloning Shrimps
48 real RF_excl_correction(const t_forcerec *fr, t_graph *g,
49 const t_mdatoms *mdatoms, const t_blocka *excl,
50 rvec x[], rvec f[], rvec *fshift, const t_pbc *pbc,
51 real lambda, real *dvdlambda)
53 /* Calculate the reaction-field energy correction for this node:
54 * epsfac q_i q_j (k_rf r_ij^2 - c_rf)
55 * and force correction for all excluded pairs, including self pairs.
57 int top, i, j, j1, j2, k, ki;
58 double q2sumA, q2sumB, ener;
59 const real *chargeA, *chargeB;
60 real ek, ec, L1, qiA, qiB, qqA, qqB, qqL, v;
64 int start = mdatoms->start;
65 int end = mdatoms->homenr+start;
67 gmx_bool bMolPBC = fr->bMolPBC;
71 /* For test particle insertion we only correct for the test molecule */
72 start = mdatoms->nr - fr->n_tpi;
75 ek = fr->epsfac*fr->k_rf;
76 ec = fr->epsfac*fr->c_rf;
77 chargeA = mdatoms->chargeA;
78 chargeB = mdatoms->chargeB;
94 if (mdatoms->nChargePerturbed == 0)
96 for (i = start; i < niat; i++)
103 /* Do the exclusions */
105 j2 = excl->index[i+1];
106 for (j = j1; j < j2; j++)
111 qqA = qiA*chargeA[k];
116 rvec_sub(x[i], x[k], dx);
117 ivec_sub(SHIFT_IVEC(g, i), SHIFT_IVEC(g, k), dt);
122 ki = pbc_dx_aiuc(pbc, x[i], x[k], dx);
126 rvec_sub(x[i], x[k], dx);
128 ener += qqA*(ek*norm2(dx) - ec);
129 svmul(-2*qqA*ek, dx, df);
132 rvec_inc(fshift[ki], df);
133 rvec_dec(fshift[CENTRAL], df);
138 ener += -0.5*ec*q2sumA;
143 for (i = start; i < niat; i++)
152 /* Do the exclusions */
154 j2 = excl->index[i+1];
155 for (j = j1; j < j2; j++)
160 qqA = qiA*chargeA[k];
161 qqB = qiB*chargeB[k];
162 if (qqA != 0 || qqB != 0)
164 qqL = L1*qqA + lambda*qqB;
167 rvec_sub(x[i], x[k], dx);
168 ivec_sub(SHIFT_IVEC(g, i), SHIFT_IVEC(g, k), dt);
173 ki = pbc_dx_aiuc(pbc, x[i], x[k], dx);
177 rvec_sub(x[i], x[k], dx);
179 v = ek*norm2(dx) - ec;
181 svmul(-2*qqL*ek, dx, df);
184 rvec_inc(fshift[ki], df);
185 rvec_dec(fshift[CENTRAL], df);
186 *dvdlambda += (qqB - qqA)*v;
191 ener += -0.5*ec*(L1*q2sumA + lambda*q2sumB);
192 *dvdlambda += -0.5*ec*(q2sumB - q2sumA);
197 fprintf(debug, "RF exclusion energy: %g\n", ener);
203 void calc_rffac(FILE *fplog, int eel, real eps_r, real eps_rf, real Rc, real Temp,
204 real zsq, matrix box,
205 real *kappa, real *krf, real *crf)
207 /* Compute constants for Generalized reaction field */
208 real k1, k2, I, vol, rmin;
215 /* Consistency check */
218 gmx_fatal(FARGS, "Temperature is %f while using"
219 " Generalized Reaction Field\n", Temp);
221 /* Ionic strength (only needed for eelGRF */
223 *kappa = sqrt(2*I/(EPSILON0*eps_rf*BOLTZ*Temp));
231 /* eps == 0 signals infinite dielectric */
234 *krf = 1/(2*Rc*Rc*Rc);
239 k2 = eps_rf*sqr((real)(*kappa*Rc));
241 *krf = ((eps_rf - eps_r)*k1 + 0.5*k2)/((2*eps_rf + eps_r)*k1 + k2)/(Rc*Rc*Rc);
243 *crf = 1/Rc + *krf*Rc*Rc;
244 rmin = pow(*krf*2.0, -1.0/3.0);
250 please_cite(fplog, "Tironi95a");
251 fprintf(fplog, "%s:\n"
252 "epsRF = %10g, I = %10g, volume = %10g, kappa = %10g\n"
253 "rc = %10g, krf = %10g, crf = %10g, epsfac = %10g\n",
254 eel_names[eel], eps_rf, I, vol, *kappa, Rc, *krf, *crf,
259 fprintf(fplog, "%s:\n"
260 "epsRF = %g, rc = %g, krf = %g, crf = %g, epsfac = %g\n",
261 eel_names[eel], eps_rf, Rc, *krf, *crf, ONE_4PI_EPS0/eps_r);
264 "The electrostatics potential has its minimum at r = %g\n",
270 void init_generalized_rf(FILE *fplog,
271 const gmx_mtop_t *mtop, const t_inputrec *ir,
275 real q, zsq, nrdf, T;
276 const gmx_moltype_t *molt;
279 if (ir->efep != efepNO && fplog)
281 fprintf(fplog, "\nWARNING: the generalized reaction field constants are determined from topology A only\n\n");
284 for (mb = 0; mb < mtop->nmolblock; mb++)
286 molt = &mtop->moltype[mtop->molblock[mb].type];
288 for (i = 0; (i < cgs->nr); i++)
291 for (j = cgs->index[i]; (j < cgs->index[i+1]); j++)
293 q += molt->atoms.atom[j].q;
295 zsq += mtop->molblock[mb].nmol*q*q;
302 for (i = 0; (i < ir->opts.ngtc); i++)
304 nrdf += ir->opts.nrdf[i];
305 T += (ir->opts.nrdf[i] * ir->opts.ref_t[i]);
309 gmx_fatal(FARGS, "No degrees of freedom!");