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, 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.
44 #include "gromacs/fileio/filenm.h"
45 #include "gromacs/legacyheaders/force.h"
46 #include "gromacs/legacyheaders/macros.h"
47 #include "gromacs/legacyheaders/nrnb.h"
48 #include "gromacs/legacyheaders/typedefs.h"
49 #include "gromacs/math/utilities.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/smalloc.h"
54 void make_wall_tables(FILE *fplog, const output_env_t oenv,
55 const t_inputrec *ir, const char *tabfn,
56 const gmx_groups_t *groups,
64 negp_pp = ir->opts.ngener - ir->nwall;
65 nm_ind = groups->grps[egcENER].nm_ind;
69 fprintf(fplog, "Reading user tables for %d energy groups with %d walls\n",
73 snew(fr->wall_tab, ir->nwall);
74 for (int w = 0; w < ir->nwall; w++)
76 snew(fr->wall_tab[w], negp_pp);
77 for (int egp = 0; egp < negp_pp; egp++)
79 /* If the energy group pair is excluded, we don't need a table */
80 if (!(fr->egp_flags[egp*ir->opts.ngener+negp_pp+w] & EGP_EXCL))
82 tab = &fr->wall_tab[w][egp];
83 sprintf(buf, "%s", tabfn);
84 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
85 *groups->grpname[nm_ind[egp]],
86 *groups->grpname[nm_ind[negp_pp+w]],
88 *tab = make_tables(fplog, oenv, fr, FALSE, buf, 0, GMX_MAKETABLES_FORCEUSER);
89 /* Since wall have no charge, we can compress the table */
90 for (int i = 0; i <= tab->n; i++)
92 for (int j = 0; j < 8; j++)
94 tab->data[8*i+j] = tab->data[12*i+4+j];
102 static void wall_error(int a, rvec *x, real r)
105 "An atom is beyond the wall: coordinates %f %f %f, distance %f\n"
106 "You might want to use the mdp option wall_r_linpot",
107 x[a][XX], x[a][YY], x[a][ZZ], r);
110 real do_walls(t_inputrec *ir, t_forcerec *fr, matrix box, t_mdatoms *md,
111 rvec x[], rvec f[], real lambda, real Vlj[], t_nrnb *nrnb)
114 int ntw[2], at, ntype, ngid, ggid, *egp_flags, *type;
115 real *nbfp, lamfac, fac_d[2], fac_r[2], Cd, Cr, Vtot;
116 real wall_z[2], r, mr, r1, r2, r4, Vd, Vr, V = 0, Fd, Fr, F = 0, dvdlambda;
119 real tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps2, Fp, VV, FF;
120 unsigned short *gid = md->cENER;
124 ngid = ir->opts.ngener;
127 egp_flags = fr->egp_flags;
129 for (int w = 0; w < nwall; w++)
131 ntw[w] = 2*ntype*ir->wall_atomtype[w];
132 switch (ir->wall_type)
135 fac_d[w] = ir->wall_density[w]*M_PI/6;
136 fac_r[w] = ir->wall_density[w]*M_PI/45;
139 fac_d[w] = ir->wall_density[w]*M_PI/2;
140 fac_r[w] = ir->wall_density[w]*M_PI/5;
147 wall_z[1] = box[ZZ][ZZ];
152 for (int lam = 0; lam < (md->nPerturbed ? 2 : 1); lam++)
172 for (int i = 0; i < md->homenr; i++)
174 for (int w = 0; w < std::min(nwall, 2); w++)
176 /* The wall energy groups are always at the end of the list */
177 ggid = gid[i]*ngid + ngid - nwall + w;
179 /* nbfp now includes the 6.0/12.0 derivative prefactors */
180 Cd = nbfp[ntw[w]+2*at]/6.0;
181 Cr = nbfp[ntw[w]+2*at+1]/12.0;
182 if (!((Cd == 0 && Cr == 0) || (egp_flags[ggid] & EGP_EXCL)))
190 r = wall_z[1] - x[i][ZZ];
192 if (r < ir->wall_r_linpot)
194 mr = ir->wall_r_linpot - r;
195 r = ir->wall_r_linpot;
201 switch (ir->wall_type)
208 tab = &(fr->wall_tab[w][gid[i]]);
209 tabscale = tab->scale;
213 n0 = static_cast<int>(rt);
216 /* Beyond the table range, set V and F to zero */
228 Geps = VFtab[nnn+2]*eps;
229 Heps2 = VFtab[nnn+3]*eps2;
230 Fp = Ft + Geps + Heps2;
232 FF = Fp + Geps + 2.0*Heps2;
239 Geps = VFtab[nnn+2]*eps;
240 Heps2 = VFtab[nnn+3]*eps2;
241 Fp = Ft + Geps + Heps2;
243 FF = Fp + Geps + 2.0*Heps2;
247 F = -lamfac*(Fd + Fr)*tabscale;
258 Vd = fac_d[w]*Cd*r2*r1;
259 Vr = fac_r[w]*Cr*r4*r4*r1;
261 F = lamfac*(9*Vr - 3*Vd)*r1;
272 Vr = fac_r[w]*Cr*r4*r4*r2;
274 F = lamfac*(10*Vr - 4*Vd)*r1;
287 F = lamfac*(12*Vr - 6*Vd)*r1;
300 Vlj[ggid] += lamfac*V;
303 /* Because of the single sum virial calculation we need
304 * to add the full virial contribution of the walls.
305 * Since the force only has a z-component, there is only
306 * a contribution to the z component of the virial tensor.
307 * We could also determine the virial contribution directly,
308 * which would be cheaper here, but that would require extra
309 * communication for f_novirsum for with virtual sites
312 xf_z[XX] -= x[i][XX]*F;
313 xf_z[YY] -= x[i][YY]*F;
314 xf_z[ZZ] -= wall_z[w]*F;
320 dvdlambda += (lam == 0 ? -1 : 1)*Vtot;
323 inc_nrnb(nrnb, eNR_WALLS, md->homenr);
326 for (int i = 0; i < DIM; i++)
328 fr->vir_wall_z[i] = -0.5*xf_z[i];