Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / 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  * 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
38 #include "gmxpre.h"
39
40 #include <math.h>
41 #include <string.h>
42
43 #include "gromacs/fileio/filenm.h"
44 #include "gromacs/legacyheaders/force.h"
45 #include "gromacs/legacyheaders/macros.h"
46 #include "gromacs/legacyheaders/nrnb.h"
47 #include "gromacs/legacyheaders/typedefs.h"
48 #include "gromacs/math/utilities.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/utility/cstringutil.h"
51 #include "gromacs/utility/smalloc.h"
52
53 void make_wall_tables(FILE *fplog, const output_env_t oenv,
54                       const t_inputrec *ir, const char *tabfn,
55                       const gmx_groups_t *groups,
56                       t_forcerec *fr)
57 {
58     int           w, negp_pp, egp, i, j;
59     int          *nm_ind;
60     char          buf[STRLEN];
61     t_forcetable *tab;
62
63     negp_pp = ir->opts.ngener - ir->nwall;
64     nm_ind  = groups->grps[egcENER].nm_ind;
65
66     if (fplog)
67     {
68         fprintf(fplog, "Reading user tables for %d energy groups with %d walls\n",
69                 negp_pp, ir->nwall);
70     }
71
72     snew(fr->wall_tab, ir->nwall);
73     for (w = 0; w < ir->nwall; w++)
74     {
75         snew(fr->wall_tab[w], negp_pp);
76         for (egp = 0; egp < negp_pp; egp++)
77         {
78             /* If the energy group pair is excluded, we don't need a table */
79             if (!(fr->egp_flags[egp*ir->opts.ngener+negp_pp+w] & EGP_EXCL))
80             {
81                 tab = &fr->wall_tab[w][egp];
82                 sprintf(buf, "%s", tabfn);
83                 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
84                         *groups->grpname[nm_ind[egp]],
85                         *groups->grpname[nm_ind[negp_pp+w]],
86                         ftp2ext(efXVG));
87                 *tab = make_tables(fplog, oenv, fr, FALSE, buf, 0, GMX_MAKETABLES_FORCEUSER);
88                 /* Since wall have no charge, we can compress the table */
89                 for (i = 0; i <= tab->n; i++)
90                 {
91                     for (j = 0; j < 8; j++)
92                     {
93                         tab->data[8*i+j] = tab->data[12*i+4+j];
94                     }
95                 }
96             }
97         }
98     }
99 }
100
101 static void wall_error(int a, rvec *x, real r)
102 {
103     gmx_fatal(FARGS,
104               "An atom is beyond the wall: coordinates %f %f %f, distance %f\n"
105               "You might want to use the mdp option wall_r_linpot",
106               x[a][XX], x[a][YY], x[a][ZZ], r);
107 }
108
109 real do_walls(t_inputrec *ir, t_forcerec *fr, matrix box, t_mdatoms *md,
110               rvec x[], rvec f[], real lambda, real Vlj[], t_nrnb *nrnb)
111 {
112     int             nwall, w, lam, i;
113     int             ntw[2], at, ntype, ngid, ggid, *egp_flags, *type;
114     real           *nbfp, lamfac, fac_d[2], fac_r[2], Cd, Cr, Vtot, Fwall[2];
115     real            wall_z[2], r, mr, r1, r2, r4, Vd, Vr, V = 0, Fd, Fr, F = 0, dvdlambda;
116     dvec            xf_z;
117     int             n0, nnn;
118     real            tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps, Heps2, Fp, VV, FF;
119     unsigned short *gid = md->cENER;
120     t_forcetable   *tab;
121
122     nwall     = ir->nwall;
123     ngid      = ir->opts.ngener;
124     ntype     = fr->ntype;
125     nbfp      = fr->nbfp;
126     egp_flags = fr->egp_flags;
127
128     for (w = 0; w < nwall; w++)
129     {
130         ntw[w] = 2*ntype*ir->wall_atomtype[w];
131         switch (ir->wall_type)
132         {
133             case ewt93:
134                 fac_d[w] = ir->wall_density[w]*M_PI/6;
135                 fac_r[w] = ir->wall_density[w]*M_PI/45;
136                 break;
137             case ewt104:
138                 fac_d[w] = ir->wall_density[w]*M_PI/2;
139                 fac_r[w] = ir->wall_density[w]*M_PI/5;
140                 break;
141             default:
142                 break;
143         }
144         Fwall[w] = 0;
145     }
146     wall_z[0] = 0;
147     wall_z[1] = box[ZZ][ZZ];
148
149     Vtot      = 0;
150     dvdlambda = 0;
151     clear_dvec(xf_z);
152     for (lam = 0; lam < (md->nPerturbed ? 2 : 1); lam++)
153     {
154         if (md->nPerturbed)
155         {
156             if (lam == 0)
157             {
158                 lamfac = 1 - lambda;
159                 type   = md->typeA;
160             }
161             else
162             {
163                 lamfac = lambda;
164                 type   = md->typeB;
165             }
166         }
167         else
168         {
169             lamfac = 1;
170             type   = md->typeA;
171         }
172         for (i = 0; i < md->homenr; i++)
173         {
174             for (w = 0; w < nwall; w++)
175             {
176                 /* The wall energy groups are always at the end of the list */
177                 ggid = gid[i]*ngid + ngid - nwall + w;
178                 at   = type[i];
179                 /* nbfp now includes the 6.0/12.0 derivative prefactors */
180                 Cd = nbfp[ntw[w]+2*at]/6.0;
181                 Cr = nbfp[ntw[w]+2*at+1]/12.0;
182                 if (!((Cd == 0 && Cr == 0) || (egp_flags[ggid] & EGP_EXCL)))
183                 {
184                     if (w == 0)
185                     {
186                         r = x[i][ZZ];
187                     }
188                     else
189                     {
190                         r = wall_z[1] - x[i][ZZ];
191                     }
192                     if (r < ir->wall_r_linpot)
193                     {
194                         mr = ir->wall_r_linpot - r;
195                         r  = ir->wall_r_linpot;
196                     }
197                     else
198                     {
199                         mr = 0;
200                     }
201                     switch (ir->wall_type)
202                     {
203                         case ewtTABLE:
204                             if (r < 0)
205                             {
206                                 wall_error(i, x, r);
207                             }
208                             tab      = &(fr->wall_tab[w][gid[i]]);
209                             tabscale = tab->scale;
210                             VFtab    = tab->data;
211
212                             rt    = r*tabscale;
213                             n0    = rt;
214                             if (n0 >= tab->n)
215                             {
216                                 /* Beyond the table range, set V and F to zero */
217                                 V     = 0;
218                                 F     = 0;
219                             }
220                             else
221                             {
222                                 eps   = rt - n0;
223                                 eps2  = eps*eps;
224                                 /* Dispersion */
225                                 nnn   = 8*n0;
226                                 Yt    = VFtab[nnn];
227                                 Ft    = VFtab[nnn+1];
228                                 Geps  = VFtab[nnn+2]*eps;
229                                 Heps2 = VFtab[nnn+3]*eps2;
230                                 Fp    = Ft + Geps + Heps2;
231                                 VV    = Yt + Fp*eps;
232                                 FF    = Fp + Geps + 2.0*Heps2;
233                                 Vd    = Cd*VV;
234                                 Fd    = Cd*FF;
235                                 /* Repulsion */
236                                 nnn   = nnn + 4;
237                                 Yt    = VFtab[nnn];
238                                 Ft    = VFtab[nnn+1];
239                                 Geps  = VFtab[nnn+2]*eps;
240                                 Heps2 = VFtab[nnn+3]*eps2;
241                                 Fp    = Ft + Geps + Heps2;
242                                 VV    = Yt + Fp*eps;
243                                 FF    = Fp + Geps + 2.0*Heps2;
244                                 Vr    = Cr*VV;
245                                 Fr    = Cr*FF;
246                                 V     = Vd + Vr;
247                                 F     = -lamfac*(Fd + Fr)*tabscale;
248                             }
249                             break;
250                         case ewt93:
251                             if (r <= 0)
252                             {
253                                 wall_error(i, x, r);
254                             }
255                             r1 = 1/r;
256                             r2 = r1*r1;
257                             r4 = r2*r2;
258                             Vd = fac_d[w]*Cd*r2*r1;
259                             Vr = fac_r[w]*Cr*r4*r4*r1;
260                             V  = Vr - Vd;
261                             F  = lamfac*(9*Vr - 3*Vd)*r1;
262                             break;
263                         case ewt104:
264                             if (r <= 0)
265                             {
266                                 wall_error(i, x, r);
267                             }
268                             r1 = 1/r;
269                             r2 = r1*r1;
270                             r4 = r2*r2;
271                             Vd = fac_d[w]*Cd*r4;
272                             Vr = fac_r[w]*Cr*r4*r4*r2;
273                             V  = Vr - Vd;
274                             F  = lamfac*(10*Vr - 4*Vd)*r1;
275                             break;
276                         case ewt126:
277                             if (r <= 0)
278                             {
279                                 wall_error(i, x, r);
280                             }
281                             r1 = 1/r;
282                             r2 = r1*r1;
283                             r4 = r2*r2;
284                             Vd = Cd*r4*r2;
285                             Vr = Cr*r4*r4*r4;
286                             V  = Vr - Vd;
287                             F  = lamfac*(12*Vr - 6*Vd)*r1;
288                             break;
289                         default:
290                             break;
291                     }
292                     if (mr > 0)
293                     {
294                         V += mr*F;
295                     }
296                     if (w == 1)
297                     {
298                         F = -F;
299                     }
300                     Vlj[ggid] += lamfac*V;
301                     Vtot      += V;
302                     f[i][ZZ]  += F;
303                     /* Because of the single sum virial calculation we need
304                      * to add  the full virial contribution of the walls.
305                      * Since the force only has a z-component, there is only
306                      * a contribution to the z component of the virial tensor.
307                      * We could also determine the virial contribution directly,
308                      * which would be cheaper here, but that would require extra
309                      * communication for f_novirsum for with virtual sites
310                      * in parallel.
311                      */
312                     xf_z[XX]  -= x[i][XX]*F;
313                     xf_z[YY]  -= x[i][YY]*F;
314                     xf_z[ZZ]  -= wall_z[w]*F;
315                 }
316             }
317         }
318         if (md->nPerturbed)
319         {
320             dvdlambda += (lam == 0 ? -1 : 1)*Vtot;
321         }
322
323         inc_nrnb(nrnb, eNR_WALLS, md->homenr);
324     }
325
326     for (i = 0; i < DIM; i++)
327     {
328         fr->vir_wall_z[i] = -0.5*xf_z[i];
329     }
330
331     return dvdlambda;
332 }