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.
47 #include "gromacs/math/utilities.h"
51 #include "gmx_fatal.h"
52 #include "gmx_fatal_collective.h"
56 #include "nonbonded.h"
65 #include "md_support.h"
66 #include "md_logging.h"
71 #include "mtop_util.h"
72 #include "nbnxn_search.h"
73 #include "nbnxn_atomdata.h"
74 #include "nbnxn_consts.h"
75 #include "gmx_omp_nthreads.h"
76 #include "gmx_detect_hardware.h"
79 /* MSVC definition for __cpuid() */
83 #include "types/nbnxn_cuda_types_ext.h"
84 #include "gpu_utils.h"
85 #include "nbnxn_cuda_data_mgmt.h"
86 #include "pmalloc_cuda.h"
88 t_forcerec *mk_forcerec(void)
98 static void pr_nbfp(FILE *fp, real *nbfp, gmx_bool bBHAM, int atnr)
102 for (i = 0; (i < atnr); i++)
104 for (j = 0; (j < atnr); j++)
106 fprintf(fp, "%2d - %2d", i, j);
109 fprintf(fp, " a=%10g, b=%10g, c=%10g\n", BHAMA(nbfp, atnr, i, j),
110 BHAMB(nbfp, atnr, i, j), BHAMC(nbfp, atnr, i, j)/6.0);
114 fprintf(fp, " c6=%10g, c12=%10g\n", C6(nbfp, atnr, i, j)/6.0,
115 C12(nbfp, atnr, i, j)/12.0);
122 static real *mk_nbfp(const gmx_ffparams_t *idef, gmx_bool bBHAM)
130 snew(nbfp, 3*atnr*atnr);
131 for (i = k = 0; (i < atnr); i++)
133 for (j = 0; (j < atnr); j++, k++)
135 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
136 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
137 /* nbfp now includes the 6.0 derivative prefactor */
138 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c*6.0;
144 snew(nbfp, 2*atnr*atnr);
145 for (i = k = 0; (i < atnr); i++)
147 for (j = 0; (j < atnr); j++, k++)
149 /* nbfp now includes the 6.0/12.0 derivative prefactors */
150 C6(nbfp, atnr, i, j) = idef->iparams[k].lj.c6*6.0;
151 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12*12.0;
159 static real *mk_nbfp_combination_rule(const gmx_ffparams_t *idef, int comb_rule)
163 real c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
167 snew(nbfp, 2*atnr*atnr);
168 for (i = 0; i < atnr; ++i)
170 for (j = 0; j < atnr; ++j)
172 c6i = idef->iparams[i*(atnr+1)].lj.c6;
173 c12i = idef->iparams[i*(atnr+1)].lj.c12;
174 c6j = idef->iparams[j*(atnr+1)].lj.c6;
175 c12j = idef->iparams[j*(atnr+1)].lj.c12;
176 c6 = sqrt(c6i * c6j);
177 c12 = sqrt(c12i * c12j);
178 if (comb_rule == eCOMB_ARITHMETIC
179 && !gmx_numzero(c6) && !gmx_numzero(c12))
181 sigmai = pow(c12i / c6i, 1.0/6.0);
182 sigmaj = pow(c12j / c6j, 1.0/6.0);
183 epsi = c6i * c6i / c12i;
184 epsj = c6j * c6j / c12j;
185 c6 = epsi * epsj * pow(0.5*(sigmai+sigmaj), 6);
186 c12 = epsi * epsj * pow(0.5*(sigmai+sigmaj), 12);
188 C6(nbfp, atnr, i, j) = c6*6.0;
189 C12(nbfp, atnr, i, j) = c12*12.0;
195 /* This routine sets fr->solvent_opt to the most common solvent in the
196 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
197 * the fr->solvent_type array with the correct type (or esolNO).
199 * Charge groups that fulfill the conditions but are not identical to the
200 * most common one will be marked as esolNO in the solvent_type array.
202 * TIP3p is identical to SPC for these purposes, so we call it
203 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
205 * NOTE: QM particle should not
206 * become an optimized solvent. Not even if there is only one charge
216 } solvent_parameters_t;
219 check_solvent_cg(const gmx_moltype_t *molt,
222 const unsigned char *qm_grpnr,
223 const t_grps *qm_grps,
225 int *n_solvent_parameters,
226 solvent_parameters_t **solvent_parameters_p,
230 const t_blocka * excl;
241 solvent_parameters_t *solvent_parameters;
243 /* We use a list with parameters for each solvent type.
244 * Every time we discover a new molecule that fulfills the basic
245 * conditions for a solvent we compare with the previous entries
246 * in these lists. If the parameters are the same we just increment
247 * the counter for that type, and otherwise we create a new type
248 * based on the current molecule.
250 * Once we've finished going through all molecules we check which
251 * solvent is most common, and mark all those molecules while we
252 * clear the flag on all others.
255 solvent_parameters = *solvent_parameters_p;
257 /* Mark the cg first as non optimized */
260 /* Check if this cg has no exclusions with atoms in other charge groups
261 * and all atoms inside the charge group excluded.
262 * We only have 3 or 4 atom solvent loops.
264 if (GET_CGINFO_EXCL_INTER(cginfo) ||
265 !GET_CGINFO_EXCL_INTRA(cginfo))
270 /* Get the indices of the first atom in this charge group */
271 j0 = molt->cgs.index[cg0];
272 j1 = molt->cgs.index[cg0+1];
274 /* Number of atoms in our molecule */
280 "Moltype '%s': there are %d atoms in this charge group\n",
284 /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
287 if (nj < 3 || nj > 4)
292 /* Check if we are doing QM on this group */
294 if (qm_grpnr != NULL)
296 for (j = j0; j < j1 && !qm; j++)
298 qm = (qm_grpnr[j] < qm_grps->nr - 1);
301 /* Cannot use solvent optimization with QM */
307 atom = molt->atoms.atom;
309 /* Still looks like a solvent, time to check parameters */
311 /* If it is perturbed (free energy) we can't use the solvent loops,
312 * so then we just skip to the next molecule.
316 for (j = j0; j < j1 && !perturbed; j++)
318 perturbed = PERTURBED(atom[j]);
326 /* Now it's only a question if the VdW and charge parameters
327 * are OK. Before doing the check we compare and see if they are
328 * identical to a possible previous solvent type.
329 * First we assign the current types and charges.
331 for (j = 0; j < nj; j++)
333 tmp_vdwtype[j] = atom[j0+j].type;
334 tmp_charge[j] = atom[j0+j].q;
337 /* Does it match any previous solvent type? */
338 for (k = 0; k < *n_solvent_parameters; k++)
343 /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
344 if ( (solvent_parameters[k].model == esolSPC && nj != 3) ||
345 (solvent_parameters[k].model == esolTIP4P && nj != 4) )
350 /* Check that types & charges match for all atoms in molecule */
351 for (j = 0; j < nj && match == TRUE; j++)
353 if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
357 if (tmp_charge[j] != solvent_parameters[k].charge[j])
364 /* Congratulations! We have a matched solvent.
365 * Flag it with this type for later processing.
368 solvent_parameters[k].count += nmol;
370 /* We are done with this charge group */
375 /* If we get here, we have a tentative new solvent type.
376 * Before we add it we must check that it fulfills the requirements
377 * of the solvent optimized loops. First determine which atoms have
380 for (j = 0; j < nj; j++)
383 tjA = tmp_vdwtype[j];
385 /* Go through all other tpes and see if any have non-zero
386 * VdW parameters when combined with this one.
388 for (k = 0; k < fr->ntype && (has_vdw[j] == FALSE); k++)
390 /* We already checked that the atoms weren't perturbed,
391 * so we only need to check state A now.
395 has_vdw[j] = (has_vdw[j] ||
396 (BHAMA(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
397 (BHAMB(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
398 (BHAMC(fr->nbfp, fr->ntype, tjA, k) != 0.0));
403 has_vdw[j] = (has_vdw[j] ||
404 (C6(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
405 (C12(fr->nbfp, fr->ntype, tjA, k) != 0.0));
410 /* Now we know all we need to make the final check and assignment. */
414 * For this we require thatn all atoms have charge,
415 * the charges on atom 2 & 3 should be the same, and only
416 * atom 1 might have VdW.
418 if (has_vdw[1] == FALSE &&
419 has_vdw[2] == FALSE &&
420 tmp_charge[0] != 0 &&
421 tmp_charge[1] != 0 &&
422 tmp_charge[2] == tmp_charge[1])
424 srenew(solvent_parameters, *n_solvent_parameters+1);
425 solvent_parameters[*n_solvent_parameters].model = esolSPC;
426 solvent_parameters[*n_solvent_parameters].count = nmol;
427 for (k = 0; k < 3; k++)
429 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
430 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
433 *cg_sp = *n_solvent_parameters;
434 (*n_solvent_parameters)++;
439 /* Or could it be a TIP4P?
440 * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
441 * Only atom 1 mght have VdW.
443 if (has_vdw[1] == FALSE &&
444 has_vdw[2] == FALSE &&
445 has_vdw[3] == FALSE &&
446 tmp_charge[0] == 0 &&
447 tmp_charge[1] != 0 &&
448 tmp_charge[2] == tmp_charge[1] &&
451 srenew(solvent_parameters, *n_solvent_parameters+1);
452 solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
453 solvent_parameters[*n_solvent_parameters].count = nmol;
454 for (k = 0; k < 4; k++)
456 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
457 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
460 *cg_sp = *n_solvent_parameters;
461 (*n_solvent_parameters)++;
465 *solvent_parameters_p = solvent_parameters;
469 check_solvent(FILE * fp,
470 const gmx_mtop_t * mtop,
472 cginfo_mb_t *cginfo_mb)
475 const t_block * mols;
476 const gmx_moltype_t *molt;
477 int mb, mol, cg_mol, at_offset, cg_offset, am, cgm, i, nmol_ch, nmol;
478 int n_solvent_parameters;
479 solvent_parameters_t *solvent_parameters;
485 fprintf(debug, "Going to determine what solvent types we have.\n");
490 n_solvent_parameters = 0;
491 solvent_parameters = NULL;
492 /* Allocate temporary array for solvent type */
493 snew(cg_sp, mtop->nmolblock);
497 for (mb = 0; mb < mtop->nmolblock; mb++)
499 molt = &mtop->moltype[mtop->molblock[mb].type];
501 /* Here we have to loop over all individual molecules
502 * because we need to check for QMMM particles.
504 snew(cg_sp[mb], cginfo_mb[mb].cg_mod);
505 nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
506 nmol = mtop->molblock[mb].nmol/nmol_ch;
507 for (mol = 0; mol < nmol_ch; mol++)
510 am = mol*cgs->index[cgs->nr];
511 for (cg_mol = 0; cg_mol < cgs->nr; cg_mol++)
513 check_solvent_cg(molt, cg_mol, nmol,
514 mtop->groups.grpnr[egcQMMM] ?
515 mtop->groups.grpnr[egcQMMM]+at_offset+am : 0,
516 &mtop->groups.grps[egcQMMM],
518 &n_solvent_parameters, &solvent_parameters,
519 cginfo_mb[mb].cginfo[cgm+cg_mol],
520 &cg_sp[mb][cgm+cg_mol]);
523 cg_offset += cgs->nr;
524 at_offset += cgs->index[cgs->nr];
527 /* Puh! We finished going through all charge groups.
528 * Now find the most common solvent model.
531 /* Most common solvent this far */
533 for (i = 0; i < n_solvent_parameters; i++)
536 solvent_parameters[i].count > solvent_parameters[bestsp].count)
544 bestsol = solvent_parameters[bestsp].model;
551 #ifdef DISABLE_WATER_NLIST
556 for (mb = 0; mb < mtop->nmolblock; mb++)
558 cgs = &mtop->moltype[mtop->molblock[mb].type].cgs;
559 nmol = (mtop->molblock[mb].nmol*cgs->nr)/cginfo_mb[mb].cg_mod;
560 for (i = 0; i < cginfo_mb[mb].cg_mod; i++)
562 if (cg_sp[mb][i] == bestsp)
564 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], bestsol);
569 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], esolNO);
576 if (bestsol != esolNO && fp != NULL)
578 fprintf(fp, "\nEnabling %s-like water optimization for %d molecules.\n\n",
580 solvent_parameters[bestsp].count);
583 sfree(solvent_parameters);
584 fr->solvent_opt = bestsol;
588 acNONE = 0, acCONSTRAINT, acSETTLE
591 static cginfo_mb_t *init_cginfo_mb(FILE *fplog, const gmx_mtop_t *mtop,
592 t_forcerec *fr, gmx_bool bNoSolvOpt,
593 gmx_bool *bExcl_IntraCGAll_InterCGNone)
596 const t_blocka *excl;
597 const gmx_moltype_t *molt;
598 const gmx_molblock_t *molb;
599 cginfo_mb_t *cginfo_mb;
602 int cg_offset, a_offset, cgm, am;
603 int mb, m, ncg_tot, cg, a0, a1, gid, ai, j, aj, excl_nalloc;
607 gmx_bool bId, *bExcl, bExclIntraAll, bExclInter, bHaveVDW, bHaveQ;
609 ncg_tot = ncg_mtop(mtop);
610 snew(cginfo_mb, mtop->nmolblock);
612 snew(type_VDW, fr->ntype);
613 for (ai = 0; ai < fr->ntype; ai++)
615 type_VDW[ai] = FALSE;
616 for (j = 0; j < fr->ntype; j++)
618 type_VDW[ai] = type_VDW[ai] ||
620 C6(fr->nbfp, fr->ntype, ai, j) != 0 ||
621 C12(fr->nbfp, fr->ntype, ai, j) != 0;
625 *bExcl_IntraCGAll_InterCGNone = TRUE;
628 snew(bExcl, excl_nalloc);
631 for (mb = 0; mb < mtop->nmolblock; mb++)
633 molb = &mtop->molblock[mb];
634 molt = &mtop->moltype[molb->type];
638 /* Check if the cginfo is identical for all molecules in this block.
639 * If so, we only need an array of the size of one molecule.
640 * Otherwise we make an array of #mol times #cgs per molecule.
644 for (m = 0; m < molb->nmol; m++)
646 am = m*cgs->index[cgs->nr];
647 for (cg = 0; cg < cgs->nr; cg++)
650 a1 = cgs->index[cg+1];
651 if (ggrpnr(&mtop->groups, egcENER, a_offset+am+a0) !=
652 ggrpnr(&mtop->groups, egcENER, a_offset +a0))
656 if (mtop->groups.grpnr[egcQMMM] != NULL)
658 for (ai = a0; ai < a1; ai++)
660 if (mtop->groups.grpnr[egcQMMM][a_offset+am+ai] !=
661 mtop->groups.grpnr[egcQMMM][a_offset +ai])
670 cginfo_mb[mb].cg_start = cg_offset;
671 cginfo_mb[mb].cg_end = cg_offset + molb->nmol*cgs->nr;
672 cginfo_mb[mb].cg_mod = (bId ? 1 : molb->nmol)*cgs->nr;
673 snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
674 cginfo = cginfo_mb[mb].cginfo;
676 /* Set constraints flags for constrained atoms */
677 snew(a_con, molt->atoms.nr);
678 for (ftype = 0; ftype < F_NRE; ftype++)
680 if (interaction_function[ftype].flags & IF_CONSTRAINT)
685 for (ia = 0; ia < molt->ilist[ftype].nr; ia += 1+nral)
689 for (a = 0; a < nral; a++)
691 a_con[molt->ilist[ftype].iatoms[ia+1+a]] =
692 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
698 for (m = 0; m < (bId ? 1 : molb->nmol); m++)
701 am = m*cgs->index[cgs->nr];
702 for (cg = 0; cg < cgs->nr; cg++)
705 a1 = cgs->index[cg+1];
707 /* Store the energy group in cginfo */
708 gid = ggrpnr(&mtop->groups, egcENER, a_offset+am+a0);
709 SET_CGINFO_GID(cginfo[cgm+cg], gid);
711 /* Check the intra/inter charge group exclusions */
712 if (a1-a0 > excl_nalloc)
714 excl_nalloc = a1 - a0;
715 srenew(bExcl, excl_nalloc);
717 /* bExclIntraAll: all intra cg interactions excluded
718 * bExclInter: any inter cg interactions excluded
720 bExclIntraAll = TRUE;
724 for (ai = a0; ai < a1; ai++)
726 /* Check VDW and electrostatic interactions */
727 bHaveVDW = bHaveVDW || (type_VDW[molt->atoms.atom[ai].type] ||
728 type_VDW[molt->atoms.atom[ai].typeB]);
729 bHaveQ = bHaveQ || (molt->atoms.atom[ai].q != 0 ||
730 molt->atoms.atom[ai].qB != 0);
732 /* Clear the exclusion list for atom ai */
733 for (aj = a0; aj < a1; aj++)
735 bExcl[aj-a0] = FALSE;
737 /* Loop over all the exclusions of atom ai */
738 for (j = excl->index[ai]; j < excl->index[ai+1]; j++)
741 if (aj < a0 || aj >= a1)
750 /* Check if ai excludes a0 to a1 */
751 for (aj = a0; aj < a1; aj++)
755 bExclIntraAll = FALSE;
762 SET_CGINFO_CONSTR(cginfo[cgm+cg]);
765 SET_CGINFO_SETTLE(cginfo[cgm+cg]);
773 SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
777 SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
779 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
781 /* The size in cginfo is currently only read with DD */
782 gmx_fatal(FARGS, "A charge group has size %d which is larger than the limit of %d atoms", a1-a0, MAX_CHARGEGROUP_SIZE);
786 SET_CGINFO_HAS_VDW(cginfo[cgm+cg]);
790 SET_CGINFO_HAS_Q(cginfo[cgm+cg]);
792 /* Store the charge group size */
793 SET_CGINFO_NATOMS(cginfo[cgm+cg], a1-a0);
795 if (!bExclIntraAll || bExclInter)
797 *bExcl_IntraCGAll_InterCGNone = FALSE;
804 cg_offset += molb->nmol*cgs->nr;
805 a_offset += molb->nmol*cgs->index[cgs->nr];
809 /* the solvent optimizer is called after the QM is initialized,
810 * because we don't want to have the QM subsystemto become an
814 check_solvent(fplog, mtop, fr, cginfo_mb);
816 if (getenv("GMX_NO_SOLV_OPT"))
820 fprintf(fplog, "Found environment variable GMX_NO_SOLV_OPT.\n"
821 "Disabling all solvent optimization\n");
823 fr->solvent_opt = esolNO;
827 fr->solvent_opt = esolNO;
829 if (!fr->solvent_opt)
831 for (mb = 0; mb < mtop->nmolblock; mb++)
833 for (cg = 0; cg < cginfo_mb[mb].cg_mod; cg++)
835 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg], esolNO);
843 static int *cginfo_expand(int nmb, cginfo_mb_t *cgi_mb)
848 ncg = cgi_mb[nmb-1].cg_end;
851 for (cg = 0; cg < ncg; cg++)
853 while (cg >= cgi_mb[mb].cg_end)
858 cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
864 static void set_chargesum(FILE *log, t_forcerec *fr, const gmx_mtop_t *mtop)
866 /*This now calculates sum for q and c6*/
867 double qsum, q2sum, q, c6sum, c6, sqrc6;
869 const t_atoms *atoms;
874 for (mb = 0; mb < mtop->nmolblock; mb++)
876 nmol = mtop->molblock[mb].nmol;
877 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
878 for (i = 0; i < atoms->nr; i++)
880 q = atoms->atom[i].q;
883 c6 = mtop->ffparams.iparams[atoms->atom[i].type*(mtop->ffparams.atnr+1)].lj.c6;
885 c6sum += nmol*sqrc6*sqrc6;
889 fr->q2sum[0] = q2sum;
890 fr->c6sum[0] = c6sum/12;
892 if (fr->efep != efepNO)
897 for (mb = 0; mb < mtop->nmolblock; mb++)
899 nmol = mtop->molblock[mb].nmol;
900 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
901 for (i = 0; i < atoms->nr; i++)
903 q = atoms->atom[i].qB;
906 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB*(mtop->ffparams.atnr+1)].lj.c6;
908 c6sum += nmol*sqrc6*sqrc6;
911 fr->q2sum[1] = q2sum;
912 fr->c6sum[1] = c6sum/12;
917 fr->qsum[1] = fr->qsum[0];
918 fr->q2sum[1] = fr->q2sum[0];
919 fr->c6sum[1] = fr->c6sum[0];
923 if (fr->efep == efepNO)
925 fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
929 fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n",
930 fr->qsum[0], fr->qsum[1]);
935 void update_forcerec(t_forcerec *fr, matrix box)
937 if (fr->eeltype == eelGRF)
939 calc_rffac(NULL, fr->eeltype, fr->epsilon_r, fr->epsilon_rf,
940 fr->rcoulomb, fr->temp, fr->zsquare, box,
941 &fr->kappa, &fr->k_rf, &fr->c_rf);
945 void set_avcsixtwelve(FILE *fplog, t_forcerec *fr, const gmx_mtop_t *mtop)
947 const t_atoms *atoms, *atoms_tpi;
948 const t_blocka *excl;
949 int mb, nmol, nmolc, i, j, tpi, tpj, j1, j2, k, n, nexcl, q;
950 gmx_int64_t npair, npair_ij, tmpi, tmpj;
951 double csix, ctwelve;
955 real *nbfp_comb = NULL;
961 /* For LJ-PME, we want to correct for the difference between the
962 * actual C6 values and the C6 values used by the LJ-PME based on
963 * combination rules. */
965 if (EVDW_PME(fr->vdwtype))
967 nbfp_comb = mk_nbfp_combination_rule(&mtop->ffparams,
968 (fr->ljpme_combination_rule == eljpmeLB) ? eCOMB_ARITHMETIC : eCOMB_GEOMETRIC);
969 for (tpi = 0; tpi < ntp; ++tpi)
971 for (tpj = 0; tpj < ntp; ++tpj)
973 C6(nbfp_comb, ntp, tpi, tpj) =
974 C6(nbfp, ntp, tpi, tpj) - C6(nbfp_comb, ntp, tpi, tpj);
975 C12(nbfp_comb, ntp, tpi, tpj) = C12(nbfp, ntp, tpi, tpj);
980 for (q = 0; q < (fr->efep == efepNO ? 1 : 2); q++)
988 /* Count the types so we avoid natoms^2 operations */
989 snew(typecount, ntp);
990 gmx_mtop_count_atomtypes(mtop, q, typecount);
992 for (tpi = 0; tpi < ntp; tpi++)
994 for (tpj = tpi; tpj < ntp; tpj++)
996 tmpi = typecount[tpi];
997 tmpj = typecount[tpj];
1000 npair_ij = tmpi*tmpj;
1004 npair_ij = tmpi*(tmpi - 1)/2;
1008 /* nbfp now includes the 6.0 derivative prefactor */
1009 csix += npair_ij*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1013 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1014 csix += npair_ij* C6(nbfp, ntp, tpi, tpj)/6.0;
1015 ctwelve += npair_ij* C12(nbfp, ntp, tpi, tpj)/12.0;
1021 /* Subtract the excluded pairs.
1022 * The main reason for substracting exclusions is that in some cases
1023 * some combinations might never occur and the parameters could have
1024 * any value. These unused values should not influence the dispersion
1027 for (mb = 0; mb < mtop->nmolblock; mb++)
1029 nmol = mtop->molblock[mb].nmol;
1030 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
1031 excl = &mtop->moltype[mtop->molblock[mb].type].excls;
1032 for (i = 0; (i < atoms->nr); i++)
1036 tpi = atoms->atom[i].type;
1040 tpi = atoms->atom[i].typeB;
1042 j1 = excl->index[i];
1043 j2 = excl->index[i+1];
1044 for (j = j1; j < j2; j++)
1051 tpj = atoms->atom[k].type;
1055 tpj = atoms->atom[k].typeB;
1059 /* nbfp now includes the 6.0 derivative prefactor */
1060 csix -= nmol*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1064 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1065 csix -= nmol*C6 (nbfp, ntp, tpi, tpj)/6.0;
1066 ctwelve -= nmol*C12(nbfp, ntp, tpi, tpj)/12.0;
1076 /* Only correct for the interaction of the test particle
1077 * with the rest of the system.
1080 &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].atoms;
1083 for (mb = 0; mb < mtop->nmolblock; mb++)
1085 nmol = mtop->molblock[mb].nmol;
1086 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
1087 for (j = 0; j < atoms->nr; j++)
1090 /* Remove the interaction of the test charge group
1093 if (mb == mtop->nmolblock-1)
1097 if (mb == 0 && nmol == 1)
1099 gmx_fatal(FARGS, "Old format tpr with TPI, please generate a new tpr file");
1104 tpj = atoms->atom[j].type;
1108 tpj = atoms->atom[j].typeB;
1110 for (i = 0; i < fr->n_tpi; i++)
1114 tpi = atoms_tpi->atom[i].type;
1118 tpi = atoms_tpi->atom[i].typeB;
1122 /* nbfp now includes the 6.0 derivative prefactor */
1123 csix += nmolc*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1127 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1128 csix += nmolc*C6 (nbfp, ntp, tpi, tpj)/6.0;
1129 ctwelve += nmolc*C12(nbfp, ntp, tpi, tpj)/12.0;
1136 if (npair - nexcl <= 0 && fplog)
1138 fprintf(fplog, "\nWARNING: There are no atom pairs for dispersion correction\n\n");
1144 csix /= npair - nexcl;
1145 ctwelve /= npair - nexcl;
1149 fprintf(debug, "Counted %d exclusions\n", nexcl);
1150 fprintf(debug, "Average C6 parameter is: %10g\n", (double)csix);
1151 fprintf(debug, "Average C12 parameter is: %10g\n", (double)ctwelve);
1153 fr->avcsix[q] = csix;
1154 fr->avctwelve[q] = ctwelve;
1157 if (EVDW_PME(fr->vdwtype))
1164 if (fr->eDispCorr == edispcAllEner ||
1165 fr->eDispCorr == edispcAllEnerPres)
1167 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1168 fr->avcsix[0], fr->avctwelve[0]);
1172 fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e\n", fr->avcsix[0]);
1178 static void set_bham_b_max(FILE *fplog, t_forcerec *fr,
1179 const gmx_mtop_t *mtop)
1181 const t_atoms *at1, *at2;
1182 int mt1, mt2, i, j, tpi, tpj, ntypes;
1188 fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
1195 for (mt1 = 0; mt1 < mtop->nmoltype; mt1++)
1197 at1 = &mtop->moltype[mt1].atoms;
1198 for (i = 0; (i < at1->nr); i++)
1200 tpi = at1->atom[i].type;
1203 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
1206 for (mt2 = mt1; mt2 < mtop->nmoltype; mt2++)
1208 at2 = &mtop->moltype[mt2].atoms;
1209 for (j = 0; (j < at2->nr); j++)
1211 tpj = at2->atom[j].type;
1214 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
1216 b = BHAMB(nbfp, ntypes, tpi, tpj);
1217 if (b > fr->bham_b_max)
1221 if ((b < bmin) || (bmin == -1))
1231 fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n",
1232 bmin, fr->bham_b_max);
1236 static void make_nbf_tables(FILE *fp, const output_env_t oenv,
1237 t_forcerec *fr, real rtab,
1238 const t_commrec *cr,
1239 const char *tabfn, char *eg1, char *eg2,
1249 fprintf(debug, "No table file name passed, can not read table, can not do non-bonded interactions\n");
1254 sprintf(buf, "%s", tabfn);
1257 /* Append the two energy group names */
1258 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
1259 eg1, eg2, ftp2ext(efXVG));
1261 nbl->table_elec_vdw = make_tables(fp, oenv, fr, MASTER(cr), buf, rtab, 0);
1262 /* Copy the contents of the table to separate coulomb and LJ tables too,
1263 * to improve cache performance.
1265 /* For performance reasons we want
1266 * the table data to be aligned to 16-byte. The pointers could be freed
1267 * but currently aren't.
1269 nbl->table_elec.interaction = GMX_TABLE_INTERACTION_ELEC;
1270 nbl->table_elec.format = nbl->table_elec_vdw.format;
1271 nbl->table_elec.r = nbl->table_elec_vdw.r;
1272 nbl->table_elec.n = nbl->table_elec_vdw.n;
1273 nbl->table_elec.scale = nbl->table_elec_vdw.scale;
1274 nbl->table_elec.scale_exp = nbl->table_elec_vdw.scale_exp;
1275 nbl->table_elec.formatsize = nbl->table_elec_vdw.formatsize;
1276 nbl->table_elec.ninteractions = 1;
1277 nbl->table_elec.stride = nbl->table_elec.formatsize * nbl->table_elec.ninteractions;
1278 snew_aligned(nbl->table_elec.data, nbl->table_elec.stride*(nbl->table_elec.n+1), 32);
1280 nbl->table_vdw.interaction = GMX_TABLE_INTERACTION_VDWREP_VDWDISP;
1281 nbl->table_vdw.format = nbl->table_elec_vdw.format;
1282 nbl->table_vdw.r = nbl->table_elec_vdw.r;
1283 nbl->table_vdw.n = nbl->table_elec_vdw.n;
1284 nbl->table_vdw.scale = nbl->table_elec_vdw.scale;
1285 nbl->table_vdw.scale_exp = nbl->table_elec_vdw.scale_exp;
1286 nbl->table_vdw.formatsize = nbl->table_elec_vdw.formatsize;
1287 nbl->table_vdw.ninteractions = 2;
1288 nbl->table_vdw.stride = nbl->table_vdw.formatsize * nbl->table_vdw.ninteractions;
1289 snew_aligned(nbl->table_vdw.data, nbl->table_vdw.stride*(nbl->table_vdw.n+1), 32);
1291 for (i = 0; i <= nbl->table_elec_vdw.n; i++)
1293 for (j = 0; j < 4; j++)
1295 nbl->table_elec.data[4*i+j] = nbl->table_elec_vdw.data[12*i+j];
1297 for (j = 0; j < 8; j++)
1299 nbl->table_vdw.data[8*i+j] = nbl->table_elec_vdw.data[12*i+4+j];
1304 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t *mtop,
1305 int *ncount, int **count)
1307 const gmx_moltype_t *molt;
1309 int mt, ftype, stride, i, j, tabnr;
1311 for (mt = 0; mt < mtop->nmoltype; mt++)
1313 molt = &mtop->moltype[mt];
1314 for (ftype = 0; ftype < F_NRE; ftype++)
1316 if (ftype == ftype1 || ftype == ftype2)
1318 il = &molt->ilist[ftype];
1319 stride = 1 + NRAL(ftype);
1320 for (i = 0; i < il->nr; i += stride)
1322 tabnr = mtop->ffparams.iparams[il->iatoms[i]].tab.table;
1325 gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
1327 if (tabnr >= *ncount)
1329 srenew(*count, tabnr+1);
1330 for (j = *ncount; j < tabnr+1; j++)
1343 static bondedtable_t *make_bonded_tables(FILE *fplog,
1344 int ftype1, int ftype2,
1345 const gmx_mtop_t *mtop,
1346 const char *basefn, const char *tabext)
1348 int i, ncount, *count;
1356 count_tables(ftype1, ftype2, mtop, &ncount, &count);
1361 for (i = 0; i < ncount; i++)
1365 sprintf(tabfn, "%s", basefn);
1366 sprintf(tabfn + strlen(basefn) - strlen(ftp2ext(efXVG)) - 1, "_%s%d.%s",
1367 tabext, i, ftp2ext(efXVG));
1368 tab[i] = make_bonded_table(fplog, tabfn, NRAL(ftype1)-2);
1377 void forcerec_set_ranges(t_forcerec *fr,
1378 int ncg_home, int ncg_force,
1380 int natoms_force_constr, int natoms_f_novirsum)
1385 /* fr->ncg_force is unused in the standard code,
1386 * but it can be useful for modified code dealing with charge groups.
1388 fr->ncg_force = ncg_force;
1389 fr->natoms_force = natoms_force;
1390 fr->natoms_force_constr = natoms_force_constr;
1392 if (fr->natoms_force_constr > fr->nalloc_force)
1394 fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1398 srenew(fr->f_twin, fr->nalloc_force);
1402 if (fr->bF_NoVirSum)
1404 fr->f_novirsum_n = natoms_f_novirsum;
1405 if (fr->f_novirsum_n > fr->f_novirsum_nalloc)
1407 fr->f_novirsum_nalloc = over_alloc_dd(fr->f_novirsum_n);
1408 srenew(fr->f_novirsum_alloc, fr->f_novirsum_nalloc);
1413 fr->f_novirsum_n = 0;
1417 static real cutoff_inf(real cutoff)
1421 cutoff = GMX_CUTOFF_INF;
1427 static void make_adress_tf_tables(FILE *fp, const output_env_t oenv,
1428 t_forcerec *fr, const t_inputrec *ir,
1429 const char *tabfn, const gmx_mtop_t *mtop,
1437 gmx_fatal(FARGS, "No thermoforce table file given. Use -tabletf to specify a file\n");
1441 snew(fr->atf_tabs, ir->adress->n_tf_grps);
1443 sprintf(buf, "%s", tabfn);
1444 for (i = 0; i < ir->adress->n_tf_grps; i++)
1446 j = ir->adress->tf_table_index[i]; /* get energy group index */
1447 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "tf_%s.%s",
1448 *(mtop->groups.grpname[mtop->groups.grps[egcENER].nm_ind[j]]), ftp2ext(efXVG));
1451 fprintf(fp, "loading tf table for energygrp index %d from %s\n", ir->adress->tf_table_index[i], buf);
1453 fr->atf_tabs[i] = make_atf_table(fp, oenv, fr, buf, box);
1458 gmx_bool can_use_allvsall(const t_inputrec *ir, gmx_bool bPrintNote, t_commrec *cr, FILE *fp)
1465 ir->rcoulomb == 0 &&
1467 ir->ePBC == epbcNONE &&
1468 ir->vdwtype == evdwCUT &&
1469 ir->coulombtype == eelCUT &&
1470 ir->efep == efepNO &&
1471 (ir->implicit_solvent == eisNO ||
1472 (ir->implicit_solvent == eisGBSA && (ir->gb_algorithm == egbSTILL ||
1473 ir->gb_algorithm == egbHCT ||
1474 ir->gb_algorithm == egbOBC))) &&
1475 getenv("GMX_NO_ALLVSALL") == NULL
1478 if (bAllvsAll && ir->opts.ngener > 1)
1480 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";
1486 fprintf(stderr, "\n%s\n", note);
1490 fprintf(fp, "\n%s\n", note);
1496 if (bAllvsAll && fp && MASTER(cr))
1498 fprintf(fp, "\nUsing accelerated all-vs-all kernels.\n\n");
1505 static void init_forcerec_f_threads(t_forcerec *fr, int nenergrp)
1509 /* These thread local data structures are used for bondeds only */
1510 fr->nthreads = gmx_omp_nthreads_get(emntBonded);
1512 if (fr->nthreads > 1)
1514 snew(fr->f_t, fr->nthreads);
1515 /* Thread 0 uses the global force and energy arrays */
1516 for (t = 1; t < fr->nthreads; t++)
1518 fr->f_t[t].f = NULL;
1519 fr->f_t[t].f_nalloc = 0;
1520 snew(fr->f_t[t].fshift, SHIFTS);
1521 fr->f_t[t].grpp.nener = nenergrp*nenergrp;
1522 for (i = 0; i < egNR; i++)
1524 snew(fr->f_t[t].grpp.ener[i], fr->f_t[t].grpp.nener);
1531 static void pick_nbnxn_kernel_cpu(const t_inputrec gmx_unused *ir,
1535 *kernel_type = nbnxnk4x4_PlainC;
1536 *ewald_excl = ewaldexclTable;
1538 #ifdef GMX_NBNXN_SIMD
1540 #ifdef GMX_NBNXN_SIMD_4XN
1541 *kernel_type = nbnxnk4xN_SIMD_4xN;
1543 #ifdef GMX_NBNXN_SIMD_2XNN
1544 /* We expect the 2xNN kernels to be faster in most cases */
1545 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1548 #if defined GMX_NBNXN_SIMD_4XN && defined GMX_X86_AVX_256
1549 if (EEL_RF(ir->coulombtype) || ir->coulombtype == eelCUT)
1551 /* The raw pair rate of the 4x8 kernel is higher than 2x(4+4),
1552 * 10% with HT, 50% without HT, but extra zeros interactions
1553 * can compensate. As we currently don't detect the actual use
1554 * of HT, switch to 4x8 to avoid a potential performance hit.
1556 *kernel_type = nbnxnk4xN_SIMD_4xN;
1559 if (getenv("GMX_NBNXN_SIMD_4XN") != NULL)
1561 #ifdef GMX_NBNXN_SIMD_4XN
1562 *kernel_type = nbnxnk4xN_SIMD_4xN;
1564 gmx_fatal(FARGS, "SIMD 4xN kernels requested, but Gromacs has been compiled without support for these kernels");
1567 if (getenv("GMX_NBNXN_SIMD_2XNN") != NULL)
1569 #ifdef GMX_NBNXN_SIMD_2XNN
1570 *kernel_type = nbnxnk4xN_SIMD_2xNN;
1572 gmx_fatal(FARGS, "SIMD 2x(N+N) kernels requested, but Gromacs has been compiled without support for these kernels");
1576 /* Analytical Ewald exclusion correction is only an option in
1577 * the SIMD kernel. On BlueGene/Q, this is faster regardless
1578 * of precision. In single precision, this is faster on
1579 * Bulldozer, and slightly faster on Sandy Bridge.
1581 #if ((defined GMX_X86_AVX_128_FMA || defined GMX_X86_AVX_256) && !defined GMX_DOUBLE) || (defined GMX_CPU_ACCELERATION_IBM_QPX)
1582 *ewald_excl = ewaldexclAnalytical;
1584 if (getenv("GMX_NBNXN_EWALD_TABLE") != NULL)
1586 *ewald_excl = ewaldexclTable;
1588 if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != NULL)
1590 *ewald_excl = ewaldexclAnalytical;
1594 #endif /* GMX_NBNXN_SIMD */
1598 const char *lookup_nbnxn_kernel_name(int kernel_type)
1600 const char *returnvalue = NULL;
1601 switch (kernel_type)
1604 returnvalue = "not set";
1606 case nbnxnk4x4_PlainC:
1607 returnvalue = "plain C";
1609 case nbnxnk4xN_SIMD_4xN:
1610 case nbnxnk4xN_SIMD_2xNN:
1611 #ifdef GMX_NBNXN_SIMD
1613 /* We have x86 SSE2 compatible SIMD */
1614 #ifdef GMX_X86_AVX_128_FMA
1615 returnvalue = "AVX-128-FMA";
1617 #if defined GMX_X86_AVX_256 || defined __AVX__
1618 /* x86 SIMD intrinsics can be converted to SSE or AVX depending
1619 * on compiler flags. As we use nearly identical intrinsics,
1620 * compiling for AVX without an AVX macros effectively results
1622 * For gcc we check for __AVX__
1623 * At least a check for icc should be added (if there is a macro)
1625 #if defined GMX_X86_AVX_256 && !defined GMX_NBNXN_HALF_WIDTH_SIMD
1626 returnvalue = "AVX-256";
1628 returnvalue = "AVX-128";
1631 #ifdef GMX_X86_SSE4_1
1632 returnvalue = "SSE4.1";
1634 returnvalue = "SSE2";
1638 #else /* GMX_X86_SSE2 */
1639 /* not GMX_X86_SSE2, but other SIMD */
1640 returnvalue = "SIMD";
1641 #endif /* GMX_X86_SSE2 */
1642 #else /* GMX_NBNXN_SIMD */
1643 returnvalue = "not available";
1644 #endif /* GMX_NBNXN_SIMD */
1646 case nbnxnk8x8x8_CUDA: returnvalue = "CUDA"; break;
1647 case nbnxnk8x8x8_PlainC: returnvalue = "plain C"; break;
1651 gmx_fatal(FARGS, "Illegal kernel type selected");
1658 static void pick_nbnxn_kernel(FILE *fp,
1659 const t_commrec *cr,
1660 gmx_bool use_cpu_acceleration,
1662 gmx_bool bEmulateGPU,
1663 const t_inputrec *ir,
1666 gmx_bool bDoNonbonded)
1668 assert(kernel_type);
1670 *kernel_type = nbnxnkNotSet;
1671 *ewald_excl = ewaldexclTable;
1675 *kernel_type = nbnxnk8x8x8_PlainC;
1679 md_print_warn(cr, fp, "Emulating a GPU run on the CPU (slow)");
1684 *kernel_type = nbnxnk8x8x8_CUDA;
1687 if (*kernel_type == nbnxnkNotSet)
1689 if (use_cpu_acceleration)
1691 pick_nbnxn_kernel_cpu(ir, kernel_type, ewald_excl);
1695 *kernel_type = nbnxnk4x4_PlainC;
1699 if (bDoNonbonded && fp != NULL)
1701 fprintf(fp, "\nUsing %s %dx%d non-bonded kernels\n\n",
1702 lookup_nbnxn_kernel_name(*kernel_type),
1703 nbnxn_kernel_pairlist_simple(*kernel_type) ? NBNXN_CPU_CLUSTER_I_SIZE : NBNXN_GPU_CLUSTER_SIZE,
1704 nbnxn_kernel_to_cj_size(*kernel_type));
1708 static void pick_nbnxn_resources(const t_commrec *cr,
1709 const gmx_hw_info_t *hwinfo,
1710 gmx_bool bDoNonbonded,
1712 gmx_bool *bEmulateGPU,
1713 const gmx_gpu_opt_t *gpu_opt)
1715 gmx_bool bEmulateGPUEnvVarSet;
1716 char gpu_err_str[STRLEN];
1720 bEmulateGPUEnvVarSet = (getenv("GMX_EMULATE_GPU") != NULL);
1722 /* Run GPU emulation mode if GMX_EMULATE_GPU is defined. Because
1723 * GPUs (currently) only handle non-bonded calculations, we will
1724 * automatically switch to emulation if non-bonded calculations are
1725 * turned off via GMX_NO_NONBONDED - this is the simple and elegant
1726 * way to turn off GPU initialization, data movement, and cleanup.
1728 * GPU emulation can be useful to assess the performance one can expect by
1729 * adding GPU(s) to the machine. The conditional below allows this even
1730 * if mdrun is compiled without GPU acceleration support.
1731 * Note that you should freezing the system as otherwise it will explode.
1733 *bEmulateGPU = (bEmulateGPUEnvVarSet ||
1735 gpu_opt->ncuda_dev_use > 0));
1737 /* Enable GPU mode when GPUs are available or no GPU emulation is requested.
1739 if (gpu_opt->ncuda_dev_use > 0 && !(*bEmulateGPU))
1741 /* Each PP node will use the intra-node id-th device from the
1742 * list of detected/selected GPUs. */
1743 if (!init_gpu(cr->rank_pp_intranode, gpu_err_str,
1744 &hwinfo->gpu_info, gpu_opt))
1746 /* At this point the init should never fail as we made sure that
1747 * we have all the GPUs we need. If it still does, we'll bail. */
1748 gmx_fatal(FARGS, "On node %d failed to initialize GPU #%d: %s",
1750 get_gpu_device_id(&hwinfo->gpu_info, gpu_opt,
1751 cr->rank_pp_intranode),
1755 /* Here we actually turn on hardware GPU acceleration */
1760 gmx_bool uses_simple_tables(int cutoff_scheme,
1761 nonbonded_verlet_t *nbv,
1764 gmx_bool bUsesSimpleTables = TRUE;
1767 switch (cutoff_scheme)
1770 bUsesSimpleTables = TRUE;
1773 assert(NULL != nbv && NULL != nbv->grp);
1774 grp_index = (group < 0) ? 0 : (nbv->ngrp - 1);
1775 bUsesSimpleTables = nbnxn_kernel_pairlist_simple(nbv->grp[grp_index].kernel_type);
1778 gmx_incons("unimplemented");
1780 return bUsesSimpleTables;
1783 static void init_ewald_f_table(interaction_const_t *ic,
1784 gmx_bool bUsesSimpleTables,
1789 if (bUsesSimpleTables)
1791 /* With a spacing of 0.0005 we are at the force summation accuracy
1792 * for the SSE kernels for "normal" atomistic simulations.
1794 ic->tabq_scale = ewald_spline3_table_scale(ic->ewaldcoeff_q,
1797 maxr = (rtab > ic->rcoulomb) ? rtab : ic->rcoulomb;
1798 ic->tabq_size = (int)(maxr*ic->tabq_scale) + 2;
1802 ic->tabq_size = GPU_EWALD_COULOMB_FORCE_TABLE_SIZE;
1803 /* Subtract 2 iso 1 to avoid access out of range due to rounding */
1804 ic->tabq_scale = (ic->tabq_size - 2)/ic->rcoulomb;
1807 sfree_aligned(ic->tabq_coul_FDV0);
1808 sfree_aligned(ic->tabq_coul_F);
1809 sfree_aligned(ic->tabq_coul_V);
1811 /* Create the original table data in FDV0 */
1812 snew_aligned(ic->tabq_coul_FDV0, ic->tabq_size*4, 32);
1813 snew_aligned(ic->tabq_coul_F, ic->tabq_size, 32);
1814 snew_aligned(ic->tabq_coul_V, ic->tabq_size, 32);
1815 table_spline3_fill_ewald_lr(ic->tabq_coul_F, ic->tabq_coul_V, ic->tabq_coul_FDV0,
1816 ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_q);
1819 void init_interaction_const_tables(FILE *fp,
1820 interaction_const_t *ic,
1821 gmx_bool bUsesSimpleTables,
1826 if (ic->eeltype == eelEWALD || EEL_PME(ic->eeltype))
1828 init_ewald_f_table(ic, bUsesSimpleTables, rtab);
1832 fprintf(fp, "Initialized non-bonded Ewald correction tables, spacing: %.2e size: %d\n\n",
1833 1/ic->tabq_scale, ic->tabq_size);
1838 static void init_interaction_const(FILE *fp,
1839 const t_commrec gmx_unused *cr,
1840 interaction_const_t **interaction_const,
1841 const t_forcerec *fr,
1844 interaction_const_t *ic;
1845 gmx_bool bUsesSimpleTables = TRUE;
1849 /* Just allocate something so we can free it */
1850 snew_aligned(ic->tabq_coul_FDV0, 16, 32);
1851 snew_aligned(ic->tabq_coul_F, 16, 32);
1852 snew_aligned(ic->tabq_coul_V, 16, 32);
1854 ic->rlist = fr->rlist;
1855 ic->rlistlong = fr->rlistlong;
1858 ic->rvdw = fr->rvdw;
1859 if (fr->vdw_modifier == eintmodPOTSHIFT)
1861 ic->sh_invrc6 = pow(ic->rvdw, -6.0);
1868 /* Electrostatics */
1869 ic->eeltype = fr->eeltype;
1870 ic->rcoulomb = fr->rcoulomb;
1871 ic->epsilon_r = fr->epsilon_r;
1872 ic->epsfac = fr->epsfac;
1875 ic->ewaldcoeff_q = fr->ewaldcoeff_q;
1876 ic->ewaldcoeff_lj = fr->ewaldcoeff_lj;
1878 if (fr->coulomb_modifier == eintmodPOTSHIFT)
1880 ic->sh_ewald = gmx_erfc(ic->ewaldcoeff_q*ic->rcoulomb);
1887 /* Reaction-field */
1888 if (EEL_RF(ic->eeltype))
1890 ic->epsilon_rf = fr->epsilon_rf;
1891 ic->k_rf = fr->k_rf;
1892 ic->c_rf = fr->c_rf;
1896 /* For plain cut-off we might use the reaction-field kernels */
1897 ic->epsilon_rf = ic->epsilon_r;
1899 if (fr->coulomb_modifier == eintmodPOTSHIFT)
1901 ic->c_rf = 1/ic->rcoulomb;
1911 fprintf(fp, "Potential shift: LJ r^-12: %.3f r^-6 %.3f",
1912 sqr(ic->sh_invrc6), ic->sh_invrc6);
1913 if (ic->eeltype == eelCUT)
1915 fprintf(fp, ", Coulomb %.3f", ic->c_rf);
1917 else if (EEL_PME(ic->eeltype))
1919 fprintf(fp, ", Ewald %.3e", ic->sh_ewald);
1924 *interaction_const = ic;
1926 if (fr->nbv != NULL && fr->nbv->bUseGPU)
1928 nbnxn_cuda_init_const(fr->nbv->cu_nbv, ic, fr->nbv->grp);
1930 /* With tMPI + GPUs some ranks may be sharing GPU(s) and therefore
1931 * also sharing texture references. To keep the code simple, we don't
1932 * treat texture references as shared resources, but this means that
1933 * the coulomb_tab and nbfp texture refs will get updated by multiple threads.
1934 * Hence, to ensure that the non-bonded kernels don't start before all
1935 * texture binding operations are finished, we need to wait for all ranks
1936 * to arrive here before continuing.
1938 * Note that we could omit this barrier if GPUs are not shared (or
1939 * texture objects are used), but as this is initialization code, there
1940 * is not point in complicating things.
1942 #ifdef GMX_THREAD_MPI
1947 #endif /* GMX_THREAD_MPI */
1950 bUsesSimpleTables = uses_simple_tables(fr->cutoff_scheme, fr->nbv, -1);
1951 init_interaction_const_tables(fp, ic, bUsesSimpleTables, rtab);
1954 static void init_nb_verlet(FILE *fp,
1955 nonbonded_verlet_t **nb_verlet,
1956 const t_inputrec *ir,
1957 const t_forcerec *fr,
1958 const t_commrec *cr,
1959 const char *nbpu_opt)
1961 nonbonded_verlet_t *nbv;
1964 gmx_bool bEmulateGPU, bHybridGPURun = FALSE;
1966 nbnxn_alloc_t *nb_alloc;
1967 nbnxn_free_t *nb_free;
1971 pick_nbnxn_resources(cr, fr->hwinfo,
1979 nbv->ngrp = (DOMAINDECOMP(cr) ? 2 : 1);
1980 for (i = 0; i < nbv->ngrp; i++)
1982 nbv->grp[i].nbl_lists.nnbl = 0;
1983 nbv->grp[i].nbat = NULL;
1984 nbv->grp[i].kernel_type = nbnxnkNotSet;
1986 if (i == 0) /* local */
1988 pick_nbnxn_kernel(fp, cr, fr->use_cpu_acceleration,
1989 nbv->bUseGPU, bEmulateGPU, ir,
1990 &nbv->grp[i].kernel_type,
1991 &nbv->grp[i].ewald_excl,
1994 else /* non-local */
1996 if (nbpu_opt != NULL && strcmp(nbpu_opt, "gpu_cpu") == 0)
1998 /* Use GPU for local, select a CPU kernel for non-local */
1999 pick_nbnxn_kernel(fp, cr, fr->use_cpu_acceleration,
2001 &nbv->grp[i].kernel_type,
2002 &nbv->grp[i].ewald_excl,
2005 bHybridGPURun = TRUE;
2009 /* Use the same kernel for local and non-local interactions */
2010 nbv->grp[i].kernel_type = nbv->grp[0].kernel_type;
2011 nbv->grp[i].ewald_excl = nbv->grp[0].ewald_excl;
2018 /* init the NxN GPU data; the last argument tells whether we'll have
2019 * both local and non-local NB calculation on GPU */
2020 nbnxn_cuda_init(fp, &nbv->cu_nbv,
2021 &fr->hwinfo->gpu_info, fr->gpu_opt,
2022 cr->rank_pp_intranode,
2023 (nbv->ngrp > 1) && !bHybridGPURun);
2025 if ((env = getenv("GMX_NB_MIN_CI")) != NULL)
2029 nbv->min_ci_balanced = strtol(env, &end, 10);
2030 if (!end || (*end != 0) || nbv->min_ci_balanced <= 0)
2032 gmx_fatal(FARGS, "Invalid value passed in GMX_NB_MIN_CI=%s, positive integer required", env);
2037 fprintf(debug, "Neighbor-list balancing parameter: %d (passed as env. var.)\n",
2038 nbv->min_ci_balanced);
2043 nbv->min_ci_balanced = nbnxn_cuda_min_ci_balanced(nbv->cu_nbv);
2046 fprintf(debug, "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU multi-processors)\n",
2047 nbv->min_ci_balanced);
2053 nbv->min_ci_balanced = 0;
2058 nbnxn_init_search(&nbv->nbs,
2059 DOMAINDECOMP(cr) ? &cr->dd->nc : NULL,
2060 DOMAINDECOMP(cr) ? domdec_zones(cr->dd) : NULL,
2061 gmx_omp_nthreads_get(emntNonbonded));
2063 for (i = 0; i < nbv->ngrp; i++)
2065 if (nbv->grp[0].kernel_type == nbnxnk8x8x8_CUDA)
2067 nb_alloc = &pmalloc;
2076 nbnxn_init_pairlist_set(&nbv->grp[i].nbl_lists,
2077 nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2078 /* 8x8x8 "non-simple" lists are ATM always combined */
2079 !nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2083 nbv->grp[0].kernel_type != nbv->grp[i].kernel_type)
2085 snew(nbv->grp[i].nbat, 1);
2086 nbnxn_atomdata_init(fp,
2088 nbv->grp[i].kernel_type,
2089 fr->ntype, fr->nbfp,
2091 nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type) ? gmx_omp_nthreads_get(emntNonbonded) : 1,
2096 nbv->grp[i].nbat = nbv->grp[0].nbat;
2101 void init_forcerec(FILE *fp,
2102 const output_env_t oenv,
2105 const t_inputrec *ir,
2106 const gmx_mtop_t *mtop,
2107 const t_commrec *cr,
2113 const char *nbpu_opt,
2114 gmx_bool bNoSolvOpt,
2117 int i, j, m, natoms, ngrp, negp_pp, negptable, egi, egj;
2122 gmx_bool bGenericKernelOnly;
2123 gmx_bool bTab, bSep14tab, bNormalnblists;
2125 int *nm_ind, egp_flags;
2127 if (fr->hwinfo == NULL)
2129 /* Detect hardware, gather information.
2130 * In mdrun, hwinfo has already been set before calling init_forcerec.
2131 * Here we ignore GPUs, as tools will not use them anyhow.
2133 fr->hwinfo = gmx_detect_hardware(fp, cr, FALSE);
2136 /* By default we turn acceleration on, but it might be turned off further down... */
2137 fr->use_cpu_acceleration = TRUE;
2139 fr->bDomDec = DOMAINDECOMP(cr);
2141 natoms = mtop->natoms;
2143 if (check_box(ir->ePBC, box))
2145 gmx_fatal(FARGS, check_box(ir->ePBC, box));
2148 /* Test particle insertion ? */
2151 /* Set to the size of the molecule to be inserted (the last one) */
2152 /* Because of old style topologies, we have to use the last cg
2153 * instead of the last molecule type.
2155 cgs = &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].cgs;
2156 fr->n_tpi = cgs->index[cgs->nr] - cgs->index[cgs->nr-1];
2157 if (fr->n_tpi != mtop->mols.index[mtop->mols.nr] - mtop->mols.index[mtop->mols.nr-1])
2159 gmx_fatal(FARGS, "The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
2167 /* Copy AdResS parameters */
2170 fr->adress_type = ir->adress->type;
2171 fr->adress_const_wf = ir->adress->const_wf;
2172 fr->adress_ex_width = ir->adress->ex_width;
2173 fr->adress_hy_width = ir->adress->hy_width;
2174 fr->adress_icor = ir->adress->icor;
2175 fr->adress_site = ir->adress->site;
2176 fr->adress_ex_forcecap = ir->adress->ex_forcecap;
2177 fr->adress_do_hybridpairs = ir->adress->do_hybridpairs;
2180 snew(fr->adress_group_explicit, ir->adress->n_energy_grps);
2181 for (i = 0; i < ir->adress->n_energy_grps; i++)
2183 fr->adress_group_explicit[i] = ir->adress->group_explicit[i];
2186 fr->n_adress_tf_grps = ir->adress->n_tf_grps;
2187 snew(fr->adress_tf_table_index, fr->n_adress_tf_grps);
2188 for (i = 0; i < fr->n_adress_tf_grps; i++)
2190 fr->adress_tf_table_index[i] = ir->adress->tf_table_index[i];
2192 copy_rvec(ir->adress->refs, fr->adress_refs);
2196 fr->adress_type = eAdressOff;
2197 fr->adress_do_hybridpairs = FALSE;
2200 /* Copy the user determined parameters */
2201 fr->userint1 = ir->userint1;
2202 fr->userint2 = ir->userint2;
2203 fr->userint3 = ir->userint3;
2204 fr->userint4 = ir->userint4;
2205 fr->userreal1 = ir->userreal1;
2206 fr->userreal2 = ir->userreal2;
2207 fr->userreal3 = ir->userreal3;
2208 fr->userreal4 = ir->userreal4;
2211 fr->fc_stepsize = ir->fc_stepsize;
2214 fr->efep = ir->efep;
2215 fr->sc_alphavdw = ir->fepvals->sc_alpha;
2216 if (ir->fepvals->bScCoul)
2218 fr->sc_alphacoul = ir->fepvals->sc_alpha;
2219 fr->sc_sigma6_min = pow(ir->fepvals->sc_sigma_min, 6);
2223 fr->sc_alphacoul = 0;
2224 fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
2226 fr->sc_power = ir->fepvals->sc_power;
2227 fr->sc_r_power = ir->fepvals->sc_r_power;
2228 fr->sc_sigma6_def = pow(ir->fepvals->sc_sigma, 6);
2230 env = getenv("GMX_SCSIGMA_MIN");
2234 sscanf(env, "%lf", &dbl);
2235 fr->sc_sigma6_min = pow(dbl, 6);
2238 fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
2242 fr->bNonbonded = TRUE;
2243 if (getenv("GMX_NO_NONBONDED") != NULL)
2245 /* turn off non-bonded calculations */
2246 fr->bNonbonded = FALSE;
2247 md_print_warn(cr, fp,
2248 "Found environment variable GMX_NO_NONBONDED.\n"
2249 "Disabling nonbonded calculations.\n");
2252 bGenericKernelOnly = FALSE;
2254 /* We now check in the NS code whether a particular combination of interactions
2255 * can be used with water optimization, and disable it if that is not the case.
2258 if (getenv("GMX_NB_GENERIC") != NULL)
2263 "Found environment variable GMX_NB_GENERIC.\n"
2264 "Disabling all interaction-specific nonbonded kernels, will only\n"
2265 "use the slow generic ones in src/gmxlib/nonbonded/nb_generic.c\n\n");
2267 bGenericKernelOnly = TRUE;
2270 if (bGenericKernelOnly == TRUE)
2275 if ( (getenv("GMX_DISABLE_CPU_ACCELERATION") != NULL) || (getenv("GMX_NOOPTIMIZEDKERNELS") != NULL) )
2277 fr->use_cpu_acceleration = FALSE;
2281 "\nFound environment variable GMX_DISABLE_CPU_ACCELERATION.\n"
2282 "Disabling all CPU architecture-specific (e.g. SSE2/SSE4/AVX) routines.\n\n");
2286 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
2288 /* Check if we can/should do all-vs-all kernels */
2289 fr->bAllvsAll = can_use_allvsall(ir, FALSE, NULL, NULL);
2290 fr->AllvsAll_work = NULL;
2291 fr->AllvsAll_workgb = NULL;
2293 /* All-vs-all kernels have not been implemented in 4.6, and
2294 * the SIMD group kernels are also buggy in this case. Non-accelerated
2295 * group kernels are OK. See Redmine #1249. */
2298 fr->bAllvsAll = FALSE;
2299 fr->use_cpu_acceleration = FALSE;
2303 "\nYour simulation settings would have triggered the efficient all-vs-all\n"
2304 "kernels in GROMACS 4.5, but these have not been implemented in GROMACS\n"
2305 "4.6. Also, we can't use the accelerated SIMD kernels here because\n"
2306 "of an unfixed bug. The reference C kernels are correct, though, so\n"
2307 "we are proceeding by disabling all CPU architecture-specific\n"
2308 "(e.g. SSE2/SSE4/AVX) routines. If performance is important, please\n"
2309 "use GROMACS 4.5.7 or try cutoff-scheme = Verlet.\n\n");
2313 /* Neighbour searching stuff */
2314 fr->cutoff_scheme = ir->cutoff_scheme;
2315 fr->bGrid = (ir->ns_type == ensGRID);
2316 fr->ePBC = ir->ePBC;
2318 /* Determine if we will do PBC for distances in bonded interactions */
2319 if (fr->ePBC == epbcNONE)
2321 fr->bMolPBC = FALSE;
2325 if (!DOMAINDECOMP(cr))
2327 /* The group cut-off scheme and SHAKE assume charge groups
2328 * are whole, but not using molpbc is faster in most cases.
2330 if (fr->cutoff_scheme == ecutsGROUP ||
2331 (ir->eConstrAlg == econtSHAKE &&
2332 (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0 ||
2333 gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0)))
2335 fr->bMolPBC = ir->bPeriodicMols;
2340 if (getenv("GMX_USE_GRAPH") != NULL)
2342 fr->bMolPBC = FALSE;
2345 fprintf(fp, "\nGMX_MOLPBC is set, using the graph for bonded interactions\n\n");
2352 fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->ePBC);
2355 fr->bGB = (ir->implicit_solvent == eisGBSA);
2357 fr->rc_scaling = ir->refcoord_scaling;
2358 copy_rvec(ir->posres_com, fr->posres_com);
2359 copy_rvec(ir->posres_comB, fr->posres_comB);
2360 fr->rlist = cutoff_inf(ir->rlist);
2361 fr->rlistlong = cutoff_inf(ir->rlistlong);
2362 fr->eeltype = ir->coulombtype;
2363 fr->vdwtype = ir->vdwtype;
2364 fr->ljpme_combination_rule = ir->ljpme_combination_rule;
2366 fr->coulomb_modifier = ir->coulomb_modifier;
2367 fr->vdw_modifier = ir->vdw_modifier;
2369 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
2370 switch (fr->eeltype)
2373 fr->nbkernel_elec_interaction = (fr->bGB) ? GMX_NBKERNEL_ELEC_GENERALIZEDBORN : GMX_NBKERNEL_ELEC_COULOMB;
2379 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2383 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2384 fr->coulomb_modifier = eintmodEXACTCUTOFF;
2393 case eelPMEUSERSWITCH:
2394 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2399 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD;
2403 gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[fr->eeltype]);
2407 /* Vdw: Translate from mdp settings to kernel format */
2408 switch (fr->vdwtype)
2414 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
2418 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
2425 case evdwENCADSHIFT:
2426 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2430 gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[fr->vdwtype]);
2434 /* These start out identical to ir, but might be altered if we e.g. tabulate the interaction in the kernel */
2435 fr->nbkernel_elec_modifier = fr->coulomb_modifier;
2436 fr->nbkernel_vdw_modifier = fr->vdw_modifier;
2438 fr->bTwinRange = fr->rlistlong > fr->rlist;
2439 fr->bEwald = (EEL_PME(fr->eeltype) || fr->eeltype == eelEWALD);
2441 fr->reppow = mtop->ffparams.reppow;
2443 if (ir->cutoff_scheme == ecutsGROUP)
2445 fr->bvdwtab = (fr->vdwtype != evdwCUT ||
2446 !gmx_within_tol(fr->reppow, 12.0, 10*GMX_DOUBLE_EPS));
2447 /* We have special kernels for standard Ewald and PME, but the pme-switch ones are tabulated above */
2448 fr->bcoultab = !(fr->eeltype == eelCUT ||
2449 fr->eeltype == eelEWALD ||
2450 fr->eeltype == eelPME ||
2451 fr->eeltype == eelRF ||
2452 fr->eeltype == eelRF_ZERO);
2454 /* If the user absolutely wants different switch/shift settings for coul/vdw, it is likely
2455 * going to be faster to tabulate the interaction than calling the generic kernel.
2457 if (fr->nbkernel_elec_modifier == eintmodPOTSWITCH && fr->nbkernel_vdw_modifier == eintmodPOTSWITCH)
2459 if ((fr->rcoulomb_switch != fr->rvdw_switch) || (fr->rcoulomb != fr->rvdw))
2461 fr->bcoultab = TRUE;
2464 else if ((fr->nbkernel_elec_modifier == eintmodPOTSHIFT && fr->nbkernel_vdw_modifier == eintmodPOTSHIFT) ||
2465 ((fr->nbkernel_elec_interaction == GMX_NBKERNEL_ELEC_REACTIONFIELD &&
2466 fr->nbkernel_elec_modifier == eintmodEXACTCUTOFF &&
2467 (fr->nbkernel_vdw_modifier == eintmodPOTSWITCH || fr->nbkernel_vdw_modifier == eintmodPOTSHIFT))))
2469 if (fr->rcoulomb != fr->rvdw)
2471 fr->bcoultab = TRUE;
2475 if (getenv("GMX_REQUIRE_TABLES"))
2478 fr->bcoultab = TRUE;
2483 fprintf(fp, "Table routines are used for coulomb: %s\n", bool_names[fr->bcoultab]);
2484 fprintf(fp, "Table routines are used for vdw: %s\n", bool_names[fr->bvdwtab ]);
2487 if (fr->bvdwtab == TRUE)
2489 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2490 fr->nbkernel_vdw_modifier = eintmodNONE;
2492 if (fr->bcoultab == TRUE)
2494 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2495 fr->nbkernel_elec_modifier = eintmodNONE;
2499 if (ir->cutoff_scheme == ecutsVERLET)
2501 if (!gmx_within_tol(fr->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2503 gmx_fatal(FARGS, "Cut-off scheme %S only supports LJ repulsion power 12", ecutscheme_names[ir->cutoff_scheme]);
2505 fr->bvdwtab = FALSE;
2506 fr->bcoultab = FALSE;
2509 /* Tables are used for direct ewald sum */
2512 if (EEL_PME(ir->coulombtype))
2516 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
2518 if (ir->coulombtype == eelP3M_AD)
2520 please_cite(fp, "Hockney1988");
2521 please_cite(fp, "Ballenegger2012");
2525 please_cite(fp, "Essmann95a");
2528 if (ir->ewald_geometry == eewg3DC)
2532 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry.\n");
2534 please_cite(fp, "In-Chul99a");
2537 fr->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
2538 init_ewald_tab(&(fr->ewald_table), ir, fp);
2541 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n",
2542 1/fr->ewaldcoeff_q);
2546 if (EVDW_PME(ir->vdwtype))
2550 fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
2552 please_cite(fp, "Essman95a");
2553 fr->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
2556 fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n",
2557 1/fr->ewaldcoeff_lj);
2561 /* Electrostatics */
2562 fr->epsilon_r = ir->epsilon_r;
2563 fr->epsilon_rf = ir->epsilon_rf;
2564 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
2565 fr->rcoulomb_switch = ir->rcoulomb_switch;
2566 fr->rcoulomb = cutoff_inf(ir->rcoulomb);
2568 /* Parameters for generalized RF */
2572 if (fr->eeltype == eelGRF)
2574 init_generalized_rf(fp, mtop, ir, fr);
2577 fr->bF_NoVirSum = (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype) ||
2578 gmx_mtop_ftype_count(mtop, F_POSRES) > 0 ||
2579 gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 ||
2580 IR_ELEC_FIELD(*ir) ||
2581 (fr->adress_icor != eAdressICOff)
2584 if (fr->cutoff_scheme == ecutsGROUP &&
2585 ncg_mtop(mtop) > fr->cg_nalloc && !DOMAINDECOMP(cr))
2587 /* Count the total number of charge groups */
2588 fr->cg_nalloc = ncg_mtop(mtop);
2589 srenew(fr->cg_cm, fr->cg_nalloc);
2591 if (fr->shift_vec == NULL)
2593 snew(fr->shift_vec, SHIFTS);
2596 if (fr->fshift == NULL)
2598 snew(fr->fshift, SHIFTS);
2601 if (fr->nbfp == NULL)
2603 fr->ntype = mtop->ffparams.atnr;
2604 fr->nbfp = mk_nbfp(&mtop->ffparams, fr->bBHAM);
2607 /* Copy the energy group exclusions */
2608 fr->egp_flags = ir->opts.egp_flags;
2610 /* Van der Waals stuff */
2611 fr->rvdw = cutoff_inf(ir->rvdw);
2612 fr->rvdw_switch = ir->rvdw_switch;
2613 if ((fr->vdwtype != evdwCUT) && (fr->vdwtype != evdwUSER) && !fr->bBHAM)
2615 if (fr->rvdw_switch >= fr->rvdw)
2617 gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)",
2618 fr->rvdw_switch, fr->rvdw);
2622 fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
2623 (fr->eeltype == eelSWITCH) ? "switched" : "shifted",
2624 fr->rvdw_switch, fr->rvdw);
2628 if (fr->bBHAM && EVDW_PME(fr->vdwtype))
2630 gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
2633 if (fr->bBHAM && (fr->vdwtype == evdwSHIFT || fr->vdwtype == evdwSWITCH))
2635 gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
2640 fprintf(fp, "Cut-off's: NS: %g Coulomb: %g %s: %g\n",
2641 fr->rlist, fr->rcoulomb, fr->bBHAM ? "BHAM" : "LJ", fr->rvdw);
2644 fr->eDispCorr = ir->eDispCorr;
2645 if (ir->eDispCorr != edispcNO)
2647 set_avcsixtwelve(fp, fr, mtop);
2652 set_bham_b_max(fp, fr, mtop);
2655 fr->gb_epsilon_solvent = ir->gb_epsilon_solvent;
2657 /* Copy the GBSA data (radius, volume and surftens for each
2658 * atomtype) from the topology atomtype section to forcerec.
2660 snew(fr->atype_radius, fr->ntype);
2661 snew(fr->atype_vol, fr->ntype);
2662 snew(fr->atype_surftens, fr->ntype);
2663 snew(fr->atype_gb_radius, fr->ntype);
2664 snew(fr->atype_S_hct, fr->ntype);
2666 if (mtop->atomtypes.nr > 0)
2668 for (i = 0; i < fr->ntype; i++)
2670 fr->atype_radius[i] = mtop->atomtypes.radius[i];
2672 for (i = 0; i < fr->ntype; i++)
2674 fr->atype_vol[i] = mtop->atomtypes.vol[i];
2676 for (i = 0; i < fr->ntype; i++)
2678 fr->atype_surftens[i] = mtop->atomtypes.surftens[i];
2680 for (i = 0; i < fr->ntype; i++)
2682 fr->atype_gb_radius[i] = mtop->atomtypes.gb_radius[i];
2684 for (i = 0; i < fr->ntype; i++)
2686 fr->atype_S_hct[i] = mtop->atomtypes.S_hct[i];
2690 /* Generate the GB table if needed */
2694 fr->gbtabscale = 2000;
2696 fr->gbtabscale = 500;
2700 fr->gbtab = make_gb_table(oenv, fr);
2702 init_gb(&fr->born, cr, fr, ir, mtop, ir->gb_algorithm);
2704 /* Copy local gb data (for dd, this is done in dd_partition_system) */
2705 if (!DOMAINDECOMP(cr))
2707 make_local_gb(cr, fr->born, ir->gb_algorithm);
2711 /* Set the charge scaling */
2712 if (fr->epsilon_r != 0)
2714 fr->epsfac = ONE_4PI_EPS0/fr->epsilon_r;
2718 /* eps = 0 is infinite dieletric: no coulomb interactions */
2722 /* Reaction field constants */
2723 if (EEL_RF(fr->eeltype))
2725 calc_rffac(fp, fr->eeltype, fr->epsilon_r, fr->epsilon_rf,
2726 fr->rcoulomb, fr->temp, fr->zsquare, box,
2727 &fr->kappa, &fr->k_rf, &fr->c_rf);
2730 /*This now calculates sum for q and c6*/
2731 set_chargesum(fp, fr, mtop);
2733 /* if we are using LR electrostatics, and they are tabulated,
2734 * the tables will contain modified coulomb interactions.
2735 * Since we want to use the non-shifted ones for 1-4
2736 * coulombic interactions, we must have an extra set of tables.
2739 /* Construct tables.
2740 * A little unnecessary to make both vdw and coul tables sometimes,
2741 * but what the heck... */
2743 bTab = fr->bcoultab || fr->bvdwtab || fr->bEwald;
2745 bSep14tab = ((!bTab || fr->eeltype != eelCUT || fr->vdwtype != evdwCUT ||
2746 fr->bBHAM || fr->bEwald) &&
2747 (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
2748 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0 ||
2749 gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0));
2751 negp_pp = ir->opts.ngener - ir->nwall;
2755 bNormalnblists = TRUE;
2760 bNormalnblists = (ir->eDispCorr != edispcNO);
2761 for (egi = 0; egi < negp_pp; egi++)
2763 for (egj = egi; egj < negp_pp; egj++)
2765 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
2766 if (!(egp_flags & EGP_EXCL))
2768 if (egp_flags & EGP_TABLE)
2774 bNormalnblists = TRUE;
2781 fr->nnblists = negptable + 1;
2785 fr->nnblists = negptable;
2787 if (fr->nnblists > 1)
2789 snew(fr->gid2nblists, ir->opts.ngener*ir->opts.ngener);
2798 snew(fr->nblists, fr->nnblists);
2800 /* This code automatically gives table length tabext without cut-off's,
2801 * in that case grompp should already have checked that we do not need
2802 * normal tables and we only generate tables for 1-4 interactions.
2804 rtab = ir->rlistlong + ir->tabext;
2808 /* make tables for ordinary interactions */
2811 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn, NULL, NULL, &fr->nblists[0]);
2814 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn, NULL, NULL, &fr->nblists[fr->nnblists/2]);
2818 fr->tab14 = fr->nblists[0].table_elec_vdw;
2828 /* Read the special tables for certain energy group pairs */
2829 nm_ind = mtop->groups.grps[egcENER].nm_ind;
2830 for (egi = 0; egi < negp_pp; egi++)
2832 for (egj = egi; egj < negp_pp; egj++)
2834 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
2835 if ((egp_flags & EGP_TABLE) && !(egp_flags & EGP_EXCL))
2837 nbl = &(fr->nblists[m]);
2838 if (fr->nnblists > 1)
2840 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = m;
2842 /* Read the table file with the two energy groups names appended */
2843 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn,
2844 *mtop->groups.grpname[nm_ind[egi]],
2845 *mtop->groups.grpname[nm_ind[egj]],
2849 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn,
2850 *mtop->groups.grpname[nm_ind[egi]],
2851 *mtop->groups.grpname[nm_ind[egj]],
2852 &fr->nblists[fr->nnblists/2+m]);
2856 else if (fr->nnblists > 1)
2858 fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = 0;
2866 /* generate extra tables with plain Coulomb for 1-4 interactions only */
2867 fr->tab14 = make_tables(fp, oenv, fr, MASTER(cr), tabpfn, rtab,
2868 GMX_MAKETABLES_14ONLY);
2871 /* Read AdResS Thermo Force table if needed */
2872 if (fr->adress_icor == eAdressICThermoForce)
2874 /* old todo replace */
2876 if (ir->adress->n_tf_grps > 0)
2878 make_adress_tf_tables(fp, oenv, fr, ir, tabfn, mtop, box);
2883 /* load the default table */
2884 snew(fr->atf_tabs, 1);
2885 fr->atf_tabs[DEFAULT_TF_TABLE] = make_atf_table(fp, oenv, fr, tabafn, box);
2890 fr->nwall = ir->nwall;
2891 if (ir->nwall && ir->wall_type == ewtTABLE)
2893 make_wall_tables(fp, oenv, ir, tabfn, &mtop->groups, fr);
2898 fcd->bondtab = make_bonded_tables(fp,
2899 F_TABBONDS, F_TABBONDSNC,
2901 fcd->angletab = make_bonded_tables(fp,
2904 fcd->dihtab = make_bonded_tables(fp,
2912 fprintf(debug, "No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
2916 /* QM/MM initialization if requested
2920 fprintf(stderr, "QM/MM calculation requested.\n");
2923 fr->bQMMM = ir->bQMMM;
2924 fr->qr = mk_QMMMrec();
2926 /* Set all the static charge group info */
2927 fr->cginfo_mb = init_cginfo_mb(fp, mtop, fr, bNoSolvOpt,
2928 &fr->bExcl_IntraCGAll_InterCGNone);
2929 if (DOMAINDECOMP(cr))
2935 fr->cginfo = cginfo_expand(mtop->nmolblock, fr->cginfo_mb);
2938 if (!DOMAINDECOMP(cr))
2940 /* When using particle decomposition, the effect of the second argument,
2941 * which sets fr->hcg, is corrected later in do_md and init_em.
2943 forcerec_set_ranges(fr, ncg_mtop(mtop), ncg_mtop(mtop),
2944 mtop->natoms, mtop->natoms, mtop->natoms);
2947 fr->print_force = print_force;
2950 /* coarse load balancing vars */
2955 /* Initialize neighbor search */
2956 init_ns(fp, cr, &fr->ns, fr, mtop);
2958 if (cr->duty & DUTY_PP)
2960 gmx_nonbonded_setup(fr, bGenericKernelOnly);
2964 gmx_setup_adress_kernels(fp,bGenericKernelOnly);
2969 /* Initialize the thread working data for bonded interactions */
2970 init_forcerec_f_threads(fr, mtop->groups.grps[egcENER].nr);
2972 snew(fr->excl_load, fr->nthreads+1);
2974 if (fr->cutoff_scheme == ecutsVERLET)
2976 if (ir->rcoulomb != ir->rvdw)
2978 gmx_fatal(FARGS, "With Verlet lists rcoulomb and rvdw should be identical");
2981 init_nb_verlet(fp, &fr->nbv, ir, fr, cr, nbpu_opt);
2984 /* fr->ic is used both by verlet and group kernels (to some extent) now */
2985 init_interaction_const(fp, cr, &fr->ic, fr, rtab);
2987 if (ir->eDispCorr != edispcNO)
2989 calc_enervirdiff(fp, ir->eDispCorr, fr);
2993 #define pr_real(fp, r) fprintf(fp, "%s: %e\n",#r, r)
2994 #define pr_int(fp, i) fprintf((fp), "%s: %d\n",#i, i)
2995 #define pr_bool(fp, b) fprintf((fp), "%s: %s\n",#b, bool_names[b])
2997 void pr_forcerec(FILE *fp, t_forcerec *fr)
3001 pr_real(fp, fr->rlist);
3002 pr_real(fp, fr->rcoulomb);
3003 pr_real(fp, fr->fudgeQQ);
3004 pr_bool(fp, fr->bGrid);
3005 pr_bool(fp, fr->bTwinRange);
3006 /*pr_int(fp,fr->cg0);
3007 pr_int(fp,fr->hcg);*/
3008 for (i = 0; i < fr->nnblists; i++)
3010 pr_int(fp, fr->nblists[i].table_elec_vdw.n);
3012 pr_real(fp, fr->rcoulomb_switch);
3013 pr_real(fp, fr->rcoulomb);
3018 void forcerec_set_excl_load(t_forcerec *fr,
3019 const gmx_localtop_t *top, const t_commrec *cr)
3022 int t, i, j, ntot, n, ntarget;
3024 if (cr != NULL && PARTDECOMP(cr))
3026 /* No OpenMP with particle decomposition */
3034 ind = top->excls.index;
3038 for (i = 0; i < top->excls.nr; i++)
3040 for (j = ind[i]; j < ind[i+1]; j++)
3049 fr->excl_load[0] = 0;
3052 for (t = 1; t <= fr->nthreads; t++)
3054 ntarget = (ntot*t)/fr->nthreads;
3055 while (i < top->excls.nr && n < ntarget)
3057 for (j = ind[i]; j < ind[i+1]; j++)
3066 fr->excl_load[t] = i;