Move mshift.* to pbcutil/
[alexxy/gromacs.git] / src / gromacs / mdlib / rf_util.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include "typedefs.h"
42 #include "force.h"
43 #include "names.h"
44 #include "copyrite.h"
45 #include "pbc.h"
46
47 #include "gromacs/pbcutil/mshift.h"
48 #include "gromacs/math/units.h"
49 #include "gromacs/math/vec.h"
50
51 real RF_excl_correction(const t_forcerec *fr, t_graph *g,
52                         const t_mdatoms *mdatoms, const t_blocka *excl,
53                         rvec x[], rvec f[], rvec *fshift, const t_pbc *pbc,
54                         real lambda, real *dvdlambda)
55 {
56     /* Calculate the reaction-field energy correction for this node:
57      * epsfac q_i q_j (k_rf r_ij^2 - c_rf)
58      * and force correction for all excluded pairs, including self pairs.
59      */
60     int         top, i, j, j1, j2, k, ki;
61     double      q2sumA, q2sumB, ener;
62     const real *chargeA, *chargeB;
63     real        ek, ec, L1, qiA, qiB, qqA, qqB, qqL, v;
64     rvec        dx, df;
65     atom_id    *AA;
66     ivec        dt;
67     int         start = 0;
68     int         end   = mdatoms->homenr;
69     int         niat;
70     gmx_bool    bMolPBC = fr->bMolPBC;
71
72     if (fr->n_tpi)
73     {
74         /* For test particle insertion we only correct for the test molecule */
75         start = mdatoms->nr - fr->n_tpi;
76     }
77
78     ek      = fr->epsfac*fr->k_rf;
79     ec      = fr->epsfac*fr->c_rf;
80     chargeA = mdatoms->chargeA;
81     chargeB = mdatoms->chargeB;
82     AA      = excl->a;
83     ki      = CENTRAL;
84
85     if (fr->bDomDec)
86     {
87         niat = excl->nr;
88     }
89     else
90     {
91         niat = end;
92     }
93
94     q2sumA = 0;
95     q2sumB = 0;
96     ener   = 0;
97     if (mdatoms->nChargePerturbed == 0)
98     {
99         for (i = start; i < niat; i++)
100         {
101             qiA = chargeA[i];
102             if (i < end)
103             {
104                 q2sumA += qiA*qiA;
105             }
106             /* Do the exclusions */
107             j1  = excl->index[i];
108             j2  = excl->index[i+1];
109             for (j = j1; j < j2; j++)
110             {
111                 k = AA[j];
112                 if (k > i)
113                 {
114                     qqA = qiA*chargeA[k];
115                     if (qqA != 0)
116                     {
117                         if (g)
118                         {
119                             rvec_sub(x[i], x[k], dx);
120                             ivec_sub(SHIFT_IVEC(g, i), SHIFT_IVEC(g, k), dt);
121                             ki = IVEC2IS(dt);
122                         }
123                         else if (bMolPBC)
124                         {
125                             ki = pbc_dx_aiuc(pbc, x[i], x[k], dx);
126                         }
127                         else
128                         {
129                             rvec_sub(x[i], x[k], dx);
130                         }
131                         ener += qqA*(ek*norm2(dx) - ec);
132                         svmul(-2*qqA*ek, dx, df);
133                         rvec_inc(f[i], df);
134                         rvec_dec(f[k], df);
135                         rvec_inc(fshift[ki], df);
136                         rvec_dec(fshift[CENTRAL], df);
137                     }
138                 }
139             }
140         }
141         ener += -0.5*ec*q2sumA;
142     }
143     else
144     {
145         L1 = 1.0 - lambda;
146         for (i = start; i < niat; i++)
147         {
148             qiA = chargeA[i];
149             qiB = chargeB[i];
150             if (i < end)
151             {
152                 q2sumA += qiA*qiA;
153                 q2sumB += qiB*qiB;
154             }
155             /* Do the exclusions */
156             j1  = excl->index[i];
157             j2  = excl->index[i+1];
158             for (j = j1; j < j2; j++)
159             {
160                 k = AA[j];
161                 if (k > i)
162                 {
163                     qqA = qiA*chargeA[k];
164                     qqB = qiB*chargeB[k];
165                     if (qqA != 0 || qqB != 0)
166                     {
167                         qqL = L1*qqA + lambda*qqB;
168                         if (g)
169                         {
170                             rvec_sub(x[i], x[k], dx);
171                             ivec_sub(SHIFT_IVEC(g, i), SHIFT_IVEC(g, k), dt);
172                             ki = IVEC2IS(dt);
173                         }
174                         else if (bMolPBC)
175                         {
176                             ki = pbc_dx_aiuc(pbc, x[i], x[k], dx);
177                         }
178                         else
179                         {
180                             rvec_sub(x[i], x[k], dx);
181                         }
182                         v     = ek*norm2(dx) - ec;
183                         ener += qqL*v;
184                         svmul(-2*qqL*ek, dx, df);
185                         rvec_inc(f[i], df);
186                         rvec_dec(f[k], df);
187                         rvec_inc(fshift[ki], df);
188                         rvec_dec(fshift[CENTRAL], df);
189                         *dvdlambda += (qqB - qqA)*v;
190                     }
191                 }
192             }
193         }
194         ener       += -0.5*ec*(L1*q2sumA + lambda*q2sumB);
195         *dvdlambda += -0.5*ec*(q2sumB - q2sumA);
196     }
197
198     if (debug)
199     {
200         fprintf(debug, "RF exclusion energy: %g\n", ener);
201     }
202
203     return ener;
204 }
205
206 void calc_rffac(FILE *fplog, int eel, real eps_r, real eps_rf, real Rc, real Temp,
207                 real zsq, matrix box,
208                 real *kappa, real *krf, real *crf)
209 {
210     /* Compute constants for Generalized reaction field */
211     real   k1, k2, I, vol, rmin;
212
213     if (EEL_RF(eel))
214     {
215         vol     = det(box);
216         if (eel == eelGRF)
217         {
218             /* Consistency check */
219             if (Temp <= 0.0)
220             {
221                 gmx_fatal(FARGS, "Temperature is %f while using"
222                           " Generalized Reaction Field\n", Temp);
223             }
224             /* Ionic strength (only needed for eelGRF */
225             I       = 0.5*zsq/vol;
226             *kappa  = sqrt(2*I/(EPSILON0*eps_rf*BOLTZ*Temp));
227         }
228         else
229         {
230             I      = 0;
231             *kappa = 0;
232         }
233
234         /* eps == 0 signals infinite dielectric */
235         if (eps_rf == 0)
236         {
237             *krf = 1/(2*Rc*Rc*Rc);
238         }
239         else
240         {
241             k1   = 1 + *kappa*Rc;
242             k2   = eps_rf*sqr((real)(*kappa*Rc));
243
244             *krf = ((eps_rf - eps_r)*k1 + 0.5*k2)/((2*eps_rf + eps_r)*k1 + k2)/(Rc*Rc*Rc);
245         }
246         *crf   = 1/Rc + *krf*Rc*Rc;
247         rmin   = pow(*krf*2.0, -1.0/3.0);
248
249         if (fplog)
250         {
251             if (eel == eelGRF)
252             {
253                 please_cite(fplog, "Tironi95a");
254                 fprintf(fplog, "%s:\n"
255                         "epsRF = %10g, I   = %10g, volume = %10g, kappa  = %10g\n"
256                         "rc    = %10g, krf = %10g, crf    = %10g, epsfac = %10g\n",
257                         eel_names[eel], eps_rf, I, vol, *kappa, Rc, *krf, *crf,
258                         ONE_4PI_EPS0/eps_r);
259             }
260             else
261             {
262                 fprintf(fplog, "%s:\n"
263                         "epsRF = %g, rc = %g, krf = %g, crf = %g, epsfac = %g\n",
264                         eel_names[eel], eps_rf, Rc, *krf, *crf, ONE_4PI_EPS0/eps_r);
265             }
266             fprintf(fplog,
267                     "The electrostatics potential has its minimum at r = %g\n",
268                     rmin);
269         }
270     }
271 }
272
273 void init_generalized_rf(FILE *fplog,
274                          const gmx_mtop_t *mtop, const t_inputrec *ir,
275                          t_forcerec *fr)
276 {
277     int                  mb, i, j;
278     real                 q, zsq, nrdf, T;
279     const gmx_moltype_t *molt;
280     const t_block       *cgs;
281
282     if (ir->efep != efepNO && fplog)
283     {
284         fprintf(fplog, "\nWARNING: the generalized reaction field constants are determined from topology A only\n\n");
285     }
286     zsq = 0.0;
287     for (mb = 0; mb < mtop->nmolblock; mb++)
288     {
289         molt = &mtop->moltype[mtop->molblock[mb].type];
290         cgs  = &molt->cgs;
291         for (i = 0; (i < cgs->nr); i++)
292         {
293             q = 0;
294             for (j = cgs->index[i]; (j < cgs->index[i+1]); j++)
295             {
296                 q += molt->atoms.atom[j].q;
297             }
298             zsq += mtop->molblock[mb].nmol*q*q;
299         }
300         fr->zsquare = zsq;
301     }
302
303     T    = 0.0;
304     nrdf = 0.0;
305     for (i = 0; (i < ir->opts.ngtc); i++)
306     {
307         nrdf += ir->opts.nrdf[i];
308         T    += (ir->opts.nrdf[i] * ir->opts.ref_t[i]);
309     }
310     if (nrdf == 0)
311     {
312         gmx_fatal(FARGS, "No degrees of freedom!");
313     }
314     fr->temp   = T/nrdf;
315 }