Remove all unnecessary HAVE_CONFIG_H
[gromacs.git] / src / gromacs / mdlib / rf_util.c
bloba6c4f11f26bd63a96ef0dc27c8b070ebbabd9e71
1 /*
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.
37 #include "config.h"
39 #include "typedefs.h"
40 #include "force.h"
41 #include "names.h"
42 #include "copyrite.h"
44 #include "gromacs/pbcutil/ishift.h"
45 #include "gromacs/pbcutil/mshift.h"
46 #include "gromacs/pbcutil/pbc.h"
47 #include "gromacs/math/units.h"
48 #include "gromacs/math/vec.h"
50 real RF_excl_correction(const t_forcerec *fr, t_graph *g,
51 const t_mdatoms *mdatoms, const t_blocka *excl,
52 rvec x[], rvec f[], rvec *fshift, const t_pbc *pbc,
53 real lambda, real *dvdlambda)
55 /* Calculate the reaction-field energy correction for this node:
56 * epsfac q_i q_j (k_rf r_ij^2 - c_rf)
57 * and force correction for all excluded pairs, including self pairs.
59 int top, i, j, j1, j2, k, ki;
60 double q2sumA, q2sumB, ener;
61 const real *chargeA, *chargeB;
62 real ek, ec, L1, qiA, qiB, qqA, qqB, qqL, v;
63 rvec dx, df;
64 atom_id *AA;
65 ivec dt;
66 int start = 0;
67 int end = mdatoms->homenr;
68 int niat;
69 gmx_bool bMolPBC = fr->bMolPBC;
71 if (fr->n_tpi)
73 /* For test particle insertion we only correct for the test molecule */
74 start = mdatoms->nr - fr->n_tpi;
77 ek = fr->epsfac*fr->k_rf;
78 ec = fr->epsfac*fr->c_rf;
79 chargeA = mdatoms->chargeA;
80 chargeB = mdatoms->chargeB;
81 AA = excl->a;
82 ki = CENTRAL;
84 if (fr->bDomDec)
86 niat = excl->nr;
88 else
90 niat = end;
93 q2sumA = 0;
94 q2sumB = 0;
95 ener = 0;
96 if (mdatoms->nChargePerturbed == 0)
98 for (i = start; i < niat; i++)
100 qiA = chargeA[i];
101 if (i < end)
103 q2sumA += qiA*qiA;
105 /* Do the exclusions */
106 j1 = excl->index[i];
107 j2 = excl->index[i+1];
108 for (j = j1; j < j2; j++)
110 k = AA[j];
111 if (k > i)
113 qqA = qiA*chargeA[k];
114 if (qqA != 0)
116 if (g)
118 rvec_sub(x[i], x[k], dx);
119 ivec_sub(SHIFT_IVEC(g, i), SHIFT_IVEC(g, k), dt);
120 ki = IVEC2IS(dt);
122 else if (bMolPBC)
124 ki = pbc_dx_aiuc(pbc, x[i], x[k], dx);
126 else
128 rvec_sub(x[i], x[k], dx);
130 ener += qqA*(ek*norm2(dx) - ec);
131 svmul(-2*qqA*ek, dx, df);
132 rvec_inc(f[i], df);
133 rvec_dec(f[k], df);
134 rvec_inc(fshift[ki], df);
135 rvec_dec(fshift[CENTRAL], df);
140 ener += -0.5*ec*q2sumA;
142 else
144 L1 = 1.0 - lambda;
145 for (i = start; i < niat; i++)
147 qiA = chargeA[i];
148 qiB = chargeB[i];
149 if (i < end)
151 q2sumA += qiA*qiA;
152 q2sumB += qiB*qiB;
154 /* Do the exclusions */
155 j1 = excl->index[i];
156 j2 = excl->index[i+1];
157 for (j = j1; j < j2; j++)
159 k = AA[j];
160 if (k > i)
162 qqA = qiA*chargeA[k];
163 qqB = qiB*chargeB[k];
164 if (qqA != 0 || qqB != 0)
166 qqL = L1*qqA + lambda*qqB;
167 if (g)
169 rvec_sub(x[i], x[k], dx);
170 ivec_sub(SHIFT_IVEC(g, i), SHIFT_IVEC(g, k), dt);
171 ki = IVEC2IS(dt);
173 else if (bMolPBC)
175 ki = pbc_dx_aiuc(pbc, x[i], x[k], dx);
177 else
179 rvec_sub(x[i], x[k], dx);
181 v = ek*norm2(dx) - ec;
182 ener += qqL*v;
183 svmul(-2*qqL*ek, dx, df);
184 rvec_inc(f[i], df);
185 rvec_dec(f[k], df);
186 rvec_inc(fshift[ki], df);
187 rvec_dec(fshift[CENTRAL], df);
188 *dvdlambda += (qqB - qqA)*v;
193 ener += -0.5*ec*(L1*q2sumA + lambda*q2sumB);
194 *dvdlambda += -0.5*ec*(q2sumB - q2sumA);
197 if (debug)
199 fprintf(debug, "RF exclusion energy: %g\n", ener);
202 return ener;
205 void calc_rffac(FILE *fplog, int eel, real eps_r, real eps_rf, real Rc, real Temp,
206 real zsq, matrix box,
207 real *kappa, real *krf, real *crf)
209 /* Compute constants for Generalized reaction field */
210 real k1, k2, I, vol, rmin;
212 if (EEL_RF(eel))
214 vol = det(box);
215 if (eel == eelGRF)
217 /* Consistency check */
218 if (Temp <= 0.0)
220 gmx_fatal(FARGS, "Temperature is %f while using"
221 " Generalized Reaction Field\n", Temp);
223 /* Ionic strength (only needed for eelGRF */
224 I = 0.5*zsq/vol;
225 *kappa = sqrt(2*I/(EPSILON0*eps_rf*BOLTZ*Temp));
227 else
229 I = 0;
230 *kappa = 0;
233 /* eps == 0 signals infinite dielectric */
234 if (eps_rf == 0)
236 *krf = 1/(2*Rc*Rc*Rc);
238 else
240 k1 = 1 + *kappa*Rc;
241 k2 = eps_rf*sqr((real)(*kappa*Rc));
243 *krf = ((eps_rf - eps_r)*k1 + 0.5*k2)/((2*eps_rf + eps_r)*k1 + k2)/(Rc*Rc*Rc);
245 *crf = 1/Rc + *krf*Rc*Rc;
246 rmin = pow(*krf*2.0, -1.0/3.0);
248 if (fplog)
250 if (eel == eelGRF)
252 please_cite(fplog, "Tironi95a");
253 fprintf(fplog, "%s:\n"
254 "epsRF = %10g, I = %10g, volume = %10g, kappa = %10g\n"
255 "rc = %10g, krf = %10g, crf = %10g, epsfac = %10g\n",
256 eel_names[eel], eps_rf, I, vol, *kappa, Rc, *krf, *crf,
257 ONE_4PI_EPS0/eps_r);
259 else
261 fprintf(fplog, "%s:\n"
262 "epsRF = %g, rc = %g, krf = %g, crf = %g, epsfac = %g\n",
263 eel_names[eel], eps_rf, Rc, *krf, *crf, ONE_4PI_EPS0/eps_r);
265 fprintf(fplog,
266 "The electrostatics potential has its minimum at r = %g\n",
267 rmin);
272 void init_generalized_rf(FILE *fplog,
273 const gmx_mtop_t *mtop, const t_inputrec *ir,
274 t_forcerec *fr)
276 int mb, i, j;
277 real q, zsq, nrdf, T;
278 const gmx_moltype_t *molt;
279 const t_block *cgs;
281 if (ir->efep != efepNO && fplog)
283 fprintf(fplog, "\nWARNING: the generalized reaction field constants are determined from topology A only\n\n");
285 zsq = 0.0;
286 for (mb = 0; mb < mtop->nmolblock; mb++)
288 molt = &mtop->moltype[mtop->molblock[mb].type];
289 cgs = &molt->cgs;
290 for (i = 0; (i < cgs->nr); i++)
292 q = 0;
293 for (j = cgs->index[i]; (j < cgs->index[i+1]); j++)
295 q += molt->atoms.atom[j].q;
297 zsq += mtop->molblock[mb].nmol*q*q;
299 fr->zsquare = zsq;
302 T = 0.0;
303 nrdf = 0.0;
304 for (i = 0; (i < ir->opts.ngtc); i++)
306 nrdf += ir->opts.nrdf[i];
307 T += (ir->opts.nrdf[i] * ir->opts.ref_t[i]);
309 if (nrdf == 0)
311 gmx_fatal(FARGS, "No degrees of freedom!");
313 fr->temp = T/nrdf;