1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2008, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
52 void make_wall_tables(FILE *fplog, const output_env_t oenv,
53 const t_inputrec *ir, const char *tabfn,
54 const gmx_groups_t *groups,
57 int w, negp_pp, egp, i, j;
62 negp_pp = ir->opts.ngener - ir->nwall;
63 nm_ind = groups->grps[egcENER].nm_ind;
67 fprintf(fplog, "Reading user tables for %d energy groups with %d walls\n",
71 snew(fr->wall_tab, ir->nwall);
72 for (w = 0; w < ir->nwall; w++)
74 snew(fr->wall_tab[w], negp_pp);
75 for (egp = 0; egp < negp_pp; egp++)
77 /* If the energy group pair is excluded, we don't need a table */
78 if (!(fr->egp_flags[egp*ir->opts.ngener+negp_pp+w] & EGP_EXCL))
80 tab = &fr->wall_tab[w][egp];
81 sprintf(buf, "%s", tabfn);
82 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
83 *groups->grpname[nm_ind[egp]],
84 *groups->grpname[nm_ind[negp_pp+w]],
86 *tab = make_tables(fplog, oenv, fr, FALSE, buf, 0, GMX_MAKETABLES_FORCEUSER);
87 /* Since wall have no charge, we can compress the table */
88 for (i = 0; i <= tab->n; i++)
90 for (j = 0; j < 8; j++)
92 tab->data[8*i+j] = tab->data[12*i+4+j];
100 static void wall_error(int a, rvec *x, real r)
103 "An atom is beyond the wall: coordinates %f %f %f, distance %f\n"
104 "You might want to use the mdp option wall_r_linpot",
105 x[a][XX], x[a][YY], x[a][ZZ], r);
108 real do_walls(t_inputrec *ir, t_forcerec *fr, matrix box, t_mdatoms *md,
109 rvec x[], rvec f[], real lambda, real Vlj[], t_nrnb *nrnb)
111 int nwall, w, lam, i;
112 int ntw[2], at, ntype, ngid, ggid, *egp_flags, *type;
113 real *nbfp, lamfac, fac_d[2], fac_r[2], Cd, Cr, Vtot, Fwall[2];
114 real wall_z[2], r, mr, r1, r2, r4, Vd, Vr, V = 0, Fd, Fr, F = 0, dvdlambda;
117 real tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps, Heps2, Fp, VV, FF;
118 unsigned short *gid = md->cENER;
122 ngid = ir->opts.ngener;
125 egp_flags = fr->egp_flags;
127 for (w = 0; w < nwall; w++)
129 ntw[w] = 2*ntype*ir->wall_atomtype[w];
130 switch (ir->wall_type)
133 fac_d[w] = ir->wall_density[w]*M_PI/6;
134 fac_r[w] = ir->wall_density[w]*M_PI/45;
137 fac_d[w] = ir->wall_density[w]*M_PI/2;
138 fac_r[w] = ir->wall_density[w]*M_PI/5;
146 wall_z[1] = box[ZZ][ZZ];
151 for (lam = 0; lam < (md->nPerturbed ? 2 : 1); lam++)
171 for (i = md->start; i < md->start+md->homenr; i++)
173 for (w = 0; w < nwall; w++)
175 /* The wall energy groups are always at the end of the list */
176 ggid = gid[i]*ngid + ngid - nwall + w;
178 /* nbfp now includes the 6.0/12.0 derivative prefactors */
179 Cd = nbfp[ntw[w]+2*at]/6.0;
180 Cr = nbfp[ntw[w]+2*at+1]/12.0;
181 if (!((Cd == 0 && Cr == 0) || (egp_flags[ggid] & EGP_EXCL)))
189 r = wall_z[1] - x[i][ZZ];
191 if (r < ir->wall_r_linpot)
193 mr = ir->wall_r_linpot - r;
194 r = ir->wall_r_linpot;
200 switch (ir->wall_type)
207 tab = &(fr->wall_tab[w][gid[i]]);
208 tabscale = tab->scale;
215 /* Beyond the table range, set V and F to zero */
227 Geps = VFtab[nnn+2]*eps;
228 Heps2 = VFtab[nnn+3]*eps2;
229 Fp = Ft + Geps + Heps2;
231 FF = Fp + Geps + 2.0*Heps2;
238 Geps = VFtab[nnn+2]*eps;
239 Heps2 = VFtab[nnn+3]*eps2;
240 Fp = Ft + Geps + Heps2;
242 FF = Fp + Geps + 2.0*Heps2;
246 F = -lamfac*(Fd + Fr)*tabscale;
257 Vd = fac_d[w]*Cd*r2*r1;
258 Vr = fac_r[w]*Cr*r4*r4*r1;
260 F = lamfac*(9*Vr - 3*Vd)*r1;
271 Vr = fac_r[w]*Cr*r4*r4*r2;
273 F = lamfac*(10*Vr - 4*Vd)*r1;
286 F = lamfac*(12*Vr - 6*Vd)*r1;
299 Vlj[ggid] += lamfac*V;
302 /* Because of the single sum virial calculation we need
303 * to add the full virial contribution of the walls.
304 * Since the force only has a z-component, there is only
305 * a contribution to the z component of the virial tensor.
306 * We could also determine the virial contribution directly,
307 * which would be cheaper here, but that would require extra
308 * communication for f_novirsum for with virtual sites
311 xf_z[XX] -= x[i][XX]*F;
312 xf_z[YY] -= x[i][YY]*F;
313 xf_z[ZZ] -= wall_z[w]*F;
319 dvdlambda += (lam == 0 ? -1 : 1)*Vtot;
322 inc_nrnb(nrnb, eNR_WALLS, md->homenr);
325 for (i = 0; i < DIM; i++)
327 fr->vir_wall_z[i] = -0.5*xf_z[i];