Enable TPI with the Verlet cut-off scheme
[gromacs.git] / src / gromacs / mdlib / forcerec.cpp
blob0b478cf1106799866f2359268347a074fee7c634
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-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, 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.
37 #include "gmxpre.h"
39 #include "forcerec.h"
41 #include "config.h"
43 #include <cassert>
44 #include <cmath>
45 #include <cstdlib>
46 #include <cstring>
48 #include <algorithm>
49 #include <memory>
51 #include "gromacs/commandline/filenm.h"
52 #include "gromacs/domdec/domdec.h"
53 #include "gromacs/domdec/domdec_struct.h"
54 #include "gromacs/ewald/ewald.h"
55 #include "gromacs/ewald/ewald_utils.h"
56 #include "gromacs/fileio/filetypes.h"
57 #include "gromacs/gmxlib/network.h"
58 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
59 #include "gromacs/gpu_utils/gpu_utils.h"
60 #include "gromacs/hardware/hw_info.h"
61 #include "gromacs/listed_forces/gpubonded.h"
62 #include "gromacs/listed_forces/manage_threading.h"
63 #include "gromacs/listed_forces/pairs.h"
64 #include "gromacs/math/functions.h"
65 #include "gromacs/math/units.h"
66 #include "gromacs/math/utilities.h"
67 #include "gromacs/math/vec.h"
68 #include "gromacs/mdlib/dispersioncorrection.h"
69 #include "gromacs/mdlib/force.h"
70 #include "gromacs/mdlib/forcerec_threading.h"
71 #include "gromacs/mdlib/gmx_omp_nthreads.h"
72 #include "gromacs/mdlib/md_support.h"
73 #include "gromacs/mdlib/qmmm.h"
74 #include "gromacs/mdlib/rf_util.h"
75 #include "gromacs/mdlib/wall.h"
76 #include "gromacs/mdtypes/commrec.h"
77 #include "gromacs/mdtypes/fcdata.h"
78 #include "gromacs/mdtypes/group.h"
79 #include "gromacs/mdtypes/iforceprovider.h"
80 #include "gromacs/mdtypes/inputrec.h"
81 #include "gromacs/mdtypes/md_enums.h"
82 #include "gromacs/nbnxm/gpu_data_mgmt.h"
83 #include "gromacs/nbnxm/nbnxm.h"
84 #include "gromacs/nbnxm/nbnxm_geometry.h"
85 #include "gromacs/pbcutil/ishift.h"
86 #include "gromacs/pbcutil/pbc.h"
87 #include "gromacs/tables/forcetable.h"
88 #include "gromacs/topology/mtop_util.h"
89 #include "gromacs/trajectory/trajectoryframe.h"
90 #include "gromacs/utility/cstringutil.h"
91 #include "gromacs/utility/exceptions.h"
92 #include "gromacs/utility/fatalerror.h"
93 #include "gromacs/utility/gmxassert.h"
94 #include "gromacs/utility/logger.h"
95 #include "gromacs/utility/physicalnodecommunicator.h"
96 #include "gromacs/utility/pleasecite.h"
97 #include "gromacs/utility/smalloc.h"
98 #include "gromacs/utility/strconvert.h"
100 static real *mk_nbfp(const gmx_ffparams_t *idef, gmx_bool bBHAM)
102 real *nbfp;
103 int i, j, k, atnr;
105 atnr = idef->atnr;
106 if (bBHAM)
108 snew(nbfp, 3*atnr*atnr);
109 for (i = k = 0; (i < atnr); i++)
111 for (j = 0; (j < atnr); j++, k++)
113 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
114 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
115 /* nbfp now includes the 6.0 derivative prefactor */
116 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c*6.0;
120 else
122 snew(nbfp, 2*atnr*atnr);
123 for (i = k = 0; (i < atnr); i++)
125 for (j = 0; (j < atnr); j++, k++)
127 /* nbfp now includes the 6.0/12.0 derivative prefactors */
128 C6(nbfp, atnr, i, j) = idef->iparams[k].lj.c6*6.0;
129 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12*12.0;
134 return nbfp;
137 static real *make_ljpme_c6grid(const gmx_ffparams_t *idef, t_forcerec *fr)
139 int i, j, k, atnr;
140 real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
141 real *grid;
143 /* For LJ-PME simulations, we correct the energies with the reciprocal space
144 * inside of the cut-off. To do this the non-bonded kernels needs to have
145 * access to the C6-values used on the reciprocal grid in pme.c
148 atnr = idef->atnr;
149 snew(grid, 2*atnr*atnr);
150 for (i = k = 0; (i < atnr); i++)
152 for (j = 0; (j < atnr); j++, k++)
154 c6i = idef->iparams[i*(atnr+1)].lj.c6;
155 c12i = idef->iparams[i*(atnr+1)].lj.c12;
156 c6j = idef->iparams[j*(atnr+1)].lj.c6;
157 c12j = idef->iparams[j*(atnr+1)].lj.c12;
158 c6 = std::sqrt(c6i * c6j);
159 if (fr->ljpme_combination_rule == eljpmeLB
160 && !gmx_numzero(c6) && !gmx_numzero(c12i) && !gmx_numzero(c12j))
162 sigmai = gmx::sixthroot(c12i / c6i);
163 sigmaj = gmx::sixthroot(c12j / c6j);
164 epsi = c6i * c6i / c12i;
165 epsj = c6j * c6j / c12j;
166 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5*(sigmai+sigmaj));
168 /* Store the elements at the same relative positions as C6 in nbfp in order
169 * to simplify access in the kernels
171 grid[2*(atnr*i+j)] = c6*6.0;
174 return grid;
177 /* This routine sets fr->solvent_opt to the most common solvent in the
178 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
179 * the fr->solvent_type array with the correct type (or esolNO).
181 * Charge groups that fulfill the conditions but are not identical to the
182 * most common one will be marked as esolNO in the solvent_type array.
184 * TIP3p is identical to SPC for these purposes, so we call it
185 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
187 * NOTE: QM particle should not
188 * become an optimized solvent. Not even if there is only one charge
189 * group in the Qm
192 typedef struct
194 int model;
195 int count;
196 int vdwtype[4];
197 real charge[4];
198 } solvent_parameters_t;
200 static void
201 check_solvent_cg(const gmx_moltype_t *molt,
202 int cg0,
203 int nmol,
204 const unsigned char *qm_grpnr,
205 const AtomGroupIndices *qm_grps,
206 t_forcerec *fr,
207 int *n_solvent_parameters,
208 solvent_parameters_t **solvent_parameters_p,
209 int cginfo,
210 int *cg_sp)
212 t_atom *atom;
213 int j, k;
214 int j0, j1, nj;
215 gmx_bool perturbed;
216 gmx_bool has_vdw[4];
217 gmx_bool match;
218 real tmp_charge[4] = { 0.0 }; /* init to zero to make gcc 7 happy */
219 int tmp_vdwtype[4] = { 0 }; /* init to zero to make gcc 7 happy */
220 int tjA;
221 gmx_bool qm;
222 solvent_parameters_t *solvent_parameters;
224 /* We use a list with parameters for each solvent type.
225 * Every time we discover a new molecule that fulfills the basic
226 * conditions for a solvent we compare with the previous entries
227 * in these lists. If the parameters are the same we just increment
228 * the counter for that type, and otherwise we create a new type
229 * based on the current molecule.
231 * Once we've finished going through all molecules we check which
232 * solvent is most common, and mark all those molecules while we
233 * clear the flag on all others.
236 solvent_parameters = *solvent_parameters_p;
238 /* Mark the cg first as non optimized */
239 *cg_sp = -1;
241 /* Check if this cg has no exclusions with atoms in other charge groups
242 * and all atoms inside the charge group excluded.
243 * We only have 3 or 4 atom solvent loops.
245 if (GET_CGINFO_EXCL_INTER(cginfo) ||
246 !GET_CGINFO_EXCL_INTRA(cginfo))
248 return;
251 /* Get the indices of the first atom in this charge group */
252 j0 = molt->cgs.index[cg0];
253 j1 = molt->cgs.index[cg0+1];
255 /* Number of atoms in our molecule */
256 nj = j1 - j0;
258 if (debug)
260 fprintf(debug,
261 "Moltype '%s': there are %d atoms in this charge group\n",
262 *molt->name, nj);
265 /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
266 * otherwise skip it.
268 if (nj < 3 || nj > 4)
270 return;
273 /* Check if we are doing QM on this group */
274 qm = FALSE;
275 if (qm_grpnr != nullptr)
277 for (j = j0; j < j1 && !qm; j++)
279 qm = (qm_grpnr[j] < qm_grps->size() - 1);
282 /* Cannot use solvent optimization with QM */
283 if (qm)
285 return;
288 atom = molt->atoms.atom;
290 /* Still looks like a solvent, time to check parameters */
292 /* If it is perturbed (free energy) we can't use the solvent loops,
293 * so then we just skip to the next molecule.
295 perturbed = FALSE;
297 for (j = j0; j < j1 && !perturbed; j++)
299 perturbed = PERTURBED(atom[j]);
302 if (perturbed)
304 return;
307 /* Now it's only a question if the VdW and charge parameters
308 * are OK. Before doing the check we compare and see if they are
309 * identical to a possible previous solvent type.
310 * First we assign the current types and charges.
312 for (j = 0; j < nj; j++)
314 tmp_vdwtype[j] = atom[j0+j].type;
315 tmp_charge[j] = atom[j0+j].q;
318 /* Does it match any previous solvent type? */
319 for (k = 0; k < *n_solvent_parameters; k++)
321 match = TRUE;
324 /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
325 if ( (solvent_parameters[k].model == esolSPC && nj != 3) ||
326 (solvent_parameters[k].model == esolTIP4P && nj != 4) )
328 match = FALSE;
331 /* Check that types & charges match for all atoms in molecule */
332 for (j = 0; j < nj && match; j++)
334 if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
336 match = FALSE;
338 if (tmp_charge[j] != solvent_parameters[k].charge[j])
340 match = FALSE;
343 if (match)
345 /* Congratulations! We have a matched solvent.
346 * Flag it with this type for later processing.
348 *cg_sp = k;
349 solvent_parameters[k].count += nmol;
351 /* We are done with this charge group */
352 return;
356 /* If we get here, we have a tentative new solvent type.
357 * Before we add it we must check that it fulfills the requirements
358 * of the solvent optimized loops. First determine which atoms have
359 * VdW interactions.
361 for (j = 0; j < nj; j++)
363 has_vdw[j] = FALSE;
364 tjA = tmp_vdwtype[j];
366 /* Go through all other tpes and see if any have non-zero
367 * VdW parameters when combined with this one.
369 for (k = 0; k < fr->ntype && (!has_vdw[j]); k++)
371 /* We already checked that the atoms weren't perturbed,
372 * so we only need to check state A now.
374 if (fr->bBHAM)
376 has_vdw[j] = (has_vdw[j] ||
377 (BHAMA(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
378 (BHAMB(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
379 (BHAMC(fr->nbfp, fr->ntype, tjA, k) != 0.0));
381 else
383 /* Standard LJ */
384 has_vdw[j] = (has_vdw[j] ||
385 (C6(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
386 (C12(fr->nbfp, fr->ntype, tjA, k) != 0.0));
391 /* Now we know all we need to make the final check and assignment. */
392 if (nj == 3)
394 /* So, is it an SPC?
395 * For this we require thatn all atoms have charge,
396 * the charges on atom 2 & 3 should be the same, and only
397 * atom 1 might have VdW.
399 if (!has_vdw[1] &&
400 !has_vdw[2] &&
401 tmp_charge[0] != 0 &&
402 tmp_charge[1] != 0 &&
403 tmp_charge[2] == tmp_charge[1])
405 srenew(solvent_parameters, *n_solvent_parameters+1);
406 solvent_parameters[*n_solvent_parameters].model = esolSPC;
407 solvent_parameters[*n_solvent_parameters].count = nmol;
408 for (k = 0; k < 3; k++)
410 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
411 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
414 *cg_sp = *n_solvent_parameters;
415 (*n_solvent_parameters)++;
418 else if (nj == 4)
420 /* Or could it be a TIP4P?
421 * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
422 * Only atom 1 mght have VdW.
424 if (!has_vdw[1] &&
425 !has_vdw[2] &&
426 !has_vdw[3] &&
427 tmp_charge[0] == 0 &&
428 tmp_charge[1] != 0 &&
429 tmp_charge[2] == tmp_charge[1] &&
430 tmp_charge[3] != 0)
432 srenew(solvent_parameters, *n_solvent_parameters+1);
433 solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
434 solvent_parameters[*n_solvent_parameters].count = nmol;
435 for (k = 0; k < 4; k++)
437 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
438 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
441 *cg_sp = *n_solvent_parameters;
442 (*n_solvent_parameters)++;
446 *solvent_parameters_p = solvent_parameters;
449 static void
450 check_solvent(FILE * fp,
451 const gmx_mtop_t * mtop,
452 t_forcerec * fr,
453 cginfo_mb_t *cginfo_mb)
455 const t_block * cgs;
456 const gmx_moltype_t *molt;
457 int mol, cg_mol, at_offset, am, cgm, i, nmol_ch, nmol;
458 int n_solvent_parameters;
459 solvent_parameters_t *solvent_parameters;
460 int **cg_sp;
461 int bestsp, bestsol;
463 if (debug)
465 fprintf(debug, "Going to determine what solvent types we have.\n");
468 n_solvent_parameters = 0;
469 solvent_parameters = nullptr;
470 /* Allocate temporary array for solvent type */
471 snew(cg_sp, mtop->molblock.size());
473 at_offset = 0;
474 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
476 molt = &mtop->moltype[mtop->molblock[mb].type];
477 cgs = &molt->cgs;
478 /* Here we have to loop over all individual molecules
479 * because we need to check for QMMM particles.
481 snew(cg_sp[mb], cginfo_mb[mb].cg_mod);
482 nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
483 nmol = mtop->molblock[mb].nmol/nmol_ch;
484 for (mol = 0; mol < nmol_ch; mol++)
486 cgm = mol*cgs->nr;
487 am = mol*cgs->index[cgs->nr];
488 for (cg_mol = 0; cg_mol < cgs->nr; cg_mol++)
490 check_solvent_cg(molt, cg_mol, nmol,
491 mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty() ?
492 nullptr : mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].data()+at_offset+am,
493 &mtop->groups.groups[SimulationAtomGroupType::QuantumMechanics],
495 &n_solvent_parameters, &solvent_parameters,
496 cginfo_mb[mb].cginfo[cgm+cg_mol],
497 &cg_sp[mb][cgm+cg_mol]);
500 at_offset += cgs->index[cgs->nr];
503 /* Puh! We finished going through all charge groups.
504 * Now find the most common solvent model.
507 /* Most common solvent this far */
508 bestsp = -2;
509 for (i = 0; i < n_solvent_parameters; i++)
511 if (bestsp == -2 ||
512 solvent_parameters[i].count > solvent_parameters[bestsp].count)
514 bestsp = i;
518 if (bestsp >= 0)
520 bestsol = solvent_parameters[bestsp].model;
522 else
524 bestsol = esolNO;
527 fr->nWatMol = 0;
528 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
530 cgs = &mtop->moltype[mtop->molblock[mb].type].cgs;
531 nmol = (mtop->molblock[mb].nmol*cgs->nr)/cginfo_mb[mb].cg_mod;
532 for (i = 0; i < cginfo_mb[mb].cg_mod; i++)
534 if (cg_sp[mb][i] == bestsp)
536 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], bestsol);
537 fr->nWatMol += nmol;
539 else
541 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], esolNO);
544 sfree(cg_sp[mb]);
546 sfree(cg_sp);
548 if (bestsol != esolNO && fp != nullptr)
550 fprintf(fp, "\nEnabling %s-like water optimization for %d molecules.\n\n",
551 esol_names[bestsol],
552 solvent_parameters[bestsp].count);
555 sfree(solvent_parameters);
556 fr->solvent_opt = bestsol;
559 enum {
560 acNONE = 0, acCONSTRAINT, acSETTLE
563 static cginfo_mb_t *init_cginfo_mb(FILE *fplog, const gmx_mtop_t *mtop,
564 t_forcerec *fr, gmx_bool bNoSolvOpt,
565 gmx_bool *bFEP_NonBonded,
566 gmx_bool *bExcl_IntraCGAll_InterCGNone)
568 const t_block *cgs;
569 const t_blocka *excl;
570 const gmx_moltype_t *molt;
571 const gmx_molblock_t *molb;
572 cginfo_mb_t *cginfo_mb;
573 gmx_bool *type_VDW;
574 int *cginfo;
575 int cg_offset, a_offset;
576 int m, cg, a0, a1, gid, ai, j, aj, excl_nalloc;
577 int *a_con;
578 int ftype;
579 int ia;
580 gmx_bool bId, *bExcl, bExclIntraAll, bExclInter, bHaveVDW, bHaveQ, bHavePerturbedAtoms;
582 snew(cginfo_mb, mtop->molblock.size());
584 snew(type_VDW, fr->ntype);
585 for (ai = 0; ai < fr->ntype; ai++)
587 type_VDW[ai] = FALSE;
588 for (j = 0; j < fr->ntype; j++)
590 type_VDW[ai] = type_VDW[ai] ||
591 fr->bBHAM ||
592 C6(fr->nbfp, fr->ntype, ai, j) != 0 ||
593 C12(fr->nbfp, fr->ntype, ai, j) != 0;
597 *bFEP_NonBonded = FALSE;
598 *bExcl_IntraCGAll_InterCGNone = TRUE;
600 excl_nalloc = 10;
601 snew(bExcl, excl_nalloc);
602 cg_offset = 0;
603 a_offset = 0;
604 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
606 molb = &mtop->molblock[mb];
607 molt = &mtop->moltype[molb->type];
608 cgs = &molt->cgs;
609 excl = &molt->excls;
611 /* Check if the cginfo is identical for all molecules in this block.
612 * If so, we only need an array of the size of one molecule.
613 * Otherwise we make an array of #mol times #cgs per molecule.
615 bId = TRUE;
616 for (m = 0; m < molb->nmol; m++)
618 int am = m*cgs->index[cgs->nr];
619 for (cg = 0; cg < cgs->nr; cg++)
621 a0 = cgs->index[cg];
622 a1 = cgs->index[cg+1];
623 if (getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset+am+a0) !=
624 getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset +a0))
626 bId = FALSE;
628 if (!mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty())
630 for (ai = a0; ai < a1; ai++)
632 if (mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset+am+ai] !=
633 mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset +ai])
635 bId = FALSE;
642 cginfo_mb[mb].cg_start = cg_offset;
643 cginfo_mb[mb].cg_end = cg_offset + molb->nmol*cgs->nr;
644 cginfo_mb[mb].cg_mod = (bId ? 1 : molb->nmol)*cgs->nr;
645 snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
646 cginfo = cginfo_mb[mb].cginfo;
648 /* Set constraints flags for constrained atoms */
649 snew(a_con, molt->atoms.nr);
650 for (ftype = 0; ftype < F_NRE; ftype++)
652 if (interaction_function[ftype].flags & IF_CONSTRAINT)
654 int nral;
656 nral = NRAL(ftype);
657 for (ia = 0; ia < molt->ilist[ftype].size(); ia += 1+nral)
659 int a;
661 for (a = 0; a < nral; a++)
663 a_con[molt->ilist[ftype].iatoms[ia+1+a]] =
664 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
670 for (m = 0; m < (bId ? 1 : molb->nmol); m++)
672 int cgm = m*cgs->nr;
673 int am = m*cgs->index[cgs->nr];
674 for (cg = 0; cg < cgs->nr; cg++)
676 a0 = cgs->index[cg];
677 a1 = cgs->index[cg+1];
679 /* Store the energy group in cginfo */
680 gid = getGroupType(mtop->groups, SimulationAtomGroupType::EnergyOutput, a_offset+am+a0);
681 SET_CGINFO_GID(cginfo[cgm+cg], gid);
683 /* Check the intra/inter charge group exclusions */
684 if (a1-a0 > excl_nalloc)
686 excl_nalloc = a1 - a0;
687 srenew(bExcl, excl_nalloc);
689 /* bExclIntraAll: all intra cg interactions excluded
690 * bExclInter: any inter cg interactions excluded
692 bExclIntraAll = TRUE;
693 bExclInter = FALSE;
694 bHaveVDW = FALSE;
695 bHaveQ = FALSE;
696 bHavePerturbedAtoms = FALSE;
697 for (ai = a0; ai < a1; ai++)
699 /* Check VDW and electrostatic interactions */
700 bHaveVDW = bHaveVDW || (type_VDW[molt->atoms.atom[ai].type] ||
701 type_VDW[molt->atoms.atom[ai].typeB]);
702 bHaveQ = bHaveQ || (molt->atoms.atom[ai].q != 0 ||
703 molt->atoms.atom[ai].qB != 0);
705 bHavePerturbedAtoms = bHavePerturbedAtoms || (PERTURBED(molt->atoms.atom[ai]) != 0);
707 /* Clear the exclusion list for atom ai */
708 for (aj = a0; aj < a1; aj++)
710 bExcl[aj-a0] = FALSE;
712 /* Loop over all the exclusions of atom ai */
713 for (j = excl->index[ai]; j < excl->index[ai+1]; j++)
715 aj = excl->a[j];
716 if (aj < a0 || aj >= a1)
718 bExclInter = TRUE;
720 else
722 bExcl[aj-a0] = TRUE;
725 /* Check if ai excludes a0 to a1 */
726 for (aj = a0; aj < a1; aj++)
728 if (!bExcl[aj-a0])
730 bExclIntraAll = FALSE;
734 switch (a_con[ai])
736 case acCONSTRAINT:
737 SET_CGINFO_CONSTR(cginfo[cgm+cg]);
738 break;
739 case acSETTLE:
740 SET_CGINFO_SETTLE(cginfo[cgm+cg]);
741 break;
742 default:
743 break;
746 if (bExclIntraAll)
748 SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
750 if (bExclInter)
752 SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
754 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
756 /* The size in cginfo is currently only read with DD */
757 gmx_fatal(FARGS, "A charge group has size %d which is larger than the limit of %d atoms", a1-a0, MAX_CHARGEGROUP_SIZE);
759 if (bHaveVDW)
761 SET_CGINFO_HAS_VDW(cginfo[cgm+cg]);
763 if (bHaveQ)
765 SET_CGINFO_HAS_Q(cginfo[cgm+cg]);
767 if (bHavePerturbedAtoms && fr->efep != efepNO)
769 SET_CGINFO_FEP(cginfo[cgm+cg]);
770 *bFEP_NonBonded = TRUE;
772 /* Store the charge group size */
773 SET_CGINFO_NATOMS(cginfo[cgm+cg], a1-a0);
775 if (!bExclIntraAll || bExclInter)
777 *bExcl_IntraCGAll_InterCGNone = FALSE;
782 sfree(a_con);
784 cg_offset += molb->nmol*cgs->nr;
785 a_offset += molb->nmol*cgs->index[cgs->nr];
787 sfree(type_VDW);
788 sfree(bExcl);
790 /* the solvent optimizer is called after the QM is initialized,
791 * because we don't want to have the QM subsystemto become an
792 * optimized solvent
795 check_solvent(fplog, mtop, fr, cginfo_mb);
797 if (getenv("GMX_NO_SOLV_OPT"))
799 if (fplog)
801 fprintf(fplog, "Found environment variable GMX_NO_SOLV_OPT.\n"
802 "Disabling all solvent optimization\n");
804 fr->solvent_opt = esolNO;
806 if (bNoSolvOpt)
808 fr->solvent_opt = esolNO;
810 if (!fr->solvent_opt)
812 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
814 for (cg = 0; cg < cginfo_mb[mb].cg_mod; cg++)
816 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg], esolNO);
821 return cginfo_mb;
824 static std::vector<int> cginfo_expand(const int nmb,
825 const cginfo_mb_t *cgi_mb)
827 const int ncg = cgi_mb[nmb - 1].cg_end;
829 std::vector<int> cginfo(ncg);
831 int mb = 0;
832 for (int cg = 0; cg < ncg; cg++)
834 while (cg >= cgi_mb[mb].cg_end)
836 mb++;
838 cginfo[cg] =
839 cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
842 return cginfo;
845 static void done_cginfo_mb(cginfo_mb_t *cginfo_mb, int numMolBlocks)
847 if (cginfo_mb == nullptr)
849 return;
851 for (int mb = 0; mb < numMolBlocks; ++mb)
853 sfree(cginfo_mb[mb].cginfo);
855 sfree(cginfo_mb);
858 /* Sets the sum of charges (squared) and C6 in the system in fr.
859 * Returns whether the system has a net charge.
861 static bool set_chargesum(FILE *log, t_forcerec *fr, const gmx_mtop_t *mtop)
863 /*This now calculates sum for q and c6*/
864 double qsum, q2sum, q, c6sum, c6;
866 qsum = 0;
867 q2sum = 0;
868 c6sum = 0;
869 for (const gmx_molblock_t &molb : mtop->molblock)
871 int nmol = molb.nmol;
872 const t_atoms *atoms = &mtop->moltype[molb.type].atoms;
873 for (int i = 0; i < atoms->nr; i++)
875 q = atoms->atom[i].q;
876 qsum += nmol*q;
877 q2sum += nmol*q*q;
878 c6 = mtop->ffparams.iparams[atoms->atom[i].type*(mtop->ffparams.atnr+1)].lj.c6;
879 c6sum += nmol*c6;
882 fr->qsum[0] = qsum;
883 fr->q2sum[0] = q2sum;
884 fr->c6sum[0] = c6sum;
886 if (fr->efep != efepNO)
888 qsum = 0;
889 q2sum = 0;
890 c6sum = 0;
891 for (const gmx_molblock_t &molb : mtop->molblock)
893 int nmol = molb.nmol;
894 const t_atoms *atoms = &mtop->moltype[molb.type].atoms;
895 for (int i = 0; i < atoms->nr; i++)
897 q = atoms->atom[i].qB;
898 qsum += nmol*q;
899 q2sum += nmol*q*q;
900 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB*(mtop->ffparams.atnr+1)].lj.c6;
901 c6sum += nmol*c6;
903 fr->qsum[1] = qsum;
904 fr->q2sum[1] = q2sum;
905 fr->c6sum[1] = c6sum;
908 else
910 fr->qsum[1] = fr->qsum[0];
911 fr->q2sum[1] = fr->q2sum[0];
912 fr->c6sum[1] = fr->c6sum[0];
914 if (log)
916 if (fr->efep == efepNO)
918 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
920 else
922 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n",
923 fr->qsum[0], fr->qsum[1]);
927 /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
928 return (std::abs(fr->qsum[0]) > 1e-4 ||
929 std::abs(fr->qsum[1]) > 1e-4);
932 static real calcBuckinghamBMax(FILE *fplog, const gmx_mtop_t *mtop)
934 const t_atoms *at1, *at2;
935 int i, j, tpi, tpj, ntypes;
936 real b, bmin;
938 if (fplog)
940 fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
942 ntypes = mtop->ffparams.atnr;
944 bmin = -1;
945 real bham_b_max = 0;
946 for (size_t mt1 = 0; mt1 < mtop->moltype.size(); mt1++)
948 at1 = &mtop->moltype[mt1].atoms;
949 for (i = 0; (i < at1->nr); i++)
951 tpi = at1->atom[i].type;
952 if (tpi >= ntypes)
954 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
957 for (size_t mt2 = mt1; mt2 < mtop->moltype.size(); mt2++)
959 at2 = &mtop->moltype[mt2].atoms;
960 for (j = 0; (j < at2->nr); j++)
962 tpj = at2->atom[j].type;
963 if (tpj >= ntypes)
965 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
967 b = mtop->ffparams.iparams[tpi*ntypes + tpj].bham.b;
968 if (b > bham_b_max)
970 bham_b_max = b;
972 if ((b < bmin) || (bmin == -1))
974 bmin = b;
980 if (fplog)
982 fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n",
983 bmin, bham_b_max);
986 return bham_b_max;
989 /*!\brief If there's bonded interactions of type \c ftype1 or \c
990 * ftype2 present in the topology, build an array of the number of
991 * interactions present for each bonded interaction index found in the
992 * topology.
994 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
995 * valid type with that parameter.
997 * \c count will be reallocated as necessary to fit the largest bonded
998 * interaction index found, and its current size will be returned in
999 * \c ncount. It will contain zero for every bonded interaction index
1000 * for which no interactions are present in the topology.
1002 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t *mtop,
1003 int *ncount, int **count)
1005 int ftype, i, j, tabnr;
1007 // Loop over all moleculetypes
1008 for (const gmx_moltype_t &molt : mtop->moltype)
1010 // Loop over all interaction types
1011 for (ftype = 0; ftype < F_NRE; ftype++)
1013 // If the current interaction type is one of the types whose tables we're trying to count...
1014 if (ftype == ftype1 || ftype == ftype2)
1016 const InteractionList &il = molt.ilist[ftype];
1017 const int stride = 1 + NRAL(ftype);
1018 // ... and there are actually some interactions for this type
1019 for (i = 0; i < il.size(); i += stride)
1021 // Find out which table index the user wanted
1022 tabnr = mtop->ffparams.iparams[il.iatoms[i]].tab.table;
1023 if (tabnr < 0)
1025 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
1027 // Make room for this index in the data structure
1028 if (tabnr >= *ncount)
1030 srenew(*count, tabnr+1);
1031 for (j = *ncount; j < tabnr+1; j++)
1033 (*count)[j] = 0;
1035 *ncount = tabnr+1;
1037 // Record that this table index is used and must have a valid file
1038 (*count)[tabnr]++;
1045 /*!\brief If there's bonded interactions of flavour \c tabext and type
1046 * \c ftype1 or \c ftype2 present in the topology, seek them in the
1047 * list of filenames passed to mdrun, and make bonded tables from
1048 * those files.
1050 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1051 * valid type with that parameter.
1053 * A fatal error occurs if no matching filename is found.
1055 static bondedtable_t *make_bonded_tables(FILE *fplog,
1056 int ftype1, int ftype2,
1057 const gmx_mtop_t *mtop,
1058 gmx::ArrayRef<const std::string> tabbfnm,
1059 const char *tabext)
1061 int ncount, *count;
1062 bondedtable_t *tab;
1064 tab = nullptr;
1066 ncount = 0;
1067 count = nullptr;
1068 count_tables(ftype1, ftype2, mtop, &ncount, &count);
1070 // Are there any relevant tabulated bond interactions?
1071 if (ncount > 0)
1073 snew(tab, ncount);
1074 for (int i = 0; i < ncount; i++)
1076 // Do any interactions exist that requires this table?
1077 if (count[i] > 0)
1079 // This pattern enforces the current requirement that
1080 // table filenames end in a characteristic sequence
1081 // before the file type extension, and avoids table 13
1082 // being recognized and used for table 1.
1083 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
1084 bool madeTable = false;
1085 for (gmx::index j = 0; j < tabbfnm.ssize() && !madeTable; ++j)
1087 if (gmx::endsWith(tabbfnm[j], patternToFind))
1089 // Finally read the table from the file found
1090 tab[i] = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1)-2);
1091 madeTable = true;
1094 if (!madeTable)
1096 bool isPlural = (ftype2 != -1);
1097 gmx_fatal(FARGS, "Tabulated interaction of type '%s%s%s' with index %d cannot be used because no table file whose name matched '%s' was passed via the gmx mdrun -tableb command-line option.",
1098 interaction_function[ftype1].longname,
1099 isPlural ? "' or '" : "",
1100 isPlural ? interaction_function[ftype2].longname : "",
1102 patternToFind.c_str());
1106 sfree(count);
1109 return tab;
1112 void forcerec_set_ranges(t_forcerec *fr,
1113 int ncg_home, int ncg_force,
1114 int natoms_force,
1115 int natoms_force_constr, int natoms_f_novirsum)
1117 fr->cg0 = 0;
1118 fr->hcg = ncg_home;
1120 /* fr->ncg_force is unused in the standard code,
1121 * but it can be useful for modified code dealing with charge groups.
1123 fr->ncg_force = ncg_force;
1124 fr->natoms_force = natoms_force;
1125 fr->natoms_force_constr = natoms_force_constr;
1127 if (fr->natoms_force_constr > fr->nalloc_force)
1129 fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1132 if (fr->haveDirectVirialContributions)
1134 fr->forceBufferForDirectVirialContributions.resize(natoms_f_novirsum);
1138 static real cutoff_inf(real cutoff)
1140 if (cutoff == 0)
1142 cutoff = GMX_CUTOFF_INF;
1145 return cutoff;
1148 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
1149 static void initCoulombEwaldParameters(FILE *fp, const t_inputrec *ir,
1150 bool systemHasNetCharge,
1151 interaction_const_t *ic)
1153 if (!EEL_PME_EWALD(ir->coulombtype))
1155 return;
1158 if (fp)
1160 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
1162 if (ir->coulombtype == eelP3M_AD)
1164 please_cite(fp, "Hockney1988");
1165 please_cite(fp, "Ballenegger2012");
1167 else
1169 please_cite(fp, "Essmann95a");
1172 if (ir->ewald_geometry == eewg3DC)
1174 if (fp)
1176 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
1177 systemHasNetCharge ? " and net charge" : "");
1179 please_cite(fp, "In-Chul99a");
1180 if (systemHasNetCharge)
1182 please_cite(fp, "Ballenegger2009");
1187 ic->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
1188 if (fp)
1190 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n",
1191 1/ic->ewaldcoeff_q);
1194 if (ic->coulomb_modifier == eintmodPOTSHIFT)
1196 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
1197 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q*ic->rcoulomb) / ic->rcoulomb;
1199 else
1201 ic->sh_ewald = 0;
1205 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
1206 static void initVdwEwaldParameters(FILE *fp, const t_inputrec *ir,
1207 interaction_const_t *ic)
1209 if (!EVDW_PME(ir->vdwtype))
1211 return;
1214 if (fp)
1216 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
1217 please_cite(fp, "Essmann95a");
1219 ic->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
1220 if (fp)
1222 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n",
1223 1/ic->ewaldcoeff_lj);
1226 if (ic->vdw_modifier == eintmodPOTSHIFT)
1228 real crc2 = gmx::square(ic->ewaldcoeff_lj*ic->rvdw);
1229 ic->sh_lj_ewald = (std::exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)/gmx::power6(ic->rvdw);
1231 else
1233 ic->sh_lj_ewald = 0;
1237 /* Generate Coulomb and/or Van der Waals Ewald long-range correction tables
1239 * Tables are generated for one or both, depending on if the pointers
1240 * are non-null. The spacing for both table sets is the same and obeys
1241 * both accuracy requirements, when relevant.
1243 static void init_ewald_f_table(const interaction_const_t &ic,
1244 EwaldCorrectionTables *coulombTables,
1245 EwaldCorrectionTables *vdwTables)
1247 const bool useCoulombTable = (EEL_PME_EWALD(ic.eeltype) && coulombTables != nullptr);
1248 const bool useVdwTable = (EVDW_PME(ic.vdwtype) && vdwTables != nullptr);
1250 /* Get the Ewald table spacing based on Coulomb and/or LJ
1251 * Ewald coefficients and rtol.
1253 const real tableScale = ewald_spline3_table_scale(ic, useCoulombTable, useVdwTable);
1255 const int tableSize = static_cast<int>(ic.rcoulomb*tableScale) + 2;
1257 if (useCoulombTable)
1259 *coulombTables = generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_q, v_q_ewald_lr);
1262 if (useVdwTable)
1264 *vdwTables = generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_lj, v_lj_ewald_lr);
1268 void init_interaction_const_tables(FILE *fp,
1269 interaction_const_t *ic)
1271 if (EEL_PME_EWALD(ic->eeltype))
1273 init_ewald_f_table(*ic, ic->coulombEwaldTables.get(), nullptr);
1274 if (fp != nullptr)
1276 fprintf(fp, "Initialized non-bonded Coulomb Ewald tables, spacing: %.2e size: %zu\n\n",
1277 1/ic->coulombEwaldTables->scale, ic->coulombEwaldTables->tableF.size());
1282 static void clear_force_switch_constants(shift_consts_t *sc)
1284 sc->c2 = 0;
1285 sc->c3 = 0;
1286 sc->cpot = 0;
1289 static void force_switch_constants(real p,
1290 real rsw, real rc,
1291 shift_consts_t *sc)
1293 /* Here we determine the coefficient for shifting the force to zero
1294 * between distance rsw and the cut-off rc.
1295 * For a potential of r^-p, we have force p*r^-(p+1).
1296 * But to save flops we absorb p in the coefficient.
1297 * Thus we get:
1298 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
1299 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
1301 sc->c2 = ((p + 1)*rsw - (p + 4)*rc)/(pow(rc, p + 2)*gmx::square(rc - rsw));
1302 sc->c3 = -((p + 1)*rsw - (p + 3)*rc)/(pow(rc, p + 2)*gmx::power3(rc - rsw));
1303 sc->cpot = -pow(rc, -p) + p*sc->c2/3*gmx::power3(rc - rsw) + p*sc->c3/4*gmx::power4(rc - rsw);
1306 static void potential_switch_constants(real rsw, real rc,
1307 switch_consts_t *sc)
1309 /* The switch function is 1 at rsw and 0 at rc.
1310 * The derivative and second derivate are zero at both ends.
1311 * rsw = max(r - r_switch, 0)
1312 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
1313 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
1314 * force = force*dsw - potential*sw
1315 * potential *= sw
1317 sc->c3 = -10/gmx::power3(rc - rsw);
1318 sc->c4 = 15/gmx::power4(rc - rsw);
1319 sc->c5 = -6/gmx::power5(rc - rsw);
1322 /*! \brief Construct interaction constants
1324 * This data is used (particularly) by search and force code for
1325 * short-range interactions. Many of these are constant for the whole
1326 * simulation; some are constant only after PME tuning completes.
1328 static void
1329 init_interaction_const(FILE *fp,
1330 interaction_const_t **interaction_const,
1331 const t_inputrec *ir,
1332 const gmx_mtop_t *mtop,
1333 bool systemHasNetCharge)
1335 interaction_const_t *ic = new interaction_const_t;
1337 ic->cutoff_scheme = ir->cutoff_scheme;
1339 ic->coulombEwaldTables = std::make_unique<EwaldCorrectionTables>();
1341 /* Lennard-Jones */
1342 ic->vdwtype = ir->vdwtype;
1343 ic->vdw_modifier = ir->vdw_modifier;
1344 ic->reppow = mtop->ffparams.reppow;
1345 ic->rvdw = cutoff_inf(ir->rvdw);
1346 ic->rvdw_switch = ir->rvdw_switch;
1347 ic->ljpme_comb_rule = ir->ljpme_combination_rule;
1348 ic->useBuckingham = (mtop->ffparams.functype[0] == F_BHAM);
1349 if (ic->useBuckingham)
1351 ic->buckinghamBMax = calcBuckinghamBMax(fp, mtop);
1354 initVdwEwaldParameters(fp, ir, ic);
1356 clear_force_switch_constants(&ic->dispersion_shift);
1357 clear_force_switch_constants(&ic->repulsion_shift);
1359 switch (ic->vdw_modifier)
1361 case eintmodPOTSHIFT:
1362 /* Only shift the potential, don't touch the force */
1363 ic->dispersion_shift.cpot = -1.0/gmx::power6(ic->rvdw);
1364 ic->repulsion_shift.cpot = -1.0/gmx::power12(ic->rvdw);
1365 break;
1366 case eintmodFORCESWITCH:
1367 /* Switch the force, switch and shift the potential */
1368 force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw,
1369 &ic->dispersion_shift);
1370 force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw,
1371 &ic->repulsion_shift);
1372 break;
1373 case eintmodPOTSWITCH:
1374 /* Switch the potential and force */
1375 potential_switch_constants(ic->rvdw_switch, ic->rvdw,
1376 &ic->vdw_switch);
1377 break;
1378 case eintmodNONE:
1379 case eintmodEXACTCUTOFF:
1380 /* Nothing to do here */
1381 break;
1382 default:
1383 gmx_incons("unimplemented potential modifier");
1386 /* Electrostatics */
1387 ic->eeltype = ir->coulombtype;
1388 ic->coulomb_modifier = ir->coulomb_modifier;
1389 ic->rcoulomb = cutoff_inf(ir->rcoulomb);
1390 ic->rcoulomb_switch = ir->rcoulomb_switch;
1391 ic->epsilon_r = ir->epsilon_r;
1393 /* Set the Coulomb energy conversion factor */
1394 if (ic->epsilon_r != 0)
1396 ic->epsfac = ONE_4PI_EPS0/ic->epsilon_r;
1398 else
1400 /* eps = 0 is infinite dieletric: no Coulomb interactions */
1401 ic->epsfac = 0;
1404 /* Reaction-field */
1405 if (EEL_RF(ic->eeltype))
1407 GMX_RELEASE_ASSERT(ic->eeltype != eelGRF_NOTUSED, "GRF is no longer supported");
1408 ic->epsilon_rf = ir->epsilon_rf;
1410 calc_rffac(fp, ic->epsilon_r, ic->epsilon_rf,
1411 ic->rcoulomb,
1412 &ic->k_rf, &ic->c_rf);
1414 else
1416 /* For plain cut-off we might use the reaction-field kernels */
1417 ic->epsilon_rf = ic->epsilon_r;
1418 ic->k_rf = 0;
1419 if (ir->coulomb_modifier == eintmodPOTSHIFT)
1421 ic->c_rf = 1/ic->rcoulomb;
1423 else
1425 ic->c_rf = 0;
1429 initCoulombEwaldParameters(fp, ir, systemHasNetCharge, ic);
1431 if (fp != nullptr)
1433 real dispersion_shift;
1435 dispersion_shift = ic->dispersion_shift.cpot;
1436 if (EVDW_PME(ic->vdwtype))
1438 dispersion_shift -= ic->sh_lj_ewald;
1440 fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
1441 ic->repulsion_shift.cpot, dispersion_shift);
1443 if (ic->eeltype == eelCUT)
1445 fprintf(fp, ", Coulomb %.e", -ic->c_rf);
1447 else if (EEL_PME(ic->eeltype))
1449 fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
1451 fprintf(fp, "\n");
1454 *interaction_const = ic;
1457 void init_forcerec(FILE *fp,
1458 const gmx::MDLogger &mdlog,
1459 t_forcerec *fr,
1460 t_fcdata *fcd,
1461 const t_inputrec *ir,
1462 const gmx_mtop_t *mtop,
1463 const t_commrec *cr,
1464 matrix box,
1465 const char *tabfn,
1466 const char *tabpfn,
1467 gmx::ArrayRef<const std::string> tabbfnm,
1468 const gmx_hw_info_t &hardwareInfo,
1469 const gmx_device_info_t *deviceInfo,
1470 const bool useGpuForBonded,
1471 gmx_bool bNoSolvOpt,
1472 real print_force,
1473 gmx_wallcycle *wcycle)
1475 real rtab;
1476 char *env;
1477 double dbl;
1478 gmx_bool bGenericKernelOnly;
1479 gmx_bool bFEP_NonBonded;
1481 /* By default we turn SIMD kernels on, but it might be turned off further down... */
1482 fr->use_simd_kernels = TRUE;
1484 if (check_box(ir->ePBC, box))
1486 gmx_fatal(FARGS, "%s", check_box(ir->ePBC, box));
1489 /* Test particle insertion ? */
1490 if (EI_TPI(ir->eI))
1492 /* Set to the size of the molecule to be inserted (the last one) */
1493 gmx::RangePartitioning molecules = gmx_mtop_molecules(*mtop);
1494 fr->n_tpi = molecules.block(molecules.numBlocks() - 1).size();
1496 else
1498 fr->n_tpi = 0;
1501 if (ir->coulombtype == eelRF_NEC_UNSUPPORTED ||
1502 ir->coulombtype == eelGRF_NOTUSED)
1504 gmx_fatal(FARGS, "%s electrostatics is no longer supported",
1505 eel_names[ir->coulombtype]);
1508 if (ir->bAdress)
1510 gmx_fatal(FARGS, "AdResS simulations are no longer supported");
1512 if (ir->useTwinRange)
1514 gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
1516 /* Copy the user determined parameters */
1517 fr->userint1 = ir->userint1;
1518 fr->userint2 = ir->userint2;
1519 fr->userint3 = ir->userint3;
1520 fr->userint4 = ir->userint4;
1521 fr->userreal1 = ir->userreal1;
1522 fr->userreal2 = ir->userreal2;
1523 fr->userreal3 = ir->userreal3;
1524 fr->userreal4 = ir->userreal4;
1526 /* Shell stuff */
1527 fr->fc_stepsize = ir->fc_stepsize;
1529 /* Free energy */
1530 fr->efep = ir->efep;
1531 fr->sc_alphavdw = ir->fepvals->sc_alpha;
1532 if (ir->fepvals->bScCoul)
1534 fr->sc_alphacoul = ir->fepvals->sc_alpha;
1535 fr->sc_sigma6_min = gmx::power6(ir->fepvals->sc_sigma_min);
1537 else
1539 fr->sc_alphacoul = 0;
1540 fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
1542 fr->sc_power = ir->fepvals->sc_power;
1543 fr->sc_r_power = ir->fepvals->sc_r_power;
1544 fr->sc_sigma6_def = gmx::power6(ir->fepvals->sc_sigma);
1546 env = getenv("GMX_SCSIGMA_MIN");
1547 if (env != nullptr)
1549 dbl = 0;
1550 sscanf(env, "%20lf", &dbl);
1551 fr->sc_sigma6_min = gmx::power6(dbl);
1552 if (fp)
1554 fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
1558 fr->bNonbonded = TRUE;
1559 if (getenv("GMX_NO_NONBONDED") != nullptr)
1561 /* turn off non-bonded calculations */
1562 fr->bNonbonded = FALSE;
1563 GMX_LOG(mdlog.warning).asParagraph().appendText(
1564 "Found environment variable GMX_NO_NONBONDED.\n"
1565 "Disabling nonbonded calculations.");
1568 bGenericKernelOnly = FALSE;
1570 /* We now check in the NS code whether a particular combination of interactions
1571 * can be used with water optimization, and disable it if that is not the case.
1574 if (getenv("GMX_NB_GENERIC") != nullptr)
1576 if (fp != nullptr)
1578 fprintf(fp,
1579 "Found environment variable GMX_NB_GENERIC.\n"
1580 "Disabling all interaction-specific nonbonded kernels, will only\n"
1581 "use the slow generic ones in src/gmxlib/nonbonded/nb_generic.c\n\n");
1583 bGenericKernelOnly = TRUE;
1586 if (bGenericKernelOnly)
1588 bNoSolvOpt = TRUE;
1591 if ( (getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr) )
1593 fr->use_simd_kernels = FALSE;
1594 if (fp != nullptr)
1596 fprintf(fp,
1597 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
1598 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
1599 "(e.g. SSE2/SSE4.1/AVX).\n\n");
1603 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
1605 /* Neighbour searching stuff */
1606 fr->cutoff_scheme = ir->cutoff_scheme;
1607 fr->bGrid = (ir->ns_type == ensGRID);
1608 fr->ePBC = ir->ePBC;
1610 /* Determine if we will do PBC for distances in bonded interactions */
1611 if (fr->ePBC == epbcNONE)
1613 fr->bMolPBC = FALSE;
1615 else
1617 if (!DOMAINDECOMP(cr))
1619 gmx_bool bSHAKE;
1621 bSHAKE = (ir->eConstrAlg == econtSHAKE &&
1622 (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0 ||
1623 gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0));
1625 /* The group cut-off scheme and SHAKE assume charge groups
1626 * are whole, but not using molpbc is faster in most cases.
1627 * With intermolecular interactions we need PBC for calculating
1628 * distances between atoms in different molecules.
1630 if (bSHAKE && !mtop->bIntermolecularInteractions)
1632 fr->bMolPBC = ir->bPeriodicMols;
1634 if (bSHAKE && fr->bMolPBC)
1636 gmx_fatal(FARGS, "SHAKE is not supported with periodic molecules");
1639 else
1641 /* Not making molecules whole is faster in most cases,
1642 * but With orientation restraints we need whole molecules.
1644 fr->bMolPBC = (fcd->orires.nr == 0);
1646 if (getenv("GMX_USE_GRAPH") != nullptr)
1648 fr->bMolPBC = FALSE;
1649 if (fp)
1651 GMX_LOG(mdlog.warning).asParagraph().appendText("GMX_USE_GRAPH is set, using the graph for bonded interactions");
1654 if (mtop->bIntermolecularInteractions)
1656 GMX_LOG(mdlog.warning).asParagraph().appendText("WARNING: Molecules linked by intermolecular interactions have to reside in the same periodic image, otherwise artifacts will occur!");
1660 GMX_RELEASE_ASSERT(fr->bMolPBC || !mtop->bIntermolecularInteractions, "We need to use PBC within molecules with inter-molecular interactions");
1662 if (bSHAKE && fr->bMolPBC)
1664 gmx_fatal(FARGS, "SHAKE is not properly supported with intermolecular interactions. For short simulations where linked molecules remain in the same periodic image, the environment variable GMX_USE_GRAPH can be used to override this check.\n");
1668 else
1670 fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->ePBC);
1674 fr->rc_scaling = ir->refcoord_scaling;
1675 copy_rvec(ir->posres_com, fr->posres_com);
1676 copy_rvec(ir->posres_comB, fr->posres_comB);
1677 fr->rlist = cutoff_inf(ir->rlist);
1678 fr->ljpme_combination_rule = ir->ljpme_combination_rule;
1680 /* This now calculates sum for q and c6*/
1681 bool systemHasNetCharge = set_chargesum(fp, fr, mtop);
1683 /* fr->ic is used both by verlet and group kernels (to some extent) now */
1684 init_interaction_const(fp, &fr->ic, ir, mtop, systemHasNetCharge);
1685 init_interaction_const_tables(fp, fr->ic);
1687 const interaction_const_t *ic = fr->ic;
1689 /* TODO: Replace this Ewald table or move it into interaction_const_t */
1690 if (ir->coulombtype == eelEWALD)
1692 init_ewald_tab(&(fr->ewald_table), ir, fp);
1695 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
1696 switch (ic->eeltype)
1698 case eelCUT:
1699 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_COULOMB;
1700 break;
1702 case eelRF:
1703 case eelRF_ZERO:
1704 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
1705 break;
1707 case eelSWITCH:
1708 case eelSHIFT:
1709 case eelUSER:
1710 case eelENCADSHIFT:
1711 case eelPMESWITCH:
1712 case eelPMEUSER:
1713 case eelPMEUSERSWITCH:
1714 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
1715 break;
1717 case eelPME:
1718 case eelP3M_AD:
1719 case eelEWALD:
1720 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD;
1721 break;
1723 default:
1724 gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[ic->eeltype]);
1726 fr->nbkernel_elec_modifier = ic->coulomb_modifier;
1728 /* Vdw: Translate from mdp settings to kernel format */
1729 switch (ic->vdwtype)
1731 case evdwCUT:
1732 if (fr->bBHAM)
1734 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
1736 else
1738 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
1740 break;
1741 case evdwPME:
1742 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD;
1743 break;
1745 case evdwSWITCH:
1746 case evdwSHIFT:
1747 case evdwUSER:
1748 case evdwENCADSHIFT:
1749 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
1750 break;
1752 default:
1753 gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[ic->vdwtype]);
1755 fr->nbkernel_vdw_modifier = ic->vdw_modifier;
1757 if (ir->cutoff_scheme == ecutsVERLET)
1759 if (!gmx_within_tol(ic->reppow, 12.0, 10*GMX_DOUBLE_EPS))
1761 gmx_fatal(FARGS, "Cut-off scheme %s only supports LJ repulsion power 12", ecutscheme_names[ir->cutoff_scheme]);
1763 /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
1764 * while mdrun does not (and never did) support this.
1766 if (EEL_USER(fr->ic->eeltype))
1768 gmx_fatal(FARGS, "Combination of %s and cutoff scheme %s is not supported",
1769 eel_names[ir->coulombtype], ecutscheme_names[ir->cutoff_scheme]);
1772 fr->bvdwtab = FALSE;
1773 fr->bcoultab = FALSE;
1776 /* 1-4 interaction electrostatics */
1777 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
1779 fr->haveDirectVirialContributions =
1780 (EEL_FULL(ic->eeltype) || EVDW_PME(ic->vdwtype) ||
1781 fr->forceProviders->hasForceProvider() ||
1782 gmx_mtop_ftype_count(mtop, F_POSRES) > 0 ||
1783 gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 ||
1784 ir->nwall > 0 ||
1785 ir->bPull ||
1786 ir->bRot ||
1787 ir->bIMD);
1789 if (fr->shift_vec == nullptr)
1791 snew(fr->shift_vec, SHIFTS);
1794 fr->shiftForces.resize(SHIFTS);
1796 if (fr->nbfp == nullptr)
1798 fr->ntype = mtop->ffparams.atnr;
1799 fr->nbfp = mk_nbfp(&mtop->ffparams, fr->bBHAM);
1800 if (EVDW_PME(ic->vdwtype))
1802 fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
1806 /* Copy the energy group exclusions */
1807 fr->egp_flags = ir->opts.egp_flags;
1809 /* Van der Waals stuff */
1810 if ((ic->vdwtype != evdwCUT) && (ic->vdwtype != evdwUSER) && !fr->bBHAM)
1812 if (ic->rvdw_switch >= ic->rvdw)
1814 gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)",
1815 ic->rvdw_switch, ic->rvdw);
1817 if (fp)
1819 fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
1820 (ic->eeltype == eelSWITCH) ? "switched" : "shifted",
1821 ic->rvdw_switch, ic->rvdw);
1825 if (fr->bBHAM && EVDW_PME(ic->vdwtype))
1827 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
1830 if (fr->bBHAM && (ic->vdwtype == evdwSHIFT || ic->vdwtype == evdwSWITCH))
1832 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
1835 if (fr->bBHAM && fr->cutoff_scheme == ecutsVERLET)
1837 gmx_fatal(FARGS, "Verlet cutoff-scheme is not supported with Buckingham");
1840 if (fp && fr->cutoff_scheme == ecutsGROUP)
1842 fprintf(fp, "Cut-off's: NS: %g Coulomb: %g %s: %g\n",
1843 fr->rlist, ic->rcoulomb, fr->bBHAM ? "BHAM" : "LJ", ic->rvdw);
1846 if (ir->implicit_solvent)
1848 gmx_fatal(FARGS, "Implict solvation is no longer supported.");
1852 /* This code automatically gives table length tabext without cut-off's,
1853 * in that case grompp should already have checked that we do not need
1854 * normal tables and we only generate tables for 1-4 interactions.
1856 rtab = ir->rlist + ir->tabext;
1858 /* We want to use unmodified tables for 1-4 coulombic
1859 * interactions, so we must in general have an extra set of
1860 * tables. */
1861 if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
1862 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0 ||
1863 gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
1865 fr->pairsTable = make_tables(fp, ic, tabpfn, rtab,
1866 GMX_MAKETABLES_14ONLY);
1869 /* Wall stuff */
1870 fr->nwall = ir->nwall;
1871 if (ir->nwall && ir->wall_type == ewtTABLE)
1873 make_wall_tables(fp, ir, tabfn, &mtop->groups, fr);
1876 if (fcd && !tabbfnm.empty())
1878 // Need to catch std::bad_alloc
1879 // TODO Don't need to catch this here, when merging with master branch
1882 fcd->bondtab = make_bonded_tables(fp,
1883 F_TABBONDS, F_TABBONDSNC,
1884 mtop, tabbfnm, "b");
1885 fcd->angletab = make_bonded_tables(fp,
1886 F_TABANGLES, -1,
1887 mtop, tabbfnm, "a");
1888 fcd->dihtab = make_bonded_tables(fp,
1889 F_TABDIHS, -1,
1890 mtop, tabbfnm, "d");
1892 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1894 else
1896 if (debug)
1898 fprintf(debug, "No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
1902 // QM/MM initialization if requested
1903 fr->bQMMM = ir->bQMMM;
1904 if (fr->bQMMM)
1906 // Initialize QM/MM if supported
1907 if (GMX_QMMM)
1909 GMX_LOG(mdlog.info).asParagraph().
1910 appendText("Large parts of the QM/MM support is deprecated, and may be removed in a future "
1911 "version. Please get in touch with the developers if you find the support useful, "
1912 "as help is needed if the functionality is to continue to be available.");
1913 fr->qr = mk_QMMMrec();
1914 init_QMMMrec(cr, mtop, ir, fr);
1916 else
1918 gmx_incons("QM/MM was requested, but is only available when GROMACS "
1919 "is configured with QM/MM support");
1923 /* Set all the static charge group info */
1924 fr->cginfo_mb = init_cginfo_mb(fp, mtop, fr, bNoSolvOpt,
1925 &bFEP_NonBonded,
1926 &fr->bExcl_IntraCGAll_InterCGNone);
1927 if (!DOMAINDECOMP(cr))
1929 fr->cginfo = cginfo_expand(mtop->molblock.size(), fr->cginfo_mb);
1932 if (!DOMAINDECOMP(cr))
1934 forcerec_set_ranges(fr, ncg_mtop(mtop), ncg_mtop(mtop),
1935 mtop->natoms, mtop->natoms, mtop->natoms);
1938 fr->print_force = print_force;
1940 /* Initialize the thread working data for bonded interactions */
1941 fr->bondedThreading =
1942 init_bonded_threading(fp, mtop->groups.groups[SimulationAtomGroupType::EnergyOutput].size());
1944 fr->nthread_ewc = gmx_omp_nthreads_get(emntBonded);
1945 snew(fr->ewc_t, fr->nthread_ewc);
1947 if (fr->cutoff_scheme == ecutsVERLET)
1949 // We checked the cut-offs in grompp, but double-check here.
1950 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
1951 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
1953 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw, "With Verlet lists and PME we should have rcoulomb>=rvdw");
1955 else
1957 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw, "With Verlet lists and no PME rcoulomb and rvdw should be identical");
1960 fr->nbv = Nbnxm::init_nb_verlet(mdlog, bFEP_NonBonded, ir, fr,
1961 cr, hardwareInfo, deviceInfo,
1962 mtop, box, wcycle);
1964 if (useGpuForBonded)
1966 auto stream = DOMAINDECOMP(cr) ?
1967 Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, Nbnxm::InteractionLocality::NonLocal) :
1968 Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, Nbnxm::InteractionLocality::Local);
1969 // TODO the heap allocation is only needed while
1970 // t_forcerec lacks a constructor.
1971 fr->gpuBonded = new gmx::GpuBonded(mtop->ffparams,
1972 stream,
1973 wcycle);
1977 if (ir->eDispCorr != edispcNO)
1979 fr->dispersionCorrection =
1980 std::make_unique<DispersionCorrection>(*mtop, *ir, fr->bBHAM,
1981 fr->ntype,
1982 gmx::arrayRefFromArray(fr->nbfp, fr->ntype*fr->ntype*2),
1983 *fr->ic, tabfn);
1984 fr->dispersionCorrection->print(mdlog);
1987 if (fp != nullptr)
1989 /* Here we switch from using mdlog, which prints the newline before
1990 * the paragraph, to our old fprintf logging, which prints the newline
1991 * after the paragraph, so we should add a newline here.
1993 fprintf(fp, "\n");
1997 /* Frees GPU memory and sets a tMPI node barrier.
1999 * Note that this function needs to be called even if GPUs are not used
2000 * in this run because the PME ranks have no knowledge of whether GPUs
2001 * are used or not, but all ranks need to enter the barrier below.
2002 * \todo Remove physical node barrier from this function after making sure
2003 * that it's not needed anymore (with a shared GPU run).
2005 void free_gpu_resources(t_forcerec *fr,
2006 const gmx::PhysicalNodeCommunicator &physicalNodeCommunicator,
2007 const gmx_gpu_info_t &gpu_info)
2009 bool isPPrankUsingGPU = (fr != nullptr) && (fr->nbv != nullptr) && fr->nbv->useGpu();
2011 /* stop the GPU profiler (only CUDA) */
2012 if (gpu_info.n_dev > 0)
2014 stopGpuProfiler();
2017 if (isPPrankUsingGPU)
2019 /* Free data in GPU memory and pinned memory before destroying the GPU context */
2020 fr->nbv.reset();
2022 delete fr->gpuBonded;
2023 fr->gpuBonded = nullptr;
2026 /* With tMPI we need to wait for all ranks to finish deallocation before
2027 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
2028 * GPU and context.
2030 * This is not a concern in OpenCL where we use one context per rank which
2031 * is freed in nbnxn_gpu_free().
2033 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
2034 * but it is easier and more futureproof to call it on the whole node.
2036 if (GMX_THREAD_MPI)
2038 physicalNodeCommunicator.barrier();
2042 void done_forcerec(t_forcerec *fr, int numMolBlocks)
2044 if (fr == nullptr)
2046 // PME-only ranks don't have a forcerec
2047 return;
2049 done_cginfo_mb(fr->cginfo_mb, numMolBlocks);
2050 sfree(fr->nbfp);
2051 delete fr->ic;
2052 sfree(fr->shift_vec);
2053 sfree(fr->ewc_t);
2054 tear_down_bonded_threading(fr->bondedThreading);
2055 GMX_RELEASE_ASSERT(fr->gpuBonded == nullptr, "Should have been deleted earlier, when used");
2056 fr->bondedThreading = nullptr;
2057 delete fr;