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, 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.
48 #include "gromacs/domdec/domdec.h"
49 #include "gromacs/ewald/ewald.h"
50 #include "gromacs/gmxlib/gpu_utils/gpu_utils.h"
51 #include "gromacs/legacyheaders/copyrite.h"
52 #include "gromacs/legacyheaders/force.h"
53 #include "gromacs/legacyheaders/gmx_detect_hardware.h"
54 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
55 #include "gromacs/legacyheaders/inputrec.h"
56 #include "gromacs/legacyheaders/macros.h"
57 #include "gromacs/legacyheaders/md_logging.h"
58 #include "gromacs/legacyheaders/md_support.h"
59 #include "gromacs/legacyheaders/names.h"
60 #include "gromacs/legacyheaders/network.h"
61 #include "gromacs/legacyheaders/nonbonded.h"
62 #include "gromacs/legacyheaders/ns.h"
63 #include "gromacs/legacyheaders/qmmm.h"
64 #include "gromacs/legacyheaders/tables.h"
65 #include "gromacs/legacyheaders/txtdump.h"
66 #include "gromacs/legacyheaders/typedefs.h"
67 #include "gromacs/legacyheaders/types/commrec.h"
68 #include "gromacs/listed-forces/manage-threading.h"
69 #include "gromacs/math/calculate-ewald-splitting-coefficient.h"
70 #include "gromacs/math/units.h"
71 #include "gromacs/math/utilities.h"
72 #include "gromacs/math/vec.h"
73 #include "gromacs/mdlib/forcerec-threading.h"
74 #include "gromacs/mdlib/nb_verlet.h"
75 #include "gromacs/mdlib/nbnxn_atomdata.h"
76 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
77 #include "gromacs/mdlib/nbnxn_search.h"
78 #include "gromacs/mdlib/nbnxn_simd.h"
79 #include "gromacs/pbcutil/ishift.h"
80 #include "gromacs/pbcutil/pbc.h"
81 #include "gromacs/simd/simd.h"
82 #include "gromacs/topology/mtop_util.h"
83 #include "gromacs/utility/cstringutil.h"
84 #include "gromacs/utility/fatalerror.h"
85 #include "gromacs/utility/smalloc.h"
87 #include "nbnxn_gpu_jit_support.h"
89 t_forcerec *mk_forcerec(void)
99 static void pr_nbfp(FILE *fp, real *nbfp, gmx_bool bBHAM, int atnr)
103 for (i = 0; (i < atnr); i++)
105 for (j = 0; (j < atnr); j++)
107 fprintf(fp, "%2d - %2d", i, j);
110 fprintf(fp, " a=%10g, b=%10g, c=%10g\n", BHAMA(nbfp, atnr, i, j),
111 BHAMB(nbfp, atnr, i, j), BHAMC(nbfp, atnr, i, j)/6.0);
115 fprintf(fp, " c6=%10g, c12=%10g\n", C6(nbfp, atnr, i, j)/6.0,
116 C12(nbfp, atnr, i, j)/12.0);
123 static real *mk_nbfp(const gmx_ffparams_t *idef, gmx_bool bBHAM)
131 snew(nbfp, 3*atnr*atnr);
132 for (i = k = 0; (i < atnr); i++)
134 for (j = 0; (j < atnr); j++, k++)
136 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
137 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
138 /* nbfp now includes the 6.0 derivative prefactor */
139 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c*6.0;
145 snew(nbfp, 2*atnr*atnr);
146 for (i = k = 0; (i < atnr); i++)
148 for (j = 0; (j < atnr); j++, k++)
150 /* nbfp now includes the 6.0/12.0 derivative prefactors */
151 C6(nbfp, atnr, i, j) = idef->iparams[k].lj.c6*6.0;
152 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12*12.0;
160 static real *make_ljpme_c6grid(const gmx_ffparams_t *idef, t_forcerec *fr)
163 real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
165 const real oneOverSix = 1.0 / 6.0;
167 /* For LJ-PME simulations, we correct the energies with the reciprocal space
168 * inside of the cut-off. To do this the non-bonded kernels needs to have
169 * access to the C6-values used on the reciprocal grid in pme.c
173 snew(grid, 2*atnr*atnr);
174 for (i = k = 0; (i < atnr); i++)
176 for (j = 0; (j < atnr); j++, k++)
178 c6i = idef->iparams[i*(atnr+1)].lj.c6;
179 c12i = idef->iparams[i*(atnr+1)].lj.c12;
180 c6j = idef->iparams[j*(atnr+1)].lj.c6;
181 c12j = idef->iparams[j*(atnr+1)].lj.c12;
182 c6 = sqrt(c6i * c6j);
183 if (fr->ljpme_combination_rule == eljpmeLB
184 && !gmx_numzero(c6) && !gmx_numzero(c12i) && !gmx_numzero(c12j))
186 sigmai = pow(c12i / c6i, oneOverSix);
187 sigmaj = pow(c12j / c6j, oneOverSix);
188 epsi = c6i * c6i / c12i;
189 epsj = c6j * c6j / c12j;
190 c6 = sqrt(epsi * epsj) * pow(0.5*(sigmai+sigmaj), 6);
192 /* Store the elements at the same relative positions as C6 in nbfp in order
193 * to simplify access in the kernels
195 grid[2*(atnr*i+j)] = c6*6.0;
201 static real *mk_nbfp_combination_rule(const gmx_ffparams_t *idef, int comb_rule)
205 real c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
207 const real oneOverSix = 1.0 / 6.0;
210 snew(nbfp, 2*atnr*atnr);
211 for (i = 0; i < atnr; ++i)
213 for (j = 0; j < atnr; ++j)
215 c6i = idef->iparams[i*(atnr+1)].lj.c6;
216 c12i = idef->iparams[i*(atnr+1)].lj.c12;
217 c6j = idef->iparams[j*(atnr+1)].lj.c6;
218 c12j = idef->iparams[j*(atnr+1)].lj.c12;
219 c6 = sqrt(c6i * c6j);
220 c12 = sqrt(c12i * c12j);
221 if (comb_rule == eCOMB_ARITHMETIC
222 && !gmx_numzero(c6) && !gmx_numzero(c12))
224 sigmai = pow(c12i / c6i, oneOverSix);
225 sigmaj = pow(c12j / c6j, oneOverSix);
226 epsi = c6i * c6i / c12i;
227 epsj = c6j * c6j / c12j;
228 c6 = sqrt(epsi * epsj) * pow(0.5*(sigmai+sigmaj), 6);
229 c12 = sqrt(epsi * epsj) * pow(0.5*(sigmai+sigmaj), 12);
231 C6(nbfp, atnr, i, j) = c6*6.0;
232 C12(nbfp, atnr, i, j) = c12*12.0;
238 /* This routine sets fr->solvent_opt to the most common solvent in the
239 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
240 * the fr->solvent_type array with the correct type (or esolNO).
242 * Charge groups that fulfill the conditions but are not identical to the
243 * most common one will be marked as esolNO in the solvent_type array.
245 * TIP3p is identical to SPC for these purposes, so we call it
246 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
248 * NOTE: QM particle should not
249 * become an optimized solvent. Not even if there is only one charge
259 } solvent_parameters_t;
262 check_solvent_cg(const gmx_moltype_t *molt,
265 const unsigned char *qm_grpnr,
266 const t_grps *qm_grps,
268 int *n_solvent_parameters,
269 solvent_parameters_t **solvent_parameters_p,
279 real tmp_charge[4] = { 0.0 }; /* init to zero to make gcc4.8 happy */
280 int tmp_vdwtype[4] = { 0 }; /* init to zero to make gcc4.8 happy */
283 solvent_parameters_t *solvent_parameters;
285 /* We use a list with parameters for each solvent type.
286 * Every time we discover a new molecule that fulfills the basic
287 * conditions for a solvent we compare with the previous entries
288 * in these lists. If the parameters are the same we just increment
289 * the counter for that type, and otherwise we create a new type
290 * based on the current molecule.
292 * Once we've finished going through all molecules we check which
293 * solvent is most common, and mark all those molecules while we
294 * clear the flag on all others.
297 solvent_parameters = *solvent_parameters_p;
299 /* Mark the cg first as non optimized */
302 /* Check if this cg has no exclusions with atoms in other charge groups
303 * and all atoms inside the charge group excluded.
304 * We only have 3 or 4 atom solvent loops.
306 if (GET_CGINFO_EXCL_INTER(cginfo) ||
307 !GET_CGINFO_EXCL_INTRA(cginfo))
312 /* Get the indices of the first atom in this charge group */
313 j0 = molt->cgs.index[cg0];
314 j1 = molt->cgs.index[cg0+1];
316 /* Number of atoms in our molecule */
322 "Moltype '%s': there are %d atoms in this charge group\n",
326 /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
329 if (nj < 3 || nj > 4)
334 /* Check if we are doing QM on this group */
336 if (qm_grpnr != NULL)
338 for (j = j0; j < j1 && !qm; j++)
340 qm = (qm_grpnr[j] < qm_grps->nr - 1);
343 /* Cannot use solvent optimization with QM */
349 atom = molt->atoms.atom;
351 /* Still looks like a solvent, time to check parameters */
353 /* If it is perturbed (free energy) we can't use the solvent loops,
354 * so then we just skip to the next molecule.
358 for (j = j0; j < j1 && !perturbed; j++)
360 perturbed = PERTURBED(atom[j]);
368 /* Now it's only a question if the VdW and charge parameters
369 * are OK. Before doing the check we compare and see if they are
370 * identical to a possible previous solvent type.
371 * First we assign the current types and charges.
373 for (j = 0; j < nj; j++)
375 tmp_vdwtype[j] = atom[j0+j].type;
376 tmp_charge[j] = atom[j0+j].q;
379 /* Does it match any previous solvent type? */
380 for (k = 0; k < *n_solvent_parameters; k++)
385 /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
386 if ( (solvent_parameters[k].model == esolSPC && nj != 3) ||
387 (solvent_parameters[k].model == esolTIP4P && nj != 4) )
392 /* Check that types & charges match for all atoms in molecule */
393 for (j = 0; j < nj && match == TRUE; j++)
395 if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
399 if (tmp_charge[j] != solvent_parameters[k].charge[j])
406 /* Congratulations! We have a matched solvent.
407 * Flag it with this type for later processing.
410 solvent_parameters[k].count += nmol;
412 /* We are done with this charge group */
417 /* If we get here, we have a tentative new solvent type.
418 * Before we add it we must check that it fulfills the requirements
419 * of the solvent optimized loops. First determine which atoms have
422 for (j = 0; j < nj; j++)
425 tjA = tmp_vdwtype[j];
427 /* Go through all other tpes and see if any have non-zero
428 * VdW parameters when combined with this one.
430 for (k = 0; k < fr->ntype && (has_vdw[j] == FALSE); k++)
432 /* We already checked that the atoms weren't perturbed,
433 * so we only need to check state A now.
437 has_vdw[j] = (has_vdw[j] ||
438 (BHAMA(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
439 (BHAMB(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
440 (BHAMC(fr->nbfp, fr->ntype, tjA, k) != 0.0));
445 has_vdw[j] = (has_vdw[j] ||
446 (C6(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
447 (C12(fr->nbfp, fr->ntype, tjA, k) != 0.0));
452 /* Now we know all we need to make the final check and assignment. */
456 * For this we require thatn all atoms have charge,
457 * the charges on atom 2 & 3 should be the same, and only
458 * atom 1 might have VdW.
460 if (has_vdw[1] == FALSE &&
461 has_vdw[2] == FALSE &&
462 tmp_charge[0] != 0 &&
463 tmp_charge[1] != 0 &&
464 tmp_charge[2] == tmp_charge[1])
466 srenew(solvent_parameters, *n_solvent_parameters+1);
467 solvent_parameters[*n_solvent_parameters].model = esolSPC;
468 solvent_parameters[*n_solvent_parameters].count = nmol;
469 for (k = 0; k < 3; k++)
471 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
472 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
475 *cg_sp = *n_solvent_parameters;
476 (*n_solvent_parameters)++;
481 /* Or could it be a TIP4P?
482 * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
483 * Only atom 1 mght have VdW.
485 if (has_vdw[1] == FALSE &&
486 has_vdw[2] == FALSE &&
487 has_vdw[3] == FALSE &&
488 tmp_charge[0] == 0 &&
489 tmp_charge[1] != 0 &&
490 tmp_charge[2] == tmp_charge[1] &&
493 srenew(solvent_parameters, *n_solvent_parameters+1);
494 solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
495 solvent_parameters[*n_solvent_parameters].count = nmol;
496 for (k = 0; k < 4; k++)
498 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
499 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
502 *cg_sp = *n_solvent_parameters;
503 (*n_solvent_parameters)++;
507 *solvent_parameters_p = solvent_parameters;
511 check_solvent(FILE * fp,
512 const gmx_mtop_t * mtop,
514 cginfo_mb_t *cginfo_mb)
517 const gmx_moltype_t *molt;
518 int mb, mol, cg_mol, at_offset, am, cgm, i, nmol_ch, nmol;
519 int n_solvent_parameters;
520 solvent_parameters_t *solvent_parameters;
526 fprintf(debug, "Going to determine what solvent types we have.\n");
529 n_solvent_parameters = 0;
530 solvent_parameters = NULL;
531 /* Allocate temporary array for solvent type */
532 snew(cg_sp, mtop->nmolblock);
535 for (mb = 0; mb < mtop->nmolblock; mb++)
537 molt = &mtop->moltype[mtop->molblock[mb].type];
539 /* Here we have to loop over all individual molecules
540 * because we need to check for QMMM particles.
542 snew(cg_sp[mb], cginfo_mb[mb].cg_mod);
543 nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
544 nmol = mtop->molblock[mb].nmol/nmol_ch;
545 for (mol = 0; mol < nmol_ch; mol++)
548 am = mol*cgs->index[cgs->nr];
549 for (cg_mol = 0; cg_mol < cgs->nr; cg_mol++)
551 check_solvent_cg(molt, cg_mol, nmol,
552 mtop->groups.grpnr[egcQMMM] ?
553 mtop->groups.grpnr[egcQMMM]+at_offset+am : 0,
554 &mtop->groups.grps[egcQMMM],
556 &n_solvent_parameters, &solvent_parameters,
557 cginfo_mb[mb].cginfo[cgm+cg_mol],
558 &cg_sp[mb][cgm+cg_mol]);
561 at_offset += cgs->index[cgs->nr];
564 /* Puh! We finished going through all charge groups.
565 * Now find the most common solvent model.
568 /* Most common solvent this far */
570 for (i = 0; i < n_solvent_parameters; i++)
573 solvent_parameters[i].count > solvent_parameters[bestsp].count)
581 bestsol = solvent_parameters[bestsp].model;
589 for (mb = 0; mb < mtop->nmolblock; mb++)
591 cgs = &mtop->moltype[mtop->molblock[mb].type].cgs;
592 nmol = (mtop->molblock[mb].nmol*cgs->nr)/cginfo_mb[mb].cg_mod;
593 for (i = 0; i < cginfo_mb[mb].cg_mod; i++)
595 if (cg_sp[mb][i] == bestsp)
597 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], bestsol);
602 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], esolNO);
609 if (bestsol != esolNO && fp != NULL)
611 fprintf(fp, "\nEnabling %s-like water optimization for %d molecules.\n\n",
613 solvent_parameters[bestsp].count);
616 sfree(solvent_parameters);
617 fr->solvent_opt = bestsol;
621 acNONE = 0, acCONSTRAINT, acSETTLE
624 static cginfo_mb_t *init_cginfo_mb(FILE *fplog, const gmx_mtop_t *mtop,
625 t_forcerec *fr, gmx_bool bNoSolvOpt,
626 gmx_bool *bFEP_NonBonded,
627 gmx_bool *bExcl_IntraCGAll_InterCGNone)
630 const t_blocka *excl;
631 const gmx_moltype_t *molt;
632 const gmx_molblock_t *molb;
633 cginfo_mb_t *cginfo_mb;
636 int cg_offset, a_offset;
637 int mb, m, cg, a0, a1, gid, ai, j, aj, excl_nalloc;
641 gmx_bool bId, *bExcl, bExclIntraAll, bExclInter, bHaveVDW, bHaveQ, bHavePerturbedAtoms;
643 snew(cginfo_mb, mtop->nmolblock);
645 snew(type_VDW, fr->ntype);
646 for (ai = 0; ai < fr->ntype; ai++)
648 type_VDW[ai] = FALSE;
649 for (j = 0; j < fr->ntype; j++)
651 type_VDW[ai] = type_VDW[ai] ||
653 C6(fr->nbfp, fr->ntype, ai, j) != 0 ||
654 C12(fr->nbfp, fr->ntype, ai, j) != 0;
658 *bFEP_NonBonded = FALSE;
659 *bExcl_IntraCGAll_InterCGNone = TRUE;
662 snew(bExcl, excl_nalloc);
665 for (mb = 0; mb < mtop->nmolblock; mb++)
667 molb = &mtop->molblock[mb];
668 molt = &mtop->moltype[molb->type];
672 /* Check if the cginfo is identical for all molecules in this block.
673 * If so, we only need an array of the size of one molecule.
674 * Otherwise we make an array of #mol times #cgs per molecule.
677 for (m = 0; m < molb->nmol; m++)
679 int am = m*cgs->index[cgs->nr];
680 for (cg = 0; cg < cgs->nr; cg++)
683 a1 = cgs->index[cg+1];
684 if (ggrpnr(&mtop->groups, egcENER, a_offset+am+a0) !=
685 ggrpnr(&mtop->groups, egcENER, a_offset +a0))
689 if (mtop->groups.grpnr[egcQMMM] != NULL)
691 for (ai = a0; ai < a1; ai++)
693 if (mtop->groups.grpnr[egcQMMM][a_offset+am+ai] !=
694 mtop->groups.grpnr[egcQMMM][a_offset +ai])
703 cginfo_mb[mb].cg_start = cg_offset;
704 cginfo_mb[mb].cg_end = cg_offset + molb->nmol*cgs->nr;
705 cginfo_mb[mb].cg_mod = (bId ? 1 : molb->nmol)*cgs->nr;
706 snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
707 cginfo = cginfo_mb[mb].cginfo;
709 /* Set constraints flags for constrained atoms */
710 snew(a_con, molt->atoms.nr);
711 for (ftype = 0; ftype < F_NRE; ftype++)
713 if (interaction_function[ftype].flags & IF_CONSTRAINT)
718 for (ia = 0; ia < molt->ilist[ftype].nr; ia += 1+nral)
722 for (a = 0; a < nral; a++)
724 a_con[molt->ilist[ftype].iatoms[ia+1+a]] =
725 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
731 for (m = 0; m < (bId ? 1 : molb->nmol); m++)
734 int am = m*cgs->index[cgs->nr];
735 for (cg = 0; cg < cgs->nr; cg++)
738 a1 = cgs->index[cg+1];
740 /* Store the energy group in cginfo */
741 gid = ggrpnr(&mtop->groups, egcENER, a_offset+am+a0);
742 SET_CGINFO_GID(cginfo[cgm+cg], gid);
744 /* Check the intra/inter charge group exclusions */
745 if (a1-a0 > excl_nalloc)
747 excl_nalloc = a1 - a0;
748 srenew(bExcl, excl_nalloc);
750 /* bExclIntraAll: all intra cg interactions excluded
751 * bExclInter: any inter cg interactions excluded
753 bExclIntraAll = TRUE;
757 bHavePerturbedAtoms = FALSE;
758 for (ai = a0; ai < a1; ai++)
760 /* Check VDW and electrostatic interactions */
761 bHaveVDW = bHaveVDW || (type_VDW[molt->atoms.atom[ai].type] ||
762 type_VDW[molt->atoms.atom[ai].typeB]);
763 bHaveQ = bHaveQ || (molt->atoms.atom[ai].q != 0 ||
764 molt->atoms.atom[ai].qB != 0);
766 bHavePerturbedAtoms = bHavePerturbedAtoms || (PERTURBED(molt->atoms.atom[ai]) != 0);
768 /* Clear the exclusion list for atom ai */
769 for (aj = a0; aj < a1; aj++)
771 bExcl[aj-a0] = FALSE;
773 /* Loop over all the exclusions of atom ai */
774 for (j = excl->index[ai]; j < excl->index[ai+1]; j++)
777 if (aj < a0 || aj >= a1)
786 /* Check if ai excludes a0 to a1 */
787 for (aj = a0; aj < a1; aj++)
791 bExclIntraAll = FALSE;
798 SET_CGINFO_CONSTR(cginfo[cgm+cg]);
801 SET_CGINFO_SETTLE(cginfo[cgm+cg]);
809 SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
813 SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
815 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
817 /* The size in cginfo is currently only read with DD */
818 gmx_fatal(FARGS, "A charge group has size %d which is larger than the limit of %d atoms", a1-a0, MAX_CHARGEGROUP_SIZE);
822 SET_CGINFO_HAS_VDW(cginfo[cgm+cg]);
826 SET_CGINFO_HAS_Q(cginfo[cgm+cg]);
828 if (bHavePerturbedAtoms && fr->efep != efepNO)
830 SET_CGINFO_FEP(cginfo[cgm+cg]);
831 *bFEP_NonBonded = TRUE;
833 /* Store the charge group size */
834 SET_CGINFO_NATOMS(cginfo[cgm+cg], a1-a0);
836 if (!bExclIntraAll || bExclInter)
838 *bExcl_IntraCGAll_InterCGNone = FALSE;
845 cg_offset += molb->nmol*cgs->nr;
846 a_offset += molb->nmol*cgs->index[cgs->nr];
850 /* the solvent optimizer is called after the QM is initialized,
851 * because we don't want to have the QM subsystemto become an
855 check_solvent(fplog, mtop, fr, cginfo_mb);
857 if (getenv("GMX_NO_SOLV_OPT"))
861 fprintf(fplog, "Found environment variable GMX_NO_SOLV_OPT.\n"
862 "Disabling all solvent optimization\n");
864 fr->solvent_opt = esolNO;
868 fr->solvent_opt = esolNO;
870 if (!fr->solvent_opt)
872 for (mb = 0; mb < mtop->nmolblock; mb++)
874 for (cg = 0; cg < cginfo_mb[mb].cg_mod; cg++)
876 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg], esolNO);
884 static int *cginfo_expand(int nmb, cginfo_mb_t *cgi_mb)
889 ncg = cgi_mb[nmb-1].cg_end;
892 for (cg = 0; cg < ncg; cg++)
894 while (cg >= cgi_mb[mb].cg_end)
899 cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
905 static void set_chargesum(FILE *log, t_forcerec *fr, const gmx_mtop_t *mtop)
907 /*This now calculates sum for q and c6*/
908 double qsum, q2sum, q, c6sum, c6;
910 const t_atoms *atoms;
915 for (mb = 0; mb < mtop->nmolblock; mb++)
917 nmol = mtop->molblock[mb].nmol;
918 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
919 for (i = 0; i < atoms->nr; i++)
921 q = atoms->atom[i].q;
924 c6 = mtop->ffparams.iparams[atoms->atom[i].type*(mtop->ffparams.atnr+1)].lj.c6;
929 fr->q2sum[0] = q2sum;
930 fr->c6sum[0] = c6sum;
932 if (fr->efep != efepNO)
937 for (mb = 0; mb < mtop->nmolblock; mb++)
939 nmol = mtop->molblock[mb].nmol;
940 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
941 for (i = 0; i < atoms->nr; i++)
943 q = atoms->atom[i].qB;
946 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB*(mtop->ffparams.atnr+1)].lj.c6;
950 fr->q2sum[1] = q2sum;
951 fr->c6sum[1] = c6sum;
956 fr->qsum[1] = fr->qsum[0];
957 fr->q2sum[1] = fr->q2sum[0];
958 fr->c6sum[1] = fr->c6sum[0];
962 if (fr->efep == efepNO)
964 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
968 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n",
969 fr->qsum[0], fr->qsum[1]);
974 void update_forcerec(t_forcerec *fr, matrix box)
976 if (fr->eeltype == eelGRF)
978 calc_rffac(NULL, fr->eeltype, fr->epsilon_r, fr->epsilon_rf,
979 fr->rcoulomb, fr->temp, fr->zsquare, box,
980 &fr->kappa, &fr->k_rf, &fr->c_rf);
984 void set_avcsixtwelve(FILE *fplog, t_forcerec *fr, const gmx_mtop_t *mtop)
986 const t_atoms *atoms, *atoms_tpi;
987 const t_blocka *excl;
988 int mb, nmol, nmolc, i, j, tpi, tpj, j1, j2, k, nexcl, q;
989 gmx_int64_t npair, npair_ij, tmpi, tmpj;
990 double csix, ctwelve;
994 real *nbfp_comb = NULL;
1000 /* For LJ-PME, we want to correct for the difference between the
1001 * actual C6 values and the C6 values used by the LJ-PME based on
1002 * combination rules. */
1004 if (EVDW_PME(fr->vdwtype))
1006 nbfp_comb = mk_nbfp_combination_rule(&mtop->ffparams,
1007 (fr->ljpme_combination_rule == eljpmeLB) ? eCOMB_ARITHMETIC : eCOMB_GEOMETRIC);
1008 for (tpi = 0; tpi < ntp; ++tpi)
1010 for (tpj = 0; tpj < ntp; ++tpj)
1012 C6(nbfp_comb, ntp, tpi, tpj) =
1013 C6(nbfp, ntp, tpi, tpj) - C6(nbfp_comb, ntp, tpi, tpj);
1014 C12(nbfp_comb, ntp, tpi, tpj) = C12(nbfp, ntp, tpi, tpj);
1019 for (q = 0; q < (fr->efep == efepNO ? 1 : 2); q++)
1027 /* Count the types so we avoid natoms^2 operations */
1028 snew(typecount, ntp);
1029 gmx_mtop_count_atomtypes(mtop, q, typecount);
1031 for (tpi = 0; tpi < ntp; tpi++)
1033 for (tpj = tpi; tpj < ntp; tpj++)
1035 tmpi = typecount[tpi];
1036 tmpj = typecount[tpj];
1039 npair_ij = tmpi*tmpj;
1043 npair_ij = tmpi*(tmpi - 1)/2;
1047 /* nbfp now includes the 6.0 derivative prefactor */
1048 csix += npair_ij*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1052 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1053 csix += npair_ij* C6(nbfp, ntp, tpi, tpj)/6.0;
1054 ctwelve += npair_ij* C12(nbfp, ntp, tpi, tpj)/12.0;
1060 /* Subtract the excluded pairs.
1061 * The main reason for substracting exclusions is that in some cases
1062 * some combinations might never occur and the parameters could have
1063 * any value. These unused values should not influence the dispersion
1066 for (mb = 0; mb < mtop->nmolblock; mb++)
1068 nmol = mtop->molblock[mb].nmol;
1069 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
1070 excl = &mtop->moltype[mtop->molblock[mb].type].excls;
1071 for (i = 0; (i < atoms->nr); i++)
1075 tpi = atoms->atom[i].type;
1079 tpi = atoms->atom[i].typeB;
1081 j1 = excl->index[i];
1082 j2 = excl->index[i+1];
1083 for (j = j1; j < j2; j++)
1090 tpj = atoms->atom[k].type;
1094 tpj = atoms->atom[k].typeB;
1098 /* nbfp now includes the 6.0 derivative prefactor */
1099 csix -= nmol*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1103 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1104 csix -= nmol*C6 (nbfp, ntp, tpi, tpj)/6.0;
1105 ctwelve -= nmol*C12(nbfp, ntp, tpi, tpj)/12.0;
1115 /* Only correct for the interaction of the test particle
1116 * with the rest of the system.
1119 &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].atoms;
1122 for (mb = 0; mb < mtop->nmolblock; mb++)
1124 nmol = mtop->molblock[mb].nmol;
1125 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
1126 for (j = 0; j < atoms->nr; j++)
1129 /* Remove the interaction of the test charge group
1132 if (mb == mtop->nmolblock-1)
1136 if (mb == 0 && nmol == 1)
1138 gmx_fatal(FARGS, "Old format tpr with TPI, please generate a new tpr file");
1143 tpj = atoms->atom[j].type;
1147 tpj = atoms->atom[j].typeB;
1149 for (i = 0; i < fr->n_tpi; i++)
1153 tpi = atoms_tpi->atom[i].type;
1157 tpi = atoms_tpi->atom[i].typeB;
1161 /* nbfp now includes the 6.0 derivative prefactor */
1162 csix += nmolc*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1166 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1167 csix += nmolc*C6 (nbfp, ntp, tpi, tpj)/6.0;
1168 ctwelve += nmolc*C12(nbfp, ntp, tpi, tpj)/12.0;
1175 if (npair - nexcl <= 0 && fplog)
1177 fprintf(fplog, "\nWARNING: There are no atom pairs for dispersion correction\n\n");
1183 csix /= npair - nexcl;
1184 ctwelve /= npair - nexcl;
1188 fprintf(debug, "Counted %d exclusions\n", nexcl);
1189 fprintf(debug, "Average C6 parameter is: %10g\n", (double)csix);
1190 fprintf(debug, "Average C12 parameter is: %10g\n", (double)ctwelve);
1192 fr->avcsix[q] = csix;
1193 fr->avctwelve[q] = ctwelve;
1196 if (EVDW_PME(fr->vdwtype))
1203 if (fr->eDispCorr == edispcAllEner ||
1204 fr->eDispCorr == edispcAllEnerPres)
1206 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1207 fr->avcsix[0], fr->avctwelve[0]);
1211 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e\n", fr->avcsix[0]);
1217 static void set_bham_b_max(FILE *fplog, t_forcerec *fr,
1218 const gmx_mtop_t *mtop)
1220 const t_atoms *at1, *at2;
1221 int mt1, mt2, i, j, tpi, tpj, ntypes;
1227 fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
1234 for (mt1 = 0; mt1 < mtop->nmoltype; mt1++)
1236 at1 = &mtop->moltype[mt1].atoms;
1237 for (i = 0; (i < at1->nr); i++)
1239 tpi = at1->atom[i].type;
1242 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
1245 for (mt2 = mt1; mt2 < mtop->nmoltype; mt2++)
1247 at2 = &mtop->moltype[mt2].atoms;
1248 for (j = 0; (j < at2->nr); j++)
1250 tpj = at2->atom[j].type;
1253 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
1255 b = BHAMB(nbfp, ntypes, tpi, tpj);
1256 if (b > fr->bham_b_max)
1260 if ((b < bmin) || (bmin == -1))
1270 fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n",
1271 bmin, fr->bham_b_max);
1275 static void make_nbf_tables(FILE *fp, const output_env_t oenv,
1276 t_forcerec *fr, real rtab,
1277 const t_commrec *cr,
1278 const char *tabfn, char *eg1, char *eg2,
1288 fprintf(debug, "No table file name passed, can not read table, can not do non-bonded interactions\n");
1293 sprintf(buf, "%s", tabfn);
1296 /* Append the two energy group names */
1297 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
1298 eg1, eg2, ftp2ext(efXVG));
1300 nbl->table_elec_vdw = make_tables(fp, oenv, fr, MASTER(cr), buf, rtab, 0);
1301 /* Copy the contents of the table to separate coulomb and LJ tables too,
1302 * to improve cache performance.
1304 /* For performance reasons we want
1305 * the table data to be aligned to 16-byte. The pointers could be freed
1306 * but currently aren't.
1308 nbl->table_elec.interaction = GMX_TABLE_INTERACTION_ELEC;
1309 nbl->table_elec.format = nbl->table_elec_vdw.format;
1310 nbl->table_elec.r = nbl->table_elec_vdw.r;
1311 nbl->table_elec.n = nbl->table_elec_vdw.n;
1312 nbl->table_elec.scale = nbl->table_elec_vdw.scale;
1313 nbl->table_elec.scale_exp = nbl->table_elec_vdw.scale_exp;
1314 nbl->table_elec.formatsize = nbl->table_elec_vdw.formatsize;
1315 nbl->table_elec.ninteractions = 1;
1316 nbl->table_elec.stride = nbl->table_elec.formatsize * nbl->table_elec.ninteractions;
1317 snew_aligned(nbl->table_elec.data, nbl->table_elec.stride*(nbl->table_elec.n+1), 32);
1319 nbl->table_vdw.interaction = GMX_TABLE_INTERACTION_VDWREP_VDWDISP;
1320 nbl->table_vdw.format = nbl->table_elec_vdw.format;
1321 nbl->table_vdw.r = nbl->table_elec_vdw.r;
1322 nbl->table_vdw.n = nbl->table_elec_vdw.n;
1323 nbl->table_vdw.scale = nbl->table_elec_vdw.scale;
1324 nbl->table_vdw.scale_exp = nbl->table_elec_vdw.scale_exp;
1325 nbl->table_vdw.formatsize = nbl->table_elec_vdw.formatsize;
1326 nbl->table_vdw.ninteractions = 2;
1327 nbl->table_vdw.stride = nbl->table_vdw.formatsize * nbl->table_vdw.ninteractions;
1328 snew_aligned(nbl->table_vdw.data, nbl->table_vdw.stride*(nbl->table_vdw.n+1), 32);
1330 for (i = 0; i <= nbl->table_elec_vdw.n; i++)
1332 for (j = 0; j < 4; j++)
1334 nbl->table_elec.data[4*i+j] = nbl->table_elec_vdw.data[12*i+j];
1336 for (j = 0; j < 8; j++)
1338 nbl->table_vdw.data[8*i+j] = nbl->table_elec_vdw.data[12*i+4+j];
1343 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t *mtop,
1344 int *ncount, int **count)
1346 const gmx_moltype_t *molt;
1348 int mt, ftype, stride, i, j, tabnr;
1350 for (mt = 0; mt < mtop->nmoltype; mt++)
1352 molt = &mtop->moltype[mt];
1353 for (ftype = 0; ftype < F_NRE; ftype++)
1355 if (ftype == ftype1 || ftype == ftype2)
1357 il = &molt->ilist[ftype];
1358 stride = 1 + NRAL(ftype);
1359 for (i = 0; i < il->nr; i += stride)
1361 tabnr = mtop->ffparams.iparams[il->iatoms[i]].tab.table;
1364 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
1366 if (tabnr >= *ncount)
1368 srenew(*count, tabnr+1);
1369 for (j = *ncount; j < tabnr+1; j++)
1382 static bondedtable_t *make_bonded_tables(FILE *fplog,
1383 int ftype1, int ftype2,
1384 const gmx_mtop_t *mtop,
1385 const char *basefn, const char *tabext)
1387 int i, ncount, *count;
1395 count_tables(ftype1, ftype2, mtop, &ncount, &count);
1400 for (i = 0; i < ncount; i++)
1404 sprintf(tabfn, "%s", basefn);
1405 sprintf(tabfn + strlen(basefn) - strlen(ftp2ext(efXVG)) - 1, "_%s%d.%s",
1406 tabext, i, ftp2ext(efXVG));
1407 tab[i] = make_bonded_table(fplog, tabfn, NRAL(ftype1)-2);
1416 void forcerec_set_ranges(t_forcerec *fr,
1417 int ncg_home, int ncg_force,
1419 int natoms_force_constr, int natoms_f_novirsum)
1424 /* fr->ncg_force is unused in the standard code,
1425 * but it can be useful for modified code dealing with charge groups.
1427 fr->ncg_force = ncg_force;
1428 fr->natoms_force = natoms_force;
1429 fr->natoms_force_constr = natoms_force_constr;
1431 if (fr->natoms_force_constr > fr->nalloc_force)
1433 fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1437 srenew(fr->f_twin, fr->nalloc_force);
1441 if (fr->bF_NoVirSum)
1443 fr->f_novirsum_n = natoms_f_novirsum;
1444 if (fr->f_novirsum_n > fr->f_novirsum_nalloc)
1446 fr->f_novirsum_nalloc = over_alloc_dd(fr->f_novirsum_n);
1447 srenew(fr->f_novirsum_alloc, fr->f_novirsum_nalloc);
1452 fr->f_novirsum_n = 0;
1456 static real cutoff_inf(real cutoff)
1460 cutoff = GMX_CUTOFF_INF;
1466 static void make_adress_tf_tables(FILE *fp, const output_env_t oenv,
1467 t_forcerec *fr, const t_inputrec *ir,
1468 const char *tabfn, const gmx_mtop_t *mtop,
1476 gmx_fatal(FARGS, "No thermoforce table file given. Use -tabletf to specify a file\n");
1480 snew(fr->atf_tabs, ir->adress->n_tf_grps);
1482 sprintf(buf, "%s", tabfn);
1483 for (i = 0; i < ir->adress->n_tf_grps; i++)
1485 j = ir->adress->tf_table_index[i]; /* get energy group index */
1486 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "tf_%s.%s",
1487 *(mtop->groups.grpname[mtop->groups.grps[egcENER].nm_ind[j]]), ftp2ext(efXVG));
1490 fprintf(fp, "loading tf table for energygrp index %d from %s\n", ir->adress->tf_table_index[i], buf);
1492 fr->atf_tabs[i] = make_atf_table(fp, oenv, fr, buf, box);
1497 gmx_bool can_use_allvsall(const t_inputrec *ir, gmx_bool bPrintNote, t_commrec *cr, FILE *fp)
1504 ir->rcoulomb == 0 &&
1506 ir->ePBC == epbcNONE &&
1507 ir->vdwtype == evdwCUT &&
1508 ir->coulombtype == eelCUT &&
1509 ir->efep == efepNO &&
1510 (ir->implicit_solvent == eisNO ||
1511 (ir->implicit_solvent == eisGBSA && (ir->gb_algorithm == egbSTILL ||
1512 ir->gb_algorithm == egbHCT ||
1513 ir->gb_algorithm == egbOBC))) &&
1514 getenv("GMX_NO_ALLVSALL") == NULL
1517 if (bAllvsAll && ir->opts.ngener > 1)
1519 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";
1525 fprintf(stderr, "\n%s\n", note);
1529 fprintf(fp, "\n%s\n", note);
1535 if (bAllvsAll && fp && MASTER(cr))
1537 fprintf(fp, "\nUsing SIMD all-vs-all kernels.\n\n");
1544 gmx_bool nbnxn_acceleration_supported(FILE *fplog,
1545 const t_commrec *cr,
1546 const t_inputrec *ir,
1549 if (!bGPU && (ir->vdwtype == evdwPME && ir->ljpme_combination_rule == eljpmeLB))
1551 md_print_warn(cr, fplog, "LJ-PME with Lorentz-Berthelot is not supported with %s, falling back to %s\n",
1552 bGPU ? "GPUs" : "SIMD kernels",
1553 bGPU ? "CPU only" : "plain-C kernels");
1561 static void pick_nbnxn_kernel_cpu(const t_inputrec gmx_unused *ir,
1565 *kernel_type = nbnxnk4x4_PlainC;
1566 *ewald_excl = ewaldexclTable;
1568 #ifdef GMX_NBNXN_SIMD
1570 #ifdef GMX_NBNXN_SIMD_4XN
1571 *kernel_type = nbnxnk4xN_SIMD_4xN;
1573 #ifdef GMX_NBNXN_SIMD_2XNN
1574 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1577 #if defined GMX_NBNXN_SIMD_2XNN && defined GMX_NBNXN_SIMD_4XN
1578 /* We need to choose if we want 2x(N+N) or 4xN kernels.
1579 * Currently this is based on the SIMD acceleration choice,
1580 * but it might be better to decide this at runtime based on CPU.
1582 * 4xN calculates more (zero) interactions, but has less pair-search
1583 * work and much better kernel instruction scheduling.
1585 * Up till now we have only seen that on Intel Sandy/Ivy Bridge,
1586 * which doesn't have FMA, both the analytical and tabulated Ewald
1587 * kernels have similar pair rates for 4x8 and 2x(4+4), so we choose
1588 * 2x(4+4) because it results in significantly fewer pairs.
1589 * For RF, the raw pair rate of the 4x8 kernel is higher than 2x(4+4),
1590 * 10% with HT, 50% without HT. As we currently don't detect the actual
1591 * use of HT, use 4x8 to avoid a potential performance hit.
1592 * On Intel Haswell 4x8 is always faster.
1594 *kernel_type = nbnxnk4xN_SIMD_4xN;
1596 #ifndef GMX_SIMD_HAVE_FMA
1597 if (EEL_PME_EWALD(ir->coulombtype) ||
1598 EVDW_PME(ir->vdwtype))
1600 /* We have Ewald kernels without FMA (Intel Sandy/Ivy Bridge).
1601 * There are enough instructions to make 2x(4+4) efficient.
1603 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1606 #endif /* GMX_NBNXN_SIMD_2XNN && GMX_NBNXN_SIMD_4XN */
1609 if (getenv("GMX_NBNXN_SIMD_4XN") != NULL)
1611 #ifdef GMX_NBNXN_SIMD_4XN
1612 *kernel_type = nbnxnk4xN_SIMD_4xN;
1614 gmx_fatal(FARGS, "SIMD 4xN kernels requested, but GROMACS has been compiled without support for these kernels");
1617 if (getenv("GMX_NBNXN_SIMD_2XNN") != NULL)
1619 #ifdef GMX_NBNXN_SIMD_2XNN
1620 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1622 gmx_fatal(FARGS, "SIMD 2x(N+N) kernels requested, but GROMACS has been compiled without support for these kernels");
1626 /* Analytical Ewald exclusion correction is only an option in
1628 * Since table lookup's don't parallelize with SIMD, analytical
1629 * will probably always be faster for a SIMD width of 8 or more.
1630 * With FMA analytical is sometimes faster for a width if 4 as well.
1631 * On BlueGene/Q, this is faster regardless of precision.
1632 * In single precision, this is faster on Bulldozer.
1634 #if GMX_SIMD_REAL_WIDTH >= 8 || \
1635 (GMX_SIMD_REAL_WIDTH >= 4 && defined GMX_SIMD_HAVE_FMA && !defined GMX_DOUBLE) || \
1636 defined GMX_SIMD_IBM_QPX
1637 *ewald_excl = ewaldexclAnalytical;
1639 if (getenv("GMX_NBNXN_EWALD_TABLE") != NULL)
1641 *ewald_excl = ewaldexclTable;
1643 if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != NULL)
1645 *ewald_excl = ewaldexclAnalytical;
1649 #endif /* GMX_NBNXN_SIMD */
1653 const char *lookup_nbnxn_kernel_name(int kernel_type)
1655 const char *returnvalue = NULL;
1656 switch (kernel_type)
1659 returnvalue = "not set";
1661 case nbnxnk4x4_PlainC:
1662 returnvalue = "plain C";
1664 case nbnxnk4xN_SIMD_4xN:
1665 case nbnxnk4xN_SIMD_2xNN:
1666 #ifdef GMX_NBNXN_SIMD
1667 #if defined GMX_SIMD_X86_SSE2
1668 returnvalue = "SSE2";
1669 #elif defined GMX_SIMD_X86_SSE4_1
1670 returnvalue = "SSE4.1";
1671 #elif defined GMX_SIMD_X86_AVX_128_FMA
1672 returnvalue = "AVX_128_FMA";
1673 #elif defined GMX_SIMD_X86_AVX_256
1674 returnvalue = "AVX_256";
1675 #elif defined GMX_SIMD_X86_AVX2_256
1676 returnvalue = "AVX2_256";
1678 returnvalue = "SIMD";
1680 #else /* GMX_NBNXN_SIMD */
1681 returnvalue = "not available";
1682 #endif /* GMX_NBNXN_SIMD */
1684 case nbnxnk8x8x8_GPU: returnvalue = "GPU"; break;
1685 case nbnxnk8x8x8_PlainC: returnvalue = "plain C"; break;
1689 gmx_fatal(FARGS, "Illegal kernel type selected");
1696 static void pick_nbnxn_kernel(FILE *fp,
1697 const t_commrec *cr,
1698 gmx_bool use_simd_kernels,
1700 gmx_bool bEmulateGPU,
1701 const t_inputrec *ir,
1704 gmx_bool bDoNonbonded)
1706 assert(kernel_type);
1708 *kernel_type = nbnxnkNotSet;
1709 *ewald_excl = ewaldexclTable;
1713 *kernel_type = nbnxnk8x8x8_PlainC;
1717 md_print_warn(cr, fp, "Emulating a GPU run on the CPU (slow)");
1722 *kernel_type = nbnxnk8x8x8_GPU;
1725 if (*kernel_type == nbnxnkNotSet)
1727 /* LJ PME with LB combination rule does 7 mesh operations.
1728 * This so slow that we don't compile SIMD non-bonded kernels for that.
1730 if (use_simd_kernels &&
1731 nbnxn_acceleration_supported(fp, cr, ir, FALSE))
1733 pick_nbnxn_kernel_cpu(ir, kernel_type, ewald_excl);
1737 *kernel_type = nbnxnk4x4_PlainC;
1741 if (bDoNonbonded && fp != NULL)
1743 fprintf(fp, "\nUsing %s %dx%d non-bonded kernels\n\n",
1744 lookup_nbnxn_kernel_name(*kernel_type),
1745 nbnxn_kernel_to_ci_size(*kernel_type),
1746 nbnxn_kernel_to_cj_size(*kernel_type));
1748 if (nbnxnk4x4_PlainC == *kernel_type ||
1749 nbnxnk8x8x8_PlainC == *kernel_type)
1751 md_print_warn(cr, fp,
1752 "WARNING: Using the slow %s kernels. This should\n"
1753 "not happen during routine usage on supported platforms.\n\n",
1754 lookup_nbnxn_kernel_name(*kernel_type));
1759 static void pick_nbnxn_resources(FILE *fp,
1760 const t_commrec *cr,
1761 const gmx_hw_info_t *hwinfo,
1762 gmx_bool bDoNonbonded,
1764 gmx_bool *bEmulateGPU,
1765 const gmx_gpu_opt_t *gpu_opt)
1767 gmx_bool bEmulateGPUEnvVarSet;
1768 char gpu_err_str[STRLEN];
1772 bEmulateGPUEnvVarSet = (getenv("GMX_EMULATE_GPU") != NULL);
1774 /* Run GPU emulation mode if GMX_EMULATE_GPU is defined. Because
1775 * GPUs (currently) only handle non-bonded calculations, we will
1776 * automatically switch to emulation if non-bonded calculations are
1777 * turned off via GMX_NO_NONBONDED - this is the simple and elegant
1778 * way to turn off GPU initialization, data movement, and cleanup.
1780 * GPU emulation can be useful to assess the performance one can expect by
1781 * adding GPU(s) to the machine. The conditional below allows this even
1782 * if mdrun is compiled without GPU acceleration support.
1783 * Note that you should freezing the system as otherwise it will explode.
1785 *bEmulateGPU = (bEmulateGPUEnvVarSet ||
1786 (!bDoNonbonded && gpu_opt->n_dev_use > 0));
1788 /* Enable GPU mode when GPUs are available or no GPU emulation is requested.
1790 if (gpu_opt->n_dev_use > 0 && !(*bEmulateGPU))
1792 /* Each PP node will use the intra-node id-th device from the
1793 * list of detected/selected GPUs. */
1794 if (!init_gpu(fp, cr->rank_pp_intranode, gpu_err_str,
1795 &hwinfo->gpu_info, gpu_opt))
1797 /* At this point the init should never fail as we made sure that
1798 * we have all the GPUs we need. If it still does, we'll bail. */
1799 /* TODO the decorating of gpu_err_str is nicer if it
1800 happens inside init_gpu. Out here, the decorating with
1801 the MPI rank makes sense. */
1802 gmx_fatal(FARGS, "On rank %d failed to initialize GPU #%d: %s",
1804 get_gpu_device_id(&hwinfo->gpu_info, gpu_opt,
1805 cr->rank_pp_intranode),
1809 /* Here we actually turn on hardware GPU acceleration */
1814 gmx_bool uses_simple_tables(int cutoff_scheme,
1815 nonbonded_verlet_t *nbv,
1818 gmx_bool bUsesSimpleTables = TRUE;
1821 switch (cutoff_scheme)
1824 bUsesSimpleTables = TRUE;
1827 assert(NULL != nbv && NULL != nbv->grp);
1828 grp_index = (group < 0) ? 0 : (nbv->ngrp - 1);
1829 bUsesSimpleTables = nbnxn_kernel_pairlist_simple(nbv->grp[grp_index].kernel_type);
1832 gmx_incons("unimplemented");
1834 return bUsesSimpleTables;
1837 static void init_ewald_f_table(interaction_const_t *ic,
1842 /* Get the Ewald table spacing based on Coulomb and/or LJ
1843 * Ewald coefficients and rtol.
1845 ic->tabq_scale = ewald_spline3_table_scale(ic);
1847 if (ic->cutoff_scheme == ecutsVERLET)
1849 maxr = ic->rcoulomb;
1853 maxr = std::max(ic->rcoulomb, rtab);
1855 ic->tabq_size = static_cast<int>(maxr*ic->tabq_scale) + 2;
1857 sfree_aligned(ic->tabq_coul_FDV0);
1858 sfree_aligned(ic->tabq_coul_F);
1859 sfree_aligned(ic->tabq_coul_V);
1861 sfree_aligned(ic->tabq_vdw_FDV0);
1862 sfree_aligned(ic->tabq_vdw_F);
1863 sfree_aligned(ic->tabq_vdw_V);
1865 if (ic->eeltype == eelEWALD || EEL_PME(ic->eeltype))
1867 /* Create the original table data in FDV0 */
1868 snew_aligned(ic->tabq_coul_FDV0, ic->tabq_size*4, 32);
1869 snew_aligned(ic->tabq_coul_F, ic->tabq_size, 32);
1870 snew_aligned(ic->tabq_coul_V, ic->tabq_size, 32);
1871 table_spline3_fill_ewald_lr(ic->tabq_coul_F, ic->tabq_coul_V, ic->tabq_coul_FDV0,
1872 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_q, v_q_ewald_lr);
1875 if (EVDW_PME(ic->vdwtype))
1877 snew_aligned(ic->tabq_vdw_FDV0, ic->tabq_size*4, 32);
1878 snew_aligned(ic->tabq_vdw_F, ic->tabq_size, 32);
1879 snew_aligned(ic->tabq_vdw_V, ic->tabq_size, 32);
1880 table_spline3_fill_ewald_lr(ic->tabq_vdw_F, ic->tabq_vdw_V, ic->tabq_vdw_FDV0,
1881 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_lj, v_lj_ewald_lr);
1885 void init_interaction_const_tables(FILE *fp,
1886 interaction_const_t *ic,
1889 if (ic->eeltype == eelEWALD || EEL_PME(ic->eeltype) || EVDW_PME(ic->vdwtype))
1891 init_ewald_f_table(ic, rtab);
1895 fprintf(fp, "Initialized non-bonded Ewald correction tables, spacing: %.2e size: %d\n\n",
1896 1/ic->tabq_scale, ic->tabq_size);
1901 static void clear_force_switch_constants(shift_consts_t *sc)
1908 static void force_switch_constants(real p,
1912 /* Here we determine the coefficient for shifting the force to zero
1913 * between distance rsw and the cut-off rc.
1914 * For a potential of r^-p, we have force p*r^-(p+1).
1915 * But to save flops we absorb p in the coefficient.
1917 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
1918 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
1920 sc->c2 = ((p + 1)*rsw - (p + 4)*rc)/(pow(rc, p + 2)*pow(rc - rsw, 2));
1921 sc->c3 = -((p + 1)*rsw - (p + 3)*rc)/(pow(rc, p + 2)*pow(rc - rsw, 3));
1922 sc->cpot = -pow(rc, -p) + p*sc->c2/3*pow(rc - rsw, 3) + p*sc->c3/4*pow(rc - rsw, 4);
1925 static void potential_switch_constants(real rsw, real rc,
1926 switch_consts_t *sc)
1928 /* The switch function is 1 at rsw and 0 at rc.
1929 * The derivative and second derivate are zero at both ends.
1930 * rsw = max(r - r_switch, 0)
1931 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
1932 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
1933 * force = force*dsw - potential*sw
1936 sc->c3 = -10*pow(rc - rsw, -3);
1937 sc->c4 = 15*pow(rc - rsw, -4);
1938 sc->c5 = -6*pow(rc - rsw, -5);
1941 /*! \brief Construct interaction constants
1943 * This data is used (particularly) by search and force code for
1944 * short-range interactions. Many of these are constant for the whole
1945 * simulation; some are constant only after PME tuning completes.
1948 init_interaction_const(FILE *fp,
1949 interaction_const_t **interaction_const,
1950 const t_forcerec *fr)
1952 interaction_const_t *ic;
1953 const real minusSix = -6.0;
1954 const real minusTwelve = -12.0;
1958 ic->cutoff_scheme = fr->cutoff_scheme;
1960 /* Just allocate something so we can free it */
1961 snew_aligned(ic->tabq_coul_FDV0, 16, 32);
1962 snew_aligned(ic->tabq_coul_F, 16, 32);
1963 snew_aligned(ic->tabq_coul_V, 16, 32);
1965 ic->rlist = fr->rlist;
1966 ic->rlistlong = fr->rlistlong;
1969 ic->vdwtype = fr->vdwtype;
1970 ic->vdw_modifier = fr->vdw_modifier;
1971 ic->rvdw = fr->rvdw;
1972 ic->rvdw_switch = fr->rvdw_switch;
1973 ic->ewaldcoeff_lj = fr->ewaldcoeff_lj;
1974 ic->ljpme_comb_rule = fr->ljpme_combination_rule;
1975 ic->sh_lj_ewald = 0;
1976 clear_force_switch_constants(&ic->dispersion_shift);
1977 clear_force_switch_constants(&ic->repulsion_shift);
1979 switch (ic->vdw_modifier)
1981 case eintmodPOTSHIFT:
1982 /* Only shift the potential, don't touch the force */
1983 ic->dispersion_shift.cpot = -pow(ic->rvdw, minusSix);
1984 ic->repulsion_shift.cpot = -pow(ic->rvdw, minusTwelve);
1985 if (EVDW_PME(ic->vdwtype))
1989 crc2 = sqr(ic->ewaldcoeff_lj*ic->rvdw);
1990 ic->sh_lj_ewald = (exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)*pow(ic->rvdw, minusSix);
1993 case eintmodFORCESWITCH:
1994 /* Switch the force, switch and shift the potential */
1995 force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw,
1996 &ic->dispersion_shift);
1997 force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw,
1998 &ic->repulsion_shift);
2000 case eintmodPOTSWITCH:
2001 /* Switch the potential and force */
2002 potential_switch_constants(ic->rvdw_switch, ic->rvdw,
2006 case eintmodEXACTCUTOFF:
2007 /* Nothing to do here */
2010 gmx_incons("unimplemented potential modifier");
2013 ic->sh_invrc6 = -ic->dispersion_shift.cpot;
2015 /* Electrostatics */
2016 ic->eeltype = fr->eeltype;
2017 ic->coulomb_modifier = fr->coulomb_modifier;
2018 ic->rcoulomb = fr->rcoulomb;
2019 ic->epsilon_r = fr->epsilon_r;
2020 ic->epsfac = fr->epsfac;
2021 ic->ewaldcoeff_q = fr->ewaldcoeff_q;
2023 if (fr->coulomb_modifier == eintmodPOTSHIFT)
2025 ic->sh_ewald = gmx_erfc(ic->ewaldcoeff_q*ic->rcoulomb);
2032 /* Reaction-field */
2033 if (EEL_RF(ic->eeltype))
2035 ic->epsilon_rf = fr->epsilon_rf;
2036 ic->k_rf = fr->k_rf;
2037 ic->c_rf = fr->c_rf;
2041 /* For plain cut-off we might use the reaction-field kernels */
2042 ic->epsilon_rf = ic->epsilon_r;
2044 if (fr->coulomb_modifier == eintmodPOTSHIFT)
2046 ic->c_rf = 1/ic->rcoulomb;
2056 real dispersion_shift;
2058 dispersion_shift = ic->dispersion_shift.cpot;
2059 if (EVDW_PME(ic->vdwtype))
2061 dispersion_shift -= ic->sh_lj_ewald;
2063 fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
2064 ic->repulsion_shift.cpot, dispersion_shift);
2066 if (ic->eeltype == eelCUT)
2068 fprintf(fp, ", Coulomb %.e", -ic->c_rf);
2070 else if (EEL_PME(ic->eeltype))
2072 fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
2077 *interaction_const = ic;
2080 static void init_nb_verlet(FILE *fp,
2081 nonbonded_verlet_t **nb_verlet,
2082 gmx_bool bFEP_NonBonded,
2083 const t_inputrec *ir,
2084 const t_forcerec *fr,
2085 const t_commrec *cr,
2086 const char *nbpu_opt)
2088 nonbonded_verlet_t *nbv;
2091 gmx_bool bEmulateGPU, bHybridGPURun = FALSE;
2093 nbnxn_alloc_t *nb_alloc;
2094 nbnxn_free_t *nb_free;
2098 pick_nbnxn_resources(fp, cr, fr->hwinfo,
2105 nbv->min_ci_balanced = 0;
2107 nbv->ngrp = (DOMAINDECOMP(cr) ? 2 : 1);
2108 for (i = 0; i < nbv->ngrp; i++)
2110 nbv->grp[i].nbl_lists.nnbl = 0;
2111 nbv->grp[i].nbat = NULL;
2112 nbv->grp[i].kernel_type = nbnxnkNotSet;
2114 if (i == 0) /* local */
2116 pick_nbnxn_kernel(fp, cr, fr->use_simd_kernels,
2117 nbv->bUseGPU, bEmulateGPU, ir,
2118 &nbv->grp[i].kernel_type,
2119 &nbv->grp[i].ewald_excl,
2122 else /* non-local */
2124 if (nbpu_opt != NULL && strcmp(nbpu_opt, "gpu_cpu") == 0)
2126 /* Use GPU for local, select a CPU kernel for non-local */
2127 pick_nbnxn_kernel(fp, cr, fr->use_simd_kernels,
2129 &nbv->grp[i].kernel_type,
2130 &nbv->grp[i].ewald_excl,
2133 bHybridGPURun = TRUE;
2137 /* Use the same kernel for local and non-local interactions */
2138 nbv->grp[i].kernel_type = nbv->grp[0].kernel_type;
2139 nbv->grp[i].ewald_excl = nbv->grp[0].ewald_excl;
2144 nbnxn_init_search(&nbv->nbs,
2145 DOMAINDECOMP(cr) ? &cr->dd->nc : NULL,
2146 DOMAINDECOMP(cr) ? domdec_zones(cr->dd) : NULL,
2148 gmx_omp_nthreads_get(emntPairsearch));
2150 for (i = 0; i < nbv->ngrp; i++)
2152 gpu_set_host_malloc_and_free(nbv->grp[0].kernel_type == nbnxnk8x8x8_GPU,
2153 &nb_alloc, &nb_free);
2155 nbnxn_init_pairlist_set(&nbv->grp[i].nbl_lists,
2156 nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2157 /* 8x8x8 "non-simple" lists are ATM always combined */
2158 !nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2162 nbv->grp[0].kernel_type != nbv->grp[i].kernel_type)
2164 gmx_bool bSimpleList;
2165 int enbnxninitcombrule;
2167 bSimpleList = nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type);
2169 if (bSimpleList && (fr->vdwtype == evdwCUT && (fr->vdw_modifier == eintmodNONE || fr->vdw_modifier == eintmodPOTSHIFT)))
2171 /* Plain LJ cut-off: we can optimize with combination rules */
2172 enbnxninitcombrule = enbnxninitcombruleDETECT;
2174 else if (fr->vdwtype == evdwPME)
2176 /* LJ-PME: we need to use a combination rule for the grid */
2177 if (fr->ljpme_combination_rule == eljpmeGEOM)
2179 enbnxninitcombrule = enbnxninitcombruleGEOM;
2183 enbnxninitcombrule = enbnxninitcombruleLB;
2188 /* We use a full combination matrix: no rule required */
2189 enbnxninitcombrule = enbnxninitcombruleNONE;
2193 snew(nbv->grp[i].nbat, 1);
2194 nbnxn_atomdata_init(fp,
2196 nbv->grp[i].kernel_type,
2198 fr->ntype, fr->nbfp,
2200 bSimpleList ? gmx_omp_nthreads_get(emntNonbonded) : 1,
2205 nbv->grp[i].nbat = nbv->grp[0].nbat;
2211 /* init the NxN GPU data; the last argument tells whether we'll have
2212 * both local and non-local NB calculation on GPU */
2213 nbnxn_gpu_init(fp, &nbv->gpu_nbv,
2214 &fr->hwinfo->gpu_info,
2218 cr->rank_pp_intranode,
2220 (nbv->ngrp > 1) && !bHybridGPURun);
2222 /* With tMPI + GPUs some ranks may be sharing GPU(s) and therefore
2223 * also sharing texture references. To keep the code simple, we don't
2224 * treat texture references as shared resources, but this means that
2225 * the coulomb_tab and nbfp texture refs will get updated by multiple threads.
2226 * Hence, to ensure that the non-bonded kernels don't start before all
2227 * texture binding operations are finished, we need to wait for all ranks
2228 * to arrive here before continuing.
2230 * Note that we could omit this barrier if GPUs are not shared (or
2231 * texture objects are used), but as this is initialization code, there
2232 * is no point in complicating things.
2234 #ifdef GMX_THREAD_MPI
2239 #endif /* GMX_THREAD_MPI */
2241 if ((env = getenv("GMX_NB_MIN_CI")) != NULL)
2245 nbv->min_ci_balanced = strtol(env, &end, 10);
2246 if (!end || (*end != 0) || nbv->min_ci_balanced <= 0)
2248 gmx_fatal(FARGS, "Invalid value passed in GMX_NB_MIN_CI=%s, positive integer required", env);
2253 fprintf(debug, "Neighbor-list balancing parameter: %d (passed as env. var.)\n",
2254 nbv->min_ci_balanced);
2259 nbv->min_ci_balanced = nbnxn_gpu_min_ci_balanced(nbv->gpu_nbv);
2262 fprintf(debug, "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU multi-processors)\n",
2263 nbv->min_ci_balanced);
2272 gmx_bool usingGpu(nonbonded_verlet_t *nbv)
2274 return nbv != NULL && nbv->bUseGPU;
2277 void init_forcerec(FILE *fp,
2278 const output_env_t oenv,
2281 const t_inputrec *ir,
2282 const gmx_mtop_t *mtop,
2283 const t_commrec *cr,
2289 const char *nbpu_opt,
2290 gmx_bool bNoSolvOpt,
2293 int i, m, negp_pp, negptable, egi, egj;
2298 gmx_bool bGenericKernelOnly;
2299 gmx_bool bMakeTables, bMakeSeparate14Table, bSomeNormalNbListsAreInUse;
2300 gmx_bool bFEP_NonBonded;
2301 int *nm_ind, egp_flags;
2303 if (fr->hwinfo == NULL)
2305 /* Detect hardware, gather information.
2306 * In mdrun, hwinfo has already been set before calling init_forcerec.
2307 * Here we ignore GPUs, as tools will not use them anyhow.
2309 fr->hwinfo = gmx_detect_hardware(fp, cr, FALSE);
2312 /* By default we turn SIMD kernels on, but it might be turned off further down... */
2313 fr->use_simd_kernels = TRUE;
2315 fr->bDomDec = DOMAINDECOMP(cr);
2317 if (check_box(ir->ePBC, box))
2319 gmx_fatal(FARGS, check_box(ir->ePBC, box));
2322 /* Test particle insertion ? */
2325 /* Set to the size of the molecule to be inserted (the last one) */
2326 /* Because of old style topologies, we have to use the last cg
2327 * instead of the last molecule type.
2329 cgs = &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].cgs;
2330 fr->n_tpi = cgs->index[cgs->nr] - cgs->index[cgs->nr-1];
2331 if (fr->n_tpi != mtop->mols.index[mtop->mols.nr] - mtop->mols.index[mtop->mols.nr-1])
2333 gmx_fatal(FARGS, "The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
2341 /* Copy AdResS parameters */
2344 fr->adress_type = ir->adress->type;
2345 fr->adress_const_wf = ir->adress->const_wf;
2346 fr->adress_ex_width = ir->adress->ex_width;
2347 fr->adress_hy_width = ir->adress->hy_width;
2348 fr->adress_icor = ir->adress->icor;
2349 fr->adress_site = ir->adress->site;
2350 fr->adress_ex_forcecap = ir->adress->ex_forcecap;
2351 fr->adress_do_hybridpairs = ir->adress->do_hybridpairs;
2354 snew(fr->adress_group_explicit, ir->adress->n_energy_grps);
2355 for (i = 0; i < ir->adress->n_energy_grps; i++)
2357 fr->adress_group_explicit[i] = ir->adress->group_explicit[i];
2360 fr->n_adress_tf_grps = ir->adress->n_tf_grps;
2361 snew(fr->adress_tf_table_index, fr->n_adress_tf_grps);
2362 for (i = 0; i < fr->n_adress_tf_grps; i++)
2364 fr->adress_tf_table_index[i] = ir->adress->tf_table_index[i];
2366 copy_rvec(ir->adress->refs, fr->adress_refs);
2370 fr->adress_type = eAdressOff;
2371 fr->adress_do_hybridpairs = FALSE;
2374 /* Copy the user determined parameters */
2375 fr->userint1 = ir->userint1;
2376 fr->userint2 = ir->userint2;
2377 fr->userint3 = ir->userint3;
2378 fr->userint4 = ir->userint4;
2379 fr->userreal1 = ir->userreal1;
2380 fr->userreal2 = ir->userreal2;
2381 fr->userreal3 = ir->userreal3;
2382 fr->userreal4 = ir->userreal4;
2385 fr->fc_stepsize = ir->fc_stepsize;
2388 fr->efep = ir->efep;
2389 fr->sc_alphavdw = ir->fepvals->sc_alpha;
2390 if (ir->fepvals->bScCoul)
2392 fr->sc_alphacoul = ir->fepvals->sc_alpha;
2393 fr->sc_sigma6_min = pow(ir->fepvals->sc_sigma_min, 6);
2397 fr->sc_alphacoul = 0;
2398 fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
2400 fr->sc_power = ir->fepvals->sc_power;
2401 fr->sc_r_power = ir->fepvals->sc_r_power;
2402 fr->sc_sigma6_def = pow(ir->fepvals->sc_sigma, 6);
2404 env = getenv("GMX_SCSIGMA_MIN");
2408 sscanf(env, "%20lf", &dbl);
2409 fr->sc_sigma6_min = pow(dbl, 6);
2412 fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
2416 fr->bNonbonded = TRUE;
2417 if (getenv("GMX_NO_NONBONDED") != NULL)
2419 /* turn off non-bonded calculations */
2420 fr->bNonbonded = FALSE;
2421 md_print_warn(cr, fp,
2422 "Found environment variable GMX_NO_NONBONDED.\n"
2423 "Disabling nonbonded calculations.\n");
2426 bGenericKernelOnly = FALSE;
2428 /* We now check in the NS code whether a particular combination of interactions
2429 * can be used with water optimization, and disable it if that is not the case.
2432 if (getenv("GMX_NB_GENERIC") != NULL)
2437 "Found environment variable GMX_NB_GENERIC.\n"
2438 "Disabling all interaction-specific nonbonded kernels, will only\n"
2439 "use the slow generic ones in src/gmxlib/nonbonded/nb_generic.c\n\n");
2441 bGenericKernelOnly = TRUE;
2444 if (bGenericKernelOnly == TRUE)
2449 if ( (getenv("GMX_DISABLE_SIMD_KERNELS") != NULL) || (getenv("GMX_NOOPTIMIZEDKERNELS") != NULL) )
2451 fr->use_simd_kernels = FALSE;
2455 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
2456 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
2457 "(e.g. SSE2/SSE4.1/AVX).\n\n");
2461 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
2463 /* Check if we can/should do all-vs-all kernels */
2464 fr->bAllvsAll = can_use_allvsall(ir, FALSE, NULL, NULL);
2465 fr->AllvsAll_work = NULL;
2466 fr->AllvsAll_workgb = NULL;
2468 /* All-vs-all kernels have not been implemented in 4.6 and later.
2469 * See Redmine #1249. */
2472 fr->bAllvsAll = FALSE;
2476 "\nYour simulation settings would have triggered the efficient all-vs-all\n"
2477 "kernels in GROMACS 4.5, but these have not been implemented in GROMACS\n"
2478 "4.6 and 5.x. If performance is important, please use GROMACS 4.5.7\n"
2479 "or try cutoff-scheme = Verlet.\n\n");
2483 /* Neighbour searching stuff */
2484 fr->cutoff_scheme = ir->cutoff_scheme;
2485 fr->bGrid = (ir->ns_type == ensGRID);
2486 fr->ePBC = ir->ePBC;
2488 if (fr->cutoff_scheme == ecutsGROUP)
2490 const char *note = "NOTE: This file uses the deprecated 'group' cutoff_scheme. This will be\n"
2491 "removed in a future release when 'verlet' supports all interaction forms.\n";
2495 fprintf(stderr, "\n%s\n", note);
2499 fprintf(fp, "\n%s\n", note);
2503 /* Determine if we will do PBC for distances in bonded interactions */
2504 if (fr->ePBC == epbcNONE)
2506 fr->bMolPBC = FALSE;
2510 if (!DOMAINDECOMP(cr))
2514 bSHAKE = (ir->eConstrAlg == econtSHAKE &&
2515 (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0 ||
2516 gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0));
2518 /* The group cut-off scheme and SHAKE assume charge groups
2519 * are whole, but not using molpbc is faster in most cases.
2520 * With intermolecular interactions we need PBC for calculating
2521 * distances between atoms in different molecules.
2523 if ((fr->cutoff_scheme == ecutsGROUP || bSHAKE) &&
2524 !mtop->bIntermolecularInteractions)
2526 fr->bMolPBC = ir->bPeriodicMols;
2528 if (bSHAKE && fr->bMolPBC)
2530 gmx_fatal(FARGS, "SHAKE is not supported with periodic molecules");
2537 if (getenv("GMX_USE_GRAPH") != NULL)
2539 fr->bMolPBC = FALSE;
2542 md_print_warn(cr, fp, "GMX_USE_GRAPH is set, using the graph for bonded interactions\n");
2545 if (mtop->bIntermolecularInteractions)
2547 md_print_warn(cr, fp, "WARNING: Molecules linked by intermolecular interactions have to reside in the same periodic image, otherwise artifacts will occur!\n");
2551 if (bSHAKE && fr->bMolPBC)
2553 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");
2559 fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->ePBC);
2562 fr->bGB = (ir->implicit_solvent == eisGBSA);
2564 fr->rc_scaling = ir->refcoord_scaling;
2565 copy_rvec(ir->posres_com, fr->posres_com);
2566 copy_rvec(ir->posres_comB, fr->posres_comB);
2567 fr->rlist = cutoff_inf(ir->rlist);
2568 fr->rlistlong = cutoff_inf(ir->rlistlong);
2569 fr->eeltype = ir->coulombtype;
2570 fr->vdwtype = ir->vdwtype;
2571 fr->ljpme_combination_rule = ir->ljpme_combination_rule;
2573 fr->coulomb_modifier = ir->coulomb_modifier;
2574 fr->vdw_modifier = ir->vdw_modifier;
2576 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
2577 switch (fr->eeltype)
2580 fr->nbkernel_elec_interaction = (fr->bGB) ? GMX_NBKERNEL_ELEC_GENERALIZEDBORN : GMX_NBKERNEL_ELEC_COULOMB;
2586 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2590 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2591 fr->coulomb_modifier = eintmodEXACTCUTOFF;
2600 case eelPMEUSERSWITCH:
2601 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2606 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD;
2610 gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[fr->eeltype]);
2614 /* Vdw: Translate from mdp settings to kernel format */
2615 switch (fr->vdwtype)
2620 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
2624 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
2628 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD;
2634 case evdwENCADSHIFT:
2635 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2639 gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[fr->vdwtype]);
2643 /* These start out identical to ir, but might be altered if we e.g. tabulate the interaction in the kernel */
2644 fr->nbkernel_elec_modifier = fr->coulomb_modifier;
2645 fr->nbkernel_vdw_modifier = fr->vdw_modifier;
2647 fr->rvdw = cutoff_inf(ir->rvdw);
2648 fr->rvdw_switch = ir->rvdw_switch;
2649 fr->rcoulomb = cutoff_inf(ir->rcoulomb);
2650 fr->rcoulomb_switch = ir->rcoulomb_switch;
2652 fr->bTwinRange = fr->rlistlong > fr->rlist;
2653 fr->bEwald = (EEL_PME(fr->eeltype) || fr->eeltype == eelEWALD);
2655 fr->reppow = mtop->ffparams.reppow;
2657 if (ir->cutoff_scheme == ecutsGROUP)
2659 fr->bvdwtab = ((fr->vdwtype != evdwCUT || !gmx_within_tol(fr->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2660 && !EVDW_PME(fr->vdwtype));
2661 /* We have special kernels for standard Ewald and PME, but the pme-switch ones are tabulated above */
2662 fr->bcoultab = !(fr->eeltype == eelCUT ||
2663 fr->eeltype == eelEWALD ||
2664 fr->eeltype == eelPME ||
2665 fr->eeltype == eelRF ||
2666 fr->eeltype == eelRF_ZERO);
2668 /* If the user absolutely wants different switch/shift settings for coul/vdw, it is likely
2669 * going to be faster to tabulate the interaction than calling the generic kernel.
2670 * However, if generic kernels have been requested we keep things analytically.
2672 if (fr->nbkernel_elec_modifier == eintmodPOTSWITCH &&
2673 fr->nbkernel_vdw_modifier == eintmodPOTSWITCH &&
2674 bGenericKernelOnly == FALSE)
2676 if ((fr->rcoulomb_switch != fr->rvdw_switch) || (fr->rcoulomb != fr->rvdw))
2678 fr->bcoultab = TRUE;
2679 /* Once we tabulate electrostatics, we can use the switch function for LJ,
2680 * which would otherwise need two tables.
2684 else if ((fr->nbkernel_elec_modifier == eintmodPOTSHIFT && fr->nbkernel_vdw_modifier == eintmodPOTSHIFT) ||
2685 ((fr->nbkernel_elec_interaction == GMX_NBKERNEL_ELEC_REACTIONFIELD &&
2686 fr->nbkernel_elec_modifier == eintmodEXACTCUTOFF &&
2687 (fr->nbkernel_vdw_modifier == eintmodPOTSWITCH || fr->nbkernel_vdw_modifier == eintmodPOTSHIFT))))
2689 if ((fr->rcoulomb != fr->rvdw) && (bGenericKernelOnly == FALSE))
2691 fr->bcoultab = TRUE;
2695 if (fr->nbkernel_elec_modifier == eintmodFORCESWITCH)
2697 fr->bcoultab = TRUE;
2699 if (fr->nbkernel_vdw_modifier == eintmodFORCESWITCH)
2704 if (getenv("GMX_REQUIRE_TABLES"))
2707 fr->bcoultab = TRUE;
2712 fprintf(fp, "Table routines are used for coulomb: %s\n", bool_names[fr->bcoultab]);
2713 fprintf(fp, "Table routines are used for vdw: %s\n", bool_names[fr->bvdwtab ]);
2716 if (fr->bvdwtab == TRUE)
2718 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2719 fr->nbkernel_vdw_modifier = eintmodNONE;
2721 if (fr->bcoultab == TRUE)
2723 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2724 fr->nbkernel_elec_modifier = eintmodNONE;
2728 if (ir->cutoff_scheme == ecutsVERLET)
2730 if (!gmx_within_tol(fr->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2732 gmx_fatal(FARGS, "Cut-off scheme %S only supports LJ repulsion power 12", ecutscheme_names[ir->cutoff_scheme]);
2734 fr->bvdwtab = FALSE;
2735 fr->bcoultab = FALSE;
2738 /* Tables are used for direct ewald sum */
2741 if (EEL_PME(ir->coulombtype))
2745 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
2747 if (ir->coulombtype == eelP3M_AD)
2749 please_cite(fp, "Hockney1988");
2750 please_cite(fp, "Ballenegger2012");
2754 please_cite(fp, "Essmann95a");
2757 if (ir->ewald_geometry == eewg3DC)
2761 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry.\n");
2763 please_cite(fp, "In-Chul99a");
2766 fr->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
2767 init_ewald_tab(&(fr->ewald_table), ir, fp);
2770 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n",
2771 1/fr->ewaldcoeff_q);
2775 if (EVDW_PME(ir->vdwtype))
2779 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
2781 please_cite(fp, "Essmann95a");
2782 fr->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
2785 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n",
2786 1/fr->ewaldcoeff_lj);
2790 /* Electrostatics */
2791 fr->epsilon_r = ir->epsilon_r;
2792 fr->epsilon_rf = ir->epsilon_rf;
2793 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
2795 /* Parameters for generalized RF */
2799 if (fr->eeltype == eelGRF)
2801 init_generalized_rf(fp, mtop, ir, fr);
2804 fr->bF_NoVirSum = (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype) ||
2805 gmx_mtop_ftype_count(mtop, F_POSRES) > 0 ||
2806 gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 ||
2807 IR_ELEC_FIELD(*ir) ||
2808 (fr->adress_icor != eAdressICOff)
2811 if (fr->cutoff_scheme == ecutsGROUP &&
2812 ncg_mtop(mtop) > fr->cg_nalloc && !DOMAINDECOMP(cr))
2814 /* Count the total number of charge groups */
2815 fr->cg_nalloc = ncg_mtop(mtop);
2816 srenew(fr->cg_cm, fr->cg_nalloc);
2818 if (fr->shift_vec == NULL)
2820 snew(fr->shift_vec, SHIFTS);
2823 if (fr->fshift == NULL)
2825 snew(fr->fshift, SHIFTS);
2828 if (fr->nbfp == NULL)
2830 fr->ntype = mtop->ffparams.atnr;
2831 fr->nbfp = mk_nbfp(&mtop->ffparams, fr->bBHAM);
2832 if (EVDW_PME(fr->vdwtype))
2834 fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
2838 /* Copy the energy group exclusions */
2839 fr->egp_flags = ir->opts.egp_flags;
2841 /* Van der Waals stuff */
2842 if ((fr->vdwtype != evdwCUT) && (fr->vdwtype != evdwUSER) && !fr->bBHAM)
2844 if (fr->rvdw_switch >= fr->rvdw)
2846 gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)",
2847 fr->rvdw_switch, fr->rvdw);
2851 fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
2852 (fr->eeltype == eelSWITCH) ? "switched" : "shifted",
2853 fr->rvdw_switch, fr->rvdw);
2857 if (fr->bBHAM && EVDW_PME(fr->vdwtype))
2859 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
2862 if (fr->bBHAM && (fr->vdwtype == evdwSHIFT || fr->vdwtype == evdwSWITCH))
2864 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
2867 if (fr->bBHAM && fr->cutoff_scheme == ecutsVERLET)
2869 gmx_fatal(FARGS, "Verlet cutoff-scheme is not supported with Buckingham");
2874 fprintf(fp, "Cut-off's: NS: %g Coulomb: %g %s: %g\n",
2875 fr->rlist, fr->rcoulomb, fr->bBHAM ? "BHAM" : "LJ", fr->rvdw);
2878 fr->eDispCorr = ir->eDispCorr;
2879 if (ir->eDispCorr != edispcNO)
2881 set_avcsixtwelve(fp, fr, mtop);
2886 set_bham_b_max(fp, fr, mtop);
2889 fr->gb_epsilon_solvent = ir->gb_epsilon_solvent;
2891 /* Copy the GBSA data (radius, volume and surftens for each
2892 * atomtype) from the topology atomtype section to forcerec.
2894 snew(fr->atype_radius, fr->ntype);
2895 snew(fr->atype_vol, fr->ntype);
2896 snew(fr->atype_surftens, fr->ntype);
2897 snew(fr->atype_gb_radius, fr->ntype);
2898 snew(fr->atype_S_hct, fr->ntype);
2900 if (mtop->atomtypes.nr > 0)
2902 for (i = 0; i < fr->ntype; i++)
2904 fr->atype_radius[i] = mtop->atomtypes.radius[i];
2906 for (i = 0; i < fr->ntype; i++)
2908 fr->atype_vol[i] = mtop->atomtypes.vol[i];
2910 for (i = 0; i < fr->ntype; i++)
2912 fr->atype_surftens[i] = mtop->atomtypes.surftens[i];
2914 for (i = 0; i < fr->ntype; i++)
2916 fr->atype_gb_radius[i] = mtop->atomtypes.gb_radius[i];
2918 for (i = 0; i < fr->ntype; i++)
2920 fr->atype_S_hct[i] = mtop->atomtypes.S_hct[i];
2924 /* Generate the GB table if needed */
2928 fr->gbtabscale = 2000;
2930 fr->gbtabscale = 500;
2934 fr->gbtab = make_gb_table(oenv, fr);
2936 init_gb(&fr->born, fr, ir, mtop, ir->gb_algorithm);
2938 /* Copy local gb data (for dd, this is done in dd_partition_system) */
2939 if (!DOMAINDECOMP(cr))
2941 make_local_gb(cr, fr->born, ir->gb_algorithm);
2945 /* Set the charge scaling */
2946 if (fr->epsilon_r != 0)
2948 fr->epsfac = ONE_4PI_EPS0/fr->epsilon_r;
2952 /* eps = 0 is infinite dieletric: no coulomb interactions */
2956 /* Reaction field constants */
2957 if (EEL_RF(fr->eeltype))
2959 calc_rffac(fp, fr->eeltype, fr->epsilon_r, fr->epsilon_rf,
2960 fr->rcoulomb, fr->temp, fr->zsquare, box,
2961 &fr->kappa, &fr->k_rf, &fr->c_rf);
2964 /*This now calculates sum for q and c6*/
2965 set_chargesum(fp, fr, mtop);
2967 /* if we are using LR electrostatics, and they are tabulated,
2968 * the tables will contain modified coulomb interactions.
2969 * Since we want to use the non-shifted ones for 1-4
2970 * coulombic interactions, we must have an extra set of tables.
2973 /* Construct tables.
2974 * A little unnecessary to make both vdw and coul tables sometimes,
2975 * but what the heck... */
2977 bMakeTables = fr->bcoultab || fr->bvdwtab || fr->bEwald ||
2978 (ir->eDispCorr != edispcNO && ir_vdw_switched(ir));
2980 bMakeSeparate14Table = ((!bMakeTables || fr->eeltype != eelCUT || fr->vdwtype != evdwCUT ||
2981 fr->coulomb_modifier != eintmodNONE ||
2982 fr->vdw_modifier != eintmodNONE ||
2983 fr->bBHAM || fr->bEwald) &&
2984 (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
2985 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0 ||
2986 gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0));
2988 negp_pp = ir->opts.ngener - ir->nwall;
2992 bSomeNormalNbListsAreInUse = TRUE;
2997 bSomeNormalNbListsAreInUse = (ir->eDispCorr != edispcNO);
2998 for (egi = 0; egi < negp_pp; egi++)
3000 for (egj = egi; egj < negp_pp; egj++)
3002 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
3003 if (!(egp_flags & EGP_EXCL))
3005 if (egp_flags & EGP_TABLE)
3011 bSomeNormalNbListsAreInUse = TRUE;
3016 if (bSomeNormalNbListsAreInUse)
3018 fr->nnblists = negptable + 1;
3022 fr->nnblists = negptable;
3024 if (fr->nnblists > 1)
3026 snew(fr->gid2nblists, ir->opts.ngener*ir->opts.ngener);
3035 snew(fr->nblists, fr->nnblists);
3037 /* This code automatically gives table length tabext without cut-off's,
3038 * in that case grompp should already have checked that we do not need
3039 * normal tables and we only generate tables for 1-4 interactions.
3041 rtab = ir->rlistlong + ir->tabext;
3045 /* make tables for ordinary interactions */
3046 if (bSomeNormalNbListsAreInUse)
3048 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn, NULL, NULL, &fr->nblists[0]);
3051 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn, NULL, NULL, &fr->nblists[fr->nnblists/2]);
3053 if (!bMakeSeparate14Table)
3055 fr->tab14 = fr->nblists[0].table_elec_vdw;
3065 /* Read the special tables for certain energy group pairs */
3066 nm_ind = mtop->groups.grps[egcENER].nm_ind;
3067 for (egi = 0; egi < negp_pp; egi++)
3069 for (egj = egi; egj < negp_pp; egj++)
3071 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
3072 if ((egp_flags & EGP_TABLE) && !(egp_flags & EGP_EXCL))
3074 if (fr->nnblists > 1)
3076 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = m;
3078 /* Read the table file with the two energy groups names appended */
3079 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn,
3080 *mtop->groups.grpname[nm_ind[egi]],
3081 *mtop->groups.grpname[nm_ind[egj]],
3085 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn,
3086 *mtop->groups.grpname[nm_ind[egi]],
3087 *mtop->groups.grpname[nm_ind[egj]],
3088 &fr->nblists[fr->nnblists/2+m]);
3092 else if (fr->nnblists > 1)
3094 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = 0;
3100 else if ((fr->eDispCorr != edispcNO) &&
3101 ((fr->vdw_modifier == eintmodPOTSWITCH) ||
3102 (fr->vdw_modifier == eintmodFORCESWITCH) ||
3103 (fr->vdw_modifier == eintmodPOTSHIFT)))
3105 /* Tables might not be used for the potential modifier interactions per se, but
3106 * we still need them to evaluate switch/shift dispersion corrections in this case.
3108 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn, NULL, NULL, &fr->nblists[0]);
3111 if (bMakeSeparate14Table)
3113 /* generate extra tables with plain Coulomb for 1-4 interactions only */
3114 fr->tab14 = make_tables(fp, oenv, fr, MASTER(cr), tabpfn, rtab,
3115 GMX_MAKETABLES_14ONLY);
3118 /* Read AdResS Thermo Force table if needed */
3119 if (fr->adress_icor == eAdressICThermoForce)
3121 /* old todo replace */
3123 if (ir->adress->n_tf_grps > 0)
3125 make_adress_tf_tables(fp, oenv, fr, ir, tabfn, mtop, box);
3130 /* load the default table */
3131 snew(fr->atf_tabs, 1);
3132 fr->atf_tabs[DEFAULT_TF_TABLE] = make_atf_table(fp, oenv, fr, tabafn, box);
3137 fr->nwall = ir->nwall;
3138 if (ir->nwall && ir->wall_type == ewtTABLE)
3140 make_wall_tables(fp, oenv, ir, tabfn, &mtop->groups, fr);
3145 fcd->bondtab = make_bonded_tables(fp,
3146 F_TABBONDS, F_TABBONDSNC,
3148 fcd->angletab = make_bonded_tables(fp,
3151 fcd->dihtab = make_bonded_tables(fp,
3159 fprintf(debug, "No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
3163 /* QM/MM initialization if requested
3167 fprintf(stderr, "QM/MM calculation requested.\n");
3170 fr->bQMMM = ir->bQMMM;
3171 fr->qr = mk_QMMMrec();
3173 /* Set all the static charge group info */
3174 fr->cginfo_mb = init_cginfo_mb(fp, mtop, fr, bNoSolvOpt,
3176 &fr->bExcl_IntraCGAll_InterCGNone);
3177 if (DOMAINDECOMP(cr))
3183 fr->cginfo = cginfo_expand(mtop->nmolblock, fr->cginfo_mb);
3186 if (!DOMAINDECOMP(cr))
3188 forcerec_set_ranges(fr, ncg_mtop(mtop), ncg_mtop(mtop),
3189 mtop->natoms, mtop->natoms, mtop->natoms);
3192 fr->print_force = print_force;
3195 /* coarse load balancing vars */
3200 /* Initialize neighbor search */
3201 init_ns(fp, cr, &fr->ns, fr, mtop);
3203 if (cr->duty & DUTY_PP)
3205 gmx_nonbonded_setup(fr, bGenericKernelOnly);
3209 gmx_setup_adress_kernels(fp,bGenericKernelOnly);
3214 /* Initialize the thread working data for bonded interactions */
3215 init_bonded_threading(fp, mtop->groups.grps[egcENER].nr,
3216 &fr->bonded_threading);
3218 fr->nthread_ewc = gmx_omp_nthreads_get(emntBonded);
3219 snew(fr->ewc_t, fr->nthread_ewc);
3220 snew(fr->excl_load, fr->nthread_ewc + 1);
3222 /* fr->ic is used both by verlet and group kernels (to some extent) now */
3223 init_interaction_const(fp, &fr->ic, fr);
3224 init_interaction_const_tables(fp, fr->ic, rtab);
3226 if (fr->cutoff_scheme == ecutsVERLET)
3228 if (ir->rcoulomb != ir->rvdw)
3230 gmx_fatal(FARGS, "With Verlet lists rcoulomb and rvdw should be identical");
3233 init_nb_verlet(fp, &fr->nbv, bFEP_NonBonded, ir, fr, cr, nbpu_opt);
3236 if (ir->eDispCorr != edispcNO)
3238 calc_enervirdiff(fp, ir->eDispCorr, fr);
3242 #define pr_real(fp, r) fprintf(fp, "%s: %e\n",#r, r)
3243 #define pr_int(fp, i) fprintf((fp), "%s: %d\n",#i, i)
3244 #define pr_bool(fp, b) fprintf((fp), "%s: %s\n",#b, bool_names[b])
3246 void pr_forcerec(FILE *fp, t_forcerec *fr)
3250 pr_real(fp, fr->rlist);
3251 pr_real(fp, fr->rcoulomb);
3252 pr_real(fp, fr->fudgeQQ);
3253 pr_bool(fp, fr->bGrid);
3254 pr_bool(fp, fr->bTwinRange);
3255 /*pr_int(fp,fr->cg0);
3256 pr_int(fp,fr->hcg);*/
3257 for (i = 0; i < fr->nnblists; i++)
3259 pr_int(fp, fr->nblists[i].table_elec_vdw.n);
3261 pr_real(fp, fr->rcoulomb_switch);
3262 pr_real(fp, fr->rcoulomb);
3267 void forcerec_set_excl_load(t_forcerec *fr,
3268 const gmx_localtop_t *top)
3271 int t, i, j, ntot, n, ntarget;
3273 ind = top->excls.index;
3277 for (i = 0; i < top->excls.nr; i++)
3279 for (j = ind[i]; j < ind[i+1]; j++)
3288 fr->excl_load[0] = 0;
3291 for (t = 1; t <= fr->nthread_ewc; t++)
3293 ntarget = (ntot*t)/fr->nthread_ewc;
3294 while (i < top->excls.nr && n < ntarget)
3296 for (j = ind[i]; j < ind[i+1]; j++)
3305 fr->excl_load[t] = i;
3309 /* Frees GPU memory and destroys the GPU context.
3311 * Note that this function needs to be called even if GPUs are not used
3312 * in this run because the PME ranks have no knowledge of whether GPUs
3313 * are used or not, but all ranks need to enter the barrier below.
3315 void free_gpu_resources(const t_forcerec *fr,
3316 const t_commrec *cr,
3317 const gmx_gpu_info_t *gpu_info,
3318 const gmx_gpu_opt_t *gpu_opt)
3320 gmx_bool bIsPPrankUsingGPU;
3321 char gpu_err_str[STRLEN];
3323 bIsPPrankUsingGPU = (cr->duty & DUTY_PP) && fr && fr->nbv && fr->nbv->bUseGPU;
3325 if (bIsPPrankUsingGPU)
3327 /* free nbnxn data in GPU memory */
3328 nbnxn_gpu_free(fr->nbv->gpu_nbv);
3330 /* With tMPI we need to wait for all ranks to finish deallocation before
3331 * destroying the context in free_gpu() as some ranks may be sharing
3333 * Note: as only PP ranks need to free GPU resources, so it is safe to
3334 * not call the barrier on PME ranks.
3336 #ifdef GMX_THREAD_MPI
3341 #endif /* GMX_THREAD_MPI */
3343 /* uninitialize GPU (by destroying the context) */
3344 if (!free_cuda_gpu(cr->rank_pp_intranode, gpu_err_str, gpu_info, gpu_opt))
3346 gmx_warning("On rank %d failed to free GPU #%d: %s",
3347 cr->nodeid, get_current_cuda_gpu_device_id(), gpu_err_str);