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
49 #include "nonbonded.h"
63 #include "mtop_util.h"
65 t_forcerec *mk_forcerec(void)
75 static void pr_nbfp(FILE *fp,real *nbfp,gmx_bool bBHAM,int atnr)
79 for(i=0; (i<atnr); i++) {
80 for(j=0; (j<atnr); j++) {
81 fprintf(fp,"%2d - %2d",i,j);
83 fprintf(fp," a=%10g, b=%10g, c=%10g\n",BHAMA(nbfp,atnr,i,j),
84 BHAMB(nbfp,atnr,i,j),BHAMC(nbfp,atnr,i,j));
86 fprintf(fp," c6=%10g, c12=%10g\n",C6(nbfp,atnr,i,j),
93 static real *mk_nbfp(const gmx_ffparams_t *idef,gmx_bool bBHAM)
100 snew(nbfp,3*atnr*atnr);
101 for(i=k=0; (i<atnr); i++) {
102 for(j=0; (j<atnr); j++,k++) {
103 BHAMA(nbfp,atnr,i,j) = idef->iparams[k].bham.a;
104 BHAMB(nbfp,atnr,i,j) = idef->iparams[k].bham.b;
105 BHAMC(nbfp,atnr,i,j) = idef->iparams[k].bham.c;
110 snew(nbfp,2*atnr*atnr);
111 for(i=k=0; (i<atnr); i++) {
112 for(j=0; (j<atnr); j++,k++) {
113 C6(nbfp,atnr,i,j) = idef->iparams[k].lj.c6;
114 C12(nbfp,atnr,i,j) = idef->iparams[k].lj.c12;
121 /* This routine sets fr->solvent_opt to the most common solvent in the
122 * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
123 * the fr->solvent_type array with the correct type (or esolNO).
125 * Charge groups that fulfill the conditions but are not identical to the
126 * most common one will be marked as esolNO in the solvent_type array.
128 * TIP3p is identical to SPC for these purposes, so we call it
129 * SPC in the arrays (Apologies to Bill Jorgensen ;-)
131 * NOTE: QM particle should not
132 * become an optimized solvent. Not even if there is only one charge
142 } solvent_parameters_t;
145 check_solvent_cg(const gmx_moltype_t *molt,
148 const unsigned char *qm_grpnr,
149 const t_grps *qm_grps,
151 int *n_solvent_parameters,
152 solvent_parameters_t **solvent_parameters_p,
156 const t_blocka * excl;
167 solvent_parameters_t *solvent_parameters;
169 /* We use a list with parameters for each solvent type.
170 * Every time we discover a new molecule that fulfills the basic
171 * conditions for a solvent we compare with the previous entries
172 * in these lists. If the parameters are the same we just increment
173 * the counter for that type, and otherwise we create a new type
174 * based on the current molecule.
176 * Once we've finished going through all molecules we check which
177 * solvent is most common, and mark all those molecules while we
178 * clear the flag on all others.
181 solvent_parameters = *solvent_parameters_p;
183 /* Mark the cg first as non optimized */
186 /* Check if this cg has no exclusions with atoms in other charge groups
187 * and all atoms inside the charge group excluded.
188 * We only have 3 or 4 atom solvent loops.
190 if (GET_CGINFO_EXCL_INTER(cginfo) ||
191 !GET_CGINFO_EXCL_INTRA(cginfo))
196 /* Get the indices of the first atom in this charge group */
197 j0 = molt->cgs.index[cg0];
198 j1 = molt->cgs.index[cg0+1];
200 /* Number of atoms in our molecule */
205 "Moltype '%s': there are %d atoms in this charge group\n",
209 /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
217 /* Check if we are doing QM on this group */
219 if (qm_grpnr != NULL)
221 for(j=j0 ; j<j1 && !qm; j++)
223 qm = (qm_grpnr[j] < qm_grps->nr - 1);
226 /* Cannot use solvent optimization with QM */
232 atom = molt->atoms.atom;
234 /* Still looks like a solvent, time to check parameters */
236 /* If it is perturbed (free energy) we can't use the solvent loops,
237 * so then we just skip to the next molecule.
241 for(j=j0; j<j1 && !perturbed; j++)
243 perturbed = PERTURBED(atom[j]);
251 /* Now it's only a question if the VdW and charge parameters
252 * are OK. Before doing the check we compare and see if they are
253 * identical to a possible previous solvent type.
254 * First we assign the current types and charges.
258 tmp_vdwtype[j] = atom[j0+j].type;
259 tmp_charge[j] = atom[j0+j].q;
262 /* Does it match any previous solvent type? */
263 for(k=0 ; k<*n_solvent_parameters; k++)
268 /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
269 if( (solvent_parameters[k].model==esolSPC && nj!=3) ||
270 (solvent_parameters[k].model==esolTIP4P && nj!=4) )
273 /* Check that types & charges match for all atoms in molecule */
274 for(j=0 ; j<nj && match==TRUE; j++)
276 if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
280 if(tmp_charge[j] != solvent_parameters[k].charge[j])
287 /* Congratulations! We have a matched solvent.
288 * Flag it with this type for later processing.
291 solvent_parameters[k].count += nmol;
293 /* We are done with this charge group */
298 /* If we get here, we have a tentative new solvent type.
299 * Before we add it we must check that it fulfills the requirements
300 * of the solvent optimized loops. First determine which atoms have
306 tjA = tmp_vdwtype[j];
308 /* Go through all other tpes and see if any have non-zero
309 * VdW parameters when combined with this one.
311 for(k=0; k<fr->ntype && (has_vdw[j]==FALSE); k++)
313 /* We already checked that the atoms weren't perturbed,
314 * so we only need to check state A now.
318 has_vdw[j] = (has_vdw[j] ||
319 (BHAMA(fr->nbfp,fr->ntype,tjA,k) != 0.0) ||
320 (BHAMB(fr->nbfp,fr->ntype,tjA,k) != 0.0) ||
321 (BHAMC(fr->nbfp,fr->ntype,tjA,k) != 0.0));
326 has_vdw[j] = (has_vdw[j] ||
327 (C6(fr->nbfp,fr->ntype,tjA,k) != 0.0) ||
328 (C12(fr->nbfp,fr->ntype,tjA,k) != 0.0));
333 /* Now we know all we need to make the final check and assignment. */
337 * For this we require thatn all atoms have charge,
338 * the charges on atom 2 & 3 should be the same, and only
339 * atom 1 should have VdW.
341 if (has_vdw[0] == TRUE &&
342 has_vdw[1] == FALSE &&
343 has_vdw[2] == FALSE &&
344 tmp_charge[0] != 0 &&
345 tmp_charge[1] != 0 &&
346 tmp_charge[2] == tmp_charge[1])
348 srenew(solvent_parameters,*n_solvent_parameters+1);
349 solvent_parameters[*n_solvent_parameters].model = esolSPC;
350 solvent_parameters[*n_solvent_parameters].count = nmol;
353 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
354 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
357 *cg_sp = *n_solvent_parameters;
358 (*n_solvent_parameters)++;
363 /* Or could it be a TIP4P?
364 * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
365 * Only atom 1 should have VdW.
367 if(has_vdw[0] == TRUE &&
368 has_vdw[1] == FALSE &&
369 has_vdw[2] == FALSE &&
370 has_vdw[3] == FALSE &&
371 tmp_charge[0] == 0 &&
372 tmp_charge[1] != 0 &&
373 tmp_charge[2] == tmp_charge[1] &&
376 srenew(solvent_parameters,*n_solvent_parameters+1);
377 solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
378 solvent_parameters[*n_solvent_parameters].count = nmol;
381 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
382 solvent_parameters[*n_solvent_parameters].charge[k] = tmp_charge[k];
385 *cg_sp = *n_solvent_parameters;
386 (*n_solvent_parameters)++;
390 *solvent_parameters_p = solvent_parameters;
394 check_solvent(FILE * fp,
395 const gmx_mtop_t * mtop,
397 cginfo_mb_t *cginfo_mb)
400 const t_block * mols;
401 const gmx_moltype_t *molt;
402 int mb,mol,cg_mol,at_offset,cg_offset,am,cgm,i,nmol_ch,nmol;
403 int n_solvent_parameters;
404 solvent_parameters_t *solvent_parameters;
410 fprintf(debug,"Going to determine what solvent types we have.\n");
415 n_solvent_parameters = 0;
416 solvent_parameters = NULL;
417 /* Allocate temporary array for solvent type */
418 snew(cg_sp,mtop->nmolblock);
422 for(mb=0; mb<mtop->nmolblock; mb++)
424 molt = &mtop->moltype[mtop->molblock[mb].type];
426 /* Here we have to loop over all individual molecules
427 * because we need to check for QMMM particles.
429 snew(cg_sp[mb],cginfo_mb[mb].cg_mod);
430 nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
431 nmol = mtop->molblock[mb].nmol/nmol_ch;
432 for(mol=0; mol<nmol_ch; mol++)
435 am = mol*cgs->index[cgs->nr];
436 for(cg_mol=0; cg_mol<cgs->nr; cg_mol++)
438 check_solvent_cg(molt,cg_mol,nmol,
439 mtop->groups.grpnr[egcQMMM] ?
440 mtop->groups.grpnr[egcQMMM]+at_offset+am : 0,
441 &mtop->groups.grps[egcQMMM],
443 &n_solvent_parameters,&solvent_parameters,
444 cginfo_mb[mb].cginfo[cgm+cg_mol],
445 &cg_sp[mb][cgm+cg_mol]);
448 cg_offset += cgs->nr;
449 at_offset += cgs->index[cgs->nr];
452 /* Puh! We finished going through all charge groups.
453 * Now find the most common solvent model.
456 /* Most common solvent this far */
458 for(i=0;i<n_solvent_parameters;i++)
461 solvent_parameters[i].count > solvent_parameters[bestsp].count)
469 bestsol = solvent_parameters[bestsp].model;
476 #ifdef DISABLE_WATER_NLIST
481 for(mb=0; mb<mtop->nmolblock; mb++)
483 cgs = &mtop->moltype[mtop->molblock[mb].type].cgs;
484 nmol = (mtop->molblock[mb].nmol*cgs->nr)/cginfo_mb[mb].cg_mod;
485 for(i=0; i<cginfo_mb[mb].cg_mod; i++)
487 if (cg_sp[mb][i] == bestsp)
489 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i],bestsol);
494 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i],esolNO);
501 if (bestsol != esolNO && fp!=NULL)
503 fprintf(fp,"\nEnabling %s-like water optimization for %d molecules.\n\n",
505 solvent_parameters[bestsp].count);
508 sfree(solvent_parameters);
509 fr->solvent_opt = bestsol;
512 static cginfo_mb_t *init_cginfo_mb(FILE *fplog,const gmx_mtop_t *mtop,
513 t_forcerec *fr,gmx_bool bNoSolvOpt,
514 gmx_bool *bExcl_IntraCGAll_InterCGNone)
517 const t_blocka *excl;
518 const gmx_moltype_t *molt;
519 const gmx_molblock_t *molb;
520 cginfo_mb_t *cginfo_mb;
522 int cg_offset,a_offset,cgm,am;
523 int mb,m,ncg_tot,cg,a0,a1,gid,ai,j,aj,excl_nalloc;
524 gmx_bool bId,*bExcl,bExclIntraAll,bExclInter;
526 ncg_tot = ncg_mtop(mtop);
527 snew(cginfo_mb,mtop->nmolblock);
529 *bExcl_IntraCGAll_InterCGNone = TRUE;
532 snew(bExcl,excl_nalloc);
535 for(mb=0; mb<mtop->nmolblock; mb++)
537 molb = &mtop->molblock[mb];
538 molt = &mtop->moltype[molb->type];
542 /* Check if the cginfo is identical for all molecules in this block.
543 * If so, we only need an array of the size of one molecule.
544 * Otherwise we make an array of #mol times #cgs per molecule.
548 for(m=0; m<molb->nmol; m++)
550 am = m*cgs->index[cgs->nr];
551 for(cg=0; cg<cgs->nr; cg++)
554 a1 = cgs->index[cg+1];
555 if (ggrpnr(&mtop->groups,egcENER,a_offset+am+a0) !=
556 ggrpnr(&mtop->groups,egcENER,a_offset +a0))
560 if (mtop->groups.grpnr[egcQMMM] != NULL)
562 for(ai=a0; ai<a1; ai++)
564 if (mtop->groups.grpnr[egcQMMM][a_offset+am+ai] !=
565 mtop->groups.grpnr[egcQMMM][a_offset +ai])
574 cginfo_mb[mb].cg_start = cg_offset;
575 cginfo_mb[mb].cg_end = cg_offset + molb->nmol*cgs->nr;
576 cginfo_mb[mb].cg_mod = (bId ? 1 : molb->nmol)*cgs->nr;
577 snew(cginfo_mb[mb].cginfo,cginfo_mb[mb].cg_mod);
578 cginfo = cginfo_mb[mb].cginfo;
580 for(m=0; m<(bId ? 1 : molb->nmol); m++)
583 am = m*cgs->index[cgs->nr];
584 for(cg=0; cg<cgs->nr; cg++)
587 a1 = cgs->index[cg+1];
589 /* Store the energy group in cginfo */
590 gid = ggrpnr(&mtop->groups,egcENER,a_offset+am+a0);
591 SET_CGINFO_GID(cginfo[cgm+cg],gid);
593 /* Check the intra/inter charge group exclusions */
594 if (a1-a0 > excl_nalloc) {
595 excl_nalloc = a1 - a0;
596 srenew(bExcl,excl_nalloc);
598 /* bExclIntraAll: all intra cg interactions excluded
599 * bExclInter: any inter cg interactions excluded
601 bExclIntraAll = TRUE;
603 for(ai=a0; ai<a1; ai++) {
604 /* Clear the exclusion list for atom ai */
605 for(aj=a0; aj<a1; aj++) {
606 bExcl[aj-a0] = FALSE;
608 /* Loop over all the exclusions of atom ai */
609 for(j=excl->index[ai]; j<excl->index[ai+1]; j++)
612 if (aj < a0 || aj >= a1)
621 /* Check if ai excludes a0 to a1 */
622 for(aj=a0; aj<a1; aj++)
626 bExclIntraAll = FALSE;
632 SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
636 SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
638 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
640 /* The size in cginfo is currently only read with DD */
641 gmx_fatal(FARGS,"A charge group has size %d which is larger than the limit of %d atoms",a1-a0,MAX_CHARGEGROUP_SIZE);
643 SET_CGINFO_NATOMS(cginfo[cgm+cg],a1-a0);
645 if (!bExclIntraAll || bExclInter)
647 *bExcl_IntraCGAll_InterCGNone = FALSE;
651 cg_offset += molb->nmol*cgs->nr;
652 a_offset += molb->nmol*cgs->index[cgs->nr];
656 /* the solvent optimizer is called after the QM is initialized,
657 * because we don't want to have the QM subsystemto become an
661 check_solvent(fplog,mtop,fr,cginfo_mb);
663 if (getenv("GMX_NO_SOLV_OPT"))
667 fprintf(fplog,"Found environment variable GMX_NO_SOLV_OPT.\n"
668 "Disabling all solvent optimization\n");
670 fr->solvent_opt = esolNO;
674 fr->solvent_opt = esolNO;
676 if (!fr->solvent_opt)
678 for(mb=0; mb<mtop->nmolblock; mb++)
680 for(cg=0; cg<cginfo_mb[mb].cg_mod; cg++)
682 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg],esolNO);
690 static int *cginfo_expand(int nmb,cginfo_mb_t *cgi_mb)
695 ncg = cgi_mb[nmb-1].cg_end;
698 for(cg=0; cg<ncg; cg++)
700 while (cg >= cgi_mb[mb].cg_end)
705 cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
711 static void set_chargesum(FILE *log,t_forcerec *fr,const gmx_mtop_t *mtop)
715 const t_atoms *atoms;
718 for(mb=0; mb<mtop->nmolblock; mb++)
720 nmol = mtop->molblock[mb].nmol;
721 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
722 for(i=0; i<atoms->nr; i++)
724 qsum += nmol*atoms->atom[i].q;
728 if (fr->efep != efepNO)
731 for(mb=0; mb<mtop->nmolblock; mb++)
733 nmol = mtop->molblock[mb].nmol;
734 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
735 for(i=0; i<atoms->nr; i++)
737 qsum += nmol*atoms->atom[i].qB;
744 fr->qsum[1] = fr->qsum[0];
747 if (fr->efep == efepNO)
748 fprintf(log,"System total charge: %.3f\n",fr->qsum[0]);
750 fprintf(log,"System total charge, top. A: %.3f top. B: %.3f\n",
751 fr->qsum[0],fr->qsum[1]);
755 void update_forcerec(FILE *log,t_forcerec *fr,matrix box)
757 if (fr->eeltype == eelGRF)
759 calc_rffac(NULL,fr->eeltype,fr->epsilon_r,fr->epsilon_rf,
760 fr->rcoulomb,fr->temp,fr->zsquare,box,
761 &fr->kappa,&fr->k_rf,&fr->c_rf);
765 void set_avcsixtwelve(FILE *fplog,t_forcerec *fr,const gmx_mtop_t *mtop)
767 const t_atoms *atoms;
768 const t_blocka *excl;
769 int mb,nmol,nmolc,i,j,tpi,tpj,j1,j2,k,n,nexcl,q;
770 #if (defined SIZEOF_LONG_LONG_INT) && (SIZEOF_LONG_LONG_INT >= 8)
771 long long int npair,npair_ij,tmpi,tmpj;
773 double npair, npair_ij,tmpi,tmpj;
784 for(q=0; q<(fr->efep==efepNO ? 1 : 2); q++) {
790 /* Count the types so we avoid natoms^2 operations */
792 for(mb=0; mb<mtop->nmolblock; mb++) {
793 nmol = mtop->molblock[mb].nmol;
794 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
795 for(i=0; i<atoms->nr; i++) {
798 tpi = atoms->atom[i].type;
802 tpi = atoms->atom[i].typeB;
804 typecount[tpi] += nmol;
807 for(tpi=0; tpi<ntp; tpi++) {
808 for(tpj=tpi; tpj<ntp; tpj++) {
809 tmpi = typecount[tpi];
810 tmpj = typecount[tpj];
813 npair_ij = tmpi*tmpj;
817 npair_ij = tmpi*(tmpi - 1)/2;
820 csix += npair_ij*BHAMC(nbfp,ntp,tpi,tpj);
822 csix += npair_ij* C6(nbfp,ntp,tpi,tpj);
823 ctwelve += npair_ij* C12(nbfp,ntp,tpi,tpj);
829 /* Subtract the excluded pairs.
830 * The main reason for substracting exclusions is that in some cases
831 * some combinations might never occur and the parameters could have
832 * any value. These unused values should not influence the dispersion
835 for(mb=0; mb<mtop->nmolblock; mb++) {
836 nmol = mtop->molblock[mb].nmol;
837 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
838 excl = &mtop->moltype[mtop->molblock[mb].type].excls;
839 for(i=0; (i<atoms->nr); i++) {
842 tpi = atoms->atom[i].type;
846 tpi = atoms->atom[i].typeB;
849 j2 = excl->index[i+1];
850 for(j=j1; j<j2; j++) {
856 tpj = atoms->atom[k].type;
860 tpj = atoms->atom[k].typeB;
863 csix -= nmol*BHAMC(nbfp,ntp,tpi,tpj);
865 csix -= nmol*C6 (nbfp,ntp,tpi,tpj);
866 ctwelve -= nmol*C12(nbfp,ntp,tpi,tpj);
874 /* Only correct for the interaction of the test particle
875 * with the rest of the system.
877 atoms = &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].atoms;
880 tpi = atoms->atom[atoms->nr-1].type;
884 tpi = atoms->atom[atoms->nr-1].typeB;
887 for(mb=0; mb<mtop->nmolblock; mb++) {
888 nmol = mtop->molblock[mb].nmol;
889 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
890 for(j=0; j<atoms->nr; j++) {
892 /* Remove the interaction of the test charge group
895 if (mb == mtop->nmolblock-1 && j >= atoms->nr - fr->n_tpi)
901 tpj = atoms->atom[j].type;
905 tpj = atoms->atom[j].typeB;
909 csix += nmolc*BHAMC(nbfp,ntp,tpi,tpj);
913 csix += nmolc*C6 (nbfp,ntp,tpi,tpj);
914 ctwelve += nmolc*C12(nbfp,ntp,tpi,tpj);
920 if (npair - nexcl <= 0 && fplog) {
921 fprintf(fplog,"\nWARNING: There are no atom pairs for dispersion correction\n\n");
925 csix /= npair - nexcl;
926 ctwelve /= npair - nexcl;
929 fprintf(debug,"Counted %d exclusions\n",nexcl);
930 fprintf(debug,"Average C6 parameter is: %10g\n",(double)csix);
931 fprintf(debug,"Average C12 parameter is: %10g\n",(double)ctwelve);
933 fr->avcsix[q] = csix;
934 fr->avctwelve[q] = ctwelve;
938 if (fr->eDispCorr == edispcAllEner ||
939 fr->eDispCorr == edispcAllEnerPres)
941 fprintf(fplog,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
942 fr->avcsix[0],fr->avctwelve[0]);
946 fprintf(fplog,"Long Range LJ corr.: <C6> %10.4e\n",fr->avcsix[0]);
952 static void set_bham_b_max(FILE *fplog,t_forcerec *fr,
953 const gmx_mtop_t *mtop)
955 const t_atoms *at1,*at2;
956 int mt1,mt2,i,j,tpi,tpj,ntypes;
962 fprintf(fplog,"Determining largest Buckingham b parameter for table\n");
969 for(mt1=0; mt1<mtop->nmoltype; mt1++)
971 at1 = &mtop->moltype[mt1].atoms;
972 for(i=0; (i<at1->nr); i++)
974 tpi = at1->atom[i].type;
976 gmx_fatal(FARGS,"Atomtype[%d] = %d, maximum = %d",i,tpi,ntypes);
978 for(mt2=mt1; mt2<mtop->nmoltype; mt2++)
980 at2 = &mtop->moltype[mt2].atoms;
981 for(j=0; (j<at2->nr); j++) {
982 tpj = at2->atom[j].type;
985 gmx_fatal(FARGS,"Atomtype[%d] = %d, maximum = %d",j,tpj,ntypes);
987 b = BHAMB(nbfp,ntypes,tpi,tpj);
988 if (b > fr->bham_b_max)
992 if ((b < bmin) || (bmin==-1))
1002 fprintf(fplog,"Buckingham b parameters, min: %g, max: %g\n",
1003 bmin,fr->bham_b_max);
1007 static void make_nbf_tables(FILE *fp,const output_env_t oenv,
1008 t_forcerec *fr,real rtab,
1009 const t_commrec *cr,
1010 const char *tabfn,char *eg1,char *eg2,
1017 if (tabfn == NULL) {
1019 fprintf(debug,"No table file name passed, can not read table, can not do non-bonded interactions\n");
1023 sprintf(buf,"%s",tabfn);
1025 /* Append the two energy group names */
1026 sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1,"_%s_%s.%s",
1027 eg1,eg2,ftp2ext(efXVG));
1028 nbl->tab = make_tables(fp,oenv,fr,MASTER(cr),buf,rtab,0);
1029 /* Copy the contents of the table to separate coulomb and LJ tables too,
1030 * to improve cache performance.
1033 /* For performance reasons we want
1034 * the table data to be aligned to 16-byte. This is accomplished
1035 * by allocating 16 bytes extra to a temporary pointer, and then
1036 * calculating an aligned pointer. This new pointer must not be
1037 * used in a free() call, but thankfully we're sloppy enough not
1041 /* 8 fp entries per vdw table point, n+1 points, and 16 bytes extra to align it. */
1042 p_tmp = malloc(8*(nbl->tab.n+1)*sizeof(real)+16);
1044 /* align it - size_t has the same same as a pointer */
1045 nbl->vdwtab = (real *) (((size_t) p_tmp + 16) & (~((size_t) 15)));
1047 /* 4 fp entries per coul table point, n+1 points, and 16 bytes extra to align it. */
1048 p_tmp = malloc(4*(nbl->tab.n+1)*sizeof(real)+16);
1050 /* align it - size_t has the same same as a pointer */
1051 nbl->coultab = (real *) (((size_t) p_tmp + 16) & (~((size_t) 15)));
1054 for(i=0; i<=nbl->tab.n; i++) {
1056 nbl->coultab[4*i+j] = nbl->tab.tab[12*i+j];
1058 nbl->vdwtab [8*i+j] = nbl->tab.tab[12*i+4+j];
1062 static void count_tables(int ftype1,int ftype2,const gmx_mtop_t *mtop,
1063 int *ncount,int **count)
1065 const gmx_moltype_t *molt;
1067 int mt,ftype,stride,i,j,tabnr;
1069 for(mt=0; mt<mtop->nmoltype; mt++)
1071 molt = &mtop->moltype[mt];
1072 for(ftype=0; ftype<F_NRE; ftype++)
1074 if (ftype == ftype1 || ftype == ftype2) {
1075 il = &molt->ilist[ftype];
1076 stride = 1 + NRAL(ftype);
1077 for(i=0; i<il->nr; i+=stride) {
1078 tabnr = mtop->ffparams.iparams[il->iatoms[i]].tab.table;
1080 gmx_fatal(FARGS,"A bonded table number is smaller than 0: %d\n",tabnr);
1081 if (tabnr >= *ncount) {
1082 srenew(*count,tabnr+1);
1083 for(j=*ncount; j<tabnr+1; j++)
1094 static bondedtable_t *make_bonded_tables(FILE *fplog,
1095 int ftype1,int ftype2,
1096 const gmx_mtop_t *mtop,
1097 const char *basefn,const char *tabext)
1099 int i,ncount,*count;
1107 count_tables(ftype1,ftype2,mtop,&ncount,&count);
1111 for(i=0; i<ncount; i++) {
1113 sprintf(tabfn,"%s",basefn);
1114 sprintf(tabfn + strlen(basefn) - strlen(ftp2ext(efXVG)) - 1,"_%s%d.%s",
1115 tabext,i,ftp2ext(efXVG));
1116 tab[i] = make_bonded_table(fplog,tabfn,NRAL(ftype1)-2);
1125 void forcerec_set_ranges(t_forcerec *fr,
1126 int ncg_home,int ncg_force,
1128 int natoms_force_constr,int natoms_f_novirsum)
1133 /* fr->ncg_force is unused in the standard code,
1134 * but it can be useful for modified code dealing with charge groups.
1136 fr->ncg_force = ncg_force;
1137 fr->natoms_force = natoms_force;
1138 fr->natoms_force_constr = natoms_force_constr;
1140 if (fr->natoms_force_constr > fr->nalloc_force)
1142 fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1146 srenew(fr->f_twin,fr->nalloc_force);
1150 if (fr->bF_NoVirSum)
1152 fr->f_novirsum_n = natoms_f_novirsum;
1153 if (fr->f_novirsum_n > fr->f_novirsum_nalloc)
1155 fr->f_novirsum_nalloc = over_alloc_dd(fr->f_novirsum_n);
1156 srenew(fr->f_novirsum_alloc,fr->f_novirsum_nalloc);
1161 fr->f_novirsum_n = 0;
1165 static real cutoff_inf(real cutoff)
1169 cutoff = GMX_CUTOFF_INF;
1175 gmx_bool can_use_allvsall(const t_inputrec *ir, const gmx_mtop_t *mtop,
1176 gmx_bool bPrintNote,t_commrec *cr,FILE *fp)
1185 ir->ePBC==epbcNONE &&
1186 ir->vdwtype==evdwCUT &&
1187 ir->coulombtype==eelCUT &&
1189 (ir->implicit_solvent == eisNO ||
1190 (ir->implicit_solvent==eisGBSA && (ir->gb_algorithm==egbSTILL ||
1191 ir->gb_algorithm==egbHCT ||
1192 ir->gb_algorithm==egbOBC))) &&
1193 getenv("GMX_NO_ALLVSALL") == NULL
1196 if (bAllvsAll && ir->opts.ngener > 1)
1198 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";
1204 fprintf(stderr,"\n%s\n",note);
1208 fprintf(fp,"\n%s\n",note);
1214 if(bAllvsAll && fp && MASTER(cr))
1216 fprintf(fp,"\nUsing accelerated all-vs-all kernels.\n\n");
1223 void init_forcerec(FILE *fp,
1224 const output_env_t oenv,
1227 const t_inputrec *ir,
1228 const gmx_mtop_t *mtop,
1229 const t_commrec *cr,
1235 gmx_bool bNoSolvOpt,
1238 int i,j,m,natoms,ngrp,negp_pp,negptable,egi,egj;
1244 gmx_bool bGenericKernelOnly;
1245 gmx_bool bTab,bSep14tab,bNormalnblists;
1247 int *nm_ind,egp_flags;
1249 fr->bDomDec = DOMAINDECOMP(cr);
1251 natoms = mtop->natoms;
1253 if (check_box(ir->ePBC,box))
1255 gmx_fatal(FARGS,check_box(ir->ePBC,box));
1258 /* Test particle insertion ? */
1259 if (EI_TPI(ir->eI)) {
1260 /* Set to the size of the molecule to be inserted (the last one) */
1261 /* Because of old style topologies, we have to use the last cg
1262 * instead of the last molecule type.
1264 cgs = &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].cgs;
1265 fr->n_tpi = cgs->index[cgs->nr] - cgs->index[cgs->nr-1];
1266 if (fr->n_tpi != mtop->mols.index[mtop->mols.nr] - mtop->mols.index[mtop->mols.nr-1]) {
1267 gmx_fatal(FARGS,"The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
1273 /* Copy the user determined parameters */
1274 fr->userint1 = ir->userint1;
1275 fr->userint2 = ir->userint2;
1276 fr->userint3 = ir->userint3;
1277 fr->userint4 = ir->userint4;
1278 fr->userreal1 = ir->userreal1;
1279 fr->userreal2 = ir->userreal2;
1280 fr->userreal3 = ir->userreal3;
1281 fr->userreal4 = ir->userreal4;
1284 fr->fc_stepsize = ir->fc_stepsize;
1287 fr->efep = ir->efep;
1288 fr->sc_alpha = ir->sc_alpha;
1289 fr->sc_power = ir->sc_power;
1290 fr->sc_sigma6_def = pow(ir->sc_sigma,6);
1291 fr->sc_sigma6_min = pow(ir->sc_sigma_min,6);
1292 env = getenv("GMX_SCSIGMA_MIN");
1296 sscanf(env,"%lf",&dbl);
1297 fr->sc_sigma6_min = pow(dbl,6);
1300 fprintf(fp,"Setting the minimum soft core sigma to %g nm\n",dbl);
1304 bGenericKernelOnly = FALSE;
1305 if (getenv("GMX_NB_GENERIC") != NULL)
1310 "Found environment variable GMX_NB_GENERIC.\n"
1311 "Disabling interaction-specific nonbonded kernels.\n\n");
1313 bGenericKernelOnly = TRUE;
1317 fr->UseOptimizedKernels = (getenv("GMX_NOOPTIMIZEDKERNELS") == NULL);
1318 if(fp && fr->UseOptimizedKernels==FALSE)
1321 "\nFound environment variable GMX_NOOPTIMIZEDKERNELS.\n"
1322 "Disabling SSE/SSE2/Altivec/ia64/Power6/Bluegene specific kernels.\n\n");
1325 /* Check if we can/should do all-vs-all kernels */
1326 fr->bAllvsAll = can_use_allvsall(ir,mtop,FALSE,NULL,NULL);
1327 fr->AllvsAll_work = NULL;
1328 fr->AllvsAll_workgb = NULL;
1332 /* Neighbour searching stuff */
1333 fr->bGrid = (ir->ns_type == ensGRID);
1334 fr->ePBC = ir->ePBC;
1335 fr->bMolPBC = ir->bPeriodicMols;
1336 fr->rc_scaling = ir->refcoord_scaling;
1337 copy_rvec(ir->posres_com,fr->posres_com);
1338 copy_rvec(ir->posres_comB,fr->posres_comB);
1339 fr->rlist = cutoff_inf(ir->rlist);
1340 fr->rlistlong = cutoff_inf(ir->rlistlong);
1341 fr->eeltype = ir->coulombtype;
1342 fr->vdwtype = ir->vdwtype;
1344 fr->bTwinRange = fr->rlistlong > fr->rlist;
1345 fr->bEwald = (EEL_PME(fr->eeltype) || fr->eeltype==eelEWALD);
1347 fr->reppow = mtop->ffparams.reppow;
1348 fr->bvdwtab = (fr->vdwtype != evdwCUT ||
1349 !gmx_within_tol(fr->reppow,12.0,10*GMX_DOUBLE_EPS));
1350 fr->bcoultab = (!(fr->eeltype == eelCUT || EEL_RF(fr->eeltype)) ||
1351 fr->eeltype == eelRF_ZERO);
1353 if (getenv("GMX_FORCE_TABLES"))
1356 fr->bcoultab = TRUE;
1360 fprintf(fp,"Table routines are used for coulomb: %s\n",bool_names[fr->bcoultab]);
1361 fprintf(fp,"Table routines are used for vdw: %s\n",bool_names[fr->bvdwtab ]);
1364 /* Tables are used for direct ewald sum */
1367 if (EEL_PME(ir->coulombtype))
1370 fprintf(fp,"Will do PME sum in reciprocal space.\n");
1371 please_cite(fp,"Essman95a");
1373 if (ir->ewald_geometry == eewg3DC)
1377 fprintf(fp,"Using the Ewald3DC correction for systems with a slab geometry.\n");
1379 please_cite(fp,"In-Chul99a");
1382 fr->ewaldcoeff=calc_ewaldcoeff(ir->rcoulomb, ir->ewald_rtol);
1383 init_ewald_tab(&(fr->ewald_table), cr, ir, fp);
1386 fprintf(fp,"Using a Gaussian width (1/beta) of %g nm for Ewald\n",
1391 /* Electrostatics */
1392 fr->epsilon_r = ir->epsilon_r;
1393 fr->epsilon_rf = ir->epsilon_rf;
1394 fr->fudgeQQ = mtop->ffparams.fudgeQQ;
1395 fr->rcoulomb_switch = ir->rcoulomb_switch;
1396 fr->rcoulomb = cutoff_inf(ir->rcoulomb);
1398 /* Parameters for generalized RF */
1402 if (fr->eeltype == eelGRF)
1404 init_generalized_rf(fp,mtop,ir,fr);
1406 else if (EEL_FULL(fr->eeltype) || (fr->eeltype == eelSHIFT) ||
1407 (fr->eeltype == eelUSER) || (fr->eeltype == eelSWITCH))
1409 /* We must use the long range cut-off for neighboursearching...
1410 * An extra range of e.g. 0.1 nm (half the size of a charge group)
1411 * is necessary for neighboursearching. This allows diffusion
1412 * into the cut-off range (between neighborlist updates),
1413 * and gives more accurate forces because all atoms within the short-range
1414 * cut-off rc must be taken into account, while the ns criterium takes
1415 * only those with the center of geometry within the cut-off.
1416 * (therefore we have to add half the size of a charge group, plus
1417 * something to account for diffusion if we have nstlist > 1)
1419 for(m=0; (m<DIM); m++)
1420 box_size[m]=box[m][m];
1422 if (fr->eeltype == eelPPPM && fr->phi == NULL)
1423 snew(fr->phi,natoms);
1425 if ((fr->eeltype==eelPPPM) || (fr->eeltype==eelPOISSON) ||
1426 (fr->eeltype == eelSHIFT && fr->rcoulomb > fr->rcoulomb_switch))
1427 set_shift_consts(fp,fr->rcoulomb_switch,fr->rcoulomb,box_size,fr);
1430 fr->bF_NoVirSum = (EEL_FULL(fr->eeltype) ||
1431 gmx_mtop_ftype_count(mtop,F_POSRES) > 0 ||
1432 IR_ELEC_FIELD(*ir));
1434 /* Mask that says whether or not this NBF list should be computed */
1435 /* if (fr->bMask == NULL) {
1436 ngrp = ir->opts.ngener*ir->opts.ngener;
1437 snew(fr->bMask,ngrp);*/
1438 /* Defaults to always */
1439 /* for(i=0; (i<ngrp); i++)
1440 fr->bMask[i] = TRUE;
1443 if (ncg_mtop(mtop) > fr->cg_nalloc && !DOMAINDECOMP(cr)) {
1444 /* Count the total number of charge groups */
1445 fr->cg_nalloc = ncg_mtop(mtop);
1446 srenew(fr->cg_cm,fr->cg_nalloc);
1448 if (fr->shift_vec == NULL)
1449 snew(fr->shift_vec,SHIFTS);
1451 if (fr->fshift == NULL)
1452 snew(fr->fshift,SHIFTS);
1454 if (fr->nbfp == NULL) {
1455 fr->ntype = mtop->ffparams.atnr;
1456 fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
1457 fr->nbfp = mk_nbfp(&mtop->ffparams,fr->bBHAM);
1460 /* Copy the energy group exclusions */
1461 fr->egp_flags = ir->opts.egp_flags;
1463 /* Van der Waals stuff */
1464 fr->rvdw = cutoff_inf(ir->rvdw);
1465 fr->rvdw_switch = ir->rvdw_switch;
1466 if ((fr->vdwtype != evdwCUT) && (fr->vdwtype != evdwUSER) && !fr->bBHAM) {
1467 if (fr->rvdw_switch >= fr->rvdw)
1468 gmx_fatal(FARGS,"rvdw_switch (%f) must be < rvdw (%f)",
1469 fr->rvdw_switch,fr->rvdw);
1471 fprintf(fp,"Using %s Lennard-Jones, switch between %g and %g nm\n",
1472 (fr->eeltype==eelSWITCH) ? "switched":"shifted",
1473 fr->rvdw_switch,fr->rvdw);
1476 if (fr->bBHAM && (fr->vdwtype == evdwSHIFT || fr->vdwtype == evdwSWITCH))
1477 gmx_fatal(FARGS,"Switch/shift interaction not supported with Buckingham");
1480 fprintf(fp,"Cut-off's: NS: %g Coulomb: %g %s: %g\n",
1481 fr->rlist,fr->rcoulomb,fr->bBHAM ? "BHAM":"LJ",fr->rvdw);
1483 fr->eDispCorr = ir->eDispCorr;
1484 if (ir->eDispCorr != edispcNO)
1486 set_avcsixtwelve(fp,fr,mtop);
1491 set_bham_b_max(fp,fr,mtop);
1494 fr->bGB = (ir->implicit_solvent == eisGBSA);
1495 fr->gb_epsilon_solvent = ir->gb_epsilon_solvent;
1497 /* Copy the GBSA data (radius, volume and surftens for each
1498 * atomtype) from the topology atomtype section to forcerec.
1500 snew(fr->atype_radius,fr->ntype);
1501 snew(fr->atype_vol,fr->ntype);
1502 snew(fr->atype_surftens,fr->ntype);
1503 snew(fr->atype_gb_radius,fr->ntype);
1504 snew(fr->atype_S_hct,fr->ntype);
1506 if (mtop->atomtypes.nr > 0)
1508 for(i=0;i<fr->ntype;i++)
1509 fr->atype_radius[i] =mtop->atomtypes.radius[i];
1510 for(i=0;i<fr->ntype;i++)
1511 fr->atype_vol[i] = mtop->atomtypes.vol[i];
1512 for(i=0;i<fr->ntype;i++)
1513 fr->atype_surftens[i] = mtop->atomtypes.surftens[i];
1514 for(i=0;i<fr->ntype;i++)
1515 fr->atype_gb_radius[i] = mtop->atomtypes.gb_radius[i];
1516 for(i=0;i<fr->ntype;i++)
1517 fr->atype_S_hct[i] = mtop->atomtypes.S_hct[i];
1520 /* Generate the GB table if needed */
1524 fr->gbtabscale=2000;
1530 fr->gbtab=make_gb_table(fp,oenv,fr,tabpfn,fr->gbtabscale);
1532 init_gb(&fr->born,cr,fr,ir,mtop,ir->rgbradii,ir->gb_algorithm);
1534 /* Copy local gb data (for dd, this is done in dd_partition_system) */
1535 if (!DOMAINDECOMP(cr))
1537 make_local_gb(cr,fr->born,ir->gb_algorithm);
1541 /* Set the charge scaling */
1542 if (fr->epsilon_r != 0)
1543 fr->epsfac = ONE_4PI_EPS0/fr->epsilon_r;
1545 /* eps = 0 is infinite dieletric: no coulomb interactions */
1548 /* Reaction field constants */
1549 if (EEL_RF(fr->eeltype))
1550 calc_rffac(fp,fr->eeltype,fr->epsilon_r,fr->epsilon_rf,
1551 fr->rcoulomb,fr->temp,fr->zsquare,box,
1552 &fr->kappa,&fr->k_rf,&fr->c_rf);
1554 set_chargesum(fp,fr,mtop);
1556 /* if we are using LR electrostatics, and they are tabulated,
1557 * the tables will contain modified coulomb interactions.
1558 * Since we want to use the non-shifted ones for 1-4
1559 * coulombic interactions, we must have an extra set of tables.
1562 /* Construct tables.
1563 * A little unnecessary to make both vdw and coul tables sometimes,
1564 * but what the heck... */
1566 bTab = fr->bcoultab || fr->bvdwtab;
1568 bSep14tab = ((!bTab || fr->eeltype!=eelCUT || fr->vdwtype!=evdwCUT ||
1570 (gmx_mtop_ftype_count(mtop,F_LJ14) > 0 ||
1571 gmx_mtop_ftype_count(mtop,F_LJC14_Q) > 0 ||
1572 gmx_mtop_ftype_count(mtop,F_LJC_PAIRS_NB) > 0));
1574 negp_pp = ir->opts.ngener - ir->nwall;
1577 bNormalnblists = TRUE;
1580 bNormalnblists = (ir->eDispCorr != edispcNO);
1581 for(egi=0; egi<negp_pp; egi++) {
1582 for(egj=egi; egj<negp_pp; egj++) {
1583 egp_flags = ir->opts.egp_flags[GID(egi,egj,ir->opts.ngener)];
1584 if (!(egp_flags & EGP_EXCL)) {
1585 if (egp_flags & EGP_TABLE) {
1588 bNormalnblists = TRUE;
1593 if (bNormalnblists) {
1594 fr->nnblists = negptable + 1;
1596 fr->nnblists = negptable;
1598 if (fr->nnblists > 1)
1599 snew(fr->gid2nblists,ir->opts.ngener*ir->opts.ngener);
1601 snew(fr->nblists,fr->nnblists);
1603 /* This code automatically gives table length tabext without cut-off's,
1604 * in that case grompp should already have checked that we do not need
1605 * normal tables and we only generate tables for 1-4 interactions.
1607 rtab = ir->rlistlong + ir->tabext;
1610 /* make tables for ordinary interactions */
1611 if (bNormalnblists) {
1612 make_nbf_tables(fp,oenv,fr,rtab,cr,tabfn,NULL,NULL,&fr->nblists[0]);
1614 fr->tab14 = fr->nblists[0].tab;
1619 if (negptable > 0) {
1620 /* Read the special tables for certain energy group pairs */
1621 nm_ind = mtop->groups.grps[egcENER].nm_ind;
1622 for(egi=0; egi<negp_pp; egi++) {
1623 for(egj=egi; egj<negp_pp; egj++) {
1624 egp_flags = ir->opts.egp_flags[GID(egi,egj,ir->opts.ngener)];
1625 if ((egp_flags & EGP_TABLE) && !(egp_flags & EGP_EXCL)) {
1626 nbl = &(fr->nblists[m]);
1627 if (fr->nnblists > 1) {
1628 fr->gid2nblists[GID(egi,egj,ir->opts.ngener)] = m;
1630 /* Read the table file with the two energy groups names appended */
1631 make_nbf_tables(fp,oenv,fr,rtab,cr,tabfn,
1632 *mtop->groups.grpname[nm_ind[egi]],
1633 *mtop->groups.grpname[nm_ind[egj]],
1636 } else if (fr->nnblists > 1) {
1637 fr->gid2nblists[GID(egi,egj,ir->opts.ngener)] = 0;
1645 /* generate extra tables with plain Coulomb for 1-4 interactions only */
1646 fr->tab14 = make_tables(fp,oenv,fr,MASTER(cr),tabpfn,rtab,
1647 GMX_MAKETABLES_14ONLY);
1651 fr->nwall = ir->nwall;
1652 if (ir->nwall && ir->wall_type==ewtTABLE)
1654 make_wall_tables(fp,oenv,ir,tabfn,&mtop->groups,fr);
1657 if (fcd && tabbfn) {
1658 fcd->bondtab = make_bonded_tables(fp,
1659 F_TABBONDS,F_TABBONDSNC,
1661 fcd->angletab = make_bonded_tables(fp,
1664 fcd->dihtab = make_bonded_tables(fp,
1669 fprintf(debug,"No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
1672 if (ir->eDispCorr != edispcNO)
1674 calc_enervirdiff(fp,ir->eDispCorr,fr);
1677 /* QM/MM initialization if requested
1681 fprintf(stderr,"QM/MM calculation requested.\n");
1684 fr->bQMMM = ir->bQMMM;
1685 fr->qr = mk_QMMMrec();
1687 /* Set all the static charge group info */
1688 fr->cginfo_mb = init_cginfo_mb(fp,mtop,fr,bNoSolvOpt,
1689 &fr->bExcl_IntraCGAll_InterCGNone);
1690 if (DOMAINDECOMP(cr)) {
1693 fr->cginfo = cginfo_expand(mtop->nmolblock,fr->cginfo_mb);
1696 if (!DOMAINDECOMP(cr))
1698 /* When using particle decomposition, the effect of the second argument,
1699 * which sets fr->hcg, is corrected later in do_md and init_em.
1701 forcerec_set_ranges(fr,ncg_mtop(mtop),ncg_mtop(mtop),
1702 mtop->natoms,mtop->natoms,mtop->natoms);
1705 fr->print_force = print_force;
1708 /* coarse load balancing vars */
1713 /* Initialize neighbor search */
1714 init_ns(fp,cr,&fr->ns,fr,mtop,box);
1716 if (cr->duty & DUTY_PP)
1717 gmx_setup_kernels(fp,bGenericKernelOnly);
1720 #define pr_real(fp,r) fprintf(fp,"%s: %e\n",#r,r)
1721 #define pr_int(fp,i) fprintf((fp),"%s: %d\n",#i,i)
1722 #define pr_bool(fp,b) fprintf((fp),"%s: %s\n",#b,bool_names[b])
1724 void pr_forcerec(FILE *fp,t_forcerec *fr,t_commrec *cr)
1728 pr_real(fp,fr->rlist);
1729 pr_real(fp,fr->rcoulomb);
1730 pr_real(fp,fr->fudgeQQ);
1731 pr_bool(fp,fr->bGrid);
1732 pr_bool(fp,fr->bTwinRange);
1733 /*pr_int(fp,fr->cg0);
1734 pr_int(fp,fr->hcg);*/
1735 for(i=0; i<fr->nnblists; i++)
1736 pr_int(fp,fr->nblists[i].tab.n);
1737 pr_real(fp,fr->rcoulomb_switch);
1738 pr_real(fp,fr->rcoulomb);