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.
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
,
65 const SimulationGroups
* groups
,
71 negp_pp
= ir
->opts
.ngener
- ir
->nwall
;
72 gmx::ArrayRef
<const int> nm_ind
= groups
->groups
[SimulationAtomGroupType::EnergyOutput
];
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
]],
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
)
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
);
125 /* Beyond the table range, set V and F to zero */
132 real eps2
= eps
* eps
;
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
;
148 Geps
= VFtab
[nnn
+ 2] * eps
;
149 Heps2
= VFtab
[nnn
+ 3] * eps2
;
150 Fp
= Ft
+ Geps
+ Heps2
;
152 FF
= Fp
+ Geps
+ 2.0 * Heps2
;
153 real Vr
= 12 * Cr
* VV
;
154 real Fr
= 12 * Cr
* FF
;
156 *F
= -(Fd
+ Fr
) * tabscale
;
160 real
do_walls(const t_inputrec
& ir
,
161 const t_forcerec
& fr
,
165 gmx::ForceWithVirial
* forceWithVirial
,
170 constexpr real sixth
= 1.0 / 6.0;
171 constexpr real twelfth
= 1.0 / 12.0;
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
)
189 fac_d
[w
] = ir
.wall_density
[w
] * M_PI
/ 6;
190 fac_r
[w
] = ir
.wall_density
[w
] * M_PI
/ 45;
193 fac_d
[w
] = ir
.wall_density
[w
] * M_PI
/ 2;
194 fac_r
[w
] = ir
.wall_density
[w
] * M_PI
/ 5;
199 const real wall_z
[2] = { 0, box
[ZZ
][ZZ
] };
201 rvec
* gmx_restrict f
= as_rvec_array(forceWithVirial
->force_
.data());
205 for (int lam
= 0; lam
< (md
.nPerturbed
? 2 : 1); lam
++)
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
)))
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
;
265 real r1
, r2
, r4
, Vd
, Vr
;
266 switch (ir
.wall_type
)
269 tableForce(r
, *fr
.wall_tab
[w
][gid
[i
]], Cd
, Cr
, &V
, &F
);
276 Vd
= fac_d
[w
] * Cd
* r2
* r1
;
277 Vr
= fac_r
[w
] * Cr
* r4
* r4
* r1
;
279 F
= lamfac
* (9 * Vr
- 3 * Vd
) * r1
;
285 Vd
= fac_d
[w
] * Cd
* r4
;
286 Vr
= fac_r
[w
] * Cr
* r4
* r4
* r2
;
288 F
= lamfac
* (10 * Vr
- 4 * Vd
) * r1
;
295 Vr
= Cr
* r4
* r4
* r4
;
297 F
= lamfac
* (12 * Vr
- 6 * Vd
) * r1
;
313 Vlj
[ggid
] += lamfac
* V
;
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
);