2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016, 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/fileio/filenm.h"
51 #include "gromacs/gmxlib/gpu_utils/gpu_utils.h"
52 #include "gromacs/legacyheaders/copyrite.h"
53 #include "gromacs/legacyheaders/force.h"
54 #include "gromacs/legacyheaders/gmx_detect_hardware.h"
55 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
56 #include "gromacs/legacyheaders/inputrec.h"
57 #include "gromacs/legacyheaders/macros.h"
58 #include "gromacs/legacyheaders/md_logging.h"
59 #include "gromacs/legacyheaders/md_support.h"
60 #include "gromacs/legacyheaders/names.h"
61 #include "gromacs/legacyheaders/network.h"
62 #include "gromacs/legacyheaders/nonbonded.h"
63 #include "gromacs/legacyheaders/ns.h"
64 #include "gromacs/legacyheaders/qmmm.h"
65 #include "gromacs/legacyheaders/tables.h"
66 #include "gromacs/legacyheaders/txtdump.h"
67 #include "gromacs/legacyheaders/typedefs.h"
68 #include "gromacs/legacyheaders/types/commrec.h"
69 #include "gromacs/listed-forces/manage-threading.h"
70 #include "gromacs/math/calculate-ewald-splitting-coefficient.h"
71 #include "gromacs/math/units.h"
72 #include "gromacs/math/utilities.h"
73 #include "gromacs/math/vec.h"
74 #include "gromacs/mdlib/forcerec-threading.h"
75 #include "gromacs/mdlib/nb_verlet.h"
76 #include "gromacs/mdlib/nbnxn_atomdata.h"
77 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
78 #include "gromacs/mdlib/nbnxn_search.h"
79 #include "gromacs/mdlib/nbnxn_simd.h"
80 #include "gromacs/pbcutil/ishift.h"
81 #include "gromacs/pbcutil/pbc.h"
82 #include "gromacs/simd/simd.h"
83 #include "gromacs/topology/mtop_util.h"
84 #include "gromacs/utility/exceptions.h"
85 #include "gromacs/utility/fatalerror.h"
86 #include "gromacs/utility/smalloc.h"
87 #include "gromacs/utility/stringutil.h"
89 #include "nbnxn_gpu_jit_support.h"
91 t_forcerec *mk_forcerec(void)
101 static void pr_nbfp(FILE *fp, real *nbfp, gmx_bool bBHAM, int atnr)
105 for (i = 0; (i < atnr); i++)
107 for (j = 0; (j < atnr); j++)
109 fprintf(fp, "%2d - %2d", i, j);
112 fprintf(fp, " a=%10g, b=%10g, c=%10g\n", BHAMA(nbfp, atnr, i, j),
113 BHAMB(nbfp, atnr, i, j), BHAMC(nbfp, atnr, i, j)/6.0);
117 fprintf(fp, " c6=%10g, c12=%10g\n", C6(nbfp, atnr, i, j)/6.0,
118 C12(nbfp, atnr, i, j)/12.0);
125 static real *mk_nbfp(const gmx_ffparams_t *idef, gmx_bool bBHAM)
133 snew(nbfp, 3*atnr*atnr);
134 for (i = k = 0; (i < atnr); i++)
136 for (j = 0; (j < atnr); j++, k++)
138 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
139 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
140 /* nbfp now includes the 6.0 derivative prefactor */
141 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c*6.0;
147 snew(nbfp, 2*atnr*atnr);
148 for (i = k = 0; (i < atnr); i++)
150 for (j = 0; (j < atnr); j++, k++)
152 /* nbfp now includes the 6.0/12.0 derivative prefactors */
153 C6(nbfp, atnr, i, j) = idef->iparams[k].lj.c6*6.0;
154 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12*12.0;
162 static real *make_ljpme_c6grid(const gmx_ffparams_t *idef, t_forcerec *fr)
165 real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
167 const real oneOverSix = 1.0 / 6.0;
169 /* For LJ-PME simulations, we correct the energies with the reciprocal space
170 * inside of the cut-off. To do this the non-bonded kernels needs to have
171 * access to the C6-values used on the reciprocal grid in pme.c
175 snew(grid, 2*atnr*atnr);
176 for (i = k = 0; (i < atnr); i++)
178 for (j = 0; (j < atnr); j++, k++)
180 c6i = idef->iparams[i*(atnr+1)].lj.c6;
181 c12i = idef->iparams[i*(atnr+1)].lj.c12;
182 c6j = idef->iparams[j*(atnr+1)].lj.c6;
183 c12j = idef->iparams[j*(atnr+1)].lj.c12;
184 c6 = sqrt(c6i * c6j);
185 if (fr->ljpme_combination_rule == eljpmeLB
186 && !gmx_numzero(c6) && !gmx_numzero(c12i) && !gmx_numzero(c12j))
188 sigmai = pow(c12i / c6i, oneOverSix);
189 sigmaj = pow(c12j / c6j, oneOverSix);
190 epsi = c6i * c6i / c12i;
191 epsj = c6j * c6j / c12j;
192 c6 = sqrt(epsi * epsj) * pow(0.5*(sigmai+sigmaj), 6);
194 /* Store the elements at the same relative positions as C6 in nbfp in order
195 * to simplify access in the kernels
197 grid[2*(atnr*i+j)] = c6*6.0;
203 static real *mk_nbfp_combination_rule(const gmx_ffparams_t *idef, int comb_rule)
207 real c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
209 const real oneOverSix = 1.0 / 6.0;
212 snew(nbfp, 2*atnr*atnr);
213 for (i = 0; i < atnr; ++i)
215 for (j = 0; j < atnr; ++j)
217 c6i = idef->iparams[i*(atnr+1)].lj.c6;
218 c12i = idef->iparams[i*(atnr+1)].lj.c12;
219 c6j = idef->iparams[j*(atnr+1)].lj.c6;
220 c12j = idef->iparams[j*(atnr+1)].lj.c12;
221 c6 = sqrt(c6i * c6j);
222 c12 = sqrt(c12i * c12j);
223 if (comb_rule == eCOMB_ARITHMETIC
224 && !gmx_numzero(c6) && !gmx_numzero(c12))
226 sigmai = pow(c12i / c6i, oneOverSix);
227 sigmaj = pow(c12j / c6j, oneOverSix);
228 epsi = c6i * c6i / c12i;
229 epsj = c6j * c6j / c12j;
230 c6 = sqrt(epsi * epsj) * pow(0.5*(sigmai+sigmaj), 6);
231 c12 = sqrt(epsi * epsj) * pow(0.5*(sigmai+sigmaj), 12);
233 C6(nbfp, atnr, i, j) = c6*6.0;
234 C12(nbfp, atnr, i, j) = c12*12.0;
240 /* This routine sets fr->solvent_opt to the most common solvent in the
241 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
242 * the fr->solvent_type array with the correct type (or esolNO).
244 * Charge groups that fulfill the conditions but are not identical to the
245 * most common one will be marked as esolNO in the solvent_type array.
247 * TIP3p is identical to SPC for these purposes, so we call it
248 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
250 * NOTE: QM particle should not
251 * become an optimized solvent. Not even if there is only one charge
261 } solvent_parameters_t;
264 check_solvent_cg(const gmx_moltype_t *molt,
267 const unsigned char *qm_grpnr,
268 const t_grps *qm_grps,
270 int *n_solvent_parameters,
271 solvent_parameters_t **solvent_parameters_p,
281 real tmp_charge[4] = { 0.0 }; /* init to zero to make gcc4.8 happy */
282 int tmp_vdwtype[4] = { 0 }; /* init to zero to make gcc4.8 happy */
285 solvent_parameters_t *solvent_parameters;
287 /* We use a list with parameters for each solvent type.
288 * Every time we discover a new molecule that fulfills the basic
289 * conditions for a solvent we compare with the previous entries
290 * in these lists. If the parameters are the same we just increment
291 * the counter for that type, and otherwise we create a new type
292 * based on the current molecule.
294 * Once we've finished going through all molecules we check which
295 * solvent is most common, and mark all those molecules while we
296 * clear the flag on all others.
299 solvent_parameters = *solvent_parameters_p;
301 /* Mark the cg first as non optimized */
304 /* Check if this cg has no exclusions with atoms in other charge groups
305 * and all atoms inside the charge group excluded.
306 * We only have 3 or 4 atom solvent loops.
308 if (GET_CGINFO_EXCL_INTER(cginfo) ||
309 !GET_CGINFO_EXCL_INTRA(cginfo))
314 /* Get the indices of the first atom in this charge group */
315 j0 = molt->cgs.index[cg0];
316 j1 = molt->cgs.index[cg0+1];
318 /* Number of atoms in our molecule */
324 "Moltype '%s': there are %d atoms in this charge group\n",
328 /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
331 if (nj < 3 || nj > 4)
336 /* Check if we are doing QM on this group */
338 if (qm_grpnr != NULL)
340 for (j = j0; j < j1 && !qm; j++)
342 qm = (qm_grpnr[j] < qm_grps->nr - 1);
345 /* Cannot use solvent optimization with QM */
351 atom = molt->atoms.atom;
353 /* Still looks like a solvent, time to check parameters */
355 /* If it is perturbed (free energy) we can't use the solvent loops,
356 * so then we just skip to the next molecule.
360 for (j = j0; j < j1 && !perturbed; j++)
362 perturbed = PERTURBED(atom[j]);
370 /* Now it's only a question if the VdW and charge parameters
371 * are OK. Before doing the check we compare and see if they are
372 * identical to a possible previous solvent type.
373 * First we assign the current types and charges.
375 for (j = 0; j < nj; j++)
377 tmp_vdwtype[j] = atom[j0+j].type;
378 tmp_charge[j] = atom[j0+j].q;
381 /* Does it match any previous solvent type? */
382 for (k = 0; k < *n_solvent_parameters; k++)
387 /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
388 if ( (solvent_parameters[k].model == esolSPC && nj != 3) ||
389 (solvent_parameters[k].model == esolTIP4P && nj != 4) )
394 /* Check that types & charges match for all atoms in molecule */
395 for (j = 0; j < nj && match == TRUE; j++)
397 if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
401 if (tmp_charge[j] != solvent_parameters[k].charge[j])
408 /* Congratulations! We have a matched solvent.
409 * Flag it with this type for later processing.
412 solvent_parameters[k].count += nmol;
414 /* We are done with this charge group */
419 /* If we get here, we have a tentative new solvent type.
420 * Before we add it we must check that it fulfills the requirements
421 * of the solvent optimized loops. First determine which atoms have
424 for (j = 0; j < nj; j++)
427 tjA = tmp_vdwtype[j];
429 /* Go through all other tpes and see if any have non-zero
430 * VdW parameters when combined with this one.
432 for (k = 0; k < fr->ntype && (has_vdw[j] == FALSE); k++)
434 /* We already checked that the atoms weren't perturbed,
435 * so we only need to check state A now.
439 has_vdw[j] = (has_vdw[j] ||
440 (BHAMA(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
441 (BHAMB(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
442 (BHAMC(fr->nbfp, fr->ntype, tjA, k) != 0.0));
447 has_vdw[j] = (has_vdw[j] ||
448 (C6(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
449 (C12(fr->nbfp, fr->ntype, tjA, k) != 0.0));
454 /* Now we know all we need to make the final check and assignment. */
458 * For this we require thatn all atoms have charge,
459 * the charges on atom 2 & 3 should be the same, and only
460 * atom 1 might have VdW.
462 if (has_vdw[1] == FALSE &&
463 has_vdw[2] == FALSE &&
464 tmp_charge[0] != 0 &&
465 tmp_charge[1] != 0 &&
466 tmp_charge[2] == tmp_charge[1])
468 srenew(solvent_parameters, *n_solvent_parameters+1);
469 solvent_parameters[*n_solvent_parameters].model = esolSPC;
470 solvent_parameters[*n_solvent_parameters].count = nmol;
471 for (k = 0; k < 3; k++)
473 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
474 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
477 *cg_sp = *n_solvent_parameters;
478 (*n_solvent_parameters)++;
483 /* Or could it be a TIP4P?
484 * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
485 * Only atom 1 mght have VdW.
487 if (has_vdw[1] == FALSE &&
488 has_vdw[2] == FALSE &&
489 has_vdw[3] == FALSE &&
490 tmp_charge[0] == 0 &&
491 tmp_charge[1] != 0 &&
492 tmp_charge[2] == tmp_charge[1] &&
495 srenew(solvent_parameters, *n_solvent_parameters+1);
496 solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
497 solvent_parameters[*n_solvent_parameters].count = nmol;
498 for (k = 0; k < 4; k++)
500 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
501 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
504 *cg_sp = *n_solvent_parameters;
505 (*n_solvent_parameters)++;
509 *solvent_parameters_p = solvent_parameters;
513 check_solvent(FILE * fp,
514 const gmx_mtop_t * mtop,
516 cginfo_mb_t *cginfo_mb)
519 const gmx_moltype_t *molt;
520 int mb, mol, cg_mol, at_offset, am, cgm, i, nmol_ch, nmol;
521 int n_solvent_parameters;
522 solvent_parameters_t *solvent_parameters;
528 fprintf(debug, "Going to determine what solvent types we have.\n");
531 n_solvent_parameters = 0;
532 solvent_parameters = NULL;
533 /* Allocate temporary array for solvent type */
534 snew(cg_sp, mtop->nmolblock);
537 for (mb = 0; mb < mtop->nmolblock; mb++)
539 molt = &mtop->moltype[mtop->molblock[mb].type];
541 /* Here we have to loop over all individual molecules
542 * because we need to check for QMMM particles.
544 snew(cg_sp[mb], cginfo_mb[mb].cg_mod);
545 nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
546 nmol = mtop->molblock[mb].nmol/nmol_ch;
547 for (mol = 0; mol < nmol_ch; mol++)
550 am = mol*cgs->index[cgs->nr];
551 for (cg_mol = 0; cg_mol < cgs->nr; cg_mol++)
553 check_solvent_cg(molt, cg_mol, nmol,
554 mtop->groups.grpnr[egcQMMM] ?
555 mtop->groups.grpnr[egcQMMM]+at_offset+am : 0,
556 &mtop->groups.grps[egcQMMM],
558 &n_solvent_parameters, &solvent_parameters,
559 cginfo_mb[mb].cginfo[cgm+cg_mol],
560 &cg_sp[mb][cgm+cg_mol]);
563 at_offset += cgs->index[cgs->nr];
566 /* Puh! We finished going through all charge groups.
567 * Now find the most common solvent model.
570 /* Most common solvent this far */
572 for (i = 0; i < n_solvent_parameters; i++)
575 solvent_parameters[i].count > solvent_parameters[bestsp].count)
583 bestsol = solvent_parameters[bestsp].model;
591 for (mb = 0; mb < mtop->nmolblock; mb++)
593 cgs = &mtop->moltype[mtop->molblock[mb].type].cgs;
594 nmol = (mtop->molblock[mb].nmol*cgs->nr)/cginfo_mb[mb].cg_mod;
595 for (i = 0; i < cginfo_mb[mb].cg_mod; i++)
597 if (cg_sp[mb][i] == bestsp)
599 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], bestsol);
604 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], esolNO);
611 if (bestsol != esolNO && fp != NULL)
613 fprintf(fp, "\nEnabling %s-like water optimization for %d molecules.\n\n",
615 solvent_parameters[bestsp].count);
618 sfree(solvent_parameters);
619 fr->solvent_opt = bestsol;
623 acNONE = 0, acCONSTRAINT, acSETTLE
626 static cginfo_mb_t *init_cginfo_mb(FILE *fplog, const gmx_mtop_t *mtop,
627 t_forcerec *fr, gmx_bool bNoSolvOpt,
628 gmx_bool *bFEP_NonBonded,
629 gmx_bool *bExcl_IntraCGAll_InterCGNone)
632 const t_blocka *excl;
633 const gmx_moltype_t *molt;
634 const gmx_molblock_t *molb;
635 cginfo_mb_t *cginfo_mb;
638 int cg_offset, a_offset;
639 int mb, m, cg, a0, a1, gid, ai, j, aj, excl_nalloc;
643 gmx_bool bId, *bExcl, bExclIntraAll, bExclInter, bHaveVDW, bHaveQ, bHavePerturbedAtoms;
645 snew(cginfo_mb, mtop->nmolblock);
647 snew(type_VDW, fr->ntype);
648 for (ai = 0; ai < fr->ntype; ai++)
650 type_VDW[ai] = FALSE;
651 for (j = 0; j < fr->ntype; j++)
653 type_VDW[ai] = type_VDW[ai] ||
655 C6(fr->nbfp, fr->ntype, ai, j) != 0 ||
656 C12(fr->nbfp, fr->ntype, ai, j) != 0;
660 *bFEP_NonBonded = FALSE;
661 *bExcl_IntraCGAll_InterCGNone = TRUE;
664 snew(bExcl, excl_nalloc);
667 for (mb = 0; mb < mtop->nmolblock; mb++)
669 molb = &mtop->molblock[mb];
670 molt = &mtop->moltype[molb->type];
674 /* Check if the cginfo is identical for all molecules in this block.
675 * If so, we only need an array of the size of one molecule.
676 * Otherwise we make an array of #mol times #cgs per molecule.
679 for (m = 0; m < molb->nmol; m++)
681 int am = m*cgs->index[cgs->nr];
682 for (cg = 0; cg < cgs->nr; cg++)
685 a1 = cgs->index[cg+1];
686 if (ggrpnr(&mtop->groups, egcENER, a_offset+am+a0) !=
687 ggrpnr(&mtop->groups, egcENER, a_offset +a0))
691 if (mtop->groups.grpnr[egcQMMM] != NULL)
693 for (ai = a0; ai < a1; ai++)
695 if (mtop->groups.grpnr[egcQMMM][a_offset+am+ai] !=
696 mtop->groups.grpnr[egcQMMM][a_offset +ai])
705 cginfo_mb[mb].cg_start = cg_offset;
706 cginfo_mb[mb].cg_end = cg_offset + molb->nmol*cgs->nr;
707 cginfo_mb[mb].cg_mod = (bId ? 1 : molb->nmol)*cgs->nr;
708 snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
709 cginfo = cginfo_mb[mb].cginfo;
711 /* Set constraints flags for constrained atoms */
712 snew(a_con, molt->atoms.nr);
713 for (ftype = 0; ftype < F_NRE; ftype++)
715 if (interaction_function[ftype].flags & IF_CONSTRAINT)
720 for (ia = 0; ia < molt->ilist[ftype].nr; ia += 1+nral)
724 for (a = 0; a < nral; a++)
726 a_con[molt->ilist[ftype].iatoms[ia+1+a]] =
727 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
733 for (m = 0; m < (bId ? 1 : molb->nmol); m++)
736 int am = m*cgs->index[cgs->nr];
737 for (cg = 0; cg < cgs->nr; cg++)
740 a1 = cgs->index[cg+1];
742 /* Store the energy group in cginfo */
743 gid = ggrpnr(&mtop->groups, egcENER, a_offset+am+a0);
744 SET_CGINFO_GID(cginfo[cgm+cg], gid);
746 /* Check the intra/inter charge group exclusions */
747 if (a1-a0 > excl_nalloc)
749 excl_nalloc = a1 - a0;
750 srenew(bExcl, excl_nalloc);
752 /* bExclIntraAll: all intra cg interactions excluded
753 * bExclInter: any inter cg interactions excluded
755 bExclIntraAll = TRUE;
759 bHavePerturbedAtoms = FALSE;
760 for (ai = a0; ai < a1; ai++)
762 /* Check VDW and electrostatic interactions */
763 bHaveVDW = bHaveVDW || (type_VDW[molt->atoms.atom[ai].type] ||
764 type_VDW[molt->atoms.atom[ai].typeB]);
765 bHaveQ = bHaveQ || (molt->atoms.atom[ai].q != 0 ||
766 molt->atoms.atom[ai].qB != 0);
768 bHavePerturbedAtoms = bHavePerturbedAtoms || (PERTURBED(molt->atoms.atom[ai]) != 0);
770 /* Clear the exclusion list for atom ai */
771 for (aj = a0; aj < a1; aj++)
773 bExcl[aj-a0] = FALSE;
775 /* Loop over all the exclusions of atom ai */
776 for (j = excl->index[ai]; j < excl->index[ai+1]; j++)
779 if (aj < a0 || aj >= a1)
788 /* Check if ai excludes a0 to a1 */
789 for (aj = a0; aj < a1; aj++)
793 bExclIntraAll = FALSE;
800 SET_CGINFO_CONSTR(cginfo[cgm+cg]);
803 SET_CGINFO_SETTLE(cginfo[cgm+cg]);
811 SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
815 SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
817 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
819 /* The size in cginfo is currently only read with DD */
820 gmx_fatal(FARGS, "A charge group has size %d which is larger than the limit of %d atoms", a1-a0, MAX_CHARGEGROUP_SIZE);
824 SET_CGINFO_HAS_VDW(cginfo[cgm+cg]);
828 SET_CGINFO_HAS_Q(cginfo[cgm+cg]);
830 if (bHavePerturbedAtoms && fr->efep != efepNO)
832 SET_CGINFO_FEP(cginfo[cgm+cg]);
833 *bFEP_NonBonded = TRUE;
835 /* Store the charge group size */
836 SET_CGINFO_NATOMS(cginfo[cgm+cg], a1-a0);
838 if (!bExclIntraAll || bExclInter)
840 *bExcl_IntraCGAll_InterCGNone = FALSE;
847 cg_offset += molb->nmol*cgs->nr;
848 a_offset += molb->nmol*cgs->index[cgs->nr];
852 /* the solvent optimizer is called after the QM is initialized,
853 * because we don't want to have the QM subsystemto become an
857 check_solvent(fplog, mtop, fr, cginfo_mb);
859 if (getenv("GMX_NO_SOLV_OPT"))
863 fprintf(fplog, "Found environment variable GMX_NO_SOLV_OPT.\n"
864 "Disabling all solvent optimization\n");
866 fr->solvent_opt = esolNO;
870 fr->solvent_opt = esolNO;
872 if (!fr->solvent_opt)
874 for (mb = 0; mb < mtop->nmolblock; mb++)
876 for (cg = 0; cg < cginfo_mb[mb].cg_mod; cg++)
878 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg], esolNO);
886 static int *cginfo_expand(int nmb, cginfo_mb_t *cgi_mb)
891 ncg = cgi_mb[nmb-1].cg_end;
894 for (cg = 0; cg < ncg; cg++)
896 while (cg >= cgi_mb[mb].cg_end)
901 cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
907 static void set_chargesum(FILE *log, t_forcerec *fr, const gmx_mtop_t *mtop)
909 /*This now calculates sum for q and c6*/
910 double qsum, q2sum, q, c6sum, c6;
912 const t_atoms *atoms;
917 for (mb = 0; mb < mtop->nmolblock; mb++)
919 nmol = mtop->molblock[mb].nmol;
920 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
921 for (i = 0; i < atoms->nr; i++)
923 q = atoms->atom[i].q;
926 c6 = mtop->ffparams.iparams[atoms->atom[i].type*(mtop->ffparams.atnr+1)].lj.c6;
931 fr->q2sum[0] = q2sum;
932 fr->c6sum[0] = c6sum;
934 if (fr->efep != efepNO)
939 for (mb = 0; mb < mtop->nmolblock; mb++)
941 nmol = mtop->molblock[mb].nmol;
942 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
943 for (i = 0; i < atoms->nr; i++)
945 q = atoms->atom[i].qB;
948 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB*(mtop->ffparams.atnr+1)].lj.c6;
952 fr->q2sum[1] = q2sum;
953 fr->c6sum[1] = c6sum;
958 fr->qsum[1] = fr->qsum[0];
959 fr->q2sum[1] = fr->q2sum[0];
960 fr->c6sum[1] = fr->c6sum[0];
964 if (fr->efep == efepNO)
966 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
970 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n",
971 fr->qsum[0], fr->qsum[1]);
976 void update_forcerec(t_forcerec *fr, matrix box)
978 if (fr->eeltype == eelGRF)
980 calc_rffac(NULL, fr->eeltype, fr->epsilon_r, fr->epsilon_rf,
981 fr->rcoulomb, fr->temp, fr->zsquare, box,
982 &fr->kappa, &fr->k_rf, &fr->c_rf);
986 void set_avcsixtwelve(FILE *fplog, t_forcerec *fr, const gmx_mtop_t *mtop)
988 const t_atoms *atoms, *atoms_tpi;
989 const t_blocka *excl;
990 int mb, nmol, nmolc, i, j, tpi, tpj, j1, j2, k, nexcl, q;
991 gmx_int64_t npair, npair_ij, tmpi, tmpj;
992 double csix, ctwelve;
996 real *nbfp_comb = NULL;
1002 /* For LJ-PME, we want to correct for the difference between the
1003 * actual C6 values and the C6 values used by the LJ-PME based on
1004 * combination rules. */
1006 if (EVDW_PME(fr->vdwtype))
1008 nbfp_comb = mk_nbfp_combination_rule(&mtop->ffparams,
1009 (fr->ljpme_combination_rule == eljpmeLB) ? eCOMB_ARITHMETIC : eCOMB_GEOMETRIC);
1010 for (tpi = 0; tpi < ntp; ++tpi)
1012 for (tpj = 0; tpj < ntp; ++tpj)
1014 C6(nbfp_comb, ntp, tpi, tpj) =
1015 C6(nbfp, ntp, tpi, tpj) - C6(nbfp_comb, ntp, tpi, tpj);
1016 C12(nbfp_comb, ntp, tpi, tpj) = C12(nbfp, ntp, tpi, tpj);
1021 for (q = 0; q < (fr->efep == efepNO ? 1 : 2); q++)
1029 /* Count the types so we avoid natoms^2 operations */
1030 snew(typecount, ntp);
1031 gmx_mtop_count_atomtypes(mtop, q, typecount);
1033 for (tpi = 0; tpi < ntp; tpi++)
1035 for (tpj = tpi; tpj < ntp; tpj++)
1037 tmpi = typecount[tpi];
1038 tmpj = typecount[tpj];
1041 npair_ij = tmpi*tmpj;
1045 npair_ij = tmpi*(tmpi - 1)/2;
1049 /* nbfp now includes the 6.0 derivative prefactor */
1050 csix += npair_ij*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1054 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1055 csix += npair_ij* C6(nbfp, ntp, tpi, tpj)/6.0;
1056 ctwelve += npair_ij* C12(nbfp, ntp, tpi, tpj)/12.0;
1062 /* Subtract the excluded pairs.
1063 * The main reason for substracting exclusions is that in some cases
1064 * some combinations might never occur and the parameters could have
1065 * any value. These unused values should not influence the dispersion
1068 for (mb = 0; mb < mtop->nmolblock; mb++)
1070 nmol = mtop->molblock[mb].nmol;
1071 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
1072 excl = &mtop->moltype[mtop->molblock[mb].type].excls;
1073 for (i = 0; (i < atoms->nr); i++)
1077 tpi = atoms->atom[i].type;
1081 tpi = atoms->atom[i].typeB;
1083 j1 = excl->index[i];
1084 j2 = excl->index[i+1];
1085 for (j = j1; j < j2; j++)
1092 tpj = atoms->atom[k].type;
1096 tpj = atoms->atom[k].typeB;
1100 /* nbfp now includes the 6.0 derivative prefactor */
1101 csix -= nmol*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1105 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1106 csix -= nmol*C6 (nbfp, ntp, tpi, tpj)/6.0;
1107 ctwelve -= nmol*C12(nbfp, ntp, tpi, tpj)/12.0;
1117 /* Only correct for the interaction of the test particle
1118 * with the rest of the system.
1121 &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].atoms;
1124 for (mb = 0; mb < mtop->nmolblock; mb++)
1126 nmol = mtop->molblock[mb].nmol;
1127 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
1128 for (j = 0; j < atoms->nr; j++)
1131 /* Remove the interaction of the test charge group
1134 if (mb == mtop->nmolblock-1)
1138 if (mb == 0 && nmol == 1)
1140 gmx_fatal(FARGS, "Old format tpr with TPI, please generate a new tpr file");
1145 tpj = atoms->atom[j].type;
1149 tpj = atoms->atom[j].typeB;
1151 for (i = 0; i < fr->n_tpi; i++)
1155 tpi = atoms_tpi->atom[i].type;
1159 tpi = atoms_tpi->atom[i].typeB;
1163 /* nbfp now includes the 6.0 derivative prefactor */
1164 csix += nmolc*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1168 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1169 csix += nmolc*C6 (nbfp, ntp, tpi, tpj)/6.0;
1170 ctwelve += nmolc*C12(nbfp, ntp, tpi, tpj)/12.0;
1177 if (npair - nexcl <= 0 && fplog)
1179 fprintf(fplog, "\nWARNING: There are no atom pairs for dispersion correction\n\n");
1185 csix /= npair - nexcl;
1186 ctwelve /= npair - nexcl;
1190 fprintf(debug, "Counted %d exclusions\n", nexcl);
1191 fprintf(debug, "Average C6 parameter is: %10g\n", (double)csix);
1192 fprintf(debug, "Average C12 parameter is: %10g\n", (double)ctwelve);
1194 fr->avcsix[q] = csix;
1195 fr->avctwelve[q] = ctwelve;
1198 if (EVDW_PME(fr->vdwtype))
1205 if (fr->eDispCorr == edispcAllEner ||
1206 fr->eDispCorr == edispcAllEnerPres)
1208 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1209 fr->avcsix[0], fr->avctwelve[0]);
1213 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e\n", fr->avcsix[0]);
1219 static void set_bham_b_max(FILE *fplog, t_forcerec *fr,
1220 const gmx_mtop_t *mtop)
1222 const t_atoms *at1, *at2;
1223 int mt1, mt2, i, j, tpi, tpj, ntypes;
1229 fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
1236 for (mt1 = 0; mt1 < mtop->nmoltype; mt1++)
1238 at1 = &mtop->moltype[mt1].atoms;
1239 for (i = 0; (i < at1->nr); i++)
1241 tpi = at1->atom[i].type;
1244 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
1247 for (mt2 = mt1; mt2 < mtop->nmoltype; mt2++)
1249 at2 = &mtop->moltype[mt2].atoms;
1250 for (j = 0; (j < at2->nr); j++)
1252 tpj = at2->atom[j].type;
1255 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
1257 b = BHAMB(nbfp, ntypes, tpi, tpj);
1258 if (b > fr->bham_b_max)
1262 if ((b < bmin) || (bmin == -1))
1272 fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n",
1273 bmin, fr->bham_b_max);
1277 static void make_nbf_tables(FILE *fp, const output_env_t oenv,
1278 t_forcerec *fr, real rtab,
1279 const t_commrec *cr,
1280 const char *tabfn, char *eg1, char *eg2,
1290 fprintf(debug, "No table file name passed, can not read table, can not do non-bonded interactions\n");
1295 sprintf(buf, "%s", tabfn);
1298 /* Append the two energy group names */
1299 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
1300 eg1, eg2, ftp2ext(efXVG));
1302 nbl->table_elec_vdw = make_tables(fp, oenv, fr, MASTER(cr), buf, rtab, 0);
1303 /* Copy the contents of the table to separate coulomb and LJ tables too,
1304 * to improve cache performance.
1306 /* For performance reasons we want
1307 * the table data to be aligned to 16-byte. The pointers could be freed
1308 * but currently aren't.
1310 nbl->table_elec.interaction = GMX_TABLE_INTERACTION_ELEC;
1311 nbl->table_elec.format = nbl->table_elec_vdw.format;
1312 nbl->table_elec.r = nbl->table_elec_vdw.r;
1313 nbl->table_elec.n = nbl->table_elec_vdw.n;
1314 nbl->table_elec.scale = nbl->table_elec_vdw.scale;
1315 nbl->table_elec.scale_exp = nbl->table_elec_vdw.scale_exp;
1316 nbl->table_elec.formatsize = nbl->table_elec_vdw.formatsize;
1317 nbl->table_elec.ninteractions = 1;
1318 nbl->table_elec.stride = nbl->table_elec.formatsize * nbl->table_elec.ninteractions;
1319 snew_aligned(nbl->table_elec.data, nbl->table_elec.stride*(nbl->table_elec.n+1), 32);
1321 nbl->table_vdw.interaction = GMX_TABLE_INTERACTION_VDWREP_VDWDISP;
1322 nbl->table_vdw.format = nbl->table_elec_vdw.format;
1323 nbl->table_vdw.r = nbl->table_elec_vdw.r;
1324 nbl->table_vdw.n = nbl->table_elec_vdw.n;
1325 nbl->table_vdw.scale = nbl->table_elec_vdw.scale;
1326 nbl->table_vdw.scale_exp = nbl->table_elec_vdw.scale_exp;
1327 nbl->table_vdw.formatsize = nbl->table_elec_vdw.formatsize;
1328 nbl->table_vdw.ninteractions = 2;
1329 nbl->table_vdw.stride = nbl->table_vdw.formatsize * nbl->table_vdw.ninteractions;
1330 snew_aligned(nbl->table_vdw.data, nbl->table_vdw.stride*(nbl->table_vdw.n+1), 32);
1332 for (i = 0; i <= nbl->table_elec_vdw.n; i++)
1334 for (j = 0; j < 4; j++)
1336 nbl->table_elec.data[4*i+j] = nbl->table_elec_vdw.data[12*i+j];
1338 for (j = 0; j < 8; j++)
1340 nbl->table_vdw.data[8*i+j] = nbl->table_elec_vdw.data[12*i+4+j];
1345 /*!\brief If there's bonded interactions of type \c ftype1 or \c
1346 * ftype2 present in the topology, build an array of the number of
1347 * interactions present for each bonded interaction index found in the
1350 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1351 * valid type with that parameter.
1353 * \c count will be reallocated as necessary to fit the largest bonded
1354 * interaction index found, and its current size will be returned in
1355 * \c ncount. It will contain zero for every bonded interaction index
1356 * for which no interactions are present in the topology.
1358 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t *mtop,
1359 int *ncount, int **count)
1361 const gmx_moltype_t *molt;
1363 int mt, ftype, stride, i, j, tabnr;
1365 // Loop over all moleculetypes
1366 for (mt = 0; mt < mtop->nmoltype; mt++)
1368 molt = &mtop->moltype[mt];
1369 // Loop over all interaction types
1370 for (ftype = 0; ftype < F_NRE; ftype++)
1372 // If the current interaction type is one of the types whose tables we're trying to count...
1373 if (ftype == ftype1 || ftype == ftype2)
1375 il = &molt->ilist[ftype];
1376 stride = 1 + NRAL(ftype);
1377 // ... and there are actually some interactions for this type
1378 for (i = 0; i < il->nr; i += stride)
1380 // Find out which table index the user wanted
1381 tabnr = mtop->ffparams.iparams[il->iatoms[i]].tab.table;
1384 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
1386 // Make room for this index in the data structure
1387 if (tabnr >= *ncount)
1389 srenew(*count, tabnr+1);
1390 for (j = *ncount; j < tabnr+1; j++)
1396 // Record that this table index is used and must have a valid file
1404 /*!\brief If there's bonded interactions of flavour \c tabext and type
1405 * \c ftype1 or \c ftype2 present in the topology, seek them in the
1406 * list of filenames passed to mdrun, and make bonded tables from
1409 * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1410 * valid type with that parameter.
1412 * A fatal error occurs if no matching filename is found.
1414 static bondedtable_t *make_bonded_tables(FILE *fplog,
1415 int ftype1, int ftype2,
1416 const gmx_mtop_t *mtop,
1417 const t_filenm *tabbfnm,
1427 count_tables(ftype1, ftype2, mtop, &ncount, &count);
1429 // Are there any relevant tabulated bond interactions?
1433 for (int i = 0; i < ncount; i++)
1435 // Do any interactions exist that requires this table?
1438 // This pattern enforces the current requirement that
1439 // table filenames end in a characteristic sequence
1440 // before the file type extension, and avoids table 13
1441 // being recognized and used for table 1.
1442 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
1443 bool madeTable = false;
1444 for (int j = 0; j < tabbfnm->nfiles && !madeTable; ++j)
1446 std::string filename(tabbfnm->fns[j]);
1447 if (gmx::endsWith(filename, patternToFind))
1449 // Finally read the table from the file found
1450 tab[i] = make_bonded_table(fplog, tabbfnm->fns[j], NRAL(ftype1)-2);
1456 bool isPlural = (ftype2 != -1);
1457 gmx_fatal(FARGS, "Tabulated interaction of type '%s%s%s' with index %d cannot be used because no table file whose name matched '%s' was passed via the gmx mdrun -tableb command-line option.",
1458 interaction_function[ftype1].longname,
1459 isPlural ? "' or '" : "",
1460 isPlural ? interaction_function[ftype2].longname : "",
1462 patternToFind.c_str());
1472 void forcerec_set_ranges(t_forcerec *fr,
1473 int ncg_home, int ncg_force,
1475 int natoms_force_constr, int natoms_f_novirsum)
1480 /* fr->ncg_force is unused in the standard code,
1481 * but it can be useful for modified code dealing with charge groups.
1483 fr->ncg_force = ncg_force;
1484 fr->natoms_force = natoms_force;
1485 fr->natoms_force_constr = natoms_force_constr;
1487 if (fr->natoms_force_constr > fr->nalloc_force)
1489 fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1493 srenew(fr->f_twin, fr->nalloc_force);
1497 if (fr->bF_NoVirSum)
1499 fr->f_novirsum_n = natoms_f_novirsum;
1500 if (fr->f_novirsum_n > fr->f_novirsum_nalloc)
1502 fr->f_novirsum_nalloc = over_alloc_dd(fr->f_novirsum_n);
1503 srenew(fr->f_novirsum_alloc, fr->f_novirsum_nalloc);
1508 fr->f_novirsum_n = 0;
1512 static real cutoff_inf(real cutoff)
1516 cutoff = GMX_CUTOFF_INF;
1522 static void make_adress_tf_tables(FILE *fp, const output_env_t oenv,
1523 t_forcerec *fr, const t_inputrec *ir,
1524 const char *tabfn, const gmx_mtop_t *mtop,
1532 gmx_fatal(FARGS, "No thermoforce table file given. Use -tabletf to specify a file\n");
1536 snew(fr->atf_tabs, ir->adress->n_tf_grps);
1538 sprintf(buf, "%s", tabfn);
1539 for (i = 0; i < ir->adress->n_tf_grps; i++)
1541 j = ir->adress->tf_table_index[i]; /* get energy group index */
1542 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "tf_%s.%s",
1543 *(mtop->groups.grpname[mtop->groups.grps[egcENER].nm_ind[j]]), ftp2ext(efXVG));
1546 fprintf(fp, "loading tf table for energygrp index %d from %s\n", ir->adress->tf_table_index[i], buf);
1548 fr->atf_tabs[i] = make_atf_table(fp, oenv, fr, buf, box);
1553 gmx_bool can_use_allvsall(const t_inputrec *ir, gmx_bool bPrintNote, t_commrec *cr, FILE *fp)
1560 ir->rcoulomb == 0 &&
1562 ir->ePBC == epbcNONE &&
1563 ir->vdwtype == evdwCUT &&
1564 ir->coulombtype == eelCUT &&
1565 ir->efep == efepNO &&
1566 (ir->implicit_solvent == eisNO ||
1567 (ir->implicit_solvent == eisGBSA && (ir->gb_algorithm == egbSTILL ||
1568 ir->gb_algorithm == egbHCT ||
1569 ir->gb_algorithm == egbOBC))) &&
1570 getenv("GMX_NO_ALLVSALL") == NULL
1573 if (bAllvsAll && ir->opts.ngener > 1)
1575 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";
1581 fprintf(stderr, "\n%s\n", note);
1585 fprintf(fp, "\n%s\n", note);
1591 if (bAllvsAll && fp && MASTER(cr))
1593 fprintf(fp, "\nUsing SIMD all-vs-all kernels.\n\n");
1600 gmx_bool nbnxn_gpu_acceleration_supported(FILE *fplog,
1601 const t_commrec *cr,
1602 const t_inputrec *ir,
1605 if (bRerunMD && ir->opts.ngener > 1)
1607 /* Rerun execution time is dominated by I/O and pair search,
1608 * so GPUs are not very useful, plus they do not support more
1609 * than one energy group. If the user requested GPUs
1610 * explicitly, a fatal error is given later. With non-reruns,
1611 * we fall back to a single whole-of system energy group
1612 * (which runs much faster than a multiple-energy-groups
1613 * implementation would), and issue a note in the .log
1614 * file. Users can re-run if they want the information. */
1615 md_print_warn(cr, fplog, "Rerun with energy groups is not implemented for GPUs, falling back to the CPU\n");
1622 gmx_bool nbnxn_simd_supported(FILE *fplog,
1623 const t_commrec *cr,
1624 const t_inputrec *ir)
1626 if (ir->vdwtype == evdwPME && ir->ljpme_combination_rule == eljpmeLB)
1628 /* LJ PME with LB combination rule does 7 mesh operations.
1629 * This so slow that we don't compile SIMD non-bonded kernels
1631 md_print_warn(cr, fplog, "LJ-PME with Lorentz-Berthelot is not supported with SIMD kernels, falling back to plain C kernels\n");
1639 static void pick_nbnxn_kernel_cpu(const t_inputrec gmx_unused *ir,
1643 *kernel_type = nbnxnk4x4_PlainC;
1644 *ewald_excl = ewaldexclTable;
1646 #ifdef GMX_NBNXN_SIMD
1648 #ifdef GMX_NBNXN_SIMD_4XN
1649 *kernel_type = nbnxnk4xN_SIMD_4xN;
1651 #ifdef GMX_NBNXN_SIMD_2XNN
1652 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1655 #if defined GMX_NBNXN_SIMD_2XNN && defined GMX_NBNXN_SIMD_4XN
1656 /* We need to choose if we want 2x(N+N) or 4xN kernels.
1657 * Currently this is based on the SIMD acceleration choice,
1658 * but it might be better to decide this at runtime based on CPU.
1660 * 4xN calculates more (zero) interactions, but has less pair-search
1661 * work and much better kernel instruction scheduling.
1663 * Up till now we have only seen that on Intel Sandy/Ivy Bridge,
1664 * which doesn't have FMA, both the analytical and tabulated Ewald
1665 * kernels have similar pair rates for 4x8 and 2x(4+4), so we choose
1666 * 2x(4+4) because it results in significantly fewer pairs.
1667 * For RF, the raw pair rate of the 4x8 kernel is higher than 2x(4+4),
1668 * 10% with HT, 50% without HT. As we currently don't detect the actual
1669 * use of HT, use 4x8 to avoid a potential performance hit.
1670 * On Intel Haswell 4x8 is always faster.
1672 *kernel_type = nbnxnk4xN_SIMD_4xN;
1674 #ifndef GMX_SIMD_HAVE_FMA
1675 if (EEL_PME_EWALD(ir->coulombtype) ||
1676 EVDW_PME(ir->vdwtype))
1678 /* We have Ewald kernels without FMA (Intel Sandy/Ivy Bridge).
1679 * There are enough instructions to make 2x(4+4) efficient.
1681 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1684 #endif /* GMX_NBNXN_SIMD_2XNN && GMX_NBNXN_SIMD_4XN */
1687 if (getenv("GMX_NBNXN_SIMD_4XN") != NULL)
1689 #ifdef GMX_NBNXN_SIMD_4XN
1690 *kernel_type = nbnxnk4xN_SIMD_4xN;
1692 gmx_fatal(FARGS, "SIMD 4xN kernels requested, but GROMACS has been compiled without support for these kernels");
1695 if (getenv("GMX_NBNXN_SIMD_2XNN") != NULL)
1697 #ifdef GMX_NBNXN_SIMD_2XNN
1698 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1700 gmx_fatal(FARGS, "SIMD 2x(N+N) kernels requested, but GROMACS has been compiled without support for these kernels");
1704 /* Analytical Ewald exclusion correction is only an option in
1706 * Since table lookup's don't parallelize with SIMD, analytical
1707 * will probably always be faster for a SIMD width of 8 or more.
1708 * With FMA analytical is sometimes faster for a width if 4 as well.
1709 * On BlueGene/Q, this is faster regardless of precision.
1710 * In single precision, this is faster on Bulldozer.
1712 #if GMX_SIMD_REAL_WIDTH >= 8 || \
1713 (GMX_SIMD_REAL_WIDTH >= 4 && defined GMX_SIMD_HAVE_FMA && !defined GMX_DOUBLE) || \
1714 defined GMX_SIMD_IBM_QPX
1715 *ewald_excl = ewaldexclAnalytical;
1717 if (getenv("GMX_NBNXN_EWALD_TABLE") != NULL)
1719 *ewald_excl = ewaldexclTable;
1721 if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != NULL)
1723 *ewald_excl = ewaldexclAnalytical;
1727 #endif /* GMX_NBNXN_SIMD */
1731 const char *lookup_nbnxn_kernel_name(int kernel_type)
1733 const char *returnvalue = NULL;
1734 switch (kernel_type)
1737 returnvalue = "not set";
1739 case nbnxnk4x4_PlainC:
1740 returnvalue = "plain C";
1742 case nbnxnk4xN_SIMD_4xN:
1743 case nbnxnk4xN_SIMD_2xNN:
1744 #ifdef GMX_NBNXN_SIMD
1745 #if defined GMX_SIMD_X86_SSE2
1746 returnvalue = "SSE2";
1747 #elif defined GMX_SIMD_X86_SSE4_1
1748 returnvalue = "SSE4.1";
1749 #elif defined GMX_SIMD_X86_AVX_128_FMA
1750 returnvalue = "AVX_128_FMA";
1751 #elif defined GMX_SIMD_X86_AVX_256
1752 returnvalue = "AVX_256";
1753 #elif defined GMX_SIMD_X86_AVX2_256
1754 returnvalue = "AVX2_256";
1756 returnvalue = "SIMD";
1758 #else /* GMX_NBNXN_SIMD */
1759 returnvalue = "not available";
1760 #endif /* GMX_NBNXN_SIMD */
1762 case nbnxnk8x8x8_GPU: returnvalue = "GPU"; break;
1763 case nbnxnk8x8x8_PlainC: returnvalue = "plain C"; break;
1767 gmx_fatal(FARGS, "Illegal kernel type selected");
1774 static void pick_nbnxn_kernel(FILE *fp,
1775 const t_commrec *cr,
1776 gmx_bool use_simd_kernels,
1778 gmx_bool bEmulateGPU,
1779 const t_inputrec *ir,
1782 gmx_bool bDoNonbonded)
1784 assert(kernel_type);
1786 *kernel_type = nbnxnkNotSet;
1787 *ewald_excl = ewaldexclTable;
1791 *kernel_type = nbnxnk8x8x8_PlainC;
1795 md_print_warn(cr, fp, "Emulating a GPU run on the CPU (slow)");
1800 *kernel_type = nbnxnk8x8x8_GPU;
1803 if (*kernel_type == nbnxnkNotSet)
1805 if (use_simd_kernels &&
1806 nbnxn_simd_supported(fp, cr, ir))
1808 pick_nbnxn_kernel_cpu(ir, kernel_type, ewald_excl);
1812 *kernel_type = nbnxnk4x4_PlainC;
1816 if (bDoNonbonded && fp != NULL)
1818 fprintf(fp, "\nUsing %s %dx%d non-bonded kernels\n\n",
1819 lookup_nbnxn_kernel_name(*kernel_type),
1820 nbnxn_kernel_to_ci_size(*kernel_type),
1821 nbnxn_kernel_to_cj_size(*kernel_type));
1823 if (nbnxnk4x4_PlainC == *kernel_type ||
1824 nbnxnk8x8x8_PlainC == *kernel_type)
1826 md_print_warn(cr, fp,
1827 "WARNING: Using the slow %s kernels. This should\n"
1828 "not happen during routine usage on supported platforms.\n\n",
1829 lookup_nbnxn_kernel_name(*kernel_type));
1834 static void pick_nbnxn_resources(FILE *fp,
1835 const t_commrec *cr,
1836 const gmx_hw_info_t *hwinfo,
1837 gmx_bool bDoNonbonded,
1839 gmx_bool *bEmulateGPU,
1840 const gmx_gpu_opt_t *gpu_opt)
1842 gmx_bool bEmulateGPUEnvVarSet;
1843 char gpu_err_str[STRLEN];
1847 bEmulateGPUEnvVarSet = (getenv("GMX_EMULATE_GPU") != NULL);
1849 /* Run GPU emulation mode if GMX_EMULATE_GPU is defined. Because
1850 * GPUs (currently) only handle non-bonded calculations, we will
1851 * automatically switch to emulation if non-bonded calculations are
1852 * turned off via GMX_NO_NONBONDED - this is the simple and elegant
1853 * way to turn off GPU initialization, data movement, and cleanup.
1855 * GPU emulation can be useful to assess the performance one can expect by
1856 * adding GPU(s) to the machine. The conditional below allows this even
1857 * if mdrun is compiled without GPU acceleration support.
1858 * Note that you should freezing the system as otherwise it will explode.
1860 *bEmulateGPU = (bEmulateGPUEnvVarSet ||
1861 (!bDoNonbonded && gpu_opt->n_dev_use > 0));
1863 /* Enable GPU mode when GPUs are available or no GPU emulation is requested.
1865 if (gpu_opt->n_dev_use > 0 && !(*bEmulateGPU))
1867 /* Each PP node will use the intra-node id-th device from the
1868 * list of detected/selected GPUs. */
1869 if (!init_gpu(fp, cr->rank_pp_intranode, gpu_err_str,
1870 &hwinfo->gpu_info, gpu_opt))
1872 /* At this point the init should never fail as we made sure that
1873 * we have all the GPUs we need. If it still does, we'll bail. */
1874 /* TODO the decorating of gpu_err_str is nicer if it
1875 happens inside init_gpu. Out here, the decorating with
1876 the MPI rank makes sense. */
1877 gmx_fatal(FARGS, "On rank %d failed to initialize GPU #%d: %s",
1879 get_gpu_device_id(&hwinfo->gpu_info, gpu_opt,
1880 cr->rank_pp_intranode),
1884 /* Here we actually turn on hardware GPU acceleration */
1889 gmx_bool uses_simple_tables(int cutoff_scheme,
1890 nonbonded_verlet_t *nbv,
1893 gmx_bool bUsesSimpleTables = TRUE;
1896 switch (cutoff_scheme)
1899 bUsesSimpleTables = TRUE;
1902 assert(NULL != nbv && NULL != nbv->grp);
1903 grp_index = (group < 0) ? 0 : (nbv->ngrp - 1);
1904 bUsesSimpleTables = nbnxn_kernel_pairlist_simple(nbv->grp[grp_index].kernel_type);
1907 gmx_incons("unimplemented");
1909 return bUsesSimpleTables;
1912 static void init_ewald_f_table(interaction_const_t *ic,
1917 /* Get the Ewald table spacing based on Coulomb and/or LJ
1918 * Ewald coefficients and rtol.
1920 ic->tabq_scale = ewald_spline3_table_scale(ic);
1922 if (ic->cutoff_scheme == ecutsVERLET)
1924 maxr = ic->rcoulomb;
1928 maxr = std::max(ic->rcoulomb, rtab);
1930 ic->tabq_size = static_cast<int>(maxr*ic->tabq_scale) + 2;
1932 sfree_aligned(ic->tabq_coul_FDV0);
1933 sfree_aligned(ic->tabq_coul_F);
1934 sfree_aligned(ic->tabq_coul_V);
1936 sfree_aligned(ic->tabq_vdw_FDV0);
1937 sfree_aligned(ic->tabq_vdw_F);
1938 sfree_aligned(ic->tabq_vdw_V);
1940 if (ic->eeltype == eelEWALD || EEL_PME(ic->eeltype))
1942 /* Create the original table data in FDV0 */
1943 snew_aligned(ic->tabq_coul_FDV0, ic->tabq_size*4, 32);
1944 snew_aligned(ic->tabq_coul_F, ic->tabq_size, 32);
1945 snew_aligned(ic->tabq_coul_V, ic->tabq_size, 32);
1946 table_spline3_fill_ewald_lr(ic->tabq_coul_F, ic->tabq_coul_V, ic->tabq_coul_FDV0,
1947 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_q, v_q_ewald_lr);
1950 if (EVDW_PME(ic->vdwtype))
1952 snew_aligned(ic->tabq_vdw_FDV0, ic->tabq_size*4, 32);
1953 snew_aligned(ic->tabq_vdw_F, ic->tabq_size, 32);
1954 snew_aligned(ic->tabq_vdw_V, ic->tabq_size, 32);
1955 table_spline3_fill_ewald_lr(ic->tabq_vdw_F, ic->tabq_vdw_V, ic->tabq_vdw_FDV0,
1956 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_lj, v_lj_ewald_lr);
1960 void init_interaction_const_tables(FILE *fp,
1961 interaction_const_t *ic,
1964 if (ic->eeltype == eelEWALD || EEL_PME(ic->eeltype) || EVDW_PME(ic->vdwtype))
1966 init_ewald_f_table(ic, rtab);
1970 fprintf(fp, "Initialized non-bonded Ewald correction tables, spacing: %.2e size: %d\n\n",
1971 1/ic->tabq_scale, ic->tabq_size);
1976 static void clear_force_switch_constants(shift_consts_t *sc)
1983 static void force_switch_constants(real p,
1987 /* Here we determine the coefficient for shifting the force to zero
1988 * between distance rsw and the cut-off rc.
1989 * For a potential of r^-p, we have force p*r^-(p+1).
1990 * But to save flops we absorb p in the coefficient.
1992 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
1993 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
1995 sc->c2 = ((p + 1)*rsw - (p + 4)*rc)/(pow(rc, p + 2)*pow(rc - rsw, 2));
1996 sc->c3 = -((p + 1)*rsw - (p + 3)*rc)/(pow(rc, p + 2)*pow(rc - rsw, 3));
1997 sc->cpot = -pow(rc, -p) + p*sc->c2/3*pow(rc - rsw, 3) + p*sc->c3/4*pow(rc - rsw, 4);
2000 static void potential_switch_constants(real rsw, real rc,
2001 switch_consts_t *sc)
2003 /* The switch function is 1 at rsw and 0 at rc.
2004 * The derivative and second derivate are zero at both ends.
2005 * rsw = max(r - r_switch, 0)
2006 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
2007 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
2008 * force = force*dsw - potential*sw
2011 sc->c3 = -10*pow(rc - rsw, -3);
2012 sc->c4 = 15*pow(rc - rsw, -4);
2013 sc->c5 = -6*pow(rc - rsw, -5);
2016 /*! \brief Construct interaction constants
2018 * This data is used (particularly) by search and force code for
2019 * short-range interactions. Many of these are constant for the whole
2020 * simulation; some are constant only after PME tuning completes.
2023 init_interaction_const(FILE *fp,
2024 interaction_const_t **interaction_const,
2025 const t_forcerec *fr)
2027 interaction_const_t *ic;
2028 const real minusSix = -6.0;
2029 const real minusTwelve = -12.0;
2033 ic->cutoff_scheme = fr->cutoff_scheme;
2035 /* Just allocate something so we can free it */
2036 snew_aligned(ic->tabq_coul_FDV0, 16, 32);
2037 snew_aligned(ic->tabq_coul_F, 16, 32);
2038 snew_aligned(ic->tabq_coul_V, 16, 32);
2040 ic->rlist = fr->rlist;
2041 ic->rlistlong = fr->rlistlong;
2044 ic->vdwtype = fr->vdwtype;
2045 ic->vdw_modifier = fr->vdw_modifier;
2046 ic->rvdw = fr->rvdw;
2047 ic->rvdw_switch = fr->rvdw_switch;
2048 ic->ewaldcoeff_lj = fr->ewaldcoeff_lj;
2049 ic->ljpme_comb_rule = fr->ljpme_combination_rule;
2050 ic->sh_lj_ewald = 0;
2051 clear_force_switch_constants(&ic->dispersion_shift);
2052 clear_force_switch_constants(&ic->repulsion_shift);
2054 switch (ic->vdw_modifier)
2056 case eintmodPOTSHIFT:
2057 /* Only shift the potential, don't touch the force */
2058 ic->dispersion_shift.cpot = -pow(ic->rvdw, minusSix);
2059 ic->repulsion_shift.cpot = -pow(ic->rvdw, minusTwelve);
2060 if (EVDW_PME(ic->vdwtype))
2064 crc2 = sqr(ic->ewaldcoeff_lj*ic->rvdw);
2065 ic->sh_lj_ewald = (exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)*pow(ic->rvdw, minusSix);
2068 case eintmodFORCESWITCH:
2069 /* Switch the force, switch and shift the potential */
2070 force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw,
2071 &ic->dispersion_shift);
2072 force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw,
2073 &ic->repulsion_shift);
2075 case eintmodPOTSWITCH:
2076 /* Switch the potential and force */
2077 potential_switch_constants(ic->rvdw_switch, ic->rvdw,
2081 case eintmodEXACTCUTOFF:
2082 /* Nothing to do here */
2085 gmx_incons("unimplemented potential modifier");
2088 ic->sh_invrc6 = -ic->dispersion_shift.cpot;
2090 /* Electrostatics */
2091 ic->eeltype = fr->eeltype;
2092 ic->coulomb_modifier = fr->coulomb_modifier;
2093 ic->rcoulomb = fr->rcoulomb;
2094 ic->epsilon_r = fr->epsilon_r;
2095 ic->epsfac = fr->epsfac;
2096 ic->ewaldcoeff_q = fr->ewaldcoeff_q;
2098 if (fr->coulomb_modifier == eintmodPOTSHIFT)
2100 ic->sh_ewald = gmx_erfc(ic->ewaldcoeff_q*ic->rcoulomb);
2107 /* Reaction-field */
2108 if (EEL_RF(ic->eeltype))
2110 ic->epsilon_rf = fr->epsilon_rf;
2111 ic->k_rf = fr->k_rf;
2112 ic->c_rf = fr->c_rf;
2116 /* For plain cut-off we might use the reaction-field kernels */
2117 ic->epsilon_rf = ic->epsilon_r;
2119 if (fr->coulomb_modifier == eintmodPOTSHIFT)
2121 ic->c_rf = 1/ic->rcoulomb;
2131 real dispersion_shift;
2133 dispersion_shift = ic->dispersion_shift.cpot;
2134 if (EVDW_PME(ic->vdwtype))
2136 dispersion_shift -= ic->sh_lj_ewald;
2138 fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
2139 ic->repulsion_shift.cpot, dispersion_shift);
2141 if (ic->eeltype == eelCUT)
2143 fprintf(fp, ", Coulomb %.e", -ic->c_rf);
2145 else if (EEL_PME(ic->eeltype))
2147 fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
2152 *interaction_const = ic;
2155 static void init_nb_verlet(FILE *fp,
2156 nonbonded_verlet_t **nb_verlet,
2157 gmx_bool bFEP_NonBonded,
2158 const t_inputrec *ir,
2159 const t_forcerec *fr,
2160 const t_commrec *cr,
2161 const char *nbpu_opt)
2163 nonbonded_verlet_t *nbv;
2166 gmx_bool bEmulateGPU, bHybridGPURun = FALSE;
2168 nbnxn_alloc_t *nb_alloc;
2169 nbnxn_free_t *nb_free;
2173 pick_nbnxn_resources(fp, cr, fr->hwinfo,
2180 nbv->min_ci_balanced = 0;
2182 nbv->ngrp = (DOMAINDECOMP(cr) ? 2 : 1);
2183 for (i = 0; i < nbv->ngrp; i++)
2185 nbv->grp[i].nbl_lists.nnbl = 0;
2186 nbv->grp[i].nbat = NULL;
2187 nbv->grp[i].kernel_type = nbnxnkNotSet;
2189 if (i == 0) /* local */
2191 pick_nbnxn_kernel(fp, cr, fr->use_simd_kernels,
2192 nbv->bUseGPU, bEmulateGPU, ir,
2193 &nbv->grp[i].kernel_type,
2194 &nbv->grp[i].ewald_excl,
2197 else /* non-local */
2199 if (nbpu_opt != NULL && strcmp(nbpu_opt, "gpu_cpu") == 0)
2201 /* Use GPU for local, select a CPU kernel for non-local */
2202 pick_nbnxn_kernel(fp, cr, fr->use_simd_kernels,
2204 &nbv->grp[i].kernel_type,
2205 &nbv->grp[i].ewald_excl,
2208 bHybridGPURun = TRUE;
2212 /* Use the same kernel for local and non-local interactions */
2213 nbv->grp[i].kernel_type = nbv->grp[0].kernel_type;
2214 nbv->grp[i].ewald_excl = nbv->grp[0].ewald_excl;
2219 nbnxn_init_search(&nbv->nbs,
2220 DOMAINDECOMP(cr) ? &cr->dd->nc : NULL,
2221 DOMAINDECOMP(cr) ? domdec_zones(cr->dd) : NULL,
2223 gmx_omp_nthreads_get(emntPairsearch));
2225 for (i = 0; i < nbv->ngrp; i++)
2227 gpu_set_host_malloc_and_free(nbv->grp[0].kernel_type == nbnxnk8x8x8_GPU,
2228 &nb_alloc, &nb_free);
2230 nbnxn_init_pairlist_set(&nbv->grp[i].nbl_lists,
2231 nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2232 /* 8x8x8 "non-simple" lists are ATM always combined */
2233 !nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2237 nbv->grp[0].kernel_type != nbv->grp[i].kernel_type)
2239 gmx_bool bSimpleList;
2240 int enbnxninitcombrule;
2242 bSimpleList = nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type);
2244 if (bSimpleList && (fr->vdwtype == evdwCUT && (fr->vdw_modifier == eintmodNONE || fr->vdw_modifier == eintmodPOTSHIFT)))
2246 /* Plain LJ cut-off: we can optimize with combination rules */
2247 enbnxninitcombrule = enbnxninitcombruleDETECT;
2249 else if (fr->vdwtype == evdwPME)
2251 /* LJ-PME: we need to use a combination rule for the grid */
2252 if (fr->ljpme_combination_rule == eljpmeGEOM)
2254 enbnxninitcombrule = enbnxninitcombruleGEOM;
2258 enbnxninitcombrule = enbnxninitcombruleLB;
2263 /* We use a full combination matrix: no rule required */
2264 enbnxninitcombrule = enbnxninitcombruleNONE;
2268 snew(nbv->grp[i].nbat, 1);
2269 nbnxn_atomdata_init(fp,
2271 nbv->grp[i].kernel_type,
2273 fr->ntype, fr->nbfp,
2275 bSimpleList ? gmx_omp_nthreads_get(emntNonbonded) : 1,
2280 nbv->grp[i].nbat = nbv->grp[0].nbat;
2286 /* init the NxN GPU data; the last argument tells whether we'll have
2287 * both local and non-local NB calculation on GPU */
2288 nbnxn_gpu_init(fp, &nbv->gpu_nbv,
2289 &fr->hwinfo->gpu_info,
2293 cr->rank_pp_intranode,
2295 (nbv->ngrp > 1) && !bHybridGPURun);
2297 /* With tMPI + GPUs some ranks may be sharing GPU(s) and therefore
2298 * also sharing texture references. To keep the code simple, we don't
2299 * treat texture references as shared resources, but this means that
2300 * the coulomb_tab and nbfp texture refs will get updated by multiple threads.
2301 * Hence, to ensure that the non-bonded kernels don't start before all
2302 * texture binding operations are finished, we need to wait for all ranks
2303 * to arrive here before continuing.
2305 * Note that we could omit this barrier if GPUs are not shared (or
2306 * texture objects are used), but as this is initialization code, there
2307 * is no point in complicating things.
2309 #ifdef GMX_THREAD_MPI
2314 #endif /* GMX_THREAD_MPI */
2316 if ((env = getenv("GMX_NB_MIN_CI")) != NULL)
2320 nbv->min_ci_balanced = strtol(env, &end, 10);
2321 if (!end || (*end != 0) || nbv->min_ci_balanced <= 0)
2323 gmx_fatal(FARGS, "Invalid value passed in GMX_NB_MIN_CI=%s, positive integer required", env);
2328 fprintf(debug, "Neighbor-list balancing parameter: %d (passed as env. var.)\n",
2329 nbv->min_ci_balanced);
2334 nbv->min_ci_balanced = nbnxn_gpu_min_ci_balanced(nbv->gpu_nbv);
2337 fprintf(debug, "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU multi-processors)\n",
2338 nbv->min_ci_balanced);
2347 gmx_bool usingGpu(nonbonded_verlet_t *nbv)
2349 return nbv != NULL && nbv->bUseGPU;
2352 void init_forcerec(FILE *fp,
2353 const output_env_t oenv,
2356 const t_inputrec *ir,
2357 const gmx_mtop_t *mtop,
2358 const t_commrec *cr,
2363 const t_filenm *tabbfnm,
2364 const char *nbpu_opt,
2365 gmx_bool bNoSolvOpt,
2368 int i, m, negp_pp, negptable, egi, egj;
2373 gmx_bool bGenericKernelOnly;
2374 gmx_bool bMakeTables, bMakeSeparate14Table, bSomeNormalNbListsAreInUse;
2375 gmx_bool bFEP_NonBonded;
2376 int *nm_ind, egp_flags;
2378 if (fr->hwinfo == NULL)
2380 /* Detect hardware, gather information.
2381 * In mdrun, hwinfo has already been set before calling init_forcerec.
2382 * Here we ignore GPUs, as tools will not use them anyhow.
2384 fr->hwinfo = gmx_detect_hardware(fp, cr, FALSE);
2387 /* By default we turn SIMD kernels on, but it might be turned off further down... */
2388 fr->use_simd_kernels = TRUE;
2390 fr->bDomDec = DOMAINDECOMP(cr);
2392 if (check_box(ir->ePBC, box))
2394 gmx_fatal(FARGS, check_box(ir->ePBC, box));
2397 /* Test particle insertion ? */
2400 /* Set to the size of the molecule to be inserted (the last one) */
2401 /* Because of old style topologies, we have to use the last cg
2402 * instead of the last molecule type.
2404 cgs = &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].cgs;
2405 fr->n_tpi = cgs->index[cgs->nr] - cgs->index[cgs->nr-1];
2406 if (fr->n_tpi != mtop->mols.index[mtop->mols.nr] - mtop->mols.index[mtop->mols.nr-1])
2408 gmx_fatal(FARGS, "The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
2416 /* Copy AdResS parameters */
2419 fr->adress_type = ir->adress->type;
2420 fr->adress_const_wf = ir->adress->const_wf;
2421 fr->adress_ex_width = ir->adress->ex_width;
2422 fr->adress_hy_width = ir->adress->hy_width;
2423 fr->adress_icor = ir->adress->icor;
2424 fr->adress_site = ir->adress->site;
2425 fr->adress_ex_forcecap = ir->adress->ex_forcecap;
2426 fr->adress_do_hybridpairs = ir->adress->do_hybridpairs;
2429 snew(fr->adress_group_explicit, ir->adress->n_energy_grps);
2430 for (i = 0; i < ir->adress->n_energy_grps; i++)
2432 fr->adress_group_explicit[i] = ir->adress->group_explicit[i];
2435 fr->n_adress_tf_grps = ir->adress->n_tf_grps;
2436 snew(fr->adress_tf_table_index, fr->n_adress_tf_grps);
2437 for (i = 0; i < fr->n_adress_tf_grps; i++)
2439 fr->adress_tf_table_index[i] = ir->adress->tf_table_index[i];
2441 copy_rvec(ir->adress->refs, fr->adress_refs);
2445 fr->adress_type = eAdressOff;
2446 fr->adress_do_hybridpairs = FALSE;
2449 /* Copy the user determined parameters */
2450 fr->userint1 = ir->userint1;
2451 fr->userint2 = ir->userint2;
2452 fr->userint3 = ir->userint3;
2453 fr->userint4 = ir->userint4;
2454 fr->userreal1 = ir->userreal1;
2455 fr->userreal2 = ir->userreal2;
2456 fr->userreal3 = ir->userreal3;
2457 fr->userreal4 = ir->userreal4;
2460 fr->fc_stepsize = ir->fc_stepsize;
2463 fr->efep = ir->efep;
2464 fr->sc_alphavdw = ir->fepvals->sc_alpha;
2465 if (ir->fepvals->bScCoul)
2467 fr->sc_alphacoul = ir->fepvals->sc_alpha;
2468 fr->sc_sigma6_min = pow(ir->fepvals->sc_sigma_min, 6);
2472 fr->sc_alphacoul = 0;
2473 fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
2475 fr->sc_power = ir->fepvals->sc_power;
2476 fr->sc_r_power = ir->fepvals->sc_r_power;
2477 fr->sc_sigma6_def = pow(ir->fepvals->sc_sigma, 6);
2479 env = getenv("GMX_SCSIGMA_MIN");
2483 sscanf(env, "%20lf", &dbl);
2484 fr->sc_sigma6_min = pow(dbl, 6);
2487 fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
2491 fr->bNonbonded = TRUE;
2492 if (getenv("GMX_NO_NONBONDED") != NULL)
2494 /* turn off non-bonded calculations */
2495 fr->bNonbonded = FALSE;
2496 md_print_warn(cr, fp,
2497 "Found environment variable GMX_NO_NONBONDED.\n"
2498 "Disabling nonbonded calculations.\n");
2501 bGenericKernelOnly = FALSE;
2503 /* We now check in the NS code whether a particular combination of interactions
2504 * can be used with water optimization, and disable it if that is not the case.
2507 if (getenv("GMX_NB_GENERIC") != NULL)
2512 "Found environment variable GMX_NB_GENERIC.\n"
2513 "Disabling all interaction-specific nonbonded kernels, will only\n"
2514 "use the slow generic ones in src/gmxlib/nonbonded/nb_generic.c\n\n");
2516 bGenericKernelOnly = TRUE;
2519 if (bGenericKernelOnly == TRUE)
2524 if ( (getenv("GMX_DISABLE_SIMD_KERNELS") != NULL) || (getenv("GMX_NOOPTIMIZEDKERNELS") != NULL) )
2526 fr->use_simd_kernels = FALSE;
2530 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
2531 "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
2532 "(e.g. SSE2/SSE4.1/AVX).\n\n");
2536 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
2538 /* Check if we can/should do all-vs-all kernels */
2539 fr->bAllvsAll = can_use_allvsall(ir, FALSE, NULL, NULL);
2540 fr->AllvsAll_work = NULL;
2541 fr->AllvsAll_workgb = NULL;
2543 /* All-vs-all kernels have not been implemented in 4.6 and later.
2544 * See Redmine #1249. */
2547 fr->bAllvsAll = FALSE;
2551 "\nYour simulation settings would have triggered the efficient all-vs-all\n"
2552 "kernels in GROMACS 4.5, but these have not been implemented in GROMACS\n"
2553 "4.6 and 5.x. If performance is important, please use GROMACS 4.5.7\n"
2554 "or try cutoff-scheme = Verlet.\n\n");
2558 /* Neighbour searching stuff */
2559 fr->cutoff_scheme = ir->cutoff_scheme;
2560 fr->bGrid = (ir->ns_type == ensGRID);
2561 fr->ePBC = ir->ePBC;
2563 if (fr->cutoff_scheme == ecutsGROUP)
2565 const char *note = "NOTE: This file uses the deprecated 'group' cutoff_scheme. This will be\n"
2566 "removed in a future release when 'verlet' supports all interaction forms.\n";
2570 fprintf(stderr, "\n%s\n", note);
2574 fprintf(fp, "\n%s\n", note);
2578 /* Determine if we will do PBC for distances in bonded interactions */
2579 if (fr->ePBC == epbcNONE)
2581 fr->bMolPBC = FALSE;
2585 if (!DOMAINDECOMP(cr))
2589 bSHAKE = (ir->eConstrAlg == econtSHAKE &&
2590 (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0 ||
2591 gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0));
2593 /* The group cut-off scheme and SHAKE assume charge groups
2594 * are whole, but not using molpbc is faster in most cases.
2595 * With intermolecular interactions we need PBC for calculating
2596 * distances between atoms in different molecules.
2598 if ((fr->cutoff_scheme == ecutsGROUP || bSHAKE) &&
2599 !mtop->bIntermolecularInteractions)
2601 fr->bMolPBC = ir->bPeriodicMols;
2603 if (bSHAKE && fr->bMolPBC)
2605 gmx_fatal(FARGS, "SHAKE is not supported with periodic molecules");
2612 if (getenv("GMX_USE_GRAPH") != NULL)
2614 fr->bMolPBC = FALSE;
2617 md_print_warn(cr, fp, "GMX_USE_GRAPH is set, using the graph for bonded interactions\n");
2620 if (mtop->bIntermolecularInteractions)
2622 md_print_warn(cr, fp, "WARNING: Molecules linked by intermolecular interactions have to reside in the same periodic image, otherwise artifacts will occur!\n");
2626 if (bSHAKE && fr->bMolPBC)
2628 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");
2634 fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->ePBC);
2637 fr->bGB = (ir->implicit_solvent == eisGBSA);
2639 fr->rc_scaling = ir->refcoord_scaling;
2640 copy_rvec(ir->posres_com, fr->posres_com);
2641 copy_rvec(ir->posres_comB, fr->posres_comB);
2642 fr->rlist = cutoff_inf(ir->rlist);
2643 fr->rlistlong = cutoff_inf(ir->rlistlong);
2644 fr->eeltype = ir->coulombtype;
2645 fr->vdwtype = ir->vdwtype;
2646 fr->ljpme_combination_rule = ir->ljpme_combination_rule;
2648 fr->coulomb_modifier = ir->coulomb_modifier;
2649 fr->vdw_modifier = ir->vdw_modifier;
2651 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
2652 switch (fr->eeltype)
2655 fr->nbkernel_elec_interaction = (fr->bGB) ? GMX_NBKERNEL_ELEC_GENERALIZEDBORN : GMX_NBKERNEL_ELEC_COULOMB;
2661 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2665 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2666 fr->coulomb_modifier = eintmodEXACTCUTOFF;
2675 case eelPMEUSERSWITCH:
2676 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2682 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD;
2686 gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[fr->eeltype]);
2690 /* Vdw: Translate from mdp settings to kernel format */
2691 switch (fr->vdwtype)
2696 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
2700 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
2704 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD;
2710 case evdwENCADSHIFT:
2711 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2715 gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[fr->vdwtype]);
2719 /* These start out identical to ir, but might be altered if we e.g. tabulate the interaction in the kernel */
2720 fr->nbkernel_elec_modifier = fr->coulomb_modifier;
2721 fr->nbkernel_vdw_modifier = fr->vdw_modifier;
2723 fr->rvdw = cutoff_inf(ir->rvdw);
2724 fr->rvdw_switch = ir->rvdw_switch;
2725 fr->rcoulomb = cutoff_inf(ir->rcoulomb);
2726 fr->rcoulomb_switch = ir->rcoulomb_switch;
2728 fr->bTwinRange = fr->rlistlong > fr->rlist;
2729 fr->bEwald = (EEL_PME(fr->eeltype) || fr->eeltype == eelEWALD);
2731 fr->reppow = mtop->ffparams.reppow;
2733 if (ir->cutoff_scheme == ecutsGROUP)
2735 fr->bvdwtab = ((fr->vdwtype != evdwCUT || !gmx_within_tol(fr->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2736 && !EVDW_PME(fr->vdwtype));
2737 /* We have special kernels for standard Ewald and PME, but the pme-switch ones are tabulated above */
2738 fr->bcoultab = !(fr->eeltype == eelCUT ||
2739 fr->eeltype == eelEWALD ||
2740 fr->eeltype == eelPME ||
2741 fr->eeltype == eelRF ||
2742 fr->eeltype == eelRF_ZERO);
2744 /* If the user absolutely wants different switch/shift settings for coul/vdw, it is likely
2745 * going to be faster to tabulate the interaction than calling the generic kernel.
2746 * However, if generic kernels have been requested we keep things analytically.
2748 if (fr->nbkernel_elec_modifier == eintmodPOTSWITCH &&
2749 fr->nbkernel_vdw_modifier == eintmodPOTSWITCH &&
2750 bGenericKernelOnly == FALSE)
2752 if ((fr->rcoulomb_switch != fr->rvdw_switch) || (fr->rcoulomb != fr->rvdw))
2754 fr->bcoultab = TRUE;
2755 /* Once we tabulate electrostatics, we can use the switch function for LJ,
2756 * which would otherwise need two tables.
2760 else if ((fr->nbkernel_elec_modifier == eintmodPOTSHIFT && fr->nbkernel_vdw_modifier == eintmodPOTSHIFT) ||
2761 ((fr->nbkernel_elec_interaction == GMX_NBKERNEL_ELEC_REACTIONFIELD &&
2762 fr->nbkernel_elec_modifier == eintmodEXACTCUTOFF &&
2763 (fr->nbkernel_vdw_modifier == eintmodPOTSWITCH || fr->nbkernel_vdw_modifier == eintmodPOTSHIFT))))
2765 if ((fr->rcoulomb != fr->rvdw) && (bGenericKernelOnly == FALSE))
2767 fr->bcoultab = TRUE;
2771 if (fr->nbkernel_elec_modifier == eintmodFORCESWITCH)
2773 fr->bcoultab = TRUE;
2775 if (fr->nbkernel_vdw_modifier == eintmodFORCESWITCH)
2780 if (getenv("GMX_REQUIRE_TABLES"))
2783 fr->bcoultab = TRUE;
2788 fprintf(fp, "Table routines are used for coulomb: %s\n", bool_names[fr->bcoultab]);
2789 fprintf(fp, "Table routines are used for vdw: %s\n", bool_names[fr->bvdwtab ]);
2792 if (fr->bvdwtab == TRUE)
2794 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2795 fr->nbkernel_vdw_modifier = eintmodNONE;
2797 if (fr->bcoultab == TRUE)
2799 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2800 fr->nbkernel_elec_modifier = eintmodNONE;
2804 if (ir->cutoff_scheme == ecutsVERLET)
2806 if (!gmx_within_tol(fr->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2808 gmx_fatal(FARGS, "Cut-off scheme %S only supports LJ repulsion power 12", ecutscheme_names[ir->cutoff_scheme]);
2810 fr->bvdwtab = FALSE;
2811 fr->bcoultab = FALSE;
2814 /* Tables are used for direct ewald sum */
2817 if (EEL_PME(ir->coulombtype))
2821 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
2823 if (ir->coulombtype == eelP3M_AD)
2825 please_cite(fp, "Hockney1988");
2826 please_cite(fp, "Ballenegger2012");
2830 please_cite(fp, "Essmann95a");
2833 if (ir->ewald_geometry == eewg3DC)
2837 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry.\n");
2839 please_cite(fp, "In-Chul99a");
2842 fr->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
2843 init_ewald_tab(&(fr->ewald_table), ir, fp);
2846 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n",
2847 1/fr->ewaldcoeff_q);
2851 if (EVDW_PME(ir->vdwtype))
2855 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
2857 please_cite(fp, "Essmann95a");
2858 fr->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
2861 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n",
2862 1/fr->ewaldcoeff_lj);
2866 /* Electrostatics */
2867 fr->epsilon_r = ir->epsilon_r;
2868 fr->epsilon_rf = ir->epsilon_rf;
2869 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
2871 /* Parameters for generalized RF */
2875 if (fr->eeltype == eelGRF)
2877 init_generalized_rf(fp, mtop, ir, fr);
2880 fr->bF_NoVirSum = (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype) ||
2881 gmx_mtop_ftype_count(mtop, F_POSRES) > 0 ||
2882 gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 ||
2883 IR_ELEC_FIELD(*ir) ||
2884 (fr->adress_icor != eAdressICOff)
2887 if (fr->cutoff_scheme == ecutsGROUP &&
2888 ncg_mtop(mtop) > fr->cg_nalloc && !DOMAINDECOMP(cr))
2890 /* Count the total number of charge groups */
2891 fr->cg_nalloc = ncg_mtop(mtop);
2892 srenew(fr->cg_cm, fr->cg_nalloc);
2894 if (fr->shift_vec == NULL)
2896 snew(fr->shift_vec, SHIFTS);
2899 if (fr->fshift == NULL)
2901 snew(fr->fshift, SHIFTS);
2904 if (fr->nbfp == NULL)
2906 fr->ntype = mtop->ffparams.atnr;
2907 fr->nbfp = mk_nbfp(&mtop->ffparams, fr->bBHAM);
2908 if (EVDW_PME(fr->vdwtype))
2910 fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
2914 /* Copy the energy group exclusions */
2915 fr->egp_flags = ir->opts.egp_flags;
2917 /* Van der Waals stuff */
2918 if ((fr->vdwtype != evdwCUT) && (fr->vdwtype != evdwUSER) && !fr->bBHAM)
2920 if (fr->rvdw_switch >= fr->rvdw)
2922 gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)",
2923 fr->rvdw_switch, fr->rvdw);
2927 fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
2928 (fr->eeltype == eelSWITCH) ? "switched" : "shifted",
2929 fr->rvdw_switch, fr->rvdw);
2933 if (fr->bBHAM && EVDW_PME(fr->vdwtype))
2935 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
2938 if (fr->bBHAM && (fr->vdwtype == evdwSHIFT || fr->vdwtype == evdwSWITCH))
2940 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
2943 if (fr->bBHAM && fr->cutoff_scheme == ecutsVERLET)
2945 gmx_fatal(FARGS, "Verlet cutoff-scheme is not supported with Buckingham");
2950 fprintf(fp, "Cut-off's: NS: %g Coulomb: %g %s: %g\n",
2951 fr->rlist, fr->rcoulomb, fr->bBHAM ? "BHAM" : "LJ", fr->rvdw);
2954 fr->eDispCorr = ir->eDispCorr;
2955 if (ir->eDispCorr != edispcNO)
2957 set_avcsixtwelve(fp, fr, mtop);
2962 set_bham_b_max(fp, fr, mtop);
2965 fr->gb_epsilon_solvent = ir->gb_epsilon_solvent;
2967 /* Copy the GBSA data (radius, volume and surftens for each
2968 * atomtype) from the topology atomtype section to forcerec.
2970 snew(fr->atype_radius, fr->ntype);
2971 snew(fr->atype_vol, fr->ntype);
2972 snew(fr->atype_surftens, fr->ntype);
2973 snew(fr->atype_gb_radius, fr->ntype);
2974 snew(fr->atype_S_hct, fr->ntype);
2976 if (mtop->atomtypes.nr > 0)
2978 for (i = 0; i < fr->ntype; i++)
2980 fr->atype_radius[i] = mtop->atomtypes.radius[i];
2982 for (i = 0; i < fr->ntype; i++)
2984 fr->atype_vol[i] = mtop->atomtypes.vol[i];
2986 for (i = 0; i < fr->ntype; i++)
2988 fr->atype_surftens[i] = mtop->atomtypes.surftens[i];
2990 for (i = 0; i < fr->ntype; i++)
2992 fr->atype_gb_radius[i] = mtop->atomtypes.gb_radius[i];
2994 for (i = 0; i < fr->ntype; i++)
2996 fr->atype_S_hct[i] = mtop->atomtypes.S_hct[i];
3000 /* Generate the GB table if needed */
3004 fr->gbtabscale = 2000;
3006 fr->gbtabscale = 500;
3010 fr->gbtab = make_gb_table(oenv, fr);
3012 init_gb(&fr->born, fr, ir, mtop, ir->gb_algorithm);
3014 /* Copy local gb data (for dd, this is done in dd_partition_system) */
3015 if (!DOMAINDECOMP(cr))
3017 make_local_gb(cr, fr->born, ir->gb_algorithm);
3021 /* Set the charge scaling */
3022 if (fr->epsilon_r != 0)
3024 fr->epsfac = ONE_4PI_EPS0/fr->epsilon_r;
3028 /* eps = 0 is infinite dieletric: no coulomb interactions */
3032 /* Reaction field constants */
3033 if (EEL_RF(fr->eeltype))
3035 calc_rffac(fp, fr->eeltype, fr->epsilon_r, fr->epsilon_rf,
3036 fr->rcoulomb, fr->temp, fr->zsquare, box,
3037 &fr->kappa, &fr->k_rf, &fr->c_rf);
3040 /*This now calculates sum for q and c6*/
3041 set_chargesum(fp, fr, mtop);
3043 /* if we are using LR electrostatics, and they are tabulated,
3044 * the tables will contain modified coulomb interactions.
3045 * Since we want to use the non-shifted ones for 1-4
3046 * coulombic interactions, we must have an extra set of tables.
3049 /* Construct tables.
3050 * A little unnecessary to make both vdw and coul tables sometimes,
3051 * but what the heck... */
3053 bMakeTables = fr->bcoultab || fr->bvdwtab || fr->bEwald ||
3054 (ir->eDispCorr != edispcNO && ir_vdw_switched(ir));
3056 bMakeSeparate14Table = ((!bMakeTables || fr->eeltype != eelCUT || fr->vdwtype != evdwCUT ||
3057 fr->coulomb_modifier != eintmodNONE ||
3058 fr->vdw_modifier != eintmodNONE ||
3059 fr->bBHAM || fr->bEwald) &&
3060 (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
3061 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0 ||
3062 gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0));
3064 negp_pp = ir->opts.ngener - ir->nwall;
3068 bSomeNormalNbListsAreInUse = TRUE;
3073 bSomeNormalNbListsAreInUse = (ir->eDispCorr != edispcNO);
3074 for (egi = 0; egi < negp_pp; egi++)
3076 for (egj = egi; egj < negp_pp; egj++)
3078 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
3079 if (!(egp_flags & EGP_EXCL))
3081 if (egp_flags & EGP_TABLE)
3087 bSomeNormalNbListsAreInUse = TRUE;
3092 if (bSomeNormalNbListsAreInUse)
3094 fr->nnblists = negptable + 1;
3098 fr->nnblists = negptable;
3100 if (fr->nnblists > 1)
3102 snew(fr->gid2nblists, ir->opts.ngener*ir->opts.ngener);
3111 snew(fr->nblists, fr->nnblists);
3113 /* This code automatically gives table length tabext without cut-off's,
3114 * in that case grompp should already have checked that we do not need
3115 * normal tables and we only generate tables for 1-4 interactions.
3117 rtab = ir->rlistlong + ir->tabext;
3121 /* make tables for ordinary interactions */
3122 if (bSomeNormalNbListsAreInUse)
3124 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn, NULL, NULL, &fr->nblists[0]);
3127 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn, NULL, NULL, &fr->nblists[fr->nnblists/2]);
3129 if (!bMakeSeparate14Table)
3131 fr->tab14 = fr->nblists[0].table_elec_vdw;
3141 /* Read the special tables for certain energy group pairs */
3142 nm_ind = mtop->groups.grps[egcENER].nm_ind;
3143 for (egi = 0; egi < negp_pp; egi++)
3145 for (egj = egi; egj < negp_pp; egj++)
3147 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
3148 if ((egp_flags & EGP_TABLE) && !(egp_flags & EGP_EXCL))
3150 if (fr->nnblists > 1)
3152 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = m;
3154 /* Read the table file with the two energy groups names appended */
3155 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn,
3156 *mtop->groups.grpname[nm_ind[egi]],
3157 *mtop->groups.grpname[nm_ind[egj]],
3161 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn,
3162 *mtop->groups.grpname[nm_ind[egi]],
3163 *mtop->groups.grpname[nm_ind[egj]],
3164 &fr->nblists[fr->nnblists/2+m]);
3168 else if (fr->nnblists > 1)
3170 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = 0;
3176 else if ((fr->eDispCorr != edispcNO) &&
3177 ((fr->vdw_modifier == eintmodPOTSWITCH) ||
3178 (fr->vdw_modifier == eintmodFORCESWITCH) ||
3179 (fr->vdw_modifier == eintmodPOTSHIFT)))
3181 /* Tables might not be used for the potential modifier interactions per se, but
3182 * we still need them to evaluate switch/shift dispersion corrections in this case.
3184 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn, NULL, NULL, &fr->nblists[0]);
3187 if (bMakeSeparate14Table)
3189 /* generate extra tables with plain Coulomb for 1-4 interactions only */
3190 fr->tab14 = make_tables(fp, oenv, fr, MASTER(cr), tabpfn, rtab,
3191 GMX_MAKETABLES_14ONLY);
3194 /* Read AdResS Thermo Force table if needed */
3195 if (fr->adress_icor == eAdressICThermoForce)
3197 /* old todo replace */
3199 if (ir->adress->n_tf_grps > 0)
3201 make_adress_tf_tables(fp, oenv, fr, ir, tabfn, mtop, box);
3206 /* load the default table */
3207 snew(fr->atf_tabs, 1);
3208 fr->atf_tabs[DEFAULT_TF_TABLE] = make_atf_table(fp, oenv, fr, tabafn, box);
3213 fr->nwall = ir->nwall;
3214 if (ir->nwall && ir->wall_type == ewtTABLE)
3216 make_wall_tables(fp, oenv, ir, tabfn, &mtop->groups, fr);
3221 // Need to catch std::bad_alloc
3222 // TODO Don't need to catch this here, when merging with master branch
3225 fcd->bondtab = make_bonded_tables(fp,
3226 F_TABBONDS, F_TABBONDSNC,
3227 mtop, tabbfnm, "b");
3228 fcd->angletab = make_bonded_tables(fp,
3230 mtop, tabbfnm, "a");
3231 fcd->dihtab = make_bonded_tables(fp,
3233 mtop, tabbfnm, "d");
3235 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
3241 fprintf(debug, "No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
3245 /* QM/MM initialization if requested
3249 fprintf(stderr, "QM/MM calculation requested.\n");
3252 fr->bQMMM = ir->bQMMM;
3253 fr->qr = mk_QMMMrec();
3255 /* Set all the static charge group info */
3256 fr->cginfo_mb = init_cginfo_mb(fp, mtop, fr, bNoSolvOpt,
3258 &fr->bExcl_IntraCGAll_InterCGNone);
3259 if (DOMAINDECOMP(cr))
3265 fr->cginfo = cginfo_expand(mtop->nmolblock, fr->cginfo_mb);
3268 if (!DOMAINDECOMP(cr))
3270 forcerec_set_ranges(fr, ncg_mtop(mtop), ncg_mtop(mtop),
3271 mtop->natoms, mtop->natoms, mtop->natoms);
3274 fr->print_force = print_force;
3277 /* coarse load balancing vars */
3282 /* Initialize neighbor search */
3283 init_ns(fp, cr, &fr->ns, fr, mtop);
3285 if (cr->duty & DUTY_PP)
3287 gmx_nonbonded_setup(fr, bGenericKernelOnly);
3291 gmx_setup_adress_kernels(fp,bGenericKernelOnly);
3296 /* Initialize the thread working data for bonded interactions */
3297 init_bonded_threading(fp, fr, mtop->groups.grps[egcENER].nr);
3299 /* fr->ic is used both by verlet and group kernels (to some extent) now */
3300 init_interaction_const(fp, &fr->ic, fr);
3301 init_interaction_const_tables(fp, fr->ic, rtab);
3303 if (fr->cutoff_scheme == ecutsVERLET)
3305 if (ir->rcoulomb != ir->rvdw)
3307 gmx_fatal(FARGS, "With Verlet lists rcoulomb and rvdw should be identical");
3310 init_nb_verlet(fp, &fr->nbv, bFEP_NonBonded, ir, fr, cr, nbpu_opt);
3313 if (ir->eDispCorr != edispcNO)
3315 calc_enervirdiff(fp, ir->eDispCorr, fr);
3319 #define pr_real(fp, r) fprintf(fp, "%s: %e\n",#r, r)
3320 #define pr_int(fp, i) fprintf((fp), "%s: %d\n",#i, i)
3321 #define pr_bool(fp, b) fprintf((fp), "%s: %s\n",#b, bool_names[b])
3323 void pr_forcerec(FILE *fp, t_forcerec *fr)
3327 pr_real(fp, fr->rlist);
3328 pr_real(fp, fr->rcoulomb);
3329 pr_real(fp, fr->fudgeQQ);
3330 pr_bool(fp, fr->bGrid);
3331 pr_bool(fp, fr->bTwinRange);
3332 /*pr_int(fp,fr->cg0);
3333 pr_int(fp,fr->hcg);*/
3334 for (i = 0; i < fr->nnblists; i++)
3336 pr_int(fp, fr->nblists[i].table_elec_vdw.n);
3338 pr_real(fp, fr->rcoulomb_switch);
3339 pr_real(fp, fr->rcoulomb);
3344 /* Frees GPU memory and destroys the GPU context.
3346 * Note that this function needs to be called even if GPUs are not used
3347 * in this run because the PME ranks have no knowledge of whether GPUs
3348 * are used or not, but all ranks need to enter the barrier below.
3350 void free_gpu_resources(const t_forcerec *fr,
3351 const t_commrec *cr,
3352 const gmx_gpu_info_t *gpu_info,
3353 const gmx_gpu_opt_t *gpu_opt)
3355 gmx_bool bIsPPrankUsingGPU;
3356 char gpu_err_str[STRLEN];
3358 bIsPPrankUsingGPU = (cr->duty & DUTY_PP) && fr && fr->nbv && fr->nbv->bUseGPU;
3360 if (bIsPPrankUsingGPU)
3362 /* free nbnxn data in GPU memory */
3363 nbnxn_gpu_free(fr->nbv->gpu_nbv);
3365 /* With tMPI we need to wait for all ranks to finish deallocation before
3366 * destroying the context in free_gpu() as some ranks may be sharing
3368 * Note: as only PP ranks need to free GPU resources, so it is safe to
3369 * not call the barrier on PME ranks.
3371 #ifdef GMX_THREAD_MPI
3376 #endif /* GMX_THREAD_MPI */
3378 /* uninitialize GPU (by destroying the context) */
3379 if (!free_cuda_gpu(cr->rank_pp_intranode, gpu_err_str, gpu_info, gpu_opt))
3381 gmx_warning("On rank %d failed to free GPU #%d: %s",
3382 cr->nodeid, get_current_cuda_gpu_device_id(), gpu_err_str);