Move listed and wall forces out of do_force_lowlevel()
[alexxy/gromacs.git] / src / gromacs / mdlib / wall.cpp
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,2015,2016,2017 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source 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 #include "gmxpre.h"
40
41 #include "wall.h"
42
43 #include <cstring>
44
45 #include <algorithm>
46
47 #include "gromacs/fileio/filetypes.h"
48 #include "gromacs/gmxlib/nrnb.h"
49 #include "gromacs/math/utilities.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/mdtypes/forceoutput.h"
52 #include "gromacs/mdtypes/forcerec.h"
53 #include "gromacs/mdtypes/inputrec.h"
54 #include "gromacs/mdtypes/md_enums.h"
55 #include "gromacs/mdtypes/mdatom.h"
56 #include "gromacs/mdtypes/nblist.h"
57 #include "gromacs/tables/forcetable.h"
58 #include "gromacs/topology/topology.h"
59 #include "gromacs/utility/arrayref.h"
60 #include "gromacs/utility/cstringutil.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/smalloc.h"
63
64 void make_wall_tables(FILE*                   fplog,
65                       const t_inputrec*       ir,
66                       const char*             tabfn,
67                       const SimulationGroups* groups,
68                       t_forcerec*             fr)
69 {
70     int  negp_pp;
71     char buf[STRLEN];
72
73     negp_pp                         = ir->opts.ngener - ir->nwall;
74     gmx::ArrayRef<const int> nm_ind = groups->groups[SimulationAtomGroupType::EnergyOutput];
75
76     if (fplog)
77     {
78         fprintf(fplog, "Reading user tables for %d energy groups with %d walls\n", negp_pp, ir->nwall);
79     }
80
81     snew(fr->wall_tab, ir->nwall);
82     for (int w = 0; w < ir->nwall; w++)
83     {
84         snew(fr->wall_tab[w], negp_pp);
85         for (int egp = 0; egp < negp_pp; egp++)
86         {
87             /* If the energy group pair is excluded, we don't need a table */
88             if (!(fr->egp_flags[egp * ir->opts.ngener + negp_pp + w] & EGP_EXCL))
89             {
90                 sprintf(buf, "%s", tabfn);
91                 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
92                         *groups->groupNames[nm_ind[egp]], *groups->groupNames[nm_ind[negp_pp + w]],
93                         ftp2ext(efXVG));
94                 fr->wall_tab[w][egp] = make_tables(fplog, fr->ic, buf, 0, GMX_MAKETABLES_FORCEUSER);
95
96                 /* Since wall have no charge, we can compress the table */
97                 for (int i = 0; i <= fr->wall_tab[w][egp]->n; i++)
98                 {
99                     for (int j = 0; j < 8; j++)
100                     {
101                         fr->wall_tab[w][egp]->data[8 * i + j] =
102                                 fr->wall_tab[w][egp]->data[12 * i + 4 + j];
103                     }
104                 }
105             }
106         }
107     }
108 }
109
110 [[noreturn]] static void wall_error(int a, gmx::ArrayRef<const gmx::RVec> x, real r)
111 {
112     gmx_fatal(FARGS,
113               "An atom is beyond the wall: coordinates %f %f %f, distance %f\n"
114               "You might want to use the mdp option wall_r_linpot",
115               x[a][XX], x[a][YY], x[a][ZZ], r);
116 }
117
118 static void tableForce(real r, const t_forcetable& tab, real Cd, real Cr, real* V, real* F)
119 {
120     const real  tabscale = tab.scale;
121     const real* VFtab    = tab.data.data();
122
123     real rt = r * tabscale;
124     int  n0 = static_cast<int>(rt);
125     if (n0 >= tab.n)
126     {
127         /* Beyond the table range, set V and F to zero */
128         *V = 0;
129         *F = 0;
130     }
131     else
132     {
133         real eps  = rt - n0;
134         real eps2 = eps * eps;
135         /* Dispersion */
136         int  nnn   = 8 * n0;
137         real Yt    = VFtab[nnn];
138         real Ft    = VFtab[nnn + 1];
139         real Geps  = VFtab[nnn + 2] * eps;
140         real Heps2 = VFtab[nnn + 3] * eps2;
141         real Fp    = Ft + Geps + Heps2;
142         real VV    = Yt + Fp * eps;
143         real FF    = Fp + Geps + 2.0 * Heps2;
144         real Vd    = 6 * Cd * VV;
145         real Fd    = 6 * Cd * FF;
146         /* Repulsion */
147         nnn     = nnn + 4;
148         Yt      = VFtab[nnn];
149         Ft      = VFtab[nnn + 1];
150         Geps    = VFtab[nnn + 2] * eps;
151         Heps2   = VFtab[nnn + 3] * eps2;
152         Fp      = Ft + Geps + Heps2;
153         VV      = Yt + Fp * eps;
154         FF      = Fp + Geps + 2.0 * Heps2;
155         real Vr = 12 * Cr * VV;
156         real Fr = 12 * Cr * FF;
157         *V      = Vd + Vr;
158         *F      = -(Fd + Fr) * tabscale;
159     }
160 }
161
162 real do_walls(const t_inputrec&              ir,
163               const t_forcerec&              fr,
164               const matrix                   box,
165               const t_mdatoms&               md,
166               gmx::ArrayRef<const gmx::RVec> x,
167               gmx::ForceWithVirial*          forceWithVirial,
168               real                           lambda,
169               real                           Vlj[],
170               t_nrnb*                        nrnb)
171 {
172     constexpr real sixth   = 1.0 / 6.0;
173     constexpr real twelfth = 1.0 / 12.0;
174
175     int                   ntw[2];
176     real                  fac_d[2], fac_r[2];
177     const unsigned short* gid = md.cENER;
178
179     const int   nwall     = ir.nwall;
180     const int   ngid      = ir.opts.ngener;
181     const int   ntype     = fr.ntype;
182     const real* nbfp      = fr.nbfp.data();
183     const int*  egp_flags = fr.egp_flags;
184
185     for (int w = 0; w < nwall; w++)
186     {
187         ntw[w] = 2 * ntype * ir.wall_atomtype[w];
188         switch (ir.wall_type)
189         {
190             case ewt93:
191                 fac_d[w] = ir.wall_density[w] * M_PI / 6;
192                 fac_r[w] = ir.wall_density[w] * M_PI / 45;
193                 break;
194             case ewt104:
195                 fac_d[w] = ir.wall_density[w] * M_PI / 2;
196                 fac_r[w] = ir.wall_density[w] * M_PI / 5;
197                 break;
198             default: break;
199         }
200     }
201     const real wall_z[2] = { 0, box[ZZ][ZZ] };
202
203     rvec* gmx_restrict f = as_rvec_array(forceWithVirial->force_.data());
204
205     real   dvdlambda = 0;
206     double sumRF     = 0;
207     for (int lam = 0; lam < (md.nPerturbed ? 2 : 1); lam++)
208     {
209         real       lamfac;
210         const int* type;
211         if (md.nPerturbed)
212         {
213             if (lam == 0)
214             {
215                 lamfac = 1 - lambda;
216                 type   = md.typeA;
217             }
218             else
219             {
220                 lamfac = lambda;
221                 type   = md.typeB;
222             }
223         }
224         else
225         {
226             lamfac = 1;
227             type   = md.typeA;
228         }
229
230         real Vlambda = 0;
231         for (int i = 0; i < md.homenr; i++)
232         {
233             for (int w = 0; w < std::min(nwall, 2); w++)
234             {
235                 /* The wall energy groups are always at the end of the list */
236                 const int ggid = gid[i] * ngid + ngid - nwall + w;
237                 const int at   = type[i];
238                 /* nbfp now includes the 6/12 derivative prefactors */
239                 const real Cd = nbfp[ntw[w] + 2 * at] * sixth;
240                 const real Cr = nbfp[ntw[w] + 2 * at + 1] * twelfth;
241                 if (!((Cd == 0 && Cr == 0) || (egp_flags[ggid] & EGP_EXCL)))
242                 {
243                     real r, mr;
244                     if (w == 0)
245                     {
246                         r = x[i][ZZ];
247                     }
248                     else
249                     {
250                         r = wall_z[1] - x[i][ZZ];
251                     }
252                     if (r < ir.wall_r_linpot)
253                     {
254                         mr = ir.wall_r_linpot - r;
255                         r  = ir.wall_r_linpot;
256                     }
257                     else
258                     {
259                         mr = 0;
260                     }
261                     if (r <= 0)
262                     {
263                         wall_error(i, x, r);
264                     }
265
266                     real V, F;
267                     real r1, r2, r4, Vd, Vr;
268                     switch (ir.wall_type)
269                     {
270                         case ewtTABLE:
271                             tableForce(r, *fr.wall_tab[w][gid[i]], Cd, Cr, &V, &F);
272                             F *= lamfac;
273                             break;
274                         case ewt93:
275                             r1 = 1 / r;
276                             r2 = r1 * r1;
277                             r4 = r2 * r2;
278                             Vd = fac_d[w] * Cd * r2 * r1;
279                             Vr = fac_r[w] * Cr * r4 * r4 * r1;
280                             V  = Vr - Vd;
281                             F  = lamfac * (9 * Vr - 3 * Vd) * r1;
282                             break;
283                         case ewt104:
284                             r1 = 1 / r;
285                             r2 = r1 * r1;
286                             r4 = r2 * r2;
287                             Vd = fac_d[w] * Cd * r4;
288                             Vr = fac_r[w] * Cr * r4 * r4 * r2;
289                             V  = Vr - Vd;
290                             F  = lamfac * (10 * Vr - 4 * Vd) * r1;
291                             break;
292                         case ewt126:
293                             r1 = 1 / r;
294                             r2 = r1 * r1;
295                             r4 = r2 * r2;
296                             Vd = Cd * r4 * r2;
297                             Vr = Cr * r4 * r4 * r4;
298                             V  = Vr - Vd;
299                             F  = lamfac * (12 * Vr - 6 * Vd) * r1;
300                             break;
301                         default:
302                             V = 0;
303                             F = 0;
304                             break;
305                     }
306                     if (mr > 0)
307                     {
308                         V += mr * F;
309                     }
310                     sumRF += r * F;
311                     if (w == 1)
312                     {
313                         F = -F;
314                     }
315                     Vlj[ggid] += lamfac * V;
316                     Vlambda += V;
317                     f[i][ZZ] += F;
318                 }
319             }
320         }
321         if (md.nPerturbed)
322         {
323             dvdlambda += (lam == 0 ? -1 : 1) * Vlambda;
324         }
325
326         inc_nrnb(nrnb, eNR_WALLS, md.homenr);
327     }
328
329     if (forceWithVirial->computeVirial_)
330     {
331         rvec virial = { 0, 0, static_cast<real>(-0.5 * sumRF) };
332         forceWithVirial->addVirialContribution(virial);
333     }
334
335     return dvdlambda;
336 }