Merge remote-tracking branch 'origin/release-4-6' into HEAD
[alexxy/gromacs.git] / src / gromacs / mdlib / rf_util.c
1 /*
2  * 
3  *                This source code is part of
4  * 
5  *                 G   R   O   M   A   C   S
6  * 
7  *          GROningen MAchine for Chemical Simulations
8  * 
9  *                        VERSION 3.2.0
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.
14
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.
19  * 
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.
26  * 
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.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  * 
32  * And Hey:
33  * GROwing Monsters And Cloning Shrimps
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include "sysstuff.h"
40 #include "typedefs.h"
41 #include "force.h"
42 #include "names.h"
43 #include "vec.h"
44 #include "physics.h"
45 #include "copyrite.h"
46 #include "pbc.h"
47
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)
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 = mdatoms->start;
66   int    end   = mdatoms->homenr+start;
67   int    niat;
68   gmx_bool   bMolPBC = fr->bMolPBC;
69
70   if (fr->n_tpi)
71     /* For test particle insertion we only correct for the test molecule */
72     start = mdatoms->nr - fr->n_tpi;
73
74   ek = fr->epsfac*fr->k_rf;
75   ec = fr->epsfac*fr->c_rf;
76   chargeA = mdatoms->chargeA;
77   chargeB = mdatoms->chargeB;
78   AA = excl->a;
79   ki = CENTRAL;
80
81   if (fr->bDomDec)
82     niat = excl->nr;
83   else
84     niat = end;
85
86   q2sumA = 0;
87   q2sumB = 0;
88   ener = 0;
89   if (mdatoms->nChargePerturbed == 0) {
90     for(i=start; i<niat; i++) {
91       qiA = chargeA[i];
92       if (i < end)
93         q2sumA += qiA*qiA;
94       /* Do the exclusions */
95       j1  = excl->index[i];
96       j2  = excl->index[i+1];
97       for(j=j1; j<j2; j++) {
98         k = AA[j];
99         if (k > i) {
100           qqA = qiA*chargeA[k];
101           if (qqA != 0) {
102             if (g) {
103               rvec_sub(x[i],x[k],dx);
104               ivec_sub(SHIFT_IVEC(g,i),SHIFT_IVEC(g,k),dt);
105               ki=IVEC2IS(dt);
106             } else if (bMolPBC) {
107               ki = pbc_dx_aiuc(pbc,x[i],x[k],dx);
108             } else
109               rvec_sub(x[i],x[k],dx);
110             ener += qqA*(ek*norm2(dx) - ec);
111             svmul(-2*qqA*ek,dx,df);
112             rvec_inc(f[i],df);
113             rvec_dec(f[k],df);
114             rvec_inc(fshift[ki],df);
115             rvec_dec(fshift[CENTRAL],df);
116           }
117         }
118       }
119     }
120     ener += -0.5*ec*q2sumA;
121   } else {
122     L1 = 1.0 - lambda;
123     for(i=start; i<niat; i++) {
124       qiA = chargeA[i];
125       qiB = chargeB[i];
126       if (i < end) {
127         q2sumA += qiA*qiA;
128         q2sumB += qiB*qiB;
129       }
130       /* Do the exclusions */
131       j1  = excl->index[i];
132       j2  = excl->index[i+1];
133       for(j=j1; j<j2; j++) {
134         k = AA[j];
135         if (k > i) {
136           qqA = qiA*chargeA[k];
137           qqB = qiB*chargeB[k];
138           if (qqA != 0 || qqB != 0) {
139             qqL = L1*qqA + lambda*qqB;
140             if (g) {
141               rvec_sub(x[i],x[k],dx);
142               ivec_sub(SHIFT_IVEC(g,i),SHIFT_IVEC(g,k),dt);
143               ki=IVEC2IS(dt);
144             } else if (bMolPBC) {
145               ki = pbc_dx_aiuc(pbc,x[i],x[k],dx);
146             } else
147               rvec_sub(x[i],x[k],dx);
148             v = ek*norm2(dx) - ec;
149             ener += qqL*v;
150             svmul(-2*qqL*ek,dx,df);
151             rvec_inc(f[i],df);
152             rvec_dec(f[k],df);
153             rvec_inc(fshift[ki],df);
154             rvec_dec(fshift[CENTRAL],df);
155             *dvdlambda += (qqB - qqA)*v;
156           }
157         }
158       }
159     }
160     ener += -0.5*ec*(L1*q2sumA + lambda*q2sumB);
161     *dvdlambda += -0.5*ec*(q2sumB - q2sumA);
162   }
163   
164   if (debug)
165     fprintf(debug,"RF exclusion energy: %g\n",ener);
166   
167   return ener;
168 }
169
170 void calc_rffac(FILE *fplog,int eel,real eps_r,real eps_rf,real Rc,real Temp,
171                 real zsq,matrix box,
172                 real *kappa,real *krf,real *crf)
173 {
174   /* Compute constants for Generalized reaction field */
175   real   k1,k2,I,vol,rmin;
176   
177   if (EEL_RF(eel)) {
178     vol     = det(box);
179     if (eel == eelGRF) {
180       /* Consistency check */
181       if (Temp <= 0.0)
182         gmx_fatal(FARGS,"Temperature is %f while using"
183                     " Generalized Reaction Field\n",Temp);
184       /* Ionic strength (only needed for eelGRF */
185       I       = 0.5*zsq/vol;
186       *kappa  = sqrt(2*I/(EPSILON0*eps_rf*BOLTZ*Temp));
187     }
188     else {
189       I      = 0;
190       *kappa = 0;
191     }
192
193     /* eps == 0 signals infinite dielectric */
194     if (eps_rf == 0) {
195       *krf = 1/(2*Rc*Rc*Rc);
196     } else {
197       k1   = 1 + *kappa*Rc;
198       k2   = eps_rf*sqr((real)(*kappa*Rc));
199       
200       *krf = ((eps_rf - eps_r)*k1 + 0.5*k2)/((2*eps_rf + eps_r)*k1 + k2)/(Rc*Rc*Rc);
201     }
202     *crf   = 1/Rc + *krf*Rc*Rc;
203     rmin   = pow(*krf*2.0,-1.0/3.0);
204     
205     if (fplog) {
206       if (eel == eelGRF) {
207         please_cite(fplog,"Tironi95a");
208         fprintf(fplog,"%s:\n"
209                 "epsRF = %10g, I   = %10g, volume = %10g, kappa  = %10g\n"
210                 "rc    = %10g, krf = %10g, crf    = %10g, epsfac = %10g\n",
211                 eel_names[eel],eps_rf,I,vol,*kappa,Rc,*krf,*crf,
212                 ONE_4PI_EPS0/eps_r);
213       } else {
214         fprintf(fplog,"%s:\n"
215                 "epsRF = %g, rc = %g, krf = %g, crf = %g, epsfac = %g\n",
216                 eel_names[eel],eps_rf,Rc,*krf,*crf,ONE_4PI_EPS0/eps_r);
217       }
218       fprintf(fplog,
219               "The electrostatics potential has its minimum at r = %g\n",
220               rmin);
221     }
222   }
223 }
224
225 void init_generalized_rf(FILE *fplog,
226                          const gmx_mtop_t *mtop,const t_inputrec *ir,
227                          t_forcerec *fr)
228 {
229   int  mb,i,j;
230   real q,zsq,nrdf,T;
231   const gmx_moltype_t *molt;
232   const t_block *cgs;
233
234   if (ir->efep != efepNO && fplog) {
235     fprintf(fplog,"\nWARNING: the generalized reaction field constants are determined from topology A only\n\n");
236   }
237   zsq = 0.0;
238   for(mb=0; mb<mtop->nmolblock; mb++) {
239     molt = &mtop->moltype[mtop->molblock[mb].type];
240     cgs  = &molt->cgs;
241     for(i=0; (i<cgs->nr); i++) {
242       q = 0;
243       for(j=cgs->index[i]; (j<cgs->index[i+1]); j++) {
244         q += molt->atoms.atom[j].q;
245       }
246       zsq += mtop->molblock[mb].nmol*q*q;
247     }
248     fr->zsquare = zsq;
249   }
250   
251   T    = 0.0;
252   nrdf = 0.0;
253   for(i=0; (i<ir->opts.ngtc); i++) {
254     nrdf += ir->opts.nrdf[i];
255     T    += (ir->opts.nrdf[i] * ir->opts.ref_t[i]);
256   }
257   if (nrdf == 0) {
258     gmx_fatal(FARGS,"No degrees of freedom!");
259   }
260   fr->temp   = T/nrdf;
261 }