Replace t_grps with std::vector
[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,2018,2019, 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 "wall.h"
41
42 #include <cstring>
43
44 #include <algorithm>
45
46 #include "gromacs/fileio/filetypes.h"
47 #include "gromacs/gmxlib/nrnb.h"
48 #include "gromacs/math/utilities.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/mdtypes/forceoutput.h"
51 #include "gromacs/mdtypes/inputrec.h"
52 #include "gromacs/mdtypes/md_enums.h"
53 #include "gromacs/mdtypes/mdatom.h"
54 #include "gromacs/mdtypes/nblist.h"
55 #include "gromacs/tables/forcetable.h"
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/utility/cstringutil.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/smalloc.h"
60
61 void make_wall_tables(FILE *fplog,
62                       const t_inputrec *ir, const char *tabfn,
63                       const SimulationGroups *groups,
64                       t_forcerec *fr)
65 {
66     int           negp_pp;
67     char          buf[STRLEN];
68
69     negp_pp = ir->opts.ngener - ir->nwall;
70     gmx::ArrayRef<const int> nm_ind  = groups->groups[SimulationAtomGroupType::EnergyOutput];
71
72     if (fplog)
73     {
74         fprintf(fplog, "Reading user tables for %d energy groups with %d walls\n",
75                 negp_pp, ir->nwall);
76     }
77
78     snew(fr->wall_tab, ir->nwall);
79     for (int w = 0; w < ir->nwall; w++)
80     {
81         snew(fr->wall_tab[w], negp_pp);
82         for (int egp = 0; egp < negp_pp; egp++)
83         {
84             /* If the energy group pair is excluded, we don't need a table */
85             if (!(fr->egp_flags[egp*ir->opts.ngener+negp_pp+w] & EGP_EXCL))
86             {
87                 sprintf(buf, "%s", tabfn);
88                 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
89                         *groups->groupNames[nm_ind[egp]],
90                         *groups->groupNames[nm_ind[negp_pp+w]],
91                         ftp2ext(efXVG));
92                 fr->wall_tab[w][egp] = make_tables(fplog, fr->ic, buf, 0,
93                                                    GMX_MAKETABLES_FORCEUSER);
94
95                 /* Since wall have no charge, we can compress the table */
96                 for (int i = 0; i <= fr->wall_tab[w][egp]->n; i++)
97                 {
98                     for (int j = 0; j < 8; j++)
99                     {
100                         fr->wall_tab[w][egp]->data[8*i+j] =
101                             fr->wall_tab[w][egp]->data[12*i+4+j];
102                     }
103                 }
104             }
105         }
106     }
107 }
108
109 [[ noreturn ]] static void wall_error(int a, const rvec *x, real r)
110 {
111     gmx_fatal(FARGS,
112               "An atom is beyond the wall: coordinates %f %f %f, distance %f\n"
113               "You might want to use the mdp option wall_r_linpot",
114               x[a][XX], x[a][YY], x[a][ZZ], r);
115 }
116
117 static void tableForce(real                r,
118                        const t_forcetable &tab,
119                        real                Cd,
120                        real                Cr,
121                        real               *V,
122                        real               *F)
123 {
124     const real  tabscale = tab.scale;
125     const real *VFtab    = tab.data;
126
127     real        rt = r*tabscale;
128     int         n0 = static_cast<int>(rt);
129     if (n0 >= tab.n)
130     {
131         /* Beyond the table range, set V and F to zero */
132         *V         = 0;
133         *F         = 0;
134     }
135     else
136     {
137         real eps   = rt - n0;
138         real eps2  = eps*eps;
139         /* Dispersion */
140         int  nnn   = 8*n0;
141         real Yt    = VFtab[nnn];
142         real Ft    = VFtab[nnn + 1];
143         real Geps  = VFtab[nnn + 2]*eps;
144         real Heps2 = VFtab[nnn + 3]*eps2;
145         real Fp    = Ft + Geps + Heps2;
146         real VV    = Yt + Fp*eps;
147         real FF    = Fp + Geps + 2.0*Heps2;
148         real Vd    = 6*Cd*VV;
149         real Fd    = 6*Cd*FF;
150         /* Repulsion */
151         nnn        = nnn + 4;
152         Yt         = VFtab[nnn];
153         Ft         = VFtab[nnn+1];
154         Geps       = VFtab[nnn+2]*eps;
155         Heps2      = VFtab[nnn+3]*eps2;
156         Fp         = Ft + Geps + Heps2;
157         VV         = Yt + Fp*eps;
158         FF         = Fp + Geps + 2.0*Heps2;
159         real Vr    = 12*Cr*VV;
160         real Fr    = 12*Cr*FF;
161         *V         = Vd + Vr;
162         *F         = -(Fd + Fr)*tabscale;
163     }
164 }
165
166 real do_walls(const t_inputrec &ir, const t_forcerec &fr,
167               const matrix box, const t_mdatoms &md,
168               const rvec *x, gmx::ForceWithVirial *forceWithVirial,
169               real lambda, real Vlj[], t_nrnb *nrnb)
170 {
171     constexpr real        sixth   = 1.0/6.0;
172     constexpr real        twelfth = 1.0/12.0;
173
174     int                   ntw[2];
175     real                  fac_d[2], fac_r[2];
176     const unsigned short *gid = md.cENER;
177
178     const int             nwall     = ir.nwall;
179     const int             ngid      = ir.opts.ngener;
180     const int             ntype     = fr.ntype;
181     const real           *nbfp      = fr.nbfp;
182     const int            *egp_flags = fr.egp_flags;
183
184     for (int w = 0; w < nwall; w++)
185     {
186         ntw[w] = 2*ntype*ir.wall_atomtype[w];
187         switch (ir.wall_type)
188         {
189             case ewt93:
190                 fac_d[w] = ir.wall_density[w]*M_PI/6;
191                 fac_r[w] = ir.wall_density[w]*M_PI/45;
192                 break;
193             case ewt104:
194                 fac_d[w] = ir.wall_density[w]*M_PI/2;
195                 fac_r[w] = ir.wall_density[w]*M_PI/5;
196                 break;
197             default:
198                 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 }