1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
50 #include "gmx_fatal.h"
51 #include "gmx_fatal_collective.h"
55 #include "nonbonded.h"
64 #include "md_support.h"
65 #include "md_logging.h"
70 #include "mtop_util.h"
71 #include "nbnxn_search.h"
72 #include "nbnxn_atomdata.h"
73 #include "nbnxn_consts.h"
75 #include "gmx_omp_nthreads.h"
78 /* MSVC definition for __cpuid() */
82 #include "types/nbnxn_cuda_types_ext.h"
83 #include "gpu_utils.h"
84 #include "nbnxn_cuda_data_mgmt.h"
85 #include "pmalloc_cuda.h"
87 t_forcerec *mk_forcerec(void)
97 static void pr_nbfp(FILE *fp,real *nbfp,gmx_bool bBHAM,int atnr)
101 for(i=0; (i<atnr); i++) {
102 for(j=0; (j<atnr); j++) {
103 fprintf(fp,"%2d - %2d",i,j);
105 fprintf(fp," a=%10g, b=%10g, c=%10g\n",BHAMA(nbfp,atnr,i,j),
106 BHAMB(nbfp,atnr,i,j),BHAMC(nbfp,atnr,i,j)/6.0);
108 fprintf(fp," c6=%10g, c12=%10g\n",C6(nbfp,atnr,i,j)/6.0,
109 C12(nbfp,atnr,i,j)/12.0);
115 static real *mk_nbfp(const gmx_ffparams_t *idef,gmx_bool bBHAM)
122 snew(nbfp,3*atnr*atnr);
123 for(i=k=0; (i<atnr); i++) {
124 for(j=0; (j<atnr); j++,k++) {
125 BHAMA(nbfp,atnr,i,j) = idef->iparams[k].bham.a;
126 BHAMB(nbfp,atnr,i,j) = idef->iparams[k].bham.b;
127 /* nbfp now includes the 6.0 derivative prefactor */
128 BHAMC(nbfp,atnr,i,j) = idef->iparams[k].bham.c*6.0;
133 snew(nbfp,2*atnr*atnr);
134 for(i=k=0; (i<atnr); i++) {
135 for(j=0; (j<atnr); j++,k++) {
136 /* nbfp now includes the 6.0/12.0 derivative prefactors */
137 C6(nbfp,atnr,i,j) = idef->iparams[k].lj.c6*6.0;
138 C12(nbfp,atnr,i,j) = idef->iparams[k].lj.c12*12.0;
146 /* This routine sets fr->solvent_opt to the most common solvent in the
147 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
148 * the fr->solvent_type array with the correct type (or esolNO).
150 * Charge groups that fulfill the conditions but are not identical to the
151 * most common one will be marked as esolNO in the solvent_type array.
153 * TIP3p is identical to SPC for these purposes, so we call it
154 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
156 * NOTE: QM particle should not
157 * become an optimized solvent. Not even if there is only one charge
167 } solvent_parameters_t;
170 check_solvent_cg(const gmx_moltype_t *molt,
173 const unsigned char *qm_grpnr,
174 const t_grps *qm_grps,
176 int *n_solvent_parameters,
177 solvent_parameters_t **solvent_parameters_p,
181 const t_blocka * excl;
192 solvent_parameters_t *solvent_parameters;
194 /* We use a list with parameters for each solvent type.
195 * Every time we discover a new molecule that fulfills the basic
196 * conditions for a solvent we compare with the previous entries
197 * in these lists. If the parameters are the same we just increment
198 * the counter for that type, and otherwise we create a new type
199 * based on the current molecule.
201 * Once we've finished going through all molecules we check which
202 * solvent is most common, and mark all those molecules while we
203 * clear the flag on all others.
206 solvent_parameters = *solvent_parameters_p;
208 /* Mark the cg first as non optimized */
211 /* Check if this cg has no exclusions with atoms in other charge groups
212 * and all atoms inside the charge group excluded.
213 * We only have 3 or 4 atom solvent loops.
215 if (GET_CGINFO_EXCL_INTER(cginfo) ||
216 !GET_CGINFO_EXCL_INTRA(cginfo))
221 /* Get the indices of the first atom in this charge group */
222 j0 = molt->cgs.index[cg0];
223 j1 = molt->cgs.index[cg0+1];
225 /* Number of atoms in our molecule */
230 "Moltype '%s': there are %d atoms in this charge group\n",
234 /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
242 /* Check if we are doing QM on this group */
244 if (qm_grpnr != NULL)
246 for(j=j0 ; j<j1 && !qm; j++)
248 qm = (qm_grpnr[j] < qm_grps->nr - 1);
251 /* Cannot use solvent optimization with QM */
257 atom = molt->atoms.atom;
259 /* Still looks like a solvent, time to check parameters */
261 /* If it is perturbed (free energy) we can't use the solvent loops,
262 * so then we just skip to the next molecule.
266 for(j=j0; j<j1 && !perturbed; j++)
268 perturbed = PERTURBED(atom[j]);
276 /* Now it's only a question if the VdW and charge parameters
277 * are OK. Before doing the check we compare and see if they are
278 * identical to a possible previous solvent type.
279 * First we assign the current types and charges.
283 tmp_vdwtype[j] = atom[j0+j].type;
284 tmp_charge[j] = atom[j0+j].q;
287 /* Does it match any previous solvent type? */
288 for(k=0 ; k<*n_solvent_parameters; k++)
293 /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
294 if( (solvent_parameters[k].model==esolSPC && nj!=3) ||
295 (solvent_parameters[k].model==esolTIP4P && nj!=4) )
298 /* Check that types & charges match for all atoms in molecule */
299 for(j=0 ; j<nj && match==TRUE; j++)
301 if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
305 if(tmp_charge[j] != solvent_parameters[k].charge[j])
312 /* Congratulations! We have a matched solvent.
313 * Flag it with this type for later processing.
316 solvent_parameters[k].count += nmol;
318 /* We are done with this charge group */
323 /* If we get here, we have a tentative new solvent type.
324 * Before we add it we must check that it fulfills the requirements
325 * of the solvent optimized loops. First determine which atoms have
331 tjA = tmp_vdwtype[j];
333 /* Go through all other tpes and see if any have non-zero
334 * VdW parameters when combined with this one.
336 for(k=0; k<fr->ntype && (has_vdw[j]==FALSE); k++)
338 /* We already checked that the atoms weren't perturbed,
339 * so we only need to check state A now.
343 has_vdw[j] = (has_vdw[j] ||
344 (BHAMA(fr->nbfp,fr->ntype,tjA,k) != 0.0) ||
345 (BHAMB(fr->nbfp,fr->ntype,tjA,k) != 0.0) ||
346 (BHAMC(fr->nbfp,fr->ntype,tjA,k) != 0.0));
351 has_vdw[j] = (has_vdw[j] ||
352 (C6(fr->nbfp,fr->ntype,tjA,k) != 0.0) ||
353 (C12(fr->nbfp,fr->ntype,tjA,k) != 0.0));
358 /* Now we know all we need to make the final check and assignment. */
362 * For this we require thatn all atoms have charge,
363 * the charges on atom 2 & 3 should be the same, and only
364 * atom 1 might have VdW.
366 if (has_vdw[1] == FALSE &&
367 has_vdw[2] == FALSE &&
368 tmp_charge[0] != 0 &&
369 tmp_charge[1] != 0 &&
370 tmp_charge[2] == tmp_charge[1])
372 srenew(solvent_parameters,*n_solvent_parameters+1);
373 solvent_parameters[*n_solvent_parameters].model = esolSPC;
374 solvent_parameters[*n_solvent_parameters].count = nmol;
377 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
378 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
381 *cg_sp = *n_solvent_parameters;
382 (*n_solvent_parameters)++;
387 /* Or could it be a TIP4P?
388 * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
389 * Only atom 1 mght have VdW.
391 if(has_vdw[1] == FALSE &&
392 has_vdw[2] == FALSE &&
393 has_vdw[3] == FALSE &&
394 tmp_charge[0] == 0 &&
395 tmp_charge[1] != 0 &&
396 tmp_charge[2] == tmp_charge[1] &&
399 srenew(solvent_parameters,*n_solvent_parameters+1);
400 solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
401 solvent_parameters[*n_solvent_parameters].count = nmol;
404 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
405 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
408 *cg_sp = *n_solvent_parameters;
409 (*n_solvent_parameters)++;
413 *solvent_parameters_p = solvent_parameters;
417 check_solvent(FILE * fp,
418 const gmx_mtop_t * mtop,
420 cginfo_mb_t *cginfo_mb)
423 const t_block * mols;
424 const gmx_moltype_t *molt;
425 int mb,mol,cg_mol,at_offset,cg_offset,am,cgm,i,nmol_ch,nmol;
426 int n_solvent_parameters;
427 solvent_parameters_t *solvent_parameters;
433 fprintf(debug,"Going to determine what solvent types we have.\n");
438 n_solvent_parameters = 0;
439 solvent_parameters = NULL;
440 /* Allocate temporary array for solvent type */
441 snew(cg_sp,mtop->nmolblock);
445 for(mb=0; mb<mtop->nmolblock; mb++)
447 molt = &mtop->moltype[mtop->molblock[mb].type];
449 /* Here we have to loop over all individual molecules
450 * because we need to check for QMMM particles.
452 snew(cg_sp[mb],cginfo_mb[mb].cg_mod);
453 nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
454 nmol = mtop->molblock[mb].nmol/nmol_ch;
455 for(mol=0; mol<nmol_ch; mol++)
458 am = mol*cgs->index[cgs->nr];
459 for(cg_mol=0; cg_mol<cgs->nr; cg_mol++)
461 check_solvent_cg(molt,cg_mol,nmol,
462 mtop->groups.grpnr[egcQMMM] ?
463 mtop->groups.grpnr[egcQMMM]+at_offset+am : 0,
464 &mtop->groups.grps[egcQMMM],
466 &n_solvent_parameters,&solvent_parameters,
467 cginfo_mb[mb].cginfo[cgm+cg_mol],
468 &cg_sp[mb][cgm+cg_mol]);
471 cg_offset += cgs->nr;
472 at_offset += cgs->index[cgs->nr];
475 /* Puh! We finished going through all charge groups.
476 * Now find the most common solvent model.
479 /* Most common solvent this far */
481 for(i=0;i<n_solvent_parameters;i++)
484 solvent_parameters[i].count > solvent_parameters[bestsp].count)
492 bestsol = solvent_parameters[bestsp].model;
499 #ifdef DISABLE_WATER_NLIST
504 for(mb=0; mb<mtop->nmolblock; mb++)
506 cgs = &mtop->moltype[mtop->molblock[mb].type].cgs;
507 nmol = (mtop->molblock[mb].nmol*cgs->nr)/cginfo_mb[mb].cg_mod;
508 for(i=0; i<cginfo_mb[mb].cg_mod; i++)
510 if (cg_sp[mb][i] == bestsp)
512 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i],bestsol);
517 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i],esolNO);
524 if (bestsol != esolNO && fp!=NULL)
526 fprintf(fp,"\nEnabling %s-like water optimization for %d molecules.\n\n",
528 solvent_parameters[bestsp].count);
531 sfree(solvent_parameters);
532 fr->solvent_opt = bestsol;
535 enum { acNONE=0, acCONSTRAINT, acSETTLE };
537 static cginfo_mb_t *init_cginfo_mb(FILE *fplog,const gmx_mtop_t *mtop,
538 t_forcerec *fr,gmx_bool bNoSolvOpt,
539 gmx_bool *bExcl_IntraCGAll_InterCGNone)
542 const t_blocka *excl;
543 const gmx_moltype_t *molt;
544 const gmx_molblock_t *molb;
545 cginfo_mb_t *cginfo_mb;
548 int cg_offset,a_offset,cgm,am;
549 int mb,m,ncg_tot,cg,a0,a1,gid,ai,j,aj,excl_nalloc;
553 gmx_bool bId,*bExcl,bExclIntraAll,bExclInter,bHaveVDW,bHaveQ;
555 ncg_tot = ncg_mtop(mtop);
556 snew(cginfo_mb,mtop->nmolblock);
558 snew(type_VDW,fr->ntype);
559 for(ai=0; ai<fr->ntype; ai++)
561 type_VDW[ai] = FALSE;
562 for(j=0; j<fr->ntype; j++)
564 type_VDW[ai] = type_VDW[ai] ||
566 C6(fr->nbfp,fr->ntype,ai,j) != 0 ||
567 C12(fr->nbfp,fr->ntype,ai,j) != 0;
571 *bExcl_IntraCGAll_InterCGNone = TRUE;
574 snew(bExcl,excl_nalloc);
577 for(mb=0; mb<mtop->nmolblock; mb++)
579 molb = &mtop->molblock[mb];
580 molt = &mtop->moltype[molb->type];
584 /* Check if the cginfo is identical for all molecules in this block.
585 * If so, we only need an array of the size of one molecule.
586 * Otherwise we make an array of #mol times #cgs per molecule.
590 for(m=0; m<molb->nmol; m++)
592 am = m*cgs->index[cgs->nr];
593 for(cg=0; cg<cgs->nr; cg++)
596 a1 = cgs->index[cg+1];
597 if (ggrpnr(&mtop->groups,egcENER,a_offset+am+a0) !=
598 ggrpnr(&mtop->groups,egcENER,a_offset +a0))
602 if (mtop->groups.grpnr[egcQMMM] != NULL)
604 for(ai=a0; ai<a1; ai++)
606 if (mtop->groups.grpnr[egcQMMM][a_offset+am+ai] !=
607 mtop->groups.grpnr[egcQMMM][a_offset +ai])
616 cginfo_mb[mb].cg_start = cg_offset;
617 cginfo_mb[mb].cg_end = cg_offset + molb->nmol*cgs->nr;
618 cginfo_mb[mb].cg_mod = (bId ? 1 : molb->nmol)*cgs->nr;
619 snew(cginfo_mb[mb].cginfo,cginfo_mb[mb].cg_mod);
620 cginfo = cginfo_mb[mb].cginfo;
622 /* Set constraints flags for constrained atoms */
623 snew(a_con,molt->atoms.nr);
624 for(ftype=0; ftype<F_NRE; ftype++)
626 if (interaction_function[ftype].flags & IF_CONSTRAINT)
631 for(ia=0; ia<molt->ilist[ftype].nr; ia+=1+nral)
635 for(a=0; a<nral; a++)
637 a_con[molt->ilist[ftype].iatoms[ia+1+a]] =
638 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
644 for(m=0; m<(bId ? 1 : molb->nmol); m++)
647 am = m*cgs->index[cgs->nr];
648 for(cg=0; cg<cgs->nr; cg++)
651 a1 = cgs->index[cg+1];
653 /* Store the energy group in cginfo */
654 gid = ggrpnr(&mtop->groups,egcENER,a_offset+am+a0);
655 SET_CGINFO_GID(cginfo[cgm+cg],gid);
657 /* Check the intra/inter charge group exclusions */
658 if (a1-a0 > excl_nalloc) {
659 excl_nalloc = a1 - a0;
660 srenew(bExcl,excl_nalloc);
662 /* bExclIntraAll: all intra cg interactions excluded
663 * bExclInter: any inter cg interactions excluded
665 bExclIntraAll = TRUE;
669 for(ai=a0; ai<a1; ai++)
671 /* Check VDW and electrostatic interactions */
672 bHaveVDW = bHaveVDW || (type_VDW[molt->atoms.atom[ai].type] ||
673 type_VDW[molt->atoms.atom[ai].typeB]);
674 bHaveQ = bHaveQ || (molt->atoms.atom[ai].q != 0 ||
675 molt->atoms.atom[ai].qB != 0);
677 /* Clear the exclusion list for atom ai */
678 for(aj=a0; aj<a1; aj++)
680 bExcl[aj-a0] = FALSE;
682 /* Loop over all the exclusions of atom ai */
683 for(j=excl->index[ai]; j<excl->index[ai+1]; j++)
686 if (aj < a0 || aj >= a1)
695 /* Check if ai excludes a0 to a1 */
696 for(aj=a0; aj<a1; aj++)
700 bExclIntraAll = FALSE;
707 SET_CGINFO_CONSTR(cginfo[cgm+cg]);
710 SET_CGINFO_SETTLE(cginfo[cgm+cg]);
718 SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
722 SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
724 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
726 /* The size in cginfo is currently only read with DD */
727 gmx_fatal(FARGS,"A charge group has size %d which is larger than the limit of %d atoms",a1-a0,MAX_CHARGEGROUP_SIZE);
731 SET_CGINFO_HAS_VDW(cginfo[cgm+cg]);
735 SET_CGINFO_HAS_Q(cginfo[cgm+cg]);
737 /* Store the charge group size */
738 SET_CGINFO_NATOMS(cginfo[cgm+cg],a1-a0);
740 if (!bExclIntraAll || bExclInter)
742 *bExcl_IntraCGAll_InterCGNone = FALSE;
749 cg_offset += molb->nmol*cgs->nr;
750 a_offset += molb->nmol*cgs->index[cgs->nr];
754 /* the solvent optimizer is called after the QM is initialized,
755 * because we don't want to have the QM subsystemto become an
759 check_solvent(fplog,mtop,fr,cginfo_mb);
761 if (getenv("GMX_NO_SOLV_OPT"))
765 fprintf(fplog,"Found environment variable GMX_NO_SOLV_OPT.\n"
766 "Disabling all solvent optimization\n");
768 fr->solvent_opt = esolNO;
772 fr->solvent_opt = esolNO;
774 if (!fr->solvent_opt)
776 for(mb=0; mb<mtop->nmolblock; mb++)
778 for(cg=0; cg<cginfo_mb[mb].cg_mod; cg++)
780 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg],esolNO);
788 static int *cginfo_expand(int nmb,cginfo_mb_t *cgi_mb)
793 ncg = cgi_mb[nmb-1].cg_end;
796 for(cg=0; cg<ncg; cg++)
798 while (cg >= cgi_mb[mb].cg_end)
803 cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
809 static void set_chargesum(FILE *log,t_forcerec *fr,const gmx_mtop_t *mtop)
813 const t_atoms *atoms;
817 for(mb=0; mb<mtop->nmolblock; mb++)
819 nmol = mtop->molblock[mb].nmol;
820 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
821 for(i=0; i<atoms->nr; i++)
823 q = atoms->atom[i].q;
829 fr->q2sum[0] = q2sum;
830 if (fr->efep != efepNO)
834 for(mb=0; mb<mtop->nmolblock; mb++)
836 nmol = mtop->molblock[mb].nmol;
837 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
838 for(i=0; i<atoms->nr; i++)
840 q = atoms->atom[i].qB;
845 fr->q2sum[1] = q2sum;
850 fr->qsum[1] = fr->qsum[0];
851 fr->q2sum[1] = fr->q2sum[0];
854 if (fr->efep == efepNO)
855 fprintf(log,"System total charge: %.3f\n",fr->qsum[0]);
857 fprintf(log,"System total charge, top. A: %.3f top. B: %.3f\n",
858 fr->qsum[0],fr->qsum[1]);
862 void update_forcerec(FILE *log,t_forcerec *fr,matrix box)
864 if (fr->eeltype == eelGRF)
866 calc_rffac(NULL,fr->eeltype,fr->epsilon_r,fr->epsilon_rf,
867 fr->rcoulomb,fr->temp,fr->zsquare,box,
868 &fr->kappa,&fr->k_rf,&fr->c_rf);
872 void set_avcsixtwelve(FILE *fplog,t_forcerec *fr,const gmx_mtop_t *mtop)
874 const t_atoms *atoms,*atoms_tpi;
875 const t_blocka *excl;
876 int mb,nmol,nmolc,i,j,tpi,tpj,j1,j2,k,n,nexcl,q;
877 #if (defined SIZEOF_LONG_LONG_INT) && (SIZEOF_LONG_LONG_INT >= 8)
878 long long int npair,npair_ij,tmpi,tmpj;
880 double npair, npair_ij,tmpi,tmpj;
891 for(q=0; q<(fr->efep==efepNO ? 1 : 2); q++) {
897 /* Count the types so we avoid natoms^2 operations */
899 for(mb=0; mb<mtop->nmolblock; mb++) {
900 nmol = mtop->molblock[mb].nmol;
901 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
902 for(i=0; i<atoms->nr; i++) {
905 tpi = atoms->atom[i].type;
909 tpi = atoms->atom[i].typeB;
911 typecount[tpi] += nmol;
914 for(tpi=0; tpi<ntp; tpi++) {
915 for(tpj=tpi; tpj<ntp; tpj++) {
916 tmpi = typecount[tpi];
917 tmpj = typecount[tpj];
920 npair_ij = tmpi*tmpj;
924 npair_ij = tmpi*(tmpi - 1)/2;
927 /* nbfp now includes the 6.0 derivative prefactor */
928 csix += npair_ij*BHAMC(nbfp,ntp,tpi,tpj)/6.0;
930 /* nbfp now includes the 6.0/12.0 derivative prefactors */
931 csix += npair_ij* C6(nbfp,ntp,tpi,tpj)/6.0;
932 ctwelve += npair_ij* C12(nbfp,ntp,tpi,tpj)/12.0;
938 /* Subtract the excluded pairs.
939 * The main reason for substracting exclusions is that in some cases
940 * some combinations might never occur and the parameters could have
941 * any value. These unused values should not influence the dispersion
944 for(mb=0; mb<mtop->nmolblock; mb++) {
945 nmol = mtop->molblock[mb].nmol;
946 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
947 excl = &mtop->moltype[mtop->molblock[mb].type].excls;
948 for(i=0; (i<atoms->nr); i++) {
951 tpi = atoms->atom[i].type;
955 tpi = atoms->atom[i].typeB;
958 j2 = excl->index[i+1];
959 for(j=j1; j<j2; j++) {
965 tpj = atoms->atom[k].type;
969 tpj = atoms->atom[k].typeB;
972 /* nbfp now includes the 6.0 derivative prefactor */
973 csix -= nmol*BHAMC(nbfp,ntp,tpi,tpj)/6.0;
975 /* nbfp now includes the 6.0/12.0 derivative prefactors */
976 csix -= nmol*C6 (nbfp,ntp,tpi,tpj)/6.0;
977 ctwelve -= nmol*C12(nbfp,ntp,tpi,tpj)/12.0;
985 /* Only correct for the interaction of the test particle
986 * with the rest of the system.
989 &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].atoms;
992 for(mb=0; mb<mtop->nmolblock; mb++) {
993 nmol = mtop->molblock[mb].nmol;
994 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
995 for(j=0; j<atoms->nr; j++) {
997 /* Remove the interaction of the test charge group
1000 if (mb == mtop->nmolblock-1)
1004 if (mb == 0 && nmol == 1)
1006 gmx_fatal(FARGS,"Old format tpr with TPI, please generate a new tpr file");
1011 tpj = atoms->atom[j].type;
1015 tpj = atoms->atom[j].typeB;
1017 for(i=0; i<fr->n_tpi; i++)
1021 tpi = atoms_tpi->atom[i].type;
1025 tpi = atoms_tpi->atom[i].typeB;
1029 /* nbfp now includes the 6.0 derivative prefactor */
1030 csix += nmolc*BHAMC(nbfp,ntp,tpi,tpj)/6.0;
1034 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1035 csix += nmolc*C6 (nbfp,ntp,tpi,tpj)/6.0;
1036 ctwelve += nmolc*C12(nbfp,ntp,tpi,tpj)/12.0;
1043 if (npair - nexcl <= 0 && fplog) {
1044 fprintf(fplog,"\nWARNING: There are no atom pairs for dispersion correction\n\n");
1048 csix /= npair - nexcl;
1049 ctwelve /= npair - nexcl;
1052 fprintf(debug,"Counted %d exclusions\n",nexcl);
1053 fprintf(debug,"Average C6 parameter is: %10g\n",(double)csix);
1054 fprintf(debug,"Average C12 parameter is: %10g\n",(double)ctwelve);
1056 fr->avcsix[q] = csix;
1057 fr->avctwelve[q] = ctwelve;
1061 if (fr->eDispCorr == edispcAllEner ||
1062 fr->eDispCorr == edispcAllEnerPres)
1064 fprintf(fplog,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1065 fr->avcsix[0],fr->avctwelve[0]);
1069 fprintf(fplog,"Long Range LJ corr.: <C6> %10.4e\n",fr->avcsix[0]);
1075 static void set_bham_b_max(FILE *fplog,t_forcerec *fr,
1076 const gmx_mtop_t *mtop)
1078 const t_atoms *at1,*at2;
1079 int mt1,mt2,i,j,tpi,tpj,ntypes;
1085 fprintf(fplog,"Determining largest Buckingham b parameter for table\n");
1092 for(mt1=0; mt1<mtop->nmoltype; mt1++)
1094 at1 = &mtop->moltype[mt1].atoms;
1095 for(i=0; (i<at1->nr); i++)
1097 tpi = at1->atom[i].type;
1099 gmx_fatal(FARGS,"Atomtype[%d] = %d, maximum = %d",i,tpi,ntypes);
1101 for(mt2=mt1; mt2<mtop->nmoltype; mt2++)
1103 at2 = &mtop->moltype[mt2].atoms;
1104 for(j=0; (j<at2->nr); j++) {
1105 tpj = at2->atom[j].type;
1108 gmx_fatal(FARGS,"Atomtype[%d] = %d, maximum = %d",j,tpj,ntypes);
1110 b = BHAMB(nbfp,ntypes,tpi,tpj);
1111 if (b > fr->bham_b_max)
1115 if ((b < bmin) || (bmin==-1))
1125 fprintf(fplog,"Buckingham b parameters, min: %g, max: %g\n",
1126 bmin,fr->bham_b_max);
1130 static void make_nbf_tables(FILE *fp,const output_env_t oenv,
1131 t_forcerec *fr,real rtab,
1132 const t_commrec *cr,
1133 const char *tabfn,char *eg1,char *eg2,
1139 if (tabfn == NULL) {
1141 fprintf(debug,"No table file name passed, can not read table, can not do non-bonded interactions\n");
1145 sprintf(buf,"%s",tabfn);
1147 /* Append the two energy group names */
1148 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1,"_%s_%s.%s",
1149 eg1,eg2,ftp2ext(efXVG));
1150 nbl->table_elec_vdw = make_tables(fp,oenv,fr,MASTER(cr),buf,rtab,0);
1151 /* Copy the contents of the table to separate coulomb and LJ tables too,
1152 * to improve cache performance.
1154 /* For performance reasons we want
1155 * the table data to be aligned to 16-byte. The pointers could be freed
1156 * but currently aren't.
1158 nbl->table_elec.interaction = GMX_TABLE_INTERACTION_ELEC;
1159 nbl->table_elec.format = nbl->table_elec_vdw.format;
1160 nbl->table_elec.r = nbl->table_elec_vdw.r;
1161 nbl->table_elec.n = nbl->table_elec_vdw.n;
1162 nbl->table_elec.scale = nbl->table_elec_vdw.scale;
1163 nbl->table_elec.scale_exp = nbl->table_elec_vdw.scale_exp;
1164 nbl->table_elec.formatsize = nbl->table_elec_vdw.formatsize;
1165 nbl->table_elec.ninteractions = 1;
1166 nbl->table_elec.stride = nbl->table_elec.formatsize * nbl->table_elec.ninteractions;
1167 snew_aligned(nbl->table_elec.data,nbl->table_elec.stride*(nbl->table_elec.n+1),16);
1169 nbl->table_vdw.interaction = GMX_TABLE_INTERACTION_VDWREP_VDWDISP;
1170 nbl->table_vdw.format = nbl->table_elec_vdw.format;
1171 nbl->table_vdw.r = nbl->table_elec_vdw.r;
1172 nbl->table_vdw.n = nbl->table_elec_vdw.n;
1173 nbl->table_vdw.scale = nbl->table_elec_vdw.scale;
1174 nbl->table_vdw.scale_exp = nbl->table_elec_vdw.scale_exp;
1175 nbl->table_vdw.formatsize = nbl->table_elec_vdw.formatsize;
1176 nbl->table_vdw.ninteractions = 2;
1177 nbl->table_vdw.stride = nbl->table_vdw.formatsize * nbl->table_vdw.ninteractions;
1178 snew_aligned(nbl->table_vdw.data,nbl->table_vdw.stride*(nbl->table_vdw.n+1),16);
1180 for(i=0; i<=nbl->table_elec_vdw.n; i++)
1183 nbl->table_elec.data[4*i+j] = nbl->table_elec_vdw.data[12*i+j];
1185 nbl->table_vdw.data[8*i+j] = nbl->table_elec_vdw.data[12*i+4+j];
1189 static void count_tables(int ftype1,int ftype2,const gmx_mtop_t *mtop,
1190 int *ncount,int **count)
1192 const gmx_moltype_t *molt;
1194 int mt,ftype,stride,i,j,tabnr;
1196 for(mt=0; mt<mtop->nmoltype; mt++)
1198 molt = &mtop->moltype[mt];
1199 for(ftype=0; ftype<F_NRE; ftype++)
1201 if (ftype == ftype1 || ftype == ftype2) {
1202 il = &molt->ilist[ftype];
1203 stride = 1 + NRAL(ftype);
1204 for(i=0; i<il->nr; i+=stride) {
1205 tabnr = mtop->ffparams.iparams[il->iatoms[i]].tab.table;
1207 gmx_fatal(FARGS,"A bonded table number is smaller than 0: %d\n",tabnr);
1208 if (tabnr >= *ncount) {
1209 srenew(*count,tabnr+1);
1210 for(j=*ncount; j<tabnr+1; j++)
1221 static bondedtable_t *make_bonded_tables(FILE *fplog,
1222 int ftype1,int ftype2,
1223 const gmx_mtop_t *mtop,
1224 const char *basefn,const char *tabext)
1226 int i,ncount,*count;
1234 count_tables(ftype1,ftype2,mtop,&ncount,&count);
1238 for(i=0; i<ncount; i++) {
1240 sprintf(tabfn,"%s",basefn);
1241 sprintf(tabfn + strlen(basefn) - strlen(ftp2ext(efXVG)) - 1,"_%s%d.%s",
1242 tabext,i,ftp2ext(efXVG));
1243 tab[i] = make_bonded_table(fplog,tabfn,NRAL(ftype1)-2);
1252 void forcerec_set_ranges(t_forcerec *fr,
1253 int ncg_home,int ncg_force,
1255 int natoms_force_constr,int natoms_f_novirsum)
1260 /* fr->ncg_force is unused in the standard code,
1261 * but it can be useful for modified code dealing with charge groups.
1263 fr->ncg_force = ncg_force;
1264 fr->natoms_force = natoms_force;
1265 fr->natoms_force_constr = natoms_force_constr;
1267 if (fr->natoms_force_constr > fr->nalloc_force)
1269 fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1273 srenew(fr->f_twin,fr->nalloc_force);
1277 if (fr->bF_NoVirSum)
1279 fr->f_novirsum_n = natoms_f_novirsum;
1280 if (fr->f_novirsum_n > fr->f_novirsum_nalloc)
1282 fr->f_novirsum_nalloc = over_alloc_dd(fr->f_novirsum_n);
1283 srenew(fr->f_novirsum_alloc,fr->f_novirsum_nalloc);
1288 fr->f_novirsum_n = 0;
1292 static real cutoff_inf(real cutoff)
1296 cutoff = GMX_CUTOFF_INF;
1302 static void make_adress_tf_tables(FILE *fp,const output_env_t oenv,
1303 t_forcerec *fr,const t_inputrec *ir,
1304 const char *tabfn, const gmx_mtop_t *mtop,
1310 if (tabfn == NULL) {
1311 gmx_fatal(FARGS,"No thermoforce table file given. Use -tabletf to specify a file\n");
1315 snew(fr->atf_tabs, ir->adress->n_tf_grps);
1317 for (i=0; i<ir->adress->n_tf_grps; i++){
1318 j = ir->adress->tf_table_index[i]; /* get energy group index */
1319 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1,"tf_%s.%s",
1320 *(mtop->groups.grpname[mtop->groups.grps[egcENER].nm_ind[j]]) ,ftp2ext(efXVG));
1321 printf("loading tf table for energygrp index %d from %s\n", ir->adress->tf_table_index[j], buf);
1322 fr->atf_tabs[i] = make_atf_table(fp,oenv,fr,buf, box);
1327 gmx_bool can_use_allvsall(const t_inputrec *ir, const gmx_mtop_t *mtop,
1328 gmx_bool bPrintNote,t_commrec *cr,FILE *fp)
1337 ir->ePBC==epbcNONE &&
1338 ir->vdwtype==evdwCUT &&
1339 ir->coulombtype==eelCUT &&
1341 (ir->implicit_solvent == eisNO ||
1342 (ir->implicit_solvent==eisGBSA && (ir->gb_algorithm==egbSTILL ||
1343 ir->gb_algorithm==egbHCT ||
1344 ir->gb_algorithm==egbOBC))) &&
1345 getenv("GMX_NO_ALLVSALL") == NULL
1348 if (bAllvsAll && ir->opts.ngener > 1)
1350 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";
1356 fprintf(stderr,"\n%s\n",note);
1360 fprintf(fp,"\n%s\n",note);
1366 if(bAllvsAll && fp && MASTER(cr))
1368 fprintf(fp,"\nUsing accelerated all-vs-all kernels.\n\n");
1375 static void init_forcerec_f_threads(t_forcerec *fr,int nenergrp)
1379 /* These thread local data structures are used for bondeds only */
1380 fr->nthreads = gmx_omp_nthreads_get(emntBonded);
1382 if (fr->nthreads > 1)
1384 snew(fr->f_t,fr->nthreads);
1385 /* Thread 0 uses the global force and energy arrays */
1386 for(t=1; t<fr->nthreads; t++)
1388 fr->f_t[t].f = NULL;
1389 fr->f_t[t].f_nalloc = 0;
1390 snew(fr->f_t[t].fshift,SHIFTS);
1391 fr->f_t[t].grpp.nener = nenergrp*nenergrp;
1392 for(i=0; i<egNR; i++)
1394 snew(fr->f_t[t].grpp.ener[i],fr->f_t[t].grpp.nener);
1401 static void pick_nbnxn_kernel_cpu(FILE *fp,
1402 const t_commrec *cr,
1403 const gmx_cpuid_t cpuid_info,
1407 *kernel_type = nbk4x4_PlainC;
1408 *ewald_excl = ewaldexclTable;
1412 /* On Intel Sandy-Bridge AVX-256 kernels are always faster.
1413 * On AMD Bulldozer AVX-256 is much slower than AVX-128.
1415 if(gmx_cpuid_feature(cpuid_info, GMX_CPUID_FEATURE_X86_AVX) == 1 &&
1416 gmx_cpuid_vendor(cpuid_info) != GMX_CPUID_VENDOR_AMD)
1418 #ifdef GMX_X86_AVX_256
1419 *kernel_type = nbk4xN_X86_SIMD256;
1421 *kernel_type = nbk4xN_X86_SIMD128;
1426 *kernel_type = nbk4xN_X86_SIMD128;
1429 if (getenv("GMX_NBNXN_AVX128") != NULL)
1431 *kernel_type = nbk4xN_X86_SIMD128;
1433 if (getenv("GMX_NBNXN_AVX256") != NULL)
1435 #ifdef GMX_X86_AVX_256
1436 *kernel_type = nbk4xN_X86_SIMD256;
1438 gmx_fatal(FARGS,"You requested AVX-256 nbnxn kernels, but GROMACS was built without AVX support");
1442 /* Analytical Ewald exclusion correction is only an option in the
1443 * x86 SIMD kernel. This is faster in single precision
1444 * on Bulldozer and slightly faster on Sandy Bridge.
1446 #if (defined GMX_X86_AVX_128_FMA || defined GMX_X86_AVX_256) && !defined GMX_DOUBLE
1447 *ewald_excl = ewaldexclAnalytical;
1449 if (getenv("GMX_NBNXN_EWALD_TABLE") != NULL)
1451 *ewald_excl = ewaldexclTable;
1453 if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != NULL)
1455 *ewald_excl = ewaldexclAnalytical;
1459 #endif /* GMX_X86_SSE2 */
1462 static void pick_nbnxn_kernel(FILE *fp,
1463 const t_commrec *cr,
1464 const gmx_hw_info_t *hwinfo,
1465 gmx_bool use_cpu_acceleration,
1470 gmx_bool bEmulateGPU, bGPU, bEmulateGPUEnvVarSet;
1471 char gpu_err_str[STRLEN];
1473 assert(kernel_type);
1475 *kernel_type = nbkNotSet;
1476 *ewald_excl = ewaldexclTable;
1478 bEmulateGPUEnvVarSet = (getenv("GMX_EMULATE_GPU") != NULL);
1480 /* if bUseGPU == NULL we don't want a GPU (e.g. hybrid mode kernel selection) */
1481 bGPU = ((bUseGPU != NULL) && hwinfo->bCanUseGPU);
1483 /* Run GPU emulation mode if GMX_EMULATE_GPU is defined. We will
1484 * automatically switch to emulation if non-bonded calculations are
1485 * turned off via GMX_NO_NONBONDED - this is the simple and elegant
1486 * way to turn off GPU initialization, data movement, and cleanup. */
1487 bEmulateGPU = (bEmulateGPUEnvVarSet || (getenv("GMX_NO_NONBONDED") != NULL && bGPU));
1489 /* Enable GPU mode when GPUs are available or GPU emulation is requested.
1490 * The latter is useful to assess the performance one can expect by adding
1491 * GPU(s) to the machine. The conditional below allows this even if mdrun
1492 * is compiled without GPU acceleration support.
1493 * Note that such a GPU acceleration performance assessment should be
1494 * carried out by setting the GMX_EMULATE_GPU and GMX_NO_NONBONDED env. vars
1495 * (and freezing the system as otherwise it would explode). */
1496 if (bGPU || bEmulateGPUEnvVarSet)
1504 /* Each PP node will use the intra-node id-th device from the
1505 * list of detected/selected GPUs. */
1506 if (!init_gpu(cr->nodeid_group_intra, gpu_err_str, &hwinfo->gpu_info))
1508 /* At this point the init should never fail as we made sure that
1509 * we have all the GPUs we need. If it still does, we'll bail. */
1510 gmx_fatal(FARGS, "On node %d failed to initialize GPU #%d: %s",
1512 get_gpu_device_id(&hwinfo->gpu_info, cr->nodeid_group_intra),
1521 *kernel_type = nbk8x8x8_PlainC;
1523 md_print_warn(cr, fp, "Emulating a GPU run on the CPU (slow)");
1527 *kernel_type = nbk8x8x8_CUDA;
1530 if (*kernel_type == nbkNotSet)
1532 if (use_cpu_acceleration)
1534 pick_nbnxn_kernel_cpu(fp,cr,hwinfo->cpuid_info,
1535 kernel_type,ewald_excl);
1539 *kernel_type = nbk4x4_PlainC;
1547 fprintf(stderr,"Using %s non-bonded kernels\n",
1548 nbk_name[*kernel_type]);
1550 fprintf(fp,"\nUsing %s non-bonded kernels\n\n",
1551 nbk_name[*kernel_type]);
1555 gmx_bool uses_simple_tables(int cutoff_scheme,
1556 nonbonded_verlet_t *nbv,
1559 gmx_bool bUsesSimpleTables = TRUE;
1562 switch(cutoff_scheme)
1565 bUsesSimpleTables = TRUE;
1568 assert(NULL != nbv && NULL != nbv->grp);
1569 grp_index = (group < 0) ? 0 : (nbv->ngrp - 1);
1570 bUsesSimpleTables = nbnxn_kernel_pairlist_simple(nbv->grp[grp_index].kernel_type);
1573 gmx_incons("unimplemented");
1575 return bUsesSimpleTables;
1578 static void init_ewald_f_table(interaction_const_t *ic,
1579 gmx_bool bUsesSimpleTables,
1584 if (bUsesSimpleTables)
1586 /* With a spacing of 0.0005 we are at the force summation accuracy
1587 * for the SSE kernels for "normal" atomistic simulations.
1589 ic->tabq_scale = ewald_spline3_table_scale(ic->ewaldcoeff,
1592 maxr = (rtab>ic->rcoulomb) ? rtab : ic->rcoulomb;
1593 ic->tabq_size = (int)(maxr*ic->tabq_scale) + 2;
1597 ic->tabq_size = GPU_EWALD_COULOMB_FORCE_TABLE_SIZE;
1598 /* Subtract 2 iso 1 to avoid access out of range due to rounding */
1599 ic->tabq_scale = (ic->tabq_size - 2)/ic->rcoulomb;
1602 sfree_aligned(ic->tabq_coul_FDV0);
1603 sfree_aligned(ic->tabq_coul_F);
1604 sfree_aligned(ic->tabq_coul_V);
1606 /* Create the original table data in FDV0 */
1607 snew_aligned(ic->tabq_coul_FDV0,ic->tabq_size*4,16);
1608 snew_aligned(ic->tabq_coul_F,ic->tabq_size,16);
1609 snew_aligned(ic->tabq_coul_V,ic->tabq_size,16);
1610 table_spline3_fill_ewald_lr(ic->tabq_coul_F,ic->tabq_coul_V,ic->tabq_coul_FDV0,
1611 ic->tabq_size,1/ic->tabq_scale,ic->ewaldcoeff);
1614 void init_interaction_const_tables(FILE *fp,
1615 interaction_const_t *ic,
1616 gmx_bool bUsesSimpleTables,
1621 if (ic->eeltype == eelEWALD || EEL_PME(ic->eeltype))
1623 init_ewald_f_table(ic,bUsesSimpleTables,rtab);
1627 fprintf(fp,"Initialized non-bonded Ewald correction tables, spacing: %.2e size: %d\n\n",
1628 1/ic->tabq_scale,ic->tabq_size);
1633 void init_interaction_const(FILE *fp,
1634 interaction_const_t **interaction_const,
1635 const t_forcerec *fr,
1638 interaction_const_t *ic;
1639 gmx_bool bUsesSimpleTables = TRUE;
1643 /* Just allocate something so we can free it */
1644 snew_aligned(ic->tabq_coul_FDV0,16,16);
1645 snew_aligned(ic->tabq_coul_F,16,16);
1646 snew_aligned(ic->tabq_coul_V,16,16);
1648 ic->rlist = fr->rlist;
1649 ic->rlistlong = fr->rlistlong;
1652 ic->rvdw = fr->rvdw;
1653 if (fr->vdw_modifier==eintmodPOTSHIFT)
1655 ic->sh_invrc6 = pow(ic->rvdw,-6.0);
1662 /* Electrostatics */
1663 ic->eeltype = fr->eeltype;
1664 ic->rcoulomb = fr->rcoulomb;
1665 ic->epsilon_r = fr->epsilon_r;
1666 ic->epsfac = fr->epsfac;
1669 ic->ewaldcoeff = fr->ewaldcoeff;
1670 if (fr->coulomb_modifier==eintmodPOTSHIFT)
1672 ic->sh_ewald = gmx_erfc(ic->ewaldcoeff*ic->rcoulomb);
1679 /* Reaction-field */
1680 if (EEL_RF(ic->eeltype))
1682 ic->epsilon_rf = fr->epsilon_rf;
1683 ic->k_rf = fr->k_rf;
1684 ic->c_rf = fr->c_rf;
1688 /* For plain cut-off we might use the reaction-field kernels */
1689 ic->epsilon_rf = ic->epsilon_r;
1691 if (fr->coulomb_modifier==eintmodPOTSHIFT)
1693 ic->c_rf = 1/ic->rcoulomb;
1703 fprintf(fp,"Potential shift: LJ r^-12: %.3f r^-6 %.3f",
1704 sqr(ic->sh_invrc6),ic->sh_invrc6);
1705 if (ic->eeltype == eelCUT)
1707 fprintf(fp,", Coulomb %.3f",ic->c_rf);
1709 else if (EEL_PME(ic->eeltype))
1711 fprintf(fp,", Ewald %.3e",ic->sh_ewald);
1716 *interaction_const = ic;
1718 if (fr->nbv != NULL && fr->nbv->bUseGPU)
1720 nbnxn_cuda_init_const(fr->nbv->cu_nbv, ic, fr->nbv);
1723 bUsesSimpleTables = uses_simple_tables(fr->cutoff_scheme, fr->nbv, -1);
1724 init_interaction_const_tables(fp,ic,bUsesSimpleTables,rtab);
1727 static void init_nb_verlet(FILE *fp,
1728 nonbonded_verlet_t **nb_verlet,
1729 const t_inputrec *ir,
1730 const t_forcerec *fr,
1731 const t_commrec *cr,
1732 const char *nbpu_opt)
1734 nonbonded_verlet_t *nbv;
1737 gmx_bool bHybridGPURun = FALSE;
1739 nbnxn_alloc_t *nb_alloc;
1740 nbnxn_free_t *nb_free;
1746 nbv->ngrp = (DOMAINDECOMP(cr) ? 2 : 1);
1747 for(i=0; i<nbv->ngrp; i++)
1749 nbv->grp[i].nbl_lists.nnbl = 0;
1750 nbv->grp[i].nbat = NULL;
1751 nbv->grp[i].kernel_type = nbkNotSet;
1753 if (i == 0) /* local */
1755 pick_nbnxn_kernel(fp, cr, fr->hwinfo, fr->use_cpu_acceleration,
1757 &nbv->grp[i].kernel_type,
1758 &nbv->grp[i].ewald_excl);
1760 else /* non-local */
1762 if (nbpu_opt != NULL && strcmp(nbpu_opt,"gpu_cpu") == 0)
1764 /* Use GPU for local, select a CPU kernel for non-local */
1765 pick_nbnxn_kernel(fp, cr, fr->hwinfo, fr->use_cpu_acceleration,
1767 &nbv->grp[i].kernel_type,
1768 &nbv->grp[i].ewald_excl);
1770 bHybridGPURun = TRUE;
1774 /* Use the same kernel for local and non-local interactions */
1775 nbv->grp[i].kernel_type = nbv->grp[0].kernel_type;
1776 nbv->grp[i].ewald_excl = nbv->grp[0].ewald_excl;
1783 /* init the NxN GPU data; the last argument tells whether we'll have
1784 * both local and non-local NB calculation on GPU */
1785 nbnxn_cuda_init(fp, &nbv->cu_nbv,
1786 &fr->hwinfo->gpu_info, cr->nodeid_group_intra,
1787 (nbv->ngrp > 1) && !bHybridGPURun);
1789 if ((env = getenv("GMX_NB_MIN_CI")) != NULL)
1793 nbv->min_ci_balanced = strtol(env, &end, 10);
1794 if (!end || (*end != 0) || nbv->min_ci_balanced <= 0)
1796 gmx_fatal(FARGS, "Invalid value passed in GMX_NB_MIN_CI=%s, positive integer required", env);
1801 fprintf(debug, "Neighbor-list balancing parameter: %d (passed as env. var.)\n",
1802 nbv->min_ci_balanced);
1807 nbv->min_ci_balanced = nbnxn_cuda_min_ci_balanced(nbv->cu_nbv);
1810 fprintf(debug, "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU multi-processors)\n",
1811 nbv->min_ci_balanced);
1817 nbv->min_ci_balanced = 0;
1822 nbnxn_init_search(&nbv->nbs,
1823 DOMAINDECOMP(cr) ? & cr->dd->nc : NULL,
1824 DOMAINDECOMP(cr) ? domdec_zones(cr->dd) : NULL,
1825 gmx_omp_nthreads_get(emntNonbonded));
1827 for(i=0; i<nbv->ngrp; i++)
1829 if (nbv->grp[0].kernel_type == nbk8x8x8_CUDA)
1831 nb_alloc = &pmalloc;
1840 nbnxn_init_pairlist_set(&nbv->grp[i].nbl_lists,
1841 nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
1842 /* 8x8x8 "non-simple" lists are ATM always combined */
1843 !nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
1847 nbv->grp[0].kernel_type != nbv->grp[i].kernel_type)
1849 snew(nbv->grp[i].nbat,1);
1850 nbnxn_atomdata_init(fp,
1852 nbv->grp[i].kernel_type,
1855 nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type) ? gmx_omp_nthreads_get(emntNonbonded) : 1,
1860 nbv->grp[i].nbat = nbv->grp[0].nbat;
1865 void init_forcerec(FILE *fp,
1866 const output_env_t oenv,
1869 const t_inputrec *ir,
1870 const gmx_mtop_t *mtop,
1871 const t_commrec *cr,
1878 const char *nbpu_opt,
1879 gmx_bool bNoSolvOpt,
1882 int i,j,m,natoms,ngrp,negp_pp,negptable,egi,egj;
1888 gmx_bool bGenericKernelOnly;
1889 gmx_bool bTab,bSep14tab,bNormalnblists;
1891 int *nm_ind,egp_flags;
1893 /* By default we turn acceleration on, but it might be turned off further down... */
1894 fr->use_cpu_acceleration = TRUE;
1896 fr->bDomDec = DOMAINDECOMP(cr);
1898 natoms = mtop->natoms;
1900 if (check_box(ir->ePBC,box))
1902 gmx_fatal(FARGS,check_box(ir->ePBC,box));
1905 /* Test particle insertion ? */
1906 if (EI_TPI(ir->eI)) {
1907 /* Set to the size of the molecule to be inserted (the last one) */
1908 /* Because of old style topologies, we have to use the last cg
1909 * instead of the last molecule type.
1911 cgs = &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].cgs;
1912 fr->n_tpi = cgs->index[cgs->nr] - cgs->index[cgs->nr-1];
1913 if (fr->n_tpi != mtop->mols.index[mtop->mols.nr] - mtop->mols.index[mtop->mols.nr-1]) {
1914 gmx_fatal(FARGS,"The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
1920 /* Copy AdResS parameters */
1922 fr->adress_type = ir->adress->type;
1923 fr->adress_const_wf = ir->adress->const_wf;
1924 fr->adress_ex_width = ir->adress->ex_width;
1925 fr->adress_hy_width = ir->adress->hy_width;
1926 fr->adress_icor = ir->adress->icor;
1927 fr->adress_site = ir->adress->site;
1928 fr->adress_ex_forcecap = ir->adress->ex_forcecap;
1929 fr->adress_do_hybridpairs = ir->adress->do_hybridpairs;
1932 snew(fr->adress_group_explicit , ir->adress->n_energy_grps);
1933 for (i=0; i< ir->adress->n_energy_grps; i++){
1934 fr->adress_group_explicit[i]= ir->adress->group_explicit[i];
1937 fr->n_adress_tf_grps = ir->adress->n_tf_grps;
1938 snew(fr->adress_tf_table_index, fr->n_adress_tf_grps);
1939 for (i=0; i< fr->n_adress_tf_grps; i++){
1940 fr->adress_tf_table_index[i]= ir->adress->tf_table_index[i];
1942 copy_rvec(ir->adress->refs,fr->adress_refs);
1944 fr->adress_type = eAdressOff;
1945 fr->adress_do_hybridpairs = FALSE;
1948 /* Copy the user determined parameters */
1949 fr->userint1 = ir->userint1;
1950 fr->userint2 = ir->userint2;
1951 fr->userint3 = ir->userint3;
1952 fr->userint4 = ir->userint4;
1953 fr->userreal1 = ir->userreal1;
1954 fr->userreal2 = ir->userreal2;
1955 fr->userreal3 = ir->userreal3;
1956 fr->userreal4 = ir->userreal4;
1959 fr->fc_stepsize = ir->fc_stepsize;
1962 fr->efep = ir->efep;
1963 fr->sc_alphavdw = ir->fepvals->sc_alpha;
1964 if (ir->fepvals->bScCoul)
1966 fr->sc_alphacoul = ir->fepvals->sc_alpha;
1967 fr->sc_sigma6_min = pow(ir->fepvals->sc_sigma_min,6);
1971 fr->sc_alphacoul = 0;
1972 fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
1974 fr->sc_power = ir->fepvals->sc_power;
1975 fr->sc_r_power = ir->fepvals->sc_r_power;
1976 fr->sc_sigma6_def = pow(ir->fepvals->sc_sigma,6);
1978 env = getenv("GMX_SCSIGMA_MIN");
1982 sscanf(env,"%lf",&dbl);
1983 fr->sc_sigma6_min = pow(dbl,6);
1986 fprintf(fp,"Setting the minimum soft core sigma to %g nm\n",dbl);
1990 fr->bNonbonded = TRUE;
1991 if (getenv("GMX_NO_NONBONDED") != NULL)
1993 /* turn off non-bonded calculations */
1994 fr->bNonbonded = FALSE;
1995 md_print_warn(cr,fp,
1996 "Found environment variable GMX_NO_NONBONDED.\n"
1997 "Disabling nonbonded calculations.\n");
2000 bGenericKernelOnly = FALSE;
2002 /* We now check in the NS code whether a particular combination of interactions
2003 * can be used with water optimization, and disable it if that is not the case.
2006 if (getenv("GMX_NB_GENERIC") != NULL)
2011 "Found environment variable GMX_NB_GENERIC.\n"
2012 "Disabling all interaction-specific nonbonded kernels, will only\n"
2013 "use the slow generic ones in src/gmxlib/nonbonded/nb_generic.c\n\n");
2015 bGenericKernelOnly = TRUE;
2018 if (bGenericKernelOnly==TRUE)
2023 if( (getenv("GMX_DISABLE_CPU_ACCELERATION") != NULL) || (getenv("GMX_NOOPTIMIZEDKERNELS") != NULL) )
2025 fr->use_cpu_acceleration = FALSE;
2029 "\nFound environment variable GMX_DISABLE_CPU_ACCELERATION.\n"
2030 "Disabling all CPU architecture-specific (e.g. SSE2/SSE4/AVX) routines.\n\n");
2034 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
2036 /* Check if we can/should do all-vs-all kernels */
2037 fr->bAllvsAll = can_use_allvsall(ir,mtop,FALSE,NULL,NULL);
2038 fr->AllvsAll_work = NULL;
2039 fr->AllvsAll_workgb = NULL;
2042 /* Neighbour searching stuff */
2043 fr->cutoff_scheme = ir->cutoff_scheme;
2044 fr->bGrid = (ir->ns_type == ensGRID);
2045 fr->ePBC = ir->ePBC;
2047 /* Determine if we will do PBC for distances in bonded interactions */
2048 if (fr->ePBC == epbcNONE)
2050 fr->bMolPBC = FALSE;
2054 if (!DOMAINDECOMP(cr))
2056 /* The group cut-off scheme and SHAKE assume charge groups
2057 * are whole, but not using molpbc is faster in most cases.
2059 if (fr->cutoff_scheme == ecutsGROUP ||
2060 (ir->eConstrAlg == econtSHAKE &&
2061 (gmx_mtop_ftype_count(mtop,F_CONSTR) > 0 ||
2062 gmx_mtop_ftype_count(mtop,F_CONSTRNC) > 0)))
2064 fr->bMolPBC = ir->bPeriodicMols;
2069 if (getenv("GMX_USE_GRAPH") != NULL)
2071 fr->bMolPBC = FALSE;
2074 fprintf(fp,"\nGMX_MOLPBC is set, using the graph for bonded interactions\n\n");
2081 fr->bMolPBC = dd_bonded_molpbc(cr->dd,fr->ePBC);
2085 fr->rc_scaling = ir->refcoord_scaling;
2086 copy_rvec(ir->posres_com,fr->posres_com);
2087 copy_rvec(ir->posres_comB,fr->posres_comB);
2088 fr->rlist = cutoff_inf(ir->rlist);
2089 fr->rlistlong = cutoff_inf(ir->rlistlong);
2090 fr->eeltype = ir->coulombtype;
2091 fr->vdwtype = ir->vdwtype;
2093 fr->coulomb_modifier = ir->coulomb_modifier;
2094 fr->vdw_modifier = ir->vdw_modifier;
2096 /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
2100 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_COULOMB;
2106 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2110 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2111 fr->coulomb_modifier = eintmodEXACTCUTOFF;
2120 case eelPMEUSERSWITCH:
2121 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2126 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD;
2130 gmx_fatal(FARGS,"Unsupported electrostatic interaction: %s",eel_names[fr->eeltype]);
2134 /* Vdw: Translate from mdp settings to kernel format */
2140 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
2144 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
2151 case evdwENCADSHIFT:
2152 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2156 gmx_fatal(FARGS,"Unsupported vdw interaction: %s",evdw_names[fr->vdwtype]);
2160 /* These start out identical to ir, but might be altered if we e.g. tabulate the interaction in the kernel */
2161 fr->nbkernel_elec_modifier = fr->coulomb_modifier;
2162 fr->nbkernel_vdw_modifier = fr->vdw_modifier;
2164 fr->bTwinRange = fr->rlistlong > fr->rlist;
2165 fr->bEwald = (EEL_PME(fr->eeltype) || fr->eeltype==eelEWALD);
2167 fr->reppow = mtop->ffparams.reppow;
2169 if (ir->cutoff_scheme == ecutsGROUP)
2171 fr->bvdwtab = (fr->vdwtype != evdwCUT ||
2172 !gmx_within_tol(fr->reppow,12.0,10*GMX_DOUBLE_EPS));
2173 /* We have special kernels for standard Ewald and PME, but the pme-switch ones are tabulated above */
2174 fr->bcoultab = !(fr->eeltype == eelCUT ||
2175 fr->eeltype == eelEWALD ||
2176 fr->eeltype == eelPME ||
2177 fr->eeltype == eelRF ||
2178 fr->eeltype == eelRF_ZERO);
2180 /* If the user absolutely wants different switch/shift settings for coul/vdw, it is likely
2181 * going to be faster to tabulate the interaction than calling the generic kernel.
2183 if(fr->nbkernel_elec_modifier==eintmodPOTSWITCH && fr->nbkernel_vdw_modifier==eintmodPOTSWITCH)
2185 if((fr->rcoulomb_switch != fr->rvdw_switch) || (fr->rcoulomb != fr->rvdw))
2187 fr->bcoultab = TRUE;
2190 else if((fr->nbkernel_elec_modifier==eintmodPOTSHIFT && fr->nbkernel_vdw_modifier==eintmodPOTSHIFT) ||
2191 ((fr->nbkernel_elec_interaction == GMX_NBKERNEL_ELEC_REACTIONFIELD &&
2192 fr->nbkernel_elec_modifier==eintmodEXACTCUTOFF &&
2193 (fr->nbkernel_vdw_modifier==eintmodPOTSWITCH || fr->nbkernel_vdw_modifier==eintmodPOTSHIFT))))
2195 if(fr->rcoulomb != fr->rvdw)
2197 fr->bcoultab = TRUE;
2201 if (getenv("GMX_REQUIRE_TABLES"))
2204 fr->bcoultab = TRUE;
2209 fprintf(fp,"Table routines are used for coulomb: %s\n",bool_names[fr->bcoultab]);
2210 fprintf(fp,"Table routines are used for vdw: %s\n",bool_names[fr->bvdwtab ]);
2213 if(fr->bvdwtab==TRUE)
2215 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2216 fr->nbkernel_vdw_modifier = eintmodNONE;
2218 if(fr->bcoultab==TRUE)
2220 fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2221 fr->nbkernel_elec_modifier = eintmodNONE;
2225 if (ir->cutoff_scheme == ecutsVERLET)
2227 if (!gmx_within_tol(fr->reppow,12.0,10*GMX_DOUBLE_EPS))
2229 gmx_fatal(FARGS,"Cut-off scheme %S only supports LJ repulsion power 12",ecutscheme_names[ir->cutoff_scheme]);
2231 fr->bvdwtab = FALSE;
2232 fr->bcoultab = FALSE;
2235 /* Tables are used for direct ewald sum */
2238 if (EEL_PME(ir->coulombtype))
2241 fprintf(fp,"Will do PME sum in reciprocal space.\n");
2242 if (ir->coulombtype == eelP3M_AD)
2244 please_cite(fp,"Hockney1988");
2245 please_cite(fp,"Ballenegger2012");
2249 please_cite(fp,"Essmann95a");
2252 if (ir->ewald_geometry == eewg3DC)
2256 fprintf(fp,"Using the Ewald3DC correction for systems with a slab geometry.\n");
2258 please_cite(fp,"In-Chul99a");
2261 fr->ewaldcoeff=calc_ewaldcoeff(ir->rcoulomb, ir->ewald_rtol);
2262 init_ewald_tab(&(fr->ewald_table), cr, ir, fp);
2265 fprintf(fp,"Using a Gaussian width (1/beta) of %g nm for Ewald\n",
2270 /* Electrostatics */
2271 fr->epsilon_r = ir->epsilon_r;
2272 fr->epsilon_rf = ir->epsilon_rf;
2273 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
2274 fr->rcoulomb_switch = ir->rcoulomb_switch;
2275 fr->rcoulomb = cutoff_inf(ir->rcoulomb);
2277 /* Parameters for generalized RF */
2281 if (fr->eeltype == eelGRF)
2283 init_generalized_rf(fp,mtop,ir,fr);
2285 else if (fr->eeltype == eelSHIFT)
2287 for(m=0; (m<DIM); m++)
2288 box_size[m]=box[m][m];
2290 if ((fr->eeltype == eelSHIFT && fr->rcoulomb > fr->rcoulomb_switch))
2291 set_shift_consts(fp,fr->rcoulomb_switch,fr->rcoulomb,box_size,fr);
2294 fr->bF_NoVirSum = (EEL_FULL(fr->eeltype) ||
2295 gmx_mtop_ftype_count(mtop,F_POSRES) > 0 ||
2296 IR_ELEC_FIELD(*ir) ||
2297 (fr->adress_icor != eAdressICOff)
2300 if (fr->cutoff_scheme == ecutsGROUP &&
2301 ncg_mtop(mtop) > fr->cg_nalloc && !DOMAINDECOMP(cr)) {
2302 /* Count the total number of charge groups */
2303 fr->cg_nalloc = ncg_mtop(mtop);
2304 srenew(fr->cg_cm,fr->cg_nalloc);
2306 if (fr->shift_vec == NULL)
2307 snew(fr->shift_vec,SHIFTS);
2309 if (fr->fshift == NULL)
2310 snew(fr->fshift,SHIFTS);
2312 if (fr->nbfp == NULL) {
2313 fr->ntype = mtop->ffparams.atnr;
2314 fr->nbfp = mk_nbfp(&mtop->ffparams,fr->bBHAM);
2317 /* Copy the energy group exclusions */
2318 fr->egp_flags = ir->opts.egp_flags;
2320 /* Van der Waals stuff */
2321 fr->rvdw = cutoff_inf(ir->rvdw);
2322 fr->rvdw_switch = ir->rvdw_switch;
2323 if ((fr->vdwtype != evdwCUT) && (fr->vdwtype != evdwUSER) && !fr->bBHAM) {
2324 if (fr->rvdw_switch >= fr->rvdw)
2325 gmx_fatal(FARGS,"rvdw_switch (%f) must be < rvdw (%f)",
2326 fr->rvdw_switch,fr->rvdw);
2328 fprintf(fp,"Using %s Lennard-Jones, switch between %g and %g nm\n",
2329 (fr->eeltype==eelSWITCH) ? "switched":"shifted",
2330 fr->rvdw_switch,fr->rvdw);
2333 if (fr->bBHAM && (fr->vdwtype == evdwSHIFT || fr->vdwtype == evdwSWITCH))
2334 gmx_fatal(FARGS,"Switch/shift interaction not supported with Buckingham");
2337 fprintf(fp,"Cut-off's: NS: %g Coulomb: %g %s: %g\n",
2338 fr->rlist,fr->rcoulomb,fr->bBHAM ? "BHAM":"LJ",fr->rvdw);
2340 fr->eDispCorr = ir->eDispCorr;
2341 if (ir->eDispCorr != edispcNO)
2343 set_avcsixtwelve(fp,fr,mtop);
2348 set_bham_b_max(fp,fr,mtop);
2351 fr->bGB = (ir->implicit_solvent == eisGBSA);
2352 fr->gb_epsilon_solvent = ir->gb_epsilon_solvent;
2354 /* Copy the GBSA data (radius, volume and surftens for each
2355 * atomtype) from the topology atomtype section to forcerec.
2357 snew(fr->atype_radius,fr->ntype);
2358 snew(fr->atype_vol,fr->ntype);
2359 snew(fr->atype_surftens,fr->ntype);
2360 snew(fr->atype_gb_radius,fr->ntype);
2361 snew(fr->atype_S_hct,fr->ntype);
2363 if (mtop->atomtypes.nr > 0)
2365 for(i=0;i<fr->ntype;i++)
2366 fr->atype_radius[i] =mtop->atomtypes.radius[i];
2367 for(i=0;i<fr->ntype;i++)
2368 fr->atype_vol[i] = mtop->atomtypes.vol[i];
2369 for(i=0;i<fr->ntype;i++)
2370 fr->atype_surftens[i] = mtop->atomtypes.surftens[i];
2371 for(i=0;i<fr->ntype;i++)
2372 fr->atype_gb_radius[i] = mtop->atomtypes.gb_radius[i];
2373 for(i=0;i<fr->ntype;i++)
2374 fr->atype_S_hct[i] = mtop->atomtypes.S_hct[i];
2377 /* Generate the GB table if needed */
2381 fr->gbtabscale=2000;
2387 fr->gbtab=make_gb_table(fp,oenv,fr,tabpfn,fr->gbtabscale);
2389 init_gb(&fr->born,cr,fr,ir,mtop,ir->rgbradii,ir->gb_algorithm);
2391 /* Copy local gb data (for dd, this is done in dd_partition_system) */
2392 if (!DOMAINDECOMP(cr))
2394 make_local_gb(cr,fr->born,ir->gb_algorithm);
2398 /* Set the charge scaling */
2399 if (fr->epsilon_r != 0)
2400 fr->epsfac = ONE_4PI_EPS0/fr->epsilon_r;
2402 /* eps = 0 is infinite dieletric: no coulomb interactions */
2405 /* Reaction field constants */
2406 if (EEL_RF(fr->eeltype))
2407 calc_rffac(fp,fr->eeltype,fr->epsilon_r,fr->epsilon_rf,
2408 fr->rcoulomb,fr->temp,fr->zsquare,box,
2409 &fr->kappa,&fr->k_rf,&fr->c_rf);
2411 set_chargesum(fp,fr,mtop);
2413 /* if we are using LR electrostatics, and they are tabulated,
2414 * the tables will contain modified coulomb interactions.
2415 * Since we want to use the non-shifted ones for 1-4
2416 * coulombic interactions, we must have an extra set of tables.
2419 /* Construct tables.
2420 * A little unnecessary to make both vdw and coul tables sometimes,
2421 * but what the heck... */
2423 bTab = fr->bcoultab || fr->bvdwtab || fr->bEwald;
2425 bSep14tab = ((!bTab || fr->eeltype!=eelCUT || fr->vdwtype!=evdwCUT ||
2426 fr->bBHAM || fr->bEwald) &&
2427 (gmx_mtop_ftype_count(mtop,F_LJ14) > 0 ||
2428 gmx_mtop_ftype_count(mtop,F_LJC14_Q) > 0 ||
2429 gmx_mtop_ftype_count(mtop,F_LJC_PAIRS_NB) > 0));
2431 negp_pp = ir->opts.ngener - ir->nwall;
2434 bNormalnblists = TRUE;
2437 bNormalnblists = (ir->eDispCorr != edispcNO);
2438 for(egi=0; egi<negp_pp; egi++) {
2439 for(egj=egi; egj<negp_pp; egj++) {
2440 egp_flags = ir->opts.egp_flags[GID(egi,egj,ir->opts.ngener)];
2441 if (!(egp_flags & EGP_EXCL)) {
2442 if (egp_flags & EGP_TABLE) {
2445 bNormalnblists = TRUE;
2450 if (bNormalnblists) {
2451 fr->nnblists = negptable + 1;
2453 fr->nnblists = negptable;
2455 if (fr->nnblists > 1)
2456 snew(fr->gid2nblists,ir->opts.ngener*ir->opts.ngener);
2458 snew(fr->nblists,fr->nnblists);
2460 /* This code automatically gives table length tabext without cut-off's,
2461 * in that case grompp should already have checked that we do not need
2462 * normal tables and we only generate tables for 1-4 interactions.
2464 rtab = ir->rlistlong + ir->tabext;
2467 /* make tables for ordinary interactions */
2468 if (bNormalnblists) {
2469 make_nbf_tables(fp,oenv,fr,rtab,cr,tabfn,NULL,NULL,&fr->nblists[0]);
2471 fr->tab14 = fr->nblists[0].table_elec_vdw;
2476 if (negptable > 0) {
2477 /* Read the special tables for certain energy group pairs */
2478 nm_ind = mtop->groups.grps[egcENER].nm_ind;
2479 for(egi=0; egi<negp_pp; egi++) {
2480 for(egj=egi; egj<negp_pp; egj++) {
2481 egp_flags = ir->opts.egp_flags[GID(egi,egj,ir->opts.ngener)];
2482 if ((egp_flags & EGP_TABLE) && !(egp_flags & EGP_EXCL)) {
2483 nbl = &(fr->nblists[m]);
2484 if (fr->nnblists > 1) {
2485 fr->gid2nblists[GID(egi,egj,ir->opts.ngener)] = m;
2487 /* Read the table file with the two energy groups names appended */
2488 make_nbf_tables(fp,oenv,fr,rtab,cr,tabfn,
2489 *mtop->groups.grpname[nm_ind[egi]],
2490 *mtop->groups.grpname[nm_ind[egj]],
2493 } else if (fr->nnblists > 1) {
2494 fr->gid2nblists[GID(egi,egj,ir->opts.ngener)] = 0;
2502 /* generate extra tables with plain Coulomb for 1-4 interactions only */
2503 fr->tab14 = make_tables(fp,oenv,fr,MASTER(cr),tabpfn,rtab,
2504 GMX_MAKETABLES_14ONLY);
2507 /* Read AdResS Thermo Force table if needed */
2508 if(fr->adress_icor == eAdressICThermoForce)
2510 /* old todo replace */
2512 if (ir->adress->n_tf_grps > 0){
2513 make_adress_tf_tables(fp,oenv,fr,ir,tabfn, mtop, box);
2516 /* load the default table */
2517 snew(fr->atf_tabs, 1);
2518 fr->atf_tabs[DEFAULT_TF_TABLE] = make_atf_table(fp,oenv,fr,tabafn, box);
2523 fr->nwall = ir->nwall;
2524 if (ir->nwall && ir->wall_type==ewtTABLE)
2526 make_wall_tables(fp,oenv,ir,tabfn,&mtop->groups,fr);
2529 if (fcd && tabbfn) {
2530 fcd->bondtab = make_bonded_tables(fp,
2531 F_TABBONDS,F_TABBONDSNC,
2533 fcd->angletab = make_bonded_tables(fp,
2536 fcd->dihtab = make_bonded_tables(fp,
2541 fprintf(debug,"No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
2544 /* QM/MM initialization if requested
2548 fprintf(stderr,"QM/MM calculation requested.\n");
2551 fr->bQMMM = ir->bQMMM;
2552 fr->qr = mk_QMMMrec();
2554 /* Set all the static charge group info */
2555 fr->cginfo_mb = init_cginfo_mb(fp,mtop,fr,bNoSolvOpt,
2556 &fr->bExcl_IntraCGAll_InterCGNone);
2557 if (DOMAINDECOMP(cr)) {
2560 fr->cginfo = cginfo_expand(mtop->nmolblock,fr->cginfo_mb);
2563 if (!DOMAINDECOMP(cr))
2565 /* When using particle decomposition, the effect of the second argument,
2566 * which sets fr->hcg, is corrected later in do_md and init_em.
2568 forcerec_set_ranges(fr,ncg_mtop(mtop),ncg_mtop(mtop),
2569 mtop->natoms,mtop->natoms,mtop->natoms);
2572 fr->print_force = print_force;
2575 /* coarse load balancing vars */
2580 /* Initialize neighbor search */
2581 init_ns(fp,cr,&fr->ns,fr,mtop,box);
2583 if (cr->duty & DUTY_PP)
2585 gmx_nonbonded_setup(fp,fr,bGenericKernelOnly);
2589 gmx_setup_adress_kernels(fp,bGenericKernelOnly);
2594 /* Initialize the thread working data for bonded interactions */
2595 init_forcerec_f_threads(fr,mtop->groups.grps[egcENER].nr);
2597 snew(fr->excl_load,fr->nthreads+1);
2599 if (fr->cutoff_scheme == ecutsVERLET)
2601 if (ir->rcoulomb != ir->rvdw)
2603 gmx_fatal(FARGS,"With Verlet lists rcoulomb and rvdw should be identical");
2606 init_nb_verlet(fp, &fr->nbv, ir, fr, cr, nbpu_opt);
2609 /* fr->ic is used both by verlet and group kernels (to some extent) now */
2610 init_interaction_const(fp, &fr->ic, fr, rtab);
2611 if (ir->eDispCorr != edispcNO)
2613 calc_enervirdiff(fp,ir->eDispCorr,fr);
2617 #define pr_real(fp,r) fprintf(fp,"%s: %e\n",#r,r)
2618 #define pr_int(fp,i) fprintf((fp),"%s: %d\n",#i,i)
2619 #define pr_bool(fp,b) fprintf((fp),"%s: %s\n",#b,bool_names[b])
2621 void pr_forcerec(FILE *fp,t_forcerec *fr,t_commrec *cr)
2625 pr_real(fp,fr->rlist);
2626 pr_real(fp,fr->rcoulomb);
2627 pr_real(fp,fr->fudgeQQ);
2628 pr_bool(fp,fr->bGrid);
2629 pr_bool(fp,fr->bTwinRange);
2630 /*pr_int(fp,fr->cg0);
2631 pr_int(fp,fr->hcg);*/
2632 for(i=0; i<fr->nnblists; i++)
2633 pr_int(fp,fr->nblists[i].table_elec_vdw.n);
2634 pr_real(fp,fr->rcoulomb_switch);
2635 pr_real(fp,fr->rcoulomb);
2640 void forcerec_set_excl_load(t_forcerec *fr,
2641 const gmx_localtop_t *top,const t_commrec *cr)
2644 int t,i,j,ntot,n,ntarget;
2646 if (cr != NULL && PARTDECOMP(cr))
2648 /* No OpenMP with particle decomposition */
2656 ind = top->excls.index;
2660 for(i=0; i<top->excls.nr; i++)
2662 for(j=ind[i]; j<ind[i+1]; j++)
2671 fr->excl_load[0] = 0;
2674 for(t=1; t<=fr->nthreads; t++)
2676 ntarget = (ntot*t)/fr->nthreads;
2677 while(i < top->excls.nr && n < ntarget)
2679 for(j=ind[i]; j<ind[i+1]; j++)
2688 fr->excl_load[t] = i;