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 static 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;
179 reallocate_nblist(nl);
184 fprintf(debug, "Initiating neighbourlist (ielec=%d, ivdw=%d, type=%d) for %s interactions,\nwith %d SR, %d LR atoms.\n",
185 nl->ielec, nl->ivdw, nl->type, gmx_nblist_geometry_names[nl->igeometry], maxsr, maxlr);
190 void init_neighbor_list(FILE *log, t_forcerec *fr, int homenr)
192 /* Make maxlr tunable! (does not seem to be a big difference though)
193 * This parameter determines the number of i particles in a long range
194 * neighbourlist. Too few means many function calls, too many means
197 int maxsr, maxsr_wat, maxlr, maxlr_wat;
198 int ielec, ielecf, ivdw, ielecmod, ielecmodf, ivdwmod, type;
200 int igeometry_def, igeometry_w, igeometry_ww;
204 /* maxsr = homenr-fr->nWatMol*3; */
209 gmx_fatal(FARGS, "%s, %d: Negative number of short range atoms.\n"
210 "Call your Gromacs dealer for assistance.", __FILE__, __LINE__);
212 /* This is just for initial allocation, so we do not reallocate
213 * all the nlist arrays many times in a row.
214 * The numbers seem very accurate, but they are uncritical.
216 maxsr_wat = min(fr->nWatMol, (homenr+2)/3);
220 maxlr_wat = min(maxsr_wat, maxlr);
224 maxlr = maxlr_wat = 0;
227 /* Determine the values for ielec/ivdw. */
228 ielec = fr->nbkernel_elec_interaction;
229 ivdw = fr->nbkernel_vdw_interaction;
230 ielecmod = fr->nbkernel_elec_modifier;
231 ivdwmod = fr->nbkernel_vdw_modifier;
232 type = GMX_NBLIST_INTERACTION_STANDARD;
234 fr->ns.bCGlist = (getenv("GMX_NBLISTCG") != 0);
237 igeometry_def = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
241 igeometry_def = GMX_NBLIST_GEOMETRY_CG_CG;
244 fprintf(log, "\nUsing charge-group - charge-group neighbor lists and kernels\n\n");
248 if (fr->solvent_opt == esolTIP4P)
250 igeometry_w = GMX_NBLIST_GEOMETRY_WATER4_PARTICLE;
251 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER4_WATER4;
255 igeometry_w = GMX_NBLIST_GEOMETRY_WATER3_PARTICLE;
256 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER3_WATER3;
259 for (i = 0; i < fr->nnblists; i++)
261 nbl = &(fr->nblists[i]);
263 if ((fr->adress_type != eAdressOff) && (i >= fr->nnblists/2))
265 type = GMX_NBLIST_INTERACTION_ADRESS;
267 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ], &nbl->nlist_lr[eNL_VDWQQ],
268 maxsr, maxlr, ivdw, ivdwmod, ielec, ielecmod, igeometry_def, type);
269 init_nblist(log, &nbl->nlist_sr[eNL_VDW], &nbl->nlist_lr[eNL_VDW],
270 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, igeometry_def, type);
271 init_nblist(log, &nbl->nlist_sr[eNL_QQ], &nbl->nlist_lr[eNL_QQ],
272 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_def, type);
273 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATER], &nbl->nlist_lr[eNL_VDWQQ_WATER],
274 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_w, type);
275 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATER], &nbl->nlist_lr[eNL_QQ_WATER],
276 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_w, type);
277 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATERWATER], &nbl->nlist_lr[eNL_VDWQQ_WATERWATER],
278 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_ww, type);
279 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATERWATER], &nbl->nlist_lr[eNL_QQ_WATERWATER],
280 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_ww, type);
282 /* Did we get the solvent loops so we can use optimized water kernels? */
283 if (nbl->nlist_sr[eNL_VDWQQ_WATER].kernelptr_vf == NULL
284 || nbl->nlist_sr[eNL_QQ_WATER].kernelptr_vf == NULL
285 #ifndef DISABLE_WATERWATER_NLIST
286 || nbl->nlist_sr[eNL_VDWQQ_WATERWATER].kernelptr_vf == NULL
287 || nbl->nlist_sr[eNL_QQ_WATERWATER].kernelptr_vf == NULL
291 fr->solvent_opt = esolNO;
292 fprintf(log, "Note: The available nonbonded kernels do not support water optimization - disabling.\n");
295 if (fr->efep != efepNO)
297 if ((fr->bEwald) && (fr->sc_alphacoul > 0)) /* need to handle long range differently if using softcore */
299 ielecf = GMX_NBKERNEL_ELEC_EWALD;
300 ielecmodf = eintmodNONE;
305 ielecmodf = ielecmod;
308 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_FREE], &nbl->nlist_lr[eNL_VDWQQ_FREE],
309 maxsr, maxlr, ivdw, ivdwmod, ielecf, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
310 init_nblist(log, &nbl->nlist_sr[eNL_VDW_FREE], &nbl->nlist_lr[eNL_VDW_FREE],
311 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
312 init_nblist(log, &nbl->nlist_sr[eNL_QQ_FREE], &nbl->nlist_lr[eNL_QQ_FREE],
313 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielecf, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
317 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
319 init_nblist(log, &fr->QMMMlist, NULL,
320 maxsr, maxlr, 0, 0, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_STANDARD);
328 fr->ns.nblist_initialized = TRUE;
331 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];
438 /* nlist length for water i molecules is treated statically
441 if (len > nlist->maxlen)
448 static inline void close_nblist(t_nblist *nlist)
450 /* Only close this nblist when it has been initialized.
451 * Avoid the creation of i-lists with no j-particles.
455 /* Some assembly kernels do not support empty lists,
456 * make sure here that we don't generate any empty lists.
457 * With the current ns code this branch is taken in two cases:
458 * No i-particles at all: nri=-1 here
459 * There are i-particles, but no j-particles; nri=0 here
465 /* Close list number nri by incrementing the count */
470 static inline void close_neighbor_lists(t_forcerec *fr, gmx_bool bMakeQMMMnblist)
476 close_nblist(&(fr->QMMMlist));
479 for (n = 0; n < fr->nnblists; n++)
481 for (i = 0; (i < eNL_NR); i++)
483 close_nblist(&(fr->nblists[n].nlist_sr[i]));
484 close_nblist(&(fr->nblists[n].nlist_lr[i]));
490 static inline void add_j_to_nblist(t_nblist *nlist, atom_id j_atom, gmx_bool bLR)
492 int nrj = nlist->nrj;
494 if (nlist->nrj >= nlist->maxnrj)
496 nlist->maxnrj = round_up_to_simd_width(over_alloc_small(nlist->nrj + 1), nlist->simd_padding_width);
500 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
501 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
504 srenew(nlist->jjnr, nlist->maxnrj);
507 nlist->jjnr[nrj] = j_atom;
511 static inline void add_j_to_nblist_cg(t_nblist *nlist,
512 atom_id j_start, int j_end,
513 t_excl *bexcl, gmx_bool i_is_j,
516 int nrj = nlist->nrj;
519 if (nlist->nrj >= nlist->maxnrj)
521 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
524 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
525 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
528 srenew(nlist->jjnr, nlist->maxnrj);
529 srenew(nlist->jjnr_end, nlist->maxnrj);
530 srenew(nlist->excl, nlist->maxnrj*MAX_CGCGSIZE);
533 nlist->jjnr[nrj] = j_start;
534 nlist->jjnr_end[nrj] = j_end;
536 if (j_end - j_start > MAX_CGCGSIZE)
538 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);
541 /* Set the exclusions */
542 for (j = j_start; j < j_end; j++)
544 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
548 /* Avoid double counting of intra-cg interactions */
549 for (j = 1; j < j_end-j_start; j++)
551 nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
559 put_in_list_t (gmx_bool bHaveVdW[],
576 put_in_list_at(gmx_bool bHaveVdW[],
592 /* The a[] index has been removed,
593 * to put it back in i_atom should be a[i0] and jj should be a[jj].
598 t_nblist * vdwc_free = NULL;
599 t_nblist * vdw_free = NULL;
600 t_nblist * coul_free = NULL;
601 t_nblist * vdwc_ww = NULL;
602 t_nblist * coul_ww = NULL;
604 int i, j, jcg, igid, gid, nbl_ind, ind_ij;
605 atom_id jj, jj0, jj1, i_atom;
610 real *charge, *chargeB;
611 real qi, qiB, qq, rlj;
612 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
613 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
617 /* Copy some pointers */
619 charge = md->chargeA;
620 chargeB = md->chargeB;
623 bPert = md->bPerturbed;
627 nicg = index[icg+1]-i0;
629 /* Get the i charge group info */
630 igid = GET_CGINFO_GID(cginfo[icg]);
632 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
637 /* Check if any of the particles involved are perturbed.
638 * If not we can do the cheaper normal put_in_list
639 * and use more solvent optimization.
641 for (i = 0; i < nicg; i++)
643 bFreeEnergy |= bPert[i0+i];
645 /* Loop over the j charge groups */
646 for (j = 0; (j < nj && !bFreeEnergy); j++)
651 /* Finally loop over the atoms in the j-charge group */
652 for (jj = jj0; jj < jj1; jj++)
654 bFreeEnergy |= bPert[jj];
659 /* Unpack pointers to neighbourlist structs */
660 if (fr->nnblists == 1)
666 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
670 nlist = fr->nblists[nbl_ind].nlist_lr;
674 nlist = fr->nblists[nbl_ind].nlist_sr;
677 if (iwater != esolNO)
679 vdwc = &nlist[eNL_VDWQQ_WATER];
680 vdw = &nlist[eNL_VDW];
681 coul = &nlist[eNL_QQ_WATER];
682 #ifndef DISABLE_WATERWATER_NLIST
683 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
684 coul_ww = &nlist[eNL_QQ_WATERWATER];
689 vdwc = &nlist[eNL_VDWQQ];
690 vdw = &nlist[eNL_VDW];
691 coul = &nlist[eNL_QQ];
696 if (iwater != esolNO)
698 /* Loop over the atoms in the i charge group */
700 gid = GID(igid, jgid, ngid);
701 /* Create new i_atom for each energy group */
702 if (bDoCoul && bDoVdW)
704 new_i_nblist(vdwc, i_atom, shift, gid);
705 #ifndef DISABLE_WATERWATER_NLIST
706 new_i_nblist(vdwc_ww, i_atom, shift, gid);
711 new_i_nblist(vdw, i_atom, shift, gid);
715 new_i_nblist(coul, i_atom, shift, gid);
716 #ifndef DISABLE_WATERWATER_NLIST
717 new_i_nblist(coul_ww, i_atom, shift, gid);
720 /* Loop over the j charge groups */
721 for (j = 0; (j < nj); j++)
731 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
733 if (iwater == esolSPC && jwater == esolSPC)
735 /* Interaction between two SPC molecules */
738 /* VdW only - only first atoms in each water interact */
739 add_j_to_nblist(vdw, jj0, bLR);
743 #ifdef DISABLE_WATERWATER_NLIST
744 /* Add entries for the three atoms - only do VdW if we need to */
747 add_j_to_nblist(coul, jj0, bLR);
751 add_j_to_nblist(vdwc, jj0, bLR);
753 add_j_to_nblist(coul, jj0+1, bLR);
754 add_j_to_nblist(coul, jj0+2, bLR);
756 /* One entry for the entire water-water interaction */
759 add_j_to_nblist(coul_ww, jj0, bLR);
763 add_j_to_nblist(vdwc_ww, jj0, bLR);
768 else if (iwater == esolTIP4P && jwater == esolTIP4P)
770 /* Interaction between two TIP4p molecules */
773 /* VdW only - only first atoms in each water interact */
774 add_j_to_nblist(vdw, jj0, bLR);
778 #ifdef DISABLE_WATERWATER_NLIST
779 /* Add entries for the four atoms - only do VdW if we need to */
782 add_j_to_nblist(vdw, jj0, bLR);
784 add_j_to_nblist(coul, jj0+1, bLR);
785 add_j_to_nblist(coul, jj0+2, bLR);
786 add_j_to_nblist(coul, jj0+3, bLR);
788 /* One entry for the entire water-water interaction */
791 add_j_to_nblist(coul_ww, jj0, bLR);
795 add_j_to_nblist(vdwc_ww, jj0, bLR);
802 /* j charge group is not water, but i is.
803 * Add entries to the water-other_atom lists; the geometry of the water
804 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
805 * so we don't care if it is SPC or TIP4P...
812 for (jj = jj0; (jj < jj1); jj++)
816 add_j_to_nblist(coul, jj, bLR);
822 for (jj = jj0; (jj < jj1); jj++)
824 if (bHaveVdW[type[jj]])
826 add_j_to_nblist(vdw, jj, bLR);
832 /* _charge_ _groups_ interact with both coulomb and LJ */
833 /* Check which atoms we should add to the lists! */
834 for (jj = jj0; (jj < jj1); jj++)
836 if (bHaveVdW[type[jj]])
840 add_j_to_nblist(vdwc, jj, bLR);
844 add_j_to_nblist(vdw, jj, bLR);
847 else if (charge[jj] != 0)
849 add_j_to_nblist(coul, jj, bLR);
856 close_i_nblist(coul);
857 close_i_nblist(vdwc);
858 #ifndef DISABLE_WATERWATER_NLIST
859 close_i_nblist(coul_ww);
860 close_i_nblist(vdwc_ww);
865 /* no solvent as i charge group */
866 /* Loop over the atoms in the i charge group */
867 for (i = 0; i < nicg; i++)
870 gid = GID(igid, jgid, ngid);
873 /* Create new i_atom for each energy group */
874 if (bDoVdW && bDoCoul)
876 new_i_nblist(vdwc, i_atom, shift, gid);
880 new_i_nblist(vdw, i_atom, shift, gid);
884 new_i_nblist(coul, i_atom, shift, gid);
886 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
887 bDoCoul_i = (bDoCoul && qi != 0);
889 if (bDoVdW_i || bDoCoul_i)
891 /* Loop over the j charge groups */
892 for (j = 0; (j < nj); j++)
896 /* Check for large charge groups */
907 /* Finally loop over the atoms in the j-charge group */
908 for (jj = jj0; jj < jj1; jj++)
910 bNotEx = NOTEXCL(bExcl, i, jj);
918 add_j_to_nblist(coul, jj, bLR);
923 if (bHaveVdW[type[jj]])
925 add_j_to_nblist(vdw, jj, bLR);
930 if (bHaveVdW[type[jj]])
934 add_j_to_nblist(vdwc, jj, bLR);
938 add_j_to_nblist(vdw, jj, bLR);
941 else if (charge[jj] != 0)
943 add_j_to_nblist(coul, jj, bLR);
951 close_i_nblist(coul);
952 close_i_nblist(vdwc);
958 /* we are doing free energy */
959 vdwc_free = &nlist[eNL_VDWQQ_FREE];
960 vdw_free = &nlist[eNL_VDW_FREE];
961 coul_free = &nlist[eNL_QQ_FREE];
962 /* Loop over the atoms in the i charge group */
963 for (i = 0; i < nicg; i++)
966 gid = GID(igid, jgid, ngid);
968 qiB = chargeB[i_atom];
970 /* Create new i_atom for each energy group */
971 if (bDoVdW && bDoCoul)
973 new_i_nblist(vdwc, i_atom, shift, gid);
977 new_i_nblist(vdw, i_atom, shift, gid);
981 new_i_nblist(coul, i_atom, shift, gid);
984 new_i_nblist(vdw_free, i_atom, shift, gid);
985 new_i_nblist(coul_free, i_atom, shift, gid);
986 new_i_nblist(vdwc_free, i_atom, shift, gid);
988 bDoVdW_i = (bDoVdW &&
989 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
990 bDoCoul_i = (bDoCoul && (qi != 0 || qiB != 0));
991 /* For TIP4P the first atom does not have a charge,
992 * but the last three do. So we should still put an atom
993 * without LJ but with charge in the water-atom neighborlist
994 * for a TIP4p i charge group.
995 * For SPC type water the first atom has LJ and charge,
996 * so there is no such problem.
998 if (iwater == esolNO)
1000 bDoCoul_i_sol = bDoCoul_i;
1004 bDoCoul_i_sol = bDoCoul;
1007 if (bDoVdW_i || bDoCoul_i_sol)
1009 /* Loop over the j charge groups */
1010 for (j = 0; (j < nj); j++)
1014 /* Check for large charge groups */
1025 /* Finally loop over the atoms in the j-charge group */
1026 bFree = bPert[i_atom];
1027 for (jj = jj0; (jj < jj1); jj++)
1029 bFreeJ = bFree || bPert[jj];
1030 /* Complicated if, because the water H's should also
1031 * see perturbed j-particles
1033 if (iwater == esolNO || i == 0 || bFreeJ)
1035 bNotEx = NOTEXCL(bExcl, i, jj);
1043 if (charge[jj] != 0 || chargeB[jj] != 0)
1045 add_j_to_nblist(coul_free, jj, bLR);
1048 else if (!bDoCoul_i)
1050 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1052 add_j_to_nblist(vdw_free, jj, bLR);
1057 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1059 if (charge[jj] != 0 || chargeB[jj] != 0)
1061 add_j_to_nblist(vdwc_free, jj, bLR);
1065 add_j_to_nblist(vdw_free, jj, bLR);
1068 else if (charge[jj] != 0 || chargeB[jj] != 0)
1070 add_j_to_nblist(coul_free, jj, bLR);
1076 /* This is done whether or not bWater is set */
1077 if (charge[jj] != 0)
1079 add_j_to_nblist(coul, jj, bLR);
1082 else if (!bDoCoul_i_sol)
1084 if (bHaveVdW[type[jj]])
1086 add_j_to_nblist(vdw, jj, bLR);
1091 if (bHaveVdW[type[jj]])
1093 if (charge[jj] != 0)
1095 add_j_to_nblist(vdwc, jj, bLR);
1099 add_j_to_nblist(vdw, jj, bLR);
1102 else if (charge[jj] != 0)
1104 add_j_to_nblist(coul, jj, bLR);
1112 close_i_nblist(vdw);
1113 close_i_nblist(coul);
1114 close_i_nblist(vdwc);
1115 close_i_nblist(vdw_free);
1116 close_i_nblist(coul_free);
1117 close_i_nblist(vdwc_free);
1123 put_in_list_adress(gmx_bool bHaveVdW[],
1139 /* The a[] index has been removed,
1140 * to put it back in i_atom should be a[i0] and jj should be a[jj].
1145 t_nblist * vdwc_adress = NULL;
1146 t_nblist * vdw_adress = NULL;
1147 t_nblist * coul_adress = NULL;
1148 t_nblist * vdwc_ww = NULL;
1149 t_nblist * coul_ww = NULL;
1151 int i, j, jcg, igid, gid, nbl_ind, nbl_ind_adress;
1152 atom_id jj, jj0, jj1, i_atom;
1157 real *charge, *chargeB;
1159 real qi, qiB, qq, rlj;
1160 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
1161 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
1163 gmx_bool j_all_atom;
1165 t_nblist *nlist, *nlist_adress;
1166 gmx_bool bEnergyGroupCG;
1168 /* Copy some pointers */
1169 cginfo = fr->cginfo;
1170 charge = md->chargeA;
1171 chargeB = md->chargeB;
1174 bPert = md->bPerturbed;
1177 /* Get atom range */
1179 nicg = index[icg+1]-i0;
1181 /* Get the i charge group info */
1182 igid = GET_CGINFO_GID(cginfo[icg]);
1184 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
1188 gmx_fatal(FARGS, "AdResS does not support free energy pertubation\n");
1191 /* Unpack pointers to neighbourlist structs */
1192 if (fr->nnblists == 2)
1199 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
1200 nbl_ind_adress = nbl_ind+fr->nnblists/2;
1204 nlist = fr->nblists[nbl_ind].nlist_lr;
1205 nlist_adress = fr->nblists[nbl_ind_adress].nlist_lr;
1209 nlist = fr->nblists[nbl_ind].nlist_sr;
1210 nlist_adress = fr->nblists[nbl_ind_adress].nlist_sr;
1214 vdwc = &nlist[eNL_VDWQQ];
1215 vdw = &nlist[eNL_VDW];
1216 coul = &nlist[eNL_QQ];
1218 vdwc_adress = &nlist_adress[eNL_VDWQQ];
1219 vdw_adress = &nlist_adress[eNL_VDW];
1220 coul_adress = &nlist_adress[eNL_QQ];
1222 /* We do not support solvent optimization with AdResS for now.
1223 For this we would need hybrid solvent-other kernels */
1225 /* no solvent as i charge group */
1226 /* Loop over the atoms in the i charge group */
1227 for (i = 0; i < nicg; i++)
1230 gid = GID(igid, jgid, ngid);
1231 qi = charge[i_atom];
1233 /* Create new i_atom for each energy group */
1234 if (bDoVdW && bDoCoul)
1236 new_i_nblist(vdwc, i_atom, shift, gid);
1237 new_i_nblist(vdwc_adress, i_atom, shift, gid);
1242 new_i_nblist(vdw, i_atom, shift, gid);
1243 new_i_nblist(vdw_adress, i_atom, shift, gid);
1248 new_i_nblist(coul, i_atom, shift, gid);
1249 new_i_nblist(coul_adress, i_atom, shift, gid);
1251 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
1252 bDoCoul_i = (bDoCoul && qi != 0);
1254 /* Here we find out whether the energy groups interaction belong to a
1255 * coarse-grained (vsite) or atomistic interaction. Note that, beacuse
1256 * interactions between coarse-grained and other (atomistic) energygroups
1257 * are excluded automatically by grompp, it is sufficient to check for
1258 * the group id of atom i (igid) */
1259 bEnergyGroupCG = !egp_explicit(fr, igid);
1261 if (bDoVdW_i || bDoCoul_i)
1263 /* Loop over the j charge groups */
1264 for (j = 0; (j < nj); j++)
1268 /* Check for large charge groups */
1279 /* Finally loop over the atoms in the j-charge group */
1280 for (jj = jj0; jj < jj1; jj++)
1282 bNotEx = NOTEXCL(bExcl, i, jj);
1284 /* Now we have to exclude interactions which will be zero
1285 * anyway due to the AdResS weights (in previous implementations
1286 * this was done in the force kernel). This is necessary as
1287 * pure interactions (those with b_hybrid=false, i.e. w_i*w_j==1 or 0)
1288 * are put into neighbour lists which will be passed to the
1289 * standard (optimized) kernels for speed. The interactions with
1290 * b_hybrid=true are placed into the _adress neighbour lists and
1291 * processed by the generic AdResS kernel.
1293 if ( (bEnergyGroupCG &&
1294 wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS ) ||
1295 ( !bEnergyGroupCG && wf[jj] <= GMX_REAL_EPS ) )
1300 b_hybrid = !((wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS) ||
1301 (wf[i_atom] <= GMX_REAL_EPS && wf[jj] <= GMX_REAL_EPS));
1307 if (charge[jj] != 0)
1311 add_j_to_nblist(coul, jj, bLR);
1315 add_j_to_nblist(coul_adress, jj, bLR);
1319 else if (!bDoCoul_i)
1321 if (bHaveVdW[type[jj]])
1325 add_j_to_nblist(vdw, jj, bLR);
1329 add_j_to_nblist(vdw_adress, jj, bLR);
1335 if (bHaveVdW[type[jj]])
1337 if (charge[jj] != 0)
1341 add_j_to_nblist(vdwc, jj, bLR);
1345 add_j_to_nblist(vdwc_adress, jj, bLR);
1352 add_j_to_nblist(vdw, jj, bLR);
1356 add_j_to_nblist(vdw_adress, jj, bLR);
1361 else if (charge[jj] != 0)
1365 add_j_to_nblist(coul, jj, bLR);
1369 add_j_to_nblist(coul_adress, jj, bLR);
1378 close_i_nblist(vdw);
1379 close_i_nblist(coul);
1380 close_i_nblist(vdwc);
1381 close_i_nblist(vdw_adress);
1382 close_i_nblist(coul_adress);
1383 close_i_nblist(vdwc_adress);
1389 put_in_list_qmmm(gmx_bool gmx_unused bHaveVdW[],
1391 t_mdatoms gmx_unused * md,
1401 gmx_bool gmx_unused bDoVdW,
1402 gmx_bool gmx_unused bDoCoul,
1403 int gmx_unused solvent_opt)
1406 int i, j, jcg, igid, gid;
1407 atom_id jj, jj0, jj1, i_atom;
1411 /* Get atom range */
1413 nicg = index[icg+1]-i0;
1415 /* Get the i charge group info */
1416 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1418 coul = &fr->QMMMlist;
1420 /* Loop over atoms in the ith charge group */
1421 for (i = 0; i < nicg; i++)
1424 gid = GID(igid, jgid, ngid);
1425 /* Create new i_atom for each energy group */
1426 new_i_nblist(coul, i_atom, shift, gid);
1428 /* Loop over the j charge groups */
1429 for (j = 0; j < nj; j++)
1433 /* Charge groups cannot have QM and MM atoms simultaneously */
1438 /* Finally loop over the atoms in the j-charge group */
1439 for (jj = jj0; jj < jj1; jj++)
1441 bNotEx = NOTEXCL(bExcl, i, jj);
1444 add_j_to_nblist(coul, jj, bLR);
1449 close_i_nblist(coul);
1454 put_in_list_cg(gmx_bool gmx_unused bHaveVdW[],
1456 t_mdatoms gmx_unused * md,
1466 gmx_bool gmx_unused bDoVdW,
1467 gmx_bool gmx_unused bDoCoul,
1468 int gmx_unused solvent_opt)
1471 int igid, gid, nbl_ind;
1475 cginfo = fr->cginfo[icg];
1477 igid = GET_CGINFO_GID(cginfo);
1478 gid = GID(igid, jgid, ngid);
1480 /* Unpack pointers to neighbourlist structs */
1481 if (fr->nnblists == 1)
1487 nbl_ind = fr->gid2nblists[gid];
1491 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1495 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1498 /* Make a new neighbor list for charge group icg.
1499 * Currently simply one neighbor list is made with LJ and Coulomb.
1500 * If required, zero interactions could be removed here
1501 * or in the force loop.
1503 new_i_nblist(vdwc, index[icg], shift, gid);
1504 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1506 for (j = 0; (j < nj); j++)
1509 /* Skip the icg-icg pairs if all self interactions are excluded */
1510 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1512 /* Here we add the j charge group jcg to the list,
1513 * exclusions are also added to the list.
1515 add_j_to_nblist_cg(vdwc, index[jcg], index[jcg+1], bExcl, icg == jcg, bLR);
1519 close_i_nblist(vdwc);
1522 static void setexcl(atom_id start, atom_id end, t_blocka *excl, gmx_bool b,
1529 for (i = start; i < end; i++)
1531 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1533 SETEXCL(bexcl, i-start, excl->a[k]);
1539 for (i = start; i < end; i++)
1541 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1543 RMEXCL(bexcl, i-start, excl->a[k]);
1549 int calc_naaj(int icg, int cgtot)
1553 if ((cgtot % 2) == 1)
1555 /* Odd number of charge groups, easy */
1556 naaj = 1 + (cgtot/2);
1558 else if ((cgtot % 4) == 0)
1560 /* Multiple of four is hard */
1597 fprintf(log, "naaj=%d\n", naaj);
1603 /************************************************
1605 * S I M P L E C O R E S T U F F
1607 ************************************************/
1609 static real calc_image_tric(rvec xi, rvec xj, matrix box,
1610 rvec b_inv, int *shift)
1612 /* This code assumes that the cut-off is smaller than
1613 * a half times the smallest diagonal element of the box.
1615 const real h25 = 2.5;
1620 /* Compute diff vector */
1621 dz = xj[ZZ] - xi[ZZ];
1622 dy = xj[YY] - xi[YY];
1623 dx = xj[XX] - xi[XX];
1625 /* Perform NINT operation, using trunc operation, therefore
1626 * we first add 2.5 then subtract 2 again
1628 tz = dz*b_inv[ZZ] + h25;
1630 dz -= tz*box[ZZ][ZZ];
1631 dy -= tz*box[ZZ][YY];
1632 dx -= tz*box[ZZ][XX];
1634 ty = dy*b_inv[YY] + h25;
1636 dy -= ty*box[YY][YY];
1637 dx -= ty*box[YY][XX];
1639 tx = dx*b_inv[XX]+h25;
1641 dx -= tx*box[XX][XX];
1643 /* Distance squared */
1644 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1646 *shift = XYZ2IS(tx, ty, tz);
1651 static real calc_image_rect(rvec xi, rvec xj, rvec box_size,
1652 rvec b_inv, int *shift)
1654 const real h15 = 1.5;
1660 /* Compute diff vector */
1661 dx = xj[XX] - xi[XX];
1662 dy = xj[YY] - xi[YY];
1663 dz = xj[ZZ] - xi[ZZ];
1665 /* Perform NINT operation, using trunc operation, therefore
1666 * we first add 1.5 then subtract 1 again
1668 tx = dx*b_inv[XX] + h15;
1669 ty = dy*b_inv[YY] + h15;
1670 tz = dz*b_inv[ZZ] + h15;
1675 /* Correct diff vector for translation */
1676 ddx = tx*box_size[XX] - dx;
1677 ddy = ty*box_size[YY] - dy;
1678 ddz = tz*box_size[ZZ] - dz;
1680 /* Distance squared */
1681 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1683 *shift = XYZ2IS(tx, ty, tz);
1688 static void add_simple(t_ns_buf *nsbuf, int nrj, atom_id cg_j,
1689 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1690 int icg, int jgid, t_block *cgs, t_excl bexcl[],
1691 int shift, t_forcerec *fr, put_in_list_t *put_in_list)
1693 if (nsbuf->nj + nrj > MAX_CG)
1695 put_in_list(bHaveVdW, ngid, md, icg, jgid, nsbuf->ncg, nsbuf->jcg,
1696 cgs->index, bexcl, shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1697 /* Reset buffer contents */
1698 nsbuf->ncg = nsbuf->nj = 0;
1700 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1704 static void ns_inner_tric(rvec x[], int icg, int *i_egp_flags,
1705 int njcg, atom_id jcg[],
1706 matrix box, rvec b_inv, real rcut2,
1707 t_block *cgs, t_ns_buf **ns_buf,
1708 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1709 t_excl bexcl[], t_forcerec *fr,
1710 put_in_list_t *put_in_list)
1714 int *cginfo = fr->cginfo;
1715 atom_id cg_j, *cgindex;
1718 cgindex = cgs->index;
1720 for (j = 0; (j < njcg); j++)
1723 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1724 if (calc_image_tric(x[icg], x[cg_j], box, b_inv, &shift) < rcut2)
1726 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1727 if (!(i_egp_flags[jgid] & EGP_EXCL))
1729 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1730 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1737 static void ns_inner_rect(rvec x[], int icg, int *i_egp_flags,
1738 int njcg, atom_id jcg[],
1739 gmx_bool bBox, rvec box_size, rvec b_inv, real rcut2,
1740 t_block *cgs, t_ns_buf **ns_buf,
1741 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1742 t_excl bexcl[], t_forcerec *fr,
1743 put_in_list_t *put_in_list)
1747 int *cginfo = fr->cginfo;
1748 atom_id cg_j, *cgindex;
1751 cgindex = cgs->index;
1755 for (j = 0; (j < njcg); j++)
1758 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1759 if (calc_image_rect(x[icg], x[cg_j], box_size, b_inv, &shift) < rcut2)
1761 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1762 if (!(i_egp_flags[jgid] & EGP_EXCL))
1764 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1765 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1773 for (j = 0; (j < njcg); j++)
1776 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1777 if ((rcut2 == 0) || (distance2(x[icg], x[cg_j]) < rcut2))
1779 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1780 if (!(i_egp_flags[jgid] & EGP_EXCL))
1782 add_simple(&ns_buf[jgid][CENTRAL], nrj, cg_j,
1783 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, CENTRAL, fr,
1791 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1793 static int ns_simple_core(t_forcerec *fr,
1794 gmx_localtop_t *top,
1796 matrix box, rvec box_size,
1797 t_excl bexcl[], atom_id *aaj,
1798 int ngid, t_ns_buf **ns_buf,
1799 put_in_list_t *put_in_list, gmx_bool bHaveVdW[])
1803 int nsearch, icg, jcg, igid, i0, nri, nn;
1806 /* atom_id *i_atoms; */
1807 t_block *cgs = &(top->cgs);
1808 t_blocka *excl = &(top->excls);
1811 gmx_bool bBox, bTriclinic;
1814 rlist2 = sqr(fr->rlist);
1816 bBox = (fr->ePBC != epbcNONE);
1819 for (m = 0; (m < DIM); m++)
1821 b_inv[m] = divide_err(1.0, box_size[m]);
1823 bTriclinic = TRICLINIC(box);
1830 cginfo = fr->cginfo;
1833 for (icg = fr->cg0; (icg < fr->hcg); icg++)
1836 i0 = cgs->index[icg];
1837 nri = cgs->index[icg+1]-i0;
1838 i_atoms = &(cgs->a[i0]);
1839 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1840 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1842 igid = GET_CGINFO_GID(cginfo[icg]);
1843 i_egp_flags = fr->egp_flags + ngid*igid;
1844 setexcl(cgs->index[icg], cgs->index[icg+1], excl, TRUE, bexcl);
1846 naaj = calc_naaj(icg, cgs->nr);
1849 ns_inner_tric(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1850 box, b_inv, rlist2, cgs, ns_buf,
1851 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1855 ns_inner_rect(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1856 bBox, box_size, b_inv, rlist2, cgs, ns_buf,
1857 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1861 for (nn = 0; (nn < ngid); nn++)
1863 for (k = 0; (k < SHIFTS); k++)
1865 nsbuf = &(ns_buf[nn][k]);
1868 put_in_list(bHaveVdW, ngid, md, icg, nn, nsbuf->ncg, nsbuf->jcg,
1869 cgs->index, bexcl, k, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1870 nsbuf->ncg = nsbuf->nj = 0;
1874 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1875 setexcl(cgs->index[icg], cgs->index[icg+1], excl, FALSE, bexcl);
1877 close_neighbor_lists(fr, FALSE);
1882 /************************************************
1884 * N S 5 G R I D S T U F F
1886 ************************************************/
1888 static inline void get_dx(int Nx, real gridx, real rc2, int xgi, real x,
1889 int *dx0, int *dx1, real *dcx2)
1917 for (i = xgi0; i >= 0; i--)
1919 dcx = (i+1)*gridx-x;
1928 for (i = xgi1; i < Nx; i++)
1941 static inline void get_dx_dd(int Nx, real gridx, real rc2, int xgi, real x,
1942 int ncpddc, int shift_min, int shift_max,
1943 int *g0, int *g1, real *dcx2)
1946 int g_min, g_max, shift_home;
1979 g_min = (shift_min == shift_home ? 0 : ncpddc);
1980 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1987 else if (shift_max < 0)
2002 /* Check one grid cell down */
2003 dcx = ((*g0 - 1) + 1)*gridx - x;
2015 /* Check one grid cell up */
2016 dcx = (*g1 + 1)*gridx - x;
2028 #define sqr(x) ((x)*(x))
2029 #define calc_dx2(XI, YI, ZI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
2030 #define calc_cyl_dx2(XI, YI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
2031 /****************************************************
2033 * F A S T N E I G H B O R S E A R C H I N G
2035 * Optimized neighboursearching routine using grid
2036 * at least 1x1x1, see GROMACS manual
2038 ****************************************************/
2041 static void get_cutoff2(t_forcerec *fr, gmx_bool bDoLongRange,
2042 real *rvdw2, real *rcoul2,
2043 real *rs2, real *rm2, real *rl2)
2045 *rs2 = sqr(fr->rlist);
2047 if (bDoLongRange && fr->bTwinRange)
2049 /* With plain cut-off or RF we need to make the list exactly
2050 * up to the cut-off and the cut-off's can be different,
2051 * so we can not simply set them to rlistlong.
2052 * To keep this code compatible with (exotic) old cases,
2053 * we also create lists up to rvdw/rcoulomb for PME and Ewald.
2054 * The interaction check should correspond to:
2055 * !ir_vdw/coulomb_might_be_zero_at_cutoff from inputrec.c.
2057 if (((fr->vdwtype == evdwCUT || fr->vdwtype == evdwPME) &&
2058 fr->vdw_modifier == eintmodNONE) ||
2059 fr->rvdw <= fr->rlist)
2061 *rvdw2 = sqr(fr->rvdw);
2065 *rvdw2 = sqr(fr->rlistlong);
2067 if (((fr->eeltype == eelCUT ||
2068 (EEL_RF(fr->eeltype) && fr->eeltype != eelRF_ZERO) ||
2069 fr->eeltype == eelPME ||
2070 fr->eeltype == eelEWALD) &&
2071 fr->coulomb_modifier == eintmodNONE) ||
2072 fr->rcoulomb <= fr->rlist)
2074 *rcoul2 = sqr(fr->rcoulomb);
2078 *rcoul2 = sqr(fr->rlistlong);
2083 /* Workaround for a gcc -O3 or -ffast-math problem */
2087 *rm2 = min(*rvdw2, *rcoul2);
2088 *rl2 = max(*rvdw2, *rcoul2);
2091 static void init_nsgrid_lists(t_forcerec *fr, int ngid, gmx_ns_t *ns)
2093 real rvdw2, rcoul2, rs2, rm2, rl2;
2096 get_cutoff2(fr, TRUE, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2098 /* Short range buffers */
2099 snew(ns->nl_sr, ngid);
2101 snew(ns->nsr, ngid);
2102 snew(ns->nlr_ljc, ngid);
2103 snew(ns->nlr_one, ngid);
2105 /* Always allocate both list types, since rcoulomb might now change with PME load balancing */
2106 /* Long range VdW and Coul buffers */
2107 snew(ns->nl_lr_ljc, ngid);
2108 /* Long range VdW or Coul only buffers */
2109 snew(ns->nl_lr_one, ngid);
2111 for (j = 0; (j < ngid); j++)
2113 snew(ns->nl_sr[j], MAX_CG);
2114 snew(ns->nl_lr_ljc[j], MAX_CG);
2115 snew(ns->nl_lr_one[j], MAX_CG);
2120 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
2125 static int nsgrid_core(t_commrec *cr, t_forcerec *fr,
2126 matrix box, int ngid,
2127 gmx_localtop_t *top,
2129 t_excl bexcl[], gmx_bool *bExcludeAlleg,
2131 put_in_list_t *put_in_list,
2132 gmx_bool bHaveVdW[],
2133 gmx_bool bDoLongRange, gmx_bool bMakeQMMMnblist)
2136 atom_id **nl_lr_ljc, **nl_lr_one, **nl_sr;
2137 int *nlr_ljc, *nlr_one, *nsr;
2138 gmx_domdec_t *dd = NULL;
2139 t_block *cgs = &(top->cgs);
2140 int *cginfo = fr->cginfo;
2141 /* atom_id *i_atoms,*cgsindex=cgs->index; */
2143 int cell_x, cell_y, cell_z;
2144 int d, tx, ty, tz, dx, dy, dz, cj;
2145 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2146 int zsh_ty, zsh_tx, ysh_tx;
2148 int dx0, dx1, dy0, dy1, dz0, dz1;
2149 int Nx, Ny, Nz, shift = -1, j, nrj, nns, nn = -1;
2150 real gridx, gridy, gridz, grid_x, grid_y, grid_z;
2151 real *dcx2, *dcy2, *dcz2;
2153 int cg0, cg1, icg = -1, cgsnr, i0, igid, nri, naaj, max_jcg;
2154 int jcg0, jcg1, jjcg, cgj0, jgid;
2155 int *grida, *gridnra, *gridind;
2156 gmx_bool rvdw_lt_rcoul, rcoul_lt_rvdw;
2157 rvec xi, *cgcm, grid_offset;
2158 real r2, rs2, rvdw2, rcoul2, rm2, rl2, XI, YI, ZI, dcx, dcy, dcz, tmp1, tmp2;
2160 gmx_bool bDomDec, bTriclinicX, bTriclinicY;
2165 bDomDec = DOMAINDECOMP(cr);
2171 bTriclinicX = ((YY < grid->npbcdim &&
2172 (!bDomDec || dd->nc[YY] == 1) && box[YY][XX] != 0) ||
2173 (ZZ < grid->npbcdim &&
2174 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][XX] != 0));
2175 bTriclinicY = (ZZ < grid->npbcdim &&
2176 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][YY] != 0);
2180 get_cutoff2(fr, bDoLongRange, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2182 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
2183 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
2185 if (bMakeQMMMnblist)
2193 nl_lr_ljc = ns->nl_lr_ljc;
2194 nl_lr_one = ns->nl_lr_one;
2195 nlr_ljc = ns->nlr_ljc;
2196 nlr_one = ns->nlr_one;
2204 gridind = grid->index;
2205 gridnra = grid->nra;
2208 gridx = grid->cell_size[XX];
2209 gridy = grid->cell_size[YY];
2210 gridz = grid->cell_size[ZZ];
2214 copy_rvec(grid->cell_offset, grid_offset);
2215 copy_ivec(grid->ncpddc, ncpddc);
2220 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2221 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
2222 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
2223 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
2224 if (zsh_tx != 0 && ysh_tx != 0)
2226 /* This could happen due to rounding, when both ratios are 0.5 */
2235 /* We only want a list for the test particle */
2244 /* Set the shift range */
2245 for (d = 0; d < DIM; d++)
2249 /* Check if we need periodicity shifts.
2250 * Without PBC or with domain decomposition we don't need them.
2252 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2259 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
2270 /* Loop over charge groups */
2271 for (icg = cg0; (icg < cg1); icg++)
2273 igid = GET_CGINFO_GID(cginfo[icg]);
2274 /* Skip this charge group if all energy groups are excluded! */
2275 if (bExcludeAlleg[igid])
2280 i0 = cgs->index[icg];
2282 if (bMakeQMMMnblist)
2284 /* Skip this charge group if it is not a QM atom while making a
2285 * QM/MM neighbourlist
2287 if (md->bQM[i0] == FALSE)
2289 continue; /* MM particle, go to next particle */
2292 /* Compute the number of charge groups that fall within the control
2295 naaj = calc_naaj(icg, cgsnr);
2302 /* make a normal neighbourlist */
2306 /* Get the j charge-group and dd cell shift ranges */
2307 dd_get_ns_ranges(cr->dd, icg, &jcg0, &jcg1, sh0, sh1);
2312 /* Compute the number of charge groups that fall within the control
2315 naaj = calc_naaj(icg, cgsnr);
2321 /* The i-particle is awlways the test particle,
2322 * so we want all j-particles
2324 max_jcg = cgsnr - 1;
2328 max_jcg = jcg1 - cgsnr;
2333 i_egp_flags = fr->egp_flags + igid*ngid;
2335 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2336 setexcl(i0, cgs->index[icg+1], &top->excls, TRUE, bexcl);
2338 ci2xyz(grid, icg, &cell_x, &cell_y, &cell_z);
2340 /* Changed iicg to icg, DvdS 990115
2341 * (but see consistency check above, DvdS 990330)
2344 fprintf(log, "icg=%5d, naaj=%5d, cell %d %d %d\n",
2345 icg, naaj, cell_x, cell_y, cell_z);
2347 /* Loop over shift vectors in three dimensions */
2348 for (tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
2350 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2351 /* Calculate range of cells in Z direction that have the shift tz */
2352 zgi = cell_z + tz*Nz;
2355 get_dx(Nz, gridz, rl2, zgi, ZI, &dz0, &dz1, dcz2);
2357 get_dx_dd(Nz, gridz, rl2, zgi, ZI-grid_offset[ZZ],
2358 ncpddc[ZZ], sh0[ZZ], sh1[ZZ], &dz0, &dz1, dcz2);
2364 for (ty = -shp[YY]; ty <= shp[YY]; ty++)
2366 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2367 /* Calculate range of cells in Y direction that have the shift ty */
2370 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2374 ygi = cell_y + ty*Ny;
2377 get_dx(Ny, gridy, rl2, ygi, YI, &dy0, &dy1, dcy2);
2379 get_dx_dd(Ny, gridy, rl2, ygi, YI-grid_offset[YY],
2380 ncpddc[YY], sh0[YY], sh1[YY], &dy0, &dy1, dcy2);
2386 for (tx = -shp[XX]; tx <= shp[XX]; tx++)
2388 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2389 /* Calculate range of cells in X direction that have the shift tx */
2392 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2396 xgi = cell_x + tx*Nx;
2399 get_dx(Nx, gridx, rl2, xgi*Nx, XI, &dx0, &dx1, dcx2);
2401 get_dx_dd(Nx, gridx, rl2, xgi, XI-grid_offset[XX],
2402 ncpddc[XX], sh0[XX], sh1[XX], &dx0, &dx1, dcx2);
2408 /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2409 * from the neigbour list as it will not interact */
2410 if (fr->adress_type != eAdressOff)
2412 if (md->wf[cgs->index[icg]] <= GMX_REAL_EPS && egp_explicit(fr, igid))
2417 /* Get shift vector */
2418 shift = XYZ2IS(tx, ty, tz);
2420 range_check(shift, 0, SHIFTS);
2422 for (nn = 0; (nn < ngid); nn++)
2429 fprintf(log, "shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2430 shift, dx0, dx1, dy0, dy1, dz0, dz1);
2431 fprintf(log, "cgcm: %8.3f %8.3f %8.3f\n", cgcm[icg][XX],
2432 cgcm[icg][YY], cgcm[icg][ZZ]);
2433 fprintf(log, "xi: %8.3f %8.3f %8.3f\n", XI, YI, ZI);
2435 for (dx = dx0; (dx <= dx1); dx++)
2437 tmp1 = rl2 - dcx2[dx];
2438 for (dy = dy0; (dy <= dy1); dy++)
2440 tmp2 = tmp1 - dcy2[dy];
2443 for (dz = dz0; (dz <= dz1); dz++)
2445 if (tmp2 > dcz2[dz])
2447 /* Find grid-cell cj in which possible neighbours are */
2448 cj = xyz2ci(Ny, Nz, dx, dy, dz);
2450 /* Check out how many cgs (nrj) there in this cell */
2453 /* Find the offset in the cg list */
2456 /* Check if all j's are out of range so we
2457 * can skip the whole cell.
2458 * Should save some time, especially with DD.
2461 (grida[cgj0] >= max_jcg &&
2462 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2468 for (j = 0; (j < nrj); j++)
2470 jjcg = grida[cgj0+j];
2472 /* check whether this guy is in range! */
2473 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2476 r2 = calc_dx2(XI, YI, ZI, cgcm[jjcg]);
2479 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2480 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2481 /* check energy group exclusions */
2482 if (!(i_egp_flags[jgid] & EGP_EXCL))
2486 if (nsr[jgid] >= MAX_CG)
2488 /* Add to short-range list */
2489 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2490 nsr[jgid], nl_sr[jgid],
2491 cgs->index, /* cgsatoms, */ bexcl,
2492 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2495 nl_sr[jgid][nsr[jgid]++] = jjcg;
2499 if (nlr_ljc[jgid] >= MAX_CG)
2501 /* Add to LJ+coulomb long-range list */
2502 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2503 nlr_ljc[jgid], nl_lr_ljc[jgid], top->cgs.index,
2504 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2507 nl_lr_ljc[jgid][nlr_ljc[jgid]++] = jjcg;
2511 if (nlr_one[jgid] >= MAX_CG)
2513 /* Add to long-range list with only coul, or only LJ */
2514 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2515 nlr_one[jgid], nl_lr_one[jgid], top->cgs.index,
2516 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2519 nl_lr_one[jgid][nlr_one[jgid]++] = jjcg;
2531 /* CHECK whether there is anything left in the buffers */
2532 for (nn = 0; (nn < ngid); nn++)
2536 put_in_list(bHaveVdW, ngid, md, icg, nn, nsr[nn], nl_sr[nn],
2537 cgs->index, /* cgsatoms, */ bexcl,
2538 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2541 if (nlr_ljc[nn] > 0)
2543 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_ljc[nn],
2544 nl_lr_ljc[nn], top->cgs.index,
2545 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2548 if (nlr_one[nn] > 0)
2550 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_one[nn],
2551 nl_lr_one[nn], top->cgs.index,
2552 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2558 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2559 setexcl(cgs->index[icg], cgs->index[icg+1], &top->excls, FALSE, bexcl);
2561 /* No need to perform any left-over force calculations anymore (as we used to do here)
2562 * since we now save the proper long-range lists for later evaluation.
2567 /* Close neighbourlists */
2568 close_neighbor_lists(fr, bMakeQMMMnblist);
2573 void ns_realloc_natoms(gmx_ns_t *ns, int natoms)
2577 if (natoms > ns->nra_alloc)
2579 ns->nra_alloc = over_alloc_dd(natoms);
2580 srenew(ns->bexcl, ns->nra_alloc);
2581 for (i = 0; i < ns->nra_alloc; i++)
2588 void init_ns(FILE *fplog, const t_commrec *cr,
2589 gmx_ns_t *ns, t_forcerec *fr,
2590 const gmx_mtop_t *mtop)
2592 int mt, icg, nr_in_cg, maxcg, i, j, jcg, ngid, ncg;
2596 /* Compute largest charge groups size (# atoms) */
2598 for (mt = 0; mt < mtop->nmoltype; mt++)
2600 cgs = &mtop->moltype[mt].cgs;
2601 for (icg = 0; (icg < cgs->nr); icg++)
2603 nr_in_cg = max(nr_in_cg, (int)(cgs->index[icg+1]-cgs->index[icg]));
2607 /* Verify whether largest charge group is <= max cg.
2608 * This is determined by the type of the local exclusion type
2609 * Exclusions are stored in bits. (If the type is not large
2610 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2612 maxcg = sizeof(t_excl)*8;
2613 if (nr_in_cg > maxcg)
2615 gmx_fatal(FARGS, "Max #atoms in a charge group: %d > %d\n",
2619 ngid = mtop->groups.grps[egcENER].nr;
2620 snew(ns->bExcludeAlleg, ngid);
2621 for (i = 0; i < ngid; i++)
2623 ns->bExcludeAlleg[i] = TRUE;
2624 for (j = 0; j < ngid; j++)
2626 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2628 ns->bExcludeAlleg[i] = FALSE;
2636 ns->grid = init_grid(fplog, fr);
2637 init_nsgrid_lists(fr, ngid, ns);
2642 snew(ns->ns_buf, ngid);
2643 for (i = 0; (i < ngid); i++)
2645 snew(ns->ns_buf[i], SHIFTS);
2647 ncg = ncg_mtop(mtop);
2648 snew(ns->simple_aaj, 2*ncg);
2649 for (jcg = 0; (jcg < ncg); jcg++)
2651 ns->simple_aaj[jcg] = jcg;
2652 ns->simple_aaj[jcg+ncg] = jcg;
2656 /* Create array that determines whether or not atoms have VdW */
2657 snew(ns->bHaveVdW, fr->ntype);
2658 for (i = 0; (i < fr->ntype); i++)
2660 for (j = 0; (j < fr->ntype); j++)
2662 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2664 ((BHAMA(fr->nbfp, fr->ntype, i, j) != 0) ||
2665 (BHAMB(fr->nbfp, fr->ntype, i, j) != 0) ||
2666 (BHAMC(fr->nbfp, fr->ntype, i, j) != 0)) :
2667 ((C6(fr->nbfp, fr->ntype, i, j) != 0) ||
2668 (C12(fr->nbfp, fr->ntype, i, j) != 0))));
2673 pr_bvec(debug, 0, "bHaveVdW", ns->bHaveVdW, fr->ntype, TRUE);
2678 if (!DOMAINDECOMP(cr))
2680 ns_realloc_natoms(ns, mtop->natoms);
2683 ns->nblist_initialized = FALSE;
2685 /* nbr list debug dump */
2687 char *ptr = getenv("GMX_DUMP_NL");
2690 ns->dump_nl = strtol(ptr, NULL, 10);
2693 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2704 int search_neighbours(FILE *log, t_forcerec *fr,
2706 gmx_localtop_t *top,
2707 gmx_groups_t *groups,
2709 t_nrnb *nrnb, t_mdatoms *md,
2711 gmx_bool bDoLongRangeNS)
2713 t_block *cgs = &(top->cgs);
2714 rvec box_size, grid_x0, grid_x1;
2716 real min_size, grid_dens;
2720 gmx_bool *i_egp_flags;
2721 int cg_start, cg_end, start, end;
2724 gmx_domdec_zones_t *dd_zones;
2725 put_in_list_t *put_in_list;
2729 /* Set some local variables */
2731 ngid = groups->grps[egcENER].nr;
2733 for (m = 0; (m < DIM); m++)
2735 box_size[m] = box[m][m];
2738 if (fr->ePBC != epbcNONE)
2740 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC, box))
2742 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.");
2746 min_size = min(box_size[XX], min(box_size[YY], box_size[ZZ]));
2747 if (2*fr->rlistlong >= min_size)
2749 gmx_fatal(FARGS, "One of the box diagonal elements has become smaller than twice the cut-off length.");
2754 if (DOMAINDECOMP(cr))
2756 ns_realloc_natoms(ns, cgs->index[cgs->nr]);
2760 /* Reset the neighbourlists */
2761 reset_neighbor_lists(fr, TRUE, TRUE);
2763 if (bGrid && bFillGrid)
2767 if (DOMAINDECOMP(cr))
2769 dd_zones = domdec_zones(cr->dd);
2775 get_nsgrid_boundaries(grid->nboundeddim, box, NULL, NULL, NULL, NULL,
2776 cgs->nr, fr->cg_cm, grid_x0, grid_x1, &grid_dens);
2778 grid_first(log, grid, NULL, NULL, box, grid_x0, grid_x1,
2779 fr->rlistlong, grid_dens);
2786 if (DOMAINDECOMP(cr))
2789 fill_grid(dd_zones, grid, end, -1, end, fr->cg_cm);
2791 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2795 fill_grid(NULL, grid, cgs->nr, fr->cg0, fr->hcg, fr->cg_cm);
2796 grid->icg0 = fr->cg0;
2797 grid->icg1 = fr->hcg;
2801 calc_elemnr(grid, start, end, cgs->nr);
2803 grid_last(grid, start, end, cgs->nr);
2808 print_grid(debug, grid);
2813 /* Set the grid cell index for the test particle only.
2814 * The cell to cg index is not corrected, but that does not matter.
2816 fill_grid(NULL, ns->grid, fr->hcg, fr->hcg-1, fr->hcg, fr->cg_cm);
2820 if (fr->adress_type == eAdressOff)
2822 if (!fr->ns.bCGlist)
2824 put_in_list = put_in_list_at;
2828 put_in_list = put_in_list_cg;
2833 put_in_list = put_in_list_adress;
2840 nsearch = nsgrid_core(cr, fr, box, ngid, top,
2841 grid, ns->bexcl, ns->bExcludeAlleg,
2842 md, put_in_list, ns->bHaveVdW,
2843 bDoLongRangeNS, FALSE);
2845 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2846 * the classical calculation. The charge-charge interaction
2847 * between QM and MM atoms is handled in the QMMM core calculation
2848 * (see QMMM.c). The VDW however, we'd like to compute classically
2849 * and the QM MM atom pairs have just been put in the
2850 * corresponding neighbourlists. in case of QMMM we still need to
2851 * fill a special QMMM neighbourlist that contains all neighbours
2852 * of the QM atoms. If bQMMM is true, this list will now be made:
2854 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
2856 nsearch += nsgrid_core(cr, fr, box, ngid, top,
2857 grid, ns->bexcl, ns->bExcludeAlleg,
2858 md, put_in_list_qmmm, ns->bHaveVdW,
2859 bDoLongRangeNS, TRUE);
2864 nsearch = ns_simple_core(fr, top, md, box, box_size,
2865 ns->bexcl, ns->simple_aaj,
2866 ngid, ns->ns_buf, put_in_list, ns->bHaveVdW);
2874 inc_nrnb(nrnb, eNR_NS, nsearch);
2875 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2880 int natoms_beyond_ns_buffer(t_inputrec *ir, t_forcerec *fr, t_block *cgs,
2881 matrix scale_tot, rvec *x)
2883 int cg0, cg1, cg, a0, a1, a, i, j;
2884 real rint, hbuf2, scale;
2886 gmx_bool bIsotropic;
2891 rint = max(ir->rcoulomb, ir->rvdw);
2892 if (ir->rlist < rint)
2894 gmx_fatal(FARGS, "The neighbor search buffer has negative size: %f nm",
2902 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2904 hbuf2 = sqr(0.5*(ir->rlist - rint));
2905 for (cg = cg0; cg < cg1; cg++)
2907 a0 = cgs->index[cg];
2908 a1 = cgs->index[cg+1];
2909 for (a = a0; a < a1; a++)
2911 if (distance2(cg_cm[cg], x[a]) > hbuf2)
2921 scale = scale_tot[0][0];
2922 for (i = 1; i < DIM; i++)
2924 /* With anisotropic scaling, the original spherical ns volumes become
2925 * ellipsoids. To avoid costly transformations we use the minimum
2926 * eigenvalue of the scaling matrix for determining the buffer size.
2927 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2929 scale = min(scale, scale_tot[i][i]);
2930 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2934 for (j = 0; j < i; j++)
2936 if (scale_tot[i][j] != 0)
2942 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2945 for (cg = cg0; cg < cg1; cg++)
2947 svmul(scale, cg_cm[cg], cgsc);
2948 a0 = cgs->index[cg];
2949 a1 = cgs->index[cg+1];
2950 for (a = a0; a < a1; a++)
2952 if (distance2(cgsc, x[a]) > hbuf2)
2961 /* Anistropic scaling */
2962 for (cg = cg0; cg < cg1; cg++)
2964 /* Since scale_tot contains the transpose of the scaling matrix,
2965 * we need to multiply with the transpose.
2967 tmvmul_ur0(scale_tot, cg_cm[cg], cgsc);
2968 a0 = cgs->index[cg];
2969 a1 = cgs->index[cg+1];
2970 for (a = a0; a < a1; a++)
2972 if (distance2(cgsc, x[a]) > hbuf2)