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.
43 #include "gromacs/legacyheaders/macros.h"
44 #include "gromacs/math/utilities.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/legacyheaders/types/commrec.h"
47 #include "gromacs/legacyheaders/network.h"
48 #include "gromacs/legacyheaders/nsgrid.h"
49 #include "gromacs/legacyheaders/force.h"
50 #include "gromacs/legacyheaders/nonbonded.h"
51 #include "gromacs/legacyheaders/ns.h"
52 #include "gromacs/legacyheaders/names.h"
53 #include "gromacs/legacyheaders/nrnb.h"
54 #include "gromacs/legacyheaders/txtdump.h"
55 #include "gromacs/topology/mtop_util.h"
57 #include "gromacs/legacyheaders/domdec.h"
60 #include "gromacs/pbcutil/ishift.h"
61 #include "gromacs/pbcutil/pbc.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/smalloc.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,
131 gmx_bool bElecAndVdwSwitchDiffers)
137 for (i = 0; (i < 2); i++)
139 nl = (i == 0) ? nl_sr : nl_lr;
140 homenr = (i == 0) ? maxsr : maxlr;
148 /* Set coul/vdw in neighborlist, and for the normal loops we determine
149 * an index of which one to call.
152 nl->ivdwmod = ivdwmod;
154 nl->ielecmod = ielecmod;
156 nl->igeometry = igeometry;
158 if (nl->type == GMX_NBLIST_INTERACTION_FREE_ENERGY)
160 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
163 /* This will also set the simd_padding_width field */
164 gmx_nonbonded_set_kernel_pointers( (i == 0) ? log : NULL, nl, bElecAndVdwSwitchDiffers);
166 /* maxnri is influenced by the number of shifts (maximum is 8)
167 * and the number of energy groups.
168 * If it is not enough, nl memory will be reallocated during the run.
169 * 4 seems to be a reasonable factor, which only causes reallocation
170 * during runs with tiny and many energygroups.
172 nl->maxnri = homenr*4;
182 reallocate_nblist(nl);
187 fprintf(debug, "Initiating neighbourlist (ielec=%d, ivdw=%d, type=%d) for %s interactions,\nwith %d SR, %d LR atoms.\n",
188 nl->ielec, nl->ivdw, nl->type, gmx_nblist_geometry_names[nl->igeometry], maxsr, maxlr);
193 void init_neighbor_list(FILE *log, t_forcerec *fr, int homenr)
195 /* Make maxlr tunable! (does not seem to be a big difference though)
196 * This parameter determines the number of i particles in a long range
197 * neighbourlist. Too few means many function calls, too many means
200 int maxsr, maxsr_wat, maxlr, maxlr_wat;
201 int ielec, ivdw, ielecmod, ivdwmod, type;
203 int igeometry_def, igeometry_w, igeometry_ww;
205 gmx_bool bElecAndVdwSwitchDiffers;
208 /* maxsr = homenr-fr->nWatMol*3; */
213 gmx_fatal(FARGS, "%s, %d: Negative number of short range atoms.\n"
214 "Call your Gromacs dealer for assistance.", __FILE__, __LINE__);
216 /* This is just for initial allocation, so we do not reallocate
217 * all the nlist arrays many times in a row.
218 * The numbers seem very accurate, but they are uncritical.
220 maxsr_wat = min(fr->nWatMol, (homenr+2)/3);
224 maxlr_wat = min(maxsr_wat, maxlr);
228 maxlr = maxlr_wat = 0;
231 /* Determine the values for ielec/ivdw. */
232 ielec = fr->nbkernel_elec_interaction;
233 ivdw = fr->nbkernel_vdw_interaction;
234 ielecmod = fr->nbkernel_elec_modifier;
235 ivdwmod = fr->nbkernel_vdw_modifier;
236 type = GMX_NBLIST_INTERACTION_STANDARD;
237 bElecAndVdwSwitchDiffers = ( (fr->rcoulomb_switch != fr->rvdw_switch) || (fr->rcoulomb != fr->rvdw));
239 fr->ns.bCGlist = (getenv("GMX_NBLISTCG") != 0);
242 igeometry_def = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
246 igeometry_def = GMX_NBLIST_GEOMETRY_CG_CG;
249 fprintf(log, "\nUsing charge-group - charge-group neighbor lists and kernels\n\n");
253 if (fr->solvent_opt == esolTIP4P)
255 igeometry_w = GMX_NBLIST_GEOMETRY_WATER4_PARTICLE;
256 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER4_WATER4;
260 igeometry_w = GMX_NBLIST_GEOMETRY_WATER3_PARTICLE;
261 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER3_WATER3;
264 for (i = 0; i < fr->nnblists; i++)
266 nbl = &(fr->nblists[i]);
268 if ((fr->adress_type != eAdressOff) && (i >= fr->nnblists/2))
270 type = GMX_NBLIST_INTERACTION_ADRESS;
272 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ], &nbl->nlist_lr[eNL_VDWQQ],
273 maxsr, maxlr, ivdw, ivdwmod, ielec, ielecmod, igeometry_def, type, bElecAndVdwSwitchDiffers);
274 init_nblist(log, &nbl->nlist_sr[eNL_VDW], &nbl->nlist_lr[eNL_VDW],
275 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, igeometry_def, type, bElecAndVdwSwitchDiffers);
276 init_nblist(log, &nbl->nlist_sr[eNL_QQ], &nbl->nlist_lr[eNL_QQ],
277 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_def, type, bElecAndVdwSwitchDiffers);
278 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATER], &nbl->nlist_lr[eNL_VDWQQ_WATER],
279 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_w, type, bElecAndVdwSwitchDiffers);
280 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATER], &nbl->nlist_lr[eNL_QQ_WATER],
281 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_w, type, bElecAndVdwSwitchDiffers);
282 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATERWATER], &nbl->nlist_lr[eNL_VDWQQ_WATERWATER],
283 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_ww, type, bElecAndVdwSwitchDiffers);
284 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATERWATER], &nbl->nlist_lr[eNL_QQ_WATERWATER],
285 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_ww, type, bElecAndVdwSwitchDiffers);
287 /* Did we get the solvent loops so we can use optimized water kernels? */
288 if (nbl->nlist_sr[eNL_VDWQQ_WATER].kernelptr_vf == NULL
289 || nbl->nlist_sr[eNL_QQ_WATER].kernelptr_vf == NULL
290 || nbl->nlist_sr[eNL_VDWQQ_WATERWATER].kernelptr_vf == NULL
291 || nbl->nlist_sr[eNL_QQ_WATERWATER].kernelptr_vf == NULL)
293 fr->solvent_opt = esolNO;
296 fprintf(log, "Note: The available nonbonded kernels do not support water optimization - disabling.\n");
300 if (fr->efep != efepNO)
302 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_FREE], &nbl->nlist_lr[eNL_VDWQQ_FREE],
303 maxsr, maxlr, ivdw, ivdwmod, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY, bElecAndVdwSwitchDiffers);
304 init_nblist(log, &nbl->nlist_sr[eNL_VDW_FREE], &nbl->nlist_lr[eNL_VDW_FREE],
305 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY, bElecAndVdwSwitchDiffers);
306 init_nblist(log, &nbl->nlist_sr[eNL_QQ_FREE], &nbl->nlist_lr[eNL_QQ_FREE],
307 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY, bElecAndVdwSwitchDiffers);
311 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
313 init_nblist(log, &fr->QMMMlist, NULL,
314 maxsr, maxlr, 0, 0, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_STANDARD, bElecAndVdwSwitchDiffers);
322 fr->ns.nblist_initialized = TRUE;
325 static void reset_nblist(t_nblist *nl)
335 static void reset_neighbor_lists(t_forcerec *fr, gmx_bool bResetSR, gmx_bool bResetLR)
341 /* only reset the short-range nblist */
342 reset_nblist(&(fr->QMMMlist));
345 for (n = 0; n < fr->nnblists; n++)
347 for (i = 0; i < eNL_NR; i++)
351 reset_nblist( &(fr->nblists[n].nlist_sr[i]) );
355 reset_nblist( &(fr->nblists[n].nlist_lr[i]) );
364 static gmx_inline void new_i_nblist(t_nblist *nlist, atom_id i_atom, int shift, int gid)
366 int i, k, nri, nshift;
370 /* Check whether we have to increase the i counter */
372 (nlist->iinr[nri] != i_atom) ||
373 (nlist->shift[nri] != shift) ||
374 (nlist->gid[nri] != gid))
376 /* This is something else. Now see if any entries have
377 * been added in the list of the previous atom.
380 ((nlist->jindex[nri+1] > nlist->jindex[nri]) &&
381 (nlist->gid[nri] != -1)))
383 /* If so increase the counter */
386 if (nlist->nri >= nlist->maxnri)
388 nlist->maxnri += over_alloc_large(nlist->nri);
389 reallocate_nblist(nlist);
392 /* Set the number of neighbours and the atom number */
393 nlist->jindex[nri+1] = nlist->jindex[nri];
394 nlist->iinr[nri] = i_atom;
395 nlist->gid[nri] = gid;
396 nlist->shift[nri] = shift;
400 /* Adding to previous list. First remove possible previous padding */
401 if (nlist->simd_padding_width > 1)
403 while (nlist->nrj > 0 && nlist->jjnr[nlist->nrj-1] < 0)
411 static gmx_inline void close_i_nblist(t_nblist *nlist)
413 int nri = nlist->nri;
418 /* Add elements up to padding. Since we allocate memory in units
419 * of the simd_padding width, we do not have to check for possible
420 * list reallocation here.
422 while ((nlist->nrj % nlist->simd_padding_width) != 0)
424 /* Use -4 here, so we can write forces for 4 atoms before real data */
425 nlist->jjnr[nlist->nrj++] = -4;
427 nlist->jindex[nri+1] = nlist->nrj;
429 len = nlist->nrj - nlist->jindex[nri];
433 static gmx_inline void close_nblist(t_nblist *nlist)
435 /* Only close this nblist when it has been initialized.
436 * Avoid the creation of i-lists with no j-particles.
440 /* Some assembly kernels do not support empty lists,
441 * make sure here that we don't generate any empty lists.
442 * With the current ns code this branch is taken in two cases:
443 * No i-particles at all: nri=-1 here
444 * There are i-particles, but no j-particles; nri=0 here
450 /* Close list number nri by incrementing the count */
455 static gmx_inline void close_neighbor_lists(t_forcerec *fr, gmx_bool bMakeQMMMnblist)
461 close_nblist(&(fr->QMMMlist));
464 for (n = 0; n < fr->nnblists; n++)
466 for (i = 0; (i < eNL_NR); i++)
468 close_nblist(&(fr->nblists[n].nlist_sr[i]));
469 close_nblist(&(fr->nblists[n].nlist_lr[i]));
475 static gmx_inline void add_j_to_nblist(t_nblist *nlist, atom_id j_atom, gmx_bool bLR)
477 int nrj = nlist->nrj;
479 if (nlist->nrj >= nlist->maxnrj)
481 nlist->maxnrj = round_up_to_simd_width(over_alloc_small(nlist->nrj + 1), nlist->simd_padding_width);
485 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
486 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
489 srenew(nlist->jjnr, nlist->maxnrj);
492 nlist->jjnr[nrj] = j_atom;
496 static gmx_inline void add_j_to_nblist_cg(t_nblist *nlist,
497 atom_id j_start, int j_end,
498 t_excl *bexcl, gmx_bool i_is_j,
501 int nrj = nlist->nrj;
504 if (nlist->nrj >= nlist->maxnrj)
506 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
509 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
510 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
513 srenew(nlist->jjnr, nlist->maxnrj);
514 srenew(nlist->jjnr_end, nlist->maxnrj);
515 srenew(nlist->excl, nlist->maxnrj*MAX_CGCGSIZE);
518 nlist->jjnr[nrj] = j_start;
519 nlist->jjnr_end[nrj] = j_end;
521 if (j_end - j_start > MAX_CGCGSIZE)
523 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);
526 /* Set the exclusions */
527 for (j = j_start; j < j_end; j++)
529 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
533 /* Avoid double counting of intra-cg interactions */
534 for (j = 1; j < j_end-j_start; j++)
536 nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
544 put_in_list_t (gmx_bool bHaveVdW[],
561 put_in_list_at(gmx_bool bHaveVdW[],
577 /* The a[] index has been removed,
578 * to put it back in i_atom should be a[i0] and jj should be a[jj].
583 t_nblist * vdwc_free = NULL;
584 t_nblist * vdw_free = NULL;
585 t_nblist * coul_free = NULL;
586 t_nblist * vdwc_ww = NULL;
587 t_nblist * coul_ww = NULL;
589 int i, j, jcg, igid, gid, nbl_ind, ind_ij;
590 atom_id jj, jj0, jj1, i_atom;
595 real *charge, *chargeB;
596 real qi, qiB, qq, rlj;
597 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
598 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
602 /* Copy some pointers */
604 charge = md->chargeA;
605 chargeB = md->chargeB;
608 bPert = md->bPerturbed;
612 nicg = index[icg+1]-i0;
614 /* Get the i charge group info */
615 igid = GET_CGINFO_GID(cginfo[icg]);
617 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
622 /* Check if any of the particles involved are perturbed.
623 * If not we can do the cheaper normal put_in_list
624 * and use more solvent optimization.
626 for (i = 0; i < nicg; i++)
628 bFreeEnergy |= bPert[i0+i];
630 /* Loop over the j charge groups */
631 for (j = 0; (j < nj && !bFreeEnergy); j++)
636 /* Finally loop over the atoms in the j-charge group */
637 for (jj = jj0; jj < jj1; jj++)
639 bFreeEnergy |= bPert[jj];
644 /* Unpack pointers to neighbourlist structs */
645 if (fr->nnblists == 1)
651 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
655 nlist = fr->nblists[nbl_ind].nlist_lr;
659 nlist = fr->nblists[nbl_ind].nlist_sr;
662 if (iwater != esolNO)
664 vdwc = &nlist[eNL_VDWQQ_WATER];
665 vdw = &nlist[eNL_VDW];
666 coul = &nlist[eNL_QQ_WATER];
667 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
668 coul_ww = &nlist[eNL_QQ_WATERWATER];
672 vdwc = &nlist[eNL_VDWQQ];
673 vdw = &nlist[eNL_VDW];
674 coul = &nlist[eNL_QQ];
679 if (iwater != esolNO)
681 /* Loop over the atoms in the i charge group */
683 gid = GID(igid, jgid, ngid);
684 /* Create new i_atom for each energy group */
685 if (bDoCoul && bDoVdW)
687 new_i_nblist(vdwc, i_atom, shift, gid);
688 new_i_nblist(vdwc_ww, i_atom, shift, gid);
692 new_i_nblist(vdw, i_atom, shift, gid);
696 new_i_nblist(coul, i_atom, shift, gid);
697 new_i_nblist(coul_ww, i_atom, shift, gid);
699 /* Loop over the j charge groups */
700 for (j = 0; (j < nj); j++)
710 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
712 if (iwater == esolSPC && jwater == esolSPC)
714 /* Interaction between two SPC molecules */
717 /* VdW only - only first atoms in each water interact */
718 add_j_to_nblist(vdw, jj0, bLR);
722 /* One entry for the entire water-water interaction */
725 add_j_to_nblist(coul_ww, jj0, bLR);
729 add_j_to_nblist(vdwc_ww, jj0, bLR);
733 else if (iwater == esolTIP4P && jwater == esolTIP4P)
735 /* Interaction between two TIP4p molecules */
738 /* VdW only - only first atoms in each water interact */
739 add_j_to_nblist(vdw, jj0, bLR);
743 /* One entry for the entire water-water interaction */
746 add_j_to_nblist(coul_ww, jj0, bLR);
750 add_j_to_nblist(vdwc_ww, jj0, bLR);
756 /* j charge group is not water, but i is.
757 * Add entries to the water-other_atom lists; the geometry of the water
758 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
759 * so we don't care if it is SPC or TIP4P...
766 for (jj = jj0; (jj < jj1); jj++)
770 add_j_to_nblist(coul, jj, bLR);
776 for (jj = jj0; (jj < jj1); jj++)
778 if (bHaveVdW[type[jj]])
780 add_j_to_nblist(vdw, jj, bLR);
786 /* _charge_ _groups_ interact with both coulomb and LJ */
787 /* Check which atoms we should add to the lists! */
788 for (jj = jj0; (jj < jj1); jj++)
790 if (bHaveVdW[type[jj]])
794 add_j_to_nblist(vdwc, jj, bLR);
798 add_j_to_nblist(vdw, jj, bLR);
801 else if (charge[jj] != 0)
803 add_j_to_nblist(coul, jj, bLR);
810 close_i_nblist(coul);
811 close_i_nblist(vdwc);
812 close_i_nblist(coul_ww);
813 close_i_nblist(vdwc_ww);
817 /* no solvent as i charge group */
818 /* Loop over the atoms in the i charge group */
819 for (i = 0; i < nicg; i++)
822 gid = GID(igid, jgid, ngid);
825 /* Create new i_atom for each energy group */
826 if (bDoVdW && bDoCoul)
828 new_i_nblist(vdwc, i_atom, shift, gid);
832 new_i_nblist(vdw, i_atom, shift, gid);
836 new_i_nblist(coul, i_atom, shift, gid);
838 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
839 bDoCoul_i = (bDoCoul && qi != 0);
841 if (bDoVdW_i || bDoCoul_i)
843 /* Loop over the j charge groups */
844 for (j = 0; (j < nj); j++)
848 /* Check for large charge groups */
859 /* Finally loop over the atoms in the j-charge group */
860 for (jj = jj0; jj < jj1; jj++)
862 bNotEx = NOTEXCL(bExcl, i, jj);
870 add_j_to_nblist(coul, jj, bLR);
875 if (bHaveVdW[type[jj]])
877 add_j_to_nblist(vdw, jj, bLR);
882 if (bHaveVdW[type[jj]])
886 add_j_to_nblist(vdwc, jj, bLR);
890 add_j_to_nblist(vdw, jj, bLR);
893 else if (charge[jj] != 0)
895 add_j_to_nblist(coul, jj, bLR);
903 close_i_nblist(coul);
904 close_i_nblist(vdwc);
910 /* we are doing free energy */
911 vdwc_free = &nlist[eNL_VDWQQ_FREE];
912 vdw_free = &nlist[eNL_VDW_FREE];
913 coul_free = &nlist[eNL_QQ_FREE];
914 /* Loop over the atoms in the i charge group */
915 for (i = 0; i < nicg; i++)
918 gid = GID(igid, jgid, ngid);
920 qiB = chargeB[i_atom];
922 /* Create new i_atom for each energy group */
923 if (bDoVdW && bDoCoul)
925 new_i_nblist(vdwc, i_atom, shift, gid);
929 new_i_nblist(vdw, i_atom, shift, gid);
933 new_i_nblist(coul, i_atom, shift, gid);
936 new_i_nblist(vdw_free, i_atom, shift, gid);
937 new_i_nblist(coul_free, i_atom, shift, gid);
938 new_i_nblist(vdwc_free, i_atom, shift, gid);
940 bDoVdW_i = (bDoVdW &&
941 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
942 bDoCoul_i = (bDoCoul && (qi != 0 || qiB != 0));
943 /* For TIP4P the first atom does not have a charge,
944 * but the last three do. So we should still put an atom
945 * without LJ but with charge in the water-atom neighborlist
946 * for a TIP4p i charge group.
947 * For SPC type water the first atom has LJ and charge,
948 * so there is no such problem.
950 if (iwater == esolNO)
952 bDoCoul_i_sol = bDoCoul_i;
956 bDoCoul_i_sol = bDoCoul;
959 if (bDoVdW_i || bDoCoul_i_sol)
961 /* Loop over the j charge groups */
962 for (j = 0; (j < nj); j++)
966 /* Check for large charge groups */
977 /* Finally loop over the atoms in the j-charge group */
978 bFree = bPert[i_atom];
979 for (jj = jj0; (jj < jj1); jj++)
981 bFreeJ = bFree || bPert[jj];
982 /* Complicated if, because the water H's should also
983 * see perturbed j-particles
985 if (iwater == esolNO || i == 0 || bFreeJ)
987 bNotEx = NOTEXCL(bExcl, i, jj);
995 if (charge[jj] != 0 || chargeB[jj] != 0)
997 add_j_to_nblist(coul_free, jj, bLR);
1000 else if (!bDoCoul_i)
1002 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1004 add_j_to_nblist(vdw_free, jj, bLR);
1009 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1011 if (charge[jj] != 0 || chargeB[jj] != 0)
1013 add_j_to_nblist(vdwc_free, jj, bLR);
1017 add_j_to_nblist(vdw_free, jj, bLR);
1020 else if (charge[jj] != 0 || chargeB[jj] != 0)
1022 add_j_to_nblist(coul_free, jj, bLR);
1028 /* This is done whether or not bWater is set */
1029 if (charge[jj] != 0)
1031 add_j_to_nblist(coul, jj, bLR);
1034 else if (!bDoCoul_i_sol)
1036 if (bHaveVdW[type[jj]])
1038 add_j_to_nblist(vdw, jj, bLR);
1043 if (bHaveVdW[type[jj]])
1045 if (charge[jj] != 0)
1047 add_j_to_nblist(vdwc, jj, bLR);
1051 add_j_to_nblist(vdw, jj, bLR);
1054 else if (charge[jj] != 0)
1056 add_j_to_nblist(coul, jj, bLR);
1064 close_i_nblist(vdw);
1065 close_i_nblist(coul);
1066 close_i_nblist(vdwc);
1067 close_i_nblist(vdw_free);
1068 close_i_nblist(coul_free);
1069 close_i_nblist(vdwc_free);
1075 put_in_list_adress(gmx_bool bHaveVdW[],
1091 /* The a[] index has been removed,
1092 * to put it back in i_atom should be a[i0] and jj should be a[jj].
1097 t_nblist * vdwc_adress = NULL;
1098 t_nblist * vdw_adress = NULL;
1099 t_nblist * coul_adress = NULL;
1100 t_nblist * vdwc_ww = NULL;
1101 t_nblist * coul_ww = NULL;
1103 int i, j, jcg, igid, gid, nbl_ind, nbl_ind_adress;
1104 atom_id jj, jj0, jj1, i_atom;
1109 real *charge, *chargeB;
1111 real qi, qiB, qq, rlj;
1112 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
1113 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
1115 gmx_bool j_all_atom;
1117 t_nblist *nlist, *nlist_adress;
1118 gmx_bool bEnergyGroupCG;
1120 /* Copy some pointers */
1121 cginfo = fr->cginfo;
1122 charge = md->chargeA;
1123 chargeB = md->chargeB;
1126 bPert = md->bPerturbed;
1129 /* Get atom range */
1131 nicg = index[icg+1]-i0;
1133 /* Get the i charge group info */
1134 igid = GET_CGINFO_GID(cginfo[icg]);
1136 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
1140 gmx_fatal(FARGS, "AdResS does not support free energy pertubation\n");
1143 /* Unpack pointers to neighbourlist structs */
1144 if (fr->nnblists == 2)
1151 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
1152 nbl_ind_adress = nbl_ind+fr->nnblists/2;
1156 nlist = fr->nblists[nbl_ind].nlist_lr;
1157 nlist_adress = fr->nblists[nbl_ind_adress].nlist_lr;
1161 nlist = fr->nblists[nbl_ind].nlist_sr;
1162 nlist_adress = fr->nblists[nbl_ind_adress].nlist_sr;
1166 vdwc = &nlist[eNL_VDWQQ];
1167 vdw = &nlist[eNL_VDW];
1168 coul = &nlist[eNL_QQ];
1170 vdwc_adress = &nlist_adress[eNL_VDWQQ];
1171 vdw_adress = &nlist_adress[eNL_VDW];
1172 coul_adress = &nlist_adress[eNL_QQ];
1174 /* We do not support solvent optimization with AdResS for now.
1175 For this we would need hybrid solvent-other kernels */
1177 /* no solvent as i charge group */
1178 /* Loop over the atoms in the i charge group */
1179 for (i = 0; i < nicg; i++)
1182 gid = GID(igid, jgid, ngid);
1183 qi = charge[i_atom];
1185 /* Create new i_atom for each energy group */
1186 if (bDoVdW && bDoCoul)
1188 new_i_nblist(vdwc, i_atom, shift, gid);
1189 new_i_nblist(vdwc_adress, i_atom, shift, gid);
1194 new_i_nblist(vdw, i_atom, shift, gid);
1195 new_i_nblist(vdw_adress, i_atom, shift, gid);
1200 new_i_nblist(coul, i_atom, shift, gid);
1201 new_i_nblist(coul_adress, i_atom, shift, gid);
1203 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
1204 bDoCoul_i = (bDoCoul && qi != 0);
1206 /* Here we find out whether the energy groups interaction belong to a
1207 * coarse-grained (vsite) or atomistic interaction. Note that, beacuse
1208 * interactions between coarse-grained and other (atomistic) energygroups
1209 * are excluded automatically by grompp, it is sufficient to check for
1210 * the group id of atom i (igid) */
1211 bEnergyGroupCG = !egp_explicit(fr, igid);
1213 if (bDoVdW_i || bDoCoul_i)
1215 /* Loop over the j charge groups */
1216 for (j = 0; (j < nj); j++)
1220 /* Check for large charge groups */
1231 /* Finally loop over the atoms in the j-charge group */
1232 for (jj = jj0; jj < jj1; jj++)
1234 bNotEx = NOTEXCL(bExcl, i, jj);
1236 /* Now we have to exclude interactions which will be zero
1237 * anyway due to the AdResS weights (in previous implementations
1238 * this was done in the force kernel). This is necessary as
1239 * pure interactions (those with b_hybrid=false, i.e. w_i*w_j==1 or 0)
1240 * are put into neighbour lists which will be passed to the
1241 * standard (optimized) kernels for speed. The interactions with
1242 * b_hybrid=true are placed into the _adress neighbour lists and
1243 * processed by the generic AdResS kernel.
1245 if ( (bEnergyGroupCG &&
1246 wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS ) ||
1247 ( !bEnergyGroupCG && wf[jj] <= GMX_REAL_EPS ) )
1252 b_hybrid = !((wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS) ||
1253 (wf[i_atom] <= GMX_REAL_EPS && wf[jj] <= GMX_REAL_EPS));
1259 if (charge[jj] != 0)
1263 add_j_to_nblist(coul, jj, bLR);
1267 add_j_to_nblist(coul_adress, jj, bLR);
1271 else if (!bDoCoul_i)
1273 if (bHaveVdW[type[jj]])
1277 add_j_to_nblist(vdw, jj, bLR);
1281 add_j_to_nblist(vdw_adress, jj, bLR);
1287 if (bHaveVdW[type[jj]])
1289 if (charge[jj] != 0)
1293 add_j_to_nblist(vdwc, jj, bLR);
1297 add_j_to_nblist(vdwc_adress, jj, bLR);
1304 add_j_to_nblist(vdw, jj, bLR);
1308 add_j_to_nblist(vdw_adress, jj, bLR);
1313 else if (charge[jj] != 0)
1317 add_j_to_nblist(coul, jj, bLR);
1321 add_j_to_nblist(coul_adress, jj, bLR);
1330 close_i_nblist(vdw);
1331 close_i_nblist(coul);
1332 close_i_nblist(vdwc);
1333 close_i_nblist(vdw_adress);
1334 close_i_nblist(coul_adress);
1335 close_i_nblist(vdwc_adress);
1341 put_in_list_qmmm(gmx_bool gmx_unused bHaveVdW[],
1343 t_mdatoms gmx_unused * md,
1353 gmx_bool gmx_unused bDoVdW,
1354 gmx_bool gmx_unused bDoCoul,
1355 int gmx_unused solvent_opt)
1358 int i, j, jcg, igid, gid;
1359 atom_id jj, jj0, jj1, i_atom;
1363 /* Get atom range */
1365 nicg = index[icg+1]-i0;
1367 /* Get the i charge group info */
1368 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1370 coul = &fr->QMMMlist;
1372 /* Loop over atoms in the ith charge group */
1373 for (i = 0; i < nicg; i++)
1376 gid = GID(igid, jgid, ngid);
1377 /* Create new i_atom for each energy group */
1378 new_i_nblist(coul, i_atom, shift, gid);
1380 /* Loop over the j charge groups */
1381 for (j = 0; j < nj; j++)
1385 /* Charge groups cannot have QM and MM atoms simultaneously */
1390 /* Finally loop over the atoms in the j-charge group */
1391 for (jj = jj0; jj < jj1; jj++)
1393 bNotEx = NOTEXCL(bExcl, i, jj);
1396 add_j_to_nblist(coul, jj, bLR);
1401 close_i_nblist(coul);
1406 put_in_list_cg(gmx_bool gmx_unused bHaveVdW[],
1408 t_mdatoms gmx_unused * md,
1418 gmx_bool gmx_unused bDoVdW,
1419 gmx_bool gmx_unused bDoCoul,
1420 int gmx_unused solvent_opt)
1423 int igid, gid, nbl_ind;
1427 cginfo = fr->cginfo[icg];
1429 igid = GET_CGINFO_GID(cginfo);
1430 gid = GID(igid, jgid, ngid);
1432 /* Unpack pointers to neighbourlist structs */
1433 if (fr->nnblists == 1)
1439 nbl_ind = fr->gid2nblists[gid];
1443 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1447 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1450 /* Make a new neighbor list for charge group icg.
1451 * Currently simply one neighbor list is made with LJ and Coulomb.
1452 * If required, zero interactions could be removed here
1453 * or in the force loop.
1455 new_i_nblist(vdwc, index[icg], shift, gid);
1456 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1458 for (j = 0; (j < nj); j++)
1461 /* Skip the icg-icg pairs if all self interactions are excluded */
1462 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1464 /* Here we add the j charge group jcg to the list,
1465 * exclusions are also added to the list.
1467 add_j_to_nblist_cg(vdwc, index[jcg], index[jcg+1], bExcl, icg == jcg, bLR);
1471 close_i_nblist(vdwc);
1474 static void setexcl(atom_id start, atom_id end, t_blocka *excl, gmx_bool b,
1481 for (i = start; i < end; i++)
1483 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1485 SETEXCL(bexcl, i-start, excl->a[k]);
1491 for (i = start; i < end; i++)
1493 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1495 RMEXCL(bexcl, i-start, excl->a[k]);
1501 int calc_naaj(int icg, int cgtot)
1505 if ((cgtot % 2) == 1)
1507 /* Odd number of charge groups, easy */
1508 naaj = 1 + (cgtot/2);
1510 else if ((cgtot % 4) == 0)
1512 /* Multiple of four is hard */
1549 fprintf(log, "naaj=%d\n", naaj);
1555 /************************************************
1557 * S I M P L E C O R E S T U F F
1559 ************************************************/
1561 static real calc_image_tric(rvec xi, rvec xj, matrix box,
1562 rvec b_inv, int *shift)
1564 /* This code assumes that the cut-off is smaller than
1565 * a half times the smallest diagonal element of the box.
1567 const real h25 = 2.5;
1572 /* Compute diff vector */
1573 dz = xj[ZZ] - xi[ZZ];
1574 dy = xj[YY] - xi[YY];
1575 dx = xj[XX] - xi[XX];
1577 /* Perform NINT operation, using trunc operation, therefore
1578 * we first add 2.5 then subtract 2 again
1580 tz = dz*b_inv[ZZ] + h25;
1582 dz -= tz*box[ZZ][ZZ];
1583 dy -= tz*box[ZZ][YY];
1584 dx -= tz*box[ZZ][XX];
1586 ty = dy*b_inv[YY] + h25;
1588 dy -= ty*box[YY][YY];
1589 dx -= ty*box[YY][XX];
1591 tx = dx*b_inv[XX]+h25;
1593 dx -= tx*box[XX][XX];
1595 /* Distance squared */
1596 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1598 *shift = XYZ2IS(tx, ty, tz);
1603 static real calc_image_rect(rvec xi, rvec xj, rvec box_size,
1604 rvec b_inv, int *shift)
1606 const real h15 = 1.5;
1612 /* Compute diff vector */
1613 dx = xj[XX] - xi[XX];
1614 dy = xj[YY] - xi[YY];
1615 dz = xj[ZZ] - xi[ZZ];
1617 /* Perform NINT operation, using trunc operation, therefore
1618 * we first add 1.5 then subtract 1 again
1620 tx = dx*b_inv[XX] + h15;
1621 ty = dy*b_inv[YY] + h15;
1622 tz = dz*b_inv[ZZ] + h15;
1627 /* Correct diff vector for translation */
1628 ddx = tx*box_size[XX] - dx;
1629 ddy = ty*box_size[YY] - dy;
1630 ddz = tz*box_size[ZZ] - dz;
1632 /* Distance squared */
1633 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1635 *shift = XYZ2IS(tx, ty, tz);
1640 static void add_simple(t_ns_buf *nsbuf, int nrj, atom_id cg_j,
1641 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1642 int icg, int jgid, t_block *cgs, t_excl bexcl[],
1643 int shift, t_forcerec *fr, put_in_list_t *put_in_list)
1645 if (nsbuf->nj + nrj > MAX_CG)
1647 put_in_list(bHaveVdW, ngid, md, icg, jgid, nsbuf->ncg, nsbuf->jcg,
1648 cgs->index, bexcl, shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1649 /* Reset buffer contents */
1650 nsbuf->ncg = nsbuf->nj = 0;
1652 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1656 static void ns_inner_tric(rvec x[], int icg, int *i_egp_flags,
1657 int njcg, atom_id jcg[],
1658 matrix box, rvec b_inv, real rcut2,
1659 t_block *cgs, t_ns_buf **ns_buf,
1660 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1661 t_excl bexcl[], t_forcerec *fr,
1662 put_in_list_t *put_in_list)
1666 int *cginfo = fr->cginfo;
1667 atom_id cg_j, *cgindex;
1670 cgindex = cgs->index;
1672 for (j = 0; (j < njcg); j++)
1675 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1676 if (calc_image_tric(x[icg], x[cg_j], box, b_inv, &shift) < rcut2)
1678 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1679 if (!(i_egp_flags[jgid] & EGP_EXCL))
1681 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1682 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1689 static void ns_inner_rect(rvec x[], int icg, int *i_egp_flags,
1690 int njcg, atom_id jcg[],
1691 gmx_bool bBox, rvec box_size, rvec b_inv, real rcut2,
1692 t_block *cgs, t_ns_buf **ns_buf,
1693 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1694 t_excl bexcl[], t_forcerec *fr,
1695 put_in_list_t *put_in_list)
1699 int *cginfo = fr->cginfo;
1700 atom_id cg_j, *cgindex;
1703 cgindex = cgs->index;
1707 for (j = 0; (j < njcg); j++)
1710 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1711 if (calc_image_rect(x[icg], x[cg_j], box_size, b_inv, &shift) < rcut2)
1713 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1714 if (!(i_egp_flags[jgid] & EGP_EXCL))
1716 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1717 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1725 for (j = 0; (j < njcg); j++)
1728 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1729 if ((rcut2 == 0) || (distance2(x[icg], x[cg_j]) < rcut2))
1731 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1732 if (!(i_egp_flags[jgid] & EGP_EXCL))
1734 add_simple(&ns_buf[jgid][CENTRAL], nrj, cg_j,
1735 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, CENTRAL, fr,
1743 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1745 static int ns_simple_core(t_forcerec *fr,
1746 gmx_localtop_t *top,
1748 matrix box, rvec box_size,
1749 t_excl bexcl[], atom_id *aaj,
1750 int ngid, t_ns_buf **ns_buf,
1751 put_in_list_t *put_in_list, gmx_bool bHaveVdW[])
1755 int nsearch, icg, jcg, igid, i0, nri, nn;
1758 /* atom_id *i_atoms; */
1759 t_block *cgs = &(top->cgs);
1760 t_blocka *excl = &(top->excls);
1763 gmx_bool bBox, bTriclinic;
1766 rlist2 = sqr(fr->rlist);
1768 bBox = (fr->ePBC != epbcNONE);
1771 for (m = 0; (m < DIM); m++)
1773 b_inv[m] = divide_err(1.0, box_size[m]);
1775 bTriclinic = TRICLINIC(box);
1782 cginfo = fr->cginfo;
1785 for (icg = fr->cg0; (icg < fr->hcg); icg++)
1788 i0 = cgs->index[icg];
1789 nri = cgs->index[icg+1]-i0;
1790 i_atoms = &(cgs->a[i0]);
1791 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1792 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1794 igid = GET_CGINFO_GID(cginfo[icg]);
1795 i_egp_flags = fr->egp_flags + ngid*igid;
1796 setexcl(cgs->index[icg], cgs->index[icg+1], excl, TRUE, bexcl);
1798 naaj = calc_naaj(icg, cgs->nr);
1801 ns_inner_tric(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1802 box, b_inv, rlist2, cgs, ns_buf,
1803 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1807 ns_inner_rect(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1808 bBox, box_size, b_inv, rlist2, cgs, ns_buf,
1809 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1813 for (nn = 0; (nn < ngid); nn++)
1815 for (k = 0; (k < SHIFTS); k++)
1817 nsbuf = &(ns_buf[nn][k]);
1820 put_in_list(bHaveVdW, ngid, md, icg, nn, nsbuf->ncg, nsbuf->jcg,
1821 cgs->index, bexcl, k, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1822 nsbuf->ncg = nsbuf->nj = 0;
1826 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1827 setexcl(cgs->index[icg], cgs->index[icg+1], excl, FALSE, bexcl);
1829 close_neighbor_lists(fr, FALSE);
1834 /************************************************
1836 * N S 5 G R I D S T U F F
1838 ************************************************/
1840 static gmx_inline void get_dx(int Nx, real gridx, real rc2, int xgi, real x,
1841 int *dx0, int *dx1, real *dcx2)
1869 for (i = xgi0; i >= 0; i--)
1871 dcx = (i+1)*gridx-x;
1880 for (i = xgi1; i < Nx; i++)
1893 static gmx_inline void get_dx_dd(int Nx, real gridx, real rc2, int xgi, real x,
1894 int ncpddc, int shift_min, int shift_max,
1895 int *g0, int *g1, real *dcx2)
1898 int g_min, g_max, shift_home;
1931 g_min = (shift_min == shift_home ? 0 : ncpddc);
1932 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1939 else if (shift_max < 0)
1954 /* Check one grid cell down */
1955 dcx = ((*g0 - 1) + 1)*gridx - x;
1967 /* Check one grid cell up */
1968 dcx = (*g1 + 1)*gridx - x;
1980 #define sqr(x) ((x)*(x))
1981 #define calc_dx2(XI, YI, ZI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
1982 #define calc_cyl_dx2(XI, YI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
1983 /****************************************************
1985 * F A S T N E I G H B O R S E A R C H I N G
1987 * Optimized neighboursearching routine using grid
1988 * at least 1x1x1, see GROMACS manual
1990 ****************************************************/
1993 static void get_cutoff2(t_forcerec *fr, gmx_bool bDoLongRange,
1994 real *rvdw2, real *rcoul2,
1995 real *rs2, real *rm2, real *rl2)
1997 *rs2 = sqr(fr->rlist);
1999 if (bDoLongRange && fr->bTwinRange)
2001 /* With plain cut-off or RF we need to make the list exactly
2002 * up to the cut-off and the cut-off's can be different,
2003 * so we can not simply set them to rlistlong.
2004 * To keep this code compatible with (exotic) old cases,
2005 * we also create lists up to rvdw/rcoulomb for PME and Ewald.
2006 * The interaction check should correspond to:
2007 * !ir_vdw/coulomb_might_be_zero_at_cutoff from inputrec.c.
2009 if (((fr->vdwtype == evdwCUT || fr->vdwtype == evdwPME) &&
2010 fr->vdw_modifier == eintmodNONE) ||
2011 fr->rvdw <= fr->rlist)
2013 *rvdw2 = sqr(fr->rvdw);
2017 *rvdw2 = sqr(fr->rlistlong);
2019 if (((fr->eeltype == eelCUT ||
2020 (EEL_RF(fr->eeltype) && fr->eeltype != eelRF_ZERO) ||
2021 fr->eeltype == eelPME ||
2022 fr->eeltype == eelEWALD) &&
2023 fr->coulomb_modifier == eintmodNONE) ||
2024 fr->rcoulomb <= fr->rlist)
2026 *rcoul2 = sqr(fr->rcoulomb);
2030 *rcoul2 = sqr(fr->rlistlong);
2035 /* Workaround for a gcc -O3 or -ffast-math problem */
2039 *rm2 = min(*rvdw2, *rcoul2);
2040 *rl2 = max(*rvdw2, *rcoul2);
2043 static void init_nsgrid_lists(t_forcerec *fr, int ngid, gmx_ns_t *ns)
2045 real rvdw2, rcoul2, rs2, rm2, rl2;
2048 get_cutoff2(fr, TRUE, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2050 /* Short range buffers */
2051 snew(ns->nl_sr, ngid);
2053 snew(ns->nsr, ngid);
2054 snew(ns->nlr_ljc, ngid);
2055 snew(ns->nlr_one, ngid);
2057 /* Always allocate both list types, since rcoulomb might now change with PME load balancing */
2058 /* Long range VdW and Coul buffers */
2059 snew(ns->nl_lr_ljc, ngid);
2060 /* Long range VdW or Coul only buffers */
2061 snew(ns->nl_lr_one, ngid);
2063 for (j = 0; (j < ngid); j++)
2065 snew(ns->nl_sr[j], MAX_CG);
2066 snew(ns->nl_lr_ljc[j], MAX_CG);
2067 snew(ns->nl_lr_one[j], MAX_CG);
2072 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
2077 static int nsgrid_core(t_commrec *cr, t_forcerec *fr,
2078 matrix box, int ngid,
2079 gmx_localtop_t *top,
2081 t_excl bexcl[], gmx_bool *bExcludeAlleg,
2083 put_in_list_t *put_in_list,
2084 gmx_bool bHaveVdW[],
2085 gmx_bool bDoLongRange, gmx_bool bMakeQMMMnblist)
2088 atom_id **nl_lr_ljc, **nl_lr_one, **nl_sr;
2089 int *nlr_ljc, *nlr_one, *nsr;
2091 t_block *cgs = &(top->cgs);
2092 int *cginfo = fr->cginfo;
2093 /* atom_id *i_atoms,*cgsindex=cgs->index; */
2095 int cell_x, cell_y, cell_z;
2096 int d, tx, ty, tz, dx, dy, dz, cj;
2097 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2098 int zsh_ty, zsh_tx, ysh_tx;
2100 int dx0, dx1, dy0, dy1, dz0, dz1;
2101 int Nx, Ny, Nz, shift = -1, j, nrj, nns, nn = -1;
2102 real gridx, gridy, gridz, grid_x, grid_y, grid_z;
2103 real *dcx2, *dcy2, *dcz2;
2105 int cg0, cg1, icg = -1, cgsnr, i0, igid, nri, naaj, max_jcg;
2106 int jcg0, jcg1, jjcg, cgj0, jgid;
2107 int *grida, *gridnra, *gridind;
2108 gmx_bool rvdw_lt_rcoul, rcoul_lt_rvdw;
2109 rvec xi, *cgcm, grid_offset;
2110 real r2, rs2, rvdw2, rcoul2, rm2, rl2, XI, YI, ZI, dcx, dcy, dcz, tmp1, tmp2;
2112 gmx_bool bDomDec, bTriclinicX, bTriclinicY;
2117 bDomDec = DOMAINDECOMP(cr);
2120 bTriclinicX = ((YY < grid->npbcdim &&
2121 (!bDomDec || dd->nc[YY] == 1) && box[YY][XX] != 0) ||
2122 (ZZ < grid->npbcdim &&
2123 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][XX] != 0));
2124 bTriclinicY = (ZZ < grid->npbcdim &&
2125 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][YY] != 0);
2129 get_cutoff2(fr, bDoLongRange, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2131 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
2132 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
2134 if (bMakeQMMMnblist)
2142 nl_lr_ljc = ns->nl_lr_ljc;
2143 nl_lr_one = ns->nl_lr_one;
2144 nlr_ljc = ns->nlr_ljc;
2145 nlr_one = ns->nlr_one;
2153 gridind = grid->index;
2154 gridnra = grid->nra;
2157 gridx = grid->cell_size[XX];
2158 gridy = grid->cell_size[YY];
2159 gridz = grid->cell_size[ZZ];
2163 copy_rvec(grid->cell_offset, grid_offset);
2164 copy_ivec(grid->ncpddc, ncpddc);
2169 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2170 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
2171 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
2172 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
2173 if (zsh_tx != 0 && ysh_tx != 0)
2175 /* This could happen due to rounding, when both ratios are 0.5 */
2184 /* We only want a list for the test particle */
2193 /* Set the shift range */
2194 for (d = 0; d < DIM; d++)
2198 /* Check if we need periodicity shifts.
2199 * Without PBC or with domain decomposition we don't need them.
2201 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2208 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
2219 /* Loop over charge groups */
2220 for (icg = cg0; (icg < cg1); icg++)
2222 igid = GET_CGINFO_GID(cginfo[icg]);
2223 /* Skip this charge group if all energy groups are excluded! */
2224 if (bExcludeAlleg[igid])
2229 i0 = cgs->index[icg];
2231 if (bMakeQMMMnblist)
2233 /* Skip this charge group if it is not a QM atom while making a
2234 * QM/MM neighbourlist
2236 if (md->bQM[i0] == FALSE)
2238 continue; /* MM particle, go to next particle */
2241 /* Compute the number of charge groups that fall within the control
2244 naaj = calc_naaj(icg, cgsnr);
2251 /* make a normal neighbourlist */
2255 /* Get the j charge-group and dd cell shift ranges */
2256 dd_get_ns_ranges(cr->dd, icg, &jcg0, &jcg1, sh0, sh1);
2261 /* Compute the number of charge groups that fall within the control
2264 naaj = calc_naaj(icg, cgsnr);
2270 /* The i-particle is awlways the test particle,
2271 * so we want all j-particles
2273 max_jcg = cgsnr - 1;
2277 max_jcg = jcg1 - cgsnr;
2282 i_egp_flags = fr->egp_flags + igid*ngid;
2284 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2285 setexcl(i0, cgs->index[icg+1], &top->excls, TRUE, bexcl);
2287 ci2xyz(grid, icg, &cell_x, &cell_y, &cell_z);
2289 /* Changed iicg to icg, DvdS 990115
2290 * (but see consistency check above, DvdS 990330)
2293 fprintf(log, "icg=%5d, naaj=%5d, cell %d %d %d\n",
2294 icg, naaj, cell_x, cell_y, cell_z);
2296 /* Loop over shift vectors in three dimensions */
2297 for (tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
2299 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2300 /* Calculate range of cells in Z direction that have the shift tz */
2301 zgi = cell_z + tz*Nz;
2304 get_dx(Nz, gridz, rl2, zgi, ZI, &dz0, &dz1, dcz2);
2306 get_dx_dd(Nz, gridz, rl2, zgi, ZI-grid_offset[ZZ],
2307 ncpddc[ZZ], sh0[ZZ], sh1[ZZ], &dz0, &dz1, dcz2);
2313 for (ty = -shp[YY]; ty <= shp[YY]; ty++)
2315 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2316 /* Calculate range of cells in Y direction that have the shift ty */
2319 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2323 ygi = cell_y + ty*Ny;
2326 get_dx(Ny, gridy, rl2, ygi, YI, &dy0, &dy1, dcy2);
2328 get_dx_dd(Ny, gridy, rl2, ygi, YI-grid_offset[YY],
2329 ncpddc[YY], sh0[YY], sh1[YY], &dy0, &dy1, dcy2);
2335 for (tx = -shp[XX]; tx <= shp[XX]; tx++)
2337 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2338 /* Calculate range of cells in X direction that have the shift tx */
2341 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2345 xgi = cell_x + tx*Nx;
2348 get_dx(Nx, gridx, rl2, xgi*Nx, XI, &dx0, &dx1, dcx2);
2350 get_dx_dd(Nx, gridx, rl2, xgi, XI-grid_offset[XX],
2351 ncpddc[XX], sh0[XX], sh1[XX], &dx0, &dx1, dcx2);
2357 /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2358 * from the neigbour list as it will not interact */
2359 if (fr->adress_type != eAdressOff)
2361 if (md->wf[cgs->index[icg]] <= GMX_REAL_EPS && egp_explicit(fr, igid))
2366 /* Get shift vector */
2367 shift = XYZ2IS(tx, ty, tz);
2369 range_check(shift, 0, SHIFTS);
2371 for (nn = 0; (nn < ngid); nn++)
2378 fprintf(log, "shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2379 shift, dx0, dx1, dy0, dy1, dz0, dz1);
2380 fprintf(log, "cgcm: %8.3f %8.3f %8.3f\n", cgcm[icg][XX],
2381 cgcm[icg][YY], cgcm[icg][ZZ]);
2382 fprintf(log, "xi: %8.3f %8.3f %8.3f\n", XI, YI, ZI);
2384 for (dx = dx0; (dx <= dx1); dx++)
2386 tmp1 = rl2 - dcx2[dx];
2387 for (dy = dy0; (dy <= dy1); dy++)
2389 tmp2 = tmp1 - dcy2[dy];
2392 for (dz = dz0; (dz <= dz1); dz++)
2394 if (tmp2 > dcz2[dz])
2396 /* Find grid-cell cj in which possible neighbours are */
2397 cj = xyz2ci(Ny, Nz, dx, dy, dz);
2399 /* Check out how many cgs (nrj) there in this cell */
2402 /* Find the offset in the cg list */
2405 /* Check if all j's are out of range so we
2406 * can skip the whole cell.
2407 * Should save some time, especially with DD.
2410 (grida[cgj0] >= max_jcg &&
2411 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2417 for (j = 0; (j < nrj); j++)
2419 jjcg = grida[cgj0+j];
2421 /* check whether this guy is in range! */
2422 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2425 r2 = calc_dx2(XI, YI, ZI, cgcm[jjcg]);
2428 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2429 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2430 /* check energy group exclusions */
2431 if (!(i_egp_flags[jgid] & EGP_EXCL))
2435 if (nsr[jgid] >= MAX_CG)
2437 /* Add to short-range list */
2438 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2439 nsr[jgid], nl_sr[jgid],
2440 cgs->index, /* cgsatoms, */ bexcl,
2441 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2444 nl_sr[jgid][nsr[jgid]++] = jjcg;
2448 if (nlr_ljc[jgid] >= MAX_CG)
2450 /* Add to LJ+coulomb long-range list */
2451 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2452 nlr_ljc[jgid], nl_lr_ljc[jgid], top->cgs.index,
2453 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2456 nl_lr_ljc[jgid][nlr_ljc[jgid]++] = jjcg;
2460 if (nlr_one[jgid] >= MAX_CG)
2462 /* Add to long-range list with only coul, or only LJ */
2463 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2464 nlr_one[jgid], nl_lr_one[jgid], top->cgs.index,
2465 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2468 nl_lr_one[jgid][nlr_one[jgid]++] = jjcg;
2480 /* CHECK whether there is anything left in the buffers */
2481 for (nn = 0; (nn < ngid); nn++)
2485 put_in_list(bHaveVdW, ngid, md, icg, nn, nsr[nn], nl_sr[nn],
2486 cgs->index, /* cgsatoms, */ bexcl,
2487 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2490 if (nlr_ljc[nn] > 0)
2492 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_ljc[nn],
2493 nl_lr_ljc[nn], top->cgs.index,
2494 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2497 if (nlr_one[nn] > 0)
2499 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_one[nn],
2500 nl_lr_one[nn], top->cgs.index,
2501 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2507 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2508 setexcl(cgs->index[icg], cgs->index[icg+1], &top->excls, FALSE, bexcl);
2510 /* No need to perform any left-over force calculations anymore (as we used to do here)
2511 * since we now save the proper long-range lists for later evaluation.
2516 /* Close neighbourlists */
2517 close_neighbor_lists(fr, bMakeQMMMnblist);
2522 void ns_realloc_natoms(gmx_ns_t *ns, int natoms)
2526 if (natoms > ns->nra_alloc)
2528 ns->nra_alloc = over_alloc_dd(natoms);
2529 srenew(ns->bexcl, ns->nra_alloc);
2530 for (i = 0; i < ns->nra_alloc; i++)
2537 void init_ns(FILE *fplog, const t_commrec *cr,
2538 gmx_ns_t *ns, t_forcerec *fr,
2539 const gmx_mtop_t *mtop)
2541 int mt, icg, nr_in_cg, maxcg, i, j, jcg, ngid, ncg;
2545 /* Compute largest charge groups size (# atoms) */
2547 for (mt = 0; mt < mtop->nmoltype; mt++)
2549 cgs = &mtop->moltype[mt].cgs;
2550 for (icg = 0; (icg < cgs->nr); icg++)
2552 nr_in_cg = max(nr_in_cg, (int)(cgs->index[icg+1]-cgs->index[icg]));
2556 /* Verify whether largest charge group is <= max cg.
2557 * This is determined by the type of the local exclusion type
2558 * Exclusions are stored in bits. (If the type is not large
2559 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2561 maxcg = sizeof(t_excl)*8;
2562 if (nr_in_cg > maxcg)
2564 gmx_fatal(FARGS, "Max #atoms in a charge group: %d > %d\n",
2568 ngid = mtop->groups.grps[egcENER].nr;
2569 snew(ns->bExcludeAlleg, ngid);
2570 for (i = 0; i < ngid; i++)
2572 ns->bExcludeAlleg[i] = TRUE;
2573 for (j = 0; j < ngid; j++)
2575 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2577 ns->bExcludeAlleg[i] = FALSE;
2585 ns->grid = init_grid(fplog, fr);
2586 init_nsgrid_lists(fr, ngid, ns);
2591 snew(ns->ns_buf, ngid);
2592 for (i = 0; (i < ngid); i++)
2594 snew(ns->ns_buf[i], SHIFTS);
2596 ncg = ncg_mtop(mtop);
2597 snew(ns->simple_aaj, 2*ncg);
2598 for (jcg = 0; (jcg < ncg); jcg++)
2600 ns->simple_aaj[jcg] = jcg;
2601 ns->simple_aaj[jcg+ncg] = jcg;
2605 /* Create array that determines whether or not atoms have VdW */
2606 snew(ns->bHaveVdW, fr->ntype);
2607 for (i = 0; (i < fr->ntype); i++)
2609 for (j = 0; (j < fr->ntype); j++)
2611 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2613 ((BHAMA(fr->nbfp, fr->ntype, i, j) != 0) ||
2614 (BHAMB(fr->nbfp, fr->ntype, i, j) != 0) ||
2615 (BHAMC(fr->nbfp, fr->ntype, i, j) != 0)) :
2616 ((C6(fr->nbfp, fr->ntype, i, j) != 0) ||
2617 (C12(fr->nbfp, fr->ntype, i, j) != 0))));
2622 pr_bvec(debug, 0, "bHaveVdW", ns->bHaveVdW, fr->ntype, TRUE);
2627 if (!DOMAINDECOMP(cr))
2629 ns_realloc_natoms(ns, mtop->natoms);
2632 ns->nblist_initialized = FALSE;
2634 /* nbr list debug dump */
2636 char *ptr = getenv("GMX_DUMP_NL");
2639 ns->dump_nl = strtol(ptr, NULL, 10);
2642 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2653 int search_neighbours(FILE *log, t_forcerec *fr,
2655 gmx_localtop_t *top,
2656 gmx_groups_t *groups,
2658 t_nrnb *nrnb, t_mdatoms *md,
2660 gmx_bool bDoLongRangeNS)
2662 t_block *cgs = &(top->cgs);
2663 rvec box_size, grid_x0, grid_x1;
2665 real min_size, grid_dens;
2669 gmx_bool *i_egp_flags;
2670 int cg_start, cg_end, start, end;
2673 gmx_domdec_zones_t *dd_zones;
2674 put_in_list_t *put_in_list;
2678 /* Set some local variables */
2680 ngid = groups->grps[egcENER].nr;
2682 for (m = 0; (m < DIM); m++)
2684 box_size[m] = box[m][m];
2687 if (fr->ePBC != epbcNONE)
2689 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC, box))
2691 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.");
2695 min_size = min(box_size[XX], min(box_size[YY], box_size[ZZ]));
2696 if (2*fr->rlistlong >= min_size)
2698 gmx_fatal(FARGS, "One of the box diagonal elements has become smaller than twice the cut-off length.");
2703 if (DOMAINDECOMP(cr))
2705 ns_realloc_natoms(ns, cgs->index[cgs->nr]);
2709 /* Reset the neighbourlists */
2710 reset_neighbor_lists(fr, TRUE, TRUE);
2712 if (bGrid && bFillGrid)
2716 if (DOMAINDECOMP(cr))
2718 dd_zones = domdec_zones(cr->dd);
2724 get_nsgrid_boundaries(grid->nboundeddim, box, NULL, NULL, NULL, NULL,
2725 cgs->nr, fr->cg_cm, grid_x0, grid_x1, &grid_dens);
2727 grid_first(log, grid, NULL, NULL, box, grid_x0, grid_x1,
2728 fr->rlistlong, grid_dens);
2735 if (DOMAINDECOMP(cr))
2738 fill_grid(dd_zones, grid, end, -1, end, fr->cg_cm);
2740 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2744 fill_grid(NULL, grid, cgs->nr, fr->cg0, fr->hcg, fr->cg_cm);
2745 grid->icg0 = fr->cg0;
2746 grid->icg1 = fr->hcg;
2750 calc_elemnr(grid, start, end, cgs->nr);
2752 grid_last(grid, start, end, cgs->nr);
2757 print_grid(debug, grid);
2762 /* Set the grid cell index for the test particle only.
2763 * The cell to cg index is not corrected, but that does not matter.
2765 fill_grid(NULL, ns->grid, fr->hcg, fr->hcg-1, fr->hcg, fr->cg_cm);
2769 if (fr->adress_type == eAdressOff)
2771 if (!fr->ns.bCGlist)
2773 put_in_list = put_in_list_at;
2777 put_in_list = put_in_list_cg;
2782 put_in_list = put_in_list_adress;
2789 nsearch = nsgrid_core(cr, fr, box, ngid, top,
2790 grid, ns->bexcl, ns->bExcludeAlleg,
2791 md, put_in_list, ns->bHaveVdW,
2792 bDoLongRangeNS, FALSE);
2794 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2795 * the classical calculation. The charge-charge interaction
2796 * between QM and MM atoms is handled in the QMMM core calculation
2797 * (see QMMM.c). The VDW however, we'd like to compute classically
2798 * and the QM MM atom pairs have just been put in the
2799 * corresponding neighbourlists. in case of QMMM we still need to
2800 * fill a special QMMM neighbourlist that contains all neighbours
2801 * of the QM atoms. If bQMMM is true, this list will now be made:
2803 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
2805 nsearch += nsgrid_core(cr, fr, box, ngid, top,
2806 grid, ns->bexcl, ns->bExcludeAlleg,
2807 md, put_in_list_qmmm, ns->bHaveVdW,
2808 bDoLongRangeNS, TRUE);
2813 nsearch = ns_simple_core(fr, top, md, box, box_size,
2814 ns->bexcl, ns->simple_aaj,
2815 ngid, ns->ns_buf, put_in_list, ns->bHaveVdW);
2823 inc_nrnb(nrnb, eNR_NS, nsearch);
2824 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2829 int natoms_beyond_ns_buffer(t_inputrec *ir, t_forcerec *fr, t_block *cgs,
2830 matrix scale_tot, rvec *x)
2832 int cg0, cg1, cg, a0, a1, a, i, j;
2833 real rint, hbuf2, scale;
2835 gmx_bool bIsotropic;
2840 rint = max(ir->rcoulomb, ir->rvdw);
2841 if (ir->rlist < rint)
2843 gmx_fatal(FARGS, "The neighbor search buffer has negative size: %f nm",
2851 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2853 hbuf2 = sqr(0.5*(ir->rlist - rint));
2854 for (cg = cg0; cg < cg1; cg++)
2856 a0 = cgs->index[cg];
2857 a1 = cgs->index[cg+1];
2858 for (a = a0; a < a1; a++)
2860 if (distance2(cg_cm[cg], x[a]) > hbuf2)
2870 scale = scale_tot[0][0];
2871 for (i = 1; i < DIM; i++)
2873 /* With anisotropic scaling, the original spherical ns volumes become
2874 * ellipsoids. To avoid costly transformations we use the minimum
2875 * eigenvalue of the scaling matrix for determining the buffer size.
2876 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2878 scale = min(scale, scale_tot[i][i]);
2879 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2883 for (j = 0; j < i; j++)
2885 if (scale_tot[i][j] != 0)
2891 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2894 for (cg = cg0; cg < cg1; cg++)
2896 svmul(scale, cg_cm[cg], cgsc);
2897 a0 = cgs->index[cg];
2898 a1 = cgs->index[cg+1];
2899 for (a = a0; a < a1; a++)
2901 if (distance2(cgsc, x[a]) > hbuf2)
2910 /* Anistropic scaling */
2911 for (cg = cg0; cg < cg1; cg++)
2913 /* Since scale_tot contains the transpose of the scaling matrix,
2914 * we need to multiply with the transpose.
2916 tmvmul_ur0(scale_tot, cg_cm[cg], cgsc);
2917 a0 = cgs->index[cg];
2918 a1 = cgs->index[cg+1];
2919 for (a = a0; a < a1; a++)
2921 if (distance2(cgsc, x[a]) > hbuf2)