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(FILE *log,
49 const t_forcerec *fr, t_graph *g,
50 const t_mdatoms *mdatoms, const t_blocka *excl,
51 rvec x[], rvec f[], rvec *fshift, const t_pbc *pbc,
52 real lambda, real *dvdlambda)
54 /* Calculate the reaction-field energy correction for this node:
55 * epsfac q_i q_j (k_rf r_ij^2 - c_rf)
56 * and force correction for all excluded pairs, including self pairs.
58 int top, i, j, j1, j2, k, ki;
59 double q2sumA, q2sumB, ener;
60 const real *chargeA, *chargeB;
61 real ek, ec, L1, qiA, qiB, qqA, qqB, qqL, v;
65 int start = mdatoms->start;
66 int end = mdatoms->homenr+start;
68 gmx_bool bMolPBC = fr->bMolPBC;
72 /* For test particle insertion we only correct for the test molecule */
73 start = mdatoms->nr - fr->n_tpi;
76 ek = fr->epsfac*fr->k_rf;
77 ec = fr->epsfac*fr->c_rf;
78 chargeA = mdatoms->chargeA;
79 chargeB = mdatoms->chargeB;
95 if (mdatoms->nChargePerturbed == 0)
97 for (i = start; i < niat; i++)
104 /* Do the exclusions */
106 j2 = excl->index[i+1];
107 for (j = j1; j < j2; j++)
112 qqA = qiA*chargeA[k];
117 rvec_sub(x[i], x[k], dx);
118 ivec_sub(SHIFT_IVEC(g, i), SHIFT_IVEC(g, k), dt);
123 ki = pbc_dx_aiuc(pbc, x[i], x[k], dx);
127 rvec_sub(x[i], x[k], dx);
129 ener += qqA*(ek*norm2(dx) - ec);
130 svmul(-2*qqA*ek, dx, df);
133 rvec_inc(fshift[ki], df);
134 rvec_dec(fshift[CENTRAL], df);
139 ener += -0.5*ec*q2sumA;
144 for (i = start; i < niat; i++)
153 /* Do the exclusions */
155 j2 = excl->index[i+1];
156 for (j = j1; j < j2; j++)
161 qqA = qiA*chargeA[k];
162 qqB = qiB*chargeB[k];
163 if (qqA != 0 || qqB != 0)
165 qqL = L1*qqA + lambda*qqB;
168 rvec_sub(x[i], x[k], dx);
169 ivec_sub(SHIFT_IVEC(g, i), SHIFT_IVEC(g, k), dt);
174 ki = pbc_dx_aiuc(pbc, x[i], x[k], dx);
178 rvec_sub(x[i], x[k], dx);
180 v = ek*norm2(dx) - ec;
182 svmul(-2*qqL*ek, dx, df);
185 rvec_inc(fshift[ki], df);
186 rvec_dec(fshift[CENTRAL], df);
187 *dvdlambda += (qqB - qqA)*v;
192 ener += -0.5*ec*(L1*q2sumA + lambda*q2sumB);
193 *dvdlambda += -0.5*ec*(q2sumB - q2sumA);
198 fprintf(debug, "RF exclusion energy: %g\n", ener);
204 void calc_rffac(FILE *fplog, int eel, real eps_r, real eps_rf, real Rc, real Temp,
205 real zsq, matrix box,
206 real *kappa, real *krf, real *crf)
208 /* Compute constants for Generalized reaction field */
209 real k1, k2, I, vol, rmin;
216 /* Consistency check */
219 gmx_fatal(FARGS, "Temperature is %f while using"
220 " Generalized Reaction Field\n", Temp);
222 /* Ionic strength (only needed for eelGRF */
224 *kappa = sqrt(2*I/(EPSILON0*eps_rf*BOLTZ*Temp));
232 /* eps == 0 signals infinite dielectric */
235 *krf = 1/(2*Rc*Rc*Rc);
240 k2 = eps_rf*sqr((real)(*kappa*Rc));
242 *krf = ((eps_rf - eps_r)*k1 + 0.5*k2)/((2*eps_rf + eps_r)*k1 + k2)/(Rc*Rc*Rc);
244 *crf = 1/Rc + *krf*Rc*Rc;
245 rmin = pow(*krf*2.0, -1.0/3.0);
251 please_cite(fplog, "Tironi95a");
252 fprintf(fplog, "%s:\n"
253 "epsRF = %10g, I = %10g, volume = %10g, kappa = %10g\n"
254 "rc = %10g, krf = %10g, crf = %10g, epsfac = %10g\n",
255 eel_names[eel], eps_rf, I, vol, *kappa, Rc, *krf, *crf,
260 fprintf(fplog, "%s:\n"
261 "epsRF = %g, rc = %g, krf = %g, crf = %g, epsfac = %g\n",
262 eel_names[eel], eps_rf, Rc, *krf, *crf, ONE_4PI_EPS0/eps_r);
265 "The electrostatics potential has its minimum at r = %g\n",
271 void init_generalized_rf(FILE *fplog,
272 const gmx_mtop_t *mtop, const t_inputrec *ir,
276 real q, zsq, nrdf, T;
277 const gmx_moltype_t *molt;
280 if (ir->efep != efepNO && fplog)
282 fprintf(fplog, "\nWARNING: the generalized reaction field constants are determined from topology A only\n\n");
285 for (mb = 0; mb < mtop->nmolblock; mb++)
287 molt = &mtop->moltype[mtop->molblock[mb].type];
289 for (i = 0; (i < cgs->nr); i++)
292 for (j = cgs->index[i]; (j < cgs->index[i+1]); j++)
294 q += molt->atoms.atom[j].q;
296 zsq += mtop->molblock[mb].nmol*q*q;
303 for (i = 0; (i < ir->opts.ngtc); i++)
305 nrdf += ir->opts.nrdf[i];
306 T += (ir->opts.nrdf[i] * ir->opts.ref_t[i]);
310 gmx_fatal(FARGS, "No degrees of freedom!");