Fixed perturbed wall interactions
[alexxy/gromacs.git] / src / mdlib / wall.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-2008, The GROMACS development team,
6  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012,2013, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
42
43 #include <math.h>
44 #include <string.h>
45 #include "maths.h"
46 #include "sysstuff.h"
47 #include "typedefs.h"
48 #include "macros.h"
49 #include "smalloc.h"
50 #include "force.h"
51 #include "filenm.h"
52 #include "nrnb.h"
53 #include "vec.h"
54
55 void make_wall_tables(FILE *fplog, const output_env_t oenv,
56                       const t_inputrec *ir, const char *tabfn,
57                       const gmx_groups_t *groups,
58                       t_forcerec *fr)
59 {
60     int           w, negp_pp, egp, i, j;
61     int          *nm_ind;
62     char          buf[STRLEN];
63     t_forcetable *tab;
64
65     negp_pp = ir->opts.ngener - ir->nwall;
66     nm_ind  = groups->grps[egcENER].nm_ind;
67
68     if (fplog)
69     {
70         fprintf(fplog, "Reading user tables for %d energy groups with %d walls\n",
71                 negp_pp, ir->nwall);
72     }
73
74     snew(fr->wall_tab, ir->nwall);
75     for (w = 0; w < ir->nwall; w++)
76     {
77         snew(fr->wall_tab[w], negp_pp);
78         for (egp = 0; egp < negp_pp; egp++)
79         {
80             /* If the energy group pair is excluded, we don't need a table */
81             if (!(fr->egp_flags[egp*ir->opts.ngener+negp_pp+w] & EGP_EXCL))
82             {
83                 tab = &fr->wall_tab[w][egp];
84                 sprintf(buf, "%s", tabfn);
85                 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
86                         *groups->grpname[nm_ind[egp]],
87                         *groups->grpname[nm_ind[negp_pp+w]],
88                         ftp2ext(efXVG));
89                 *tab = make_tables(fplog, oenv, fr, FALSE, buf, 0, GMX_MAKETABLES_FORCEUSER);
90                 /* Since wall have no charge, we can compress the table */
91                 for (i = 0; i <= tab->n; i++)
92                 {
93                     for (j = 0; j < 8; j++)
94                     {
95                         tab->data[8*i+j] = tab->data[12*i+4+j];
96                     }
97                 }
98             }
99         }
100     }
101 }
102
103 static void wall_error(int a, rvec *x, real r)
104 {
105     gmx_fatal(FARGS,
106               "An atom is beyond the wall: coordinates %f %f %f, distance %f\n"
107               "You might want to use the mdp option wall_r_linpot",
108               x[a][XX], x[a][YY], x[a][ZZ], r);
109 }
110
111 real do_walls(t_inputrec *ir, t_forcerec *fr, matrix box, t_mdatoms *md,
112               rvec x[], rvec f[], real lambda, real Vlj[], t_nrnb *nrnb)
113 {
114     int             nwall, w, lam, i;
115     int             ntw[2], at, ntype, ngid, ggid, *egp_flags, *type;
116     real           *nbfp, lamfac, fac_d[2], fac_r[2], Cd, Cr, Vtot, Fwall[2];
117     real            wall_z[2], r, mr, r1, r2, r4, Vd, Vr, V = 0, Fd, Fr, F = 0, dvdlambda;
118     dvec            xf_z;
119     int             n0, nnn;
120     real            tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps, Heps2, Fp, VV, FF;
121     unsigned short *gid = md->cENER;
122     t_forcetable   *tab;
123
124     nwall     = ir->nwall;
125     ngid      = ir->opts.ngener;
126     ntype     = fr->ntype;
127     nbfp      = fr->nbfp;
128     egp_flags = fr->egp_flags;
129
130     for (w = 0; w < nwall; w++)
131     {
132         ntw[w] = 2*ntype*ir->wall_atomtype[w];
133         switch (ir->wall_type)
134         {
135             case ewt93:
136                 fac_d[w] = ir->wall_density[w]*M_PI/6;
137                 fac_r[w] = ir->wall_density[w]*M_PI/45;
138                 break;
139             case ewt104:
140                 fac_d[w] = ir->wall_density[w]*M_PI/2;
141                 fac_r[w] = ir->wall_density[w]*M_PI/5;
142                 break;
143             default:
144                 break;
145         }
146         Fwall[w] = 0;
147     }
148     wall_z[0] = 0;
149     wall_z[1] = box[ZZ][ZZ];
150
151     Vtot      = 0;
152     dvdlambda = 0;
153     clear_dvec(xf_z);
154     for (lam = 0; lam < (md->nPerturbed ? 2 : 1); lam++)
155     {
156         if (md->nPerturbed)
157         {
158             if (lam == 0)
159             {
160                 lamfac = 1 - lambda;
161                 type   = md->typeA;
162             }
163             else
164             {
165                 lamfac = lambda;
166                 type   = md->typeB;
167             }
168         }
169         else
170         {
171             lamfac = 1;
172             type   = md->typeA;
173         }
174         for (i = md->start; i < md->start+md->homenr; i++)
175         {
176             for (w = 0; w < nwall; w++)
177             {
178                 /* The wall energy groups are always at the end of the list */
179                 ggid = gid[i]*ngid + ngid - nwall + w;
180                 at   = type[i];
181                 /* nbfp now includes the 6.0/12.0 derivative prefactors */
182                 Cd = nbfp[ntw[w]+2*at]/6.0;
183                 Cr = nbfp[ntw[w]+2*at+1]/12.0;
184                 if (!((Cd == 0 && Cr == 0) || (egp_flags[ggid] & EGP_EXCL)))
185                 {
186                     if (w == 0)
187                     {
188                         r = x[i][ZZ];
189                     }
190                     else
191                     {
192                         r = wall_z[1] - x[i][ZZ];
193                     }
194                     if (r < ir->wall_r_linpot)
195                     {
196                         mr = ir->wall_r_linpot - r;
197                         r  = ir->wall_r_linpot;
198                     }
199                     else
200                     {
201                         mr = 0;
202                     }
203                     switch (ir->wall_type)
204                     {
205                         case ewtTABLE:
206                             if (r < 0)
207                             {
208                                 wall_error(i, x, r);
209                             }
210                             tab      = &(fr->wall_tab[w][gid[i]]);
211                             tabscale = tab->scale;
212                             VFtab    = tab->data;
213
214                             rt    = r*tabscale;
215                             n0    = rt;
216                             if (n0 >= tab->n)
217                             {
218                                 /* Beyond the table range, set V and F to zero */
219                                 V     = 0;
220                                 F     = 0;
221                             }
222                             else
223                             {
224                                 eps   = rt - n0;
225                                 eps2  = eps*eps;
226                                 /* Dispersion */
227                                 nnn   = 8*n0;
228                                 Yt    = VFtab[nnn];
229                                 Ft    = VFtab[nnn+1];
230                                 Geps  = VFtab[nnn+2]*eps;
231                                 Heps2 = VFtab[nnn+3]*eps2;
232                                 Fp    = Ft + Geps + Heps2;
233                                 VV    = Yt + Fp*eps;
234                                 FF    = Fp + Geps + 2.0*Heps2;
235                                 Vd    = Cd*VV;
236                                 Fd    = Cd*FF;
237                                 /* Repulsion */
238                                 nnn   = nnn + 4;
239                                 Yt    = VFtab[nnn];
240                                 Ft    = VFtab[nnn+1];
241                                 Geps  = VFtab[nnn+2]*eps;
242                                 Heps2 = VFtab[nnn+3]*eps2;
243                                 Fp    = Ft + Geps + Heps2;
244                                 VV    = Yt + Fp*eps;
245                                 FF    = Fp + Geps + 2.0*Heps2;
246                                 Vr    = Cr*VV;
247                                 Fr    = Cr*FF;
248                                 V     = Vd + Vr;
249                                 F     = -lamfac*(Fd + Fr)*tabscale;
250                             }
251                             break;
252                         case ewt93:
253                             if (r <= 0)
254                             {
255                                 wall_error(i, x, r);
256                             }
257                             r1 = 1/r;
258                             r2 = r1*r1;
259                             r4 = r2*r2;
260                             Vd = fac_d[w]*Cd*r2*r1;
261                             Vr = fac_r[w]*Cr*r4*r4*r1;
262                             V  = Vr - Vd;
263                             F  = lamfac*(9*Vr - 3*Vd)*r1;
264                             break;
265                         case ewt104:
266                             if (r <= 0)
267                             {
268                                 wall_error(i, x, r);
269                             }
270                             r1 = 1/r;
271                             r2 = r1*r1;
272                             r4 = r2*r2;
273                             Vd = fac_d[w]*Cd*r4;
274                             Vr = fac_r[w]*Cr*r4*r4*r2;
275                             V  = Vr - Vd;
276                             F  = lamfac*(10*Vr - 4*Vd)*r1;
277                             break;
278                         case ewt126:
279                             if (r <= 0)
280                             {
281                                 wall_error(i, x, r);
282                             }
283                             r1 = 1/r;
284                             r2 = r1*r1;
285                             r4 = r2*r2;
286                             Vd = Cd*r4*r2;
287                             Vr = Cr*r4*r4*r4;
288                             V  = Vr - Vd;
289                             F  = lamfac*(12*Vr - 6*Vd)*r1;
290                             break;
291                         default:
292                             break;
293                     }
294                     if (mr > 0)
295                     {
296                         V += mr*F;
297                     }
298                     if (w == 1)
299                     {
300                         F = -F;
301                     }
302                     Vlj[ggid] += lamfac*V;
303                     Vtot      += V;
304                     f[i][ZZ]  += F;
305                     /* Because of the single sum virial calculation we need
306                      * to add  the full virial contribution of the walls.
307                      * Since the force only has a z-component, there is only
308                      * a contribution to the z component of the virial tensor.
309                      * We could also determine the virial contribution directly,
310                      * which would be cheaper here, but that would require extra
311                      * communication for f_novirsum for with virtual sites
312                      * in parallel.
313                      */
314                     xf_z[XX]  -= x[i][XX]*F;
315                     xf_z[YY]  -= x[i][YY]*F;
316                     xf_z[ZZ]  -= wall_z[w]*F;
317                 }
318             }
319         }
320         if (md->nPerturbed)
321         {
322             dvdlambda += (lam == 0 ? -1 : 1)*Vtot;
323         }
324
325         inc_nrnb(nrnb, eNR_WALLS, md->homenr);
326     }
327
328     for (i = 0; i < DIM; i++)
329     {
330         fr->vir_wall_z[i] = -0.5*xf_z[i];
331     }
332
333     return dvdlambda;
334 }