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