2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
44 #include "gromacs/utility/smalloc.h"
46 #include "gromacs/math/utilities.h"
48 #include "types/commrec.h"
52 #include "nonbonded.h"
56 #include "gmx_fatal.h"
59 #include "mtop_util.h"
66 * E X C L U S I O N H A N D L I N G
70 static void SETEXCL_(t_excl e[], atom_id i, atom_id j)
74 static void RMEXCL_(t_excl e[], atom_id i, atom_id j)
76 e[j] = e[j] & ~(1<<i);
78 static gmx_bool ISEXCL_(t_excl e[], atom_id i, atom_id j)
80 return (gmx_bool)(e[j] & (1<<i));
82 static gmx_bool NOTEXCL_(t_excl e[], atom_id i, atom_id j)
84 return !(ISEXCL(e, i, j));
87 #define SETEXCL(e, i, j) (e)[((atom_id) (j))] |= (1<<((atom_id) (i)))
88 #define RMEXCL(e, i, j) (e)[((atom_id) (j))] &= (~(1<<((atom_id) (i))))
89 #define ISEXCL(e, i, j) (gmx_bool) ((e)[((atom_id) (j))] & (1<<((atom_id) (i))))
90 #define NOTEXCL(e, i, j) !(ISEXCL(e, i, j))
94 round_up_to_simd_width(int length, int simd_width)
96 int offset, newlength;
98 offset = (simd_width > 0) ? length % simd_width : 0;
100 return (offset == 0) ? length : length-offset+simd_width;
102 /************************************************
104 * U T I L I T I E S F O R N S
106 ************************************************/
108 void reallocate_nblist(t_nblist *nl)
112 fprintf(debug, "reallocating neigborlist (ielec=%d, ivdw=%d, igeometry=%d, type=%d), maxnri=%d\n",
113 nl->ielec, nl->ivdw, nl->igeometry, nl->type, nl->maxnri);
115 srenew(nl->iinr, nl->maxnri);
116 if (nl->igeometry == GMX_NBLIST_GEOMETRY_CG_CG)
118 srenew(nl->iinr_end, nl->maxnri);
120 srenew(nl->gid, nl->maxnri);
121 srenew(nl->shift, nl->maxnri);
122 srenew(nl->jindex, nl->maxnri+1);
126 static void init_nblist(FILE *log, t_nblist *nl_sr, t_nblist *nl_lr,
127 int maxsr, int maxlr,
128 int ivdw, int ivdwmod,
129 int ielec, int ielecmod,
130 int igeometry, int type)
136 for (i = 0; (i < 2); i++)
138 nl = (i == 0) ? nl_sr : nl_lr;
139 homenr = (i == 0) ? maxsr : maxlr;
147 /* Set coul/vdw in neighborlist, and for the normal loops we determine
148 * an index of which one to call.
151 nl->ivdwmod = ivdwmod;
153 nl->ielecmod = ielecmod;
155 nl->igeometry = igeometry;
157 if (nl->type == GMX_NBLIST_INTERACTION_FREE_ENERGY)
159 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
162 /* This will also set the simd_padding_width field */
163 gmx_nonbonded_set_kernel_pointers( (i == 0) ? log : NULL, nl);
165 /* maxnri is influenced by the number of shifts (maximum is 8)
166 * and the number of energy groups.
167 * If it is not enough, nl memory will be reallocated during the run.
168 * 4 seems to be a reasonable factor, which only causes reallocation
169 * during runs with tiny and many energygroups.
171 nl->maxnri = homenr*4;
181 reallocate_nblist(nl);
186 fprintf(debug, "Initiating neighbourlist (ielec=%d, ivdw=%d, type=%d) for %s interactions,\nwith %d SR, %d LR atoms.\n",
187 nl->ielec, nl->ivdw, nl->type, gmx_nblist_geometry_names[nl->igeometry], maxsr, maxlr);
192 void init_neighbor_list(FILE *log, t_forcerec *fr, int homenr)
194 /* Make maxlr tunable! (does not seem to be a big difference though)
195 * This parameter determines the number of i particles in a long range
196 * neighbourlist. Too few means many function calls, too many means
199 int maxsr, maxsr_wat, maxlr, maxlr_wat;
200 int ielec, ielecf, ivdw, ielecmod, ielecmodf, ivdwmod, type;
202 int igeometry_def, igeometry_w, igeometry_ww;
206 /* maxsr = homenr-fr->nWatMol*3; */
211 gmx_fatal(FARGS, "%s, %d: Negative number of short range atoms.\n"
212 "Call your Gromacs dealer for assistance.", __FILE__, __LINE__);
214 /* This is just for initial allocation, so we do not reallocate
215 * all the nlist arrays many times in a row.
216 * The numbers seem very accurate, but they are uncritical.
218 maxsr_wat = min(fr->nWatMol, (homenr+2)/3);
222 maxlr_wat = min(maxsr_wat, maxlr);
226 maxlr = maxlr_wat = 0;
229 /* Determine the values for ielec/ivdw. */
230 ielec = fr->nbkernel_elec_interaction;
231 ivdw = fr->nbkernel_vdw_interaction;
232 ielecmod = fr->nbkernel_elec_modifier;
233 ivdwmod = fr->nbkernel_vdw_modifier;
234 type = GMX_NBLIST_INTERACTION_STANDARD;
236 fr->ns.bCGlist = (getenv("GMX_NBLISTCG") != 0);
239 igeometry_def = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
243 igeometry_def = GMX_NBLIST_GEOMETRY_CG_CG;
246 fprintf(log, "\nUsing charge-group - charge-group neighbor lists and kernels\n\n");
250 if (fr->solvent_opt == esolTIP4P)
252 igeometry_w = GMX_NBLIST_GEOMETRY_WATER4_PARTICLE;
253 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER4_WATER4;
257 igeometry_w = GMX_NBLIST_GEOMETRY_WATER3_PARTICLE;
258 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER3_WATER3;
261 for (i = 0; i < fr->nnblists; i++)
263 nbl = &(fr->nblists[i]);
265 if ((fr->adress_type != eAdressOff) && (i >= fr->nnblists/2))
267 type = GMX_NBLIST_INTERACTION_ADRESS;
269 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ], &nbl->nlist_lr[eNL_VDWQQ],
270 maxsr, maxlr, ivdw, ivdwmod, ielec, ielecmod, igeometry_def, type);
271 init_nblist(log, &nbl->nlist_sr[eNL_VDW], &nbl->nlist_lr[eNL_VDW],
272 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, igeometry_def, type);
273 init_nblist(log, &nbl->nlist_sr[eNL_QQ], &nbl->nlist_lr[eNL_QQ],
274 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_def, type);
275 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATER], &nbl->nlist_lr[eNL_VDWQQ_WATER],
276 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_w, type);
277 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATER], &nbl->nlist_lr[eNL_QQ_WATER],
278 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_w, type);
279 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATERWATER], &nbl->nlist_lr[eNL_VDWQQ_WATERWATER],
280 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_ww, type);
281 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATERWATER], &nbl->nlist_lr[eNL_QQ_WATERWATER],
282 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_ww, type);
284 /* Did we get the solvent loops so we can use optimized water kernels? */
285 if (nbl->nlist_sr[eNL_VDWQQ_WATER].kernelptr_vf == NULL
286 || nbl->nlist_sr[eNL_QQ_WATER].kernelptr_vf == NULL
287 #ifndef DISABLE_WATERWATER_NLIST
288 || nbl->nlist_sr[eNL_VDWQQ_WATERWATER].kernelptr_vf == NULL
289 || nbl->nlist_sr[eNL_QQ_WATERWATER].kernelptr_vf == NULL
293 fr->solvent_opt = esolNO;
294 fprintf(log, "Note: The available nonbonded kernels do not support water optimization - disabling.\n");
297 if (fr->efep != efepNO)
299 if ((fr->bEwald) && (fr->sc_alphacoul > 0)) /* need to handle long range differently if using softcore */
301 ielecf = GMX_NBKERNEL_ELEC_EWALD;
302 ielecmodf = eintmodNONE;
307 ielecmodf = ielecmod;
310 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_FREE], &nbl->nlist_lr[eNL_VDWQQ_FREE],
311 maxsr, maxlr, ivdw, ivdwmod, ielecf, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
312 init_nblist(log, &nbl->nlist_sr[eNL_VDW_FREE], &nbl->nlist_lr[eNL_VDW_FREE],
313 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
314 init_nblist(log, &nbl->nlist_sr[eNL_QQ_FREE], &nbl->nlist_lr[eNL_QQ_FREE],
315 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielecf, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
319 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
321 init_nblist(log, &fr->QMMMlist, NULL,
322 maxsr, maxlr, 0, 0, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_STANDARD);
330 fr->ns.nblist_initialized = TRUE;
333 static void reset_nblist(t_nblist *nl)
343 static void reset_neighbor_lists(t_forcerec *fr, gmx_bool bResetSR, gmx_bool bResetLR)
349 /* only reset the short-range nblist */
350 reset_nblist(&(fr->QMMMlist));
353 for (n = 0; n < fr->nnblists; n++)
355 for (i = 0; i < eNL_NR; i++)
359 reset_nblist( &(fr->nblists[n].nlist_sr[i]) );
363 reset_nblist( &(fr->nblists[n].nlist_lr[i]) );
372 static gmx_inline void new_i_nblist(t_nblist *nlist, atom_id i_atom, int shift, int gid)
374 int i, k, nri, nshift;
378 /* Check whether we have to increase the i counter */
380 (nlist->iinr[nri] != i_atom) ||
381 (nlist->shift[nri] != shift) ||
382 (nlist->gid[nri] != gid))
384 /* This is something else. Now see if any entries have
385 * been added in the list of the previous atom.
388 ((nlist->jindex[nri+1] > nlist->jindex[nri]) &&
389 (nlist->gid[nri] != -1)))
391 /* If so increase the counter */
394 if (nlist->nri >= nlist->maxnri)
396 nlist->maxnri += over_alloc_large(nlist->nri);
397 reallocate_nblist(nlist);
400 /* Set the number of neighbours and the atom number */
401 nlist->jindex[nri+1] = nlist->jindex[nri];
402 nlist->iinr[nri] = i_atom;
403 nlist->gid[nri] = gid;
404 nlist->shift[nri] = shift;
408 /* Adding to previous list. First remove possible previous padding */
409 if (nlist->simd_padding_width > 1)
411 while (nlist->nrj > 0 && nlist->jjnr[nlist->nrj-1] < 0)
419 static gmx_inline void close_i_nblist(t_nblist *nlist)
421 int nri = nlist->nri;
426 /* Add elements up to padding. Since we allocate memory in units
427 * of the simd_padding width, we do not have to check for possible
428 * list reallocation here.
430 while ((nlist->nrj % nlist->simd_padding_width) != 0)
432 /* Use -4 here, so we can write forces for 4 atoms before real data */
433 nlist->jjnr[nlist->nrj++] = -4;
435 nlist->jindex[nri+1] = nlist->nrj;
437 len = nlist->nrj - nlist->jindex[nri];
441 static gmx_inline void close_nblist(t_nblist *nlist)
443 /* Only close this nblist when it has been initialized.
444 * Avoid the creation of i-lists with no j-particles.
448 /* Some assembly kernels do not support empty lists,
449 * make sure here that we don't generate any empty lists.
450 * With the current ns code this branch is taken in two cases:
451 * No i-particles at all: nri=-1 here
452 * There are i-particles, but no j-particles; nri=0 here
458 /* Close list number nri by incrementing the count */
463 static gmx_inline void close_neighbor_lists(t_forcerec *fr, gmx_bool bMakeQMMMnblist)
469 close_nblist(&(fr->QMMMlist));
472 for (n = 0; n < fr->nnblists; n++)
474 for (i = 0; (i < eNL_NR); i++)
476 close_nblist(&(fr->nblists[n].nlist_sr[i]));
477 close_nblist(&(fr->nblists[n].nlist_lr[i]));
483 static gmx_inline void add_j_to_nblist(t_nblist *nlist, atom_id j_atom, gmx_bool bLR)
485 int nrj = nlist->nrj;
487 if (nlist->nrj >= nlist->maxnrj)
489 nlist->maxnrj = round_up_to_simd_width(over_alloc_small(nlist->nrj + 1), nlist->simd_padding_width);
493 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
494 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
497 srenew(nlist->jjnr, nlist->maxnrj);
500 nlist->jjnr[nrj] = j_atom;
504 static gmx_inline void add_j_to_nblist_cg(t_nblist *nlist,
505 atom_id j_start, int j_end,
506 t_excl *bexcl, gmx_bool i_is_j,
509 int nrj = nlist->nrj;
512 if (nlist->nrj >= nlist->maxnrj)
514 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
517 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
518 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
521 srenew(nlist->jjnr, nlist->maxnrj);
522 srenew(nlist->jjnr_end, nlist->maxnrj);
523 srenew(nlist->excl, nlist->maxnrj*MAX_CGCGSIZE);
526 nlist->jjnr[nrj] = j_start;
527 nlist->jjnr_end[nrj] = j_end;
529 if (j_end - j_start > MAX_CGCGSIZE)
531 gmx_fatal(FARGS, "The charge-group - charge-group neighborlist do not support charge groups larger than %d, found a charge group of size %d", MAX_CGCGSIZE, j_end-j_start);
534 /* Set the exclusions */
535 for (j = j_start; j < j_end; j++)
537 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
541 /* Avoid double counting of intra-cg interactions */
542 for (j = 1; j < j_end-j_start; j++)
544 nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
552 put_in_list_t (gmx_bool bHaveVdW[],
569 put_in_list_at(gmx_bool bHaveVdW[],
585 /* The a[] index has been removed,
586 * to put it back in i_atom should be a[i0] and jj should be a[jj].
591 t_nblist * vdwc_free = NULL;
592 t_nblist * vdw_free = NULL;
593 t_nblist * coul_free = NULL;
594 t_nblist * vdwc_ww = NULL;
595 t_nblist * coul_ww = NULL;
597 int i, j, jcg, igid, gid, nbl_ind, ind_ij;
598 atom_id jj, jj0, jj1, i_atom;
603 real *charge, *chargeB;
604 real qi, qiB, qq, rlj;
605 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
606 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
610 /* Copy some pointers */
612 charge = md->chargeA;
613 chargeB = md->chargeB;
616 bPert = md->bPerturbed;
620 nicg = index[icg+1]-i0;
622 /* Get the i charge group info */
623 igid = GET_CGINFO_GID(cginfo[icg]);
625 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
630 /* Check if any of the particles involved are perturbed.
631 * If not we can do the cheaper normal put_in_list
632 * and use more solvent optimization.
634 for (i = 0; i < nicg; i++)
636 bFreeEnergy |= bPert[i0+i];
638 /* Loop over the j charge groups */
639 for (j = 0; (j < nj && !bFreeEnergy); j++)
644 /* Finally loop over the atoms in the j-charge group */
645 for (jj = jj0; jj < jj1; jj++)
647 bFreeEnergy |= bPert[jj];
652 /* Unpack pointers to neighbourlist structs */
653 if (fr->nnblists == 1)
659 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
663 nlist = fr->nblists[nbl_ind].nlist_lr;
667 nlist = fr->nblists[nbl_ind].nlist_sr;
670 if (iwater != esolNO)
672 vdwc = &nlist[eNL_VDWQQ_WATER];
673 vdw = &nlist[eNL_VDW];
674 coul = &nlist[eNL_QQ_WATER];
675 #ifndef DISABLE_WATERWATER_NLIST
676 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
677 coul_ww = &nlist[eNL_QQ_WATERWATER];
682 vdwc = &nlist[eNL_VDWQQ];
683 vdw = &nlist[eNL_VDW];
684 coul = &nlist[eNL_QQ];
689 if (iwater != esolNO)
691 /* Loop over the atoms in the i charge group */
693 gid = GID(igid, jgid, ngid);
694 /* Create new i_atom for each energy group */
695 if (bDoCoul && bDoVdW)
697 new_i_nblist(vdwc, i_atom, shift, gid);
698 #ifndef DISABLE_WATERWATER_NLIST
699 new_i_nblist(vdwc_ww, i_atom, shift, gid);
704 new_i_nblist(vdw, i_atom, shift, gid);
708 new_i_nblist(coul, i_atom, shift, gid);
709 #ifndef DISABLE_WATERWATER_NLIST
710 new_i_nblist(coul_ww, i_atom, shift, gid);
713 /* Loop over the j charge groups */
714 for (j = 0; (j < nj); j++)
724 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
726 if (iwater == esolSPC && jwater == esolSPC)
728 /* Interaction between two SPC molecules */
731 /* VdW only - only first atoms in each water interact */
732 add_j_to_nblist(vdw, jj0, bLR);
736 #ifdef DISABLE_WATERWATER_NLIST
737 /* Add entries for the three atoms - only do VdW if we need to */
740 add_j_to_nblist(coul, jj0, bLR);
744 add_j_to_nblist(vdwc, jj0, bLR);
746 add_j_to_nblist(coul, jj0+1, bLR);
747 add_j_to_nblist(coul, jj0+2, bLR);
749 /* One entry for the entire water-water interaction */
752 add_j_to_nblist(coul_ww, jj0, bLR);
756 add_j_to_nblist(vdwc_ww, jj0, bLR);
761 else if (iwater == esolTIP4P && jwater == esolTIP4P)
763 /* Interaction between two TIP4p molecules */
766 /* VdW only - only first atoms in each water interact */
767 add_j_to_nblist(vdw, jj0, bLR);
771 #ifdef DISABLE_WATERWATER_NLIST
772 /* Add entries for the four atoms - only do VdW if we need to */
775 add_j_to_nblist(vdw, jj0, bLR);
777 add_j_to_nblist(coul, jj0+1, bLR);
778 add_j_to_nblist(coul, jj0+2, bLR);
779 add_j_to_nblist(coul, jj0+3, bLR);
781 /* One entry for the entire water-water interaction */
784 add_j_to_nblist(coul_ww, jj0, bLR);
788 add_j_to_nblist(vdwc_ww, jj0, bLR);
795 /* j charge group is not water, but i is.
796 * Add entries to the water-other_atom lists; the geometry of the water
797 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
798 * so we don't care if it is SPC or TIP4P...
805 for (jj = jj0; (jj < jj1); jj++)
809 add_j_to_nblist(coul, jj, bLR);
815 for (jj = jj0; (jj < jj1); jj++)
817 if (bHaveVdW[type[jj]])
819 add_j_to_nblist(vdw, jj, bLR);
825 /* _charge_ _groups_ interact with both coulomb and LJ */
826 /* Check which atoms we should add to the lists! */
827 for (jj = jj0; (jj < jj1); jj++)
829 if (bHaveVdW[type[jj]])
833 add_j_to_nblist(vdwc, jj, bLR);
837 add_j_to_nblist(vdw, jj, bLR);
840 else if (charge[jj] != 0)
842 add_j_to_nblist(coul, jj, bLR);
849 close_i_nblist(coul);
850 close_i_nblist(vdwc);
851 #ifndef DISABLE_WATERWATER_NLIST
852 close_i_nblist(coul_ww);
853 close_i_nblist(vdwc_ww);
858 /* no solvent as i charge group */
859 /* Loop over the atoms in the i charge group */
860 for (i = 0; i < nicg; i++)
863 gid = GID(igid, jgid, ngid);
866 /* Create new i_atom for each energy group */
867 if (bDoVdW && bDoCoul)
869 new_i_nblist(vdwc, i_atom, shift, gid);
873 new_i_nblist(vdw, i_atom, shift, gid);
877 new_i_nblist(coul, i_atom, shift, gid);
879 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
880 bDoCoul_i = (bDoCoul && qi != 0);
882 if (bDoVdW_i || bDoCoul_i)
884 /* Loop over the j charge groups */
885 for (j = 0; (j < nj); j++)
889 /* Check for large charge groups */
900 /* Finally loop over the atoms in the j-charge group */
901 for (jj = jj0; jj < jj1; jj++)
903 bNotEx = NOTEXCL(bExcl, i, jj);
911 add_j_to_nblist(coul, jj, bLR);
916 if (bHaveVdW[type[jj]])
918 add_j_to_nblist(vdw, jj, bLR);
923 if (bHaveVdW[type[jj]])
927 add_j_to_nblist(vdwc, jj, bLR);
931 add_j_to_nblist(vdw, jj, bLR);
934 else if (charge[jj] != 0)
936 add_j_to_nblist(coul, jj, bLR);
944 close_i_nblist(coul);
945 close_i_nblist(vdwc);
951 /* we are doing free energy */
952 vdwc_free = &nlist[eNL_VDWQQ_FREE];
953 vdw_free = &nlist[eNL_VDW_FREE];
954 coul_free = &nlist[eNL_QQ_FREE];
955 /* Loop over the atoms in the i charge group */
956 for (i = 0; i < nicg; i++)
959 gid = GID(igid, jgid, ngid);
961 qiB = chargeB[i_atom];
963 /* Create new i_atom for each energy group */
964 if (bDoVdW && bDoCoul)
966 new_i_nblist(vdwc, i_atom, shift, gid);
970 new_i_nblist(vdw, i_atom, shift, gid);
974 new_i_nblist(coul, i_atom, shift, gid);
977 new_i_nblist(vdw_free, i_atom, shift, gid);
978 new_i_nblist(coul_free, i_atom, shift, gid);
979 new_i_nblist(vdwc_free, i_atom, shift, gid);
981 bDoVdW_i = (bDoVdW &&
982 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
983 bDoCoul_i = (bDoCoul && (qi != 0 || qiB != 0));
984 /* For TIP4P the first atom does not have a charge,
985 * but the last three do. So we should still put an atom
986 * without LJ but with charge in the water-atom neighborlist
987 * for a TIP4p i charge group.
988 * For SPC type water the first atom has LJ and charge,
989 * so there is no such problem.
991 if (iwater == esolNO)
993 bDoCoul_i_sol = bDoCoul_i;
997 bDoCoul_i_sol = bDoCoul;
1000 if (bDoVdW_i || bDoCoul_i_sol)
1002 /* Loop over the j charge groups */
1003 for (j = 0; (j < nj); j++)
1007 /* Check for large charge groups */
1018 /* Finally loop over the atoms in the j-charge group */
1019 bFree = bPert[i_atom];
1020 for (jj = jj0; (jj < jj1); jj++)
1022 bFreeJ = bFree || bPert[jj];
1023 /* Complicated if, because the water H's should also
1024 * see perturbed j-particles
1026 if (iwater == esolNO || i == 0 || bFreeJ)
1028 bNotEx = NOTEXCL(bExcl, i, jj);
1036 if (charge[jj] != 0 || chargeB[jj] != 0)
1038 add_j_to_nblist(coul_free, jj, bLR);
1041 else if (!bDoCoul_i)
1043 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1045 add_j_to_nblist(vdw_free, jj, bLR);
1050 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1052 if (charge[jj] != 0 || chargeB[jj] != 0)
1054 add_j_to_nblist(vdwc_free, jj, bLR);
1058 add_j_to_nblist(vdw_free, jj, bLR);
1061 else if (charge[jj] != 0 || chargeB[jj] != 0)
1063 add_j_to_nblist(coul_free, jj, bLR);
1069 /* This is done whether or not bWater is set */
1070 if (charge[jj] != 0)
1072 add_j_to_nblist(coul, jj, bLR);
1075 else if (!bDoCoul_i_sol)
1077 if (bHaveVdW[type[jj]])
1079 add_j_to_nblist(vdw, jj, bLR);
1084 if (bHaveVdW[type[jj]])
1086 if (charge[jj] != 0)
1088 add_j_to_nblist(vdwc, jj, bLR);
1092 add_j_to_nblist(vdw, jj, bLR);
1095 else if (charge[jj] != 0)
1097 add_j_to_nblist(coul, jj, bLR);
1105 close_i_nblist(vdw);
1106 close_i_nblist(coul);
1107 close_i_nblist(vdwc);
1108 close_i_nblist(vdw_free);
1109 close_i_nblist(coul_free);
1110 close_i_nblist(vdwc_free);
1116 put_in_list_adress(gmx_bool bHaveVdW[],
1132 /* The a[] index has been removed,
1133 * to put it back in i_atom should be a[i0] and jj should be a[jj].
1138 t_nblist * vdwc_adress = NULL;
1139 t_nblist * vdw_adress = NULL;
1140 t_nblist * coul_adress = NULL;
1141 t_nblist * vdwc_ww = NULL;
1142 t_nblist * coul_ww = NULL;
1144 int i, j, jcg, igid, gid, nbl_ind, nbl_ind_adress;
1145 atom_id jj, jj0, jj1, i_atom;
1150 real *charge, *chargeB;
1152 real qi, qiB, qq, rlj;
1153 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
1154 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
1156 gmx_bool j_all_atom;
1158 t_nblist *nlist, *nlist_adress;
1159 gmx_bool bEnergyGroupCG;
1161 /* Copy some pointers */
1162 cginfo = fr->cginfo;
1163 charge = md->chargeA;
1164 chargeB = md->chargeB;
1167 bPert = md->bPerturbed;
1170 /* Get atom range */
1172 nicg = index[icg+1]-i0;
1174 /* Get the i charge group info */
1175 igid = GET_CGINFO_GID(cginfo[icg]);
1177 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
1181 gmx_fatal(FARGS, "AdResS does not support free energy pertubation\n");
1184 /* Unpack pointers to neighbourlist structs */
1185 if (fr->nnblists == 2)
1192 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
1193 nbl_ind_adress = nbl_ind+fr->nnblists/2;
1197 nlist = fr->nblists[nbl_ind].nlist_lr;
1198 nlist_adress = fr->nblists[nbl_ind_adress].nlist_lr;
1202 nlist = fr->nblists[nbl_ind].nlist_sr;
1203 nlist_adress = fr->nblists[nbl_ind_adress].nlist_sr;
1207 vdwc = &nlist[eNL_VDWQQ];
1208 vdw = &nlist[eNL_VDW];
1209 coul = &nlist[eNL_QQ];
1211 vdwc_adress = &nlist_adress[eNL_VDWQQ];
1212 vdw_adress = &nlist_adress[eNL_VDW];
1213 coul_adress = &nlist_adress[eNL_QQ];
1215 /* We do not support solvent optimization with AdResS for now.
1216 For this we would need hybrid solvent-other kernels */
1218 /* no solvent as i charge group */
1219 /* Loop over the atoms in the i charge group */
1220 for (i = 0; i < nicg; i++)
1223 gid = GID(igid, jgid, ngid);
1224 qi = charge[i_atom];
1226 /* Create new i_atom for each energy group */
1227 if (bDoVdW && bDoCoul)
1229 new_i_nblist(vdwc, i_atom, shift, gid);
1230 new_i_nblist(vdwc_adress, i_atom, shift, gid);
1235 new_i_nblist(vdw, i_atom, shift, gid);
1236 new_i_nblist(vdw_adress, i_atom, shift, gid);
1241 new_i_nblist(coul, i_atom, shift, gid);
1242 new_i_nblist(coul_adress, i_atom, shift, gid);
1244 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
1245 bDoCoul_i = (bDoCoul && qi != 0);
1247 /* Here we find out whether the energy groups interaction belong to a
1248 * coarse-grained (vsite) or atomistic interaction. Note that, beacuse
1249 * interactions between coarse-grained and other (atomistic) energygroups
1250 * are excluded automatically by grompp, it is sufficient to check for
1251 * the group id of atom i (igid) */
1252 bEnergyGroupCG = !egp_explicit(fr, igid);
1254 if (bDoVdW_i || bDoCoul_i)
1256 /* Loop over the j charge groups */
1257 for (j = 0; (j < nj); j++)
1261 /* Check for large charge groups */
1272 /* Finally loop over the atoms in the j-charge group */
1273 for (jj = jj0; jj < jj1; jj++)
1275 bNotEx = NOTEXCL(bExcl, i, jj);
1277 /* Now we have to exclude interactions which will be zero
1278 * anyway due to the AdResS weights (in previous implementations
1279 * this was done in the force kernel). This is necessary as
1280 * pure interactions (those with b_hybrid=false, i.e. w_i*w_j==1 or 0)
1281 * are put into neighbour lists which will be passed to the
1282 * standard (optimized) kernels for speed. The interactions with
1283 * b_hybrid=true are placed into the _adress neighbour lists and
1284 * processed by the generic AdResS kernel.
1286 if ( (bEnergyGroupCG &&
1287 wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS ) ||
1288 ( !bEnergyGroupCG && wf[jj] <= GMX_REAL_EPS ) )
1293 b_hybrid = !((wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS) ||
1294 (wf[i_atom] <= GMX_REAL_EPS && wf[jj] <= GMX_REAL_EPS));
1300 if (charge[jj] != 0)
1304 add_j_to_nblist(coul, jj, bLR);
1308 add_j_to_nblist(coul_adress, jj, bLR);
1312 else if (!bDoCoul_i)
1314 if (bHaveVdW[type[jj]])
1318 add_j_to_nblist(vdw, jj, bLR);
1322 add_j_to_nblist(vdw_adress, jj, bLR);
1328 if (bHaveVdW[type[jj]])
1330 if (charge[jj] != 0)
1334 add_j_to_nblist(vdwc, jj, bLR);
1338 add_j_to_nblist(vdwc_adress, jj, bLR);
1345 add_j_to_nblist(vdw, jj, bLR);
1349 add_j_to_nblist(vdw_adress, jj, bLR);
1354 else if (charge[jj] != 0)
1358 add_j_to_nblist(coul, jj, bLR);
1362 add_j_to_nblist(coul_adress, jj, bLR);
1371 close_i_nblist(vdw);
1372 close_i_nblist(coul);
1373 close_i_nblist(vdwc);
1374 close_i_nblist(vdw_adress);
1375 close_i_nblist(coul_adress);
1376 close_i_nblist(vdwc_adress);
1382 put_in_list_qmmm(gmx_bool gmx_unused bHaveVdW[],
1384 t_mdatoms gmx_unused * md,
1394 gmx_bool gmx_unused bDoVdW,
1395 gmx_bool gmx_unused bDoCoul,
1396 int gmx_unused solvent_opt)
1399 int i, j, jcg, igid, gid;
1400 atom_id jj, jj0, jj1, i_atom;
1404 /* Get atom range */
1406 nicg = index[icg+1]-i0;
1408 /* Get the i charge group info */
1409 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1411 coul = &fr->QMMMlist;
1413 /* Loop over atoms in the ith charge group */
1414 for (i = 0; i < nicg; i++)
1417 gid = GID(igid, jgid, ngid);
1418 /* Create new i_atom for each energy group */
1419 new_i_nblist(coul, i_atom, shift, gid);
1421 /* Loop over the j charge groups */
1422 for (j = 0; j < nj; j++)
1426 /* Charge groups cannot have QM and MM atoms simultaneously */
1431 /* Finally loop over the atoms in the j-charge group */
1432 for (jj = jj0; jj < jj1; jj++)
1434 bNotEx = NOTEXCL(bExcl, i, jj);
1437 add_j_to_nblist(coul, jj, bLR);
1442 close_i_nblist(coul);
1447 put_in_list_cg(gmx_bool gmx_unused bHaveVdW[],
1449 t_mdatoms gmx_unused * md,
1459 gmx_bool gmx_unused bDoVdW,
1460 gmx_bool gmx_unused bDoCoul,
1461 int gmx_unused solvent_opt)
1464 int igid, gid, nbl_ind;
1468 cginfo = fr->cginfo[icg];
1470 igid = GET_CGINFO_GID(cginfo);
1471 gid = GID(igid, jgid, ngid);
1473 /* Unpack pointers to neighbourlist structs */
1474 if (fr->nnblists == 1)
1480 nbl_ind = fr->gid2nblists[gid];
1484 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1488 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1491 /* Make a new neighbor list for charge group icg.
1492 * Currently simply one neighbor list is made with LJ and Coulomb.
1493 * If required, zero interactions could be removed here
1494 * or in the force loop.
1496 new_i_nblist(vdwc, index[icg], shift, gid);
1497 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1499 for (j = 0; (j < nj); j++)
1502 /* Skip the icg-icg pairs if all self interactions are excluded */
1503 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1505 /* Here we add the j charge group jcg to the list,
1506 * exclusions are also added to the list.
1508 add_j_to_nblist_cg(vdwc, index[jcg], index[jcg+1], bExcl, icg == jcg, bLR);
1512 close_i_nblist(vdwc);
1515 static void setexcl(atom_id start, atom_id end, t_blocka *excl, gmx_bool b,
1522 for (i = start; i < end; i++)
1524 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1526 SETEXCL(bexcl, i-start, excl->a[k]);
1532 for (i = start; i < end; i++)
1534 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1536 RMEXCL(bexcl, i-start, excl->a[k]);
1542 int calc_naaj(int icg, int cgtot)
1546 if ((cgtot % 2) == 1)
1548 /* Odd number of charge groups, easy */
1549 naaj = 1 + (cgtot/2);
1551 else if ((cgtot % 4) == 0)
1553 /* Multiple of four is hard */
1590 fprintf(log, "naaj=%d\n", naaj);
1596 /************************************************
1598 * S I M P L E C O R E S T U F F
1600 ************************************************/
1602 static real calc_image_tric(rvec xi, rvec xj, matrix box,
1603 rvec b_inv, int *shift)
1605 /* This code assumes that the cut-off is smaller than
1606 * a half times the smallest diagonal element of the box.
1608 const real h25 = 2.5;
1613 /* Compute diff vector */
1614 dz = xj[ZZ] - xi[ZZ];
1615 dy = xj[YY] - xi[YY];
1616 dx = xj[XX] - xi[XX];
1618 /* Perform NINT operation, using trunc operation, therefore
1619 * we first add 2.5 then subtract 2 again
1621 tz = dz*b_inv[ZZ] + h25;
1623 dz -= tz*box[ZZ][ZZ];
1624 dy -= tz*box[ZZ][YY];
1625 dx -= tz*box[ZZ][XX];
1627 ty = dy*b_inv[YY] + h25;
1629 dy -= ty*box[YY][YY];
1630 dx -= ty*box[YY][XX];
1632 tx = dx*b_inv[XX]+h25;
1634 dx -= tx*box[XX][XX];
1636 /* Distance squared */
1637 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1639 *shift = XYZ2IS(tx, ty, tz);
1644 static real calc_image_rect(rvec xi, rvec xj, rvec box_size,
1645 rvec b_inv, int *shift)
1647 const real h15 = 1.5;
1653 /* Compute diff vector */
1654 dx = xj[XX] - xi[XX];
1655 dy = xj[YY] - xi[YY];
1656 dz = xj[ZZ] - xi[ZZ];
1658 /* Perform NINT operation, using trunc operation, therefore
1659 * we first add 1.5 then subtract 1 again
1661 tx = dx*b_inv[XX] + h15;
1662 ty = dy*b_inv[YY] + h15;
1663 tz = dz*b_inv[ZZ] + h15;
1668 /* Correct diff vector for translation */
1669 ddx = tx*box_size[XX] - dx;
1670 ddy = ty*box_size[YY] - dy;
1671 ddz = tz*box_size[ZZ] - dz;
1673 /* Distance squared */
1674 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1676 *shift = XYZ2IS(tx, ty, tz);
1681 static void add_simple(t_ns_buf *nsbuf, int nrj, atom_id cg_j,
1682 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1683 int icg, int jgid, t_block *cgs, t_excl bexcl[],
1684 int shift, t_forcerec *fr, put_in_list_t *put_in_list)
1686 if (nsbuf->nj + nrj > MAX_CG)
1688 put_in_list(bHaveVdW, ngid, md, icg, jgid, nsbuf->ncg, nsbuf->jcg,
1689 cgs->index, bexcl, shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1690 /* Reset buffer contents */
1691 nsbuf->ncg = nsbuf->nj = 0;
1693 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1697 static void ns_inner_tric(rvec x[], int icg, int *i_egp_flags,
1698 int njcg, atom_id jcg[],
1699 matrix box, rvec b_inv, real rcut2,
1700 t_block *cgs, t_ns_buf **ns_buf,
1701 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1702 t_excl bexcl[], t_forcerec *fr,
1703 put_in_list_t *put_in_list)
1707 int *cginfo = fr->cginfo;
1708 atom_id cg_j, *cgindex;
1711 cgindex = cgs->index;
1713 for (j = 0; (j < njcg); j++)
1716 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1717 if (calc_image_tric(x[icg], x[cg_j], box, b_inv, &shift) < rcut2)
1719 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1720 if (!(i_egp_flags[jgid] & EGP_EXCL))
1722 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1723 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1730 static void ns_inner_rect(rvec x[], int icg, int *i_egp_flags,
1731 int njcg, atom_id jcg[],
1732 gmx_bool bBox, rvec box_size, rvec b_inv, real rcut2,
1733 t_block *cgs, t_ns_buf **ns_buf,
1734 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1735 t_excl bexcl[], t_forcerec *fr,
1736 put_in_list_t *put_in_list)
1740 int *cginfo = fr->cginfo;
1741 atom_id cg_j, *cgindex;
1744 cgindex = cgs->index;
1748 for (j = 0; (j < njcg); j++)
1751 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1752 if (calc_image_rect(x[icg], x[cg_j], box_size, b_inv, &shift) < rcut2)
1754 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1755 if (!(i_egp_flags[jgid] & EGP_EXCL))
1757 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1758 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1766 for (j = 0; (j < njcg); j++)
1769 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1770 if ((rcut2 == 0) || (distance2(x[icg], x[cg_j]) < rcut2))
1772 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1773 if (!(i_egp_flags[jgid] & EGP_EXCL))
1775 add_simple(&ns_buf[jgid][CENTRAL], nrj, cg_j,
1776 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, CENTRAL, fr,
1784 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1786 static int ns_simple_core(t_forcerec *fr,
1787 gmx_localtop_t *top,
1789 matrix box, rvec box_size,
1790 t_excl bexcl[], atom_id *aaj,
1791 int ngid, t_ns_buf **ns_buf,
1792 put_in_list_t *put_in_list, gmx_bool bHaveVdW[])
1796 int nsearch, icg, jcg, igid, i0, nri, nn;
1799 /* atom_id *i_atoms; */
1800 t_block *cgs = &(top->cgs);
1801 t_blocka *excl = &(top->excls);
1804 gmx_bool bBox, bTriclinic;
1807 rlist2 = sqr(fr->rlist);
1809 bBox = (fr->ePBC != epbcNONE);
1812 for (m = 0; (m < DIM); m++)
1814 b_inv[m] = divide_err(1.0, box_size[m]);
1816 bTriclinic = TRICLINIC(box);
1823 cginfo = fr->cginfo;
1826 for (icg = fr->cg0; (icg < fr->hcg); icg++)
1829 i0 = cgs->index[icg];
1830 nri = cgs->index[icg+1]-i0;
1831 i_atoms = &(cgs->a[i0]);
1832 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1833 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1835 igid = GET_CGINFO_GID(cginfo[icg]);
1836 i_egp_flags = fr->egp_flags + ngid*igid;
1837 setexcl(cgs->index[icg], cgs->index[icg+1], excl, TRUE, bexcl);
1839 naaj = calc_naaj(icg, cgs->nr);
1842 ns_inner_tric(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1843 box, b_inv, rlist2, cgs, ns_buf,
1844 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1848 ns_inner_rect(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1849 bBox, box_size, b_inv, rlist2, cgs, ns_buf,
1850 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1854 for (nn = 0; (nn < ngid); nn++)
1856 for (k = 0; (k < SHIFTS); k++)
1858 nsbuf = &(ns_buf[nn][k]);
1861 put_in_list(bHaveVdW, ngid, md, icg, nn, nsbuf->ncg, nsbuf->jcg,
1862 cgs->index, bexcl, k, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1863 nsbuf->ncg = nsbuf->nj = 0;
1867 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1868 setexcl(cgs->index[icg], cgs->index[icg+1], excl, FALSE, bexcl);
1870 close_neighbor_lists(fr, FALSE);
1875 /************************************************
1877 * N S 5 G R I D S T U F F
1879 ************************************************/
1881 static gmx_inline void get_dx(int Nx, real gridx, real rc2, int xgi, real x,
1882 int *dx0, int *dx1, real *dcx2)
1910 for (i = xgi0; i >= 0; i--)
1912 dcx = (i+1)*gridx-x;
1921 for (i = xgi1; i < Nx; i++)
1934 static gmx_inline void get_dx_dd(int Nx, real gridx, real rc2, int xgi, real x,
1935 int ncpddc, int shift_min, int shift_max,
1936 int *g0, int *g1, real *dcx2)
1939 int g_min, g_max, shift_home;
1972 g_min = (shift_min == shift_home ? 0 : ncpddc);
1973 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1980 else if (shift_max < 0)
1995 /* Check one grid cell down */
1996 dcx = ((*g0 - 1) + 1)*gridx - x;
2008 /* Check one grid cell up */
2009 dcx = (*g1 + 1)*gridx - x;
2021 #define sqr(x) ((x)*(x))
2022 #define calc_dx2(XI, YI, ZI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
2023 #define calc_cyl_dx2(XI, YI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
2024 /****************************************************
2026 * F A S T N E I G H B O R S E A R C H I N G
2028 * Optimized neighboursearching routine using grid
2029 * at least 1x1x1, see GROMACS manual
2031 ****************************************************/
2034 static void get_cutoff2(t_forcerec *fr, gmx_bool bDoLongRange,
2035 real *rvdw2, real *rcoul2,
2036 real *rs2, real *rm2, real *rl2)
2038 *rs2 = sqr(fr->rlist);
2040 if (bDoLongRange && fr->bTwinRange)
2042 /* With plain cut-off or RF we need to make the list exactly
2043 * up to the cut-off and the cut-off's can be different,
2044 * so we can not simply set them to rlistlong.
2045 * To keep this code compatible with (exotic) old cases,
2046 * we also create lists up to rvdw/rcoulomb for PME and Ewald.
2047 * The interaction check should correspond to:
2048 * !ir_vdw/coulomb_might_be_zero_at_cutoff from inputrec.c.
2050 if (((fr->vdwtype == evdwCUT || fr->vdwtype == evdwPME) &&
2051 fr->vdw_modifier == eintmodNONE) ||
2052 fr->rvdw <= fr->rlist)
2054 *rvdw2 = sqr(fr->rvdw);
2058 *rvdw2 = sqr(fr->rlistlong);
2060 if (((fr->eeltype == eelCUT ||
2061 (EEL_RF(fr->eeltype) && fr->eeltype != eelRF_ZERO) ||
2062 fr->eeltype == eelPME ||
2063 fr->eeltype == eelEWALD) &&
2064 fr->coulomb_modifier == eintmodNONE) ||
2065 fr->rcoulomb <= fr->rlist)
2067 *rcoul2 = sqr(fr->rcoulomb);
2071 *rcoul2 = sqr(fr->rlistlong);
2076 /* Workaround for a gcc -O3 or -ffast-math problem */
2080 *rm2 = min(*rvdw2, *rcoul2);
2081 *rl2 = max(*rvdw2, *rcoul2);
2084 static void init_nsgrid_lists(t_forcerec *fr, int ngid, gmx_ns_t *ns)
2086 real rvdw2, rcoul2, rs2, rm2, rl2;
2089 get_cutoff2(fr, TRUE, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2091 /* Short range buffers */
2092 snew(ns->nl_sr, ngid);
2094 snew(ns->nsr, ngid);
2095 snew(ns->nlr_ljc, ngid);
2096 snew(ns->nlr_one, ngid);
2098 /* Always allocate both list types, since rcoulomb might now change with PME load balancing */
2099 /* Long range VdW and Coul buffers */
2100 snew(ns->nl_lr_ljc, ngid);
2101 /* Long range VdW or Coul only buffers */
2102 snew(ns->nl_lr_one, ngid);
2104 for (j = 0; (j < ngid); j++)
2106 snew(ns->nl_sr[j], MAX_CG);
2107 snew(ns->nl_lr_ljc[j], MAX_CG);
2108 snew(ns->nl_lr_one[j], MAX_CG);
2113 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
2118 static int nsgrid_core(t_commrec *cr, t_forcerec *fr,
2119 matrix box, int ngid,
2120 gmx_localtop_t *top,
2122 t_excl bexcl[], gmx_bool *bExcludeAlleg,
2124 put_in_list_t *put_in_list,
2125 gmx_bool bHaveVdW[],
2126 gmx_bool bDoLongRange, gmx_bool bMakeQMMMnblist)
2129 atom_id **nl_lr_ljc, **nl_lr_one, **nl_sr;
2130 int *nlr_ljc, *nlr_one, *nsr;
2131 gmx_domdec_t *dd = NULL;
2132 t_block *cgs = &(top->cgs);
2133 int *cginfo = fr->cginfo;
2134 /* atom_id *i_atoms,*cgsindex=cgs->index; */
2136 int cell_x, cell_y, cell_z;
2137 int d, tx, ty, tz, dx, dy, dz, cj;
2138 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2139 int zsh_ty, zsh_tx, ysh_tx;
2141 int dx0, dx1, dy0, dy1, dz0, dz1;
2142 int Nx, Ny, Nz, shift = -1, j, nrj, nns, nn = -1;
2143 real gridx, gridy, gridz, grid_x, grid_y, grid_z;
2144 real *dcx2, *dcy2, *dcz2;
2146 int cg0, cg1, icg = -1, cgsnr, i0, igid, nri, naaj, max_jcg;
2147 int jcg0, jcg1, jjcg, cgj0, jgid;
2148 int *grida, *gridnra, *gridind;
2149 gmx_bool rvdw_lt_rcoul, rcoul_lt_rvdw;
2150 rvec xi, *cgcm, grid_offset;
2151 real r2, rs2, rvdw2, rcoul2, rm2, rl2, XI, YI, ZI, dcx, dcy, dcz, tmp1, tmp2;
2153 gmx_bool bDomDec, bTriclinicX, bTriclinicY;
2158 bDomDec = DOMAINDECOMP(cr);
2164 bTriclinicX = ((YY < grid->npbcdim &&
2165 (!bDomDec || dd->nc[YY] == 1) && box[YY][XX] != 0) ||
2166 (ZZ < grid->npbcdim &&
2167 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][XX] != 0));
2168 bTriclinicY = (ZZ < grid->npbcdim &&
2169 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][YY] != 0);
2173 get_cutoff2(fr, bDoLongRange, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2175 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
2176 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
2178 if (bMakeQMMMnblist)
2186 nl_lr_ljc = ns->nl_lr_ljc;
2187 nl_lr_one = ns->nl_lr_one;
2188 nlr_ljc = ns->nlr_ljc;
2189 nlr_one = ns->nlr_one;
2197 gridind = grid->index;
2198 gridnra = grid->nra;
2201 gridx = grid->cell_size[XX];
2202 gridy = grid->cell_size[YY];
2203 gridz = grid->cell_size[ZZ];
2207 copy_rvec(grid->cell_offset, grid_offset);
2208 copy_ivec(grid->ncpddc, ncpddc);
2213 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2214 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
2215 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
2216 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
2217 if (zsh_tx != 0 && ysh_tx != 0)
2219 /* This could happen due to rounding, when both ratios are 0.5 */
2228 /* We only want a list for the test particle */
2237 /* Set the shift range */
2238 for (d = 0; d < DIM; d++)
2242 /* Check if we need periodicity shifts.
2243 * Without PBC or with domain decomposition we don't need them.
2245 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2252 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
2263 /* Loop over charge groups */
2264 for (icg = cg0; (icg < cg1); icg++)
2266 igid = GET_CGINFO_GID(cginfo[icg]);
2267 /* Skip this charge group if all energy groups are excluded! */
2268 if (bExcludeAlleg[igid])
2273 i0 = cgs->index[icg];
2275 if (bMakeQMMMnblist)
2277 /* Skip this charge group if it is not a QM atom while making a
2278 * QM/MM neighbourlist
2280 if (md->bQM[i0] == FALSE)
2282 continue; /* MM particle, go to next particle */
2285 /* Compute the number of charge groups that fall within the control
2288 naaj = calc_naaj(icg, cgsnr);
2295 /* make a normal neighbourlist */
2299 /* Get the j charge-group and dd cell shift ranges */
2300 dd_get_ns_ranges(cr->dd, icg, &jcg0, &jcg1, sh0, sh1);
2305 /* Compute the number of charge groups that fall within the control
2308 naaj = calc_naaj(icg, cgsnr);
2314 /* The i-particle is awlways the test particle,
2315 * so we want all j-particles
2317 max_jcg = cgsnr - 1;
2321 max_jcg = jcg1 - cgsnr;
2326 i_egp_flags = fr->egp_flags + igid*ngid;
2328 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2329 setexcl(i0, cgs->index[icg+1], &top->excls, TRUE, bexcl);
2331 ci2xyz(grid, icg, &cell_x, &cell_y, &cell_z);
2333 /* Changed iicg to icg, DvdS 990115
2334 * (but see consistency check above, DvdS 990330)
2337 fprintf(log, "icg=%5d, naaj=%5d, cell %d %d %d\n",
2338 icg, naaj, cell_x, cell_y, cell_z);
2340 /* Loop over shift vectors in three dimensions */
2341 for (tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
2343 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2344 /* Calculate range of cells in Z direction that have the shift tz */
2345 zgi = cell_z + tz*Nz;
2348 get_dx(Nz, gridz, rl2, zgi, ZI, &dz0, &dz1, dcz2);
2350 get_dx_dd(Nz, gridz, rl2, zgi, ZI-grid_offset[ZZ],
2351 ncpddc[ZZ], sh0[ZZ], sh1[ZZ], &dz0, &dz1, dcz2);
2357 for (ty = -shp[YY]; ty <= shp[YY]; ty++)
2359 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2360 /* Calculate range of cells in Y direction that have the shift ty */
2363 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2367 ygi = cell_y + ty*Ny;
2370 get_dx(Ny, gridy, rl2, ygi, YI, &dy0, &dy1, dcy2);
2372 get_dx_dd(Ny, gridy, rl2, ygi, YI-grid_offset[YY],
2373 ncpddc[YY], sh0[YY], sh1[YY], &dy0, &dy1, dcy2);
2379 for (tx = -shp[XX]; tx <= shp[XX]; tx++)
2381 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2382 /* Calculate range of cells in X direction that have the shift tx */
2385 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2389 xgi = cell_x + tx*Nx;
2392 get_dx(Nx, gridx, rl2, xgi*Nx, XI, &dx0, &dx1, dcx2);
2394 get_dx_dd(Nx, gridx, rl2, xgi, XI-grid_offset[XX],
2395 ncpddc[XX], sh0[XX], sh1[XX], &dx0, &dx1, dcx2);
2401 /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2402 * from the neigbour list as it will not interact */
2403 if (fr->adress_type != eAdressOff)
2405 if (md->wf[cgs->index[icg]] <= GMX_REAL_EPS && egp_explicit(fr, igid))
2410 /* Get shift vector */
2411 shift = XYZ2IS(tx, ty, tz);
2413 range_check(shift, 0, SHIFTS);
2415 for (nn = 0; (nn < ngid); nn++)
2422 fprintf(log, "shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2423 shift, dx0, dx1, dy0, dy1, dz0, dz1);
2424 fprintf(log, "cgcm: %8.3f %8.3f %8.3f\n", cgcm[icg][XX],
2425 cgcm[icg][YY], cgcm[icg][ZZ]);
2426 fprintf(log, "xi: %8.3f %8.3f %8.3f\n", XI, YI, ZI);
2428 for (dx = dx0; (dx <= dx1); dx++)
2430 tmp1 = rl2 - dcx2[dx];
2431 for (dy = dy0; (dy <= dy1); dy++)
2433 tmp2 = tmp1 - dcy2[dy];
2436 for (dz = dz0; (dz <= dz1); dz++)
2438 if (tmp2 > dcz2[dz])
2440 /* Find grid-cell cj in which possible neighbours are */
2441 cj = xyz2ci(Ny, Nz, dx, dy, dz);
2443 /* Check out how many cgs (nrj) there in this cell */
2446 /* Find the offset in the cg list */
2449 /* Check if all j's are out of range so we
2450 * can skip the whole cell.
2451 * Should save some time, especially with DD.
2454 (grida[cgj0] >= max_jcg &&
2455 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2461 for (j = 0; (j < nrj); j++)
2463 jjcg = grida[cgj0+j];
2465 /* check whether this guy is in range! */
2466 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2469 r2 = calc_dx2(XI, YI, ZI, cgcm[jjcg]);
2472 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2473 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2474 /* check energy group exclusions */
2475 if (!(i_egp_flags[jgid] & EGP_EXCL))
2479 if (nsr[jgid] >= MAX_CG)
2481 /* Add to short-range list */
2482 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2483 nsr[jgid], nl_sr[jgid],
2484 cgs->index, /* cgsatoms, */ bexcl,
2485 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2488 nl_sr[jgid][nsr[jgid]++] = jjcg;
2492 if (nlr_ljc[jgid] >= MAX_CG)
2494 /* Add to LJ+coulomb long-range list */
2495 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2496 nlr_ljc[jgid], nl_lr_ljc[jgid], top->cgs.index,
2497 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2500 nl_lr_ljc[jgid][nlr_ljc[jgid]++] = jjcg;
2504 if (nlr_one[jgid] >= MAX_CG)
2506 /* Add to long-range list with only coul, or only LJ */
2507 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2508 nlr_one[jgid], nl_lr_one[jgid], top->cgs.index,
2509 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2512 nl_lr_one[jgid][nlr_one[jgid]++] = jjcg;
2524 /* CHECK whether there is anything left in the buffers */
2525 for (nn = 0; (nn < ngid); nn++)
2529 put_in_list(bHaveVdW, ngid, md, icg, nn, nsr[nn], nl_sr[nn],
2530 cgs->index, /* cgsatoms, */ bexcl,
2531 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2534 if (nlr_ljc[nn] > 0)
2536 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_ljc[nn],
2537 nl_lr_ljc[nn], top->cgs.index,
2538 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2541 if (nlr_one[nn] > 0)
2543 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_one[nn],
2544 nl_lr_one[nn], top->cgs.index,
2545 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2551 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2552 setexcl(cgs->index[icg], cgs->index[icg+1], &top->excls, FALSE, bexcl);
2554 /* No need to perform any left-over force calculations anymore (as we used to do here)
2555 * since we now save the proper long-range lists for later evaluation.
2560 /* Close neighbourlists */
2561 close_neighbor_lists(fr, bMakeQMMMnblist);
2566 void ns_realloc_natoms(gmx_ns_t *ns, int natoms)
2570 if (natoms > ns->nra_alloc)
2572 ns->nra_alloc = over_alloc_dd(natoms);
2573 srenew(ns->bexcl, ns->nra_alloc);
2574 for (i = 0; i < ns->nra_alloc; i++)
2581 void init_ns(FILE *fplog, const t_commrec *cr,
2582 gmx_ns_t *ns, t_forcerec *fr,
2583 const gmx_mtop_t *mtop)
2585 int mt, icg, nr_in_cg, maxcg, i, j, jcg, ngid, ncg;
2589 /* Compute largest charge groups size (# atoms) */
2591 for (mt = 0; mt < mtop->nmoltype; mt++)
2593 cgs = &mtop->moltype[mt].cgs;
2594 for (icg = 0; (icg < cgs->nr); icg++)
2596 nr_in_cg = max(nr_in_cg, (int)(cgs->index[icg+1]-cgs->index[icg]));
2600 /* Verify whether largest charge group is <= max cg.
2601 * This is determined by the type of the local exclusion type
2602 * Exclusions are stored in bits. (If the type is not large
2603 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2605 maxcg = sizeof(t_excl)*8;
2606 if (nr_in_cg > maxcg)
2608 gmx_fatal(FARGS, "Max #atoms in a charge group: %d > %d\n",
2612 ngid = mtop->groups.grps[egcENER].nr;
2613 snew(ns->bExcludeAlleg, ngid);
2614 for (i = 0; i < ngid; i++)
2616 ns->bExcludeAlleg[i] = TRUE;
2617 for (j = 0; j < ngid; j++)
2619 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2621 ns->bExcludeAlleg[i] = FALSE;
2629 ns->grid = init_grid(fplog, fr);
2630 init_nsgrid_lists(fr, ngid, ns);
2635 snew(ns->ns_buf, ngid);
2636 for (i = 0; (i < ngid); i++)
2638 snew(ns->ns_buf[i], SHIFTS);
2640 ncg = ncg_mtop(mtop);
2641 snew(ns->simple_aaj, 2*ncg);
2642 for (jcg = 0; (jcg < ncg); jcg++)
2644 ns->simple_aaj[jcg] = jcg;
2645 ns->simple_aaj[jcg+ncg] = jcg;
2649 /* Create array that determines whether or not atoms have VdW */
2650 snew(ns->bHaveVdW, fr->ntype);
2651 for (i = 0; (i < fr->ntype); i++)
2653 for (j = 0; (j < fr->ntype); j++)
2655 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2657 ((BHAMA(fr->nbfp, fr->ntype, i, j) != 0) ||
2658 (BHAMB(fr->nbfp, fr->ntype, i, j) != 0) ||
2659 (BHAMC(fr->nbfp, fr->ntype, i, j) != 0)) :
2660 ((C6(fr->nbfp, fr->ntype, i, j) != 0) ||
2661 (C12(fr->nbfp, fr->ntype, i, j) != 0))));
2666 pr_bvec(debug, 0, "bHaveVdW", ns->bHaveVdW, fr->ntype, TRUE);
2671 if (!DOMAINDECOMP(cr))
2673 ns_realloc_natoms(ns, mtop->natoms);
2676 ns->nblist_initialized = FALSE;
2678 /* nbr list debug dump */
2680 char *ptr = getenv("GMX_DUMP_NL");
2683 ns->dump_nl = strtol(ptr, NULL, 10);
2686 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2697 int search_neighbours(FILE *log, t_forcerec *fr,
2699 gmx_localtop_t *top,
2700 gmx_groups_t *groups,
2702 t_nrnb *nrnb, t_mdatoms *md,
2704 gmx_bool bDoLongRangeNS)
2706 t_block *cgs = &(top->cgs);
2707 rvec box_size, grid_x0, grid_x1;
2709 real min_size, grid_dens;
2713 gmx_bool *i_egp_flags;
2714 int cg_start, cg_end, start, end;
2717 gmx_domdec_zones_t *dd_zones;
2718 put_in_list_t *put_in_list;
2722 /* Set some local variables */
2724 ngid = groups->grps[egcENER].nr;
2726 for (m = 0; (m < DIM); m++)
2728 box_size[m] = box[m][m];
2731 if (fr->ePBC != epbcNONE)
2733 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC, box))
2735 gmx_fatal(FARGS, "One of the box vectors has become shorter than twice the cut-off length or box_yy-|box_zy| or box_zz has become smaller than the cut-off.");
2739 min_size = min(box_size[XX], min(box_size[YY], box_size[ZZ]));
2740 if (2*fr->rlistlong >= min_size)
2742 gmx_fatal(FARGS, "One of the box diagonal elements has become smaller than twice the cut-off length.");
2747 if (DOMAINDECOMP(cr))
2749 ns_realloc_natoms(ns, cgs->index[cgs->nr]);
2753 /* Reset the neighbourlists */
2754 reset_neighbor_lists(fr, TRUE, TRUE);
2756 if (bGrid && bFillGrid)
2760 if (DOMAINDECOMP(cr))
2762 dd_zones = domdec_zones(cr->dd);
2768 get_nsgrid_boundaries(grid->nboundeddim, box, NULL, NULL, NULL, NULL,
2769 cgs->nr, fr->cg_cm, grid_x0, grid_x1, &grid_dens);
2771 grid_first(log, grid, NULL, NULL, box, grid_x0, grid_x1,
2772 fr->rlistlong, grid_dens);
2779 if (DOMAINDECOMP(cr))
2782 fill_grid(dd_zones, grid, end, -1, end, fr->cg_cm);
2784 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2788 fill_grid(NULL, grid, cgs->nr, fr->cg0, fr->hcg, fr->cg_cm);
2789 grid->icg0 = fr->cg0;
2790 grid->icg1 = fr->hcg;
2794 calc_elemnr(grid, start, end, cgs->nr);
2796 grid_last(grid, start, end, cgs->nr);
2801 print_grid(debug, grid);
2806 /* Set the grid cell index for the test particle only.
2807 * The cell to cg index is not corrected, but that does not matter.
2809 fill_grid(NULL, ns->grid, fr->hcg, fr->hcg-1, fr->hcg, fr->cg_cm);
2813 if (fr->adress_type == eAdressOff)
2815 if (!fr->ns.bCGlist)
2817 put_in_list = put_in_list_at;
2821 put_in_list = put_in_list_cg;
2826 put_in_list = put_in_list_adress;
2833 nsearch = nsgrid_core(cr, fr, box, ngid, top,
2834 grid, ns->bexcl, ns->bExcludeAlleg,
2835 md, put_in_list, ns->bHaveVdW,
2836 bDoLongRangeNS, FALSE);
2838 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2839 * the classical calculation. The charge-charge interaction
2840 * between QM and MM atoms is handled in the QMMM core calculation
2841 * (see QMMM.c). The VDW however, we'd like to compute classically
2842 * and the QM MM atom pairs have just been put in the
2843 * corresponding neighbourlists. in case of QMMM we still need to
2844 * fill a special QMMM neighbourlist that contains all neighbours
2845 * of the QM atoms. If bQMMM is true, this list will now be made:
2847 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
2849 nsearch += nsgrid_core(cr, fr, box, ngid, top,
2850 grid, ns->bexcl, ns->bExcludeAlleg,
2851 md, put_in_list_qmmm, ns->bHaveVdW,
2852 bDoLongRangeNS, TRUE);
2857 nsearch = ns_simple_core(fr, top, md, box, box_size,
2858 ns->bexcl, ns->simple_aaj,
2859 ngid, ns->ns_buf, put_in_list, ns->bHaveVdW);
2867 inc_nrnb(nrnb, eNR_NS, nsearch);
2868 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2873 int natoms_beyond_ns_buffer(t_inputrec *ir, t_forcerec *fr, t_block *cgs,
2874 matrix scale_tot, rvec *x)
2876 int cg0, cg1, cg, a0, a1, a, i, j;
2877 real rint, hbuf2, scale;
2879 gmx_bool bIsotropic;
2884 rint = max(ir->rcoulomb, ir->rvdw);
2885 if (ir->rlist < rint)
2887 gmx_fatal(FARGS, "The neighbor search buffer has negative size: %f nm",
2895 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2897 hbuf2 = sqr(0.5*(ir->rlist - rint));
2898 for (cg = cg0; cg < cg1; cg++)
2900 a0 = cgs->index[cg];
2901 a1 = cgs->index[cg+1];
2902 for (a = a0; a < a1; a++)
2904 if (distance2(cg_cm[cg], x[a]) > hbuf2)
2914 scale = scale_tot[0][0];
2915 for (i = 1; i < DIM; i++)
2917 /* With anisotropic scaling, the original spherical ns volumes become
2918 * ellipsoids. To avoid costly transformations we use the minimum
2919 * eigenvalue of the scaling matrix for determining the buffer size.
2920 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2922 scale = min(scale, scale_tot[i][i]);
2923 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2927 for (j = 0; j < i; j++)
2929 if (scale_tot[i][j] != 0)
2935 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2938 for (cg = cg0; cg < cg1; cg++)
2940 svmul(scale, cg_cm[cg], cgsc);
2941 a0 = cgs->index[cg];
2942 a1 = cgs->index[cg+1];
2943 for (a = a0; a < a1; a++)
2945 if (distance2(cgsc, x[a]) > hbuf2)
2954 /* Anistropic scaling */
2955 for (cg = cg0; cg < cg1; cg++)
2957 /* Since scale_tot contains the transpose of the scaling matrix,
2958 * we need to multiply with the transpose.
2960 tmvmul_ur0(scale_tot, cg_cm[cg], cgsc);
2961 a0 = cgs->index[cg];
2962 a1 = cgs->index[cg+1];
2963 for (a = a0; a < a1; a++)
2965 if (distance2(cgsc, x[a]) > hbuf2)