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.
46 #include "gromacs/math/utilities.h"
51 #include "nonbonded.h"
55 #include "gmx_fatal.h"
58 #include "mtop_util.h"
65 * E X C L U S I O N H A N D L I N G
69 static void SETEXCL_(t_excl e[], atom_id i, atom_id j)
73 static void RMEXCL_(t_excl e[], atom_id i, atom_id j)
75 e[j] = e[j] & ~(1<<i);
77 static gmx_bool ISEXCL_(t_excl e[], atom_id i, atom_id j)
79 return (gmx_bool)(e[j] & (1<<i));
81 static gmx_bool NOTEXCL_(t_excl e[], atom_id i, atom_id j)
83 return !(ISEXCL(e, i, j));
86 #define SETEXCL(e, i, j) (e)[((atom_id) (j))] |= (1<<((atom_id) (i)))
87 #define RMEXCL(e, i, j) (e)[((atom_id) (j))] &= (~(1<<((atom_id) (i))))
88 #define ISEXCL(e, i, j) (gmx_bool) ((e)[((atom_id) (j))] & (1<<((atom_id) (i))))
89 #define NOTEXCL(e, i, j) !(ISEXCL(e, i, j))
93 round_up_to_simd_width(int length, int simd_width)
95 int offset, newlength;
97 offset = (simd_width > 0) ? length % simd_width : 0;
99 return (offset == 0) ? length : length-offset+simd_width;
101 /************************************************
103 * U T I L I T I E S F O R N S
105 ************************************************/
107 void reallocate_nblist(t_nblist *nl)
111 fprintf(debug, "reallocating neigborlist (ielec=%d, ivdw=%d, igeometry=%d, type=%d), maxnri=%d\n",
112 nl->ielec, nl->ivdw, nl->igeometry, nl->type, nl->maxnri);
114 srenew(nl->iinr, nl->maxnri);
115 if (nl->igeometry == GMX_NBLIST_GEOMETRY_CG_CG)
117 srenew(nl->iinr_end, nl->maxnri);
119 srenew(nl->gid, nl->maxnri);
120 srenew(nl->shift, nl->maxnri);
121 srenew(nl->jindex, nl->maxnri+1);
125 static void init_nblist(FILE *log, t_nblist *nl_sr, t_nblist *nl_lr,
126 int maxsr, int maxlr,
127 int ivdw, int ivdwmod,
128 int ielec, int ielecmod,
129 int igeometry, int type)
135 for (i = 0; (i < 2); i++)
137 nl = (i == 0) ? nl_sr : nl_lr;
138 homenr = (i == 0) ? maxsr : maxlr;
146 /* Set coul/vdw in neighborlist, and for the normal loops we determine
147 * an index of which one to call.
150 nl->ivdwmod = ivdwmod;
152 nl->ielecmod = ielecmod;
154 nl->igeometry = igeometry;
156 if (nl->type == GMX_NBLIST_INTERACTION_FREE_ENERGY)
158 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
161 /* This will also set the simd_padding_width field */
162 gmx_nonbonded_set_kernel_pointers( (i == 0) ? log : NULL, nl);
164 /* maxnri is influenced by the number of shifts (maximum is 8)
165 * and the number of energy groups.
166 * If it is not enough, nl memory will be reallocated during the run.
167 * 4 seems to be a reasonable factor, which only causes reallocation
168 * during runs with tiny and many energygroups.
170 nl->maxnri = homenr*4;
180 reallocate_nblist(nl);
185 fprintf(debug, "Initiating neighbourlist (ielec=%d, ivdw=%d, type=%d) for %s interactions,\nwith %d SR, %d LR atoms.\n",
186 nl->ielec, nl->ivdw, nl->type, gmx_nblist_geometry_names[nl->igeometry], maxsr, maxlr);
191 void init_neighbor_list(FILE *log, t_forcerec *fr, int homenr)
193 /* Make maxlr tunable! (does not seem to be a big difference though)
194 * This parameter determines the number of i particles in a long range
195 * neighbourlist. Too few means many function calls, too many means
198 int maxsr, maxsr_wat, maxlr, maxlr_wat;
199 int ielec, ielecf, ivdw, ielecmod, ielecmodf, ivdwmod, type;
201 int igeometry_def, igeometry_w, igeometry_ww;
205 /* maxsr = homenr-fr->nWatMol*3; */
210 gmx_fatal(FARGS, "%s, %d: Negative number of short range atoms.\n"
211 "Call your Gromacs dealer for assistance.", __FILE__, __LINE__);
213 /* This is just for initial allocation, so we do not reallocate
214 * all the nlist arrays many times in a row.
215 * The numbers seem very accurate, but they are uncritical.
217 maxsr_wat = min(fr->nWatMol, (homenr+2)/3);
221 maxlr_wat = min(maxsr_wat, maxlr);
225 maxlr = maxlr_wat = 0;
228 /* Determine the values for ielec/ivdw. */
229 ielec = fr->nbkernel_elec_interaction;
230 ivdw = fr->nbkernel_vdw_interaction;
231 ielecmod = fr->nbkernel_elec_modifier;
232 ivdwmod = fr->nbkernel_vdw_modifier;
233 type = GMX_NBLIST_INTERACTION_STANDARD;
235 fr->ns.bCGlist = (getenv("GMX_NBLISTCG") != 0);
238 igeometry_def = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
242 igeometry_def = GMX_NBLIST_GEOMETRY_CG_CG;
245 fprintf(log, "\nUsing charge-group - charge-group neighbor lists and kernels\n\n");
249 if (fr->solvent_opt == esolTIP4P)
251 igeometry_w = GMX_NBLIST_GEOMETRY_WATER4_PARTICLE;
252 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER4_WATER4;
256 igeometry_w = GMX_NBLIST_GEOMETRY_WATER3_PARTICLE;
257 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER3_WATER3;
260 for (i = 0; i < fr->nnblists; i++)
262 nbl = &(fr->nblists[i]);
264 if ((fr->adress_type != eAdressOff) && (i >= fr->nnblists/2))
266 type = GMX_NBLIST_INTERACTION_ADRESS;
268 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ], &nbl->nlist_lr[eNL_VDWQQ],
269 maxsr, maxlr, ivdw, ivdwmod, ielec, ielecmod, igeometry_def, type);
270 init_nblist(log, &nbl->nlist_sr[eNL_VDW], &nbl->nlist_lr[eNL_VDW],
271 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, igeometry_def, type);
272 init_nblist(log, &nbl->nlist_sr[eNL_QQ], &nbl->nlist_lr[eNL_QQ],
273 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_def, type);
274 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATER], &nbl->nlist_lr[eNL_VDWQQ_WATER],
275 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_w, type);
276 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATER], &nbl->nlist_lr[eNL_QQ_WATER],
277 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_w, type);
278 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATERWATER], &nbl->nlist_lr[eNL_VDWQQ_WATERWATER],
279 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_ww, type);
280 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATERWATER], &nbl->nlist_lr[eNL_QQ_WATERWATER],
281 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_ww, type);
283 /* Did we get the solvent loops so we can use optimized water kernels? */
284 if (nbl->nlist_sr[eNL_VDWQQ_WATER].kernelptr_vf == NULL
285 || nbl->nlist_sr[eNL_QQ_WATER].kernelptr_vf == NULL
286 #ifndef DISABLE_WATERWATER_NLIST
287 || nbl->nlist_sr[eNL_VDWQQ_WATERWATER].kernelptr_vf == NULL
288 || nbl->nlist_sr[eNL_QQ_WATERWATER].kernelptr_vf == NULL
292 fr->solvent_opt = esolNO;
293 fprintf(log, "Note: The available nonbonded kernels do not support water optimization - disabling.\n");
296 if (fr->efep != efepNO)
298 if ((fr->bEwald) && (fr->sc_alphacoul > 0)) /* need to handle long range differently if using softcore */
300 ielecf = GMX_NBKERNEL_ELEC_EWALD;
301 ielecmodf = eintmodNONE;
306 ielecmodf = ielecmod;
309 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_FREE], &nbl->nlist_lr[eNL_VDWQQ_FREE],
310 maxsr, maxlr, ivdw, ivdwmod, ielecf, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
311 init_nblist(log, &nbl->nlist_sr[eNL_VDW_FREE], &nbl->nlist_lr[eNL_VDW_FREE],
312 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
313 init_nblist(log, &nbl->nlist_sr[eNL_QQ_FREE], &nbl->nlist_lr[eNL_QQ_FREE],
314 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielecf, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
318 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
320 init_nblist(log, &fr->QMMMlist, NULL,
321 maxsr, maxlr, 0, 0, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_STANDARD);
329 fr->ns.nblist_initialized = TRUE;
332 static void reset_nblist(t_nblist *nl)
342 static void reset_neighbor_lists(t_forcerec *fr, gmx_bool bResetSR, gmx_bool bResetLR)
348 /* only reset the short-range nblist */
349 reset_nblist(&(fr->QMMMlist));
352 for (n = 0; n < fr->nnblists; n++)
354 for (i = 0; i < eNL_NR; i++)
358 reset_nblist( &(fr->nblists[n].nlist_sr[i]) );
362 reset_nblist( &(fr->nblists[n].nlist_lr[i]) );
371 static inline void new_i_nblist(t_nblist *nlist, atom_id i_atom, int shift, int gid)
373 int i, k, nri, nshift;
377 /* Check whether we have to increase the i counter */
379 (nlist->iinr[nri] != i_atom) ||
380 (nlist->shift[nri] != shift) ||
381 (nlist->gid[nri] != gid))
383 /* This is something else. Now see if any entries have
384 * been added in the list of the previous atom.
387 ((nlist->jindex[nri+1] > nlist->jindex[nri]) &&
388 (nlist->gid[nri] != -1)))
390 /* If so increase the counter */
393 if (nlist->nri >= nlist->maxnri)
395 nlist->maxnri += over_alloc_large(nlist->nri);
396 reallocate_nblist(nlist);
399 /* Set the number of neighbours and the atom number */
400 nlist->jindex[nri+1] = nlist->jindex[nri];
401 nlist->iinr[nri] = i_atom;
402 nlist->gid[nri] = gid;
403 nlist->shift[nri] = shift;
407 /* Adding to previous list. First remove possible previous padding */
408 if (nlist->simd_padding_width > 1)
410 while (nlist->nrj > 0 && nlist->jjnr[nlist->nrj-1] < 0)
418 static inline void close_i_nblist(t_nblist *nlist)
420 int nri = nlist->nri;
425 /* Add elements up to padding. Since we allocate memory in units
426 * of the simd_padding width, we do not have to check for possible
427 * list reallocation here.
429 while ((nlist->nrj % nlist->simd_padding_width) != 0)
431 /* Use -4 here, so we can write forces for 4 atoms before real data */
432 nlist->jjnr[nlist->nrj++] = -4;
434 nlist->jindex[nri+1] = nlist->nrj;
436 len = nlist->nrj - nlist->jindex[nri];
440 static inline void close_nblist(t_nblist *nlist)
442 /* Only close this nblist when it has been initialized.
443 * Avoid the creation of i-lists with no j-particles.
447 /* Some assembly kernels do not support empty lists,
448 * make sure here that we don't generate any empty lists.
449 * With the current ns code this branch is taken in two cases:
450 * No i-particles at all: nri=-1 here
451 * There are i-particles, but no j-particles; nri=0 here
457 /* Close list number nri by incrementing the count */
462 static inline void close_neighbor_lists(t_forcerec *fr, gmx_bool bMakeQMMMnblist)
468 close_nblist(&(fr->QMMMlist));
471 for (n = 0; n < fr->nnblists; n++)
473 for (i = 0; (i < eNL_NR); i++)
475 close_nblist(&(fr->nblists[n].nlist_sr[i]));
476 close_nblist(&(fr->nblists[n].nlist_lr[i]));
482 static inline void add_j_to_nblist(t_nblist *nlist, atom_id j_atom, gmx_bool bLR)
484 int nrj = nlist->nrj;
486 if (nlist->nrj >= nlist->maxnrj)
488 nlist->maxnrj = round_up_to_simd_width(over_alloc_small(nlist->nrj + 1), nlist->simd_padding_width);
492 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
493 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
496 srenew(nlist->jjnr, nlist->maxnrj);
499 nlist->jjnr[nrj] = j_atom;
503 static inline void add_j_to_nblist_cg(t_nblist *nlist,
504 atom_id j_start, int j_end,
505 t_excl *bexcl, gmx_bool i_is_j,
508 int nrj = nlist->nrj;
511 if (nlist->nrj >= nlist->maxnrj)
513 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
516 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
517 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
520 srenew(nlist->jjnr, nlist->maxnrj);
521 srenew(nlist->jjnr_end, nlist->maxnrj);
522 srenew(nlist->excl, nlist->maxnrj*MAX_CGCGSIZE);
525 nlist->jjnr[nrj] = j_start;
526 nlist->jjnr_end[nrj] = j_end;
528 if (j_end - j_start > MAX_CGCGSIZE)
530 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);
533 /* Set the exclusions */
534 for (j = j_start; j < j_end; j++)
536 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
540 /* Avoid double counting of intra-cg interactions */
541 for (j = 1; j < j_end-j_start; j++)
543 nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
551 put_in_list_t (gmx_bool bHaveVdW[],
568 put_in_list_at(gmx_bool bHaveVdW[],
584 /* The a[] index has been removed,
585 * to put it back in i_atom should be a[i0] and jj should be a[jj].
590 t_nblist * vdwc_free = NULL;
591 t_nblist * vdw_free = NULL;
592 t_nblist * coul_free = NULL;
593 t_nblist * vdwc_ww = NULL;
594 t_nblist * coul_ww = NULL;
596 int i, j, jcg, igid, gid, nbl_ind, ind_ij;
597 atom_id jj, jj0, jj1, i_atom;
602 real *charge, *chargeB;
603 real qi, qiB, qq, rlj;
604 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
605 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
609 /* Copy some pointers */
611 charge = md->chargeA;
612 chargeB = md->chargeB;
615 bPert = md->bPerturbed;
619 nicg = index[icg+1]-i0;
621 /* Get the i charge group info */
622 igid = GET_CGINFO_GID(cginfo[icg]);
624 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
629 /* Check if any of the particles involved are perturbed.
630 * If not we can do the cheaper normal put_in_list
631 * and use more solvent optimization.
633 for (i = 0; i < nicg; i++)
635 bFreeEnergy |= bPert[i0+i];
637 /* Loop over the j charge groups */
638 for (j = 0; (j < nj && !bFreeEnergy); j++)
643 /* Finally loop over the atoms in the j-charge group */
644 for (jj = jj0; jj < jj1; jj++)
646 bFreeEnergy |= bPert[jj];
651 /* Unpack pointers to neighbourlist structs */
652 if (fr->nnblists == 1)
658 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
662 nlist = fr->nblists[nbl_ind].nlist_lr;
666 nlist = fr->nblists[nbl_ind].nlist_sr;
669 if (iwater != esolNO)
671 vdwc = &nlist[eNL_VDWQQ_WATER];
672 vdw = &nlist[eNL_VDW];
673 coul = &nlist[eNL_QQ_WATER];
674 #ifndef DISABLE_WATERWATER_NLIST
675 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
676 coul_ww = &nlist[eNL_QQ_WATERWATER];
681 vdwc = &nlist[eNL_VDWQQ];
682 vdw = &nlist[eNL_VDW];
683 coul = &nlist[eNL_QQ];
688 if (iwater != esolNO)
690 /* Loop over the atoms in the i charge group */
692 gid = GID(igid, jgid, ngid);
693 /* Create new i_atom for each energy group */
694 if (bDoCoul && bDoVdW)
696 new_i_nblist(vdwc, i_atom, shift, gid);
697 #ifndef DISABLE_WATERWATER_NLIST
698 new_i_nblist(vdwc_ww, i_atom, shift, gid);
703 new_i_nblist(vdw, i_atom, shift, gid);
707 new_i_nblist(coul, i_atom, shift, gid);
708 #ifndef DISABLE_WATERWATER_NLIST
709 new_i_nblist(coul_ww, i_atom, shift, gid);
712 /* Loop over the j charge groups */
713 for (j = 0; (j < nj); j++)
723 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
725 if (iwater == esolSPC && jwater == esolSPC)
727 /* Interaction between two SPC molecules */
730 /* VdW only - only first atoms in each water interact */
731 add_j_to_nblist(vdw, jj0, bLR);
735 #ifdef DISABLE_WATERWATER_NLIST
736 /* Add entries for the three atoms - only do VdW if we need to */
739 add_j_to_nblist(coul, jj0, bLR);
743 add_j_to_nblist(vdwc, jj0, bLR);
745 add_j_to_nblist(coul, jj0+1, bLR);
746 add_j_to_nblist(coul, jj0+2, bLR);
748 /* One entry for the entire water-water interaction */
751 add_j_to_nblist(coul_ww, jj0, bLR);
755 add_j_to_nblist(vdwc_ww, jj0, bLR);
760 else if (iwater == esolTIP4P && jwater == esolTIP4P)
762 /* Interaction between two TIP4p molecules */
765 /* VdW only - only first atoms in each water interact */
766 add_j_to_nblist(vdw, jj0, bLR);
770 #ifdef DISABLE_WATERWATER_NLIST
771 /* Add entries for the four atoms - only do VdW if we need to */
774 add_j_to_nblist(vdw, jj0, bLR);
776 add_j_to_nblist(coul, jj0+1, bLR);
777 add_j_to_nblist(coul, jj0+2, bLR);
778 add_j_to_nblist(coul, jj0+3, bLR);
780 /* One entry for the entire water-water interaction */
783 add_j_to_nblist(coul_ww, jj0, bLR);
787 add_j_to_nblist(vdwc_ww, jj0, bLR);
794 /* j charge group is not water, but i is.
795 * Add entries to the water-other_atom lists; the geometry of the water
796 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
797 * so we don't care if it is SPC or TIP4P...
804 for (jj = jj0; (jj < jj1); jj++)
808 add_j_to_nblist(coul, jj, bLR);
814 for (jj = jj0; (jj < jj1); jj++)
816 if (bHaveVdW[type[jj]])
818 add_j_to_nblist(vdw, jj, bLR);
824 /* _charge_ _groups_ interact with both coulomb and LJ */
825 /* Check which atoms we should add to the lists! */
826 for (jj = jj0; (jj < jj1); jj++)
828 if (bHaveVdW[type[jj]])
832 add_j_to_nblist(vdwc, jj, bLR);
836 add_j_to_nblist(vdw, jj, bLR);
839 else if (charge[jj] != 0)
841 add_j_to_nblist(coul, jj, bLR);
848 close_i_nblist(coul);
849 close_i_nblist(vdwc);
850 #ifndef DISABLE_WATERWATER_NLIST
851 close_i_nblist(coul_ww);
852 close_i_nblist(vdwc_ww);
857 /* no solvent as i charge group */
858 /* Loop over the atoms in the i charge group */
859 for (i = 0; i < nicg; i++)
862 gid = GID(igid, jgid, ngid);
865 /* Create new i_atom for each energy group */
866 if (bDoVdW && bDoCoul)
868 new_i_nblist(vdwc, i_atom, shift, gid);
872 new_i_nblist(vdw, i_atom, shift, gid);
876 new_i_nblist(coul, i_atom, shift, gid);
878 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
879 bDoCoul_i = (bDoCoul && qi != 0);
881 if (bDoVdW_i || bDoCoul_i)
883 /* Loop over the j charge groups */
884 for (j = 0; (j < nj); j++)
888 /* Check for large charge groups */
899 /* Finally loop over the atoms in the j-charge group */
900 for (jj = jj0; jj < jj1; jj++)
902 bNotEx = NOTEXCL(bExcl, i, jj);
910 add_j_to_nblist(coul, jj, bLR);
915 if (bHaveVdW[type[jj]])
917 add_j_to_nblist(vdw, jj, bLR);
922 if (bHaveVdW[type[jj]])
926 add_j_to_nblist(vdwc, jj, bLR);
930 add_j_to_nblist(vdw, jj, bLR);
933 else if (charge[jj] != 0)
935 add_j_to_nblist(coul, jj, bLR);
943 close_i_nblist(coul);
944 close_i_nblist(vdwc);
950 /* we are doing free energy */
951 vdwc_free = &nlist[eNL_VDWQQ_FREE];
952 vdw_free = &nlist[eNL_VDW_FREE];
953 coul_free = &nlist[eNL_QQ_FREE];
954 /* Loop over the atoms in the i charge group */
955 for (i = 0; i < nicg; i++)
958 gid = GID(igid, jgid, ngid);
960 qiB = chargeB[i_atom];
962 /* Create new i_atom for each energy group */
963 if (bDoVdW && bDoCoul)
965 new_i_nblist(vdwc, i_atom, shift, gid);
969 new_i_nblist(vdw, i_atom, shift, gid);
973 new_i_nblist(coul, i_atom, shift, gid);
976 new_i_nblist(vdw_free, i_atom, shift, gid);
977 new_i_nblist(coul_free, i_atom, shift, gid);
978 new_i_nblist(vdwc_free, i_atom, shift, gid);
980 bDoVdW_i = (bDoVdW &&
981 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
982 bDoCoul_i = (bDoCoul && (qi != 0 || qiB != 0));
983 /* For TIP4P the first atom does not have a charge,
984 * but the last three do. So we should still put an atom
985 * without LJ but with charge in the water-atom neighborlist
986 * for a TIP4p i charge group.
987 * For SPC type water the first atom has LJ and charge,
988 * so there is no such problem.
990 if (iwater == esolNO)
992 bDoCoul_i_sol = bDoCoul_i;
996 bDoCoul_i_sol = bDoCoul;
999 if (bDoVdW_i || bDoCoul_i_sol)
1001 /* Loop over the j charge groups */
1002 for (j = 0; (j < nj); j++)
1006 /* Check for large charge groups */
1017 /* Finally loop over the atoms in the j-charge group */
1018 bFree = bPert[i_atom];
1019 for (jj = jj0; (jj < jj1); jj++)
1021 bFreeJ = bFree || bPert[jj];
1022 /* Complicated if, because the water H's should also
1023 * see perturbed j-particles
1025 if (iwater == esolNO || i == 0 || bFreeJ)
1027 bNotEx = NOTEXCL(bExcl, i, jj);
1035 if (charge[jj] != 0 || chargeB[jj] != 0)
1037 add_j_to_nblist(coul_free, jj, bLR);
1040 else if (!bDoCoul_i)
1042 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1044 add_j_to_nblist(vdw_free, jj, bLR);
1049 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1051 if (charge[jj] != 0 || chargeB[jj] != 0)
1053 add_j_to_nblist(vdwc_free, jj, bLR);
1057 add_j_to_nblist(vdw_free, jj, bLR);
1060 else if (charge[jj] != 0 || chargeB[jj] != 0)
1062 add_j_to_nblist(coul_free, jj, bLR);
1068 /* This is done whether or not bWater is set */
1069 if (charge[jj] != 0)
1071 add_j_to_nblist(coul, jj, bLR);
1074 else if (!bDoCoul_i_sol)
1076 if (bHaveVdW[type[jj]])
1078 add_j_to_nblist(vdw, jj, bLR);
1083 if (bHaveVdW[type[jj]])
1085 if (charge[jj] != 0)
1087 add_j_to_nblist(vdwc, jj, bLR);
1091 add_j_to_nblist(vdw, jj, bLR);
1094 else if (charge[jj] != 0)
1096 add_j_to_nblist(coul, jj, bLR);
1104 close_i_nblist(vdw);
1105 close_i_nblist(coul);
1106 close_i_nblist(vdwc);
1107 close_i_nblist(vdw_free);
1108 close_i_nblist(coul_free);
1109 close_i_nblist(vdwc_free);
1115 put_in_list_adress(gmx_bool bHaveVdW[],
1131 /* The a[] index has been removed,
1132 * to put it back in i_atom should be a[i0] and jj should be a[jj].
1137 t_nblist * vdwc_adress = NULL;
1138 t_nblist * vdw_adress = NULL;
1139 t_nblist * coul_adress = NULL;
1140 t_nblist * vdwc_ww = NULL;
1141 t_nblist * coul_ww = NULL;
1143 int i, j, jcg, igid, gid, nbl_ind, nbl_ind_adress;
1144 atom_id jj, jj0, jj1, i_atom;
1149 real *charge, *chargeB;
1151 real qi, qiB, qq, rlj;
1152 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
1153 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
1155 gmx_bool j_all_atom;
1157 t_nblist *nlist, *nlist_adress;
1158 gmx_bool bEnergyGroupCG;
1160 /* Copy some pointers */
1161 cginfo = fr->cginfo;
1162 charge = md->chargeA;
1163 chargeB = md->chargeB;
1166 bPert = md->bPerturbed;
1169 /* Get atom range */
1171 nicg = index[icg+1]-i0;
1173 /* Get the i charge group info */
1174 igid = GET_CGINFO_GID(cginfo[icg]);
1176 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
1180 gmx_fatal(FARGS, "AdResS does not support free energy pertubation\n");
1183 /* Unpack pointers to neighbourlist structs */
1184 if (fr->nnblists == 2)
1191 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
1192 nbl_ind_adress = nbl_ind+fr->nnblists/2;
1196 nlist = fr->nblists[nbl_ind].nlist_lr;
1197 nlist_adress = fr->nblists[nbl_ind_adress].nlist_lr;
1201 nlist = fr->nblists[nbl_ind].nlist_sr;
1202 nlist_adress = fr->nblists[nbl_ind_adress].nlist_sr;
1206 vdwc = &nlist[eNL_VDWQQ];
1207 vdw = &nlist[eNL_VDW];
1208 coul = &nlist[eNL_QQ];
1210 vdwc_adress = &nlist_adress[eNL_VDWQQ];
1211 vdw_adress = &nlist_adress[eNL_VDW];
1212 coul_adress = &nlist_adress[eNL_QQ];
1214 /* We do not support solvent optimization with AdResS for now.
1215 For this we would need hybrid solvent-other kernels */
1217 /* no solvent as i charge group */
1218 /* Loop over the atoms in the i charge group */
1219 for (i = 0; i < nicg; i++)
1222 gid = GID(igid, jgid, ngid);
1223 qi = charge[i_atom];
1225 /* Create new i_atom for each energy group */
1226 if (bDoVdW && bDoCoul)
1228 new_i_nblist(vdwc, i_atom, shift, gid);
1229 new_i_nblist(vdwc_adress, i_atom, shift, gid);
1234 new_i_nblist(vdw, i_atom, shift, gid);
1235 new_i_nblist(vdw_adress, i_atom, shift, gid);
1240 new_i_nblist(coul, i_atom, shift, gid);
1241 new_i_nblist(coul_adress, i_atom, shift, gid);
1243 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
1244 bDoCoul_i = (bDoCoul && qi != 0);
1246 /* Here we find out whether the energy groups interaction belong to a
1247 * coarse-grained (vsite) or atomistic interaction. Note that, beacuse
1248 * interactions between coarse-grained and other (atomistic) energygroups
1249 * are excluded automatically by grompp, it is sufficient to check for
1250 * the group id of atom i (igid) */
1251 bEnergyGroupCG = !egp_explicit(fr, igid);
1253 if (bDoVdW_i || bDoCoul_i)
1255 /* Loop over the j charge groups */
1256 for (j = 0; (j < nj); j++)
1260 /* Check for large charge groups */
1271 /* Finally loop over the atoms in the j-charge group */
1272 for (jj = jj0; jj < jj1; jj++)
1274 bNotEx = NOTEXCL(bExcl, i, jj);
1276 /* Now we have to exclude interactions which will be zero
1277 * anyway due to the AdResS weights (in previous implementations
1278 * this was done in the force kernel). This is necessary as
1279 * pure interactions (those with b_hybrid=false, i.e. w_i*w_j==1 or 0)
1280 * are put into neighbour lists which will be passed to the
1281 * standard (optimized) kernels for speed. The interactions with
1282 * b_hybrid=true are placed into the _adress neighbour lists and
1283 * processed by the generic AdResS kernel.
1285 if ( (bEnergyGroupCG &&
1286 wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS ) ||
1287 ( !bEnergyGroupCG && wf[jj] <= GMX_REAL_EPS ) )
1292 b_hybrid = !((wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS) ||
1293 (wf[i_atom] <= GMX_REAL_EPS && wf[jj] <= GMX_REAL_EPS));
1299 if (charge[jj] != 0)
1303 add_j_to_nblist(coul, jj, bLR);
1307 add_j_to_nblist(coul_adress, jj, bLR);
1311 else if (!bDoCoul_i)
1313 if (bHaveVdW[type[jj]])
1317 add_j_to_nblist(vdw, jj, bLR);
1321 add_j_to_nblist(vdw_adress, jj, bLR);
1327 if (bHaveVdW[type[jj]])
1329 if (charge[jj] != 0)
1333 add_j_to_nblist(vdwc, jj, bLR);
1337 add_j_to_nblist(vdwc_adress, jj, bLR);
1344 add_j_to_nblist(vdw, jj, bLR);
1348 add_j_to_nblist(vdw_adress, jj, bLR);
1353 else if (charge[jj] != 0)
1357 add_j_to_nblist(coul, jj, bLR);
1361 add_j_to_nblist(coul_adress, jj, bLR);
1370 close_i_nblist(vdw);
1371 close_i_nblist(coul);
1372 close_i_nblist(vdwc);
1373 close_i_nblist(vdw_adress);
1374 close_i_nblist(coul_adress);
1375 close_i_nblist(vdwc_adress);
1381 put_in_list_qmmm(gmx_bool gmx_unused bHaveVdW[],
1383 t_mdatoms gmx_unused * md,
1393 gmx_bool gmx_unused bDoVdW,
1394 gmx_bool gmx_unused bDoCoul,
1395 int gmx_unused solvent_opt)
1398 int i, j, jcg, igid, gid;
1399 atom_id jj, jj0, jj1, i_atom;
1403 /* Get atom range */
1405 nicg = index[icg+1]-i0;
1407 /* Get the i charge group info */
1408 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1410 coul = &fr->QMMMlist;
1412 /* Loop over atoms in the ith charge group */
1413 for (i = 0; i < nicg; i++)
1416 gid = GID(igid, jgid, ngid);
1417 /* Create new i_atom for each energy group */
1418 new_i_nblist(coul, i_atom, shift, gid);
1420 /* Loop over the j charge groups */
1421 for (j = 0; j < nj; j++)
1425 /* Charge groups cannot have QM and MM atoms simultaneously */
1430 /* Finally loop over the atoms in the j-charge group */
1431 for (jj = jj0; jj < jj1; jj++)
1433 bNotEx = NOTEXCL(bExcl, i, jj);
1436 add_j_to_nblist(coul, jj, bLR);
1441 close_i_nblist(coul);
1446 put_in_list_cg(gmx_bool gmx_unused bHaveVdW[],
1448 t_mdatoms gmx_unused * md,
1458 gmx_bool gmx_unused bDoVdW,
1459 gmx_bool gmx_unused bDoCoul,
1460 int gmx_unused solvent_opt)
1463 int igid, gid, nbl_ind;
1467 cginfo = fr->cginfo[icg];
1469 igid = GET_CGINFO_GID(cginfo);
1470 gid = GID(igid, jgid, ngid);
1472 /* Unpack pointers to neighbourlist structs */
1473 if (fr->nnblists == 1)
1479 nbl_ind = fr->gid2nblists[gid];
1483 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1487 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1490 /* Make a new neighbor list for charge group icg.
1491 * Currently simply one neighbor list is made with LJ and Coulomb.
1492 * If required, zero interactions could be removed here
1493 * or in the force loop.
1495 new_i_nblist(vdwc, index[icg], shift, gid);
1496 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1498 for (j = 0; (j < nj); j++)
1501 /* Skip the icg-icg pairs if all self interactions are excluded */
1502 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1504 /* Here we add the j charge group jcg to the list,
1505 * exclusions are also added to the list.
1507 add_j_to_nblist_cg(vdwc, index[jcg], index[jcg+1], bExcl, icg == jcg, bLR);
1511 close_i_nblist(vdwc);
1514 static void setexcl(atom_id start, atom_id end, t_blocka *excl, gmx_bool b,
1521 for (i = start; i < end; i++)
1523 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1525 SETEXCL(bexcl, i-start, excl->a[k]);
1531 for (i = start; i < end; i++)
1533 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1535 RMEXCL(bexcl, i-start, excl->a[k]);
1541 int calc_naaj(int icg, int cgtot)
1545 if ((cgtot % 2) == 1)
1547 /* Odd number of charge groups, easy */
1548 naaj = 1 + (cgtot/2);
1550 else if ((cgtot % 4) == 0)
1552 /* Multiple of four is hard */
1589 fprintf(log, "naaj=%d\n", naaj);
1595 /************************************************
1597 * S I M P L E C O R E S T U F F
1599 ************************************************/
1601 static real calc_image_tric(rvec xi, rvec xj, matrix box,
1602 rvec b_inv, int *shift)
1604 /* This code assumes that the cut-off is smaller than
1605 * a half times the smallest diagonal element of the box.
1607 const real h25 = 2.5;
1612 /* Compute diff vector */
1613 dz = xj[ZZ] - xi[ZZ];
1614 dy = xj[YY] - xi[YY];
1615 dx = xj[XX] - xi[XX];
1617 /* Perform NINT operation, using trunc operation, therefore
1618 * we first add 2.5 then subtract 2 again
1620 tz = dz*b_inv[ZZ] + h25;
1622 dz -= tz*box[ZZ][ZZ];
1623 dy -= tz*box[ZZ][YY];
1624 dx -= tz*box[ZZ][XX];
1626 ty = dy*b_inv[YY] + h25;
1628 dy -= ty*box[YY][YY];
1629 dx -= ty*box[YY][XX];
1631 tx = dx*b_inv[XX]+h25;
1633 dx -= tx*box[XX][XX];
1635 /* Distance squared */
1636 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1638 *shift = XYZ2IS(tx, ty, tz);
1643 static real calc_image_rect(rvec xi, rvec xj, rvec box_size,
1644 rvec b_inv, int *shift)
1646 const real h15 = 1.5;
1652 /* Compute diff vector */
1653 dx = xj[XX] - xi[XX];
1654 dy = xj[YY] - xi[YY];
1655 dz = xj[ZZ] - xi[ZZ];
1657 /* Perform NINT operation, using trunc operation, therefore
1658 * we first add 1.5 then subtract 1 again
1660 tx = dx*b_inv[XX] + h15;
1661 ty = dy*b_inv[YY] + h15;
1662 tz = dz*b_inv[ZZ] + h15;
1667 /* Correct diff vector for translation */
1668 ddx = tx*box_size[XX] - dx;
1669 ddy = ty*box_size[YY] - dy;
1670 ddz = tz*box_size[ZZ] - dz;
1672 /* Distance squared */
1673 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1675 *shift = XYZ2IS(tx, ty, tz);
1680 static void add_simple(t_ns_buf *nsbuf, int nrj, atom_id cg_j,
1681 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1682 int icg, int jgid, t_block *cgs, t_excl bexcl[],
1683 int shift, t_forcerec *fr, put_in_list_t *put_in_list)
1685 if (nsbuf->nj + nrj > MAX_CG)
1687 put_in_list(bHaveVdW, ngid, md, icg, jgid, nsbuf->ncg, nsbuf->jcg,
1688 cgs->index, bexcl, shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1689 /* Reset buffer contents */
1690 nsbuf->ncg = nsbuf->nj = 0;
1692 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1696 static void ns_inner_tric(rvec x[], int icg, int *i_egp_flags,
1697 int njcg, atom_id jcg[],
1698 matrix box, rvec b_inv, real rcut2,
1699 t_block *cgs, t_ns_buf **ns_buf,
1700 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1701 t_excl bexcl[], t_forcerec *fr,
1702 put_in_list_t *put_in_list)
1706 int *cginfo = fr->cginfo;
1707 atom_id cg_j, *cgindex;
1710 cgindex = cgs->index;
1712 for (j = 0; (j < njcg); j++)
1715 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1716 if (calc_image_tric(x[icg], x[cg_j], box, b_inv, &shift) < rcut2)
1718 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1719 if (!(i_egp_flags[jgid] & EGP_EXCL))
1721 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1722 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1729 static void ns_inner_rect(rvec x[], int icg, int *i_egp_flags,
1730 int njcg, atom_id jcg[],
1731 gmx_bool bBox, rvec box_size, rvec b_inv, real rcut2,
1732 t_block *cgs, t_ns_buf **ns_buf,
1733 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1734 t_excl bexcl[], t_forcerec *fr,
1735 put_in_list_t *put_in_list)
1739 int *cginfo = fr->cginfo;
1740 atom_id cg_j, *cgindex;
1743 cgindex = cgs->index;
1747 for (j = 0; (j < njcg); j++)
1750 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1751 if (calc_image_rect(x[icg], x[cg_j], box_size, b_inv, &shift) < rcut2)
1753 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1754 if (!(i_egp_flags[jgid] & EGP_EXCL))
1756 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1757 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1765 for (j = 0; (j < njcg); j++)
1768 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1769 if ((rcut2 == 0) || (distance2(x[icg], x[cg_j]) < rcut2))
1771 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1772 if (!(i_egp_flags[jgid] & EGP_EXCL))
1774 add_simple(&ns_buf[jgid][CENTRAL], nrj, cg_j,
1775 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, CENTRAL, fr,
1783 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1785 static int ns_simple_core(t_forcerec *fr,
1786 gmx_localtop_t *top,
1788 matrix box, rvec box_size,
1789 t_excl bexcl[], atom_id *aaj,
1790 int ngid, t_ns_buf **ns_buf,
1791 put_in_list_t *put_in_list, gmx_bool bHaveVdW[])
1795 int nsearch, icg, jcg, igid, i0, nri, nn;
1798 /* atom_id *i_atoms; */
1799 t_block *cgs = &(top->cgs);
1800 t_blocka *excl = &(top->excls);
1803 gmx_bool bBox, bTriclinic;
1806 rlist2 = sqr(fr->rlist);
1808 bBox = (fr->ePBC != epbcNONE);
1811 for (m = 0; (m < DIM); m++)
1813 b_inv[m] = divide_err(1.0, box_size[m]);
1815 bTriclinic = TRICLINIC(box);
1822 cginfo = fr->cginfo;
1825 for (icg = fr->cg0; (icg < fr->hcg); icg++)
1828 i0 = cgs->index[icg];
1829 nri = cgs->index[icg+1]-i0;
1830 i_atoms = &(cgs->a[i0]);
1831 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1832 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1834 igid = GET_CGINFO_GID(cginfo[icg]);
1835 i_egp_flags = fr->egp_flags + ngid*igid;
1836 setexcl(cgs->index[icg], cgs->index[icg+1], excl, TRUE, bexcl);
1838 naaj = calc_naaj(icg, cgs->nr);
1841 ns_inner_tric(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1842 box, b_inv, rlist2, cgs, ns_buf,
1843 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1847 ns_inner_rect(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1848 bBox, box_size, b_inv, rlist2, cgs, ns_buf,
1849 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1853 for (nn = 0; (nn < ngid); nn++)
1855 for (k = 0; (k < SHIFTS); k++)
1857 nsbuf = &(ns_buf[nn][k]);
1860 put_in_list(bHaveVdW, ngid, md, icg, nn, nsbuf->ncg, nsbuf->jcg,
1861 cgs->index, bexcl, k, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1862 nsbuf->ncg = nsbuf->nj = 0;
1866 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1867 setexcl(cgs->index[icg], cgs->index[icg+1], excl, FALSE, bexcl);
1869 close_neighbor_lists(fr, FALSE);
1874 /************************************************
1876 * N S 5 G R I D S T U F F
1878 ************************************************/
1880 static inline void get_dx(int Nx, real gridx, real rc2, int xgi, real x,
1881 int *dx0, int *dx1, real *dcx2)
1909 for (i = xgi0; i >= 0; i--)
1911 dcx = (i+1)*gridx-x;
1920 for (i = xgi1; i < Nx; i++)
1933 static inline void get_dx_dd(int Nx, real gridx, real rc2, int xgi, real x,
1934 int ncpddc, int shift_min, int shift_max,
1935 int *g0, int *g1, real *dcx2)
1938 int g_min, g_max, shift_home;
1971 g_min = (shift_min == shift_home ? 0 : ncpddc);
1972 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1979 else if (shift_max < 0)
1994 /* Check one grid cell down */
1995 dcx = ((*g0 - 1) + 1)*gridx - x;
2007 /* Check one grid cell up */
2008 dcx = (*g1 + 1)*gridx - x;
2020 #define sqr(x) ((x)*(x))
2021 #define calc_dx2(XI, YI, ZI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
2022 #define calc_cyl_dx2(XI, YI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
2023 /****************************************************
2025 * F A S T N E I G H B O R S E A R C H I N G
2027 * Optimized neighboursearching routine using grid
2028 * at least 1x1x1, see GROMACS manual
2030 ****************************************************/
2033 static void get_cutoff2(t_forcerec *fr, gmx_bool bDoLongRange,
2034 real *rvdw2, real *rcoul2,
2035 real *rs2, real *rm2, real *rl2)
2037 *rs2 = sqr(fr->rlist);
2039 if (bDoLongRange && fr->bTwinRange)
2041 /* With plain cut-off or RF we need to make the list exactly
2042 * up to the cut-off and the cut-off's can be different,
2043 * so we can not simply set them to rlistlong.
2044 * To keep this code compatible with (exotic) old cases,
2045 * we also create lists up to rvdw/rcoulomb for PME and Ewald.
2046 * The interaction check should correspond to:
2047 * !ir_vdw/coulomb_might_be_zero_at_cutoff from inputrec.c.
2049 if (((fr->vdwtype == evdwCUT || fr->vdwtype == evdwPME) &&
2050 fr->vdw_modifier == eintmodNONE) ||
2051 fr->rvdw <= fr->rlist)
2053 *rvdw2 = sqr(fr->rvdw);
2057 *rvdw2 = sqr(fr->rlistlong);
2059 if (((fr->eeltype == eelCUT ||
2060 (EEL_RF(fr->eeltype) && fr->eeltype != eelRF_ZERO) ||
2061 fr->eeltype == eelPME ||
2062 fr->eeltype == eelEWALD) &&
2063 fr->coulomb_modifier == eintmodNONE) ||
2064 fr->rcoulomb <= fr->rlist)
2066 *rcoul2 = sqr(fr->rcoulomb);
2070 *rcoul2 = sqr(fr->rlistlong);
2075 /* Workaround for a gcc -O3 or -ffast-math problem */
2079 *rm2 = min(*rvdw2, *rcoul2);
2080 *rl2 = max(*rvdw2, *rcoul2);
2083 static void init_nsgrid_lists(t_forcerec *fr, int ngid, gmx_ns_t *ns)
2085 real rvdw2, rcoul2, rs2, rm2, rl2;
2088 get_cutoff2(fr, TRUE, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2090 /* Short range buffers */
2091 snew(ns->nl_sr, ngid);
2093 snew(ns->nsr, ngid);
2094 snew(ns->nlr_ljc, ngid);
2095 snew(ns->nlr_one, ngid);
2097 /* Always allocate both list types, since rcoulomb might now change with PME load balancing */
2098 /* Long range VdW and Coul buffers */
2099 snew(ns->nl_lr_ljc, ngid);
2100 /* Long range VdW or Coul only buffers */
2101 snew(ns->nl_lr_one, ngid);
2103 for (j = 0; (j < ngid); j++)
2105 snew(ns->nl_sr[j], MAX_CG);
2106 snew(ns->nl_lr_ljc[j], MAX_CG);
2107 snew(ns->nl_lr_one[j], MAX_CG);
2112 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
2117 static int nsgrid_core(t_commrec *cr, t_forcerec *fr,
2118 matrix box, int ngid,
2119 gmx_localtop_t *top,
2121 t_excl bexcl[], gmx_bool *bExcludeAlleg,
2123 put_in_list_t *put_in_list,
2124 gmx_bool bHaveVdW[],
2125 gmx_bool bDoLongRange, gmx_bool bMakeQMMMnblist)
2128 atom_id **nl_lr_ljc, **nl_lr_one, **nl_sr;
2129 int *nlr_ljc, *nlr_one, *nsr;
2130 gmx_domdec_t *dd = NULL;
2131 t_block *cgs = &(top->cgs);
2132 int *cginfo = fr->cginfo;
2133 /* atom_id *i_atoms,*cgsindex=cgs->index; */
2135 int cell_x, cell_y, cell_z;
2136 int d, tx, ty, tz, dx, dy, dz, cj;
2137 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2138 int zsh_ty, zsh_tx, ysh_tx;
2140 int dx0, dx1, dy0, dy1, dz0, dz1;
2141 int Nx, Ny, Nz, shift = -1, j, nrj, nns, nn = -1;
2142 real gridx, gridy, gridz, grid_x, grid_y, grid_z;
2143 real *dcx2, *dcy2, *dcz2;
2145 int cg0, cg1, icg = -1, cgsnr, i0, igid, nri, naaj, max_jcg;
2146 int jcg0, jcg1, jjcg, cgj0, jgid;
2147 int *grida, *gridnra, *gridind;
2148 gmx_bool rvdw_lt_rcoul, rcoul_lt_rvdw;
2149 rvec xi, *cgcm, grid_offset;
2150 real r2, rs2, rvdw2, rcoul2, rm2, rl2, XI, YI, ZI, dcx, dcy, dcz, tmp1, tmp2;
2152 gmx_bool bDomDec, bTriclinicX, bTriclinicY;
2157 bDomDec = DOMAINDECOMP(cr);
2163 bTriclinicX = ((YY < grid->npbcdim &&
2164 (!bDomDec || dd->nc[YY] == 1) && box[YY][XX] != 0) ||
2165 (ZZ < grid->npbcdim &&
2166 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][XX] != 0));
2167 bTriclinicY = (ZZ < grid->npbcdim &&
2168 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][YY] != 0);
2172 get_cutoff2(fr, bDoLongRange, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2174 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
2175 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
2177 if (bMakeQMMMnblist)
2185 nl_lr_ljc = ns->nl_lr_ljc;
2186 nl_lr_one = ns->nl_lr_one;
2187 nlr_ljc = ns->nlr_ljc;
2188 nlr_one = ns->nlr_one;
2196 gridind = grid->index;
2197 gridnra = grid->nra;
2200 gridx = grid->cell_size[XX];
2201 gridy = grid->cell_size[YY];
2202 gridz = grid->cell_size[ZZ];
2206 copy_rvec(grid->cell_offset, grid_offset);
2207 copy_ivec(grid->ncpddc, ncpddc);
2212 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2213 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
2214 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
2215 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
2216 if (zsh_tx != 0 && ysh_tx != 0)
2218 /* This could happen due to rounding, when both ratios are 0.5 */
2227 /* We only want a list for the test particle */
2236 /* Set the shift range */
2237 for (d = 0; d < DIM; d++)
2241 /* Check if we need periodicity shifts.
2242 * Without PBC or with domain decomposition we don't need them.
2244 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2251 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
2262 /* Loop over charge groups */
2263 for (icg = cg0; (icg < cg1); icg++)
2265 igid = GET_CGINFO_GID(cginfo[icg]);
2266 /* Skip this charge group if all energy groups are excluded! */
2267 if (bExcludeAlleg[igid])
2272 i0 = cgs->index[icg];
2274 if (bMakeQMMMnblist)
2276 /* Skip this charge group if it is not a QM atom while making a
2277 * QM/MM neighbourlist
2279 if (md->bQM[i0] == FALSE)
2281 continue; /* MM particle, go to next particle */
2284 /* Compute the number of charge groups that fall within the control
2287 naaj = calc_naaj(icg, cgsnr);
2294 /* make a normal neighbourlist */
2298 /* Get the j charge-group and dd cell shift ranges */
2299 dd_get_ns_ranges(cr->dd, icg, &jcg0, &jcg1, sh0, sh1);
2304 /* Compute the number of charge groups that fall within the control
2307 naaj = calc_naaj(icg, cgsnr);
2313 /* The i-particle is awlways the test particle,
2314 * so we want all j-particles
2316 max_jcg = cgsnr - 1;
2320 max_jcg = jcg1 - cgsnr;
2325 i_egp_flags = fr->egp_flags + igid*ngid;
2327 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2328 setexcl(i0, cgs->index[icg+1], &top->excls, TRUE, bexcl);
2330 ci2xyz(grid, icg, &cell_x, &cell_y, &cell_z);
2332 /* Changed iicg to icg, DvdS 990115
2333 * (but see consistency check above, DvdS 990330)
2336 fprintf(log, "icg=%5d, naaj=%5d, cell %d %d %d\n",
2337 icg, naaj, cell_x, cell_y, cell_z);
2339 /* Loop over shift vectors in three dimensions */
2340 for (tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
2342 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2343 /* Calculate range of cells in Z direction that have the shift tz */
2344 zgi = cell_z + tz*Nz;
2347 get_dx(Nz, gridz, rl2, zgi, ZI, &dz0, &dz1, dcz2);
2349 get_dx_dd(Nz, gridz, rl2, zgi, ZI-grid_offset[ZZ],
2350 ncpddc[ZZ], sh0[ZZ], sh1[ZZ], &dz0, &dz1, dcz2);
2356 for (ty = -shp[YY]; ty <= shp[YY]; ty++)
2358 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2359 /* Calculate range of cells in Y direction that have the shift ty */
2362 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2366 ygi = cell_y + ty*Ny;
2369 get_dx(Ny, gridy, rl2, ygi, YI, &dy0, &dy1, dcy2);
2371 get_dx_dd(Ny, gridy, rl2, ygi, YI-grid_offset[YY],
2372 ncpddc[YY], sh0[YY], sh1[YY], &dy0, &dy1, dcy2);
2378 for (tx = -shp[XX]; tx <= shp[XX]; tx++)
2380 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2381 /* Calculate range of cells in X direction that have the shift tx */
2384 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2388 xgi = cell_x + tx*Nx;
2391 get_dx(Nx, gridx, rl2, xgi*Nx, XI, &dx0, &dx1, dcx2);
2393 get_dx_dd(Nx, gridx, rl2, xgi, XI-grid_offset[XX],
2394 ncpddc[XX], sh0[XX], sh1[XX], &dx0, &dx1, dcx2);
2400 /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2401 * from the neigbour list as it will not interact */
2402 if (fr->adress_type != eAdressOff)
2404 if (md->wf[cgs->index[icg]] <= GMX_REAL_EPS && egp_explicit(fr, igid))
2409 /* Get shift vector */
2410 shift = XYZ2IS(tx, ty, tz);
2412 range_check(shift, 0, SHIFTS);
2414 for (nn = 0; (nn < ngid); nn++)
2421 fprintf(log, "shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2422 shift, dx0, dx1, dy0, dy1, dz0, dz1);
2423 fprintf(log, "cgcm: %8.3f %8.3f %8.3f\n", cgcm[icg][XX],
2424 cgcm[icg][YY], cgcm[icg][ZZ]);
2425 fprintf(log, "xi: %8.3f %8.3f %8.3f\n", XI, YI, ZI);
2427 for (dx = dx0; (dx <= dx1); dx++)
2429 tmp1 = rl2 - dcx2[dx];
2430 for (dy = dy0; (dy <= dy1); dy++)
2432 tmp2 = tmp1 - dcy2[dy];
2435 for (dz = dz0; (dz <= dz1); dz++)
2437 if (tmp2 > dcz2[dz])
2439 /* Find grid-cell cj in which possible neighbours are */
2440 cj = xyz2ci(Ny, Nz, dx, dy, dz);
2442 /* Check out how many cgs (nrj) there in this cell */
2445 /* Find the offset in the cg list */
2448 /* Check if all j's are out of range so we
2449 * can skip the whole cell.
2450 * Should save some time, especially with DD.
2453 (grida[cgj0] >= max_jcg &&
2454 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2460 for (j = 0; (j < nrj); j++)
2462 jjcg = grida[cgj0+j];
2464 /* check whether this guy is in range! */
2465 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2468 r2 = calc_dx2(XI, YI, ZI, cgcm[jjcg]);
2471 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2472 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2473 /* check energy group exclusions */
2474 if (!(i_egp_flags[jgid] & EGP_EXCL))
2478 if (nsr[jgid] >= MAX_CG)
2480 /* Add to short-range list */
2481 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2482 nsr[jgid], nl_sr[jgid],
2483 cgs->index, /* cgsatoms, */ bexcl,
2484 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2487 nl_sr[jgid][nsr[jgid]++] = jjcg;
2491 if (nlr_ljc[jgid] >= MAX_CG)
2493 /* Add to LJ+coulomb long-range list */
2494 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2495 nlr_ljc[jgid], nl_lr_ljc[jgid], top->cgs.index,
2496 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2499 nl_lr_ljc[jgid][nlr_ljc[jgid]++] = jjcg;
2503 if (nlr_one[jgid] >= MAX_CG)
2505 /* Add to long-range list with only coul, or only LJ */
2506 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2507 nlr_one[jgid], nl_lr_one[jgid], top->cgs.index,
2508 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2511 nl_lr_one[jgid][nlr_one[jgid]++] = jjcg;
2523 /* CHECK whether there is anything left in the buffers */
2524 for (nn = 0; (nn < ngid); nn++)
2528 put_in_list(bHaveVdW, ngid, md, icg, nn, nsr[nn], nl_sr[nn],
2529 cgs->index, /* cgsatoms, */ bexcl,
2530 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2533 if (nlr_ljc[nn] > 0)
2535 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_ljc[nn],
2536 nl_lr_ljc[nn], top->cgs.index,
2537 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2540 if (nlr_one[nn] > 0)
2542 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_one[nn],
2543 nl_lr_one[nn], top->cgs.index,
2544 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2550 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2551 setexcl(cgs->index[icg], cgs->index[icg+1], &top->excls, FALSE, bexcl);
2553 /* No need to perform any left-over force calculations anymore (as we used to do here)
2554 * since we now save the proper long-range lists for later evaluation.
2559 /* Close neighbourlists */
2560 close_neighbor_lists(fr, bMakeQMMMnblist);
2565 void ns_realloc_natoms(gmx_ns_t *ns, int natoms)
2569 if (natoms > ns->nra_alloc)
2571 ns->nra_alloc = over_alloc_dd(natoms);
2572 srenew(ns->bexcl, ns->nra_alloc);
2573 for (i = 0; i < ns->nra_alloc; i++)
2580 void init_ns(FILE *fplog, const t_commrec *cr,
2581 gmx_ns_t *ns, t_forcerec *fr,
2582 const gmx_mtop_t *mtop)
2584 int mt, icg, nr_in_cg, maxcg, i, j, jcg, ngid, ncg;
2588 /* Compute largest charge groups size (# atoms) */
2590 for (mt = 0; mt < mtop->nmoltype; mt++)
2592 cgs = &mtop->moltype[mt].cgs;
2593 for (icg = 0; (icg < cgs->nr); icg++)
2595 nr_in_cg = max(nr_in_cg, (int)(cgs->index[icg+1]-cgs->index[icg]));
2599 /* Verify whether largest charge group is <= max cg.
2600 * This is determined by the type of the local exclusion type
2601 * Exclusions are stored in bits. (If the type is not large
2602 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2604 maxcg = sizeof(t_excl)*8;
2605 if (nr_in_cg > maxcg)
2607 gmx_fatal(FARGS, "Max #atoms in a charge group: %d > %d\n",
2611 ngid = mtop->groups.grps[egcENER].nr;
2612 snew(ns->bExcludeAlleg, ngid);
2613 for (i = 0; i < ngid; i++)
2615 ns->bExcludeAlleg[i] = TRUE;
2616 for (j = 0; j < ngid; j++)
2618 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2620 ns->bExcludeAlleg[i] = FALSE;
2628 ns->grid = init_grid(fplog, fr);
2629 init_nsgrid_lists(fr, ngid, ns);
2634 snew(ns->ns_buf, ngid);
2635 for (i = 0; (i < ngid); i++)
2637 snew(ns->ns_buf[i], SHIFTS);
2639 ncg = ncg_mtop(mtop);
2640 snew(ns->simple_aaj, 2*ncg);
2641 for (jcg = 0; (jcg < ncg); jcg++)
2643 ns->simple_aaj[jcg] = jcg;
2644 ns->simple_aaj[jcg+ncg] = jcg;
2648 /* Create array that determines whether or not atoms have VdW */
2649 snew(ns->bHaveVdW, fr->ntype);
2650 for (i = 0; (i < fr->ntype); i++)
2652 for (j = 0; (j < fr->ntype); j++)
2654 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2656 ((BHAMA(fr->nbfp, fr->ntype, i, j) != 0) ||
2657 (BHAMB(fr->nbfp, fr->ntype, i, j) != 0) ||
2658 (BHAMC(fr->nbfp, fr->ntype, i, j) != 0)) :
2659 ((C6(fr->nbfp, fr->ntype, i, j) != 0) ||
2660 (C12(fr->nbfp, fr->ntype, i, j) != 0))));
2665 pr_bvec(debug, 0, "bHaveVdW", ns->bHaveVdW, fr->ntype, TRUE);
2670 if (!DOMAINDECOMP(cr))
2672 ns_realloc_natoms(ns, mtop->natoms);
2675 ns->nblist_initialized = FALSE;
2677 /* nbr list debug dump */
2679 char *ptr = getenv("GMX_DUMP_NL");
2682 ns->dump_nl = strtol(ptr, NULL, 10);
2685 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2696 int search_neighbours(FILE *log, t_forcerec *fr,
2698 gmx_localtop_t *top,
2699 gmx_groups_t *groups,
2701 t_nrnb *nrnb, t_mdatoms *md,
2703 gmx_bool bDoLongRangeNS)
2705 t_block *cgs = &(top->cgs);
2706 rvec box_size, grid_x0, grid_x1;
2708 real min_size, grid_dens;
2712 gmx_bool *i_egp_flags;
2713 int cg_start, cg_end, start, end;
2716 gmx_domdec_zones_t *dd_zones;
2717 put_in_list_t *put_in_list;
2721 /* Set some local variables */
2723 ngid = groups->grps[egcENER].nr;
2725 for (m = 0; (m < DIM); m++)
2727 box_size[m] = box[m][m];
2730 if (fr->ePBC != epbcNONE)
2732 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC, box))
2734 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.");
2738 min_size = min(box_size[XX], min(box_size[YY], box_size[ZZ]));
2739 if (2*fr->rlistlong >= min_size)
2741 gmx_fatal(FARGS, "One of the box diagonal elements has become smaller than twice the cut-off length.");
2746 if (DOMAINDECOMP(cr))
2748 ns_realloc_natoms(ns, cgs->index[cgs->nr]);
2752 /* Reset the neighbourlists */
2753 reset_neighbor_lists(fr, TRUE, TRUE);
2755 if (bGrid && bFillGrid)
2759 if (DOMAINDECOMP(cr))
2761 dd_zones = domdec_zones(cr->dd);
2767 get_nsgrid_boundaries(grid->nboundeddim, box, NULL, NULL, NULL, NULL,
2768 cgs->nr, fr->cg_cm, grid_x0, grid_x1, &grid_dens);
2770 grid_first(log, grid, NULL, NULL, box, grid_x0, grid_x1,
2771 fr->rlistlong, grid_dens);
2778 if (DOMAINDECOMP(cr))
2781 fill_grid(dd_zones, grid, end, -1, end, fr->cg_cm);
2783 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2787 fill_grid(NULL, grid, cgs->nr, fr->cg0, fr->hcg, fr->cg_cm);
2788 grid->icg0 = fr->cg0;
2789 grid->icg1 = fr->hcg;
2793 calc_elemnr(grid, start, end, cgs->nr);
2795 grid_last(grid, start, end, cgs->nr);
2800 print_grid(debug, grid);
2805 /* Set the grid cell index for the test particle only.
2806 * The cell to cg index is not corrected, but that does not matter.
2808 fill_grid(NULL, ns->grid, fr->hcg, fr->hcg-1, fr->hcg, fr->cg_cm);
2812 if (fr->adress_type == eAdressOff)
2814 if (!fr->ns.bCGlist)
2816 put_in_list = put_in_list_at;
2820 put_in_list = put_in_list_cg;
2825 put_in_list = put_in_list_adress;
2832 nsearch = nsgrid_core(cr, fr, box, ngid, top,
2833 grid, ns->bexcl, ns->bExcludeAlleg,
2834 md, put_in_list, ns->bHaveVdW,
2835 bDoLongRangeNS, FALSE);
2837 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2838 * the classical calculation. The charge-charge interaction
2839 * between QM and MM atoms is handled in the QMMM core calculation
2840 * (see QMMM.c). The VDW however, we'd like to compute classically
2841 * and the QM MM atom pairs have just been put in the
2842 * corresponding neighbourlists. in case of QMMM we still need to
2843 * fill a special QMMM neighbourlist that contains all neighbours
2844 * of the QM atoms. If bQMMM is true, this list will now be made:
2846 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
2848 nsearch += nsgrid_core(cr, fr, box, ngid, top,
2849 grid, ns->bexcl, ns->bExcludeAlleg,
2850 md, put_in_list_qmmm, ns->bHaveVdW,
2851 bDoLongRangeNS, TRUE);
2856 nsearch = ns_simple_core(fr, top, md, box, box_size,
2857 ns->bexcl, ns->simple_aaj,
2858 ngid, ns->ns_buf, put_in_list, ns->bHaveVdW);
2866 inc_nrnb(nrnb, eNR_NS, nsearch);
2867 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2872 int natoms_beyond_ns_buffer(t_inputrec *ir, t_forcerec *fr, t_block *cgs,
2873 matrix scale_tot, rvec *x)
2875 int cg0, cg1, cg, a0, a1, a, i, j;
2876 real rint, hbuf2, scale;
2878 gmx_bool bIsotropic;
2883 rint = max(ir->rcoulomb, ir->rvdw);
2884 if (ir->rlist < rint)
2886 gmx_fatal(FARGS, "The neighbor search buffer has negative size: %f nm",
2894 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2896 hbuf2 = sqr(0.5*(ir->rlist - rint));
2897 for (cg = cg0; cg < cg1; cg++)
2899 a0 = cgs->index[cg];
2900 a1 = cgs->index[cg+1];
2901 for (a = a0; a < a1; a++)
2903 if (distance2(cg_cm[cg], x[a]) > hbuf2)
2913 scale = scale_tot[0][0];
2914 for (i = 1; i < DIM; i++)
2916 /* With anisotropic scaling, the original spherical ns volumes become
2917 * ellipsoids. To avoid costly transformations we use the minimum
2918 * eigenvalue of the scaling matrix for determining the buffer size.
2919 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2921 scale = min(scale, scale_tot[i][i]);
2922 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2926 for (j = 0; j < i; j++)
2928 if (scale_tot[i][j] != 0)
2934 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2937 for (cg = cg0; cg < cg1; cg++)
2939 svmul(scale, cg_cm[cg], cgsc);
2940 a0 = cgs->index[cg];
2941 a1 = cgs->index[cg+1];
2942 for (a = a0; a < a1; a++)
2944 if (distance2(cgsc, x[a]) > hbuf2)
2953 /* Anistropic scaling */
2954 for (cg = cg0; cg < cg1; cg++)
2956 /* Since scale_tot contains the transpose of the scaling matrix,
2957 * we need to multiply with the transpose.
2959 tmvmul_ur0(scale_tot, cg_cm[cg], cgsc);
2960 a0 = cgs->index[cg];
2961 a1 = cgs->index[cg+1];
2962 for (a = a0; a < a1; a++)
2964 if (distance2(cgsc, x[a]) > hbuf2)