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, 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.
45 #include "types/commrec.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/math/utilities.h"
49 #include "gromacs/utility/smalloc.h"
51 #include "gromacs/utility/fatalerror.h"
55 #include "nonbonded.h"
64 #include "md_support.h"
65 #include "md_logging.h"
69 #include "mtop_util.h"
70 #include "nbnxn_simd.h"
71 #include "nbnxn_search.h"
72 #include "nbnxn_atomdata.h"
73 #include "nbnxn_consts.h"
74 #include "gmx_omp_nthreads.h"
75 #include "gmx_detect_hardware.h"
78 #include "types/nbnxn_cuda_types_ext.h"
79 #include "gpu_utils.h"
80 #include "nbnxn_cuda_data_mgmt.h"
81 #include "pmalloc_cuda.h"
83 t_forcerec *mk_forcerec(void)
93 static void pr_nbfp(FILE *fp, real *nbfp, gmx_bool bBHAM, int atnr)
97 for (i = 0; (i < atnr); i++)
99 for (j = 0; (j < atnr); j++)
101 fprintf(fp, "%2d - %2d", i, j);
104 fprintf(fp, " a=%10g, b=%10g, c=%10g\n", BHAMA(nbfp, atnr, i, j),
105 BHAMB(nbfp, atnr, i, j), BHAMC(nbfp, atnr, i, j)/6.0);
109 fprintf(fp, " c6=%10g, c12=%10g\n", C6(nbfp, atnr, i, j)/6.0,
110 C12(nbfp, atnr, i, j)/12.0);
117 static real *mk_nbfp(const gmx_ffparams_t *idef, gmx_bool bBHAM)
125 snew(nbfp, 3*atnr*atnr);
126 for (i = k = 0; (i < atnr); i++)
128 for (j = 0; (j < atnr); j++, k++)
130 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
131 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
132 /* nbfp now includes the 6.0 derivative prefactor */
133 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c*6.0;
139 snew(nbfp, 2*atnr*atnr);
140 for (i = k = 0; (i < atnr); i++)
142 for (j = 0; (j < atnr); j++, k++)
144 /* nbfp now includes the 6.0/12.0 derivative prefactors */
145 C6(nbfp, atnr, i, j) = idef->iparams[k].lj.c6*6.0;
146 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12*12.0;
154 static real *make_ljpme_c6grid(const gmx_ffparams_t *idef, t_forcerec *fr)
157 real c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
160 /* For LJ-PME simulations, we correct the energies with the reciprocal space
161 * inside of the cut-off. To do this the non-bonded kernels needs to have
162 * access to the C6-values used on the reciprocal grid in pme.c
166 snew(grid, 2*atnr*atnr);
167 for (i = k = 0; (i < atnr); i++)
169 for (j = 0; (j < atnr); j++, k++)
171 c6i = idef->iparams[i*(atnr+1)].lj.c6;
172 c12i = idef->iparams[i*(atnr+1)].lj.c12;
173 c6j = idef->iparams[j*(atnr+1)].lj.c6;
174 c12j = idef->iparams[j*(atnr+1)].lj.c12;
175 c6 = sqrt(c6i * c6j);
176 if (fr->ljpme_combination_rule == eljpmeLB
177 && !gmx_numzero(c6) && !gmx_numzero(c12i) && !gmx_numzero(c12j))
179 sigmai = pow(c12i / c6i, 1.0/6.0);
180 sigmaj = pow(c12j / c6j, 1.0/6.0);
181 epsi = c6i * c6i / c12i;
182 epsj = c6j * c6j / c12j;
183 c6 = sqrt(epsi * epsj) * pow(0.5*(sigmai+sigmaj), 6);
185 /* Store the elements at the same relative positions as C6 in nbfp in order
186 * to simplify access in the kernels
188 grid[2*(atnr*i+j)] = c6*6.0;
194 static real *mk_nbfp_combination_rule(const gmx_ffparams_t *idef, int comb_rule)
198 real c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
202 snew(nbfp, 2*atnr*atnr);
203 for (i = 0; i < atnr; ++i)
205 for (j = 0; j < atnr; ++j)
207 c6i = idef->iparams[i*(atnr+1)].lj.c6;
208 c12i = idef->iparams[i*(atnr+1)].lj.c12;
209 c6j = idef->iparams[j*(atnr+1)].lj.c6;
210 c12j = idef->iparams[j*(atnr+1)].lj.c12;
211 c6 = sqrt(c6i * c6j);
212 c12 = sqrt(c12i * c12j);
213 if (comb_rule == eCOMB_ARITHMETIC
214 && !gmx_numzero(c6) && !gmx_numzero(c12))
216 sigmai = pow(c12i / c6i, 1.0/6.0);
217 sigmaj = pow(c12j / c6j, 1.0/6.0);
218 epsi = c6i * c6i / c12i;
219 epsj = c6j * c6j / c12j;
220 c6 = epsi * epsj * pow(0.5*(sigmai+sigmaj), 6);
221 c12 = epsi * epsj * pow(0.5*(sigmai+sigmaj), 12);
223 C6(nbfp, atnr, i, j) = c6*6.0;
224 C12(nbfp, atnr, i, j) = c12*12.0;
230 /* This routine sets fr->solvent_opt to the most common solvent in the
231 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
232 * the fr->solvent_type array with the correct type (or esolNO).
234 * Charge groups that fulfill the conditions but are not identical to the
235 * most common one will be marked as esolNO in the solvent_type array.
237 * TIP3p is identical to SPC for these purposes, so we call it
238 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
240 * NOTE: QM particle should not
241 * become an optimized solvent. Not even if there is only one charge
251 } solvent_parameters_t;
254 check_solvent_cg(const gmx_moltype_t *molt,
257 const unsigned char *qm_grpnr,
258 const t_grps *qm_grps,
260 int *n_solvent_parameters,
261 solvent_parameters_t **solvent_parameters_p,
265 const t_blocka *excl;
272 real tmp_charge[4] = { 0.0 }; /* init to zero to make gcc4.8 happy */
273 int tmp_vdwtype[4] = { 0 }; /* init to zero to make gcc4.8 happy */
276 solvent_parameters_t *solvent_parameters;
278 /* We use a list with parameters for each solvent type.
279 * Every time we discover a new molecule that fulfills the basic
280 * conditions for a solvent we compare with the previous entries
281 * in these lists. If the parameters are the same we just increment
282 * the counter for that type, and otherwise we create a new type
283 * based on the current molecule.
285 * Once we've finished going through all molecules we check which
286 * solvent is most common, and mark all those molecules while we
287 * clear the flag on all others.
290 solvent_parameters = *solvent_parameters_p;
292 /* Mark the cg first as non optimized */
295 /* Check if this cg has no exclusions with atoms in other charge groups
296 * and all atoms inside the charge group excluded.
297 * We only have 3 or 4 atom solvent loops.
299 if (GET_CGINFO_EXCL_INTER(cginfo) ||
300 !GET_CGINFO_EXCL_INTRA(cginfo))
305 /* Get the indices of the first atom in this charge group */
306 j0 = molt->cgs.index[cg0];
307 j1 = molt->cgs.index[cg0+1];
309 /* Number of atoms in our molecule */
315 "Moltype '%s': there are %d atoms in this charge group\n",
319 /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
322 if (nj < 3 || nj > 4)
327 /* Check if we are doing QM on this group */
329 if (qm_grpnr != NULL)
331 for (j = j0; j < j1 && !qm; j++)
333 qm = (qm_grpnr[j] < qm_grps->nr - 1);
336 /* Cannot use solvent optimization with QM */
342 atom = molt->atoms.atom;
344 /* Still looks like a solvent, time to check parameters */
346 /* If it is perturbed (free energy) we can't use the solvent loops,
347 * so then we just skip to the next molecule.
351 for (j = j0; j < j1 && !perturbed; j++)
353 perturbed = PERTURBED(atom[j]);
361 /* Now it's only a question if the VdW and charge parameters
362 * are OK. Before doing the check we compare and see if they are
363 * identical to a possible previous solvent type.
364 * First we assign the current types and charges.
366 for (j = 0; j < nj; j++)
368 tmp_vdwtype[j] = atom[j0+j].type;
369 tmp_charge[j] = atom[j0+j].q;
372 /* Does it match any previous solvent type? */
373 for (k = 0; k < *n_solvent_parameters; k++)
378 /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
379 if ( (solvent_parameters[k].model == esolSPC && nj != 3) ||
380 (solvent_parameters[k].model == esolTIP4P && nj != 4) )
385 /* Check that types & charges match for all atoms in molecule */
386 for (j = 0; j < nj && match == TRUE; j++)
388 if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
392 if (tmp_charge[j] != solvent_parameters[k].charge[j])
399 /* Congratulations! We have a matched solvent.
400 * Flag it with this type for later processing.
403 solvent_parameters[k].count += nmol;
405 /* We are done with this charge group */
410 /* If we get here, we have a tentative new solvent type.
411 * Before we add it we must check that it fulfills the requirements
412 * of the solvent optimized loops. First determine which atoms have
415 for (j = 0; j < nj; j++)
418 tjA = tmp_vdwtype[j];
420 /* Go through all other tpes and see if any have non-zero
421 * VdW parameters when combined with this one.
423 for (k = 0; k < fr->ntype && (has_vdw[j] == FALSE); k++)
425 /* We already checked that the atoms weren't perturbed,
426 * so we only need to check state A now.
430 has_vdw[j] = (has_vdw[j] ||
431 (BHAMA(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
432 (BHAMB(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
433 (BHAMC(fr->nbfp, fr->ntype, tjA, k) != 0.0));
438 has_vdw[j] = (has_vdw[j] ||
439 (C6(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
440 (C12(fr->nbfp, fr->ntype, tjA, k) != 0.0));
445 /* Now we know all we need to make the final check and assignment. */
449 * For this we require thatn all atoms have charge,
450 * the charges on atom 2 & 3 should be the same, and only
451 * atom 1 might have VdW.
453 if (has_vdw[1] == FALSE &&
454 has_vdw[2] == FALSE &&
455 tmp_charge[0] != 0 &&
456 tmp_charge[1] != 0 &&
457 tmp_charge[2] == tmp_charge[1])
459 srenew(solvent_parameters, *n_solvent_parameters+1);
460 solvent_parameters[*n_solvent_parameters].model = esolSPC;
461 solvent_parameters[*n_solvent_parameters].count = nmol;
462 for (k = 0; k < 3; k++)
464 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
465 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
468 *cg_sp = *n_solvent_parameters;
469 (*n_solvent_parameters)++;
474 /* Or could it be a TIP4P?
475 * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
476 * Only atom 1 mght have VdW.
478 if (has_vdw[1] == FALSE &&
479 has_vdw[2] == FALSE &&
480 has_vdw[3] == FALSE &&
481 tmp_charge[0] == 0 &&
482 tmp_charge[1] != 0 &&
483 tmp_charge[2] == tmp_charge[1] &&
486 srenew(solvent_parameters, *n_solvent_parameters+1);
487 solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
488 solvent_parameters[*n_solvent_parameters].count = nmol;
489 for (k = 0; k < 4; k++)
491 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
492 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
495 *cg_sp = *n_solvent_parameters;
496 (*n_solvent_parameters)++;
500 *solvent_parameters_p = solvent_parameters;
504 check_solvent(FILE * fp,
505 const gmx_mtop_t * mtop,
507 cginfo_mb_t *cginfo_mb)
510 const t_block * mols;
511 const gmx_moltype_t *molt;
512 int mb, mol, cg_mol, at_offset, cg_offset, am, cgm, i, nmol_ch, nmol;
513 int n_solvent_parameters;
514 solvent_parameters_t *solvent_parameters;
520 fprintf(debug, "Going to determine what solvent types we have.\n");
525 n_solvent_parameters = 0;
526 solvent_parameters = NULL;
527 /* Allocate temporary array for solvent type */
528 snew(cg_sp, mtop->nmolblock);
532 for (mb = 0; mb < mtop->nmolblock; mb++)
534 molt = &mtop->moltype[mtop->molblock[mb].type];
536 /* Here we have to loop over all individual molecules
537 * because we need to check for QMMM particles.
539 snew(cg_sp[mb], cginfo_mb[mb].cg_mod);
540 nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
541 nmol = mtop->molblock[mb].nmol/nmol_ch;
542 for (mol = 0; mol < nmol_ch; mol++)
545 am = mol*cgs->index[cgs->nr];
546 for (cg_mol = 0; cg_mol < cgs->nr; cg_mol++)
548 check_solvent_cg(molt, cg_mol, nmol,
549 mtop->groups.grpnr[egcQMMM] ?
550 mtop->groups.grpnr[egcQMMM]+at_offset+am : 0,
551 &mtop->groups.grps[egcQMMM],
553 &n_solvent_parameters, &solvent_parameters,
554 cginfo_mb[mb].cginfo[cgm+cg_mol],
555 &cg_sp[mb][cgm+cg_mol]);
558 cg_offset += cgs->nr;
559 at_offset += cgs->index[cgs->nr];
562 /* Puh! We finished going through all charge groups.
563 * Now find the most common solvent model.
566 /* Most common solvent this far */
568 for (i = 0; i < n_solvent_parameters; i++)
571 solvent_parameters[i].count > solvent_parameters[bestsp].count)
579 bestsol = solvent_parameters[bestsp].model;
586 #ifdef DISABLE_WATER_NLIST
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, cgm, am;
639 int mb, m, ncg_tot, cg, a0, a1, gid, ai, j, aj, excl_nalloc;
643 gmx_bool bId, *bExcl, bExclIntraAll, bExclInter, bHaveVDW, bHaveQ, bHavePerturbedAtoms;
645 ncg_tot = ncg_mtop(mtop);
646 snew(cginfo_mb, mtop->nmolblock);
648 snew(type_VDW, fr->ntype);
649 for (ai = 0; ai < fr->ntype; ai++)
651 type_VDW[ai] = FALSE;
652 for (j = 0; j < fr->ntype; j++)
654 type_VDW[ai] = type_VDW[ai] ||
656 C6(fr->nbfp, fr->ntype, ai, j) != 0 ||
657 C12(fr->nbfp, fr->ntype, ai, j) != 0;
661 *bFEP_NonBonded = FALSE;
662 *bExcl_IntraCGAll_InterCGNone = TRUE;
665 snew(bExcl, excl_nalloc);
668 for (mb = 0; mb < mtop->nmolblock; mb++)
670 molb = &mtop->molblock[mb];
671 molt = &mtop->moltype[molb->type];
675 /* Check if the cginfo is identical for all molecules in this block.
676 * If so, we only need an array of the size of one molecule.
677 * Otherwise we make an array of #mol times #cgs per molecule.
681 for (m = 0; m < molb->nmol; m++)
683 am = m*cgs->index[cgs->nr];
684 for (cg = 0; cg < cgs->nr; cg++)
687 a1 = cgs->index[cg+1];
688 if (ggrpnr(&mtop->groups, egcENER, a_offset+am+a0) !=
689 ggrpnr(&mtop->groups, egcENER, a_offset +a0))
693 if (mtop->groups.grpnr[egcQMMM] != NULL)
695 for (ai = a0; ai < a1; ai++)
697 if (mtop->groups.grpnr[egcQMMM][a_offset+am+ai] !=
698 mtop->groups.grpnr[egcQMMM][a_offset +ai])
707 cginfo_mb[mb].cg_start = cg_offset;
708 cginfo_mb[mb].cg_end = cg_offset + molb->nmol*cgs->nr;
709 cginfo_mb[mb].cg_mod = (bId ? 1 : molb->nmol)*cgs->nr;
710 snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
711 cginfo = cginfo_mb[mb].cginfo;
713 /* Set constraints flags for constrained atoms */
714 snew(a_con, molt->atoms.nr);
715 for (ftype = 0; ftype < F_NRE; ftype++)
717 if (interaction_function[ftype].flags & IF_CONSTRAINT)
722 for (ia = 0; ia < molt->ilist[ftype].nr; ia += 1+nral)
726 for (a = 0; a < nral; a++)
728 a_con[molt->ilist[ftype].iatoms[ia+1+a]] =
729 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
735 for (m = 0; m < (bId ? 1 : molb->nmol); m++)
738 am = m*cgs->index[cgs->nr];
739 for (cg = 0; cg < cgs->nr; cg++)
742 a1 = cgs->index[cg+1];
744 /* Store the energy group in cginfo */
745 gid = ggrpnr(&mtop->groups, egcENER, a_offset+am+a0);
746 SET_CGINFO_GID(cginfo[cgm+cg], gid);
748 /* Check the intra/inter charge group exclusions */
749 if (a1-a0 > excl_nalloc)
751 excl_nalloc = a1 - a0;
752 srenew(bExcl, excl_nalloc);
754 /* bExclIntraAll: all intra cg interactions excluded
755 * bExclInter: any inter cg interactions excluded
757 bExclIntraAll = TRUE;
761 bHavePerturbedAtoms = FALSE;
762 for (ai = a0; ai < a1; ai++)
764 /* Check VDW and electrostatic interactions */
765 bHaveVDW = bHaveVDW || (type_VDW[molt->atoms.atom[ai].type] ||
766 type_VDW[molt->atoms.atom[ai].typeB]);
767 bHaveQ = bHaveQ || (molt->atoms.atom[ai].q != 0 ||
768 molt->atoms.atom[ai].qB != 0);
770 bHavePerturbedAtoms = bHavePerturbedAtoms || (PERTURBED(molt->atoms.atom[ai]) != 0);
772 /* Clear the exclusion list for atom ai */
773 for (aj = a0; aj < a1; aj++)
775 bExcl[aj-a0] = FALSE;
777 /* Loop over all the exclusions of atom ai */
778 for (j = excl->index[ai]; j < excl->index[ai+1]; j++)
781 if (aj < a0 || aj >= a1)
790 /* Check if ai excludes a0 to a1 */
791 for (aj = a0; aj < a1; aj++)
795 bExclIntraAll = FALSE;
802 SET_CGINFO_CONSTR(cginfo[cgm+cg]);
805 SET_CGINFO_SETTLE(cginfo[cgm+cg]);
813 SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
817 SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
819 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
821 /* The size in cginfo is currently only read with DD */
822 gmx_fatal(FARGS, "A charge group has size %d which is larger than the limit of %d atoms", a1-a0, MAX_CHARGEGROUP_SIZE);
826 SET_CGINFO_HAS_VDW(cginfo[cgm+cg]);
830 SET_CGINFO_HAS_Q(cginfo[cgm+cg]);
832 if (bHavePerturbedAtoms && fr->efep != efepNO)
834 SET_CGINFO_FEP(cginfo[cgm+cg]);
835 *bFEP_NonBonded = TRUE;
837 /* Store the charge group size */
838 SET_CGINFO_NATOMS(cginfo[cgm+cg], a1-a0);
840 if (!bExclIntraAll || bExclInter)
842 *bExcl_IntraCGAll_InterCGNone = FALSE;
849 cg_offset += molb->nmol*cgs->nr;
850 a_offset += molb->nmol*cgs->index[cgs->nr];
854 /* the solvent optimizer is called after the QM is initialized,
855 * because we don't want to have the QM subsystemto become an
859 check_solvent(fplog, mtop, fr, cginfo_mb);
861 if (getenv("GMX_NO_SOLV_OPT"))
865 fprintf(fplog, "Found environment variable GMX_NO_SOLV_OPT.\n"
866 "Disabling all solvent optimization\n");
868 fr->solvent_opt = esolNO;
872 fr->solvent_opt = esolNO;
874 if (!fr->solvent_opt)
876 for (mb = 0; mb < mtop->nmolblock; mb++)
878 for (cg = 0; cg < cginfo_mb[mb].cg_mod; cg++)
880 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg], esolNO);
888 static int *cginfo_expand(int nmb, cginfo_mb_t *cgi_mb)
893 ncg = cgi_mb[nmb-1].cg_end;
896 for (cg = 0; cg < ncg; cg++)
898 while (cg >= cgi_mb[mb].cg_end)
903 cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
909 static void set_chargesum(FILE *log, t_forcerec *fr, const gmx_mtop_t *mtop)
911 /*This now calculates sum for q and c6*/
912 double qsum, q2sum, q, c6sum, c6;
914 const t_atoms *atoms;
919 for (mb = 0; mb < mtop->nmolblock; mb++)
921 nmol = mtop->molblock[mb].nmol;
922 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
923 for (i = 0; i < atoms->nr; i++)
925 q = atoms->atom[i].q;
928 c6 = mtop->ffparams.iparams[atoms->atom[i].type*(mtop->ffparams.atnr+1)].lj.c6;
933 fr->q2sum[0] = q2sum;
934 fr->c6sum[0] = c6sum;
936 if (fr->efep != efepNO)
941 for (mb = 0; mb < mtop->nmolblock; mb++)
943 nmol = mtop->molblock[mb].nmol;
944 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
945 for (i = 0; i < atoms->nr; i++)
947 q = atoms->atom[i].qB;
950 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB*(mtop->ffparams.atnr+1)].lj.c6;
954 fr->q2sum[1] = q2sum;
955 fr->c6sum[1] = c6sum;
960 fr->qsum[1] = fr->qsum[0];
961 fr->q2sum[1] = fr->q2sum[0];
962 fr->c6sum[1] = fr->c6sum[0];
966 if (fr->efep == efepNO)
968 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
972 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n",
973 fr->qsum[0], fr->qsum[1]);
978 void update_forcerec(t_forcerec *fr, matrix box)
980 if (fr->eeltype == eelGRF)
982 calc_rffac(NULL, fr->eeltype, fr->epsilon_r, fr->epsilon_rf,
983 fr->rcoulomb, fr->temp, fr->zsquare, box,
984 &fr->kappa, &fr->k_rf, &fr->c_rf);
988 void set_avcsixtwelve(FILE *fplog, t_forcerec *fr, const gmx_mtop_t *mtop)
990 const t_atoms *atoms, *atoms_tpi;
991 const t_blocka *excl;
992 int mb, nmol, nmolc, i, j, tpi, tpj, j1, j2, k, n, nexcl, q;
993 gmx_int64_t npair, npair_ij, tmpi, tmpj;
994 double csix, ctwelve;
998 real *nbfp_comb = NULL;
1004 /* For LJ-PME, we want to correct for the difference between the
1005 * actual C6 values and the C6 values used by the LJ-PME based on
1006 * combination rules. */
1008 if (EVDW_PME(fr->vdwtype))
1010 nbfp_comb = mk_nbfp_combination_rule(&mtop->ffparams,
1011 (fr->ljpme_combination_rule == eljpmeLB) ? eCOMB_ARITHMETIC : eCOMB_GEOMETRIC);
1012 for (tpi = 0; tpi < ntp; ++tpi)
1014 for (tpj = 0; tpj < ntp; ++tpj)
1016 C6(nbfp_comb, ntp, tpi, tpj) =
1017 C6(nbfp, ntp, tpi, tpj) - C6(nbfp_comb, ntp, tpi, tpj);
1018 C12(nbfp_comb, ntp, tpi, tpj) = C12(nbfp, ntp, tpi, tpj);
1023 for (q = 0; q < (fr->efep == efepNO ? 1 : 2); q++)
1031 /* Count the types so we avoid natoms^2 operations */
1032 snew(typecount, ntp);
1033 gmx_mtop_count_atomtypes(mtop, q, typecount);
1035 for (tpi = 0; tpi < ntp; tpi++)
1037 for (tpj = tpi; tpj < ntp; tpj++)
1039 tmpi = typecount[tpi];
1040 tmpj = typecount[tpj];
1043 npair_ij = tmpi*tmpj;
1047 npair_ij = tmpi*(tmpi - 1)/2;
1051 /* nbfp now includes the 6.0 derivative prefactor */
1052 csix += npair_ij*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1056 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1057 csix += npair_ij* C6(nbfp, ntp, tpi, tpj)/6.0;
1058 ctwelve += npair_ij* C12(nbfp, ntp, tpi, tpj)/12.0;
1064 /* Subtract the excluded pairs.
1065 * The main reason for substracting exclusions is that in some cases
1066 * some combinations might never occur and the parameters could have
1067 * any value. These unused values should not influence the dispersion
1070 for (mb = 0; mb < mtop->nmolblock; mb++)
1072 nmol = mtop->molblock[mb].nmol;
1073 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
1074 excl = &mtop->moltype[mtop->molblock[mb].type].excls;
1075 for (i = 0; (i < atoms->nr); i++)
1079 tpi = atoms->atom[i].type;
1083 tpi = atoms->atom[i].typeB;
1085 j1 = excl->index[i];
1086 j2 = excl->index[i+1];
1087 for (j = j1; j < j2; j++)
1094 tpj = atoms->atom[k].type;
1098 tpj = atoms->atom[k].typeB;
1102 /* nbfp now includes the 6.0 derivative prefactor */
1103 csix -= nmol*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1107 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1108 csix -= nmol*C6 (nbfp, ntp, tpi, tpj)/6.0;
1109 ctwelve -= nmol*C12(nbfp, ntp, tpi, tpj)/12.0;
1119 /* Only correct for the interaction of the test particle
1120 * with the rest of the system.
1123 &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].atoms;
1126 for (mb = 0; mb < mtop->nmolblock; mb++)
1128 nmol = mtop->molblock[mb].nmol;
1129 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
1130 for (j = 0; j < atoms->nr; j++)
1133 /* Remove the interaction of the test charge group
1136 if (mb == mtop->nmolblock-1)
1140 if (mb == 0 && nmol == 1)
1142 gmx_fatal(FARGS, "Old format tpr with TPI, please generate a new tpr file");
1147 tpj = atoms->atom[j].type;
1151 tpj = atoms->atom[j].typeB;
1153 for (i = 0; i < fr->n_tpi; i++)
1157 tpi = atoms_tpi->atom[i].type;
1161 tpi = atoms_tpi->atom[i].typeB;
1165 /* nbfp now includes the 6.0 derivative prefactor */
1166 csix += nmolc*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1170 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1171 csix += nmolc*C6 (nbfp, ntp, tpi, tpj)/6.0;
1172 ctwelve += nmolc*C12(nbfp, ntp, tpi, tpj)/12.0;
1179 if (npair - nexcl <= 0 && fplog)
1181 fprintf(fplog, "\nWARNING: There are no atom pairs for dispersion correction\n\n");
1187 csix /= npair - nexcl;
1188 ctwelve /= npair - nexcl;
1192 fprintf(debug, "Counted %d exclusions\n", nexcl);
1193 fprintf(debug, "Average C6 parameter is: %10g\n", (double)csix);
1194 fprintf(debug, "Average C12 parameter is: %10g\n", (double)ctwelve);
1196 fr->avcsix[q] = csix;
1197 fr->avctwelve[q] = ctwelve;
1200 if (EVDW_PME(fr->vdwtype))
1207 if (fr->eDispCorr == edispcAllEner ||
1208 fr->eDispCorr == edispcAllEnerPres)
1210 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1211 fr->avcsix[0], fr->avctwelve[0]);
1215 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e\n", fr->avcsix[0]);
1221 static void set_bham_b_max(FILE *fplog, t_forcerec *fr,
1222 const gmx_mtop_t *mtop)
1224 const t_atoms *at1, *at2;
1225 int mt1, mt2, i, j, tpi, tpj, ntypes;
1231 fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
1238 for (mt1 = 0; mt1 < mtop->nmoltype; mt1++)
1240 at1 = &mtop->moltype[mt1].atoms;
1241 for (i = 0; (i < at1->nr); i++)
1243 tpi = at1->atom[i].type;
1246 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
1249 for (mt2 = mt1; mt2 < mtop->nmoltype; mt2++)
1251 at2 = &mtop->moltype[mt2].atoms;
1252 for (j = 0; (j < at2->nr); j++)
1254 tpj = at2->atom[j].type;
1257 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
1259 b = BHAMB(nbfp, ntypes, tpi, tpj);
1260 if (b > fr->bham_b_max)
1264 if ((b < bmin) || (bmin == -1))
1274 fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n",
1275 bmin, fr->bham_b_max);
1279 static void make_nbf_tables(FILE *fp, const output_env_t oenv,
1280 t_forcerec *fr, real rtab,
1281 const t_commrec *cr,
1282 const char *tabfn, char *eg1, char *eg2,
1292 fprintf(debug, "No table file name passed, can not read table, can not do non-bonded interactions\n");
1297 sprintf(buf, "%s", tabfn);
1300 /* Append the two energy group names */
1301 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
1302 eg1, eg2, ftp2ext(efXVG));
1304 nbl->table_elec_vdw = make_tables(fp, oenv, fr, MASTER(cr), buf, rtab, 0);
1305 /* Copy the contents of the table to separate coulomb and LJ tables too,
1306 * to improve cache performance.
1308 /* For performance reasons we want
1309 * the table data to be aligned to 16-byte. The pointers could be freed
1310 * but currently aren't.
1312 nbl->table_elec.interaction = GMX_TABLE_INTERACTION_ELEC;
1313 nbl->table_elec.format = nbl->table_elec_vdw.format;
1314 nbl->table_elec.r = nbl->table_elec_vdw.r;
1315 nbl->table_elec.n = nbl->table_elec_vdw.n;
1316 nbl->table_elec.scale = nbl->table_elec_vdw.scale;
1317 nbl->table_elec.scale_exp = nbl->table_elec_vdw.scale_exp;
1318 nbl->table_elec.formatsize = nbl->table_elec_vdw.formatsize;
1319 nbl->table_elec.ninteractions = 1;
1320 nbl->table_elec.stride = nbl->table_elec.formatsize * nbl->table_elec.ninteractions;
1321 snew_aligned(nbl->table_elec.data, nbl->table_elec.stride*(nbl->table_elec.n+1), 32);
1323 nbl->table_vdw.interaction = GMX_TABLE_INTERACTION_VDWREP_VDWDISP;
1324 nbl->table_vdw.format = nbl->table_elec_vdw.format;
1325 nbl->table_vdw.r = nbl->table_elec_vdw.r;
1326 nbl->table_vdw.n = nbl->table_elec_vdw.n;
1327 nbl->table_vdw.scale = nbl->table_elec_vdw.scale;
1328 nbl->table_vdw.scale_exp = nbl->table_elec_vdw.scale_exp;
1329 nbl->table_vdw.formatsize = nbl->table_elec_vdw.formatsize;
1330 nbl->table_vdw.ninteractions = 2;
1331 nbl->table_vdw.stride = nbl->table_vdw.formatsize * nbl->table_vdw.ninteractions;
1332 snew_aligned(nbl->table_vdw.data, nbl->table_vdw.stride*(nbl->table_vdw.n+1), 32);
1334 for (i = 0; i <= nbl->table_elec_vdw.n; i++)
1336 for (j = 0; j < 4; j++)
1338 nbl->table_elec.data[4*i+j] = nbl->table_elec_vdw.data[12*i+j];
1340 for (j = 0; j < 8; j++)
1342 nbl->table_vdw.data[8*i+j] = nbl->table_elec_vdw.data[12*i+4+j];
1347 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t *mtop,
1348 int *ncount, int **count)
1350 const gmx_moltype_t *molt;
1352 int mt, ftype, stride, i, j, tabnr;
1354 for (mt = 0; mt < mtop->nmoltype; mt++)
1356 molt = &mtop->moltype[mt];
1357 for (ftype = 0; ftype < F_NRE; ftype++)
1359 if (ftype == ftype1 || ftype == ftype2)
1361 il = &molt->ilist[ftype];
1362 stride = 1 + NRAL(ftype);
1363 for (i = 0; i < il->nr; i += stride)
1365 tabnr = mtop->ffparams.iparams[il->iatoms[i]].tab.table;
1368 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
1370 if (tabnr >= *ncount)
1372 srenew(*count, tabnr+1);
1373 for (j = *ncount; j < tabnr+1; j++)
1386 static bondedtable_t *make_bonded_tables(FILE *fplog,
1387 int ftype1, int ftype2,
1388 const gmx_mtop_t *mtop,
1389 const char *basefn, const char *tabext)
1391 int i, ncount, *count;
1399 count_tables(ftype1, ftype2, mtop, &ncount, &count);
1404 for (i = 0; i < ncount; i++)
1408 sprintf(tabfn, "%s", basefn);
1409 sprintf(tabfn + strlen(basefn) - strlen(ftp2ext(efXVG)) - 1, "_%s%d.%s",
1410 tabext, i, ftp2ext(efXVG));
1411 tab[i] = make_bonded_table(fplog, tabfn, NRAL(ftype1)-2);
1420 void forcerec_set_ranges(t_forcerec *fr,
1421 int ncg_home, int ncg_force,
1423 int natoms_force_constr, int natoms_f_novirsum)
1428 /* fr->ncg_force is unused in the standard code,
1429 * but it can be useful for modified code dealing with charge groups.
1431 fr->ncg_force = ncg_force;
1432 fr->natoms_force = natoms_force;
1433 fr->natoms_force_constr = natoms_force_constr;
1435 if (fr->natoms_force_constr > fr->nalloc_force)
1437 fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1441 srenew(fr->f_twin, fr->nalloc_force);
1445 if (fr->bF_NoVirSum)
1447 fr->f_novirsum_n = natoms_f_novirsum;
1448 if (fr->f_novirsum_n > fr->f_novirsum_nalloc)
1450 fr->f_novirsum_nalloc = over_alloc_dd(fr->f_novirsum_n);
1451 srenew(fr->f_novirsum_alloc, fr->f_novirsum_nalloc);
1456 fr->f_novirsum_n = 0;
1460 static real cutoff_inf(real cutoff)
1464 cutoff = GMX_CUTOFF_INF;
1470 static void make_adress_tf_tables(FILE *fp, const output_env_t oenv,
1471 t_forcerec *fr, const t_inputrec *ir,
1472 const char *tabfn, const gmx_mtop_t *mtop,
1480 gmx_fatal(FARGS, "No thermoforce table file given. Use -tabletf to specify a file\n");
1484 snew(fr->atf_tabs, ir->adress->n_tf_grps);
1486 sprintf(buf, "%s", tabfn);
1487 for (i = 0; i < ir->adress->n_tf_grps; i++)
1489 j = ir->adress->tf_table_index[i]; /* get energy group index */
1490 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "tf_%s.%s",
1491 *(mtop->groups.grpname[mtop->groups.grps[egcENER].nm_ind[j]]), ftp2ext(efXVG));
1494 fprintf(fp, "loading tf table for energygrp index %d from %s\n", ir->adress->tf_table_index[i], buf);
1496 fr->atf_tabs[i] = make_atf_table(fp, oenv, fr, buf, box);
1501 gmx_bool can_use_allvsall(const t_inputrec *ir, gmx_bool bPrintNote, t_commrec *cr, FILE *fp)
1508 ir->rcoulomb == 0 &&
1510 ir->ePBC == epbcNONE &&
1511 ir->vdwtype == evdwCUT &&
1512 ir->coulombtype == eelCUT &&
1513 ir->efep == efepNO &&
1514 (ir->implicit_solvent == eisNO ||
1515 (ir->implicit_solvent == eisGBSA && (ir->gb_algorithm == egbSTILL ||
1516 ir->gb_algorithm == egbHCT ||
1517 ir->gb_algorithm == egbOBC))) &&
1518 getenv("GMX_NO_ALLVSALL") == NULL
1521 if (bAllvsAll && ir->opts.ngener > 1)
1523 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";
1529 fprintf(stderr, "\n%s\n", note);
1533 fprintf(fp, "\n%s\n", note);
1539 if (bAllvsAll && fp && MASTER(cr))
1541 fprintf(fp, "\nUsing SIMD all-vs-all kernels.\n\n");
1548 static void init_forcerec_f_threads(t_forcerec *fr, int nenergrp)
1552 /* These thread local data structures are used for bondeds only */
1553 fr->nthreads = gmx_omp_nthreads_get(emntBonded);
1555 if (fr->nthreads > 1)
1557 snew(fr->f_t, fr->nthreads);
1558 /* Thread 0 uses the global force and energy arrays */
1559 for (t = 1; t < fr->nthreads; t++)
1561 fr->f_t[t].f = NULL;
1562 fr->f_t[t].f_nalloc = 0;
1563 snew(fr->f_t[t].fshift, SHIFTS);
1564 fr->f_t[t].grpp.nener = nenergrp*nenergrp;
1565 for (i = 0; i < egNR; i++)
1567 snew(fr->f_t[t].grpp.ener[i], fr->f_t[t].grpp.nener);
1574 gmx_bool nbnxn_acceleration_supported(FILE *fplog,
1575 const t_commrec *cr,
1576 const t_inputrec *ir,
1579 if (!bGPU && (ir->vdwtype == evdwPME && ir->ljpme_combination_rule == eljpmeLB))
1581 md_print_warn(cr, fplog, "LJ-PME with Lorentz-Berthelot is not supported with %s, falling back to %s\n",
1582 bGPU ? "GPUs" : "SIMD kernels",
1583 bGPU ? "CPU only" : "plain-C kernels");
1591 static void pick_nbnxn_kernel_cpu(const t_inputrec gmx_unused *ir,
1595 *kernel_type = nbnxnk4x4_PlainC;
1596 *ewald_excl = ewaldexclTable;
1598 #ifdef GMX_NBNXN_SIMD
1600 #ifdef GMX_NBNXN_SIMD_4XN
1601 *kernel_type = nbnxnk4xN_SIMD_4xN;
1603 #ifdef GMX_NBNXN_SIMD_2XNN
1604 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1607 #if defined GMX_NBNXN_SIMD_2XNN && defined GMX_NBNXN_SIMD_4XN
1608 /* We need to choose if we want 2x(N+N) or 4xN kernels.
1609 * Currently this is based on the SIMD acceleration choice,
1610 * but it might be better to decide this at runtime based on CPU.
1612 * 4xN calculates more (zero) interactions, but has less pair-search
1613 * work and much better kernel instruction scheduling.
1615 * Up till now we have only seen that on Intel Sandy/Ivy Bridge,
1616 * which doesn't have FMA, both the analytical and tabulated Ewald
1617 * kernels have similar pair rates for 4x8 and 2x(4+4), so we choose
1618 * 2x(4+4) because it results in significantly fewer pairs.
1619 * For RF, the raw pair rate of the 4x8 kernel is higher than 2x(4+4),
1620 * 10% with HT, 50% without HT. As we currently don't detect the actual
1621 * use of HT, use 4x8 to avoid a potential performance hit.
1622 * On Intel Haswell 4x8 is always faster.
1624 *kernel_type = nbnxnk4xN_SIMD_4xN;
1626 #ifndef GMX_SIMD_HAVE_FMA
1627 if (EEL_PME_EWALD(ir->coulombtype) ||
1628 EVDW_PME(ir->vdwtype))
1630 /* We have Ewald kernels without FMA (Intel Sandy/Ivy Bridge).
1631 * There are enough instructions to make 2x(4+4) efficient.
1633 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1636 #endif /* GMX_NBNXN_SIMD_2XNN && GMX_NBNXN_SIMD_4XN */
1639 if (getenv("GMX_NBNXN_SIMD_4XN") != NULL)
1641 #ifdef GMX_NBNXN_SIMD_4XN
1642 *kernel_type = nbnxnk4xN_SIMD_4xN;
1644 gmx_fatal(FARGS, "SIMD 4xN kernels requested, but Gromacs has been compiled without support for these kernels");
1647 if (getenv("GMX_NBNXN_SIMD_2XNN") != NULL)
1649 #ifdef GMX_NBNXN_SIMD_2XNN
1650 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1652 gmx_fatal(FARGS, "SIMD 2x(N+N) kernels requested, but Gromacs has been compiled without support for these kernels");
1656 /* Analytical Ewald exclusion correction is only an option in
1658 * Since table lookup's don't parallelize with SIMD, analytical
1659 * will probably always be faster for a SIMD width of 8 or more.
1660 * With FMA analytical is sometimes faster for a width if 4 as well.
1661 * On BlueGene/Q, this is faster regardless of precision.
1662 * In single precision, this is faster on Bulldozer.
1664 #if GMX_SIMD_REAL_WIDTH >= 8 || \
1665 (GMX_SIMD_REAL_WIDTH >= 4 && defined GMX_SIMD_HAVE_FMA && !defined GMX_DOUBLE) || \
1666 defined GMX_SIMD_IBM_QPX
1667 *ewald_excl = ewaldexclAnalytical;
1669 if (getenv("GMX_NBNXN_EWALD_TABLE") != NULL)
1671 *ewald_excl = ewaldexclTable;
1673 if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != NULL)
1675 *ewald_excl = ewaldexclAnalytical;
1679 #endif /* GMX_NBNXN_SIMD */
1683 const char *lookup_nbnxn_kernel_name(int kernel_type)
1685 const char *returnvalue = NULL;
1686 switch (kernel_type)
1689 returnvalue = "not set";
1691 case nbnxnk4x4_PlainC:
1692 returnvalue = "plain C";
1694 case nbnxnk4xN_SIMD_4xN:
1695 case nbnxnk4xN_SIMD_2xNN:
1696 #ifdef GMX_NBNXN_SIMD
1697 #if defined GMX_SIMD_X86_SSE2
1698 returnvalue = "SSE2";
1699 #elif defined GMX_SIMD_X86_SSE4_1
1700 returnvalue = "SSE4.1";
1701 #elif defined GMX_SIMD_X86_AVX_128_FMA
1702 returnvalue = "AVX_128_FMA";
1703 #elif defined GMX_SIMD_X86_AVX_256
1704 returnvalue = "AVX_256";
1705 #elif defined GMX_SIMD_X86_AVX2_256
1706 returnvalue = "AVX2_256";
1708 returnvalue = "SIMD";
1710 #else /* GMX_NBNXN_SIMD */
1711 returnvalue = "not available";
1712 #endif /* GMX_NBNXN_SIMD */
1714 case nbnxnk8x8x8_CUDA: returnvalue = "CUDA"; break;
1715 case nbnxnk8x8x8_PlainC: returnvalue = "plain C"; break;
1719 gmx_fatal(FARGS, "Illegal kernel type selected");
1726 static void pick_nbnxn_kernel(FILE *fp,
1727 const t_commrec *cr,
1728 gmx_bool use_simd_kernels,
1730 gmx_bool bEmulateGPU,
1731 const t_inputrec *ir,
1734 gmx_bool bDoNonbonded)
1736 assert(kernel_type);
1738 *kernel_type = nbnxnkNotSet;
1739 *ewald_excl = ewaldexclTable;
1743 *kernel_type = nbnxnk8x8x8_PlainC;
1747 md_print_warn(cr, fp, "Emulating a GPU run on the CPU (slow)");
1752 *kernel_type = nbnxnk8x8x8_CUDA;
1755 if (*kernel_type == nbnxnkNotSet)
1757 /* LJ PME with LB combination rule does 7 mesh operations.
1758 * This so slow that we don't compile SIMD non-bonded kernels for that.
1760 if (use_simd_kernels &&
1761 nbnxn_acceleration_supported(fp, cr, ir, FALSE))
1763 pick_nbnxn_kernel_cpu(ir, kernel_type, ewald_excl);
1767 *kernel_type = nbnxnk4x4_PlainC;
1771 if (bDoNonbonded && fp != NULL)
1773 fprintf(fp, "\nUsing %s %dx%d non-bonded kernels\n\n",
1774 lookup_nbnxn_kernel_name(*kernel_type),
1775 nbnxn_kernel_pairlist_simple(*kernel_type) ? NBNXN_CPU_CLUSTER_I_SIZE : NBNXN_GPU_CLUSTER_SIZE,
1776 nbnxn_kernel_to_cj_size(*kernel_type));
1780 static void pick_nbnxn_resources(const t_commrec *cr,
1781 const gmx_hw_info_t *hwinfo,
1782 gmx_bool bDoNonbonded,
1784 gmx_bool *bEmulateGPU,
1785 const gmx_gpu_opt_t *gpu_opt)
1787 gmx_bool bEmulateGPUEnvVarSet;
1788 char gpu_err_str[STRLEN];
1792 bEmulateGPUEnvVarSet = (getenv("GMX_EMULATE_GPU") != NULL);
1794 /* Run GPU emulation mode if GMX_EMULATE_GPU is defined. Because
1795 * GPUs (currently) only handle non-bonded calculations, we will
1796 * automatically switch to emulation if non-bonded calculations are
1797 * turned off via GMX_NO_NONBONDED - this is the simple and elegant
1798 * way to turn off GPU initialization, data movement, and cleanup.
1800 * GPU emulation can be useful to assess the performance one can expect by
1801 * adding GPU(s) to the machine. The conditional below allows this even
1802 * if mdrun is compiled without GPU acceleration support.
1803 * Note that you should freezing the system as otherwise it will explode.
1805 *bEmulateGPU = (bEmulateGPUEnvVarSet ||
1807 gpu_opt->ncuda_dev_use > 0));
1809 /* Enable GPU mode when GPUs are available or no GPU emulation is requested.
1811 if (gpu_opt->ncuda_dev_use > 0 && !(*bEmulateGPU))
1813 /* Each PP node will use the intra-node id-th device from the
1814 * list of detected/selected GPUs. */
1815 if (!init_gpu(cr->rank_pp_intranode, gpu_err_str,
1816 &hwinfo->gpu_info, gpu_opt))
1818 /* At this point the init should never fail as we made sure that
1819 * we have all the GPUs we need. If it still does, we'll bail. */
1820 gmx_fatal(FARGS, "On node %d failed to initialize GPU #%d: %s",
1822 get_gpu_device_id(&hwinfo->gpu_info, gpu_opt,
1823 cr->rank_pp_intranode),
1827 /* Here we actually turn on hardware GPU acceleration */
1832 gmx_bool uses_simple_tables(int cutoff_scheme,
1833 nonbonded_verlet_t *nbv,
1836 gmx_bool bUsesSimpleTables = TRUE;
1839 switch (cutoff_scheme)
1842 bUsesSimpleTables = TRUE;
1845 assert(NULL != nbv && NULL != nbv->grp);
1846 grp_index = (group < 0) ? 0 : (nbv->ngrp - 1);
1847 bUsesSimpleTables = nbnxn_kernel_pairlist_simple(nbv->grp[grp_index].kernel_type);
1850 gmx_incons("unimplemented");
1852 return bUsesSimpleTables;
1855 static void init_ewald_f_table(interaction_const_t *ic,
1856 gmx_bool bUsesSimpleTables,
1861 if (bUsesSimpleTables)
1863 /* With a spacing of 0.0005 we are at the force summation accuracy
1864 * for the SSE kernels for "normal" atomistic simulations.
1866 ic->tabq_scale = ewald_spline3_table_scale(ic->ewaldcoeff_q,
1869 maxr = (rtab > ic->rcoulomb) ? rtab : ic->rcoulomb;
1870 ic->tabq_size = (int)(maxr*ic->tabq_scale) + 2;
1874 ic->tabq_size = GPU_EWALD_COULOMB_FORCE_TABLE_SIZE;
1875 /* Subtract 2 iso 1 to avoid access out of range due to rounding */
1876 ic->tabq_scale = (ic->tabq_size - 2)/ic->rcoulomb;
1879 sfree_aligned(ic->tabq_coul_FDV0);
1880 sfree_aligned(ic->tabq_coul_F);
1881 sfree_aligned(ic->tabq_coul_V);
1883 sfree_aligned(ic->tabq_vdw_FDV0);
1884 sfree_aligned(ic->tabq_vdw_F);
1885 sfree_aligned(ic->tabq_vdw_V);
1887 if (ic->eeltype == eelEWALD || EEL_PME(ic->eeltype))
1889 /* Create the original table data in FDV0 */
1890 snew_aligned(ic->tabq_coul_FDV0, ic->tabq_size*4, 32);
1891 snew_aligned(ic->tabq_coul_F, ic->tabq_size, 32);
1892 snew_aligned(ic->tabq_coul_V, ic->tabq_size, 32);
1893 table_spline3_fill_ewald_lr(ic->tabq_coul_F, ic->tabq_coul_V, ic->tabq_coul_FDV0,
1894 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_q, v_q_ewald_lr);
1897 if (EVDW_PME(ic->vdwtype))
1899 snew_aligned(ic->tabq_vdw_FDV0, ic->tabq_size*4, 32);
1900 snew_aligned(ic->tabq_vdw_F, ic->tabq_size, 32);
1901 snew_aligned(ic->tabq_vdw_V, ic->tabq_size, 32);
1902 table_spline3_fill_ewald_lr(ic->tabq_vdw_F, ic->tabq_vdw_V, ic->tabq_vdw_FDV0,
1903 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_lj, v_lj_ewald_lr);
1907 void init_interaction_const_tables(FILE *fp,
1908 interaction_const_t *ic,
1909 gmx_bool bUsesSimpleTables,
1914 if (ic->eeltype == eelEWALD || EEL_PME(ic->eeltype) || EVDW_PME(ic->vdwtype))
1916 init_ewald_f_table(ic, bUsesSimpleTables, rtab);
1920 fprintf(fp, "Initialized non-bonded Ewald correction tables, spacing: %.2e size: %d\n\n",
1921 1/ic->tabq_scale, ic->tabq_size);
1926 static void clear_force_switch_constants(shift_consts_t *sc)
1933 static void force_switch_constants(real p,
1937 /* Here we determine the coefficient for shifting the force to zero
1938 * between distance rsw and the cut-off rc.
1939 * For a potential of r^-p, we have force p*r^-(p+1).
1940 * But to save flops we absorb p in the coefficient.
1942 * force/p = r^-(p+1) + c2*r^2 + c3*r^3
1943 * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
1945 sc->c2 = ((p + 1)*rsw - (p + 4)*rc)/(pow(rc, p + 2)*pow(rc - rsw, 2));
1946 sc->c3 = -((p + 1)*rsw - (p + 3)*rc)/(pow(rc, p + 2)*pow(rc - rsw, 3));
1947 sc->cpot = -pow(rc, -p) + p*sc->c2/3*pow(rc - rsw, 3) + p*sc->c3/4*pow(rc - rsw, 4);
1950 static void potential_switch_constants(real rsw, real rc,
1951 switch_consts_t *sc)
1953 /* The switch function is 1 at rsw and 0 at rc.
1954 * The derivative and second derivate are zero at both ends.
1955 * rsw = max(r - r_switch, 0)
1956 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
1957 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
1958 * force = force*dsw - potential*sw
1961 sc->c3 = -10*pow(rc - rsw, -3);
1962 sc->c4 = 15*pow(rc - rsw, -4);
1963 sc->c5 = -6*pow(rc - rsw, -5);
1967 init_interaction_const(FILE *fp,
1968 const t_commrec gmx_unused *cr,
1969 interaction_const_t **interaction_const,
1970 const t_forcerec *fr,
1973 interaction_const_t *ic;
1974 gmx_bool bUsesSimpleTables = TRUE;
1978 /* Just allocate something so we can free it */
1979 snew_aligned(ic->tabq_coul_FDV0, 16, 32);
1980 snew_aligned(ic->tabq_coul_F, 16, 32);
1981 snew_aligned(ic->tabq_coul_V, 16, 32);
1983 ic->rlist = fr->rlist;
1984 ic->rlistlong = fr->rlistlong;
1987 ic->vdwtype = fr->vdwtype;
1988 ic->vdw_modifier = fr->vdw_modifier;
1989 ic->rvdw = fr->rvdw;
1990 ic->rvdw_switch = fr->rvdw_switch;
1991 ic->ewaldcoeff_lj = fr->ewaldcoeff_lj;
1992 ic->ljpme_comb_rule = fr->ljpme_combination_rule;
1993 ic->sh_lj_ewald = 0;
1994 clear_force_switch_constants(&ic->dispersion_shift);
1995 clear_force_switch_constants(&ic->repulsion_shift);
1997 switch (ic->vdw_modifier)
1999 case eintmodPOTSHIFT:
2000 /* Only shift the potential, don't touch the force */
2001 ic->dispersion_shift.cpot = -pow(ic->rvdw, -6.0);
2002 ic->repulsion_shift.cpot = -pow(ic->rvdw, -12.0);
2003 if (EVDW_PME(ic->vdwtype))
2007 crc2 = sqr(ic->ewaldcoeff_lj*ic->rvdw);
2008 ic->sh_lj_ewald = (exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)*pow(ic->rvdw, -6.0);
2011 case eintmodFORCESWITCH:
2012 /* Switch the force, switch and shift the potential */
2013 force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw,
2014 &ic->dispersion_shift);
2015 force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw,
2016 &ic->repulsion_shift);
2018 case eintmodPOTSWITCH:
2019 /* Switch the potential and force */
2020 potential_switch_constants(ic->rvdw_switch, ic->rvdw,
2024 case eintmodEXACTCUTOFF:
2025 /* Nothing to do here */
2028 gmx_incons("unimplemented potential modifier");
2031 ic->sh_invrc6 = -ic->dispersion_shift.cpot;
2033 /* Electrostatics */
2034 ic->eeltype = fr->eeltype;
2035 ic->coulomb_modifier = fr->coulomb_modifier;
2036 ic->rcoulomb = fr->rcoulomb;
2037 ic->epsilon_r = fr->epsilon_r;
2038 ic->epsfac = fr->epsfac;
2039 ic->ewaldcoeff_q = fr->ewaldcoeff_q;
2041 if (fr->coulomb_modifier == eintmodPOTSHIFT)
2043 ic->sh_ewald = gmx_erfc(ic->ewaldcoeff_q*ic->rcoulomb);
2050 /* Reaction-field */
2051 if (EEL_RF(ic->eeltype))
2053 ic->epsilon_rf = fr->epsilon_rf;
2054 ic->k_rf = fr->k_rf;
2055 ic->c_rf = fr->c_rf;
2059 /* For plain cut-off we might use the reaction-field kernels */
2060 ic->epsilon_rf = ic->epsilon_r;
2062 if (fr->coulomb_modifier == eintmodPOTSHIFT)
2064 ic->c_rf = 1/ic->rcoulomb;
2074 real dispersion_shift;
2076 dispersion_shift = ic->dispersion_shift.cpot;
2077 if (EVDW_PME(ic->vdwtype))
2079 dispersion_shift -= ic->sh_lj_ewald;
2081 fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
2082 ic->repulsion_shift.cpot, dispersion_shift);
2084 if (ic->eeltype == eelCUT)
2086 fprintf(fp, ", Coulomb %.e", -ic->c_rf);
2088 else if (EEL_PME(ic->eeltype))
2090 fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
2095 *interaction_const = ic;
2097 if (fr->nbv != NULL && fr->nbv->bUseGPU)
2099 nbnxn_cuda_init_const(fr->nbv->cu_nbv, ic, fr->nbv->grp);
2101 /* With tMPI + GPUs some ranks may be sharing GPU(s) and therefore
2102 * also sharing texture references. To keep the code simple, we don't
2103 * treat texture references as shared resources, but this means that
2104 * the coulomb_tab and nbfp texture refs will get updated by multiple threads.
2105 * Hence, to ensure that the non-bonded kernels don't start before all
2106 * texture binding operations are finished, we need to wait for all ranks
2107 * to arrive here before continuing.
2109 * Note that we could omit this barrier if GPUs are not shared (or
2110 * texture objects are used), but as this is initialization code, there
2111 * is not point in complicating things.
2113 #ifdef GMX_THREAD_MPI
2118 #endif /* GMX_THREAD_MPI */
2121 bUsesSimpleTables = uses_simple_tables(fr->cutoff_scheme, fr->nbv, -1);
2122 init_interaction_const_tables(fp, ic, bUsesSimpleTables, rtab);
2125 static void init_nb_verlet(FILE *fp,
2126 nonbonded_verlet_t **nb_verlet,
2127 gmx_bool bFEP_NonBonded,
2128 const t_inputrec *ir,
2129 const t_forcerec *fr,
2130 const t_commrec *cr,
2131 const char *nbpu_opt)
2133 nonbonded_verlet_t *nbv;
2136 gmx_bool bEmulateGPU, bHybridGPURun = FALSE;
2138 nbnxn_alloc_t *nb_alloc;
2139 nbnxn_free_t *nb_free;
2143 pick_nbnxn_resources(cr, fr->hwinfo,
2151 nbv->ngrp = (DOMAINDECOMP(cr) ? 2 : 1);
2152 for (i = 0; i < nbv->ngrp; i++)
2154 nbv->grp[i].nbl_lists.nnbl = 0;
2155 nbv->grp[i].nbat = NULL;
2156 nbv->grp[i].kernel_type = nbnxnkNotSet;
2158 if (i == 0) /* local */
2160 pick_nbnxn_kernel(fp, cr, fr->use_simd_kernels,
2161 nbv->bUseGPU, bEmulateGPU, ir,
2162 &nbv->grp[i].kernel_type,
2163 &nbv->grp[i].ewald_excl,
2166 else /* non-local */
2168 if (nbpu_opt != NULL && strcmp(nbpu_opt, "gpu_cpu") == 0)
2170 /* Use GPU for local, select a CPU kernel for non-local */
2171 pick_nbnxn_kernel(fp, cr, fr->use_simd_kernels,
2173 &nbv->grp[i].kernel_type,
2174 &nbv->grp[i].ewald_excl,
2177 bHybridGPURun = TRUE;
2181 /* Use the same kernel for local and non-local interactions */
2182 nbv->grp[i].kernel_type = nbv->grp[0].kernel_type;
2183 nbv->grp[i].ewald_excl = nbv->grp[0].ewald_excl;
2190 /* init the NxN GPU data; the last argument tells whether we'll have
2191 * both local and non-local NB calculation on GPU */
2192 nbnxn_cuda_init(fp, &nbv->cu_nbv,
2193 &fr->hwinfo->gpu_info, fr->gpu_opt,
2194 cr->rank_pp_intranode,
2195 (nbv->ngrp > 1) && !bHybridGPURun);
2197 if ((env = getenv("GMX_NB_MIN_CI")) != NULL)
2201 nbv->min_ci_balanced = strtol(env, &end, 10);
2202 if (!end || (*end != 0) || nbv->min_ci_balanced <= 0)
2204 gmx_fatal(FARGS, "Invalid value passed in GMX_NB_MIN_CI=%s, positive integer required", env);
2209 fprintf(debug, "Neighbor-list balancing parameter: %d (passed as env. var.)\n",
2210 nbv->min_ci_balanced);
2215 nbv->min_ci_balanced = nbnxn_cuda_min_ci_balanced(nbv->cu_nbv);
2218 fprintf(debug, "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU multi-processors)\n",
2219 nbv->min_ci_balanced);
2225 nbv->min_ci_balanced = 0;
2230 nbnxn_init_search(&nbv->nbs,
2231 DOMAINDECOMP(cr) ? &cr->dd->nc : NULL,
2232 DOMAINDECOMP(cr) ? domdec_zones(cr->dd) : NULL,
2234 gmx_omp_nthreads_get(emntNonbonded));
2236 for (i = 0; i < nbv->ngrp; i++)
2238 if (nbv->grp[0].kernel_type == nbnxnk8x8x8_CUDA)
2240 nb_alloc = &pmalloc;
2249 nbnxn_init_pairlist_set(&nbv->grp[i].nbl_lists,
2250 nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2251 /* 8x8x8 "non-simple" lists are ATM always combined */
2252 !nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2256 nbv->grp[0].kernel_type != nbv->grp[i].kernel_type)
2258 gmx_bool bSimpleList;
2259 int enbnxninitcombrule;
2261 bSimpleList = nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type);
2263 if (bSimpleList && (fr->vdwtype == evdwCUT && (fr->vdw_modifier == eintmodNONE || fr->vdw_modifier == eintmodPOTSHIFT)))
2265 /* Plain LJ cut-off: we can optimize with combination rules */
2266 enbnxninitcombrule = enbnxninitcombruleDETECT;
2268 else if (fr->vdwtype == evdwPME)
2270 /* LJ-PME: we need to use a combination rule for the grid */
2271 if (fr->ljpme_combination_rule == eljpmeGEOM)
2273 enbnxninitcombrule = enbnxninitcombruleGEOM;
2277 enbnxninitcombrule = enbnxninitcombruleLB;
2282 /* We use a full combination matrix: no rule required */
2283 enbnxninitcombrule = enbnxninitcombruleNONE;
2287 snew(nbv->grp[i].nbat, 1);
2288 nbnxn_atomdata_init(fp,
2290 nbv->grp[i].kernel_type,
2292 fr->ntype, fr->nbfp,
2294 bSimpleList ? gmx_omp_nthreads_get(emntNonbonded) : 1,
2299 nbv->grp[i].nbat = nbv->grp[0].nbat;
2304 void init_forcerec(FILE *fp,
2305 const output_env_t oenv,
2308 const t_inputrec *ir,
2309 const gmx_mtop_t *mtop,
2310 const t_commrec *cr,
2316 const char *nbpu_opt,
2317 gmx_bool bNoSolvOpt,
2320 int i, j, m, natoms, ngrp, negp_pp, negptable, egi, egj;
2325 gmx_bool bGenericKernelOnly;
2326 gmx_bool bMakeTables, bMakeSeparate14Table, bSomeNormalNbListsAreInUse;
2327 gmx_bool bFEP_NonBonded;
2329 int *nm_ind, egp_flags;
2331 if (fr->hwinfo == NULL)
2333 /* Detect hardware, gather information.
2334 * In mdrun, hwinfo has already been set before calling init_forcerec.
2335 * Here we ignore GPUs, as tools will not use them anyhow.
2337 fr->hwinfo = gmx_detect_hardware(fp, cr, FALSE);
2340 /* By default we turn SIMD kernels on, but it might be turned off further down... */
2341 fr->use_simd_kernels = TRUE;
2343 fr->bDomDec = DOMAINDECOMP(cr);
2345 natoms = mtop->natoms;
2347 if (check_box(ir->ePBC, box))
2349 gmx_fatal(FARGS, check_box(ir->ePBC, box));
2352 /* Test particle insertion ? */
2355 /* Set to the size of the molecule to be inserted (the last one) */
2356 /* Because of old style topologies, we have to use the last cg
2357 * instead of the last molecule type.
2359 cgs = &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].cgs;
2360 fr->n_tpi = cgs->index[cgs->nr] - cgs->index[cgs->nr-1];
2361 if (fr->n_tpi != mtop->mols.index[mtop->mols.nr] - mtop->mols.index[mtop->mols.nr-1])
2363 gmx_fatal(FARGS, "The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
2371 /* Copy AdResS parameters */
2374 fr->adress_type = ir->adress->type;
2375 fr->adress_const_wf = ir->adress->const_wf;
2376 fr->adress_ex_width = ir->adress->ex_width;
2377 fr->adress_hy_width = ir->adress->hy_width;
2378 fr->adress_icor = ir->adress->icor;
2379 fr->adress_site = ir->adress->site;
2380 fr->adress_ex_forcecap = ir->adress->ex_forcecap;
2381 fr->adress_do_hybridpairs = ir->adress->do_hybridpairs;
2384 snew(fr->adress_group_explicit, ir->adress->n_energy_grps);
2385 for (i = 0; i < ir->adress->n_energy_grps; i++)
2387 fr->adress_group_explicit[i] = ir->adress->group_explicit[i];
2390 fr->n_adress_tf_grps = ir->adress->n_tf_grps;
2391 snew(fr->adress_tf_table_index, fr->n_adress_tf_grps);
2392 for (i = 0; i < fr->n_adress_tf_grps; i++)
2394 fr->adress_tf_table_index[i] = ir->adress->tf_table_index[i];
2396 copy_rvec(ir->adress->refs, fr->adress_refs);
2400 fr->adress_type = eAdressOff;
2401 fr->adress_do_hybridpairs = FALSE;
2404 /* Copy the user determined parameters */
2405 fr->userint1 = ir->userint1;
2406 fr->userint2 = ir->userint2;
2407 fr->userint3 = ir->userint3;
2408 fr->userint4 = ir->userint4;
2409 fr->userreal1 = ir->userreal1;
2410 fr->userreal2 = ir->userreal2;
2411 fr->userreal3 = ir->userreal3;
2412 fr->userreal4 = ir->userreal4;
2415 fr->fc_stepsize = ir->fc_stepsize;
2418 fr->efep = ir->efep;
2419 fr->sc_alphavdw = ir->fepvals->sc_alpha;
2420 if (ir->fepvals->bScCoul)
2422 fr->sc_alphacoul = ir->fepvals->sc_alpha;
2423 fr->sc_sigma6_min = pow(ir->fepvals->sc_sigma_min, 6);
2427 fr->sc_alphacoul = 0;
2428 fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
2430 fr->sc_power = ir->fepvals->sc_power;
2431 fr->sc_r_power = ir->fepvals->sc_r_power;
2432 fr->sc_sigma6_def = pow(ir->fepvals->sc_sigma, 6);
2434 env = getenv("GMX_SCSIGMA_MIN");
2438 sscanf(env, "%lf", &dbl);
2439 fr->sc_sigma6_min = pow(dbl, 6);
2442 fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
2446 fr->bNonbonded = TRUE;
2447 if (getenv("GMX_NO_NONBONDED") != NULL)
2449 /* turn off non-bonded calculations */
2450 fr->bNonbonded = FALSE;
2451 md_print_warn(cr, fp,
2452 "Found environment variable GMX_NO_NONBONDED.\n"
2453 "Disabling nonbonded calculations.\n");
2456 bGenericKernelOnly = FALSE;
2458 /* We now check in the NS code whether a particular combination of interactions
2459 * can be used with water optimization, and disable it if that is not the case.
2462 if (getenv("GMX_NB_GENERIC") != NULL)
2467 "Found environment variable GMX_NB_GENERIC.\n"
2468 "Disabling all interaction-specific nonbonded kernels, will only\n"
2469 "use the slow generic ones in src/gmxlib/nonbonded/nb_generic.c\n\n");
2471 bGenericKernelOnly = TRUE;
2474 if (bGenericKernelOnly == TRUE)
2479 if ( (getenv("GMX_DISABLE_SIMD_KERNELS") != NULL) || (getenv("GMX_NOOPTIMIZEDKERNELS") != NULL) )
2481 fr->use_simd_kernels = FALSE;
2485 "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
2486 "Disabling the usage of any SIMD-specific kernel routines (e.g. SSE2/SSE4.1/AVX).\n\n");
2490 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
2492 /* Check if we can/should do all-vs-all kernels */
2493 fr->bAllvsAll = can_use_allvsall(ir, FALSE, NULL, NULL);
2494 fr->AllvsAll_work = NULL;
2495 fr->AllvsAll_workgb = NULL;
2497 /* All-vs-all kernels have not been implemented in 4.6, and
2498 * the SIMD group kernels are also buggy in this case. Non-SIMD
2499 * group kernels are OK. See Redmine #1249. */
2502 fr->bAllvsAll = FALSE;
2503 fr->use_simd_kernels = FALSE;
2507 "\nYour simulation settings would have triggered the efficient all-vs-all\n"
2508 "kernels in GROMACS 4.5, but these have not been implemented in GROMACS\n"
2509 "4.6. Also, we can't use the accelerated SIMD kernels here because\n"
2510 "of an unfixed bug. The reference C kernels are correct, though, so\n"
2511 "we are proceeding by disabling all CPU architecture-specific\n"
2512 "(e.g. SSE2/SSE4/AVX) routines. If performance is important, please\n"
2513 "use GROMACS 4.5.7 or try cutoff-scheme = Verlet.\n\n");
2517 /* Neighbour searching stuff */
2518 fr->cutoff_scheme = ir->cutoff_scheme;
2519 fr->bGrid = (ir->ns_type == ensGRID);
2520 fr->ePBC = ir->ePBC;
2522 if (fr->cutoff_scheme == ecutsGROUP)
2524 const char *note = "NOTE: This file uses the deprecated 'group' cutoff_scheme. This will be\n"
2525 "removed in a future release when 'verlet' supports all interaction forms.\n";
2529 fprintf(stderr, "\n%s\n", note);
2533 fprintf(fp, "\n%s\n", note);
2537 /* Determine if we will do PBC for distances in bonded interactions */
2538 if (fr->ePBC == epbcNONE)
2540 fr->bMolPBC = FALSE;
2544 if (!DOMAINDECOMP(cr))
2546 /* The group cut-off scheme and SHAKE assume charge groups
2547 * are whole, but not using molpbc is faster in most cases.
2549 if (fr->cutoff_scheme == ecutsGROUP ||
2550 (ir->eConstrAlg == econtSHAKE &&
2551 (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0 ||
2552 gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0)))
2554 fr->bMolPBC = ir->bPeriodicMols;
2559 if (getenv("GMX_USE_GRAPH") != NULL)
2561 fr->bMolPBC = FALSE;
2564 fprintf(fp, "\nGMX_MOLPBC is set, using the graph for bonded interactions\n\n");
2571 fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->ePBC);
2574 fr->bGB = (ir->implicit_solvent == eisGBSA);
2576 fr->rc_scaling = ir->refcoord_scaling;
2577 copy_rvec(ir->posres_com, fr->posres_com);
2578 copy_rvec(ir->posres_comB, fr->posres_comB);
2579 fr->rlist = cutoff_inf(ir->rlist);
2580 fr->rlistlong = cutoff_inf(ir->rlistlong);
2581 fr->eeltype = ir->coulombtype;
2582 fr->vdwtype = ir->vdwtype;
2583 fr->ljpme_combination_rule = ir->ljpme_combination_rule;
2585 fr->coulomb_modifier = ir->coulomb_modifier;
2586 fr->vdw_modifier = ir->vdw_modifier;
2588 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
2589 switch (fr->eeltype)
2592 fr->nbkernel_elec_interaction = (fr->bGB) ? GMX_NBKERNEL_ELEC_GENERALIZEDBORN : GMX_NBKERNEL_ELEC_COULOMB;
2598 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2602 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2603 fr->coulomb_modifier = eintmodEXACTCUTOFF;
2612 case eelPMEUSERSWITCH:
2613 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2618 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD;
2622 gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[fr->eeltype]);
2626 /* Vdw: Translate from mdp settings to kernel format */
2627 switch (fr->vdwtype)
2632 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
2636 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
2640 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD;
2646 case evdwENCADSHIFT:
2647 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2651 gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[fr->vdwtype]);
2655 /* These start out identical to ir, but might be altered if we e.g. tabulate the interaction in the kernel */
2656 fr->nbkernel_elec_modifier = fr->coulomb_modifier;
2657 fr->nbkernel_vdw_modifier = fr->vdw_modifier;
2659 fr->bTwinRange = fr->rlistlong > fr->rlist;
2660 fr->bEwald = (EEL_PME(fr->eeltype) || fr->eeltype == eelEWALD);
2662 fr->reppow = mtop->ffparams.reppow;
2664 if (ir->cutoff_scheme == ecutsGROUP)
2666 fr->bvdwtab = ((fr->vdwtype != evdwCUT || !gmx_within_tol(fr->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2667 && !EVDW_PME(fr->vdwtype));
2668 /* We have special kernels for standard Ewald and PME, but the pme-switch ones are tabulated above */
2669 fr->bcoultab = !(fr->eeltype == eelCUT ||
2670 fr->eeltype == eelEWALD ||
2671 fr->eeltype == eelPME ||
2672 fr->eeltype == eelRF ||
2673 fr->eeltype == eelRF_ZERO);
2675 /* If the user absolutely wants different switch/shift settings for coul/vdw, it is likely
2676 * going to be faster to tabulate the interaction than calling the generic kernel.
2678 if (fr->nbkernel_elec_modifier == eintmodPOTSWITCH && fr->nbkernel_vdw_modifier == eintmodPOTSWITCH)
2680 if ((fr->rcoulomb_switch != fr->rvdw_switch) || (fr->rcoulomb != fr->rvdw))
2682 fr->bcoultab = TRUE;
2685 else if ((fr->nbkernel_elec_modifier == eintmodPOTSHIFT && fr->nbkernel_vdw_modifier == eintmodPOTSHIFT) ||
2686 ((fr->nbkernel_elec_interaction == GMX_NBKERNEL_ELEC_REACTIONFIELD &&
2687 fr->nbkernel_elec_modifier == eintmodEXACTCUTOFF &&
2688 (fr->nbkernel_vdw_modifier == eintmodPOTSWITCH || fr->nbkernel_vdw_modifier == eintmodPOTSHIFT))))
2690 if (fr->rcoulomb != fr->rvdw)
2692 fr->bcoultab = TRUE;
2696 if (getenv("GMX_REQUIRE_TABLES"))
2699 fr->bcoultab = TRUE;
2704 fprintf(fp, "Table routines are used for coulomb: %s\n", bool_names[fr->bcoultab]);
2705 fprintf(fp, "Table routines are used for vdw: %s\n", bool_names[fr->bvdwtab ]);
2708 if (fr->bvdwtab == TRUE)
2710 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2711 fr->nbkernel_vdw_modifier = eintmodNONE;
2713 if (fr->bcoultab == TRUE)
2715 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2716 fr->nbkernel_elec_modifier = eintmodNONE;
2720 if (ir->cutoff_scheme == ecutsVERLET)
2722 if (!gmx_within_tol(fr->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2724 gmx_fatal(FARGS, "Cut-off scheme %S only supports LJ repulsion power 12", ecutscheme_names[ir->cutoff_scheme]);
2726 fr->bvdwtab = FALSE;
2727 fr->bcoultab = FALSE;
2730 /* Tables are used for direct ewald sum */
2733 if (EEL_PME(ir->coulombtype))
2737 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
2739 if (ir->coulombtype == eelP3M_AD)
2741 please_cite(fp, "Hockney1988");
2742 please_cite(fp, "Ballenegger2012");
2746 please_cite(fp, "Essmann95a");
2749 if (ir->ewald_geometry == eewg3DC)
2753 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry.\n");
2755 please_cite(fp, "In-Chul99a");
2758 fr->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
2759 init_ewald_tab(&(fr->ewald_table), ir, fp);
2762 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n",
2763 1/fr->ewaldcoeff_q);
2767 if (EVDW_PME(ir->vdwtype))
2771 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
2773 please_cite(fp, "Essmann95a");
2774 fr->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
2777 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n",
2778 1/fr->ewaldcoeff_lj);
2782 /* Electrostatics */
2783 fr->epsilon_r = ir->epsilon_r;
2784 fr->epsilon_rf = ir->epsilon_rf;
2785 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
2786 fr->rcoulomb_switch = ir->rcoulomb_switch;
2787 fr->rcoulomb = cutoff_inf(ir->rcoulomb);
2789 /* Parameters for generalized RF */
2793 if (fr->eeltype == eelGRF)
2795 init_generalized_rf(fp, mtop, ir, fr);
2798 fr->bF_NoVirSum = (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype) ||
2799 gmx_mtop_ftype_count(mtop, F_POSRES) > 0 ||
2800 gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 ||
2801 IR_ELEC_FIELD(*ir) ||
2802 (fr->adress_icor != eAdressICOff)
2805 if (fr->cutoff_scheme == ecutsGROUP &&
2806 ncg_mtop(mtop) > fr->cg_nalloc && !DOMAINDECOMP(cr))
2808 /* Count the total number of charge groups */
2809 fr->cg_nalloc = ncg_mtop(mtop);
2810 srenew(fr->cg_cm, fr->cg_nalloc);
2812 if (fr->shift_vec == NULL)
2814 snew(fr->shift_vec, SHIFTS);
2817 if (fr->fshift == NULL)
2819 snew(fr->fshift, SHIFTS);
2822 if (fr->nbfp == NULL)
2824 fr->ntype = mtop->ffparams.atnr;
2825 fr->nbfp = mk_nbfp(&mtop->ffparams, fr->bBHAM);
2826 if (EVDW_PME(fr->vdwtype))
2828 fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
2832 /* Copy the energy group exclusions */
2833 fr->egp_flags = ir->opts.egp_flags;
2835 /* Van der Waals stuff */
2836 fr->rvdw = cutoff_inf(ir->rvdw);
2837 fr->rvdw_switch = ir->rvdw_switch;
2838 if ((fr->vdwtype != evdwCUT) && (fr->vdwtype != evdwUSER) && !fr->bBHAM)
2840 if (fr->rvdw_switch >= fr->rvdw)
2842 gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)",
2843 fr->rvdw_switch, fr->rvdw);
2847 fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
2848 (fr->eeltype == eelSWITCH) ? "switched" : "shifted",
2849 fr->rvdw_switch, fr->rvdw);
2853 if (fr->bBHAM && EVDW_PME(fr->vdwtype))
2855 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
2858 if (fr->bBHAM && (fr->vdwtype == evdwSHIFT || fr->vdwtype == evdwSWITCH))
2860 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
2865 fprintf(fp, "Cut-off's: NS: %g Coulomb: %g %s: %g\n",
2866 fr->rlist, fr->rcoulomb, fr->bBHAM ? "BHAM" : "LJ", fr->rvdw);
2869 fr->eDispCorr = ir->eDispCorr;
2870 if (ir->eDispCorr != edispcNO)
2872 set_avcsixtwelve(fp, fr, mtop);
2877 set_bham_b_max(fp, fr, mtop);
2880 fr->gb_epsilon_solvent = ir->gb_epsilon_solvent;
2882 /* Copy the GBSA data (radius, volume and surftens for each
2883 * atomtype) from the topology atomtype section to forcerec.
2885 snew(fr->atype_radius, fr->ntype);
2886 snew(fr->atype_vol, fr->ntype);
2887 snew(fr->atype_surftens, fr->ntype);
2888 snew(fr->atype_gb_radius, fr->ntype);
2889 snew(fr->atype_S_hct, fr->ntype);
2891 if (mtop->atomtypes.nr > 0)
2893 for (i = 0; i < fr->ntype; i++)
2895 fr->atype_radius[i] = mtop->atomtypes.radius[i];
2897 for (i = 0; i < fr->ntype; i++)
2899 fr->atype_vol[i] = mtop->atomtypes.vol[i];
2901 for (i = 0; i < fr->ntype; i++)
2903 fr->atype_surftens[i] = mtop->atomtypes.surftens[i];
2905 for (i = 0; i < fr->ntype; i++)
2907 fr->atype_gb_radius[i] = mtop->atomtypes.gb_radius[i];
2909 for (i = 0; i < fr->ntype; i++)
2911 fr->atype_S_hct[i] = mtop->atomtypes.S_hct[i];
2915 /* Generate the GB table if needed */
2919 fr->gbtabscale = 2000;
2921 fr->gbtabscale = 500;
2925 fr->gbtab = make_gb_table(oenv, fr);
2927 init_gb(&fr->born, fr, ir, mtop, ir->gb_algorithm);
2929 /* Copy local gb data (for dd, this is done in dd_partition_system) */
2930 if (!DOMAINDECOMP(cr))
2932 make_local_gb(cr, fr->born, ir->gb_algorithm);
2936 /* Set the charge scaling */
2937 if (fr->epsilon_r != 0)
2939 fr->epsfac = ONE_4PI_EPS0/fr->epsilon_r;
2943 /* eps = 0 is infinite dieletric: no coulomb interactions */
2947 /* Reaction field constants */
2948 if (EEL_RF(fr->eeltype))
2950 calc_rffac(fp, fr->eeltype, fr->epsilon_r, fr->epsilon_rf,
2951 fr->rcoulomb, fr->temp, fr->zsquare, box,
2952 &fr->kappa, &fr->k_rf, &fr->c_rf);
2955 /*This now calculates sum for q and c6*/
2956 set_chargesum(fp, fr, mtop);
2958 /* if we are using LR electrostatics, and they are tabulated,
2959 * the tables will contain modified coulomb interactions.
2960 * Since we want to use the non-shifted ones for 1-4
2961 * coulombic interactions, we must have an extra set of tables.
2964 /* Construct tables.
2965 * A little unnecessary to make both vdw and coul tables sometimes,
2966 * but what the heck... */
2968 bMakeTables = fr->bcoultab || fr->bvdwtab || fr->bEwald ||
2969 (ir->eDispCorr != edispcNO && ir_vdw_switched(ir));
2971 bMakeSeparate14Table = ((!bMakeTables || fr->eeltype != eelCUT || fr->vdwtype != evdwCUT ||
2972 fr->bBHAM || fr->bEwald) &&
2973 (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
2974 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0 ||
2975 gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0));
2977 negp_pp = ir->opts.ngener - ir->nwall;
2981 bSomeNormalNbListsAreInUse = TRUE;
2986 bSomeNormalNbListsAreInUse = (ir->eDispCorr != edispcNO);
2987 for (egi = 0; egi < negp_pp; egi++)
2989 for (egj = egi; egj < negp_pp; egj++)
2991 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
2992 if (!(egp_flags & EGP_EXCL))
2994 if (egp_flags & EGP_TABLE)
3000 bSomeNormalNbListsAreInUse = TRUE;
3005 if (bSomeNormalNbListsAreInUse)
3007 fr->nnblists = negptable + 1;
3011 fr->nnblists = negptable;
3013 if (fr->nnblists > 1)
3015 snew(fr->gid2nblists, ir->opts.ngener*ir->opts.ngener);
3024 snew(fr->nblists, fr->nnblists);
3026 /* This code automatically gives table length tabext without cut-off's,
3027 * in that case grompp should already have checked that we do not need
3028 * normal tables and we only generate tables for 1-4 interactions.
3030 rtab = ir->rlistlong + ir->tabext;
3034 /* make tables for ordinary interactions */
3035 if (bSomeNormalNbListsAreInUse)
3037 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn, NULL, NULL, &fr->nblists[0]);
3040 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn, NULL, NULL, &fr->nblists[fr->nnblists/2]);
3042 if (!bMakeSeparate14Table)
3044 fr->tab14 = fr->nblists[0].table_elec_vdw;
3054 /* Read the special tables for certain energy group pairs */
3055 nm_ind = mtop->groups.grps[egcENER].nm_ind;
3056 for (egi = 0; egi < negp_pp; egi++)
3058 for (egj = egi; egj < negp_pp; egj++)
3060 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
3061 if ((egp_flags & EGP_TABLE) && !(egp_flags & EGP_EXCL))
3063 nbl = &(fr->nblists[m]);
3064 if (fr->nnblists > 1)
3066 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = m;
3068 /* Read the table file with the two energy groups names appended */
3069 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn,
3070 *mtop->groups.grpname[nm_ind[egi]],
3071 *mtop->groups.grpname[nm_ind[egj]],
3075 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn,
3076 *mtop->groups.grpname[nm_ind[egi]],
3077 *mtop->groups.grpname[nm_ind[egj]],
3078 &fr->nblists[fr->nnblists/2+m]);
3082 else if (fr->nnblists > 1)
3084 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = 0;
3090 if (bMakeSeparate14Table)
3092 /* generate extra tables with plain Coulomb for 1-4 interactions only */
3093 fr->tab14 = make_tables(fp, oenv, fr, MASTER(cr), tabpfn, rtab,
3094 GMX_MAKETABLES_14ONLY);
3097 /* Read AdResS Thermo Force table if needed */
3098 if (fr->adress_icor == eAdressICThermoForce)
3100 /* old todo replace */
3102 if (ir->adress->n_tf_grps > 0)
3104 make_adress_tf_tables(fp, oenv, fr, ir, tabfn, mtop, box);
3109 /* load the default table */
3110 snew(fr->atf_tabs, 1);
3111 fr->atf_tabs[DEFAULT_TF_TABLE] = make_atf_table(fp, oenv, fr, tabafn, box);
3116 fr->nwall = ir->nwall;
3117 if (ir->nwall && ir->wall_type == ewtTABLE)
3119 make_wall_tables(fp, oenv, ir, tabfn, &mtop->groups, fr);
3124 fcd->bondtab = make_bonded_tables(fp,
3125 F_TABBONDS, F_TABBONDSNC,
3127 fcd->angletab = make_bonded_tables(fp,
3130 fcd->dihtab = make_bonded_tables(fp,
3138 fprintf(debug, "No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
3142 /* QM/MM initialization if requested
3146 fprintf(stderr, "QM/MM calculation requested.\n");
3149 fr->bQMMM = ir->bQMMM;
3150 fr->qr = mk_QMMMrec();
3152 /* Set all the static charge group info */
3153 fr->cginfo_mb = init_cginfo_mb(fp, mtop, fr, bNoSolvOpt,
3155 &fr->bExcl_IntraCGAll_InterCGNone);
3156 if (DOMAINDECOMP(cr))
3162 fr->cginfo = cginfo_expand(mtop->nmolblock, fr->cginfo_mb);
3165 if (!DOMAINDECOMP(cr))
3167 forcerec_set_ranges(fr, ncg_mtop(mtop), ncg_mtop(mtop),
3168 mtop->natoms, mtop->natoms, mtop->natoms);
3171 fr->print_force = print_force;
3174 /* coarse load balancing vars */
3179 /* Initialize neighbor search */
3180 init_ns(fp, cr, &fr->ns, fr, mtop);
3182 if (cr->duty & DUTY_PP)
3184 gmx_nonbonded_setup(fr, bGenericKernelOnly);
3188 gmx_setup_adress_kernels(fp,bGenericKernelOnly);
3193 /* Initialize the thread working data for bonded interactions */
3194 init_forcerec_f_threads(fr, mtop->groups.grps[egcENER].nr);
3196 snew(fr->excl_load, fr->nthreads+1);
3198 if (fr->cutoff_scheme == ecutsVERLET)
3200 if (ir->rcoulomb != ir->rvdw)
3202 gmx_fatal(FARGS, "With Verlet lists rcoulomb and rvdw should be identical");
3205 init_nb_verlet(fp, &fr->nbv, bFEP_NonBonded, ir, fr, cr, nbpu_opt);
3208 /* fr->ic is used both by verlet and group kernels (to some extent) now */
3209 init_interaction_const(fp, cr, &fr->ic, fr, rtab);
3211 if (ir->eDispCorr != edispcNO)
3213 calc_enervirdiff(fp, ir->eDispCorr, fr);
3217 #define pr_real(fp, r) fprintf(fp, "%s: %e\n",#r, r)
3218 #define pr_int(fp, i) fprintf((fp), "%s: %d\n",#i, i)
3219 #define pr_bool(fp, b) fprintf((fp), "%s: %s\n",#b, bool_names[b])
3221 void pr_forcerec(FILE *fp, t_forcerec *fr)
3225 pr_real(fp, fr->rlist);
3226 pr_real(fp, fr->rcoulomb);
3227 pr_real(fp, fr->fudgeQQ);
3228 pr_bool(fp, fr->bGrid);
3229 pr_bool(fp, fr->bTwinRange);
3230 /*pr_int(fp,fr->cg0);
3231 pr_int(fp,fr->hcg);*/
3232 for (i = 0; i < fr->nnblists; i++)
3234 pr_int(fp, fr->nblists[i].table_elec_vdw.n);
3236 pr_real(fp, fr->rcoulomb_switch);
3237 pr_real(fp, fr->rcoulomb);
3242 void forcerec_set_excl_load(t_forcerec *fr,
3243 const gmx_localtop_t *top)
3246 int t, i, j, ntot, n, ntarget;
3248 ind = top->excls.index;
3252 for (i = 0; i < top->excls.nr; i++)
3254 for (j = ind[i]; j < ind[i+1]; j++)
3263 fr->excl_load[0] = 0;
3266 for (t = 1; t <= fr->nthreads; t++)
3268 ntarget = (ntot*t)/fr->nthreads;
3269 while (i < top->excls.nr && n < ntarget)
3271 for (j = ind[i]; j < ind[i+1]; j++)
3280 fr->excl_load[t] = i;