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