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 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.
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.
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.
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.
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.
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.
47 #include "gromacs/fileio/filetypes.h"
48 #include "gromacs/gmxlib/nrnb.h"
49 #include "gromacs/math/utilities.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/mdtypes/forceoutput.h"
52 #include "gromacs/mdtypes/forcerec.h"
53 #include "gromacs/mdtypes/inputrec.h"
54 #include "gromacs/mdtypes/md_enums.h"
55 #include "gromacs/mdtypes/mdatom.h"
56 #include "gromacs/mdtypes/nblist.h"
57 #include "gromacs/tables/forcetable.h"
58 #include "gromacs/topology/topology.h"
59 #include "gromacs/utility/arrayref.h"
60 #include "gromacs/utility/cstringutil.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/smalloc.h"
64 void make_wall_tables(FILE* fplog,
67 const SimulationGroups* groups,
73 negp_pp = ir->opts.ngener - ir->nwall;
74 gmx::ArrayRef<const int> nm_ind = groups->groups[SimulationAtomGroupType::EnergyOutput];
78 fprintf(fplog, "Reading user tables for %d energy groups with %d walls\n", negp_pp, ir->nwall);
81 fr->wall_tab.resize(ir->nwall);
82 for (int w = 0; w < ir->nwall; w++)
84 fr->wall_tab[w].resize(negp_pp);
85 for (int egp = 0; egp < negp_pp; egp++)
87 /* If the energy group pair is excluded, we don't need a table */
88 if (!(fr->egp_flags[egp * ir->opts.ngener + negp_pp + w] & EGP_EXCL))
90 sprintf(buf, "%s", tabfn);
91 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1,
93 *groups->groupNames[nm_ind[egp]],
94 *groups->groupNames[nm_ind[negp_pp + w]],
96 fr->wall_tab[w][egp] = make_tables(fplog, fr->ic, buf, 0, GMX_MAKETABLES_FORCEUSER);
98 /* Since wall have no charge, we can compress the table */
99 for (int i = 0; i <= fr->wall_tab[w][egp]->n; i++)
101 for (int j = 0; j < 8; j++)
103 fr->wall_tab[w][egp]->data[8 * i + j] =
104 fr->wall_tab[w][egp]->data[12 * i + 4 + j];
112 [[noreturn]] static void wall_error(int a, gmx::ArrayRef<const gmx::RVec> x, real r)
115 "An atom is beyond the wall: coordinates %f %f %f, distance %f\n"
116 "You might want to use the mdp option wall_r_linpot",
123 static void tableForce(real r, const t_forcetable& tab, real Cd, real Cr, real* V, real* F)
125 const real tabscale = tab.scale;
126 const real* VFtab = tab.data.data();
128 real rt = r * tabscale;
129 int n0 = static_cast<int>(rt);
132 /* Beyond the table range, set V and F to zero */
139 real eps2 = eps * eps;
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;
155 Geps = VFtab[nnn + 2] * eps;
156 Heps2 = VFtab[nnn + 3] * eps2;
157 Fp = Ft + Geps + Heps2;
159 FF = Fp + Geps + 2.0 * Heps2;
160 real Vr = 12 * Cr * VV;
161 real Fr = 12 * Cr * FF;
163 *F = -(Fd + Fr) * tabscale;
167 real do_walls(const t_inputrec& ir,
168 const t_forcerec& fr,
171 gmx::ArrayRef<const gmx::RVec> x,
172 gmx::ForceWithVirial* forceWithVirial,
177 constexpr real sixth = 1.0 / 6.0;
178 constexpr real twelfth = 1.0 / 12.0;
181 real fac_d[2], fac_r[2];
182 const unsigned short* gid = md.cENER;
184 const int nwall = ir.nwall;
185 const int ngid = ir.opts.ngener;
186 const int ntype = fr.ntype;
187 const real* nbfp = fr.nbfp.data();
188 const int* egp_flags = fr.egp_flags;
190 for (int w = 0; w < nwall; w++)
192 ntw[w] = 2 * ntype * ir.wall_atomtype[w];
193 switch (ir.wall_type)
196 fac_d[w] = ir.wall_density[w] * M_PI / 6;
197 fac_r[w] = ir.wall_density[w] * M_PI / 45;
200 fac_d[w] = ir.wall_density[w] * M_PI / 2;
201 fac_r[w] = ir.wall_density[w] * M_PI / 5;
206 const real wall_z[2] = { 0, box[ZZ][ZZ] };
208 rvec* gmx_restrict f = as_rvec_array(forceWithVirial->force_.data());
212 for (int lam = 0; lam < (md.nPerturbed ? 2 : 1); lam++)
236 for (int i = 0; i < md.homenr; i++)
238 for (int w = 0; w < std::min(nwall, 2); w++)
240 /* The wall energy groups are always at the end of the list */
241 const int ggid = gid[i] * ngid + ngid - nwall + w;
242 const int at = type[i];
243 /* nbfp now includes the 6/12 derivative prefactors */
244 const real Cd = nbfp[ntw[w] + 2 * at] * sixth;
245 const real Cr = nbfp[ntw[w] + 2 * at + 1] * twelfth;
246 if (!((Cd == 0 && Cr == 0) || (egp_flags[ggid] & EGP_EXCL)))
255 r = wall_z[1] - x[i][ZZ];
257 if (r < ir.wall_r_linpot)
259 mr = ir.wall_r_linpot - r;
260 r = ir.wall_r_linpot;
272 real r1, r2, r4, Vd, Vr;
273 switch (ir.wall_type)
276 tableForce(r, *fr.wall_tab[w][gid[i]], Cd, Cr, &V, &F);
283 Vd = fac_d[w] * Cd * r2 * r1;
284 Vr = fac_r[w] * Cr * r4 * r4 * r1;
286 F = lamfac * (9 * Vr - 3 * Vd) * r1;
292 Vd = fac_d[w] * Cd * r4;
293 Vr = fac_r[w] * Cr * r4 * r4 * r2;
295 F = lamfac * (10 * Vr - 4 * Vd) * r1;
302 Vr = Cr * r4 * r4 * r4;
304 F = lamfac * (12 * Vr - 6 * Vd) * r1;
320 Vlj[ggid] += lamfac * V;
328 dvdlambda += (lam == 0 ? -1 : 1) * Vlambda;
331 inc_nrnb(nrnb, eNR_WALLS, md.homenr);
334 if (forceWithVirial->computeVirial_)
336 rvec virial = { 0, 0, static_cast<real>(-0.5 * sumRF) };
337 forceWithVirial->addVirialContribution(virial);