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/fatalerror.h"
84 #include "gromacs/utility/smalloc.h"
86 #include "nbnxn_gpu_jit_support.h"
88 t_forcerec *mk_forcerec(void)
98 static void pr_nbfp(FILE *fp, real *nbfp, gmx_bool bBHAM, int atnr)
102 for (i = 0; (i < atnr); i++)
104 for (j = 0; (j < atnr); j++)
106 fprintf(fp, "%2d - %2d", i, j);
109 fprintf(fp, " a=%10g, b=%10g, c=%10g\n", BHAMA(nbfp, atnr, i, j),
110 BHAMB(nbfp, atnr, i, j), BHAMC(nbfp, atnr, i, j)/6.0);
114 fprintf(fp, " c6=%10g, c12=%10g\n", C6(nbfp, atnr, i, j)/6.0,
115 C12(nbfp, atnr, i, j)/12.0);
122 static real *mk_nbfp(const gmx_ffparams_t *idef, gmx_bool bBHAM)
130 snew(nbfp, 3*atnr*atnr);
131 for (i = k = 0; (i < atnr); i++)
133 for (j = 0; (j < atnr); j++, k++)
135 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
136 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
137 /* nbfp now includes the 6.0 derivative prefactor */
138 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c*6.0;
144 snew(nbfp, 2*atnr*atnr);
145 for (i = k = 0; (i < atnr); i++)
147 for (j = 0; (j < atnr); j++, k++)
149 /* nbfp now includes the 6.0/12.0 derivative prefactors */
150 C6(nbfp, atnr, i, j) = idef->iparams[k].lj.c6*6.0;
151 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12*12.0;
159 static real *make_ljpme_c6grid(const gmx_ffparams_t *idef, t_forcerec *fr)
162 real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
164 const real oneOverSix = 1.0 / 6.0;
166 /* For LJ-PME simulations, we correct the energies with the reciprocal space
167 * inside of the cut-off. To do this the non-bonded kernels needs to have
168 * access to the C6-values used on the reciprocal grid in pme.c
172 snew(grid, 2*atnr*atnr);
173 for (i = k = 0; (i < atnr); i++)
175 for (j = 0; (j < atnr); j++, k++)
177 c6i = idef->iparams[i*(atnr+1)].lj.c6;
178 c12i = idef->iparams[i*(atnr+1)].lj.c12;
179 c6j = idef->iparams[j*(atnr+1)].lj.c6;
180 c12j = idef->iparams[j*(atnr+1)].lj.c12;
181 c6 = sqrt(c6i * c6j);
182 if (fr->ljpme_combination_rule == eljpmeLB
183 && !gmx_numzero(c6) && !gmx_numzero(c12i) && !gmx_numzero(c12j))
185 sigmai = pow(c12i / c6i, oneOverSix);
186 sigmaj = pow(c12j / c6j, oneOverSix);
187 epsi = c6i * c6i / c12i;
188 epsj = c6j * c6j / c12j;
189 c6 = sqrt(epsi * epsj) * pow(0.5*(sigmai+sigmaj), 6);
191 /* Store the elements at the same relative positions as C6 in nbfp in order
192 * to simplify access in the kernels
194 grid[2*(atnr*i+j)] = c6*6.0;
200 static real *mk_nbfp_combination_rule(const gmx_ffparams_t *idef, int comb_rule)
204 real c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
206 const real oneOverSix = 1.0 / 6.0;
209 snew(nbfp, 2*atnr*atnr);
210 for (i = 0; i < atnr; ++i)
212 for (j = 0; j < atnr; ++j)
214 c6i = idef->iparams[i*(atnr+1)].lj.c6;
215 c12i = idef->iparams[i*(atnr+1)].lj.c12;
216 c6j = idef->iparams[j*(atnr+1)].lj.c6;
217 c12j = idef->iparams[j*(atnr+1)].lj.c12;
218 c6 = sqrt(c6i * c6j);
219 c12 = sqrt(c12i * c12j);
220 if (comb_rule == eCOMB_ARITHMETIC
221 && !gmx_numzero(c6) && !gmx_numzero(c12))
223 sigmai = pow(c12i / c6i, oneOverSix);
224 sigmaj = pow(c12j / c6j, oneOverSix);
225 epsi = c6i * c6i / c12i;
226 epsj = c6j * c6j / c12j;
227 c6 = sqrt(epsi * epsj) * pow(0.5*(sigmai+sigmaj), 6);
228 c12 = sqrt(epsi * epsj) * pow(0.5*(sigmai+sigmaj), 12);
230 C6(nbfp, atnr, i, j) = c6*6.0;
231 C12(nbfp, atnr, i, j) = c12*12.0;
237 /* This routine sets fr->solvent_opt to the most common solvent in the
238 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
239 * the fr->solvent_type array with the correct type (or esolNO).
241 * Charge groups that fulfill the conditions but are not identical to the
242 * most common one will be marked as esolNO in the solvent_type array.
244 * TIP3p is identical to SPC for these purposes, so we call it
245 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
247 * NOTE: QM particle should not
248 * become an optimized solvent. Not even if there is only one charge
258 } solvent_parameters_t;
261 check_solvent_cg(const gmx_moltype_t *molt,
264 const unsigned char *qm_grpnr,
265 const t_grps *qm_grps,
267 int *n_solvent_parameters,
268 solvent_parameters_t **solvent_parameters_p,
278 real tmp_charge[4] = { 0.0 }; /* init to zero to make gcc4.8 happy */
279 int tmp_vdwtype[4] = { 0 }; /* init to zero to make gcc4.8 happy */
282 solvent_parameters_t *solvent_parameters;
284 /* We use a list with parameters for each solvent type.
285 * Every time we discover a new molecule that fulfills the basic
286 * conditions for a solvent we compare with the previous entries
287 * in these lists. If the parameters are the same we just increment
288 * the counter for that type, and otherwise we create a new type
289 * based on the current molecule.
291 * Once we've finished going through all molecules we check which
292 * solvent is most common, and mark all those molecules while we
293 * clear the flag on all others.
296 solvent_parameters = *solvent_parameters_p;
298 /* Mark the cg first as non optimized */
301 /* Check if this cg has no exclusions with atoms in other charge groups
302 * and all atoms inside the charge group excluded.
303 * We only have 3 or 4 atom solvent loops.
305 if (GET_CGINFO_EXCL_INTER(cginfo) ||
306 !GET_CGINFO_EXCL_INTRA(cginfo))
311 /* Get the indices of the first atom in this charge group */
312 j0 = molt->cgs.index[cg0];
313 j1 = molt->cgs.index[cg0+1];
315 /* Number of atoms in our molecule */
321 "Moltype '%s': there are %d atoms in this charge group\n",
325 /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
328 if (nj < 3 || nj > 4)
333 /* Check if we are doing QM on this group */
335 if (qm_grpnr != NULL)
337 for (j = j0; j < j1 && !qm; j++)
339 qm = (qm_grpnr[j] < qm_grps->nr - 1);
342 /* Cannot use solvent optimization with QM */
348 atom = molt->atoms.atom;
350 /* Still looks like a solvent, time to check parameters */
352 /* If it is perturbed (free energy) we can't use the solvent loops,
353 * so then we just skip to the next molecule.
357 for (j = j0; j < j1 && !perturbed; j++)
359 perturbed = PERTURBED(atom[j]);
367 /* Now it's only a question if the VdW and charge parameters
368 * are OK. Before doing the check we compare and see if they are
369 * identical to a possible previous solvent type.
370 * First we assign the current types and charges.
372 for (j = 0; j < nj; j++)
374 tmp_vdwtype[j] = atom[j0+j].type;
375 tmp_charge[j] = atom[j0+j].q;
378 /* Does it match any previous solvent type? */
379 for (k = 0; k < *n_solvent_parameters; k++)
384 /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
385 if ( (solvent_parameters[k].model == esolSPC && nj != 3) ||
386 (solvent_parameters[k].model == esolTIP4P && nj != 4) )
391 /* Check that types & charges match for all atoms in molecule */
392 for (j = 0; j < nj && match == TRUE; j++)
394 if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
398 if (tmp_charge[j] != solvent_parameters[k].charge[j])
405 /* Congratulations! We have a matched solvent.
406 * Flag it with this type for later processing.
409 solvent_parameters[k].count += nmol;
411 /* We are done with this charge group */
416 /* If we get here, we have a tentative new solvent type.
417 * Before we add it we must check that it fulfills the requirements
418 * of the solvent optimized loops. First determine which atoms have
421 for (j = 0; j < nj; j++)
424 tjA = tmp_vdwtype[j];
426 /* Go through all other tpes and see if any have non-zero
427 * VdW parameters when combined with this one.
429 for (k = 0; k < fr->ntype && (has_vdw[j] == FALSE); k++)
431 /* We already checked that the atoms weren't perturbed,
432 * so we only need to check state A now.
436 has_vdw[j] = (has_vdw[j] ||
437 (BHAMA(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
438 (BHAMB(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
439 (BHAMC(fr->nbfp, fr->ntype, tjA, k) != 0.0));
444 has_vdw[j] = (has_vdw[j] ||
445 (C6(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
446 (C12(fr->nbfp, fr->ntype, tjA, k) != 0.0));
451 /* Now we know all we need to make the final check and assignment. */
455 * For this we require thatn all atoms have charge,
456 * the charges on atom 2 & 3 should be the same, and only
457 * atom 1 might have VdW.
459 if (has_vdw[1] == FALSE &&
460 has_vdw[2] == FALSE &&
461 tmp_charge[0] != 0 &&
462 tmp_charge[1] != 0 &&
463 tmp_charge[2] == tmp_charge[1])
465 srenew(solvent_parameters, *n_solvent_parameters+1);
466 solvent_parameters[*n_solvent_parameters].model = esolSPC;
467 solvent_parameters[*n_solvent_parameters].count = nmol;
468 for (k = 0; k < 3; k++)
470 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
471 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
474 *cg_sp = *n_solvent_parameters;
475 (*n_solvent_parameters)++;
480 /* Or could it be a TIP4P?
481 * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
482 * Only atom 1 mght have VdW.
484 if (has_vdw[1] == FALSE &&
485 has_vdw[2] == FALSE &&
486 has_vdw[3] == FALSE &&
487 tmp_charge[0] == 0 &&
488 tmp_charge[1] != 0 &&
489 tmp_charge[2] == tmp_charge[1] &&
492 srenew(solvent_parameters, *n_solvent_parameters+1);
493 solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
494 solvent_parameters[*n_solvent_parameters].count = nmol;
495 for (k = 0; k < 4; k++)
497 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
498 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
501 *cg_sp = *n_solvent_parameters;
502 (*n_solvent_parameters)++;
506 *solvent_parameters_p = solvent_parameters;
510 check_solvent(FILE * fp,
511 const gmx_mtop_t * mtop,
513 cginfo_mb_t *cginfo_mb)
516 const gmx_moltype_t *molt;
517 int mb, mol, cg_mol, at_offset, am, cgm, i, nmol_ch, nmol;
518 int n_solvent_parameters;
519 solvent_parameters_t *solvent_parameters;
525 fprintf(debug, "Going to determine what solvent types we have.\n");
528 n_solvent_parameters = 0;
529 solvent_parameters = NULL;
530 /* Allocate temporary array for solvent type */
531 snew(cg_sp, mtop->nmolblock);
534 for (mb = 0; mb < mtop->nmolblock; mb++)
536 molt = &mtop->moltype[mtop->molblock[mb].type];
538 /* Here we have to loop over all individual molecules
539 * because we need to check for QMMM particles.
541 snew(cg_sp[mb], cginfo_mb[mb].cg_mod);
542 nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
543 nmol = mtop->molblock[mb].nmol/nmol_ch;
544 for (mol = 0; mol < nmol_ch; mol++)
547 am = mol*cgs->index[cgs->nr];
548 for (cg_mol = 0; cg_mol < cgs->nr; cg_mol++)
550 check_solvent_cg(molt, cg_mol, nmol,
551 mtop->groups.grpnr[egcQMMM] ?
552 mtop->groups.grpnr[egcQMMM]+at_offset+am : 0,
553 &mtop->groups.grps[egcQMMM],
555 &n_solvent_parameters, &solvent_parameters,
556 cginfo_mb[mb].cginfo[cgm+cg_mol],
557 &cg_sp[mb][cgm+cg_mol]);
560 at_offset += cgs->index[cgs->nr];
563 /* Puh! We finished going through all charge groups.
564 * Now find the most common solvent model.
567 /* Most common solvent this far */
569 for (i = 0; i < n_solvent_parameters; i++)
572 solvent_parameters[i].count > solvent_parameters[bestsp].count)
580 bestsol = solvent_parameters[bestsp].model;
588 for (mb = 0; mb < mtop->nmolblock; mb++)
590 cgs = &mtop->moltype[mtop->molblock[mb].type].cgs;
591 nmol = (mtop->molblock[mb].nmol*cgs->nr)/cginfo_mb[mb].cg_mod;
592 for (i = 0; i < cginfo_mb[mb].cg_mod; i++)
594 if (cg_sp[mb][i] == bestsp)
596 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], bestsol);
601 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], esolNO);
608 if (bestsol != esolNO && fp != NULL)
610 fprintf(fp, "\nEnabling %s-like water optimization for %d molecules.\n\n",
612 solvent_parameters[bestsp].count);
615 sfree(solvent_parameters);
616 fr->solvent_opt = bestsol;
620 acNONE = 0, acCONSTRAINT, acSETTLE
623 static cginfo_mb_t *init_cginfo_mb(FILE *fplog, const gmx_mtop_t *mtop,
624 t_forcerec *fr, gmx_bool bNoSolvOpt,
625 gmx_bool *bFEP_NonBonded,
626 gmx_bool *bExcl_IntraCGAll_InterCGNone)
629 const t_blocka *excl;
630 const gmx_moltype_t *molt;
631 const gmx_molblock_t *molb;
632 cginfo_mb_t *cginfo_mb;
635 int cg_offset, a_offset;
636 int mb, m, cg, a0, a1, gid, ai, j, aj, excl_nalloc;
640 gmx_bool bId, *bExcl, bExclIntraAll, bExclInter, bHaveVDW, bHaveQ, bHavePerturbedAtoms;
642 snew(cginfo_mb, mtop->nmolblock);
644 snew(type_VDW, fr->ntype);
645 for (ai = 0; ai < fr->ntype; ai++)
647 type_VDW[ai] = FALSE;
648 for (j = 0; j < fr->ntype; j++)
650 type_VDW[ai] = type_VDW[ai] ||
652 C6(fr->nbfp, fr->ntype, ai, j) != 0 ||
653 C12(fr->nbfp, fr->ntype, ai, j) != 0;
657 *bFEP_NonBonded = FALSE;
658 *bExcl_IntraCGAll_InterCGNone = TRUE;
661 snew(bExcl, excl_nalloc);
664 for (mb = 0; mb < mtop->nmolblock; mb++)
666 molb = &mtop->molblock[mb];
667 molt = &mtop->moltype[molb->type];
671 /* Check if the cginfo is identical for all molecules in this block.
672 * If so, we only need an array of the size of one molecule.
673 * Otherwise we make an array of #mol times #cgs per molecule.
676 for (m = 0; m < molb->nmol; m++)
678 int am = m*cgs->index[cgs->nr];
679 for (cg = 0; cg < cgs->nr; cg++)
682 a1 = cgs->index[cg+1];
683 if (ggrpnr(&mtop->groups, egcENER, a_offset+am+a0) !=
684 ggrpnr(&mtop->groups, egcENER, a_offset +a0))
688 if (mtop->groups.grpnr[egcQMMM] != NULL)
690 for (ai = a0; ai < a1; ai++)
692 if (mtop->groups.grpnr[egcQMMM][a_offset+am+ai] !=
693 mtop->groups.grpnr[egcQMMM][a_offset +ai])
702 cginfo_mb[mb].cg_start = cg_offset;
703 cginfo_mb[mb].cg_end = cg_offset + molb->nmol*cgs->nr;
704 cginfo_mb[mb].cg_mod = (bId ? 1 : molb->nmol)*cgs->nr;
705 snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
706 cginfo = cginfo_mb[mb].cginfo;
708 /* Set constraints flags for constrained atoms */
709 snew(a_con, molt->atoms.nr);
710 for (ftype = 0; ftype < F_NRE; ftype++)
712 if (interaction_function[ftype].flags & IF_CONSTRAINT)
717 for (ia = 0; ia < molt->ilist[ftype].nr; ia += 1+nral)
721 for (a = 0; a < nral; a++)
723 a_con[molt->ilist[ftype].iatoms[ia+1+a]] =
724 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
730 for (m = 0; m < (bId ? 1 : molb->nmol); m++)
733 int am = m*cgs->index[cgs->nr];
734 for (cg = 0; cg < cgs->nr; cg++)
737 a1 = cgs->index[cg+1];
739 /* Store the energy group in cginfo */
740 gid = ggrpnr(&mtop->groups, egcENER, a_offset+am+a0);
741 SET_CGINFO_GID(cginfo[cgm+cg], gid);
743 /* Check the intra/inter charge group exclusions */
744 if (a1-a0 > excl_nalloc)
746 excl_nalloc = a1 - a0;
747 srenew(bExcl, excl_nalloc);
749 /* bExclIntraAll: all intra cg interactions excluded
750 * bExclInter: any inter cg interactions excluded
752 bExclIntraAll = TRUE;
756 bHavePerturbedAtoms = FALSE;
757 for (ai = a0; ai < a1; ai++)
759 /* Check VDW and electrostatic interactions */
760 bHaveVDW = bHaveVDW || (type_VDW[molt->atoms.atom[ai].type] ||
761 type_VDW[molt->atoms.atom[ai].typeB]);
762 bHaveQ = bHaveQ || (molt->atoms.atom[ai].q != 0 ||
763 molt->atoms.atom[ai].qB != 0);
765 bHavePerturbedAtoms = bHavePerturbedAtoms || (PERTURBED(molt->atoms.atom[ai]) != 0);
767 /* Clear the exclusion list for atom ai */
768 for (aj = a0; aj < a1; aj++)
770 bExcl[aj-a0] = FALSE;
772 /* Loop over all the exclusions of atom ai */
773 for (j = excl->index[ai]; j < excl->index[ai+1]; j++)
776 if (aj < a0 || aj >= a1)
785 /* Check if ai excludes a0 to a1 */
786 for (aj = a0; aj < a1; aj++)
790 bExclIntraAll = FALSE;
797 SET_CGINFO_CONSTR(cginfo[cgm+cg]);
800 SET_CGINFO_SETTLE(cginfo[cgm+cg]);
808 SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
812 SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
814 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
816 /* The size in cginfo is currently only read with DD */
817 gmx_fatal(FARGS, "A charge group has size %d which is larger than the limit of %d atoms", a1-a0, MAX_CHARGEGROUP_SIZE);
821 SET_CGINFO_HAS_VDW(cginfo[cgm+cg]);
825 SET_CGINFO_HAS_Q(cginfo[cgm+cg]);
827 if (bHavePerturbedAtoms && fr->efep != efepNO)
829 SET_CGINFO_FEP(cginfo[cgm+cg]);
830 *bFEP_NonBonded = TRUE;
832 /* Store the charge group size */
833 SET_CGINFO_NATOMS(cginfo[cgm+cg], a1-a0);
835 if (!bExclIntraAll || bExclInter)
837 *bExcl_IntraCGAll_InterCGNone = FALSE;
844 cg_offset += molb->nmol*cgs->nr;
845 a_offset += molb->nmol*cgs->index[cgs->nr];
849 /* the solvent optimizer is called after the QM is initialized,
850 * because we don't want to have the QM subsystemto become an
854 check_solvent(fplog, mtop, fr, cginfo_mb);
856 if (getenv("GMX_NO_SOLV_OPT"))
860 fprintf(fplog, "Found environment variable GMX_NO_SOLV_OPT.\n"
861 "Disabling all solvent optimization\n");
863 fr->solvent_opt = esolNO;
867 fr->solvent_opt = esolNO;
869 if (!fr->solvent_opt)
871 for (mb = 0; mb < mtop->nmolblock; mb++)
873 for (cg = 0; cg < cginfo_mb[mb].cg_mod; cg++)
875 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg], esolNO);
883 static int *cginfo_expand(int nmb, cginfo_mb_t *cgi_mb)
888 ncg = cgi_mb[nmb-1].cg_end;
891 for (cg = 0; cg < ncg; cg++)
893 while (cg >= cgi_mb[mb].cg_end)
898 cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
904 static void set_chargesum(FILE *log, t_forcerec *fr, const gmx_mtop_t *mtop)
906 /*This now calculates sum for q and c6*/
907 double qsum, q2sum, q, c6sum, c6;
909 const t_atoms *atoms;
914 for (mb = 0; mb < mtop->nmolblock; mb++)
916 nmol = mtop->molblock[mb].nmol;
917 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
918 for (i = 0; i < atoms->nr; i++)
920 q = atoms->atom[i].q;
923 c6 = mtop->ffparams.iparams[atoms->atom[i].type*(mtop->ffparams.atnr+1)].lj.c6;
928 fr->q2sum[0] = q2sum;
929 fr->c6sum[0] = c6sum;
931 if (fr->efep != efepNO)
936 for (mb = 0; mb < mtop->nmolblock; mb++)
938 nmol = mtop->molblock[mb].nmol;
939 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
940 for (i = 0; i < atoms->nr; i++)
942 q = atoms->atom[i].qB;
945 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB*(mtop->ffparams.atnr+1)].lj.c6;
949 fr->q2sum[1] = q2sum;
950 fr->c6sum[1] = c6sum;
955 fr->qsum[1] = fr->qsum[0];
956 fr->q2sum[1] = fr->q2sum[0];
957 fr->c6sum[1] = fr->c6sum[0];
961 if (fr->efep == efepNO)
963 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
967 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n",
968 fr->qsum[0], fr->qsum[1]);
973 void update_forcerec(t_forcerec *fr, matrix box)
975 if (fr->eeltype == eelGRF)
977 calc_rffac(NULL, fr->eeltype, fr->epsilon_r, fr->epsilon_rf,
978 fr->rcoulomb, fr->temp, fr->zsquare, box,
979 &fr->kappa, &fr->k_rf, &fr->c_rf);
983 void set_avcsixtwelve(FILE *fplog, t_forcerec *fr, const gmx_mtop_t *mtop)
985 const t_atoms *atoms, *atoms_tpi;
986 const t_blocka *excl;
987 int mb, nmol, nmolc, i, j, tpi, tpj, j1, j2, k, nexcl, q;
988 gmx_int64_t npair, npair_ij, tmpi, tmpj;
989 double csix, ctwelve;
993 real *nbfp_comb = NULL;
999 /* For LJ-PME, we want to correct for the difference between the
1000 * actual C6 values and the C6 values used by the LJ-PME based on
1001 * combination rules. */
1003 if (EVDW_PME(fr->vdwtype))
1005 nbfp_comb = mk_nbfp_combination_rule(&mtop->ffparams,
1006 (fr->ljpme_combination_rule == eljpmeLB) ? eCOMB_ARITHMETIC : eCOMB_GEOMETRIC);
1007 for (tpi = 0; tpi < ntp; ++tpi)
1009 for (tpj = 0; tpj < ntp; ++tpj)
1011 C6(nbfp_comb, ntp, tpi, tpj) =
1012 C6(nbfp, ntp, tpi, tpj) - C6(nbfp_comb, ntp, tpi, tpj);
1013 C12(nbfp_comb, ntp, tpi, tpj) = C12(nbfp, ntp, tpi, tpj);
1018 for (q = 0; q < (fr->efep == efepNO ? 1 : 2); q++)
1026 /* Count the types so we avoid natoms^2 operations */
1027 snew(typecount, ntp);
1028 gmx_mtop_count_atomtypes(mtop, q, typecount);
1030 for (tpi = 0; tpi < ntp; tpi++)
1032 for (tpj = tpi; tpj < ntp; tpj++)
1034 tmpi = typecount[tpi];
1035 tmpj = typecount[tpj];
1038 npair_ij = tmpi*tmpj;
1042 npair_ij = tmpi*(tmpi - 1)/2;
1046 /* nbfp now includes the 6.0 derivative prefactor */
1047 csix += npair_ij*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1051 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1052 csix += npair_ij* C6(nbfp, ntp, tpi, tpj)/6.0;
1053 ctwelve += npair_ij* C12(nbfp, ntp, tpi, tpj)/12.0;
1059 /* Subtract the excluded pairs.
1060 * The main reason for substracting exclusions is that in some cases
1061 * some combinations might never occur and the parameters could have
1062 * any value. These unused values should not influence the dispersion
1065 for (mb = 0; mb < mtop->nmolblock; mb++)
1067 nmol = mtop->molblock[mb].nmol;
1068 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
1069 excl = &mtop->moltype[mtop->molblock[mb].type].excls;
1070 for (i = 0; (i < atoms->nr); i++)
1074 tpi = atoms->atom[i].type;
1078 tpi = atoms->atom[i].typeB;
1080 j1 = excl->index[i];
1081 j2 = excl->index[i+1];
1082 for (j = j1; j < j2; j++)
1089 tpj = atoms->atom[k].type;
1093 tpj = atoms->atom[k].typeB;
1097 /* nbfp now includes the 6.0 derivative prefactor */
1098 csix -= nmol*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1102 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1103 csix -= nmol*C6 (nbfp, ntp, tpi, tpj)/6.0;
1104 ctwelve -= nmol*C12(nbfp, ntp, tpi, tpj)/12.0;
1114 /* Only correct for the interaction of the test particle
1115 * with the rest of the system.
1118 &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].atoms;
1121 for (mb = 0; mb < mtop->nmolblock; mb++)
1123 nmol = mtop->molblock[mb].nmol;
1124 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
1125 for (j = 0; j < atoms->nr; j++)
1128 /* Remove the interaction of the test charge group
1131 if (mb == mtop->nmolblock-1)
1135 if (mb == 0 && nmol == 1)
1137 gmx_fatal(FARGS, "Old format tpr with TPI, please generate a new tpr file");
1142 tpj = atoms->atom[j].type;
1146 tpj = atoms->atom[j].typeB;
1148 for (i = 0; i < fr->n_tpi; i++)
1152 tpi = atoms_tpi->atom[i].type;
1156 tpi = atoms_tpi->atom[i].typeB;
1160 /* nbfp now includes the 6.0 derivative prefactor */
1161 csix += nmolc*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1165 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1166 csix += nmolc*C6 (nbfp, ntp, tpi, tpj)/6.0;
1167 ctwelve += nmolc*C12(nbfp, ntp, tpi, tpj)/12.0;
1174 if (npair - nexcl <= 0 && fplog)
1176 fprintf(fplog, "\nWARNING: There are no atom pairs for dispersion correction\n\n");
1182 csix /= npair - nexcl;
1183 ctwelve /= npair - nexcl;
1187 fprintf(debug, "Counted %d exclusions\n", nexcl);
1188 fprintf(debug, "Average C6 parameter is: %10g\n", (double)csix);
1189 fprintf(debug, "Average C12 parameter is: %10g\n", (double)ctwelve);
1191 fr->avcsix[q] = csix;
1192 fr->avctwelve[q] = ctwelve;
1195 if (EVDW_PME(fr->vdwtype))
1202 if (fr->eDispCorr == edispcAllEner ||
1203 fr->eDispCorr == edispcAllEnerPres)
1205 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1206 fr->avcsix[0], fr->avctwelve[0]);
1210 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e\n", fr->avcsix[0]);
1216 static void set_bham_b_max(FILE *fplog, t_forcerec *fr,
1217 const gmx_mtop_t *mtop)
1219 const t_atoms *at1, *at2;
1220 int mt1, mt2, i, j, tpi, tpj, ntypes;
1226 fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
1233 for (mt1 = 0; mt1 < mtop->nmoltype; mt1++)
1235 at1 = &mtop->moltype[mt1].atoms;
1236 for (i = 0; (i < at1->nr); i++)
1238 tpi = at1->atom[i].type;
1241 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
1244 for (mt2 = mt1; mt2 < mtop->nmoltype; mt2++)
1246 at2 = &mtop->moltype[mt2].atoms;
1247 for (j = 0; (j < at2->nr); j++)
1249 tpj = at2->atom[j].type;
1252 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
1254 b = BHAMB(nbfp, ntypes, tpi, tpj);
1255 if (b > fr->bham_b_max)
1259 if ((b < bmin) || (bmin == -1))
1269 fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n",
1270 bmin, fr->bham_b_max);
1274 static void make_nbf_tables(FILE *fp, const output_env_t oenv,
1275 t_forcerec *fr, real rtab,
1276 const t_commrec *cr,
1277 const char *tabfn, char *eg1, char *eg2,
1287 fprintf(debug, "No table file name passed, can not read table, can not do non-bonded interactions\n");
1292 sprintf(buf, "%s", tabfn);
1295 /* Append the two energy group names */
1296 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
1297 eg1, eg2, ftp2ext(efXVG));
1299 nbl->table_elec_vdw = make_tables(fp, oenv, fr, MASTER(cr), buf, rtab, 0);
1300 /* Copy the contents of the table to separate coulomb and LJ tables too,
1301 * to improve cache performance.
1303 /* For performance reasons we want
1304 * the table data to be aligned to 16-byte. The pointers could be freed
1305 * but currently aren't.
1307 nbl->table_elec.interaction = GMX_TABLE_INTERACTION_ELEC;
1308 nbl->table_elec.format = nbl->table_elec_vdw.format;
1309 nbl->table_elec.r = nbl->table_elec_vdw.r;
1310 nbl->table_elec.n = nbl->table_elec_vdw.n;
1311 nbl->table_elec.scale = nbl->table_elec_vdw.scale;
1312 nbl->table_elec.scale_exp = nbl->table_elec_vdw.scale_exp;
1313 nbl->table_elec.formatsize = nbl->table_elec_vdw.formatsize;
1314 nbl->table_elec.ninteractions = 1;
1315 nbl->table_elec.stride = nbl->table_elec.formatsize * nbl->table_elec.ninteractions;
1316 snew_aligned(nbl->table_elec.data, nbl->table_elec.stride*(nbl->table_elec.n+1), 32);
1318 nbl->table_vdw.interaction = GMX_TABLE_INTERACTION_VDWREP_VDWDISP;
1319 nbl->table_vdw.format = nbl->table_elec_vdw.format;
1320 nbl->table_vdw.r = nbl->table_elec_vdw.r;
1321 nbl->table_vdw.n = nbl->table_elec_vdw.n;
1322 nbl->table_vdw.scale = nbl->table_elec_vdw.scale;
1323 nbl->table_vdw.scale_exp = nbl->table_elec_vdw.scale_exp;
1324 nbl->table_vdw.formatsize = nbl->table_elec_vdw.formatsize;
1325 nbl->table_vdw.ninteractions = 2;
1326 nbl->table_vdw.stride = nbl->table_vdw.formatsize * nbl->table_vdw.ninteractions;
1327 snew_aligned(nbl->table_vdw.data, nbl->table_vdw.stride*(nbl->table_vdw.n+1), 32);
1329 for (i = 0; i <= nbl->table_elec_vdw.n; i++)
1331 for (j = 0; j < 4; j++)
1333 nbl->table_elec.data[4*i+j] = nbl->table_elec_vdw.data[12*i+j];
1335 for (j = 0; j < 8; j++)
1337 nbl->table_vdw.data[8*i+j] = nbl->table_elec_vdw.data[12*i+4+j];
1342 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t *mtop,
1343 int *ncount, int **count)
1345 const gmx_moltype_t *molt;
1347 int mt, ftype, stride, i, j, tabnr;
1349 for (mt = 0; mt < mtop->nmoltype; mt++)
1351 molt = &mtop->moltype[mt];
1352 for (ftype = 0; ftype < F_NRE; ftype++)
1354 if (ftype == ftype1 || ftype == ftype2)
1356 il = &molt->ilist[ftype];
1357 stride = 1 + NRAL(ftype);
1358 for (i = 0; i < il->nr; i += stride)
1360 tabnr = mtop->ffparams.iparams[il->iatoms[i]].tab.table;
1363 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
1365 if (tabnr >= *ncount)
1367 srenew(*count, tabnr+1);
1368 for (j = *ncount; j < tabnr+1; j++)
1381 static bondedtable_t *make_bonded_tables(FILE *fplog,
1382 int ftype1, int ftype2,
1383 const gmx_mtop_t *mtop,
1384 const char *basefn, const char *tabext)
1386 int i, ncount, *count;
1394 count_tables(ftype1, ftype2, mtop, &ncount, &count);
1399 for (i = 0; i < ncount; i++)
1403 sprintf(tabfn, "%s", basefn);
1404 sprintf(tabfn + strlen(basefn) - strlen(ftp2ext(efXVG)) - 1, "_%s%d.%s",
1405 tabext, i, ftp2ext(efXVG));
1406 tab[i] = make_bonded_table(fplog, tabfn, NRAL(ftype1)-2);
1415 void forcerec_set_ranges(t_forcerec *fr,
1416 int ncg_home, int ncg_force,
1418 int natoms_force_constr, int natoms_f_novirsum)
1423 /* fr->ncg_force is unused in the standard code,
1424 * but it can be useful for modified code dealing with charge groups.
1426 fr->ncg_force = ncg_force;
1427 fr->natoms_force = natoms_force;
1428 fr->natoms_force_constr = natoms_force_constr;
1430 if (fr->natoms_force_constr > fr->nalloc_force)
1432 fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1436 srenew(fr->f_twin, fr->nalloc_force);
1440 if (fr->bF_NoVirSum)
1442 fr->f_novirsum_n = natoms_f_novirsum;
1443 if (fr->f_novirsum_n > fr->f_novirsum_nalloc)
1445 fr->f_novirsum_nalloc = over_alloc_dd(fr->f_novirsum_n);
1446 srenew(fr->f_novirsum_alloc, fr->f_novirsum_nalloc);
1451 fr->f_novirsum_n = 0;
1455 static real cutoff_inf(real cutoff)
1459 cutoff = GMX_CUTOFF_INF;
1465 static void make_adress_tf_tables(FILE *fp, const output_env_t oenv,
1466 t_forcerec *fr, const t_inputrec *ir,
1467 const char *tabfn, const gmx_mtop_t *mtop,
1475 gmx_fatal(FARGS, "No thermoforce table file given. Use -tabletf to specify a file\n");
1479 snew(fr->atf_tabs, ir->adress->n_tf_grps);
1481 sprintf(buf, "%s", tabfn);
1482 for (i = 0; i < ir->adress->n_tf_grps; i++)
1484 j = ir->adress->tf_table_index[i]; /* get energy group index */
1485 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "tf_%s.%s",
1486 *(mtop->groups.grpname[mtop->groups.grps[egcENER].nm_ind[j]]), ftp2ext(efXVG));
1489 fprintf(fp, "loading tf table for energygrp index %d from %s\n", ir->adress->tf_table_index[i], buf);
1491 fr->atf_tabs[i] = make_atf_table(fp, oenv, fr, buf, box);
1496 gmx_bool can_use_allvsall(const t_inputrec *ir, gmx_bool bPrintNote, t_commrec *cr, FILE *fp)
1503 ir->rcoulomb == 0 &&
1505 ir->ePBC == epbcNONE &&
1506 ir->vdwtype == evdwCUT &&
1507 ir->coulombtype == eelCUT &&
1508 ir->efep == efepNO &&
1509 (ir->implicit_solvent == eisNO ||
1510 (ir->implicit_solvent == eisGBSA && (ir->gb_algorithm == egbSTILL ||
1511 ir->gb_algorithm == egbHCT ||
1512 ir->gb_algorithm == egbOBC))) &&
1513 getenv("GMX_NO_ALLVSALL") == NULL
1516 if (bAllvsAll && ir->opts.ngener > 1)
1518 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";
1524 fprintf(stderr, "\n%s\n", note);
1528 fprintf(fp, "\n%s\n", note);
1534 if (bAllvsAll && fp && MASTER(cr))
1536 fprintf(fp, "\nUsing SIMD all-vs-all kernels.\n\n");
1543 gmx_bool nbnxn_acceleration_supported(FILE *fplog,
1544 const t_commrec *cr,
1545 const t_inputrec *ir,
1548 if (!bGPU && (ir->vdwtype == evdwPME && ir->ljpme_combination_rule == eljpmeLB))
1550 md_print_warn(cr, fplog, "LJ-PME with Lorentz-Berthelot is not supported with %s, falling back to %s\n",
1551 bGPU ? "GPUs" : "SIMD kernels",
1552 bGPU ? "CPU only" : "plain-C kernels");
1560 static void pick_nbnxn_kernel_cpu(const t_inputrec gmx_unused *ir,
1564 *kernel_type = nbnxnk4x4_PlainC;
1565 *ewald_excl = ewaldexclTable;
1567 #ifdef GMX_NBNXN_SIMD
1569 #ifdef GMX_NBNXN_SIMD_4XN
1570 *kernel_type = nbnxnk4xN_SIMD_4xN;
1572 #ifdef GMX_NBNXN_SIMD_2XNN
1573 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1576 #if defined GMX_NBNXN_SIMD_2XNN && defined GMX_NBNXN_SIMD_4XN
1577 /* We need to choose if we want 2x(N+N) or 4xN kernels.
1578 * Currently this is based on the SIMD acceleration choice,
1579 * but it might be better to decide this at runtime based on CPU.
1581 * 4xN calculates more (zero) interactions, but has less pair-search
1582 * work and much better kernel instruction scheduling.
1584 * Up till now we have only seen that on Intel Sandy/Ivy Bridge,
1585 * which doesn't have FMA, both the analytical and tabulated Ewald
1586 * kernels have similar pair rates for 4x8 and 2x(4+4), so we choose
1587 * 2x(4+4) because it results in significantly fewer pairs.
1588 * For RF, the raw pair rate of the 4x8 kernel is higher than 2x(4+4),
1589 * 10% with HT, 50% without HT. As we currently don't detect the actual
1590 * use of HT, use 4x8 to avoid a potential performance hit.
1591 * On Intel Haswell 4x8 is always faster.
1593 *kernel_type = nbnxnk4xN_SIMD_4xN;
1595 #ifndef GMX_SIMD_HAVE_FMA
1596 if (EEL_PME_EWALD(ir->coulombtype) ||
1597 EVDW_PME(ir->vdwtype))
1599 /* We have Ewald kernels without FMA (Intel Sandy/Ivy Bridge).
1600 * There are enough instructions to make 2x(4+4) efficient.
1602 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1605 #endif /* GMX_NBNXN_SIMD_2XNN && GMX_NBNXN_SIMD_4XN */
1608 if (getenv("GMX_NBNXN_SIMD_4XN") != NULL)
1610 #ifdef GMX_NBNXN_SIMD_4XN
1611 *kernel_type = nbnxnk4xN_SIMD_4xN;
1613 gmx_fatal(FARGS, "SIMD 4xN kernels requested, but GROMACS has been compiled without support for these kernels");
1616 if (getenv("GMX_NBNXN_SIMD_2XNN") != NULL)
1618 #ifdef GMX_NBNXN_SIMD_2XNN
1619 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1621 gmx_fatal(FARGS, "SIMD 2x(N+N) kernels requested, but GROMACS has been compiled without support for these kernels");
1625 /* Analytical Ewald exclusion correction is only an option in
1627 * Since table lookup's don't parallelize with SIMD, analytical
1628 * will probably always be faster for a SIMD width of 8 or more.
1629 * With FMA analytical is sometimes faster for a width if 4 as well.
1630 * On BlueGene/Q, this is faster regardless of precision.
1631 * In single precision, this is faster on Bulldozer.
1633 #if GMX_SIMD_REAL_WIDTH >= 8 || \
1634 (GMX_SIMD_REAL_WIDTH >= 4 && defined GMX_SIMD_HAVE_FMA && !defined GMX_DOUBLE) || \
1635 defined GMX_SIMD_IBM_QPX
1636 *ewald_excl = ewaldexclAnalytical;
1638 if (getenv("GMX_NBNXN_EWALD_TABLE") != NULL)
1640 *ewald_excl = ewaldexclTable;
1642 if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != NULL)
1644 *ewald_excl = ewaldexclAnalytical;
1648 #endif /* GMX_NBNXN_SIMD */
1652 const char *lookup_nbnxn_kernel_name(int kernel_type)
1654 const char *returnvalue = NULL;
1655 switch (kernel_type)
1658 returnvalue = "not set";
1660 case nbnxnk4x4_PlainC:
1661 returnvalue = "plain C";
1663 case nbnxnk4xN_SIMD_4xN:
1664 case nbnxnk4xN_SIMD_2xNN:
1665 #ifdef GMX_NBNXN_SIMD
1666 #if defined GMX_SIMD_X86_SSE2
1667 returnvalue = "SSE2";
1668 #elif defined GMX_SIMD_X86_SSE4_1
1669 returnvalue = "SSE4.1";
1670 #elif defined GMX_SIMD_X86_AVX_128_FMA
1671 returnvalue = "AVX_128_FMA";
1672 #elif defined GMX_SIMD_X86_AVX_256
1673 returnvalue = "AVX_256";
1674 #elif defined GMX_SIMD_X86_AVX2_256
1675 returnvalue = "AVX2_256";
1677 returnvalue = "SIMD";
1679 #else /* GMX_NBNXN_SIMD */
1680 returnvalue = "not available";
1681 #endif /* GMX_NBNXN_SIMD */
1683 case nbnxnk8x8x8_GPU: returnvalue = "GPU"; break;
1684 case nbnxnk8x8x8_PlainC: returnvalue = "plain C"; break;
1688 gmx_fatal(FARGS, "Illegal kernel type selected");
1695 static void pick_nbnxn_kernel(FILE *fp,
1696 const t_commrec *cr,
1697 gmx_bool use_simd_kernels,
1699 gmx_bool bEmulateGPU,
1700 const t_inputrec *ir,
1703 gmx_bool bDoNonbonded)
1705 assert(kernel_type);
1707 *kernel_type = nbnxnkNotSet;
1708 *ewald_excl = ewaldexclTable;
1712 *kernel_type = nbnxnk8x8x8_PlainC;
1716 md_print_warn(cr, fp, "Emulating a GPU run on the CPU (slow)");
1721 *kernel_type = nbnxnk8x8x8_GPU;
1724 if (*kernel_type == nbnxnkNotSet)
1726 /* LJ PME with LB combination rule does 7 mesh operations.
1727 * This so slow that we don't compile SIMD non-bonded kernels for that.
1729 if (use_simd_kernels &&
1730 nbnxn_acceleration_supported(fp, cr, ir, FALSE))
1732 pick_nbnxn_kernel_cpu(ir, kernel_type, ewald_excl);
1736 *kernel_type = nbnxnk4x4_PlainC;
1740 if (bDoNonbonded && fp != NULL)
1742 fprintf(fp, "\nUsing %s %dx%d non-bonded kernels\n\n",
1743 lookup_nbnxn_kernel_name(*kernel_type),
1744 nbnxn_kernel_to_ci_size(*kernel_type),
1745 nbnxn_kernel_to_cj_size(*kernel_type));
1747 if (nbnxnk4x4_PlainC == *kernel_type ||
1748 nbnxnk8x8x8_PlainC == *kernel_type)
1750 md_print_warn(cr, fp,
1751 "WARNING: Using the slow %s kernels. This should\n"
1752 "not happen during routine usage on supported platforms.\n\n",
1753 lookup_nbnxn_kernel_name(*kernel_type));
1758 static void pick_nbnxn_resources(FILE *fp,
1759 const t_commrec *cr,
1760 const gmx_hw_info_t *hwinfo,
1761 gmx_bool bDoNonbonded,
1763 gmx_bool *bEmulateGPU,
1764 const gmx_gpu_opt_t *gpu_opt)
1766 gmx_bool bEmulateGPUEnvVarSet;
1767 char gpu_err_str[STRLEN];
1771 bEmulateGPUEnvVarSet = (getenv("GMX_EMULATE_GPU") != NULL);
1773 /* Run GPU emulation mode if GMX_EMULATE_GPU is defined. Because
1774 * GPUs (currently) only handle non-bonded calculations, we will
1775 * automatically switch to emulation if non-bonded calculations are
1776 * turned off via GMX_NO_NONBONDED - this is the simple and elegant
1777 * way to turn off GPU initialization, data movement, and cleanup.
1779 * GPU emulation can be useful to assess the performance one can expect by
1780 * adding GPU(s) to the machine. The conditional below allows this even
1781 * if mdrun is compiled without GPU acceleration support.
1782 * Note that you should freezing the system as otherwise it will explode.
1784 *bEmulateGPU = (bEmulateGPUEnvVarSet ||
1785 (!bDoNonbonded && gpu_opt->n_dev_use > 0));
1787 /* Enable GPU mode when GPUs are available or no GPU emulation is requested.
1789 if (gpu_opt->n_dev_use > 0 && !(*bEmulateGPU))
1791 /* Each PP node will use the intra-node id-th device from the
1792 * list of detected/selected GPUs. */
1793 if (!init_gpu(fp, cr->rank_pp_intranode, gpu_err_str,
1794 &hwinfo->gpu_info, gpu_opt))
1796 /* At this point the init should never fail as we made sure that
1797 * we have all the GPUs we need. If it still does, we'll bail. */
1798 /* TODO the decorating of gpu_err_str is nicer if it
1799 happens inside init_gpu. Out here, the decorating with
1800 the MPI rank makes sense. */
1801 gmx_fatal(FARGS, "On rank %d failed to initialize GPU #%d: %s",
1803 get_cuda_gpu_device_id(&hwinfo->gpu_info, gpu_opt,
1804 cr->rank_pp_intranode),
1808 /* Here we actually turn on hardware GPU acceleration */
1813 gmx_bool uses_simple_tables(int cutoff_scheme,
1814 nonbonded_verlet_t *nbv,
1817 gmx_bool bUsesSimpleTables = TRUE;
1820 switch (cutoff_scheme)
1823 bUsesSimpleTables = TRUE;
1826 assert(NULL != nbv && NULL != nbv->grp);
1827 grp_index = (group < 0) ? 0 : (nbv->ngrp - 1);
1828 bUsesSimpleTables = nbnxn_kernel_pairlist_simple(nbv->grp[grp_index].kernel_type);
1831 gmx_incons("unimplemented");
1833 return bUsesSimpleTables;
1836 static void init_ewald_f_table(interaction_const_t *ic,
1841 /* Get the Ewald table spacing based on Coulomb and/or LJ
1842 * Ewald coefficients and rtol.
1844 ic->tabq_scale = ewald_spline3_table_scale(ic);
1846 if (ic->cutoff_scheme == ecutsVERLET)
1848 maxr = ic->rcoulomb;
1852 maxr = std::max(ic->rcoulomb, rtab);
1854 ic->tabq_size = static_cast<int>(maxr*ic->tabq_scale) + 2;
1856 sfree_aligned(ic->tabq_coul_FDV0);
1857 sfree_aligned(ic->tabq_coul_F);
1858 sfree_aligned(ic->tabq_coul_V);
1860 sfree_aligned(ic->tabq_vdw_FDV0);
1861 sfree_aligned(ic->tabq_vdw_F);
1862 sfree_aligned(ic->tabq_vdw_V);
1864 if (ic->eeltype == eelEWALD || EEL_PME(ic->eeltype))
1866 /* Create the original table data in FDV0 */
1867 snew_aligned(ic->tabq_coul_FDV0, ic->tabq_size*4, 32);
1868 snew_aligned(ic->tabq_coul_F, ic->tabq_size, 32);
1869 snew_aligned(ic->tabq_coul_V, ic->tabq_size, 32);
1870 table_spline3_fill_ewald_lr(ic->tabq_coul_F, ic->tabq_coul_V, ic->tabq_coul_FDV0,
1871 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_q, v_q_ewald_lr);
1874 if (EVDW_PME(ic->vdwtype))
1876 snew_aligned(ic->tabq_vdw_FDV0, ic->tabq_size*4, 32);
1877 snew_aligned(ic->tabq_vdw_F, ic->tabq_size, 32);
1878 snew_aligned(ic->tabq_vdw_V, ic->tabq_size, 32);
1879 table_spline3_fill_ewald_lr(ic->tabq_vdw_F, ic->tabq_vdw_V, ic->tabq_vdw_FDV0,
1880 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_lj, v_lj_ewald_lr);
1884 void init_interaction_const_tables(FILE *fp,
1885 interaction_const_t *ic,
1888 if (ic->eeltype == eelEWALD || EEL_PME(ic->eeltype) || EVDW_PME(ic->vdwtype))
1890 init_ewald_f_table(ic, rtab);
1894 fprintf(fp, "Initialized non-bonded Ewald correction tables, spacing: %.2e size: %d\n\n",
1895 1/ic->tabq_scale, ic->tabq_size);
1900 static void clear_force_switch_constants(shift_consts_t *sc)
1907 static void force_switch_constants(real p,
1911 /* Here we determine the coefficient for shifting the force to zero
1912 * between distance rsw and the cut-off rc.
1913 * For a potential of r^-p, we have force p*r^-(p+1).
1914 * But to save flops we absorb p in the coefficient.
1916 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
1917 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
1919 sc->c2 = ((p + 1)*rsw - (p + 4)*rc)/(pow(rc, p + 2)*pow(rc - rsw, 2));
1920 sc->c3 = -((p + 1)*rsw - (p + 3)*rc)/(pow(rc, p + 2)*pow(rc - rsw, 3));
1921 sc->cpot = -pow(rc, -p) + p*sc->c2/3*pow(rc - rsw, 3) + p*sc->c3/4*pow(rc - rsw, 4);
1924 static void potential_switch_constants(real rsw, real rc,
1925 switch_consts_t *sc)
1927 /* The switch function is 1 at rsw and 0 at rc.
1928 * The derivative and second derivate are zero at both ends.
1929 * rsw = max(r - r_switch, 0)
1930 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
1931 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
1932 * force = force*dsw - potential*sw
1935 sc->c3 = -10*pow(rc - rsw, -3);
1936 sc->c4 = 15*pow(rc - rsw, -4);
1937 sc->c5 = -6*pow(rc - rsw, -5);
1940 /*! \brief Construct interaction constants
1942 * This data is used (particularly) by search and force code for
1943 * short-range interactions. Many of these are constant for the whole
1944 * simulation; some are constant only after PME tuning completes.
1947 init_interaction_const(FILE *fp,
1948 interaction_const_t **interaction_const,
1949 const t_forcerec *fr)
1951 interaction_const_t *ic;
1952 const real minusSix = -6.0;
1953 const real minusTwelve = -12.0;
1957 ic->cutoff_scheme = fr->cutoff_scheme;
1959 /* Just allocate something so we can free it */
1960 snew_aligned(ic->tabq_coul_FDV0, 16, 32);
1961 snew_aligned(ic->tabq_coul_F, 16, 32);
1962 snew_aligned(ic->tabq_coul_V, 16, 32);
1964 ic->rlist = fr->rlist;
1965 ic->rlistlong = fr->rlistlong;
1968 ic->vdwtype = fr->vdwtype;
1969 ic->vdw_modifier = fr->vdw_modifier;
1970 ic->rvdw = fr->rvdw;
1971 ic->rvdw_switch = fr->rvdw_switch;
1972 ic->ewaldcoeff_lj = fr->ewaldcoeff_lj;
1973 ic->ljpme_comb_rule = fr->ljpme_combination_rule;
1974 ic->sh_lj_ewald = 0;
1975 clear_force_switch_constants(&ic->dispersion_shift);
1976 clear_force_switch_constants(&ic->repulsion_shift);
1978 switch (ic->vdw_modifier)
1980 case eintmodPOTSHIFT:
1981 /* Only shift the potential, don't touch the force */
1982 ic->dispersion_shift.cpot = -pow(ic->rvdw, minusSix);
1983 ic->repulsion_shift.cpot = -pow(ic->rvdw, minusTwelve);
1984 if (EVDW_PME(ic->vdwtype))
1988 crc2 = sqr(ic->ewaldcoeff_lj*ic->rvdw);
1989 ic->sh_lj_ewald = (exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)*pow(ic->rvdw, minusSix);
1992 case eintmodFORCESWITCH:
1993 /* Switch the force, switch and shift the potential */
1994 force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw,
1995 &ic->dispersion_shift);
1996 force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw,
1997 &ic->repulsion_shift);
1999 case eintmodPOTSWITCH:
2000 /* Switch the potential and force */
2001 potential_switch_constants(ic->rvdw_switch, ic->rvdw,
2005 case eintmodEXACTCUTOFF:
2006 /* Nothing to do here */
2009 gmx_incons("unimplemented potential modifier");
2012 ic->sh_invrc6 = -ic->dispersion_shift.cpot;
2014 /* Electrostatics */
2015 ic->eeltype = fr->eeltype;
2016 ic->coulomb_modifier = fr->coulomb_modifier;
2017 ic->rcoulomb = fr->rcoulomb;
2018 ic->epsilon_r = fr->epsilon_r;
2019 ic->epsfac = fr->epsfac;
2020 ic->ewaldcoeff_q = fr->ewaldcoeff_q;
2022 if (fr->coulomb_modifier == eintmodPOTSHIFT)
2024 ic->sh_ewald = gmx_erfc(ic->ewaldcoeff_q*ic->rcoulomb);
2031 /* Reaction-field */
2032 if (EEL_RF(ic->eeltype))
2034 ic->epsilon_rf = fr->epsilon_rf;
2035 ic->k_rf = fr->k_rf;
2036 ic->c_rf = fr->c_rf;
2040 /* For plain cut-off we might use the reaction-field kernels */
2041 ic->epsilon_rf = ic->epsilon_r;
2043 if (fr->coulomb_modifier == eintmodPOTSHIFT)
2045 ic->c_rf = 1/ic->rcoulomb;
2055 real dispersion_shift;
2057 dispersion_shift = ic->dispersion_shift.cpot;
2058 if (EVDW_PME(ic->vdwtype))
2060 dispersion_shift -= ic->sh_lj_ewald;
2062 fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
2063 ic->repulsion_shift.cpot, dispersion_shift);
2065 if (ic->eeltype == eelCUT)
2067 fprintf(fp, ", Coulomb %.e", -ic->c_rf);
2069 else if (EEL_PME(ic->eeltype))
2071 fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
2076 *interaction_const = ic;
2079 /*! \brief Manage initialization within the NBNXN module of
2080 * run-time constants.
2083 initialize_gpu_constants(const t_commrec gmx_unused *cr,
2084 interaction_const_t *interaction_const,
2085 const struct nonbonded_verlet_t *nbv)
2087 if (nbv != NULL && nbv->bUseGPU)
2089 nbnxn_gpu_init_const(nbv->gpu_nbv, interaction_const, nbv->grp);
2091 /* With tMPI + GPUs some ranks may be sharing GPU(s) and therefore
2092 * also sharing texture references. To keep the code simple, we don't
2093 * treat texture references as shared resources, but this means that
2094 * the coulomb_tab and nbfp texture refs will get updated by multiple threads.
2095 * Hence, to ensure that the non-bonded kernels don't start before all
2096 * texture binding operations are finished, we need to wait for all ranks
2097 * to arrive here before continuing.
2099 * Note that we could omit this barrier if GPUs are not shared (or
2100 * texture objects are used), but as this is initialization code, there
2101 * is no point in complicating things.
2103 #ifdef GMX_THREAD_MPI
2108 #endif /* GMX_THREAD_MPI */
2113 static void init_nb_verlet(FILE *fp,
2114 nonbonded_verlet_t **nb_verlet,
2115 gmx_bool bFEP_NonBonded,
2116 const t_inputrec *ir,
2117 const t_forcerec *fr,
2118 const t_commrec *cr,
2119 const char *nbpu_opt)
2121 nonbonded_verlet_t *nbv;
2124 gmx_bool bEmulateGPU, bHybridGPURun = FALSE;
2126 nbnxn_alloc_t *nb_alloc;
2127 nbnxn_free_t *nb_free;
2131 pick_nbnxn_resources(fp, cr, fr->hwinfo,
2139 nbv->ngrp = (DOMAINDECOMP(cr) ? 2 : 1);
2140 for (i = 0; i < nbv->ngrp; i++)
2142 nbv->grp[i].nbl_lists.nnbl = 0;
2143 nbv->grp[i].nbat = NULL;
2144 nbv->grp[i].kernel_type = nbnxnkNotSet;
2146 if (i == 0) /* local */
2148 pick_nbnxn_kernel(fp, cr, fr->use_simd_kernels,
2149 nbv->bUseGPU, bEmulateGPU, ir,
2150 &nbv->grp[i].kernel_type,
2151 &nbv->grp[i].ewald_excl,
2154 else /* non-local */
2156 if (nbpu_opt != NULL && strcmp(nbpu_opt, "gpu_cpu") == 0)
2158 /* Use GPU for local, select a CPU kernel for non-local */
2159 pick_nbnxn_kernel(fp, cr, fr->use_simd_kernels,
2161 &nbv->grp[i].kernel_type,
2162 &nbv->grp[i].ewald_excl,
2165 bHybridGPURun = TRUE;
2169 /* Use the same kernel for local and non-local interactions */
2170 nbv->grp[i].kernel_type = nbv->grp[0].kernel_type;
2171 nbv->grp[i].ewald_excl = nbv->grp[0].ewald_excl;
2178 nbnxn_gpu_compile_kernels(cr->rank_pp_intranode, cr->nodeid, &fr->hwinfo->gpu_info, fr->gpu_opt, fr->ic);
2180 /* init the NxN GPU data; the last argument tells whether we'll have
2181 * both local and non-local NB calculation on GPU */
2182 nbnxn_gpu_init(fp, &nbv->gpu_nbv,
2183 &fr->hwinfo->gpu_info, fr->gpu_opt,
2184 cr->rank_pp_intranode,
2185 (nbv->ngrp > 1) && !bHybridGPURun);
2187 if ((env = getenv("GMX_NB_MIN_CI")) != NULL)
2191 nbv->min_ci_balanced = strtol(env, &end, 10);
2192 if (!end || (*end != 0) || nbv->min_ci_balanced <= 0)
2194 gmx_fatal(FARGS, "Invalid value passed in GMX_NB_MIN_CI=%s, positive integer required", env);
2199 fprintf(debug, "Neighbor-list balancing parameter: %d (passed as env. var.)\n",
2200 nbv->min_ci_balanced);
2205 nbv->min_ci_balanced = nbnxn_gpu_min_ci_balanced(nbv->gpu_nbv);
2208 fprintf(debug, "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU multi-processors)\n",
2209 nbv->min_ci_balanced);
2215 nbv->min_ci_balanced = 0;
2220 nbnxn_init_search(&nbv->nbs,
2221 DOMAINDECOMP(cr) ? &cr->dd->nc : NULL,
2222 DOMAINDECOMP(cr) ? domdec_zones(cr->dd) : NULL,
2224 gmx_omp_nthreads_get(emntPairsearch));
2226 for (i = 0; i < nbv->ngrp; i++)
2228 gpu_set_host_malloc_and_free(nbv->grp[0].kernel_type == nbnxnk8x8x8_GPU,
2229 &nb_alloc, &nb_free);
2231 nbnxn_init_pairlist_set(&nbv->grp[i].nbl_lists,
2232 nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2233 /* 8x8x8 "non-simple" lists are ATM always combined */
2234 !nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2238 nbv->grp[0].kernel_type != nbv->grp[i].kernel_type)
2240 gmx_bool bSimpleList;
2241 int enbnxninitcombrule;
2243 bSimpleList = nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type);
2245 if (bSimpleList && (fr->vdwtype == evdwCUT && (fr->vdw_modifier == eintmodNONE || fr->vdw_modifier == eintmodPOTSHIFT)))
2247 /* Plain LJ cut-off: we can optimize with combination rules */
2248 enbnxninitcombrule = enbnxninitcombruleDETECT;
2250 else if (fr->vdwtype == evdwPME)
2252 /* LJ-PME: we need to use a combination rule for the grid */
2253 if (fr->ljpme_combination_rule == eljpmeGEOM)
2255 enbnxninitcombrule = enbnxninitcombruleGEOM;
2259 enbnxninitcombrule = enbnxninitcombruleLB;
2264 /* We use a full combination matrix: no rule required */
2265 enbnxninitcombrule = enbnxninitcombruleNONE;
2269 snew(nbv->grp[i].nbat, 1);
2270 nbnxn_atomdata_init(fp,
2272 nbv->grp[i].kernel_type,
2274 fr->ntype, fr->nbfp,
2276 bSimpleList ? gmx_omp_nthreads_get(emntNonbonded) : 1,
2281 nbv->grp[i].nbat = nbv->grp[0].nbat;
2286 gmx_bool usingGpu(nonbonded_verlet_t *nbv)
2288 return nbv != NULL && nbv->bUseGPU;
2291 void init_forcerec(FILE *fp,
2292 const output_env_t oenv,
2295 const t_inputrec *ir,
2296 const gmx_mtop_t *mtop,
2297 const t_commrec *cr,
2303 const char *nbpu_opt,
2304 gmx_bool bNoSolvOpt,
2307 int i, m, negp_pp, negptable, egi, egj;
2312 gmx_bool bGenericKernelOnly;
2313 gmx_bool bMakeTables, bMakeSeparate14Table, bSomeNormalNbListsAreInUse;
2314 gmx_bool bFEP_NonBonded;
2315 int *nm_ind, egp_flags;
2317 if (fr->hwinfo == NULL)
2319 /* Detect hardware, gather information.
2320 * In mdrun, hwinfo has already been set before calling init_forcerec.
2321 * Here we ignore GPUs, as tools will not use them anyhow.
2323 fr->hwinfo = gmx_detect_hardware(fp, cr, FALSE);
2326 /* By default we turn SIMD kernels on, but it might be turned off further down... */
2327 fr->use_simd_kernels = TRUE;
2329 fr->bDomDec = DOMAINDECOMP(cr);
2331 if (check_box(ir->ePBC, box))
2333 gmx_fatal(FARGS, check_box(ir->ePBC, box));
2336 /* Test particle insertion ? */
2339 /* Set to the size of the molecule to be inserted (the last one) */
2340 /* Because of old style topologies, we have to use the last cg
2341 * instead of the last molecule type.
2343 cgs = &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].cgs;
2344 fr->n_tpi = cgs->index[cgs->nr] - cgs->index[cgs->nr-1];
2345 if (fr->n_tpi != mtop->mols.index[mtop->mols.nr] - mtop->mols.index[mtop->mols.nr-1])
2347 gmx_fatal(FARGS, "The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
2355 /* Copy AdResS parameters */
2358 fr->adress_type = ir->adress->type;
2359 fr->adress_const_wf = ir->adress->const_wf;
2360 fr->adress_ex_width = ir->adress->ex_width;
2361 fr->adress_hy_width = ir->adress->hy_width;
2362 fr->adress_icor = ir->adress->icor;
2363 fr->adress_site = ir->adress->site;
2364 fr->adress_ex_forcecap = ir->adress->ex_forcecap;
2365 fr->adress_do_hybridpairs = ir->adress->do_hybridpairs;
2368 snew(fr->adress_group_explicit, ir->adress->n_energy_grps);
2369 for (i = 0; i < ir->adress->n_energy_grps; i++)
2371 fr->adress_group_explicit[i] = ir->adress->group_explicit[i];
2374 fr->n_adress_tf_grps = ir->adress->n_tf_grps;
2375 snew(fr->adress_tf_table_index, fr->n_adress_tf_grps);
2376 for (i = 0; i < fr->n_adress_tf_grps; i++)
2378 fr->adress_tf_table_index[i] = ir->adress->tf_table_index[i];
2380 copy_rvec(ir->adress->refs, fr->adress_refs);
2384 fr->adress_type = eAdressOff;
2385 fr->adress_do_hybridpairs = FALSE;
2388 /* Copy the user determined parameters */
2389 fr->userint1 = ir->userint1;
2390 fr->userint2 = ir->userint2;
2391 fr->userint3 = ir->userint3;
2392 fr->userint4 = ir->userint4;
2393 fr->userreal1 = ir->userreal1;
2394 fr->userreal2 = ir->userreal2;
2395 fr->userreal3 = ir->userreal3;
2396 fr->userreal4 = ir->userreal4;
2399 fr->fc_stepsize = ir->fc_stepsize;
2402 fr->efep = ir->efep;
2403 fr->sc_alphavdw = ir->fepvals->sc_alpha;
2404 if (ir->fepvals->bScCoul)
2406 fr->sc_alphacoul = ir->fepvals->sc_alpha;
2407 fr->sc_sigma6_min = pow(ir->fepvals->sc_sigma_min, 6);
2411 fr->sc_alphacoul = 0;
2412 fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
2414 fr->sc_power = ir->fepvals->sc_power;
2415 fr->sc_r_power = ir->fepvals->sc_r_power;
2416 fr->sc_sigma6_def = pow(ir->fepvals->sc_sigma, 6);
2418 env = getenv("GMX_SCSIGMA_MIN");
2422 sscanf(env, "%20lf", &dbl);
2423 fr->sc_sigma6_min = pow(dbl, 6);
2426 fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
2430 fr->bNonbonded = TRUE;
2431 if (getenv("GMX_NO_NONBONDED") != NULL)
2433 /* turn off non-bonded calculations */
2434 fr->bNonbonded = FALSE;
2435 md_print_warn(cr, fp,
2436 "Found environment variable GMX_NO_NONBONDED.\n"
2437 "Disabling nonbonded calculations.\n");
2440 bGenericKernelOnly = FALSE;
2442 /* We now check in the NS code whether a particular combination of interactions
2443 * can be used with water optimization, and disable it if that is not the case.
2446 if (getenv("GMX_NB_GENERIC") != NULL)
2451 "Found environment variable GMX_NB_GENERIC.\n"
2452 "Disabling all interaction-specific nonbonded kernels, will only\n"
2453 "use the slow generic ones in src/gmxlib/nonbonded/nb_generic.c\n\n");
2455 bGenericKernelOnly = TRUE;
2458 if (bGenericKernelOnly == TRUE)
2463 if ( (getenv("GMX_DISABLE_SIMD_KERNELS") != NULL) || (getenv("GMX_NOOPTIMIZEDKERNELS") != NULL) )
2465 fr->use_simd_kernels = FALSE;
2469 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
2470 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
2471 "(e.g. SSE2/SSE4.1/AVX).\n\n");
2475 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
2477 /* Check if we can/should do all-vs-all kernels */
2478 fr->bAllvsAll = can_use_allvsall(ir, FALSE, NULL, NULL);
2479 fr->AllvsAll_work = NULL;
2480 fr->AllvsAll_workgb = NULL;
2482 /* All-vs-all kernels have not been implemented in 4.6 and later.
2483 * See Redmine #1249. */
2486 fr->bAllvsAll = FALSE;
2490 "\nYour simulation settings would have triggered the efficient all-vs-all\n"
2491 "kernels in GROMACS 4.5, but these have not been implemented in GROMACS\n"
2492 "4.6 and 5.x. If performance is important, please use GROMACS 4.5.7\n"
2493 "or try cutoff-scheme = Verlet.\n\n");
2497 /* Neighbour searching stuff */
2498 fr->cutoff_scheme = ir->cutoff_scheme;
2499 fr->bGrid = (ir->ns_type == ensGRID);
2500 fr->ePBC = ir->ePBC;
2502 if (fr->cutoff_scheme == ecutsGROUP)
2504 const char *note = "NOTE: This file uses the deprecated 'group' cutoff_scheme. This will be\n"
2505 "removed in a future release when 'verlet' supports all interaction forms.\n";
2509 fprintf(stderr, "\n%s\n", note);
2513 fprintf(fp, "\n%s\n", note);
2517 /* Determine if we will do PBC for distances in bonded interactions */
2518 if (fr->ePBC == epbcNONE)
2520 fr->bMolPBC = FALSE;
2524 if (!DOMAINDECOMP(cr))
2528 bSHAKE = (ir->eConstrAlg == econtSHAKE &&
2529 (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0 ||
2530 gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0));
2532 /* The group cut-off scheme and SHAKE assume charge groups
2533 * are whole, but not using molpbc is faster in most cases.
2534 * With intermolecular interactions we need PBC for calculating
2535 * distances between atoms in different molecules.
2537 if ((fr->cutoff_scheme == ecutsGROUP || bSHAKE) &&
2538 !mtop->bIntermolecularInteractions)
2540 fr->bMolPBC = ir->bPeriodicMols;
2542 if (bSHAKE && fr->bMolPBC)
2544 gmx_fatal(FARGS, "SHAKE is not supported with periodic molecules");
2551 if (getenv("GMX_USE_GRAPH") != NULL)
2553 fr->bMolPBC = FALSE;
2556 md_print_warn(cr, fp, "GMX_USE_GRAPH is set, using the graph for bonded interactions\n");
2559 if (mtop->bIntermolecularInteractions)
2561 md_print_warn(cr, fp, "WARNING: Molecules linked by intermolecular interactions have to reside in the same periodic image, otherwise artifacts will occur!\n");
2565 if (bSHAKE && fr->bMolPBC)
2567 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");
2573 fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->ePBC);
2576 fr->bGB = (ir->implicit_solvent == eisGBSA);
2578 fr->rc_scaling = ir->refcoord_scaling;
2579 copy_rvec(ir->posres_com, fr->posres_com);
2580 copy_rvec(ir->posres_comB, fr->posres_comB);
2581 fr->rlist = cutoff_inf(ir->rlist);
2582 fr->rlistlong = cutoff_inf(ir->rlistlong);
2583 fr->eeltype = ir->coulombtype;
2584 fr->vdwtype = ir->vdwtype;
2585 fr->ljpme_combination_rule = ir->ljpme_combination_rule;
2587 fr->coulomb_modifier = ir->coulomb_modifier;
2588 fr->vdw_modifier = ir->vdw_modifier;
2590 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
2591 switch (fr->eeltype)
2594 fr->nbkernel_elec_interaction = (fr->bGB) ? GMX_NBKERNEL_ELEC_GENERALIZEDBORN : GMX_NBKERNEL_ELEC_COULOMB;
2600 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2604 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2605 fr->coulomb_modifier = eintmodEXACTCUTOFF;
2614 case eelPMEUSERSWITCH:
2615 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2620 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD;
2624 gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[fr->eeltype]);
2628 /* Vdw: Translate from mdp settings to kernel format */
2629 switch (fr->vdwtype)
2634 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
2638 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
2642 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD;
2648 case evdwENCADSHIFT:
2649 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2653 gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[fr->vdwtype]);
2657 /* These start out identical to ir, but might be altered if we e.g. tabulate the interaction in the kernel */
2658 fr->nbkernel_elec_modifier = fr->coulomb_modifier;
2659 fr->nbkernel_vdw_modifier = fr->vdw_modifier;
2661 fr->rvdw = cutoff_inf(ir->rvdw);
2662 fr->rvdw_switch = ir->rvdw_switch;
2663 fr->rcoulomb = cutoff_inf(ir->rcoulomb);
2664 fr->rcoulomb_switch = ir->rcoulomb_switch;
2666 fr->bTwinRange = fr->rlistlong > fr->rlist;
2667 fr->bEwald = (EEL_PME(fr->eeltype) || fr->eeltype == eelEWALD);
2669 fr->reppow = mtop->ffparams.reppow;
2671 if (ir->cutoff_scheme == ecutsGROUP)
2673 fr->bvdwtab = ((fr->vdwtype != evdwCUT || !gmx_within_tol(fr->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2674 && !EVDW_PME(fr->vdwtype));
2675 /* We have special kernels for standard Ewald and PME, but the pme-switch ones are tabulated above */
2676 fr->bcoultab = !(fr->eeltype == eelCUT ||
2677 fr->eeltype == eelEWALD ||
2678 fr->eeltype == eelPME ||
2679 fr->eeltype == eelRF ||
2680 fr->eeltype == eelRF_ZERO);
2682 /* If the user absolutely wants different switch/shift settings for coul/vdw, it is likely
2683 * going to be faster to tabulate the interaction than calling the generic kernel.
2684 * However, if generic kernels have been requested we keep things analytically.
2686 if (fr->nbkernel_elec_modifier == eintmodPOTSWITCH &&
2687 fr->nbkernel_vdw_modifier == eintmodPOTSWITCH &&
2688 bGenericKernelOnly == FALSE)
2690 if ((fr->rcoulomb_switch != fr->rvdw_switch) || (fr->rcoulomb != fr->rvdw))
2692 fr->bcoultab = TRUE;
2693 /* Once we tabulate electrostatics, we can use the switch function for LJ,
2694 * which would otherwise need two tables.
2698 else if ((fr->nbkernel_elec_modifier == eintmodPOTSHIFT && fr->nbkernel_vdw_modifier == eintmodPOTSHIFT) ||
2699 ((fr->nbkernel_elec_interaction == GMX_NBKERNEL_ELEC_REACTIONFIELD &&
2700 fr->nbkernel_elec_modifier == eintmodEXACTCUTOFF &&
2701 (fr->nbkernel_vdw_modifier == eintmodPOTSWITCH || fr->nbkernel_vdw_modifier == eintmodPOTSHIFT))))
2703 if ((fr->rcoulomb != fr->rvdw) && (bGenericKernelOnly == FALSE))
2705 fr->bcoultab = TRUE;
2709 if (fr->nbkernel_elec_modifier == eintmodFORCESWITCH)
2711 fr->bcoultab = TRUE;
2713 if (fr->nbkernel_vdw_modifier == eintmodFORCESWITCH)
2718 if (getenv("GMX_REQUIRE_TABLES"))
2721 fr->bcoultab = TRUE;
2726 fprintf(fp, "Table routines are used for coulomb: %s\n", bool_names[fr->bcoultab]);
2727 fprintf(fp, "Table routines are used for vdw: %s\n", bool_names[fr->bvdwtab ]);
2730 if (fr->bvdwtab == TRUE)
2732 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2733 fr->nbkernel_vdw_modifier = eintmodNONE;
2735 if (fr->bcoultab == TRUE)
2737 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2738 fr->nbkernel_elec_modifier = eintmodNONE;
2742 if (ir->cutoff_scheme == ecutsVERLET)
2744 if (!gmx_within_tol(fr->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2746 gmx_fatal(FARGS, "Cut-off scheme %S only supports LJ repulsion power 12", ecutscheme_names[ir->cutoff_scheme]);
2748 fr->bvdwtab = FALSE;
2749 fr->bcoultab = FALSE;
2752 /* Tables are used for direct ewald sum */
2755 if (EEL_PME(ir->coulombtype))
2759 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
2761 if (ir->coulombtype == eelP3M_AD)
2763 please_cite(fp, "Hockney1988");
2764 please_cite(fp, "Ballenegger2012");
2768 please_cite(fp, "Essmann95a");
2771 if (ir->ewald_geometry == eewg3DC)
2775 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry.\n");
2777 please_cite(fp, "In-Chul99a");
2780 fr->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
2781 init_ewald_tab(&(fr->ewald_table), ir, fp);
2784 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n",
2785 1/fr->ewaldcoeff_q);
2789 if (EVDW_PME(ir->vdwtype))
2793 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
2795 please_cite(fp, "Essmann95a");
2796 fr->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
2799 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n",
2800 1/fr->ewaldcoeff_lj);
2804 /* Electrostatics */
2805 fr->epsilon_r = ir->epsilon_r;
2806 fr->epsilon_rf = ir->epsilon_rf;
2807 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
2809 /* Parameters for generalized RF */
2813 if (fr->eeltype == eelGRF)
2815 init_generalized_rf(fp, mtop, ir, fr);
2818 fr->bF_NoVirSum = (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype) ||
2819 gmx_mtop_ftype_count(mtop, F_POSRES) > 0 ||
2820 gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 ||
2821 IR_ELEC_FIELD(*ir) ||
2822 (fr->adress_icor != eAdressICOff)
2825 if (fr->cutoff_scheme == ecutsGROUP &&
2826 ncg_mtop(mtop) > fr->cg_nalloc && !DOMAINDECOMP(cr))
2828 /* Count the total number of charge groups */
2829 fr->cg_nalloc = ncg_mtop(mtop);
2830 srenew(fr->cg_cm, fr->cg_nalloc);
2832 if (fr->shift_vec == NULL)
2834 snew(fr->shift_vec, SHIFTS);
2837 if (fr->fshift == NULL)
2839 snew(fr->fshift, SHIFTS);
2842 if (fr->nbfp == NULL)
2844 fr->ntype = mtop->ffparams.atnr;
2845 fr->nbfp = mk_nbfp(&mtop->ffparams, fr->bBHAM);
2846 if (EVDW_PME(fr->vdwtype))
2848 fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
2852 /* Copy the energy group exclusions */
2853 fr->egp_flags = ir->opts.egp_flags;
2855 /* Van der Waals stuff */
2856 if ((fr->vdwtype != evdwCUT) && (fr->vdwtype != evdwUSER) && !fr->bBHAM)
2858 if (fr->rvdw_switch >= fr->rvdw)
2860 gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)",
2861 fr->rvdw_switch, fr->rvdw);
2865 fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
2866 (fr->eeltype == eelSWITCH) ? "switched" : "shifted",
2867 fr->rvdw_switch, fr->rvdw);
2871 if (fr->bBHAM && EVDW_PME(fr->vdwtype))
2873 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
2876 if (fr->bBHAM && (fr->vdwtype == evdwSHIFT || fr->vdwtype == evdwSWITCH))
2878 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
2881 if (fr->bBHAM && fr->cutoff_scheme == ecutsVERLET)
2883 gmx_fatal(FARGS, "Verlet cutoff-scheme is not supported with Buckingham");
2888 fprintf(fp, "Cut-off's: NS: %g Coulomb: %g %s: %g\n",
2889 fr->rlist, fr->rcoulomb, fr->bBHAM ? "BHAM" : "LJ", fr->rvdw);
2892 fr->eDispCorr = ir->eDispCorr;
2893 if (ir->eDispCorr != edispcNO)
2895 set_avcsixtwelve(fp, fr, mtop);
2900 set_bham_b_max(fp, fr, mtop);
2903 fr->gb_epsilon_solvent = ir->gb_epsilon_solvent;
2905 /* Copy the GBSA data (radius, volume and surftens for each
2906 * atomtype) from the topology atomtype section to forcerec.
2908 snew(fr->atype_radius, fr->ntype);
2909 snew(fr->atype_vol, fr->ntype);
2910 snew(fr->atype_surftens, fr->ntype);
2911 snew(fr->atype_gb_radius, fr->ntype);
2912 snew(fr->atype_S_hct, fr->ntype);
2914 if (mtop->atomtypes.nr > 0)
2916 for (i = 0; i < fr->ntype; i++)
2918 fr->atype_radius[i] = mtop->atomtypes.radius[i];
2920 for (i = 0; i < fr->ntype; i++)
2922 fr->atype_vol[i] = mtop->atomtypes.vol[i];
2924 for (i = 0; i < fr->ntype; i++)
2926 fr->atype_surftens[i] = mtop->atomtypes.surftens[i];
2928 for (i = 0; i < fr->ntype; i++)
2930 fr->atype_gb_radius[i] = mtop->atomtypes.gb_radius[i];
2932 for (i = 0; i < fr->ntype; i++)
2934 fr->atype_S_hct[i] = mtop->atomtypes.S_hct[i];
2938 /* Generate the GB table if needed */
2942 fr->gbtabscale = 2000;
2944 fr->gbtabscale = 500;
2948 fr->gbtab = make_gb_table(oenv, fr);
2950 init_gb(&fr->born, fr, ir, mtop, ir->gb_algorithm);
2952 /* Copy local gb data (for dd, this is done in dd_partition_system) */
2953 if (!DOMAINDECOMP(cr))
2955 make_local_gb(cr, fr->born, ir->gb_algorithm);
2959 /* Set the charge scaling */
2960 if (fr->epsilon_r != 0)
2962 fr->epsfac = ONE_4PI_EPS0/fr->epsilon_r;
2966 /* eps = 0 is infinite dieletric: no coulomb interactions */
2970 /* Reaction field constants */
2971 if (EEL_RF(fr->eeltype))
2973 calc_rffac(fp, fr->eeltype, fr->epsilon_r, fr->epsilon_rf,
2974 fr->rcoulomb, fr->temp, fr->zsquare, box,
2975 &fr->kappa, &fr->k_rf, &fr->c_rf);
2978 /*This now calculates sum for q and c6*/
2979 set_chargesum(fp, fr, mtop);
2981 /* if we are using LR electrostatics, and they are tabulated,
2982 * the tables will contain modified coulomb interactions.
2983 * Since we want to use the non-shifted ones for 1-4
2984 * coulombic interactions, we must have an extra set of tables.
2987 /* Construct tables.
2988 * A little unnecessary to make both vdw and coul tables sometimes,
2989 * but what the heck... */
2991 bMakeTables = fr->bcoultab || fr->bvdwtab || fr->bEwald ||
2992 (ir->eDispCorr != edispcNO && ir_vdw_switched(ir));
2994 bMakeSeparate14Table = ((!bMakeTables || fr->eeltype != eelCUT || fr->vdwtype != evdwCUT ||
2995 fr->coulomb_modifier != eintmodNONE ||
2996 fr->vdw_modifier != eintmodNONE ||
2997 fr->bBHAM || fr->bEwald) &&
2998 (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
2999 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0 ||
3000 gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0));
3002 negp_pp = ir->opts.ngener - ir->nwall;
3006 bSomeNormalNbListsAreInUse = TRUE;
3011 bSomeNormalNbListsAreInUse = (ir->eDispCorr != edispcNO);
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);
3049 snew(fr->nblists, fr->nnblists);
3051 /* This code automatically gives table length tabext without cut-off's,
3052 * in that case grompp should already have checked that we do not need
3053 * normal tables and we only generate tables for 1-4 interactions.
3055 rtab = ir->rlistlong + ir->tabext;
3059 /* make tables for ordinary interactions */
3060 if (bSomeNormalNbListsAreInUse)
3062 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn, NULL, NULL, &fr->nblists[0]);
3065 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn, NULL, NULL, &fr->nblists[fr->nnblists/2]);
3067 if (!bMakeSeparate14Table)
3069 fr->tab14 = fr->nblists[0].table_elec_vdw;
3079 /* Read the special tables for certain energy group pairs */
3080 nm_ind = mtop->groups.grps[egcENER].nm_ind;
3081 for (egi = 0; egi < negp_pp; egi++)
3083 for (egj = egi; egj < negp_pp; egj++)
3085 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
3086 if ((egp_flags & EGP_TABLE) && !(egp_flags & EGP_EXCL))
3088 if (fr->nnblists > 1)
3090 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = m;
3092 /* Read the table file with the two energy groups names appended */
3093 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn,
3094 *mtop->groups.grpname[nm_ind[egi]],
3095 *mtop->groups.grpname[nm_ind[egj]],
3099 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn,
3100 *mtop->groups.grpname[nm_ind[egi]],
3101 *mtop->groups.grpname[nm_ind[egj]],
3102 &fr->nblists[fr->nnblists/2+m]);
3106 else if (fr->nnblists > 1)
3108 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = 0;
3114 else if ((fr->eDispCorr != edispcNO) &&
3115 ((fr->vdw_modifier == eintmodPOTSWITCH) ||
3116 (fr->vdw_modifier == eintmodFORCESWITCH) ||
3117 (fr->vdw_modifier == eintmodPOTSHIFT)))
3119 /* Tables might not be used for the potential modifier interactions per se, but
3120 * we still need them to evaluate switch/shift dispersion corrections in this case.
3122 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn, NULL, NULL, &fr->nblists[0]);
3125 if (bMakeSeparate14Table)
3127 /* generate extra tables with plain Coulomb for 1-4 interactions only */
3128 fr->tab14 = make_tables(fp, oenv, fr, MASTER(cr), tabpfn, rtab,
3129 GMX_MAKETABLES_14ONLY);
3132 /* Read AdResS Thermo Force table if needed */
3133 if (fr->adress_icor == eAdressICThermoForce)
3135 /* old todo replace */
3137 if (ir->adress->n_tf_grps > 0)
3139 make_adress_tf_tables(fp, oenv, fr, ir, tabfn, mtop, box);
3144 /* load the default table */
3145 snew(fr->atf_tabs, 1);
3146 fr->atf_tabs[DEFAULT_TF_TABLE] = make_atf_table(fp, oenv, fr, tabafn, box);
3151 fr->nwall = ir->nwall;
3152 if (ir->nwall && ir->wall_type == ewtTABLE)
3154 make_wall_tables(fp, oenv, ir, tabfn, &mtop->groups, fr);
3159 fcd->bondtab = make_bonded_tables(fp,
3160 F_TABBONDS, F_TABBONDSNC,
3162 fcd->angletab = make_bonded_tables(fp,
3165 fcd->dihtab = make_bonded_tables(fp,
3173 fprintf(debug, "No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
3177 /* QM/MM initialization if requested
3181 fprintf(stderr, "QM/MM calculation requested.\n");
3184 fr->bQMMM = ir->bQMMM;
3185 fr->qr = mk_QMMMrec();
3187 /* Set all the static charge group info */
3188 fr->cginfo_mb = init_cginfo_mb(fp, mtop, fr, bNoSolvOpt,
3190 &fr->bExcl_IntraCGAll_InterCGNone);
3191 if (DOMAINDECOMP(cr))
3197 fr->cginfo = cginfo_expand(mtop->nmolblock, fr->cginfo_mb);
3200 if (!DOMAINDECOMP(cr))
3202 forcerec_set_ranges(fr, ncg_mtop(mtop), ncg_mtop(mtop),
3203 mtop->natoms, mtop->natoms, mtop->natoms);
3206 fr->print_force = print_force;
3209 /* coarse load balancing vars */
3214 /* Initialize neighbor search */
3215 init_ns(fp, cr, &fr->ns, fr, mtop);
3217 if (cr->duty & DUTY_PP)
3219 gmx_nonbonded_setup(fr, bGenericKernelOnly);
3223 gmx_setup_adress_kernels(fp,bGenericKernelOnly);
3228 /* Initialize the thread working data for bonded interactions */
3229 init_bonded_threading(fp, fr, mtop->groups.grps[egcENER].nr);
3231 snew(fr->excl_load, fr->nthreads+1);
3233 /* fr->ic is used both by verlet and group kernels (to some extent) now */
3234 init_interaction_const(fp, &fr->ic, fr);
3236 if (fr->cutoff_scheme == ecutsVERLET)
3238 if (ir->rcoulomb != ir->rvdw)
3240 gmx_fatal(FARGS, "With Verlet lists rcoulomb and rvdw should be identical");
3243 init_nb_verlet(fp, &fr->nbv, bFEP_NonBonded, ir, fr, cr, nbpu_opt);
3246 init_interaction_const_tables(fp, fr->ic, rtab);
3248 initialize_gpu_constants(cr, fr->ic, fr->nbv);
3250 if (ir->eDispCorr != edispcNO)
3252 calc_enervirdiff(fp, ir->eDispCorr, fr);
3256 #define pr_real(fp, r) fprintf(fp, "%s: %e\n",#r, r)
3257 #define pr_int(fp, i) fprintf((fp), "%s: %d\n",#i, i)
3258 #define pr_bool(fp, b) fprintf((fp), "%s: %s\n",#b, bool_names[b])
3260 void pr_forcerec(FILE *fp, t_forcerec *fr)
3264 pr_real(fp, fr->rlist);
3265 pr_real(fp, fr->rcoulomb);
3266 pr_real(fp, fr->fudgeQQ);
3267 pr_bool(fp, fr->bGrid);
3268 pr_bool(fp, fr->bTwinRange);
3269 /*pr_int(fp,fr->cg0);
3270 pr_int(fp,fr->hcg);*/
3271 for (i = 0; i < fr->nnblists; i++)
3273 pr_int(fp, fr->nblists[i].table_elec_vdw.n);
3275 pr_real(fp, fr->rcoulomb_switch);
3276 pr_real(fp, fr->rcoulomb);
3281 void forcerec_set_excl_load(t_forcerec *fr,
3282 const gmx_localtop_t *top)
3285 int t, i, j, ntot, n, ntarget;
3287 ind = top->excls.index;
3291 for (i = 0; i < top->excls.nr; i++)
3293 for (j = ind[i]; j < ind[i+1]; j++)
3302 fr->excl_load[0] = 0;
3305 for (t = 1; t <= fr->nthreads; t++)
3307 ntarget = (ntot*t)/fr->nthreads;
3308 while (i < top->excls.nr && n < ntarget)
3310 for (j = ind[i]; j < ind[i+1]; j++)
3319 fr->excl_load[t] = i;
3323 /* Frees GPU memory and destroys the GPU context.
3325 * Note that this function needs to be called even if GPUs are not used
3326 * in this run because the PME ranks have no knowledge of whether GPUs
3327 * are used or not, but all ranks need to enter the barrier below.
3329 void free_gpu_resources(const t_forcerec *fr,
3330 const t_commrec *cr,
3331 const gmx_gpu_info_t *gpu_info,
3332 const gmx_gpu_opt_t *gpu_opt)
3334 gmx_bool bIsPPrankUsingGPU;
3335 char gpu_err_str[STRLEN];
3337 bIsPPrankUsingGPU = (cr->duty & DUTY_PP) && fr && fr->nbv && fr->nbv->bUseGPU;
3339 if (bIsPPrankUsingGPU)
3341 /* free nbnxn data in GPU memory */
3342 nbnxn_gpu_free(fr->nbv->gpu_nbv);
3344 /* With tMPI we need to wait for all ranks to finish deallocation before
3345 * destroying the context in free_gpu() as some ranks may be sharing
3347 * Note: as only PP ranks need to free GPU resources, so it is safe to
3348 * not call the barrier on PME ranks.
3350 #ifdef GMX_THREAD_MPI
3355 #endif /* GMX_THREAD_MPI */
3357 /* uninitialize GPU (by destroying the context) */
3358 if (!free_cuda_gpu(cr->rank_pp_intranode, gpu_err_str, gpu_info, gpu_opt))
3360 gmx_warning("On rank %d failed to free GPU #%d: %s",
3361 cr->nodeid, get_current_cuda_gpu_device_id(), gpu_err_str);