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, 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/fileio/filetypes.h"
56 #include "gromacs/gmxlib/md_logging.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/detecthardware.h"
61 #include "gromacs/listed-forces/manage-threading.h"
62 #include "gromacs/listed-forces/pairs.h"
63 #include "gromacs/math/calculate-ewald-splitting-coefficient.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/force.h"
69 #include "gromacs/mdlib/forcerec-threading.h"
70 #include "gromacs/mdlib/gmx_omp_nthreads.h"
71 #include "gromacs/mdlib/md_support.h"
72 #include "gromacs/mdlib/nb_verlet.h"
73 #include "gromacs/mdlib/nbnxn_atomdata.h"
74 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
75 #include "gromacs/mdlib/nbnxn_search.h"
76 #include "gromacs/mdlib/nbnxn_simd.h"
77 #include "gromacs/mdlib/nbnxn_util.h"
78 #include "gromacs/mdlib/ns.h"
79 #include "gromacs/mdlib/qmmm.h"
80 #include "gromacs/mdlib/sim_util.h"
81 #include "gromacs/mdtypes/commrec.h"
82 #include "gromacs/mdtypes/fcdata.h"
83 #include "gromacs/mdtypes/group.h"
84 #include "gromacs/mdtypes/inputrec.h"
85 #include "gromacs/mdtypes/md_enums.h"
86 #include "gromacs/pbcutil/ishift.h"
87 #include "gromacs/pbcutil/pbc.h"
88 #include "gromacs/simd/simd.h"
89 #include "gromacs/tables/forcetable.h"
90 #include "gromacs/topology/mtop_util.h"
91 #include "gromacs/trajectory/trajectoryframe.h"
92 #include "gromacs/utility/cstringutil.h"
93 #include "gromacs/utility/exceptions.h"
94 #include "gromacs/utility/fatalerror.h"
95 #include "gromacs/utility/gmxassert.h"
96 #include "gromacs/utility/pleasecite.h"
97 #include "gromacs/utility/smalloc.h"
98 #include "gromacs/utility/stringutil.h"
100 #include "nbnxn_gpu_jit_support.h"
102 const char *egrp_nm[egNR+1] = {
103 "Coul-SR", "LJ-SR", "Buck-SR",
104 "Coul-14", "LJ-14", NULL
107 t_forcerec *mk_forcerec(void)
117 static void pr_nbfp(FILE *fp, real *nbfp, gmx_bool bBHAM, int atnr)
121 for (i = 0; (i < atnr); i++)
123 for (j = 0; (j < atnr); j++)
125 fprintf(fp, "%2d - %2d", i, j);
128 fprintf(fp, " a=%10g, b=%10g, c=%10g\n", BHAMA(nbfp, atnr, i, j),
129 BHAMB(nbfp, atnr, i, j), BHAMC(nbfp, atnr, i, j)/6.0);
133 fprintf(fp, " c6=%10g, c12=%10g\n", C6(nbfp, atnr, i, j)/6.0,
134 C12(nbfp, atnr, i, j)/12.0);
141 static real *mk_nbfp(const gmx_ffparams_t *idef, gmx_bool bBHAM)
149 snew(nbfp, 3*atnr*atnr);
150 for (i = k = 0; (i < atnr); i++)
152 for (j = 0; (j < atnr); j++, k++)
154 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
155 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
156 /* nbfp now includes the 6.0 derivative prefactor */
157 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c*6.0;
163 snew(nbfp, 2*atnr*atnr);
164 for (i = k = 0; (i < atnr); i++)
166 for (j = 0; (j < atnr); j++, k++)
168 /* nbfp now includes the 6.0/12.0 derivative prefactors */
169 C6(nbfp, atnr, i, j) = idef->iparams[k].lj.c6*6.0;
170 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12*12.0;
178 static real *make_ljpme_c6grid(const gmx_ffparams_t *idef, t_forcerec *fr)
181 real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
184 /* For LJ-PME simulations, we correct the energies with the reciprocal space
185 * inside of the cut-off. To do this the non-bonded kernels needs to have
186 * access to the C6-values used on the reciprocal grid in pme.c
190 snew(grid, 2*atnr*atnr);
191 for (i = k = 0; (i < atnr); i++)
193 for (j = 0; (j < atnr); j++, k++)
195 c6i = idef->iparams[i*(atnr+1)].lj.c6;
196 c12i = idef->iparams[i*(atnr+1)].lj.c12;
197 c6j = idef->iparams[j*(atnr+1)].lj.c6;
198 c12j = idef->iparams[j*(atnr+1)].lj.c12;
199 c6 = std::sqrt(c6i * c6j);
200 if (fr->ljpme_combination_rule == eljpmeLB
201 && !gmx_numzero(c6) && !gmx_numzero(c12i) && !gmx_numzero(c12j))
203 sigmai = gmx::sixthroot(c12i / c6i);
204 sigmaj = gmx::sixthroot(c12j / c6j);
205 epsi = c6i * c6i / c12i;
206 epsj = c6j * c6j / c12j;
207 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5*(sigmai+sigmaj));
209 /* Store the elements at the same relative positions as C6 in nbfp in order
210 * to simplify access in the kernels
212 grid[2*(atnr*i+j)] = c6*6.0;
218 static real *mk_nbfp_combination_rule(const gmx_ffparams_t *idef, int comb_rule)
222 real c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
226 snew(nbfp, 2*atnr*atnr);
227 for (i = 0; i < atnr; ++i)
229 for (j = 0; j < atnr; ++j)
231 c6i = idef->iparams[i*(atnr+1)].lj.c6;
232 c12i = idef->iparams[i*(atnr+1)].lj.c12;
233 c6j = idef->iparams[j*(atnr+1)].lj.c6;
234 c12j = idef->iparams[j*(atnr+1)].lj.c12;
235 c6 = std::sqrt(c6i * c6j);
236 c12 = std::sqrt(c12i * c12j);
237 if (comb_rule == eCOMB_ARITHMETIC
238 && !gmx_numzero(c6) && !gmx_numzero(c12))
240 sigmai = gmx::sixthroot(c12i / c6i);
241 sigmaj = gmx::sixthroot(c12j / c6j);
242 epsi = c6i * c6i / c12i;
243 epsj = c6j * c6j / c12j;
244 c6 = std::sqrt(epsi * epsj) * gmx::power6(0.5*(sigmai+sigmaj));
245 c12 = std::sqrt(epsi * epsj) * gmx::power12(0.5*(sigmai+sigmaj));
247 C6(nbfp, atnr, i, j) = c6*6.0;
248 C12(nbfp, atnr, i, j) = c12*12.0;
254 /* This routine sets fr->solvent_opt to the most common solvent in the
255 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
256 * the fr->solvent_type array with the correct type (or esolNO).
258 * Charge groups that fulfill the conditions but are not identical to the
259 * most common one will be marked as esolNO in the solvent_type array.
261 * TIP3p is identical to SPC for these purposes, so we call it
262 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
264 * NOTE: QM particle should not
265 * become an optimized solvent. Not even if there is only one charge
275 } solvent_parameters_t;
278 check_solvent_cg(const gmx_moltype_t *molt,
281 const unsigned char *qm_grpnr,
282 const t_grps *qm_grps,
284 int *n_solvent_parameters,
285 solvent_parameters_t **solvent_parameters_p,
295 real tmp_charge[4] = { 0.0 }; /* init to zero to make gcc4.8 happy */
296 int tmp_vdwtype[4] = { 0 }; /* init to zero to make gcc4.8 happy */
299 solvent_parameters_t *solvent_parameters;
301 /* We use a list with parameters for each solvent type.
302 * Every time we discover a new molecule that fulfills the basic
303 * conditions for a solvent we compare with the previous entries
304 * in these lists. If the parameters are the same we just increment
305 * the counter for that type, and otherwise we create a new type
306 * based on the current molecule.
308 * Once we've finished going through all molecules we check which
309 * solvent is most common, and mark all those molecules while we
310 * clear the flag on all others.
313 solvent_parameters = *solvent_parameters_p;
315 /* Mark the cg first as non optimized */
318 /* Check if this cg has no exclusions with atoms in other charge groups
319 * and all atoms inside the charge group excluded.
320 * We only have 3 or 4 atom solvent loops.
322 if (GET_CGINFO_EXCL_INTER(cginfo) ||
323 !GET_CGINFO_EXCL_INTRA(cginfo))
328 /* Get the indices of the first atom in this charge group */
329 j0 = molt->cgs.index[cg0];
330 j1 = molt->cgs.index[cg0+1];
332 /* Number of atoms in our molecule */
338 "Moltype '%s': there are %d atoms in this charge group\n",
342 /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
345 if (nj < 3 || nj > 4)
350 /* Check if we are doing QM on this group */
352 if (qm_grpnr != NULL)
354 for (j = j0; j < j1 && !qm; j++)
356 qm = (qm_grpnr[j] < qm_grps->nr - 1);
359 /* Cannot use solvent optimization with QM */
365 atom = molt->atoms.atom;
367 /* Still looks like a solvent, time to check parameters */
369 /* If it is perturbed (free energy) we can't use the solvent loops,
370 * so then we just skip to the next molecule.
374 for (j = j0; j < j1 && !perturbed; j++)
376 perturbed = PERTURBED(atom[j]);
384 /* Now it's only a question if the VdW and charge parameters
385 * are OK. Before doing the check we compare and see if they are
386 * identical to a possible previous solvent type.
387 * First we assign the current types and charges.
389 for (j = 0; j < nj; j++)
391 tmp_vdwtype[j] = atom[j0+j].type;
392 tmp_charge[j] = atom[j0+j].q;
395 /* Does it match any previous solvent type? */
396 for (k = 0; k < *n_solvent_parameters; k++)
401 /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
402 if ( (solvent_parameters[k].model == esolSPC && nj != 3) ||
403 (solvent_parameters[k].model == esolTIP4P && nj != 4) )
408 /* Check that types & charges match for all atoms in molecule */
409 for (j = 0; j < nj && match == TRUE; j++)
411 if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
415 if (tmp_charge[j] != solvent_parameters[k].charge[j])
422 /* Congratulations! We have a matched solvent.
423 * Flag it with this type for later processing.
426 solvent_parameters[k].count += nmol;
428 /* We are done with this charge group */
433 /* If we get here, we have a tentative new solvent type.
434 * Before we add it we must check that it fulfills the requirements
435 * of the solvent optimized loops. First determine which atoms have
438 for (j = 0; j < nj; j++)
441 tjA = tmp_vdwtype[j];
443 /* Go through all other tpes and see if any have non-zero
444 * VdW parameters when combined with this one.
446 for (k = 0; k < fr->ntype && (has_vdw[j] == FALSE); k++)
448 /* We already checked that the atoms weren't perturbed,
449 * so we only need to check state A now.
453 has_vdw[j] = (has_vdw[j] ||
454 (BHAMA(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
455 (BHAMB(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
456 (BHAMC(fr->nbfp, fr->ntype, tjA, k) != 0.0));
461 has_vdw[j] = (has_vdw[j] ||
462 (C6(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
463 (C12(fr->nbfp, fr->ntype, tjA, k) != 0.0));
468 /* Now we know all we need to make the final check and assignment. */
472 * For this we require thatn all atoms have charge,
473 * the charges on atom 2 & 3 should be the same, and only
474 * atom 1 might have VdW.
476 if (has_vdw[1] == FALSE &&
477 has_vdw[2] == FALSE &&
478 tmp_charge[0] != 0 &&
479 tmp_charge[1] != 0 &&
480 tmp_charge[2] == tmp_charge[1])
482 srenew(solvent_parameters, *n_solvent_parameters+1);
483 solvent_parameters[*n_solvent_parameters].model = esolSPC;
484 solvent_parameters[*n_solvent_parameters].count = nmol;
485 for (k = 0; k < 3; k++)
487 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
488 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
491 *cg_sp = *n_solvent_parameters;
492 (*n_solvent_parameters)++;
497 /* Or could it be a TIP4P?
498 * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
499 * Only atom 1 mght have VdW.
501 if (has_vdw[1] == FALSE &&
502 has_vdw[2] == FALSE &&
503 has_vdw[3] == FALSE &&
504 tmp_charge[0] == 0 &&
505 tmp_charge[1] != 0 &&
506 tmp_charge[2] == tmp_charge[1] &&
509 srenew(solvent_parameters, *n_solvent_parameters+1);
510 solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
511 solvent_parameters[*n_solvent_parameters].count = nmol;
512 for (k = 0; k < 4; k++)
514 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
515 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
518 *cg_sp = *n_solvent_parameters;
519 (*n_solvent_parameters)++;
523 *solvent_parameters_p = solvent_parameters;
527 check_solvent(FILE * fp,
528 const gmx_mtop_t * mtop,
530 cginfo_mb_t *cginfo_mb)
533 const gmx_moltype_t *molt;
534 int mb, mol, cg_mol, at_offset, am, cgm, i, nmol_ch, nmol;
535 int n_solvent_parameters;
536 solvent_parameters_t *solvent_parameters;
542 fprintf(debug, "Going to determine what solvent types we have.\n");
545 n_solvent_parameters = 0;
546 solvent_parameters = NULL;
547 /* Allocate temporary array for solvent type */
548 snew(cg_sp, mtop->nmolblock);
551 for (mb = 0; mb < mtop->nmolblock; mb++)
553 molt = &mtop->moltype[mtop->molblock[mb].type];
555 /* Here we have to loop over all individual molecules
556 * because we need to check for QMMM particles.
558 snew(cg_sp[mb], cginfo_mb[mb].cg_mod);
559 nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
560 nmol = mtop->molblock[mb].nmol/nmol_ch;
561 for (mol = 0; mol < nmol_ch; mol++)
564 am = mol*cgs->index[cgs->nr];
565 for (cg_mol = 0; cg_mol < cgs->nr; cg_mol++)
567 check_solvent_cg(molt, cg_mol, nmol,
568 mtop->groups.grpnr[egcQMMM] ?
569 mtop->groups.grpnr[egcQMMM]+at_offset+am : 0,
570 &mtop->groups.grps[egcQMMM],
572 &n_solvent_parameters, &solvent_parameters,
573 cginfo_mb[mb].cginfo[cgm+cg_mol],
574 &cg_sp[mb][cgm+cg_mol]);
577 at_offset += cgs->index[cgs->nr];
580 /* Puh! We finished going through all charge groups.
581 * Now find the most common solvent model.
584 /* Most common solvent this far */
586 for (i = 0; i < n_solvent_parameters; i++)
589 solvent_parameters[i].count > solvent_parameters[bestsp].count)
597 bestsol = solvent_parameters[bestsp].model;
605 for (mb = 0; mb < mtop->nmolblock; mb++)
607 cgs = &mtop->moltype[mtop->molblock[mb].type].cgs;
608 nmol = (mtop->molblock[mb].nmol*cgs->nr)/cginfo_mb[mb].cg_mod;
609 for (i = 0; i < cginfo_mb[mb].cg_mod; i++)
611 if (cg_sp[mb][i] == bestsp)
613 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], bestsol);
618 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], esolNO);
625 if (bestsol != esolNO && fp != NULL)
627 fprintf(fp, "\nEnabling %s-like water optimization for %d molecules.\n\n",
629 solvent_parameters[bestsp].count);
632 sfree(solvent_parameters);
633 fr->solvent_opt = bestsol;
637 acNONE = 0, acCONSTRAINT, acSETTLE
640 static cginfo_mb_t *init_cginfo_mb(FILE *fplog, const gmx_mtop_t *mtop,
641 t_forcerec *fr, gmx_bool bNoSolvOpt,
642 gmx_bool *bFEP_NonBonded,
643 gmx_bool *bExcl_IntraCGAll_InterCGNone)
646 const t_blocka *excl;
647 const gmx_moltype_t *molt;
648 const gmx_molblock_t *molb;
649 cginfo_mb_t *cginfo_mb;
652 int cg_offset, a_offset;
653 int mb, m, cg, a0, a1, gid, ai, j, aj, excl_nalloc;
657 gmx_bool bId, *bExcl, bExclIntraAll, bExclInter, bHaveVDW, bHaveQ, bHavePerturbedAtoms;
659 snew(cginfo_mb, mtop->nmolblock);
661 snew(type_VDW, fr->ntype);
662 for (ai = 0; ai < fr->ntype; ai++)
664 type_VDW[ai] = FALSE;
665 for (j = 0; j < fr->ntype; j++)
667 type_VDW[ai] = type_VDW[ai] ||
669 C6(fr->nbfp, fr->ntype, ai, j) != 0 ||
670 C12(fr->nbfp, fr->ntype, ai, j) != 0;
674 *bFEP_NonBonded = FALSE;
675 *bExcl_IntraCGAll_InterCGNone = TRUE;
678 snew(bExcl, excl_nalloc);
681 for (mb = 0; mb < mtop->nmolblock; mb++)
683 molb = &mtop->molblock[mb];
684 molt = &mtop->moltype[molb->type];
688 /* Check if the cginfo is identical for all molecules in this block.
689 * If so, we only need an array of the size of one molecule.
690 * Otherwise we make an array of #mol times #cgs per molecule.
693 for (m = 0; m < molb->nmol; m++)
695 int am = m*cgs->index[cgs->nr];
696 for (cg = 0; cg < cgs->nr; cg++)
699 a1 = cgs->index[cg+1];
700 if (ggrpnr(&mtop->groups, egcENER, a_offset+am+a0) !=
701 ggrpnr(&mtop->groups, egcENER, a_offset +a0))
705 if (mtop->groups.grpnr[egcQMMM] != NULL)
707 for (ai = a0; ai < a1; ai++)
709 if (mtop->groups.grpnr[egcQMMM][a_offset+am+ai] !=
710 mtop->groups.grpnr[egcQMMM][a_offset +ai])
719 cginfo_mb[mb].cg_start = cg_offset;
720 cginfo_mb[mb].cg_end = cg_offset + molb->nmol*cgs->nr;
721 cginfo_mb[mb].cg_mod = (bId ? 1 : molb->nmol)*cgs->nr;
722 snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
723 cginfo = cginfo_mb[mb].cginfo;
725 /* Set constraints flags for constrained atoms */
726 snew(a_con, molt->atoms.nr);
727 for (ftype = 0; ftype < F_NRE; ftype++)
729 if (interaction_function[ftype].flags & IF_CONSTRAINT)
734 for (ia = 0; ia < molt->ilist[ftype].nr; ia += 1+nral)
738 for (a = 0; a < nral; a++)
740 a_con[molt->ilist[ftype].iatoms[ia+1+a]] =
741 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
747 for (m = 0; m < (bId ? 1 : molb->nmol); m++)
750 int am = m*cgs->index[cgs->nr];
751 for (cg = 0; cg < cgs->nr; cg++)
754 a1 = cgs->index[cg+1];
756 /* Store the energy group in cginfo */
757 gid = ggrpnr(&mtop->groups, egcENER, a_offset+am+a0);
758 SET_CGINFO_GID(cginfo[cgm+cg], gid);
760 /* Check the intra/inter charge group exclusions */
761 if (a1-a0 > excl_nalloc)
763 excl_nalloc = a1 - a0;
764 srenew(bExcl, excl_nalloc);
766 /* bExclIntraAll: all intra cg interactions excluded
767 * bExclInter: any inter cg interactions excluded
769 bExclIntraAll = TRUE;
773 bHavePerturbedAtoms = FALSE;
774 for (ai = a0; ai < a1; ai++)
776 /* Check VDW and electrostatic interactions */
777 bHaveVDW = bHaveVDW || (type_VDW[molt->atoms.atom[ai].type] ||
778 type_VDW[molt->atoms.atom[ai].typeB]);
779 bHaveQ = bHaveQ || (molt->atoms.atom[ai].q != 0 ||
780 molt->atoms.atom[ai].qB != 0);
782 bHavePerturbedAtoms = bHavePerturbedAtoms || (PERTURBED(molt->atoms.atom[ai]) != 0);
784 /* Clear the exclusion list for atom ai */
785 for (aj = a0; aj < a1; aj++)
787 bExcl[aj-a0] = FALSE;
789 /* Loop over all the exclusions of atom ai */
790 for (j = excl->index[ai]; j < excl->index[ai+1]; j++)
793 if (aj < a0 || aj >= a1)
802 /* Check if ai excludes a0 to a1 */
803 for (aj = a0; aj < a1; aj++)
807 bExclIntraAll = FALSE;
814 SET_CGINFO_CONSTR(cginfo[cgm+cg]);
817 SET_CGINFO_SETTLE(cginfo[cgm+cg]);
825 SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
829 SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
831 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
833 /* The size in cginfo is currently only read with DD */
834 gmx_fatal(FARGS, "A charge group has size %d which is larger than the limit of %d atoms", a1-a0, MAX_CHARGEGROUP_SIZE);
838 SET_CGINFO_HAS_VDW(cginfo[cgm+cg]);
842 SET_CGINFO_HAS_Q(cginfo[cgm+cg]);
844 if (bHavePerturbedAtoms && fr->efep != efepNO)
846 SET_CGINFO_FEP(cginfo[cgm+cg]);
847 *bFEP_NonBonded = TRUE;
849 /* Store the charge group size */
850 SET_CGINFO_NATOMS(cginfo[cgm+cg], a1-a0);
852 if (!bExclIntraAll || bExclInter)
854 *bExcl_IntraCGAll_InterCGNone = FALSE;
861 cg_offset += molb->nmol*cgs->nr;
862 a_offset += molb->nmol*cgs->index[cgs->nr];
866 /* the solvent optimizer is called after the QM is initialized,
867 * because we don't want to have the QM subsystemto become an
871 check_solvent(fplog, mtop, fr, cginfo_mb);
873 if (getenv("GMX_NO_SOLV_OPT"))
877 fprintf(fplog, "Found environment variable GMX_NO_SOLV_OPT.\n"
878 "Disabling all solvent optimization\n");
880 fr->solvent_opt = esolNO;
884 fr->solvent_opt = esolNO;
886 if (!fr->solvent_opt)
888 for (mb = 0; mb < mtop->nmolblock; mb++)
890 for (cg = 0; cg < cginfo_mb[mb].cg_mod; cg++)
892 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg], esolNO);
900 static int *cginfo_expand(int nmb, cginfo_mb_t *cgi_mb)
905 ncg = cgi_mb[nmb-1].cg_end;
908 for (cg = 0; cg < ncg; cg++)
910 while (cg >= cgi_mb[mb].cg_end)
915 cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
921 static void set_chargesum(FILE *log, t_forcerec *fr, const gmx_mtop_t *mtop)
923 /*This now calculates sum for q and c6*/
924 double qsum, q2sum, q, c6sum, c6;
926 const t_atoms *atoms;
931 for (mb = 0; mb < mtop->nmolblock; mb++)
933 nmol = mtop->molblock[mb].nmol;
934 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
935 for (i = 0; i < atoms->nr; i++)
937 q = atoms->atom[i].q;
940 c6 = mtop->ffparams.iparams[atoms->atom[i].type*(mtop->ffparams.atnr+1)].lj.c6;
945 fr->q2sum[0] = q2sum;
946 fr->c6sum[0] = c6sum;
948 if (fr->efep != efepNO)
953 for (mb = 0; mb < mtop->nmolblock; mb++)
955 nmol = mtop->molblock[mb].nmol;
956 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
957 for (i = 0; i < atoms->nr; i++)
959 q = atoms->atom[i].qB;
962 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB*(mtop->ffparams.atnr+1)].lj.c6;
966 fr->q2sum[1] = q2sum;
967 fr->c6sum[1] = c6sum;
972 fr->qsum[1] = fr->qsum[0];
973 fr->q2sum[1] = fr->q2sum[0];
974 fr->c6sum[1] = fr->c6sum[0];
978 if (fr->efep == efepNO)
980 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
984 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n",
985 fr->qsum[0], fr->qsum[1]);
990 void update_forcerec(t_forcerec *fr, matrix box)
992 if (fr->eeltype == eelGRF)
994 calc_rffac(NULL, fr->eeltype, fr->epsilon_r, fr->epsilon_rf,
995 fr->rcoulomb, fr->temp, fr->zsquare, box,
996 &fr->kappa, &fr->k_rf, &fr->c_rf);
1000 void set_avcsixtwelve(FILE *fplog, t_forcerec *fr, const gmx_mtop_t *mtop)
1002 const t_atoms *atoms, *atoms_tpi;
1003 const t_blocka *excl;
1004 int mb, nmol, nmolc, i, j, tpi, tpj, j1, j2, k, nexcl, q;
1005 gmx_int64_t npair, npair_ij, tmpi, tmpj;
1006 double csix, ctwelve;
1007 int ntp, *typecount;
1010 real *nbfp_comb = NULL;
1016 /* For LJ-PME, we want to correct for the difference between the
1017 * actual C6 values and the C6 values used by the LJ-PME based on
1018 * combination rules. */
1020 if (EVDW_PME(fr->vdwtype))
1022 nbfp_comb = mk_nbfp_combination_rule(&mtop->ffparams,
1023 (fr->ljpme_combination_rule == eljpmeLB) ? eCOMB_ARITHMETIC : eCOMB_GEOMETRIC);
1024 for (tpi = 0; tpi < ntp; ++tpi)
1026 for (tpj = 0; tpj < ntp; ++tpj)
1028 C6(nbfp_comb, ntp, tpi, tpj) =
1029 C6(nbfp, ntp, tpi, tpj) - C6(nbfp_comb, ntp, tpi, tpj);
1030 C12(nbfp_comb, ntp, tpi, tpj) = C12(nbfp, ntp, tpi, tpj);
1035 for (q = 0; q < (fr->efep == efepNO ? 1 : 2); q++)
1043 /* Count the types so we avoid natoms^2 operations */
1044 snew(typecount, ntp);
1045 gmx_mtop_count_atomtypes(mtop, q, typecount);
1047 for (tpi = 0; tpi < ntp; tpi++)
1049 for (tpj = tpi; tpj < ntp; tpj++)
1051 tmpi = typecount[tpi];
1052 tmpj = typecount[tpj];
1055 npair_ij = tmpi*tmpj;
1059 npair_ij = tmpi*(tmpi - 1)/2;
1063 /* nbfp now includes the 6.0 derivative prefactor */
1064 csix += npair_ij*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1068 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1069 csix += npair_ij* C6(nbfp, ntp, tpi, tpj)/6.0;
1070 ctwelve += npair_ij* C12(nbfp, ntp, tpi, tpj)/12.0;
1076 /* Subtract the excluded pairs.
1077 * The main reason for substracting exclusions is that in some cases
1078 * some combinations might never occur and the parameters could have
1079 * any value. These unused values should not influence the dispersion
1082 for (mb = 0; mb < mtop->nmolblock; mb++)
1084 nmol = mtop->molblock[mb].nmol;
1085 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
1086 excl = &mtop->moltype[mtop->molblock[mb].type].excls;
1087 for (i = 0; (i < atoms->nr); i++)
1091 tpi = atoms->atom[i].type;
1095 tpi = atoms->atom[i].typeB;
1097 j1 = excl->index[i];
1098 j2 = excl->index[i+1];
1099 for (j = j1; j < j2; j++)
1106 tpj = atoms->atom[k].type;
1110 tpj = atoms->atom[k].typeB;
1114 /* nbfp now includes the 6.0 derivative prefactor */
1115 csix -= nmol*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1119 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1120 csix -= nmol*C6 (nbfp, ntp, tpi, tpj)/6.0;
1121 ctwelve -= nmol*C12(nbfp, ntp, tpi, tpj)/12.0;
1131 /* Only correct for the interaction of the test particle
1132 * with the rest of the system.
1135 &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].atoms;
1138 for (mb = 0; mb < mtop->nmolblock; mb++)
1140 nmol = mtop->molblock[mb].nmol;
1141 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
1142 for (j = 0; j < atoms->nr; j++)
1145 /* Remove the interaction of the test charge group
1148 if (mb == mtop->nmolblock-1)
1152 if (mb == 0 && nmol == 1)
1154 gmx_fatal(FARGS, "Old format tpr with TPI, please generate a new tpr file");
1159 tpj = atoms->atom[j].type;
1163 tpj = atoms->atom[j].typeB;
1165 for (i = 0; i < fr->n_tpi; i++)
1169 tpi = atoms_tpi->atom[i].type;
1173 tpi = atoms_tpi->atom[i].typeB;
1177 /* nbfp now includes the 6.0 derivative prefactor */
1178 csix += nmolc*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1182 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1183 csix += nmolc*C6 (nbfp, ntp, tpi, tpj)/6.0;
1184 ctwelve += nmolc*C12(nbfp, ntp, tpi, tpj)/12.0;
1191 if (npair - nexcl <= 0 && fplog)
1193 fprintf(fplog, "\nWARNING: There are no atom pairs for dispersion correction\n\n");
1199 csix /= npair - nexcl;
1200 ctwelve /= npair - nexcl;
1204 fprintf(debug, "Counted %d exclusions\n", nexcl);
1205 fprintf(debug, "Average C6 parameter is: %10g\n", (double)csix);
1206 fprintf(debug, "Average C12 parameter is: %10g\n", (double)ctwelve);
1208 fr->avcsix[q] = csix;
1209 fr->avctwelve[q] = ctwelve;
1212 if (EVDW_PME(fr->vdwtype))
1219 if (fr->eDispCorr == edispcAllEner ||
1220 fr->eDispCorr == edispcAllEnerPres)
1222 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1223 fr->avcsix[0], fr->avctwelve[0]);
1227 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e\n", fr->avcsix[0]);
1233 static void set_bham_b_max(FILE *fplog, t_forcerec *fr,
1234 const gmx_mtop_t *mtop)
1236 const t_atoms *at1, *at2;
1237 int mt1, mt2, i, j, tpi, tpj, ntypes;
1243 fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
1250 for (mt1 = 0; mt1 < mtop->nmoltype; mt1++)
1252 at1 = &mtop->moltype[mt1].atoms;
1253 for (i = 0; (i < at1->nr); i++)
1255 tpi = at1->atom[i].type;
1258 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
1261 for (mt2 = mt1; mt2 < mtop->nmoltype; mt2++)
1263 at2 = &mtop->moltype[mt2].atoms;
1264 for (j = 0; (j < at2->nr); j++)
1266 tpj = at2->atom[j].type;
1269 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
1271 b = BHAMB(nbfp, ntypes, tpi, tpj);
1272 if (b > fr->bham_b_max)
1276 if ((b < bmin) || (bmin == -1))
1286 fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n",
1287 bmin, fr->bham_b_max);
1291 static void make_nbf_tables(FILE *fp,
1292 t_forcerec *fr, real rtab,
1293 const char *tabfn, char *eg1, char *eg2,
1303 fprintf(debug, "No table file name passed, can not read table, can not do non-bonded interactions\n");
1308 sprintf(buf, "%s", tabfn);
1311 /* Append the two energy group names */
1312 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
1313 eg1, eg2, ftp2ext(efXVG));
1315 nbl->table_elec_vdw = make_tables(fp, fr, buf, rtab, 0);
1316 /* Copy the contents of the table to separate coulomb and LJ tables too,
1317 * to improve cache performance.
1319 /* For performance reasons we want
1320 * the table data to be aligned to 16-byte. The pointers could be freed
1321 * but currently aren't.
1323 snew(nbl->table_elec, 1);
1324 nbl->table_elec->interaction = GMX_TABLE_INTERACTION_ELEC;
1325 nbl->table_elec->format = nbl->table_elec_vdw->format;
1326 nbl->table_elec->r = nbl->table_elec_vdw->r;
1327 nbl->table_elec->n = nbl->table_elec_vdw->n;
1328 nbl->table_elec->scale = nbl->table_elec_vdw->scale;
1329 nbl->table_elec->formatsize = nbl->table_elec_vdw->formatsize;
1330 nbl->table_elec->ninteractions = 1;
1331 nbl->table_elec->stride = nbl->table_elec->formatsize * nbl->table_elec->ninteractions;
1332 snew_aligned(nbl->table_elec->data, nbl->table_elec->stride*(nbl->table_elec->n+1), 32);
1334 snew(nbl->table_vdw, 1);
1335 nbl->table_vdw->interaction = GMX_TABLE_INTERACTION_VDWREP_VDWDISP;
1336 nbl->table_vdw->format = nbl->table_elec_vdw->format;
1337 nbl->table_vdw->r = nbl->table_elec_vdw->r;
1338 nbl->table_vdw->n = nbl->table_elec_vdw->n;
1339 nbl->table_vdw->scale = nbl->table_elec_vdw->scale;
1340 nbl->table_vdw->formatsize = nbl->table_elec_vdw->formatsize;
1341 nbl->table_vdw->ninteractions = 2;
1342 nbl->table_vdw->stride = nbl->table_vdw->formatsize * nbl->table_vdw->ninteractions;
1343 snew_aligned(nbl->table_vdw->data, nbl->table_vdw->stride*(nbl->table_vdw->n+1), 32);
1345 for (i = 0; i <= nbl->table_elec_vdw->n; i++)
1347 for (j = 0; j < 4; j++)
1349 nbl->table_elec->data[4*i+j] = nbl->table_elec_vdw->data[12*i+j];
1351 for (j = 0; j < 8; j++)
1353 nbl->table_vdw->data[8*i+j] = nbl->table_elec_vdw->data[12*i+4+j];
1358 /*!\brief If there's bonded interactions of type \c ftype1 or \c
1359 * ftype2 present in the topology, build an array of the number of
1360 * interactions present for each bonded interaction index found in the
1363 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1364 * valid type with that parameter.
1366 * \c count will be reallocated as necessary to fit the largest bonded
1367 * interaction index found, and its current size will be returned in
1368 * \c ncount. It will contain zero for every bonded interaction index
1369 * for which no interactions are present in the topology.
1371 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t *mtop,
1372 int *ncount, int **count)
1374 const gmx_moltype_t *molt;
1376 int mt, ftype, stride, i, j, tabnr;
1378 // Loop over all moleculetypes
1379 for (mt = 0; mt < mtop->nmoltype; mt++)
1381 molt = &mtop->moltype[mt];
1382 // Loop over all interaction types
1383 for (ftype = 0; ftype < F_NRE; ftype++)
1385 // If the current interaction type is one of the types whose tables we're trying to count...
1386 if (ftype == ftype1 || ftype == ftype2)
1388 il = &molt->ilist[ftype];
1389 stride = 1 + NRAL(ftype);
1390 // ... and there are actually some interactions for this type
1391 for (i = 0; i < il->nr; i += stride)
1393 // Find out which table index the user wanted
1394 tabnr = mtop->ffparams.iparams[il->iatoms[i]].tab.table;
1397 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
1399 // Make room for this index in the data structure
1400 if (tabnr >= *ncount)
1402 srenew(*count, tabnr+1);
1403 for (j = *ncount; j < tabnr+1; j++)
1409 // Record that this table index is used and must have a valid file
1417 /*!\brief If there's bonded interactions of flavour \c tabext and type
1418 * \c ftype1 or \c ftype2 present in the topology, seek them in the
1419 * list of filenames passed to mdrun, and make bonded tables from
1422 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1423 * valid type with that parameter.
1425 * A fatal error occurs if no matching filename is found.
1427 static bondedtable_t *make_bonded_tables(FILE *fplog,
1428 int ftype1, int ftype2,
1429 const gmx_mtop_t *mtop,
1430 const t_filenm *tabbfnm,
1440 count_tables(ftype1, ftype2, mtop, &ncount, &count);
1442 // Are there any relevant tabulated bond interactions?
1446 for (int i = 0; i < ncount; i++)
1448 // Do any interactions exist that requires this table?
1451 // This pattern enforces the current requirement that
1452 // table filenames end in a characteristic sequence
1453 // before the file type extension, and avoids table 13
1454 // being recognized and used for table 1.
1455 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
1456 bool madeTable = false;
1457 for (int j = 0; j < tabbfnm->nfiles && !madeTable; ++j)
1459 std::string filename(tabbfnm->fns[j]);
1460 if (gmx::endsWith(filename, patternToFind))
1462 // Finally read the table from the file found
1463 tab[i] = make_bonded_table(fplog, tabbfnm->fns[j], NRAL(ftype1)-2);
1469 bool isPlural = (ftype2 != -1);
1470 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.",
1471 interaction_function[ftype1].longname,
1472 isPlural ? "' or '" : "",
1473 isPlural ? interaction_function[ftype2].longname : "",
1475 patternToFind.c_str());
1485 void forcerec_set_ranges(t_forcerec *fr,
1486 int ncg_home, int ncg_force,
1488 int natoms_force_constr, int natoms_f_novirsum)
1493 /* fr->ncg_force is unused in the standard code,
1494 * but it can be useful for modified code dealing with charge groups.
1496 fr->ncg_force = ncg_force;
1497 fr->natoms_force = natoms_force;
1498 fr->natoms_force_constr = natoms_force_constr;
1500 if (fr->natoms_force_constr > fr->nalloc_force)
1502 fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1505 if (fr->bF_NoVirSum)
1507 fr->f_novirsum_n = natoms_f_novirsum;
1508 if (fr->f_novirsum_n > fr->f_novirsum_nalloc)
1510 fr->f_novirsum_nalloc = over_alloc_dd(fr->f_novirsum_n);
1511 srenew(fr->f_novirsum_alloc, fr->f_novirsum_nalloc);
1516 fr->f_novirsum_n = 0;
1520 static real cutoff_inf(real cutoff)
1524 cutoff = GMX_CUTOFF_INF;
1530 gmx_bool can_use_allvsall(const t_inputrec *ir, gmx_bool bPrintNote, t_commrec *cr, FILE *fp)
1537 ir->rcoulomb == 0 &&
1539 ir->ePBC == epbcNONE &&
1540 ir->vdwtype == evdwCUT &&
1541 ir->coulombtype == eelCUT &&
1542 ir->efep == efepNO &&
1543 (ir->implicit_solvent == eisNO ||
1544 (ir->implicit_solvent == eisGBSA && (ir->gb_algorithm == egbSTILL ||
1545 ir->gb_algorithm == egbHCT ||
1546 ir->gb_algorithm == egbOBC))) &&
1547 getenv("GMX_NO_ALLVSALL") == NULL
1550 if (bAllvsAll && ir->opts.ngener > 1)
1552 const char *note = "NOTE: Can not use all-vs-all force loops, because there are multiple energy monitor groups; you might get significantly higher performance when using only a single energy monitor group.\n";
1558 fprintf(stderr, "\n%s\n", note);
1562 fprintf(fp, "\n%s\n", note);
1568 if (bAllvsAll && fp && MASTER(cr))
1570 fprintf(fp, "\nUsing SIMD all-vs-all kernels.\n\n");
1577 gmx_bool nbnxn_gpu_acceleration_supported(FILE *fplog,
1578 const t_commrec *cr,
1579 const t_inputrec *ir,
1582 if (bRerunMD && ir->opts.ngener > 1)
1584 /* Rerun execution time is dominated by I/O and pair search,
1585 * so GPUs are not very useful, plus they do not support more
1586 * than one energy group. If the user requested GPUs
1587 * explicitly, a fatal error is given later. With non-reruns,
1588 * we fall back to a single whole-of system energy group
1589 * (which runs much faster than a multiple-energy-groups
1590 * implementation would), and issue a note in the .log
1591 * file. Users can re-run if they want the information. */
1592 md_print_warn(cr, fplog, "Rerun with energy groups is not implemented for GPUs, falling back to the CPU\n");
1599 gmx_bool nbnxn_simd_supported(FILE *fplog,
1600 const t_commrec *cr,
1601 const t_inputrec *ir)
1603 if (ir->vdwtype == evdwPME && ir->ljpme_combination_rule == eljpmeLB)
1605 /* LJ PME with LB combination rule does 7 mesh operations.
1606 * This so slow that we don't compile SIMD non-bonded kernels
1608 md_print_warn(cr, fplog, "LJ-PME with Lorentz-Berthelot is not supported with SIMD kernels, falling back to plain C kernels\n");
1616 static void pick_nbnxn_kernel_cpu(const t_inputrec gmx_unused *ir,
1620 *kernel_type = nbnxnk4x4_PlainC;
1621 *ewald_excl = ewaldexclTable;
1625 #ifdef GMX_NBNXN_SIMD_4XN
1626 *kernel_type = nbnxnk4xN_SIMD_4xN;
1628 #ifdef GMX_NBNXN_SIMD_2XNN
1629 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1632 #if defined GMX_NBNXN_SIMD_2XNN && defined GMX_NBNXN_SIMD_4XN
1633 /* We need to choose if we want 2x(N+N) or 4xN kernels.
1634 * Currently this is based on the SIMD acceleration choice,
1635 * but it might be better to decide this at runtime based on CPU.
1637 * 4xN calculates more (zero) interactions, but has less pair-search
1638 * work and much better kernel instruction scheduling.
1640 * Up till now we have only seen that on Intel Sandy/Ivy Bridge,
1641 * which doesn't have FMA, both the analytical and tabulated Ewald
1642 * kernels have similar pair rates for 4x8 and 2x(4+4), so we choose
1643 * 2x(4+4) because it results in significantly fewer pairs.
1644 * For RF, the raw pair rate of the 4x8 kernel is higher than 2x(4+4),
1645 * 10% with HT, 50% without HT. As we currently don't detect the actual
1646 * use of HT, use 4x8 to avoid a potential performance hit.
1647 * On Intel Haswell 4x8 is always faster.
1649 *kernel_type = nbnxnk4xN_SIMD_4xN;
1651 #if !GMX_SIMD_HAVE_FMA
1652 if (EEL_PME_EWALD(ir->coulombtype) ||
1653 EVDW_PME(ir->vdwtype))
1655 /* We have Ewald kernels without FMA (Intel Sandy/Ivy Bridge).
1656 * There are enough instructions to make 2x(4+4) efficient.
1658 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1661 #endif /* GMX_NBNXN_SIMD_2XNN && GMX_NBNXN_SIMD_4XN */
1664 if (getenv("GMX_NBNXN_SIMD_4XN") != NULL)
1666 #ifdef GMX_NBNXN_SIMD_4XN
1667 *kernel_type = nbnxnk4xN_SIMD_4xN;
1669 gmx_fatal(FARGS, "SIMD 4xN kernels requested, but GROMACS has been compiled without support for these kernels");
1672 if (getenv("GMX_NBNXN_SIMD_2XNN") != NULL)
1674 #ifdef GMX_NBNXN_SIMD_2XNN
1675 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1677 gmx_fatal(FARGS, "SIMD 2x(N+N) kernels requested, but GROMACS has been compiled without support for these kernels");
1681 /* Analytical Ewald exclusion correction is only an option in
1683 * Since table lookup's don't parallelize with SIMD, analytical
1684 * will probably always be faster for a SIMD width of 8 or more.
1685 * With FMA analytical is sometimes faster for a width if 4 as well.
1686 * On BlueGene/Q, this is faster regardless of precision.
1687 * In single precision, this is faster on Bulldozer.
1689 #if GMX_SIMD_REAL_WIDTH >= 8 || \
1690 (GMX_SIMD_REAL_WIDTH >= 4 && GMX_SIMD_HAVE_FMA && !GMX_DOUBLE) || GMX_SIMD_IBM_QPX
1691 *ewald_excl = ewaldexclAnalytical;
1693 if (getenv("GMX_NBNXN_EWALD_TABLE") != NULL)
1695 *ewald_excl = ewaldexclTable;
1697 if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != NULL)
1699 *ewald_excl = ewaldexclAnalytical;
1707 const char *lookup_nbnxn_kernel_name(int kernel_type)
1709 const char *returnvalue = NULL;
1710 switch (kernel_type)
1713 returnvalue = "not set";
1715 case nbnxnk4x4_PlainC:
1716 returnvalue = "plain C";
1718 case nbnxnk4xN_SIMD_4xN:
1719 case nbnxnk4xN_SIMD_2xNN:
1721 returnvalue = "SIMD";
1723 returnvalue = "not available";
1726 case nbnxnk8x8x8_GPU: returnvalue = "GPU"; break;
1727 case nbnxnk8x8x8_PlainC: returnvalue = "plain C"; break;
1731 gmx_fatal(FARGS, "Illegal kernel type selected");
1738 static void pick_nbnxn_kernel(FILE *fp,
1739 const t_commrec *cr,
1740 gmx_bool use_simd_kernels,
1742 gmx_bool bEmulateGPU,
1743 const t_inputrec *ir,
1746 gmx_bool bDoNonbonded)
1748 assert(kernel_type);
1750 *kernel_type = nbnxnkNotSet;
1751 *ewald_excl = ewaldexclTable;
1755 *kernel_type = nbnxnk8x8x8_PlainC;
1759 md_print_warn(cr, fp, "Emulating a GPU run on the CPU (slow)");
1764 *kernel_type = nbnxnk8x8x8_GPU;
1767 if (*kernel_type == nbnxnkNotSet)
1769 if (use_simd_kernels &&
1770 nbnxn_simd_supported(fp, cr, ir))
1772 pick_nbnxn_kernel_cpu(ir, kernel_type, ewald_excl);
1776 *kernel_type = nbnxnk4x4_PlainC;
1780 if (bDoNonbonded && fp != NULL)
1782 fprintf(fp, "\nUsing %s %dx%d non-bonded kernels\n\n",
1783 lookup_nbnxn_kernel_name(*kernel_type),
1784 nbnxn_kernel_to_cluster_i_size(*kernel_type),
1785 nbnxn_kernel_to_cluster_j_size(*kernel_type));
1787 if (nbnxnk4x4_PlainC == *kernel_type ||
1788 nbnxnk8x8x8_PlainC == *kernel_type)
1790 md_print_warn(cr, fp,
1791 "WARNING: Using the slow %s kernels. This should\n"
1792 "not happen during routine usage on supported platforms.\n\n",
1793 lookup_nbnxn_kernel_name(*kernel_type));
1798 static void pick_nbnxn_resources(FILE *fp,
1799 const t_commrec *cr,
1800 const gmx_hw_info_t *hwinfo,
1801 gmx_bool bDoNonbonded,
1803 gmx_bool *bEmulateGPU,
1804 const gmx_gpu_opt_t *gpu_opt)
1806 gmx_bool bEmulateGPUEnvVarSet;
1807 char gpu_err_str[STRLEN];
1811 bEmulateGPUEnvVarSet = (getenv("GMX_EMULATE_GPU") != NULL);
1813 /* Run GPU emulation mode if GMX_EMULATE_GPU is defined. Because
1814 * GPUs (currently) only handle non-bonded calculations, we will
1815 * automatically switch to emulation if non-bonded calculations are
1816 * turned off via GMX_NO_NONBONDED - this is the simple and elegant
1817 * way to turn off GPU initialization, data movement, and cleanup.
1819 * GPU emulation can be useful to assess the performance one can expect by
1820 * adding GPU(s) to the machine. The conditional below allows this even
1821 * if mdrun is compiled without GPU acceleration support.
1822 * Note that you should freezing the system as otherwise it will explode.
1824 *bEmulateGPU = (bEmulateGPUEnvVarSet ||
1825 (!bDoNonbonded && gpu_opt->n_dev_use > 0));
1827 /* Enable GPU mode when GPUs are available or no GPU emulation is requested.
1829 if (gpu_opt->n_dev_use > 0 && !(*bEmulateGPU))
1831 /* Each PP node will use the intra-node id-th device from the
1832 * list of detected/selected GPUs. */
1833 if (!init_gpu(fp, cr->rank_pp_intranode, gpu_err_str,
1834 &hwinfo->gpu_info, gpu_opt))
1836 /* At this point the init should never fail as we made sure that
1837 * we have all the GPUs we need. If it still does, we'll bail. */
1838 /* TODO the decorating of gpu_err_str is nicer if it
1839 happens inside init_gpu. Out here, the decorating with
1840 the MPI rank makes sense. */
1841 gmx_fatal(FARGS, "On rank %d failed to initialize GPU #%d: %s",
1843 get_gpu_device_id(&hwinfo->gpu_info, gpu_opt,
1844 cr->rank_pp_intranode),
1848 /* Here we actually turn on hardware GPU acceleration */
1853 gmx_bool uses_simple_tables(int cutoff_scheme,
1854 nonbonded_verlet_t *nbv,
1857 gmx_bool bUsesSimpleTables = TRUE;
1860 switch (cutoff_scheme)
1863 bUsesSimpleTables = TRUE;
1866 assert(NULL != nbv && NULL != nbv->grp);
1867 grp_index = (group < 0) ? 0 : (nbv->ngrp - 1);
1868 bUsesSimpleTables = nbnxn_kernel_pairlist_simple(nbv->grp[grp_index].kernel_type);
1871 gmx_incons("unimplemented");
1873 return bUsesSimpleTables;
1876 static void init_ewald_f_table(interaction_const_t *ic,
1881 /* Get the Ewald table spacing based on Coulomb and/or LJ
1882 * Ewald coefficients and rtol.
1884 ic->tabq_scale = ewald_spline3_table_scale(ic);
1886 if (ic->cutoff_scheme == ecutsVERLET)
1888 maxr = ic->rcoulomb;
1892 maxr = std::max(ic->rcoulomb, rtab);
1894 ic->tabq_size = static_cast<int>(maxr*ic->tabq_scale) + 2;
1896 sfree_aligned(ic->tabq_coul_FDV0);
1897 sfree_aligned(ic->tabq_coul_F);
1898 sfree_aligned(ic->tabq_coul_V);
1900 sfree_aligned(ic->tabq_vdw_FDV0);
1901 sfree_aligned(ic->tabq_vdw_F);
1902 sfree_aligned(ic->tabq_vdw_V);
1904 if (EEL_PME_EWALD(ic->eeltype))
1906 /* Create the original table data in FDV0 */
1907 snew_aligned(ic->tabq_coul_FDV0, ic->tabq_size*4, 32);
1908 snew_aligned(ic->tabq_coul_F, ic->tabq_size, 32);
1909 snew_aligned(ic->tabq_coul_V, ic->tabq_size, 32);
1910 table_spline3_fill_ewald_lr(ic->tabq_coul_F, ic->tabq_coul_V, ic->tabq_coul_FDV0,
1911 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_q, v_q_ewald_lr);
1914 if (EVDW_PME(ic->vdwtype))
1916 snew_aligned(ic->tabq_vdw_FDV0, ic->tabq_size*4, 32);
1917 snew_aligned(ic->tabq_vdw_F, ic->tabq_size, 32);
1918 snew_aligned(ic->tabq_vdw_V, ic->tabq_size, 32);
1919 table_spline3_fill_ewald_lr(ic->tabq_vdw_F, ic->tabq_vdw_V, ic->tabq_vdw_FDV0,
1920 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_lj, v_lj_ewald_lr);
1924 void init_interaction_const_tables(FILE *fp,
1925 interaction_const_t *ic,
1928 if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
1930 init_ewald_f_table(ic, rtab);
1934 fprintf(fp, "Initialized non-bonded Ewald correction tables, spacing: %.2e size: %d\n\n",
1935 1/ic->tabq_scale, ic->tabq_size);
1940 static void clear_force_switch_constants(shift_consts_t *sc)
1947 static void force_switch_constants(real p,
1951 /* Here we determine the coefficient for shifting the force to zero
1952 * between distance rsw and the cut-off rc.
1953 * For a potential of r^-p, we have force p*r^-(p+1).
1954 * But to save flops we absorb p in the coefficient.
1956 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
1957 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
1959 sc->c2 = ((p + 1)*rsw - (p + 4)*rc)/(pow(rc, p + 2)*gmx::square(rc - rsw));
1960 sc->c3 = -((p + 1)*rsw - (p + 3)*rc)/(pow(rc, p + 2)*gmx::power3(rc - rsw));
1961 sc->cpot = -pow(rc, -p) + p*sc->c2/3*gmx::power3(rc - rsw) + p*sc->c3/4*gmx::power4(rc - rsw);
1964 static void potential_switch_constants(real rsw, real rc,
1965 switch_consts_t *sc)
1967 /* The switch function is 1 at rsw and 0 at rc.
1968 * The derivative and second derivate are zero at both ends.
1969 * rsw = max(r - r_switch, 0)
1970 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
1971 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
1972 * force = force*dsw - potential*sw
1975 sc->c3 = -10/gmx::power3(rc - rsw);
1976 sc->c4 = 15/gmx::power4(rc - rsw);
1977 sc->c5 = -6/gmx::power5(rc - rsw);
1980 /*! \brief Construct interaction constants
1982 * This data is used (particularly) by search and force code for
1983 * short-range interactions. Many of these are constant for the whole
1984 * simulation; some are constant only after PME tuning completes.
1987 init_interaction_const(FILE *fp,
1988 interaction_const_t **interaction_const,
1989 const t_forcerec *fr)
1991 interaction_const_t *ic;
1995 ic->cutoff_scheme = fr->cutoff_scheme;
1997 /* Just allocate something so we can free it */
1998 snew_aligned(ic->tabq_coul_FDV0, 16, 32);
1999 snew_aligned(ic->tabq_coul_F, 16, 32);
2000 snew_aligned(ic->tabq_coul_V, 16, 32);
2002 ic->rlist = fr->rlist;
2005 ic->vdwtype = fr->vdwtype;
2006 ic->vdw_modifier = fr->vdw_modifier;
2007 ic->rvdw = fr->rvdw;
2008 ic->rvdw_switch = fr->rvdw_switch;
2009 ic->ewaldcoeff_lj = fr->ewaldcoeff_lj;
2010 ic->ljpme_comb_rule = fr->ljpme_combination_rule;
2011 ic->sh_lj_ewald = 0;
2012 clear_force_switch_constants(&ic->dispersion_shift);
2013 clear_force_switch_constants(&ic->repulsion_shift);
2015 switch (ic->vdw_modifier)
2017 case eintmodPOTSHIFT:
2018 /* Only shift the potential, don't touch the force */
2019 ic->dispersion_shift.cpot = -1.0/gmx::power6(ic->rvdw);
2020 ic->repulsion_shift.cpot = -1.0/gmx::power12(ic->rvdw);
2021 if (EVDW_PME(ic->vdwtype))
2025 crc2 = gmx::square(ic->ewaldcoeff_lj*ic->rvdw);
2026 ic->sh_lj_ewald = (std::exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)/gmx::power6(ic->rvdw);
2029 case eintmodFORCESWITCH:
2030 /* Switch the force, switch and shift the potential */
2031 force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw,
2032 &ic->dispersion_shift);
2033 force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw,
2034 &ic->repulsion_shift);
2036 case eintmodPOTSWITCH:
2037 /* Switch the potential and force */
2038 potential_switch_constants(ic->rvdw_switch, ic->rvdw,
2042 case eintmodEXACTCUTOFF:
2043 /* Nothing to do here */
2046 gmx_incons("unimplemented potential modifier");
2049 ic->sh_invrc6 = -ic->dispersion_shift.cpot;
2051 /* Electrostatics */
2052 ic->eeltype = fr->eeltype;
2053 ic->coulomb_modifier = fr->coulomb_modifier;
2054 ic->rcoulomb = fr->rcoulomb;
2055 ic->epsilon_r = fr->epsilon_r;
2056 ic->epsfac = fr->epsfac;
2057 ic->ewaldcoeff_q = fr->ewaldcoeff_q;
2059 if (EEL_PME_EWALD(ic->eeltype) && ic->coulomb_modifier == eintmodPOTSHIFT)
2061 GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
2062 ic->sh_ewald = std::erfc(ic->ewaldcoeff_q*ic->rcoulomb) / ic->rcoulomb;
2069 /* Reaction-field */
2070 if (EEL_RF(ic->eeltype))
2072 ic->epsilon_rf = fr->epsilon_rf;
2073 ic->k_rf = fr->k_rf;
2074 ic->c_rf = fr->c_rf;
2078 /* For plain cut-off we might use the reaction-field kernels */
2079 ic->epsilon_rf = ic->epsilon_r;
2081 if (fr->coulomb_modifier == eintmodPOTSHIFT)
2083 ic->c_rf = 1/ic->rcoulomb;
2093 real dispersion_shift;
2095 dispersion_shift = ic->dispersion_shift.cpot;
2096 if (EVDW_PME(ic->vdwtype))
2098 dispersion_shift -= ic->sh_lj_ewald;
2100 fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
2101 ic->repulsion_shift.cpot, dispersion_shift);
2103 if (ic->eeltype == eelCUT)
2105 fprintf(fp, ", Coulomb %.e", -ic->c_rf);
2107 else if (EEL_PME(ic->eeltype))
2109 fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
2114 *interaction_const = ic;
2117 static void init_nb_verlet(FILE *fp,
2118 nonbonded_verlet_t **nb_verlet,
2119 gmx_bool bFEP_NonBonded,
2120 const t_inputrec *ir,
2121 const t_forcerec *fr,
2122 const t_commrec *cr,
2123 const char *nbpu_opt)
2125 nonbonded_verlet_t *nbv;
2128 gmx_bool bEmulateGPU, bHybridGPURun = FALSE;
2130 nbnxn_alloc_t *nb_alloc;
2131 nbnxn_free_t *nb_free;
2135 pick_nbnxn_resources(fp, cr, fr->hwinfo,
2142 nbv->min_ci_balanced = 0;
2144 nbv->ngrp = (DOMAINDECOMP(cr) ? 2 : 1);
2145 for (i = 0; i < nbv->ngrp; i++)
2147 nbv->grp[i].nbl_lists.nnbl = 0;
2148 nbv->grp[i].nbat = NULL;
2149 nbv->grp[i].kernel_type = nbnxnkNotSet;
2151 if (i == 0) /* local */
2153 pick_nbnxn_kernel(fp, cr, fr->use_simd_kernels,
2154 nbv->bUseGPU, bEmulateGPU, ir,
2155 &nbv->grp[i].kernel_type,
2156 &nbv->grp[i].ewald_excl,
2159 else /* non-local */
2161 if (nbpu_opt != NULL && strcmp(nbpu_opt, "gpu_cpu") == 0)
2163 /* Use GPU for local, select a CPU kernel for non-local */
2164 pick_nbnxn_kernel(fp, cr, fr->use_simd_kernels,
2166 &nbv->grp[i].kernel_type,
2167 &nbv->grp[i].ewald_excl,
2170 bHybridGPURun = TRUE;
2174 /* Use the same kernel for local and non-local interactions */
2175 nbv->grp[i].kernel_type = nbv->grp[0].kernel_type;
2176 nbv->grp[i].ewald_excl = nbv->grp[0].ewald_excl;
2181 nbnxn_init_search(&nbv->nbs,
2182 DOMAINDECOMP(cr) ? &cr->dd->nc : NULL,
2183 DOMAINDECOMP(cr) ? domdec_zones(cr->dd) : NULL,
2185 gmx_omp_nthreads_get(emntPairsearch));
2187 for (i = 0; i < nbv->ngrp; i++)
2189 gpu_set_host_malloc_and_free(nbv->grp[0].kernel_type == nbnxnk8x8x8_GPU,
2190 &nb_alloc, &nb_free);
2192 nbnxn_init_pairlist_set(&nbv->grp[i].nbl_lists,
2193 nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2194 /* 8x8x8 "non-simple" lists are ATM always combined */
2195 !nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2199 nbv->grp[0].kernel_type != nbv->grp[i].kernel_type)
2201 gmx_bool bSimpleList;
2202 int enbnxninitcombrule;
2204 bSimpleList = nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type);
2206 if (fr->vdwtype == evdwCUT &&
2207 (fr->vdw_modifier == eintmodNONE ||
2208 fr->vdw_modifier == eintmodPOTSHIFT) &&
2209 getenv("GMX_NO_LJ_COMB_RULE") == NULL)
2211 /* Plain LJ cut-off: we can optimize with combination rules */
2212 enbnxninitcombrule = enbnxninitcombruleDETECT;
2214 else if (fr->vdwtype == evdwPME)
2216 /* LJ-PME: we need to use a combination rule for the grid */
2217 if (fr->ljpme_combination_rule == eljpmeGEOM)
2219 enbnxninitcombrule = enbnxninitcombruleGEOM;
2223 enbnxninitcombrule = enbnxninitcombruleLB;
2228 /* We use a full combination matrix: no rule required */
2229 enbnxninitcombrule = enbnxninitcombruleNONE;
2233 snew(nbv->grp[i].nbat, 1);
2234 nbnxn_atomdata_init(fp,
2236 nbv->grp[i].kernel_type,
2238 fr->ntype, fr->nbfp,
2240 bSimpleList ? gmx_omp_nthreads_get(emntNonbonded) : 1,
2245 nbv->grp[i].nbat = nbv->grp[0].nbat;
2251 /* init the NxN GPU data; the last argument tells whether we'll have
2252 * both local and non-local NB calculation on GPU */
2253 nbnxn_gpu_init(&nbv->gpu_nbv,
2254 &fr->hwinfo->gpu_info,
2258 cr->rank_pp_intranode,
2260 (nbv->ngrp > 1) && !bHybridGPURun);
2262 /* With tMPI + GPUs some ranks may be sharing GPU(s) and therefore
2263 * also sharing texture references. To keep the code simple, we don't
2264 * treat texture references as shared resources, but this means that
2265 * the coulomb_tab and nbfp texture refs will get updated by multiple threads.
2266 * Hence, to ensure that the non-bonded kernels don't start before all
2267 * texture binding operations are finished, we need to wait for all ranks
2268 * to arrive here before continuing.
2270 * Note that we could omit this barrier if GPUs are not shared (or
2271 * texture objects are used), but as this is initialization code, there
2272 * is no point in complicating things.
2279 #endif /* GMX_THREAD_MPI */
2281 if ((env = getenv("GMX_NB_MIN_CI")) != NULL)
2285 nbv->min_ci_balanced = strtol(env, &end, 10);
2286 if (!end || (*end != 0) || nbv->min_ci_balanced < 0)
2288 gmx_fatal(FARGS, "Invalid value passed in GMX_NB_MIN_CI=%s, non-negative integer required", env);
2293 fprintf(debug, "Neighbor-list balancing parameter: %d (passed as env. var.)\n",
2294 nbv->min_ci_balanced);
2299 nbv->min_ci_balanced = nbnxn_gpu_min_ci_balanced(nbv->gpu_nbv);
2302 fprintf(debug, "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU multi-processors)\n",
2303 nbv->min_ci_balanced);
2312 gmx_bool usingGpu(nonbonded_verlet_t *nbv)
2314 return nbv != NULL && nbv->bUseGPU;
2317 void init_forcerec(FILE *fp,
2320 const t_inputrec *ir,
2321 const gmx_mtop_t *mtop,
2322 const t_commrec *cr,
2326 const t_filenm *tabbfnm,
2327 const char *nbpu_opt,
2328 gmx_bool bNoSolvOpt,
2331 int i, m, negp_pp, negptable, egi, egj;
2336 gmx_bool bGenericKernelOnly;
2337 gmx_bool needGroupSchemeTables, bSomeNormalNbListsAreInUse;
2338 gmx_bool bFEP_NonBonded;
2339 int *nm_ind, egp_flags;
2341 if (fr->hwinfo == NULL)
2343 /* Detect hardware, gather information.
2344 * In mdrun, hwinfo has already been set before calling init_forcerec.
2345 * Here we ignore GPUs, as tools will not use them anyhow.
2347 fr->hwinfo = gmx_detect_hardware(fp, cr, FALSE);
2350 /* By default we turn SIMD kernels on, but it might be turned off further down... */
2351 fr->use_simd_kernels = TRUE;
2353 fr->bDomDec = DOMAINDECOMP(cr);
2355 if (check_box(ir->ePBC, box))
2357 gmx_fatal(FARGS, check_box(ir->ePBC, box));
2360 /* Test particle insertion ? */
2363 /* Set to the size of the molecule to be inserted (the last one) */
2364 /* Because of old style topologies, we have to use the last cg
2365 * instead of the last molecule type.
2367 cgs = &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].cgs;
2368 fr->n_tpi = cgs->index[cgs->nr] - cgs->index[cgs->nr-1];
2369 if (fr->n_tpi != mtop->mols.index[mtop->mols.nr] - mtop->mols.index[mtop->mols.nr-1])
2371 gmx_fatal(FARGS, "The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
2379 if (ir->coulombtype == eelRF_NEC_UNSUPPORTED)
2381 gmx_fatal(FARGS, "%s electrostatics is no longer supported",
2382 eel_names[ir->coulombtype]);
2387 gmx_fatal(FARGS, "AdResS simulations are no longer supported");
2389 if (ir->useTwinRange)
2391 gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
2393 /* Copy the user determined parameters */
2394 fr->userint1 = ir->userint1;
2395 fr->userint2 = ir->userint2;
2396 fr->userint3 = ir->userint3;
2397 fr->userint4 = ir->userint4;
2398 fr->userreal1 = ir->userreal1;
2399 fr->userreal2 = ir->userreal2;
2400 fr->userreal3 = ir->userreal3;
2401 fr->userreal4 = ir->userreal4;
2404 fr->fc_stepsize = ir->fc_stepsize;
2407 fr->efep = ir->efep;
2408 fr->sc_alphavdw = ir->fepvals->sc_alpha;
2409 if (ir->fepvals->bScCoul)
2411 fr->sc_alphacoul = ir->fepvals->sc_alpha;
2412 fr->sc_sigma6_min = gmx::power6(ir->fepvals->sc_sigma_min);
2416 fr->sc_alphacoul = 0;
2417 fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
2419 fr->sc_power = ir->fepvals->sc_power;
2420 fr->sc_r_power = ir->fepvals->sc_r_power;
2421 fr->sc_sigma6_def = gmx::power6(ir->fepvals->sc_sigma);
2423 env = getenv("GMX_SCSIGMA_MIN");
2427 sscanf(env, "%20lf", &dbl);
2428 fr->sc_sigma6_min = gmx::power6(dbl);
2431 fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
2435 fr->bNonbonded = TRUE;
2436 if (getenv("GMX_NO_NONBONDED") != NULL)
2438 /* turn off non-bonded calculations */
2439 fr->bNonbonded = FALSE;
2440 md_print_warn(cr, fp,
2441 "Found environment variable GMX_NO_NONBONDED.\n"
2442 "Disabling nonbonded calculations.\n");
2445 bGenericKernelOnly = FALSE;
2447 /* We now check in the NS code whether a particular combination of interactions
2448 * can be used with water optimization, and disable it if that is not the case.
2451 if (getenv("GMX_NB_GENERIC") != NULL)
2456 "Found environment variable GMX_NB_GENERIC.\n"
2457 "Disabling all interaction-specific nonbonded kernels, will only\n"
2458 "use the slow generic ones in src/gmxlib/nonbonded/nb_generic.c\n\n");
2460 bGenericKernelOnly = TRUE;
2463 if (bGenericKernelOnly == TRUE)
2468 if ( (getenv("GMX_DISABLE_SIMD_KERNELS") != NULL) || (getenv("GMX_NOOPTIMIZEDKERNELS") != NULL) )
2470 fr->use_simd_kernels = FALSE;
2474 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
2475 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
2476 "(e.g. SSE2/SSE4.1/AVX).\n\n");
2480 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
2482 /* Check if we can/should do all-vs-all kernels */
2483 fr->bAllvsAll = can_use_allvsall(ir, FALSE, NULL, NULL);
2484 fr->AllvsAll_work = NULL;
2485 fr->AllvsAll_workgb = NULL;
2487 /* All-vs-all kernels have not been implemented in 4.6 and later.
2488 * See Redmine #1249. */
2491 fr->bAllvsAll = FALSE;
2495 "\nYour simulation settings would have triggered the efficient all-vs-all\n"
2496 "kernels in GROMACS 4.5, but these have not been implemented in GROMACS\n"
2497 "4.6 and 5.x. If performance is important, please use GROMACS 4.5.7\n"
2498 "or try cutoff-scheme = Verlet.\n\n");
2502 /* Neighbour searching stuff */
2503 fr->cutoff_scheme = ir->cutoff_scheme;
2504 fr->bGrid = (ir->ns_type == ensGRID);
2505 fr->ePBC = ir->ePBC;
2507 if (fr->cutoff_scheme == ecutsGROUP)
2509 const char *note = "NOTE: This file uses the deprecated 'group' cutoff_scheme. This will be\n"
2510 "removed in a future release when 'verlet' supports all interaction forms.\n";
2514 fprintf(stderr, "\n%s\n", note);
2518 fprintf(fp, "\n%s\n", note);
2522 /* Determine if we will do PBC for distances in bonded interactions */
2523 if (fr->ePBC == epbcNONE)
2525 fr->bMolPBC = FALSE;
2529 if (!DOMAINDECOMP(cr))
2533 bSHAKE = (ir->eConstrAlg == econtSHAKE &&
2534 (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0 ||
2535 gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0));
2537 /* The group cut-off scheme and SHAKE assume charge groups
2538 * are whole, but not using molpbc is faster in most cases.
2539 * With intermolecular interactions we need PBC for calculating
2540 * distances between atoms in different molecules.
2542 if ((fr->cutoff_scheme == ecutsGROUP || bSHAKE) &&
2543 !mtop->bIntermolecularInteractions)
2545 fr->bMolPBC = ir->bPeriodicMols;
2547 if (bSHAKE && fr->bMolPBC)
2549 gmx_fatal(FARGS, "SHAKE is not supported with periodic molecules");
2556 if (getenv("GMX_USE_GRAPH") != NULL)
2558 fr->bMolPBC = FALSE;
2561 md_print_warn(cr, fp, "GMX_USE_GRAPH is set, using the graph for bonded interactions\n");
2564 if (mtop->bIntermolecularInteractions)
2566 md_print_warn(cr, fp, "WARNING: Molecules linked by intermolecular interactions have to reside in the same periodic image, otherwise artifacts will occur!\n");
2570 if (bSHAKE && fr->bMolPBC)
2572 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");
2578 fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->ePBC);
2581 fr->bGB = (ir->implicit_solvent == eisGBSA);
2583 fr->rc_scaling = ir->refcoord_scaling;
2584 copy_rvec(ir->posres_com, fr->posres_com);
2585 copy_rvec(ir->posres_comB, fr->posres_comB);
2586 fr->rlist = cutoff_inf(ir->rlist);
2587 fr->eeltype = ir->coulombtype;
2588 fr->vdwtype = ir->vdwtype;
2589 fr->ljpme_combination_rule = ir->ljpme_combination_rule;
2591 fr->coulomb_modifier = ir->coulomb_modifier;
2592 fr->vdw_modifier = ir->vdw_modifier;
2594 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
2595 switch (fr->eeltype)
2598 fr->nbkernel_elec_interaction = (fr->bGB) ? GMX_NBKERNEL_ELEC_GENERALIZEDBORN : GMX_NBKERNEL_ELEC_COULOMB;
2603 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2607 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2608 fr->coulomb_modifier = eintmodEXACTCUTOFF;
2617 case eelPMEUSERSWITCH:
2618 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2624 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD;
2628 gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[fr->eeltype]);
2632 /* Vdw: Translate from mdp settings to kernel format */
2633 switch (fr->vdwtype)
2638 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
2642 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
2646 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD;
2652 case evdwENCADSHIFT:
2653 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2657 gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[fr->vdwtype]);
2661 /* These start out identical to ir, but might be altered if we e.g. tabulate the interaction in the kernel */
2662 fr->nbkernel_elec_modifier = fr->coulomb_modifier;
2663 fr->nbkernel_vdw_modifier = fr->vdw_modifier;
2665 fr->rvdw = cutoff_inf(ir->rvdw);
2666 fr->rvdw_switch = ir->rvdw_switch;
2667 fr->rcoulomb = cutoff_inf(ir->rcoulomb);
2668 fr->rcoulomb_switch = ir->rcoulomb_switch;
2670 fr->bEwald = EEL_PME_EWALD(fr->eeltype);
2672 fr->reppow = mtop->ffparams.reppow;
2674 if (ir->cutoff_scheme == ecutsGROUP)
2676 fr->bvdwtab = ((fr->vdwtype != evdwCUT || !gmx_within_tol(fr->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2677 && !EVDW_PME(fr->vdwtype));
2678 /* We have special kernels for standard Ewald and PME, but the pme-switch ones are tabulated above */
2679 fr->bcoultab = !(fr->eeltype == eelCUT ||
2680 fr->eeltype == eelEWALD ||
2681 fr->eeltype == eelPME ||
2682 fr->eeltype == eelRF ||
2683 fr->eeltype == eelRF_ZERO);
2685 /* If the user absolutely wants different switch/shift settings for coul/vdw, it is likely
2686 * going to be faster to tabulate the interaction than calling the generic kernel.
2687 * However, if generic kernels have been requested we keep things analytically.
2689 if (fr->nbkernel_elec_modifier == eintmodPOTSWITCH &&
2690 fr->nbkernel_vdw_modifier == eintmodPOTSWITCH &&
2691 bGenericKernelOnly == FALSE)
2693 if ((fr->rcoulomb_switch != fr->rvdw_switch) || (fr->rcoulomb != fr->rvdw))
2695 fr->bcoultab = TRUE;
2696 /* Once we tabulate electrostatics, we can use the switch function for LJ,
2697 * which would otherwise need two tables.
2701 else if ((fr->nbkernel_elec_modifier == eintmodPOTSHIFT && fr->nbkernel_vdw_modifier == eintmodPOTSHIFT) ||
2702 ((fr->nbkernel_elec_interaction == GMX_NBKERNEL_ELEC_REACTIONFIELD &&
2703 fr->nbkernel_elec_modifier == eintmodEXACTCUTOFF &&
2704 (fr->nbkernel_vdw_modifier == eintmodPOTSWITCH || fr->nbkernel_vdw_modifier == eintmodPOTSHIFT))))
2706 if ((fr->rcoulomb != fr->rvdw) && (bGenericKernelOnly == FALSE))
2708 fr->bcoultab = TRUE;
2712 if (fr->nbkernel_elec_modifier == eintmodFORCESWITCH)
2714 fr->bcoultab = TRUE;
2716 if (fr->nbkernel_vdw_modifier == eintmodFORCESWITCH)
2721 if (getenv("GMX_REQUIRE_TABLES"))
2724 fr->bcoultab = TRUE;
2729 fprintf(fp, "Table routines are used for coulomb: %s\n",
2730 gmx::boolToString(fr->bcoultab));
2731 fprintf(fp, "Table routines are used for vdw: %s\n",
2732 gmx::boolToString(fr->bvdwtab));
2735 if (fr->bvdwtab == TRUE)
2737 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2738 fr->nbkernel_vdw_modifier = eintmodNONE;
2740 if (fr->bcoultab == TRUE)
2742 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2743 fr->nbkernel_elec_modifier = eintmodNONE;
2747 if (ir->cutoff_scheme == ecutsVERLET)
2749 if (!gmx_within_tol(fr->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2751 gmx_fatal(FARGS, "Cut-off scheme %S only supports LJ repulsion power 12", ecutscheme_names[ir->cutoff_scheme]);
2753 /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
2754 * while mdrun does not (and never did) support this.
2756 if (EEL_USER(fr->eeltype))
2758 gmx_fatal(FARGS, "Combination of %s and cutoff scheme %s is not supported",
2759 eel_names[ir->coulombtype], ecutscheme_names[ir->cutoff_scheme]);
2762 fr->bvdwtab = FALSE;
2763 fr->bcoultab = FALSE;
2766 /* Tables are used for direct ewald sum */
2769 if (EEL_PME(ir->coulombtype))
2773 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
2775 if (ir->coulombtype == eelP3M_AD)
2777 please_cite(fp, "Hockney1988");
2778 please_cite(fp, "Ballenegger2012");
2782 please_cite(fp, "Essmann95a");
2785 if (ir->ewald_geometry == eewg3DC)
2789 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry.\n");
2791 please_cite(fp, "In-Chul99a");
2794 fr->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
2795 init_ewald_tab(&(fr->ewald_table), ir, fp);
2798 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n",
2799 1/fr->ewaldcoeff_q);
2803 if (EVDW_PME(ir->vdwtype))
2807 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
2809 please_cite(fp, "Essmann95a");
2810 fr->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
2813 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n",
2814 1/fr->ewaldcoeff_lj);
2818 /* Electrostatics */
2819 fr->epsilon_r = ir->epsilon_r;
2820 fr->epsilon_rf = ir->epsilon_rf;
2821 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
2823 /* Parameters for generalized RF */
2827 if (fr->eeltype == eelGRF)
2829 init_generalized_rf(fp, mtop, ir, fr);
2832 fr->bF_NoVirSum = (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype) ||
2833 gmx_mtop_ftype_count(mtop, F_POSRES) > 0 ||
2834 gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 ||
2835 inputrecElecField(ir)
2838 if (fr->cutoff_scheme == ecutsGROUP &&
2839 ncg_mtop(mtop) > fr->cg_nalloc && !DOMAINDECOMP(cr))
2841 /* Count the total number of charge groups */
2842 fr->cg_nalloc = ncg_mtop(mtop);
2843 srenew(fr->cg_cm, fr->cg_nalloc);
2845 if (fr->shift_vec == NULL)
2847 snew(fr->shift_vec, SHIFTS);
2850 if (fr->fshift == NULL)
2852 snew(fr->fshift, SHIFTS);
2855 if (fr->nbfp == NULL)
2857 fr->ntype = mtop->ffparams.atnr;
2858 fr->nbfp = mk_nbfp(&mtop->ffparams, fr->bBHAM);
2859 if (EVDW_PME(fr->vdwtype))
2861 fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
2865 /* Copy the energy group exclusions */
2866 fr->egp_flags = ir->opts.egp_flags;
2868 /* Van der Waals stuff */
2869 if ((fr->vdwtype != evdwCUT) && (fr->vdwtype != evdwUSER) && !fr->bBHAM)
2871 if (fr->rvdw_switch >= fr->rvdw)
2873 gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)",
2874 fr->rvdw_switch, fr->rvdw);
2878 fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
2879 (fr->eeltype == eelSWITCH) ? "switched" : "shifted",
2880 fr->rvdw_switch, fr->rvdw);
2884 if (fr->bBHAM && EVDW_PME(fr->vdwtype))
2886 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
2889 if (fr->bBHAM && (fr->vdwtype == evdwSHIFT || fr->vdwtype == evdwSWITCH))
2891 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
2894 if (fr->bBHAM && fr->cutoff_scheme == ecutsVERLET)
2896 gmx_fatal(FARGS, "Verlet cutoff-scheme is not supported with Buckingham");
2901 fprintf(fp, "Cut-off's: NS: %g Coulomb: %g %s: %g\n",
2902 fr->rlist, fr->rcoulomb, fr->bBHAM ? "BHAM" : "LJ", fr->rvdw);
2905 fr->eDispCorr = ir->eDispCorr;
2906 fr->numAtomsForDispersionCorrection = mtop->natoms;
2907 if (ir->eDispCorr != edispcNO)
2909 set_avcsixtwelve(fp, fr, mtop);
2914 set_bham_b_max(fp, fr, mtop);
2917 fr->gb_epsilon_solvent = ir->gb_epsilon_solvent;
2919 /* Copy the GBSA data (radius, volume and surftens for each
2920 * atomtype) from the topology atomtype section to forcerec.
2922 snew(fr->atype_radius, fr->ntype);
2923 snew(fr->atype_vol, fr->ntype);
2924 snew(fr->atype_surftens, fr->ntype);
2925 snew(fr->atype_gb_radius, fr->ntype);
2926 snew(fr->atype_S_hct, fr->ntype);
2928 if (mtop->atomtypes.nr > 0)
2930 for (i = 0; i < fr->ntype; i++)
2932 fr->atype_radius[i] = mtop->atomtypes.radius[i];
2934 for (i = 0; i < fr->ntype; i++)
2936 fr->atype_vol[i] = mtop->atomtypes.vol[i];
2938 for (i = 0; i < fr->ntype; i++)
2940 fr->atype_surftens[i] = mtop->atomtypes.surftens[i];
2942 for (i = 0; i < fr->ntype; i++)
2944 fr->atype_gb_radius[i] = mtop->atomtypes.gb_radius[i];
2946 for (i = 0; i < fr->ntype; i++)
2948 fr->atype_S_hct[i] = mtop->atomtypes.S_hct[i];
2952 /* Generate the GB table if needed */
2956 fr->gbtabscale = 2000;
2958 fr->gbtabscale = 500;
2962 fr->gbtab = make_gb_table(fr);
2964 init_gb(&fr->born, fr, ir, mtop, ir->gb_algorithm);
2966 /* Copy local gb data (for dd, this is done in dd_partition_system) */
2967 if (!DOMAINDECOMP(cr))
2969 make_local_gb(cr, fr->born, ir->gb_algorithm);
2973 /* Set the charge scaling */
2974 if (fr->epsilon_r != 0)
2976 fr->epsfac = ONE_4PI_EPS0/fr->epsilon_r;
2980 /* eps = 0 is infinite dieletric: no coulomb interactions */
2984 /* Reaction field constants */
2985 if (EEL_RF(fr->eeltype))
2987 calc_rffac(fp, fr->eeltype, fr->epsilon_r, fr->epsilon_rf,
2988 fr->rcoulomb, fr->temp, fr->zsquare, box,
2989 &fr->kappa, &fr->k_rf, &fr->c_rf);
2992 /*This now calculates sum for q and c6*/
2993 set_chargesum(fp, fr, mtop);
2995 /* Construct tables for the group scheme. A little unnecessary to
2996 * make both vdw and coul tables sometimes, but what the
2997 * heck. Note that both cutoff schemes construct Ewald tables in
2998 * init_interaction_const_tables. */
2999 needGroupSchemeTables = (ir->cutoff_scheme == ecutsGROUP &&
3000 (fr->bcoultab || fr->bvdwtab));
3002 negp_pp = ir->opts.ngener - ir->nwall;
3004 if (!needGroupSchemeTables)
3006 bSomeNormalNbListsAreInUse = TRUE;
3011 bSomeNormalNbListsAreInUse = FALSE;
3012 for (egi = 0; egi < negp_pp; egi++)
3014 for (egj = egi; egj < negp_pp; egj++)
3016 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
3017 if (!(egp_flags & EGP_EXCL))
3019 if (egp_flags & EGP_TABLE)
3025 bSomeNormalNbListsAreInUse = TRUE;
3030 if (bSomeNormalNbListsAreInUse)
3032 fr->nnblists = negptable + 1;
3036 fr->nnblists = negptable;
3038 if (fr->nnblists > 1)
3040 snew(fr->gid2nblists, ir->opts.ngener*ir->opts.ngener);
3044 snew(fr->nblists, fr->nnblists);
3046 /* This code automatically gives table length tabext without cut-off's,
3047 * in that case grompp should already have checked that we do not need
3048 * normal tables and we only generate tables for 1-4 interactions.
3050 rtab = ir->rlist + ir->tabext;
3052 if (needGroupSchemeTables)
3054 /* make tables for ordinary interactions */
3055 if (bSomeNormalNbListsAreInUse)
3057 make_nbf_tables(fp, fr, rtab, tabfn, NULL, NULL, &fr->nblists[0]);
3066 /* Read the special tables for certain energy group pairs */
3067 nm_ind = mtop->groups.grps[egcENER].nm_ind;
3068 for (egi = 0; egi < negp_pp; egi++)
3070 for (egj = egi; egj < negp_pp; egj++)
3072 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
3073 if ((egp_flags & EGP_TABLE) && !(egp_flags & EGP_EXCL))
3075 if (fr->nnblists > 1)
3077 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = m;
3079 /* Read the table file with the two energy groups names appended */
3080 make_nbf_tables(fp, fr, rtab, tabfn,
3081 *mtop->groups.grpname[nm_ind[egi]],
3082 *mtop->groups.grpname[nm_ind[egj]],
3086 else if (fr->nnblists > 1)
3088 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = 0;
3095 /* Tables might not be used for the potential modifier
3096 * interactions per se, but we still need them to evaluate
3097 * switch/shift dispersion corrections in this case. */
3098 if (fr->eDispCorr != edispcNO)
3100 fr->dispersionCorrectionTable = makeDispersionCorrectionTable(fp, fr, rtab, tabfn);
3103 /* We want to use unmodified tables for 1-4 coulombic
3104 * interactions, so we must in general have an extra set of
3106 if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
3107 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0 ||
3108 gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
3110 fr->pairsTable = make_tables(fp, fr, tabpfn, rtab,
3111 GMX_MAKETABLES_14ONLY);
3115 fr->nwall = ir->nwall;
3116 if (ir->nwall && ir->wall_type == ewtTABLE)
3118 make_wall_tables(fp, ir, tabfn, &mtop->groups, fr);
3123 // Need to catch std::bad_alloc
3124 // TODO Don't need to catch this here, when merging with master branch
3127 fcd->bondtab = make_bonded_tables(fp,
3128 F_TABBONDS, F_TABBONDSNC,
3129 mtop, tabbfnm, "b");
3130 fcd->angletab = make_bonded_tables(fp,
3132 mtop, tabbfnm, "a");
3133 fcd->dihtab = make_bonded_tables(fp,
3135 mtop, tabbfnm, "d");
3137 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
3143 fprintf(debug, "No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
3147 /* QM/MM initialization if requested
3151 fprintf(stderr, "QM/MM calculation requested.\n");
3154 fr->bQMMM = ir->bQMMM;
3155 fr->qr = mk_QMMMrec();
3157 /* Set all the static charge group info */
3158 fr->cginfo_mb = init_cginfo_mb(fp, mtop, fr, bNoSolvOpt,
3160 &fr->bExcl_IntraCGAll_InterCGNone);
3161 if (DOMAINDECOMP(cr))
3167 fr->cginfo = cginfo_expand(mtop->nmolblock, fr->cginfo_mb);
3170 if (!DOMAINDECOMP(cr))
3172 forcerec_set_ranges(fr, ncg_mtop(mtop), ncg_mtop(mtop),
3173 mtop->natoms, mtop->natoms, mtop->natoms);
3176 fr->print_force = print_force;
3179 /* coarse load balancing vars */
3184 /* Initialize neighbor search */
3186 init_ns(fp, cr, fr->ns, fr, mtop);
3188 if (cr->duty & DUTY_PP)
3190 gmx_nonbonded_setup(fr, bGenericKernelOnly);
3193 /* Initialize the thread working data for bonded interactions */
3194 init_bonded_threading(fp, mtop->groups.grps[egcENER].nr,
3195 &fr->bonded_threading);
3197 fr->nthread_ewc = gmx_omp_nthreads_get(emntBonded);
3198 snew(fr->ewc_t, fr->nthread_ewc);
3200 /* fr->ic is used both by verlet and group kernels (to some extent) now */
3201 init_interaction_const(fp, &fr->ic, fr);
3202 init_interaction_const_tables(fp, fr->ic, rtab);
3204 if (fr->cutoff_scheme == ecutsVERLET)
3206 // We checked the cut-offs in grompp, but double-check here.
3207 // We have PME+LJcutoff kernels for rcoulomb>rvdw.
3208 if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
3210 GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw, "With Verlet lists and PME we should have rcoulomb>=rvdw");
3214 GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw, "With Verlet lists and no PME rcoulomb and rvdw should be identical");
3217 init_nb_verlet(fp, &fr->nbv, bFEP_NonBonded, ir, fr, cr, nbpu_opt);
3220 if (ir->eDispCorr != edispcNO)
3222 calc_enervirdiff(fp, ir->eDispCorr, fr);
3226 #define pr_real(fp, r) fprintf(fp, "%s: %e\n",#r, r)
3227 #define pr_int(fp, i) fprintf((fp), "%s: %d\n",#i, i)
3228 #define pr_bool(fp, b) fprintf((fp), "%s: %s\n",#b, gmx::boolToString(b))
3230 void pr_forcerec(FILE *fp, t_forcerec *fr)
3234 pr_real(fp, fr->rlist);
3235 pr_real(fp, fr->rcoulomb);
3236 pr_real(fp, fr->fudgeQQ);
3237 pr_bool(fp, fr->bGrid);
3238 /*pr_int(fp,fr->cg0);
3239 pr_int(fp,fr->hcg);*/
3240 for (i = 0; i < fr->nnblists; i++)
3242 pr_int(fp, fr->nblists[i].table_elec_vdw->n);
3244 pr_real(fp, fr->rcoulomb_switch);
3245 pr_real(fp, fr->rcoulomb);
3250 /* Frees GPU memory and destroys the GPU context.
3252 * Note that this function needs to be called even if GPUs are not used
3253 * in this run because the PME ranks have no knowledge of whether GPUs
3254 * are used or not, but all ranks need to enter the barrier below.
3256 void free_gpu_resources(const t_forcerec *fr,
3257 const t_commrec *cr,
3258 const gmx_gpu_info_t *gpu_info,
3259 const gmx_gpu_opt_t *gpu_opt)
3261 gmx_bool bIsPPrankUsingGPU;
3262 char gpu_err_str[STRLEN];
3264 bIsPPrankUsingGPU = (cr->duty & DUTY_PP) && fr && fr->nbv && fr->nbv->bUseGPU;
3266 if (bIsPPrankUsingGPU)
3268 /* free nbnxn data in GPU memory */
3269 nbnxn_gpu_free(fr->nbv->gpu_nbv);
3270 /* stop the GPU profiler (only CUDA) */
3273 /* With tMPI we need to wait for all ranks to finish deallocation before
3274 * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
3277 * This is not a concern in OpenCL where we use one context per rank which
3278 * is freed in nbnxn_gpu_free().
3280 * Note: as only PP ranks need to free GPU resources, so it is safe to
3281 * not call the barrier on PME ranks.
3288 #endif /* GMX_THREAD_MPI */
3290 /* uninitialize GPU (by destroying the context) */
3291 if (!free_cuda_gpu(cr->rank_pp_intranode, gpu_err_str, gpu_info, gpu_opt))
3293 gmx_warning("On rank %d failed to free GPU #%d: %s",
3294 cr->nodeid, get_current_cuda_gpu_device_id(), gpu_err_str);