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.
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)
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;
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;
137 static real *make_ljpme_c6grid(const gmx_ffparams_t *idef, t_forcerec *fr)
140 real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
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
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;
178 acNONE = 0, acCONSTRAINT, acSETTLE
181 static cginfo_mb_t *init_cginfo_mb(const gmx_mtop_t *mtop,
182 const t_forcerec *fr,
183 gmx_bool *bFEP_NonBonded,
184 gmx_bool *bExcl_IntraCGAll_InterCGNone)
187 const t_blocka *excl;
188 const gmx_moltype_t *molt;
189 const gmx_molblock_t *molb;
190 cginfo_mb_t *cginfo_mb;
193 int cg_offset, a_offset;
194 int m, cg, a0, a1, gid, ai, j, aj, excl_nalloc;
198 gmx_bool bId, *bExcl, bExclIntraAll, bExclInter, bHaveVDW, bHaveQ, bHavePerturbedAtoms;
200 snew(cginfo_mb, mtop->molblock.size());
202 snew(type_VDW, fr->ntype);
203 for (ai = 0; ai < fr->ntype; ai++)
205 type_VDW[ai] = FALSE;
206 for (j = 0; j < fr->ntype; j++)
208 type_VDW[ai] = type_VDW[ai] ||
210 C6(fr->nbfp, fr->ntype, ai, j) != 0 ||
211 C12(fr->nbfp, fr->ntype, ai, j) != 0;
215 *bFEP_NonBonded = FALSE;
216 *bExcl_IntraCGAll_InterCGNone = TRUE;
219 snew(bExcl, excl_nalloc);
222 for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
224 molb = &mtop->molblock[mb];
225 molt = &mtop->moltype[molb->type];
229 /* Check if the cginfo is identical for all molecules in this block.
230 * If so, we only need an array of the size of one molecule.
231 * Otherwise we make an array of #mol times #cgs per molecule.
234 for (m = 0; m < molb->nmol; m++)
236 int am = m*cgs->index[cgs->nr];
237 for (cg = 0; cg < cgs->nr; cg++)
240 a1 = cgs->index[cg+1];
241 if (getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset+am+a0) !=
242 getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset +a0))
246 if (!mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty())
248 for (ai = a0; ai < a1; ai++)
250 if (mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset+am+ai] !=
251 mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset +ai])
260 cginfo_mb[mb].cg_start = cg_offset;
261 cginfo_mb[mb].cg_end = cg_offset + molb->nmol*cgs->nr;
262 cginfo_mb[mb].cg_mod = (bId ? 1 : molb->nmol)*cgs->nr;
263 snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
264 cginfo = cginfo_mb[mb].cginfo;
266 /* Set constraints flags for constrained atoms */
267 snew(a_con, molt->atoms.nr);
268 for (ftype = 0; ftype < F_NRE; ftype++)
270 if (interaction_function[ftype].flags & IF_CONSTRAINT)
275 for (ia = 0; ia < molt->ilist[ftype].size(); ia += 1+nral)
279 for (a = 0; a < nral; a++)
281 a_con[molt->ilist[ftype].iatoms[ia+1+a]] =
282 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
288 for (m = 0; m < (bId ? 1 : molb->nmol); m++)
291 int am = m*cgs->index[cgs->nr];
292 for (cg = 0; cg < cgs->nr; cg++)
295 a1 = cgs->index[cg+1];
297 /* Store the energy group in cginfo */
298 gid = getGroupType(mtop->groups, SimulationAtomGroupType::EnergyOutput, a_offset+am+a0);
299 SET_CGINFO_GID(cginfo[cgm+cg], gid);
301 /* Check the intra/inter charge group exclusions */
302 if (a1-a0 > excl_nalloc)
304 excl_nalloc = a1 - a0;
305 srenew(bExcl, excl_nalloc);
307 /* bExclIntraAll: all intra cg interactions excluded
308 * bExclInter: any inter cg interactions excluded
310 bExclIntraAll = TRUE;
314 bHavePerturbedAtoms = FALSE;
315 for (ai = a0; ai < a1; ai++)
317 /* Check VDW and electrostatic interactions */
318 bHaveVDW = bHaveVDW || (type_VDW[molt->atoms.atom[ai].type] ||
319 type_VDW[molt->atoms.atom[ai].typeB]);
320 bHaveQ = bHaveQ || (molt->atoms.atom[ai].q != 0 ||
321 molt->atoms.atom[ai].qB != 0);
323 bHavePerturbedAtoms = bHavePerturbedAtoms || (PERTURBED(molt->atoms.atom[ai]) != 0);
325 /* Clear the exclusion list for atom ai */
326 for (aj = a0; aj < a1; aj++)
328 bExcl[aj-a0] = FALSE;
330 /* Loop over all the exclusions of atom ai */
331 for (j = excl->index[ai]; j < excl->index[ai+1]; j++)
334 if (aj < a0 || aj >= a1)
343 /* Check if ai excludes a0 to a1 */
344 for (aj = a0; aj < a1; aj++)
348 bExclIntraAll = FALSE;
355 SET_CGINFO_CONSTR(cginfo[cgm+cg]);
358 SET_CGINFO_SETTLE(cginfo[cgm+cg]);
366 SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
370 SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
372 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
374 /* The size in cginfo is currently only read with DD */
375 gmx_fatal(FARGS, "A charge group has size %d which is larger than the limit of %d atoms", a1-a0, MAX_CHARGEGROUP_SIZE);
379 SET_CGINFO_HAS_VDW(cginfo[cgm+cg]);
383 SET_CGINFO_HAS_Q(cginfo[cgm+cg]);
385 if (bHavePerturbedAtoms && fr->efep != efepNO)
387 SET_CGINFO_FEP(cginfo[cgm+cg]);
388 *bFEP_NonBonded = TRUE;
391 if (!bExclIntraAll || bExclInter)
393 *bExcl_IntraCGAll_InterCGNone = FALSE;
400 cg_offset += molb->nmol*cgs->nr;
401 a_offset += molb->nmol*cgs->index[cgs->nr];
409 static std::vector<int> cginfo_expand(const int nmb,
410 const cginfo_mb_t *cgi_mb)
412 const int ncg = cgi_mb[nmb - 1].cg_end;
414 std::vector<int> cginfo(ncg);
417 for (int cg = 0; cg < ncg; cg++)
419 while (cg >= cgi_mb[mb].cg_end)
424 cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
430 static void done_cginfo_mb(cginfo_mb_t *cginfo_mb, int numMolBlocks)
432 if (cginfo_mb == nullptr)
436 for (int mb = 0; mb < numMolBlocks; ++mb)
438 sfree(cginfo_mb[mb].cginfo);
443 /* Sets the sum of charges (squared) and C6 in the system in fr.
444 * Returns whether the system has a net charge.
446 static bool set_chargesum(FILE *log, t_forcerec *fr, const gmx_mtop_t *mtop)
448 /*This now calculates sum for q and c6*/
449 double qsum, q2sum, q, c6sum, c6;
454 for (const gmx_molblock_t &molb : mtop->molblock)
456 int nmol = molb.nmol;
457 const t_atoms *atoms = &mtop->moltype[molb.type].atoms;
458 for (int i = 0; i < atoms->nr; i++)
460 q = atoms->atom[i].q;
463 c6 = mtop->ffparams.iparams[atoms->atom[i].type*(mtop->ffparams.atnr+1)].lj.c6;
468 fr->q2sum[0] = q2sum;
469 fr->c6sum[0] = c6sum;
471 if (fr->efep != efepNO)
476 for (const gmx_molblock_t &molb : mtop->molblock)
478 int nmol = molb.nmol;
479 const t_atoms *atoms = &mtop->moltype[molb.type].atoms;
480 for (int i = 0; i < atoms->nr; i++)
482 q = atoms->atom[i].qB;
485 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB*(mtop->ffparams.atnr+1)].lj.c6;
489 fr->q2sum[1] = q2sum;
490 fr->c6sum[1] = c6sum;
495 fr->qsum[1] = fr->qsum[0];
496 fr->q2sum[1] = fr->q2sum[0];
497 fr->c6sum[1] = fr->c6sum[0];
501 if (fr->efep == efepNO)
503 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
507 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n",
508 fr->qsum[0], fr->qsum[1]);
512 /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
513 return (std::abs(fr->qsum[0]) > 1e-4 ||
514 std::abs(fr->qsum[1]) > 1e-4);
517 static real calcBuckinghamBMax(FILE *fplog, const gmx_mtop_t *mtop)
519 const t_atoms *at1, *at2;
520 int i, j, tpi, tpj, ntypes;
525 fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
527 ntypes = mtop->ffparams.atnr;
531 for (size_t mt1 = 0; mt1 < mtop->moltype.size(); mt1++)
533 at1 = &mtop->moltype[mt1].atoms;
534 for (i = 0; (i < at1->nr); i++)
536 tpi = at1->atom[i].type;
539 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
542 for (size_t mt2 = mt1; mt2 < mtop->moltype.size(); mt2++)
544 at2 = &mtop->moltype[mt2].atoms;
545 for (j = 0; (j < at2->nr); j++)
547 tpj = at2->atom[j].type;
550 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
552 b = mtop->ffparams.iparams[tpi*ntypes + tpj].bham.b;
557 if ((b < bmin) || (bmin == -1))
567 fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n",
574 /*!\brief If there's bonded interactions of type \c ftype1 or \c
575 * ftype2 present in the topology, build an array of the number of
576 * interactions present for each bonded interaction index found in the
579 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
580 * valid type with that parameter.
582 * \c count will be reallocated as necessary to fit the largest bonded
583 * interaction index found, and its current size will be returned in
584 * \c ncount. It will contain zero for every bonded interaction index
585 * for which no interactions are present in the topology.
587 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t *mtop,
588 int *ncount, int **count)
590 int ftype, i, j, tabnr;
592 // Loop over all moleculetypes
593 for (const gmx_moltype_t &molt : mtop->moltype)
595 // Loop over all interaction types
596 for (ftype = 0; ftype < F_NRE; ftype++)
598 // If the current interaction type is one of the types whose tables we're trying to count...
599 if (ftype == ftype1 || ftype == ftype2)
601 const InteractionList &il = molt.ilist[ftype];
602 const int stride = 1 + NRAL(ftype);
603 // ... and there are actually some interactions for this type
604 for (i = 0; i < il.size(); i += stride)
606 // Find out which table index the user wanted
607 tabnr = mtop->ffparams.iparams[il.iatoms[i]].tab.table;
610 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
612 // Make room for this index in the data structure
613 if (tabnr >= *ncount)
615 srenew(*count, tabnr+1);
616 for (j = *ncount; j < tabnr+1; j++)
622 // Record that this table index is used and must have a valid file
630 /*!\brief If there's bonded interactions of flavour \c tabext and type
631 * \c ftype1 or \c ftype2 present in the topology, seek them in the
632 * list of filenames passed to mdrun, and make bonded tables from
635 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
636 * valid type with that parameter.
638 * A fatal error occurs if no matching filename is found.
640 static bondedtable_t *make_bonded_tables(FILE *fplog,
641 int ftype1, int ftype2,
642 const gmx_mtop_t *mtop,
643 gmx::ArrayRef<const std::string> tabbfnm,
653 count_tables(ftype1, ftype2, mtop, &ncount, &count);
655 // Are there any relevant tabulated bond interactions?
659 for (int i = 0; i < ncount; i++)
661 // Do any interactions exist that requires this table?
664 // This pattern enforces the current requirement that
665 // table filenames end in a characteristic sequence
666 // before the file type extension, and avoids table 13
667 // being recognized and used for table 1.
668 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
669 bool madeTable = false;
670 for (gmx::index j = 0; j < tabbfnm.ssize() && !madeTable; ++j)
672 if (gmx::endsWith(tabbfnm[j], patternToFind))
674 // Finally read the table from the file found
675 tab[i] = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1)-2);
681 bool isPlural = (ftype2 != -1);
682 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.",
683 interaction_function[ftype1].longname,
684 isPlural ? "' or '" : "",
685 isPlural ? interaction_function[ftype2].longname : "",
687 patternToFind.c_str());
697 void forcerec_set_ranges(t_forcerec *fr,
698 int ncg_home, int ncg_force,
700 int natoms_force_constr, int natoms_f_novirsum)
705 /* fr->ncg_force is unused in the standard code,
706 * but it can be useful for modified code dealing with charge groups.
708 fr->ncg_force = ncg_force;
709 fr->natoms_force = natoms_force;
710 fr->natoms_force_constr = natoms_force_constr;
712 if (fr->natoms_force_constr > fr->nalloc_force)
714 fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
717 if (fr->haveDirectVirialContributions)
719 fr->forceBufferForDirectVirialContributions.resize(natoms_f_novirsum);
723 static real cutoff_inf(real cutoff)
727 cutoff = GMX_CUTOFF_INF;
733 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
734 static void initCoulombEwaldParameters(FILE *fp, const t_inputrec *ir,
735 bool systemHasNetCharge,
736 interaction_const_t *ic)
738 if (!EEL_PME_EWALD(ir->coulombtype))
745 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
747 if (ir->coulombtype == eelP3M_AD)
749 please_cite(fp, "Hockney1988");
750 please_cite(fp, "Ballenegger2012");
754 please_cite(fp, "Essmann95a");
757 if (ir->ewald_geometry == eewg3DC)
761 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
762 systemHasNetCharge ? " and net charge" : "");
764 please_cite(fp, "In-Chul99a");
765 if (systemHasNetCharge)
767 please_cite(fp, "Ballenegger2009");
772 ic->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
775 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n",
779 if (ic->coulomb_modifier == eintmodPOTSHIFT)
781 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
782 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q*ic->rcoulomb) / ic->rcoulomb;
790 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
791 static void initVdwEwaldParameters(FILE *fp, const t_inputrec *ir,
792 interaction_const_t *ic)
794 if (!EVDW_PME(ir->vdwtype))
801 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
802 please_cite(fp, "Essmann95a");
804 ic->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
807 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n",
808 1/ic->ewaldcoeff_lj);
811 if (ic->vdw_modifier == eintmodPOTSHIFT)
813 real crc2 = gmx::square(ic->ewaldcoeff_lj*ic->rvdw);
814 ic->sh_lj_ewald = (std::exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)/gmx::power6(ic->rvdw);
822 /* Generate Coulomb and/or Van der Waals Ewald long-range correction tables
824 * Tables are generated for one or both, depending on if the pointers
825 * are non-null. The spacing for both table sets is the same and obeys
826 * both accuracy requirements, when relevant.
828 static void init_ewald_f_table(const interaction_const_t &ic,
829 EwaldCorrectionTables *coulombTables,
830 EwaldCorrectionTables *vdwTables)
832 const bool useCoulombTable = (EEL_PME_EWALD(ic.eeltype) && coulombTables != nullptr);
833 const bool useVdwTable = (EVDW_PME(ic.vdwtype) && vdwTables != nullptr);
835 /* Get the Ewald table spacing based on Coulomb and/or LJ
836 * Ewald coefficients and rtol.
838 const real tableScale = ewald_spline3_table_scale(ic, useCoulombTable, useVdwTable);
840 const int tableSize = static_cast<int>(ic.rcoulomb*tableScale) + 2;
844 *coulombTables = generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_q, v_q_ewald_lr);
849 *vdwTables = generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_lj, v_lj_ewald_lr);
853 void init_interaction_const_tables(FILE *fp,
854 interaction_const_t *ic)
856 if (EEL_PME_EWALD(ic->eeltype))
858 init_ewald_f_table(*ic, ic->coulombEwaldTables.get(), nullptr);
861 fprintf(fp, "Initialized non-bonded Coulomb Ewald tables, spacing: %.2e size: %zu\n\n",
862 1/ic->coulombEwaldTables->scale, ic->coulombEwaldTables->tableF.size());
867 static void clear_force_switch_constants(shift_consts_t *sc)
874 static void force_switch_constants(real p,
878 /* Here we determine the coefficient for shifting the force to zero
879 * between distance rsw and the cut-off rc.
880 * For a potential of r^-p, we have force p*r^-(p+1).
881 * But to save flops we absorb p in the coefficient.
883 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
884 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
886 sc->c2 = ((p + 1)*rsw - (p + 4)*rc)/(pow(rc, p + 2)*gmx::square(rc - rsw));
887 sc->c3 = -((p + 1)*rsw - (p + 3)*rc)/(pow(rc, p + 2)*gmx::power3(rc - rsw));
888 sc->cpot = -pow(rc, -p) + p*sc->c2/3*gmx::power3(rc - rsw) + p*sc->c3/4*gmx::power4(rc - rsw);
891 static void potential_switch_constants(real rsw, real rc,
894 /* The switch function is 1 at rsw and 0 at rc.
895 * The derivative and second derivate are zero at both ends.
896 * rsw = max(r - r_switch, 0)
897 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
898 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
899 * force = force*dsw - potential*sw
902 sc->c3 = -10/gmx::power3(rc - rsw);
903 sc->c4 = 15/gmx::power4(rc - rsw);
904 sc->c5 = -6/gmx::power5(rc - rsw);
907 /*! \brief Construct interaction constants
909 * This data is used (particularly) by search and force code for
910 * short-range interactions. Many of these are constant for the whole
911 * simulation; some are constant only after PME tuning completes.
914 init_interaction_const(FILE *fp,
915 interaction_const_t **interaction_const,
916 const t_inputrec *ir,
917 const gmx_mtop_t *mtop,
918 bool systemHasNetCharge)
920 interaction_const_t *ic = new interaction_const_t;
922 ic->cutoff_scheme = ir->cutoff_scheme;
924 ic->coulombEwaldTables = std::make_unique<EwaldCorrectionTables>();
927 ic->vdwtype = ir->vdwtype;
928 ic->vdw_modifier = ir->vdw_modifier;
929 ic->reppow = mtop->ffparams.reppow;
930 ic->rvdw = cutoff_inf(ir->rvdw);
931 ic->rvdw_switch = ir->rvdw_switch;
932 ic->ljpme_comb_rule = ir->ljpme_combination_rule;
933 ic->useBuckingham = (mtop->ffparams.functype[0] == F_BHAM);
934 if (ic->useBuckingham)
936 ic->buckinghamBMax = calcBuckinghamBMax(fp, mtop);
939 initVdwEwaldParameters(fp, ir, ic);
941 clear_force_switch_constants(&ic->dispersion_shift);
942 clear_force_switch_constants(&ic->repulsion_shift);
944 switch (ic->vdw_modifier)
946 case eintmodPOTSHIFT:
947 /* Only shift the potential, don't touch the force */
948 ic->dispersion_shift.cpot = -1.0/gmx::power6(ic->rvdw);
949 ic->repulsion_shift.cpot = -1.0/gmx::power12(ic->rvdw);
951 case eintmodFORCESWITCH:
952 /* Switch the force, switch and shift the potential */
953 force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw,
954 &ic->dispersion_shift);
955 force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw,
956 &ic->repulsion_shift);
958 case eintmodPOTSWITCH:
959 /* Switch the potential and force */
960 potential_switch_constants(ic->rvdw_switch, ic->rvdw,
964 case eintmodEXACTCUTOFF:
965 /* Nothing to do here */
968 gmx_incons("unimplemented potential modifier");
972 ic->eeltype = ir->coulombtype;
973 ic->coulomb_modifier = ir->coulomb_modifier;
974 ic->rcoulomb = cutoff_inf(ir->rcoulomb);
975 ic->rcoulomb_switch = ir->rcoulomb_switch;
976 ic->epsilon_r = ir->epsilon_r;
978 /* Set the Coulomb energy conversion factor */
979 if (ic->epsilon_r != 0)
981 ic->epsfac = ONE_4PI_EPS0/ic->epsilon_r;
985 /* eps = 0 is infinite dieletric: no Coulomb interactions */
990 if (EEL_RF(ic->eeltype))
992 GMX_RELEASE_ASSERT(ic->eeltype != eelGRF_NOTUSED, "GRF is no longer supported");
993 ic->epsilon_rf = ir->epsilon_rf;
995 calc_rffac(fp, ic->epsilon_r, ic->epsilon_rf,
997 &ic->k_rf, &ic->c_rf);
1001 /* For plain cut-off we might use the reaction-field kernels */
1002 ic->epsilon_rf = ic->epsilon_r;
1004 if (ir->coulomb_modifier == eintmodPOTSHIFT)
1006 ic->c_rf = 1/ic->rcoulomb;
1014 initCoulombEwaldParameters(fp, ir, systemHasNetCharge, ic);
1018 real dispersion_shift;
1020 dispersion_shift = ic->dispersion_shift.cpot;
1021 if (EVDW_PME(ic->vdwtype))
1023 dispersion_shift -= ic->sh_lj_ewald;
1025 fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
1026 ic->repulsion_shift.cpot, dispersion_shift);
1028 if (ic->eeltype == eelCUT)
1030 fprintf(fp, ", Coulomb %.e", -ic->c_rf);
1032 else if (EEL_PME(ic->eeltype))
1034 fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
1039 *interaction_const = ic;
1042 void init_forcerec(FILE *fp,
1043 const gmx::MDLogger &mdlog,
1046 const t_inputrec *ir,
1047 const gmx_mtop_t *mtop,
1048 const t_commrec *cr,
1052 gmx::ArrayRef<const std::string> tabbfnm,
1053 const gmx_hw_info_t &hardwareInfo,
1054 const gmx_device_info_t *deviceInfo,
1055 const bool useGpuForBonded,
1057 gmx_wallcycle *wcycle)
1062 gmx_bool bFEP_NonBonded;
1064 /* By default we turn SIMD kernels on, but it might be turned off further down... */
1065 fr->use_simd_kernels = TRUE;
1067 if (check_box(ir->ePBC, box))
1069 gmx_fatal(FARGS, "%s", check_box(ir->ePBC, box));
1072 /* Test particle insertion ? */
1075 /* Set to the size of the molecule to be inserted (the last one) */
1076 gmx::RangePartitioning molecules = gmx_mtop_molecules(*mtop);
1077 fr->n_tpi = molecules.block(molecules.numBlocks() - 1).size();
1084 if (ir->coulombtype == eelRF_NEC_UNSUPPORTED ||
1085 ir->coulombtype == eelGRF_NOTUSED)
1087 gmx_fatal(FARGS, "%s electrostatics is no longer supported",
1088 eel_names[ir->coulombtype]);
1093 gmx_fatal(FARGS, "AdResS simulations are no longer supported");
1095 if (ir->useTwinRange)
1097 gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
1099 /* Copy the user determined parameters */
1100 fr->userint1 = ir->userint1;
1101 fr->userint2 = ir->userint2;
1102 fr->userint3 = ir->userint3;
1103 fr->userint4 = ir->userint4;
1104 fr->userreal1 = ir->userreal1;
1105 fr->userreal2 = ir->userreal2;
1106 fr->userreal3 = ir->userreal3;
1107 fr->userreal4 = ir->userreal4;
1110 fr->fc_stepsize = ir->fc_stepsize;
1113 fr->efep = ir->efep;
1114 fr->sc_alphavdw = ir->fepvals->sc_alpha;
1115 if (ir->fepvals->bScCoul)
1117 fr->sc_alphacoul = ir->fepvals->sc_alpha;
1118 fr->sc_sigma6_min = gmx::power6(ir->fepvals->sc_sigma_min);
1122 fr->sc_alphacoul = 0;
1123 fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
1125 fr->sc_power = ir->fepvals->sc_power;
1126 fr->sc_r_power = ir->fepvals->sc_r_power;
1127 fr->sc_sigma6_def = gmx::power6(ir->fepvals->sc_sigma);
1129 env = getenv("GMX_SCSIGMA_MIN");
1133 sscanf(env, "%20lf", &dbl);
1134 fr->sc_sigma6_min = gmx::power6(dbl);
1137 fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
1141 fr->bNonbonded = TRUE;
1142 if (getenv("GMX_NO_NONBONDED") != nullptr)
1144 /* turn off non-bonded calculations */
1145 fr->bNonbonded = FALSE;
1146 GMX_LOG(mdlog.warning).asParagraph().appendText(
1147 "Found environment variable GMX_NO_NONBONDED.\n"
1148 "Disabling nonbonded calculations.");
1151 if ( (getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr) )
1153 fr->use_simd_kernels = FALSE;
1157 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
1158 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
1159 "(e.g. SSE2/SSE4.1/AVX).\n\n");
1163 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
1165 /* Neighbour searching stuff */
1166 fr->cutoff_scheme = ir->cutoff_scheme;
1167 fr->ePBC = ir->ePBC;
1169 /* Determine if we will do PBC for distances in bonded interactions */
1170 if (fr->ePBC == epbcNONE)
1172 fr->bMolPBC = FALSE;
1176 const bool useEwaldSurfaceCorrection =
1177 (EEL_PME_EWALD(ir->coulombtype) && ir->epsilon_surface != 0);
1178 if (!DOMAINDECOMP(cr))
1182 bSHAKE = (ir->eConstrAlg == econtSHAKE &&
1183 (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0 ||
1184 gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0));
1186 /* The group cut-off scheme and SHAKE assume charge groups
1187 * are whole, but not using molpbc is faster in most cases.
1188 * With intermolecular interactions we need PBC for calculating
1189 * distances between atoms in different molecules.
1191 if (bSHAKE && !mtop->bIntermolecularInteractions)
1193 fr->bMolPBC = ir->bPeriodicMols;
1195 if (bSHAKE && fr->bMolPBC)
1197 gmx_fatal(FARGS, "SHAKE is not supported with periodic molecules");
1202 /* Not making molecules whole is faster in most cases,
1203 * but with orientation restraints or non-tinfoil boundary
1204 * conditions we need whole molecules.
1206 fr->bMolPBC = (fcd->orires.nr == 0 && !useEwaldSurfaceCorrection);
1208 if (getenv("GMX_USE_GRAPH") != nullptr)
1210 fr->bMolPBC = FALSE;
1213 GMX_LOG(mdlog.warning).asParagraph().appendText("GMX_USE_GRAPH is set, using the graph for bonded interactions");
1216 if (mtop->bIntermolecularInteractions)
1218 GMX_LOG(mdlog.warning).asParagraph().appendText("WARNING: Molecules linked by intermolecular interactions have to reside in the same periodic image, otherwise artifacts will occur!");
1222 GMX_RELEASE_ASSERT(fr->bMolPBC || !mtop->bIntermolecularInteractions, "We need to use PBC within molecules with inter-molecular interactions");
1224 if (bSHAKE && fr->bMolPBC)
1226 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");
1232 fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->ePBC);
1234 if (useEwaldSurfaceCorrection && !dd_moleculesAreAlwaysWhole(*cr->dd))
1237 "You requested dipole correction (epsilon_surface > 0), but molecules are broken "
1238 "over periodic boundary conditions by the domain decomposition. "
1239 "Run without domain decomposition instead.");
1243 if (useEwaldSurfaceCorrection)
1245 GMX_RELEASE_ASSERT((!DOMAINDECOMP(cr) && !fr->bMolPBC) ||
1246 (DOMAINDECOMP(cr) && dd_moleculesAreAlwaysWhole(*cr->dd)),
1247 "Molecules can not be broken by PBC with epsilon_surface > 0");
1251 fr->rc_scaling = ir->refcoord_scaling;
1252 copy_rvec(ir->posres_com, fr->posres_com);
1253 copy_rvec(ir->posres_comB, fr->posres_comB);
1254 fr->rlist = cutoff_inf(ir->rlist);
1255 fr->ljpme_combination_rule = ir->ljpme_combination_rule;
1257 /* This now calculates sum for q and c6*/
1258 bool systemHasNetCharge = set_chargesum(fp, fr, mtop);
1260 /* fr->ic is used both by verlet and group kernels (to some extent) now */
1261 init_interaction_const(fp, &fr->ic, ir, mtop, systemHasNetCharge);
1262 init_interaction_const_tables(fp, fr->ic);
1264 const interaction_const_t *ic = fr->ic;
1266 /* TODO: Replace this Ewald table or move it into interaction_const_t */
1267 if (ir->coulombtype == eelEWALD)
1269 init_ewald_tab(&(fr->ewald_table), ir, fp);
1272 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
1273 switch (ic->eeltype)
1276 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_COULOMB;
1281 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
1290 case eelPMEUSERSWITCH:
1291 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
1297 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD;
1301 gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[ic->eeltype]);
1303 fr->nbkernel_elec_modifier = ic->coulomb_modifier;
1305 /* Vdw: Translate from mdp settings to kernel format */
1306 switch (ic->vdwtype)
1311 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
1315 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
1319 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD;
1325 case evdwENCADSHIFT:
1326 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
1330 gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[ic->vdwtype]);
1332 fr->nbkernel_vdw_modifier = ic->vdw_modifier;
1334 if (ir->cutoff_scheme == ecutsVERLET)
1336 if (!gmx_within_tol(ic->reppow, 12.0, 10*GMX_DOUBLE_EPS))
1338 gmx_fatal(FARGS, "Cut-off scheme %s only supports LJ repulsion power 12", ecutscheme_names[ir->cutoff_scheme]);
1340 /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
1341 * while mdrun does not (and never did) support this.
1343 if (EEL_USER(fr->ic->eeltype))
1345 gmx_fatal(FARGS, "Combination of %s and cutoff scheme %s is not supported",
1346 eel_names[ir->coulombtype], ecutscheme_names[ir->cutoff_scheme]);
1349 fr->bvdwtab = FALSE;
1350 fr->bcoultab = FALSE;
1353 /* 1-4 interaction electrostatics */
1354 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
1356 fr->haveDirectVirialContributions =
1357 (EEL_FULL(ic->eeltype) || EVDW_PME(ic->vdwtype) ||
1358 fr->forceProviders->hasForceProvider() ||
1359 gmx_mtop_ftype_count(mtop, F_POSRES) > 0 ||
1360 gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 ||
1366 if (fr->shift_vec == nullptr)
1368 snew(fr->shift_vec, SHIFTS);
1371 fr->shiftForces.resize(SHIFTS);
1373 if (fr->nbfp == nullptr)
1375 fr->ntype = mtop->ffparams.atnr;
1376 fr->nbfp = mk_nbfp(&mtop->ffparams, fr->bBHAM);
1377 if (EVDW_PME(ic->vdwtype))
1379 fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
1383 /* Copy the energy group exclusions */
1384 fr->egp_flags = ir->opts.egp_flags;
1386 /* Van der Waals stuff */
1387 if ((ic->vdwtype != evdwCUT) && (ic->vdwtype != evdwUSER) && !fr->bBHAM)
1389 if (ic->rvdw_switch >= ic->rvdw)
1391 gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)",
1392 ic->rvdw_switch, ic->rvdw);
1396 fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
1397 (ic->eeltype == eelSWITCH) ? "switched" : "shifted",
1398 ic->rvdw_switch, ic->rvdw);
1402 if (fr->bBHAM && EVDW_PME(ic->vdwtype))
1404 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
1407 if (fr->bBHAM && (ic->vdwtype == evdwSHIFT || ic->vdwtype == evdwSWITCH))
1409 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
1412 if (fr->bBHAM && fr->cutoff_scheme == ecutsVERLET)
1414 gmx_fatal(FARGS, "Verlet cutoff-scheme is not supported with Buckingham");
1417 if (fp && fr->cutoff_scheme == ecutsGROUP)
1419 fprintf(fp, "Cut-off's: NS: %g Coulomb: %g %s: %g\n",
1420 fr->rlist, ic->rcoulomb, fr->bBHAM ? "BHAM" : "LJ", ic->rvdw);
1423 if (ir->implicit_solvent)
1425 gmx_fatal(FARGS, "Implict solvation is no longer supported.");
1429 /* This code automatically gives table length tabext without cut-off's,
1430 * in that case grompp should already have checked that we do not need
1431 * normal tables and we only generate tables for 1-4 interactions.
1433 rtab = ir->rlist + ir->tabext;
1435 /* We want to use unmodified tables for 1-4 coulombic
1436 * interactions, so we must in general have an extra set of
1438 if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
1439 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0 ||
1440 gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
1442 fr->pairsTable = make_tables(fp, ic, tabpfn, rtab,
1443 GMX_MAKETABLES_14ONLY);
1447 fr->nwall = ir->nwall;
1448 if (ir->nwall && ir->wall_type == ewtTABLE)
1450 make_wall_tables(fp, ir, tabfn, &mtop->groups, fr);
1453 if (fcd && !tabbfnm.empty())
1455 // Need to catch std::bad_alloc
1456 // TODO Don't need to catch this here, when merging with master branch
1459 fcd->bondtab = make_bonded_tables(fp,
1460 F_TABBONDS, F_TABBONDSNC,
1461 mtop, tabbfnm, "b");
1462 fcd->angletab = make_bonded_tables(fp,
1464 mtop, tabbfnm, "a");
1465 fcd->dihtab = make_bonded_tables(fp,
1467 mtop, tabbfnm, "d");
1469 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1475 fprintf(debug, "No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
1479 // QM/MM initialization if requested
1480 fr->bQMMM = ir->bQMMM;
1483 // Initialize QM/MM if supported
1486 GMX_LOG(mdlog.info).asParagraph().
1487 appendText("Large parts of the QM/MM support is deprecated, and may be removed in a future "
1488 "version. Please get in touch with the developers if you find the support useful, "
1489 "as help is needed if the functionality is to continue to be available.");
1490 fr->qr = mk_QMMMrec();
1491 init_QMMMrec(cr, mtop, ir, fr);
1495 gmx_incons("QM/MM was requested, but is only available when GROMACS "
1496 "is configured with QM/MM support");
1500 /* Set all the static charge group info */
1501 fr->cginfo_mb = init_cginfo_mb(mtop, fr,
1503 &fr->bExcl_IntraCGAll_InterCGNone);
1504 if (!DOMAINDECOMP(cr))
1506 fr->cginfo = cginfo_expand(mtop->molblock.size(), fr->cginfo_mb);
1509 if (!DOMAINDECOMP(cr))
1511 forcerec_set_ranges(fr, ncg_mtop(mtop), ncg_mtop(mtop),
1512 mtop->natoms, mtop->natoms, mtop->natoms);
1515 fr->print_force = print_force;
1517 /* Initialize the thread working data for bonded interactions */
1518 fr->bondedThreading =
1519 init_bonded_threading(fp, mtop->groups.groups[SimulationAtomGroupType::EnergyOutput].size());
1521 fr->nthread_ewc = gmx_omp_nthreads_get(emntBonded);
1522 snew(fr->ewc_t, fr->nthread_ewc);
1524 if (fr->cutoff_scheme == ecutsVERLET)
1526 // We checked the cut-offs in grompp, but double-check here.
1527 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
1528 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
1530 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw, "With Verlet lists and PME we should have rcoulomb>=rvdw");
1534 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw, "With Verlet lists and no PME rcoulomb and rvdw should be identical");
1537 fr->nbv = Nbnxm::init_nb_verlet(mdlog, bFEP_NonBonded, ir, fr,
1538 cr, hardwareInfo, deviceInfo,
1541 if (useGpuForBonded)
1543 auto stream = DOMAINDECOMP(cr) ?
1544 Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, Nbnxm::InteractionLocality::NonLocal) :
1545 Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, Nbnxm::InteractionLocality::Local);
1546 // TODO the heap allocation is only needed while
1547 // t_forcerec lacks a constructor.
1548 fr->gpuBonded = new gmx::GpuBonded(mtop->ffparams,
1554 if (ir->eDispCorr != edispcNO)
1556 fr->dispersionCorrection =
1557 std::make_unique<DispersionCorrection>(*mtop, *ir, fr->bBHAM,
1559 gmx::arrayRefFromArray(fr->nbfp, fr->ntype*fr->ntype*2),
1561 fr->dispersionCorrection->print(mdlog);
1566 /* Here we switch from using mdlog, which prints the newline before
1567 * the paragraph, to our old fprintf logging, which prints the newline
1568 * after the paragraph, so we should add a newline here.
1574 /* Frees GPU memory and sets a tMPI node barrier.
1576 * Note that this function needs to be called even if GPUs are not used
1577 * in this run because the PME ranks have no knowledge of whether GPUs
1578 * are used or not, but all ranks need to enter the barrier below.
1579 * \todo Remove physical node barrier from this function after making sure
1580 * that it's not needed anymore (with a shared GPU run).
1582 void free_gpu_resources(t_forcerec *fr,
1583 const gmx::PhysicalNodeCommunicator &physicalNodeCommunicator,
1584 const gmx_gpu_info_t &gpu_info)
1586 bool isPPrankUsingGPU = (fr != nullptr) && (fr->nbv != nullptr) && fr->nbv->useGpu();
1588 /* stop the GPU profiler (only CUDA) */
1589 if (gpu_info.n_dev > 0)
1594 if (isPPrankUsingGPU)
1596 /* Free data in GPU memory and pinned memory before destroying the GPU context */
1599 delete fr->gpuBonded;
1600 fr->gpuBonded = nullptr;
1603 /* With tMPI we need to wait for all ranks to finish deallocation before
1604 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
1607 * This is not a concern in OpenCL where we use one context per rank which
1608 * is freed in nbnxn_gpu_free().
1610 * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1611 * but it is easier and more futureproof to call it on the whole node.
1615 physicalNodeCommunicator.barrier();
1619 void done_forcerec(t_forcerec *fr, int numMolBlocks)
1623 // PME-only ranks don't have a forcerec
1626 done_cginfo_mb(fr->cginfo_mb, numMolBlocks);
1629 sfree(fr->shift_vec);
1631 tear_down_bonded_threading(fr->bondedThreading);
1632 GMX_RELEASE_ASSERT(fr->gpuBonded == nullptr, "Should have been deleted earlier, when used");
1633 fr->bondedThreading = nullptr;