Split lines with many copyright years
[gromacs.git] / src / gromacs / mdlib / wall.cpp
blob434f8d471eec67dc7e69e90a6b999ce3eb62d540
1 /*
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, 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.
39 #include "gmxpre.h"
41 #include "wall.h"
43 #include <cstring>
45 #include <algorithm>
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/inputrec.h"
53 #include "gromacs/mdtypes/md_enums.h"
54 #include "gromacs/mdtypes/mdatom.h"
55 #include "gromacs/mdtypes/nblist.h"
56 #include "gromacs/tables/forcetable.h"
57 #include "gromacs/topology/topology.h"
58 #include "gromacs/utility/cstringutil.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/smalloc.h"
62 void make_wall_tables(FILE* fplog,
63 const t_inputrec* ir,
64 const char* tabfn,
65 const SimulationGroups* groups,
66 t_forcerec* fr)
68 int negp_pp;
69 char buf[STRLEN];
71 negp_pp = ir->opts.ngener - ir->nwall;
72 gmx::ArrayRef<const int> nm_ind = groups->groups[SimulationAtomGroupType::EnergyOutput];
74 if (fplog)
76 fprintf(fplog, "Reading user tables for %d energy groups with %d walls\n", negp_pp, ir->nwall);
79 snew(fr->wall_tab, ir->nwall);
80 for (int w = 0; w < ir->nwall; w++)
82 snew(fr->wall_tab[w], negp_pp);
83 for (int egp = 0; egp < negp_pp; egp++)
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))
88 sprintf(buf, "%s", tabfn);
89 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
90 *groups->groupNames[nm_ind[egp]], *groups->groupNames[nm_ind[negp_pp + w]],
91 ftp2ext(efXVG));
92 fr->wall_tab[w][egp] = make_tables(fplog, fr->ic, buf, 0, GMX_MAKETABLES_FORCEUSER);
94 /* Since wall have no charge, we can compress the table */
95 for (int i = 0; i <= fr->wall_tab[w][egp]->n; i++)
97 for (int j = 0; j < 8; j++)
99 fr->wall_tab[w][egp]->data[8 * i + j] =
100 fr->wall_tab[w][egp]->data[12 * i + 4 + j];
108 [[noreturn]] static void wall_error(int a, const rvec* x, real r)
110 gmx_fatal(FARGS,
111 "An atom is beyond the wall: coordinates %f %f %f, distance %f\n"
112 "You might want to use the mdp option wall_r_linpot",
113 x[a][XX], x[a][YY], x[a][ZZ], r);
116 static void tableForce(real r, const t_forcetable& tab, real Cd, real Cr, real* V, real* F)
118 const real tabscale = tab.scale;
119 const real* VFtab = tab.data;
121 real rt = r * tabscale;
122 int n0 = static_cast<int>(rt);
123 if (n0 >= tab.n)
125 /* Beyond the table range, set V and F to zero */
126 *V = 0;
127 *F = 0;
129 else
131 real eps = rt - n0;
132 real eps2 = eps * eps;
133 /* Dispersion */
134 int nnn = 8 * n0;
135 real Yt = VFtab[nnn];
136 real Ft = VFtab[nnn + 1];
137 real Geps = VFtab[nnn + 2] * eps;
138 real Heps2 = VFtab[nnn + 3] * eps2;
139 real Fp = Ft + Geps + Heps2;
140 real VV = Yt + Fp * eps;
141 real FF = Fp + Geps + 2.0 * Heps2;
142 real Vd = 6 * Cd * VV;
143 real Fd = 6 * Cd * FF;
144 /* Repulsion */
145 nnn = nnn + 4;
146 Yt = VFtab[nnn];
147 Ft = VFtab[nnn + 1];
148 Geps = VFtab[nnn + 2] * eps;
149 Heps2 = VFtab[nnn + 3] * eps2;
150 Fp = Ft + Geps + Heps2;
151 VV = Yt + Fp * eps;
152 FF = Fp + Geps + 2.0 * Heps2;
153 real Vr = 12 * Cr * VV;
154 real Fr = 12 * Cr * FF;
155 *V = Vd + Vr;
156 *F = -(Fd + Fr) * tabscale;
160 real do_walls(const t_inputrec& ir,
161 const t_forcerec& fr,
162 const matrix box,
163 const t_mdatoms& md,
164 const rvec* x,
165 gmx::ForceWithVirial* forceWithVirial,
166 real lambda,
167 real Vlj[],
168 t_nrnb* nrnb)
170 constexpr real sixth = 1.0 / 6.0;
171 constexpr real twelfth = 1.0 / 12.0;
173 int ntw[2];
174 real fac_d[2], fac_r[2];
175 const unsigned short* gid = md.cENER;
177 const int nwall = ir.nwall;
178 const int ngid = ir.opts.ngener;
179 const int ntype = fr.ntype;
180 const real* nbfp = fr.nbfp.data();
181 const int* egp_flags = fr.egp_flags;
183 for (int w = 0; w < nwall; w++)
185 ntw[w] = 2 * ntype * ir.wall_atomtype[w];
186 switch (ir.wall_type)
188 case ewt93:
189 fac_d[w] = ir.wall_density[w] * M_PI / 6;
190 fac_r[w] = ir.wall_density[w] * M_PI / 45;
191 break;
192 case ewt104:
193 fac_d[w] = ir.wall_density[w] * M_PI / 2;
194 fac_r[w] = ir.wall_density[w] * M_PI / 5;
195 break;
196 default: break;
199 const real wall_z[2] = { 0, box[ZZ][ZZ] };
201 rvec* gmx_restrict f = as_rvec_array(forceWithVirial->force_.data());
203 real dvdlambda = 0;
204 double sumRF = 0;
205 for (int lam = 0; lam < (md.nPerturbed ? 2 : 1); lam++)
207 real lamfac;
208 const int* type;
209 if (md.nPerturbed)
211 if (lam == 0)
213 lamfac = 1 - lambda;
214 type = md.typeA;
216 else
218 lamfac = lambda;
219 type = md.typeB;
222 else
224 lamfac = 1;
225 type = md.typeA;
228 real Vlambda = 0;
229 for (int i = 0; i < md.homenr; i++)
231 for (int w = 0; w < std::min(nwall, 2); w++)
233 /* The wall energy groups are always at the end of the list */
234 const int ggid = gid[i] * ngid + ngid - nwall + w;
235 const int at = type[i];
236 /* nbfp now includes the 6/12 derivative prefactors */
237 const real Cd = nbfp[ntw[w] + 2 * at] * sixth;
238 const real Cr = nbfp[ntw[w] + 2 * at + 1] * twelfth;
239 if (!((Cd == 0 && Cr == 0) || (egp_flags[ggid] & EGP_EXCL)))
241 real r, mr;
242 if (w == 0)
244 r = x[i][ZZ];
246 else
248 r = wall_z[1] - x[i][ZZ];
250 if (r < ir.wall_r_linpot)
252 mr = ir.wall_r_linpot - r;
253 r = ir.wall_r_linpot;
255 else
257 mr = 0;
259 if (r <= 0)
261 wall_error(i, x, r);
264 real V, F;
265 real r1, r2, r4, Vd, Vr;
266 switch (ir.wall_type)
268 case ewtTABLE:
269 tableForce(r, *fr.wall_tab[w][gid[i]], Cd, Cr, &V, &F);
270 F *= lamfac;
271 break;
272 case ewt93:
273 r1 = 1 / r;
274 r2 = r1 * r1;
275 r4 = r2 * r2;
276 Vd = fac_d[w] * Cd * r2 * r1;
277 Vr = fac_r[w] * Cr * r4 * r4 * r1;
278 V = Vr - Vd;
279 F = lamfac * (9 * Vr - 3 * Vd) * r1;
280 break;
281 case ewt104:
282 r1 = 1 / r;
283 r2 = r1 * r1;
284 r4 = r2 * r2;
285 Vd = fac_d[w] * Cd * r4;
286 Vr = fac_r[w] * Cr * r4 * r4 * r2;
287 V = Vr - Vd;
288 F = lamfac * (10 * Vr - 4 * Vd) * r1;
289 break;
290 case ewt126:
291 r1 = 1 / r;
292 r2 = r1 * r1;
293 r4 = r2 * r2;
294 Vd = Cd * r4 * r2;
295 Vr = Cr * r4 * r4 * r4;
296 V = Vr - Vd;
297 F = lamfac * (12 * Vr - 6 * Vd) * r1;
298 break;
299 default:
300 V = 0;
301 F = 0;
302 break;
304 if (mr > 0)
306 V += mr * F;
308 sumRF += r * F;
309 if (w == 1)
311 F = -F;
313 Vlj[ggid] += lamfac * V;
314 Vlambda += V;
315 f[i][ZZ] += F;
319 if (md.nPerturbed)
321 dvdlambda += (lam == 0 ? -1 : 1) * Vlambda;
324 inc_nrnb(nrnb, eNR_WALLS, md.homenr);
327 if (forceWithVirial->computeVirial_)
329 rvec virial = { 0, 0, static_cast<real>(-0.5 * sumRF) };
330 forceWithVirial->addVirialContribution(virial);
333 return dvdlambda;