3f086e6898a9367e5b7e32bedeb04b655294d79a
[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     int          *nm_ind;
68     char          buf[STRLEN];
69
70     negp_pp = ir->opts.ngener - ir->nwall;
71     nm_ind  = groups->groups[SimulationAtomGroupType::EnergyOutput].nm_ind;
72
73     if (fplog)
74     {
75         fprintf(fplog, "Reading user tables for %d energy groups with %d walls\n",
76                 negp_pp, ir->nwall);
77     }
78
79     snew(fr->wall_tab, ir->nwall);
80     for (int w = 0; w < ir->nwall; w++)
81     {
82         snew(fr->wall_tab[w], negp_pp);
83         for (int egp = 0; egp < negp_pp; egp++)
84         {
85             /* If the energy group pair is excluded, we don't need a table */
86             if (!(fr->egp_flags[egp*ir->opts.ngener+negp_pp+w] & EGP_EXCL))
87             {
88                 sprintf(buf, "%s", tabfn);
89                 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
90                         *groups->groupNames[nm_ind[egp]],
91                         *groups->groupNames[nm_ind[negp_pp+w]],
92                         ftp2ext(efXVG));
93                 fr->wall_tab[w][egp] = make_tables(fplog, fr->ic, buf, 0,
94                                                    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, const 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,
119                        const t_forcetable &tab,
120                        real                Cd,
121                        real                Cr,
122                        real               *V,
123                        real               *F)
124 {
125     const real  tabscale = tab.scale;
126     const real *VFtab    = tab.data;
127
128     real        rt = r*tabscale;
129     int         n0 = static_cast<int>(rt);
130     if (n0 >= tab.n)
131     {
132         /* Beyond the table range, set V and F to zero */
133         *V         = 0;
134         *F         = 0;
135     }
136     else
137     {
138         real eps   = rt - n0;
139         real eps2  = eps*eps;
140         /* Dispersion */
141         int  nnn   = 8*n0;
142         real Yt    = VFtab[nnn];
143         real Ft    = VFtab[nnn + 1];
144         real Geps  = VFtab[nnn + 2]*eps;
145         real Heps2 = VFtab[nnn + 3]*eps2;
146         real Fp    = Ft + Geps + Heps2;
147         real VV    = Yt + Fp*eps;
148         real FF    = Fp + Geps + 2.0*Heps2;
149         real Vd    = 6*Cd*VV;
150         real Fd    = 6*Cd*FF;
151         /* Repulsion */
152         nnn        = nnn + 4;
153         Yt         = VFtab[nnn];
154         Ft         = VFtab[nnn+1];
155         Geps       = VFtab[nnn+2]*eps;
156         Heps2      = VFtab[nnn+3]*eps2;
157         Fp         = Ft + Geps + Heps2;
158         VV         = Yt + Fp*eps;
159         FF         = Fp + Geps + 2.0*Heps2;
160         real Vr    = 12*Cr*VV;
161         real Fr    = 12*Cr*FF;
162         *V         = Vd + Vr;
163         *F         = -(Fd + Fr)*tabscale;
164     }
165 }
166
167 real do_walls(const t_inputrec &ir, const t_forcerec &fr,
168               const matrix box, const t_mdatoms &md,
169               const rvec *x, gmx::ForceWithVirial *forceWithVirial,
170               real lambda, real Vlj[], 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;
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:
199                 break;
200         }
201     }
202     const real          wall_z[2] = { 0, box[ZZ][ZZ] };
203
204     rvec * gmx_restrict f = as_rvec_array(forceWithVirial->force_.data());
205
206     real                dvdlambda = 0;
207     double              sumRF     = 0;
208     for (int lam = 0; lam < (md.nPerturbed ? 2 : 1); lam++)
209     {
210         real       lamfac;
211         const int *type;
212         if (md.nPerturbed)
213         {
214             if (lam == 0)
215             {
216                 lamfac = 1 - lambda;
217                 type   = md.typeA;
218             }
219             else
220             {
221                 lamfac = lambda;
222                 type   = md.typeB;
223             }
224         }
225         else
226         {
227             lamfac = 1;
228             type   = md.typeA;
229         }
230
231         real Vlambda = 0;
232         for (int i = 0; i < md.homenr; i++)
233         {
234             for (int w = 0; w < std::min(nwall, 2); w++)
235             {
236                 /* The wall energy groups are always at the end of the list */
237                 const int  ggid = gid[i]*ngid + ngid - nwall + w;
238                 const int  at   = type[i];
239                 /* nbfp now includes the 6/12 derivative prefactors */
240                 const real Cd = nbfp[ntw[w] + 2*at]*sixth;
241                 const real Cr = nbfp[ntw[w] + 2*at + 1]*twelfth;
242                 if (!((Cd == 0 && Cr == 0) || (egp_flags[ggid] & EGP_EXCL)))
243                 {
244                     real r, mr;
245                     if (w == 0)
246                     {
247                         r = x[i][ZZ];
248                     }
249                     else
250                     {
251                         r = wall_z[1] - x[i][ZZ];
252                     }
253                     if (r < ir.wall_r_linpot)
254                     {
255                         mr = ir.wall_r_linpot - r;
256                         r  = ir.wall_r_linpot;
257                     }
258                     else
259                     {
260                         mr = 0;
261                     }
262                     if (r <= 0)
263                     {
264                         wall_error(i, x, r);
265                     }
266
267                     real V, F;
268                     real r1, r2, r4, Vd, Vr;
269                     switch (ir.wall_type)
270                     {
271                         case ewtTABLE:
272                             tableForce(r, *fr.wall_tab[w][gid[i]], Cd, Cr, &V, &F);
273                             F *= lamfac;
274                             break;
275                         case ewt93:
276                             r1 = 1/r;
277                             r2 = r1*r1;
278                             r4 = r2*r2;
279                             Vd = fac_d[w]*Cd*r2*r1;
280                             Vr = fac_r[w]*Cr*r4*r4*r1;
281                             V  = Vr - Vd;
282                             F  = lamfac*(9*Vr - 3*Vd)*r1;
283                             break;
284                         case ewt104:
285                             r1 = 1/r;
286                             r2 = r1*r1;
287                             r4 = r2*r2;
288                             Vd = fac_d[w]*Cd*r4;
289                             Vr = fac_r[w]*Cr*r4*r4*r2;
290                             V  = Vr - Vd;
291                             F  = lamfac*(10*Vr - 4*Vd)*r1;
292                             break;
293                         case ewt126:
294                             r1 = 1/r;
295                             r2 = r1*r1;
296                             r4 = r2*r2;
297                             Vd = Cd*r4*r2;
298                             Vr = Cr*r4*r4*r4;
299                             V  = Vr - Vd;
300                             F  = lamfac*(12*Vr - 6*Vd)*r1;
301                             break;
302                         default:
303                             V  = 0;
304                             F  = 0;
305                             break;
306                     }
307                     if (mr > 0)
308                     {
309                         V     += mr*F;
310                     }
311                     sumRF     += r*F;
312                     if (w == 1)
313                     {
314                         F      = -F;
315                     }
316                     Vlj[ggid] += lamfac*V;
317                     Vlambda   += V;
318                     f[i][ZZ]  += F;
319                 }
320             }
321         }
322         if (md.nPerturbed)
323         {
324             dvdlambda += (lam == 0 ? -1 : 1)*Vlambda;
325         }
326
327         inc_nrnb(nrnb, eNR_WALLS, md.homenr);
328     }
329
330     if (forceWithVirial->computeVirial_)
331     {
332         rvec virial = { 0, 0, static_cast<real>(-0.5*sumRF) };
333         forceWithVirial->addVirialContribution(virial);
334     }
335
336     return dvdlambda;
337 }