2 * This file is part of the GROMACS molecular simulation package.
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, 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.
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.
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.
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.
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.
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.
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/inputrec.h"
51 #include "gromacs/mdtypes/md_enums.h"
52 #include "gromacs/mdtypes/mdatom.h"
53 #include "gromacs/mdtypes/nblist.h"
54 #include "gromacs/tables/forcetable.h"
55 #include "gromacs/topology/topology.h"
56 #include "gromacs/utility/cstringutil.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/smalloc.h"
60 void make_wall_tables(FILE *fplog,
61 const t_inputrec *ir, const char *tabfn,
62 const gmx_groups_t *groups,
69 negp_pp = ir->opts.ngener - ir->nwall;
70 nm_ind = groups->grps[egcENER].nm_ind;
74 fprintf(fplog, "Reading user tables for %d energy groups with %d walls\n",
78 snew(fr->wall_tab, ir->nwall);
79 for (int w = 0; w < ir->nwall; w++)
81 snew(fr->wall_tab[w], negp_pp);
82 for (int egp = 0; egp < negp_pp; egp++)
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))
87 sprintf(buf, "%s", tabfn);
88 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
89 *groups->grpname[nm_ind[egp]],
90 *groups->grpname[nm_ind[negp_pp+w]],
92 fr->wall_tab[w][egp] = make_tables(fplog, fr->ic, buf, 0,
93 GMX_MAKETABLES_FORCEUSER);
95 /* Since wall have no charge, we can compress the table */
96 for (int i = 0; i <= fr->wall_tab[w][egp]->n; i++)
98 for (int j = 0; j < 8; j++)
100 fr->wall_tab[w][egp]->data[8*i+j] =
101 fr->wall_tab[w][egp]->data[12*i+4+j];
109 static void wall_error(int a, const rvec *x, real r)
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);
117 real do_walls(const t_inputrec *ir, t_forcerec *fr, matrix box, const t_mdatoms *md,
118 const rvec x[], rvec f[], real lambda, real Vlj[], t_nrnb *nrnb)
121 int ntw[2], at, ntype, ngid, ggid, *egp_flags, *type;
122 real *nbfp, lamfac, fac_d[2], fac_r[2], Cd, Cr;
123 real wall_z[2], r, mr, r1, r2, r4, Vd, Vr, V = 0, Fd, Fr, F = 0, dvdlambda;
126 real tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps2, Fp, VV, FF;
127 unsigned short *gid = md->cENER;
131 ngid = ir->opts.ngener;
134 egp_flags = fr->egp_flags;
136 for (int w = 0; w < nwall; w++)
138 ntw[w] = 2*ntype*ir->wall_atomtype[w];
139 switch (ir->wall_type)
142 fac_d[w] = ir->wall_density[w]*M_PI/6;
143 fac_r[w] = ir->wall_density[w]*M_PI/45;
146 fac_d[w] = ir->wall_density[w]*M_PI/2;
147 fac_r[w] = ir->wall_density[w]*M_PI/5;
154 wall_z[1] = box[ZZ][ZZ];
158 for (int lam = 0; lam < (md->nPerturbed ? 2 : 1); lam++)
180 for (int i = 0; i < md->homenr; i++)
182 for (int w = 0; w < std::min(nwall, 2); w++)
184 /* The wall energy groups are always at the end of the list */
185 ggid = gid[i]*ngid + ngid - nwall + w;
187 /* nbfp now includes the 6/12 derivative prefactors */
188 Cd = nbfp[ntw[w]+2*at]/6;
189 Cr = nbfp[ntw[w]+2*at+1]/12;
190 if (!((Cd == 0 && Cr == 0) || (egp_flags[ggid] & EGP_EXCL)))
198 r = wall_z[1] - x[i][ZZ];
200 if (r < ir->wall_r_linpot)
202 mr = ir->wall_r_linpot - r;
203 r = ir->wall_r_linpot;
209 switch (ir->wall_type)
216 tab = fr->wall_tab[w][gid[i]];
217 tabscale = tab->scale;
221 n0 = static_cast<int>(rt);
224 /* Beyond the table range, set V and F to zero */
236 Geps = VFtab[nnn+2]*eps;
237 Heps2 = VFtab[nnn+3]*eps2;
238 Fp = Ft + Geps + Heps2;
240 FF = Fp + Geps + 2.0*Heps2;
247 Geps = VFtab[nnn+2]*eps;
248 Heps2 = VFtab[nnn+3]*eps2;
249 Fp = Ft + Geps + Heps2;
251 FF = Fp + Geps + 2.0*Heps2;
255 F = -lamfac*(Fd + Fr)*tabscale;
266 Vd = fac_d[w]*Cd*r2*r1;
267 Vr = fac_r[w]*Cr*r4*r4*r1;
269 F = lamfac*(9*Vr - 3*Vd)*r1;
280 Vr = fac_r[w]*Cr*r4*r4*r2;
282 F = lamfac*(10*Vr - 4*Vd)*r1;
295 F = lamfac*(12*Vr - 6*Vd)*r1;
308 Vlj[ggid] += lamfac*V;
311 /* Because of the single sum virial calculation we need
312 * to add the full virial contribution of the walls.
313 * Since the force only has a z-component, there is only
314 * a contribution to the z component of the virial tensor.
315 * We could also determine the virial contribution directly,
316 * which would be cheaper here, but that would require extra
317 * communication for f_novirsum for with virtual sites
320 xf_z[XX] -= x[i][XX]*F;
321 xf_z[YY] -= x[i][YY]*F;
322 xf_z[ZZ] -= wall_z[w]*F;
328 dvdlambda += (lam == 0 ? -1 : 1)*Vlambda;
331 inc_nrnb(nrnb, eNR_WALLS, md->homenr);
334 for (int i = 0; i < DIM; i++)
336 fr->vir_wall_z[i] = -0.5*xf_z[i];