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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
52 #include "nonbonded.h"
56 #include "gmx_fatal.h"
59 #include "mtop_util.h"
66 * E X C L U S I O N H A N D L I N G
70 static void SETEXCL_(t_excl e[], atom_id i, atom_id j)
74 static void RMEXCL_(t_excl e[], atom_id i, atom_id j)
76 e[j] = e[j] & ~(1<<i);
78 static gmx_bool ISEXCL_(t_excl e[], atom_id i, atom_id j)
80 return (gmx_bool)(e[j] & (1<<i));
82 static gmx_bool NOTEXCL_(t_excl e[], atom_id i, atom_id j)
84 return !(ISEXCL(e, i, j));
87 #define SETEXCL(e, i, j) (e)[((atom_id) (j))] |= (1<<((atom_id) (i)))
88 #define RMEXCL(e, i, j) (e)[((atom_id) (j))] &= (~(1<<((atom_id) (i))))
89 #define ISEXCL(e, i, j) (gmx_bool) ((e)[((atom_id) (j))] & (1<<((atom_id) (i))))
90 #define NOTEXCL(e, i, j) !(ISEXCL(e, i, j))
94 round_up_to_simd_width(int length, int simd_width)
96 int offset, newlength;
98 offset = (simd_width > 0) ? length % simd_width : 0;
100 return (offset == 0) ? length : length-offset+simd_width;
102 /************************************************
104 * U T I L I T I E S F O R N S
106 ************************************************/
108 static void reallocate_nblist(t_nblist *nl)
112 fprintf(debug, "reallocating neigborlist (ielec=%d, ivdw=%d, igeometry=%d, type=%d), maxnri=%d\n",
113 nl->ielec, nl->ivdw, nl->igeometry, nl->type, nl->maxnri);
115 srenew(nl->iinr, nl->maxnri);
116 if (nl->igeometry == GMX_NBLIST_GEOMETRY_CG_CG)
118 srenew(nl->iinr_end, nl->maxnri);
120 srenew(nl->gid, nl->maxnri);
121 srenew(nl->shift, nl->maxnri);
122 srenew(nl->jindex, nl->maxnri+1);
126 static void init_nblist(FILE *log, t_nblist *nl_sr, t_nblist *nl_lr,
127 int maxsr, int maxlr,
128 int ivdw, int ivdwmod,
129 int ielec, int ielecmod,
130 int igeometry, int type)
136 for (i = 0; (i < 2); i++)
138 nl = (i == 0) ? nl_sr : nl_lr;
139 homenr = (i == 0) ? maxsr : maxlr;
147 /* Set coul/vdw in neighborlist, and for the normal loops we determine
148 * an index of which one to call.
151 nl->ivdwmod = ivdwmod;
153 nl->ielecmod = ielecmod;
155 nl->igeometry = igeometry;
157 if (nl->type == GMX_NBLIST_INTERACTION_FREE_ENERGY)
159 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
162 /* This will also set the simd_padding_width field */
163 gmx_nonbonded_set_kernel_pointers( (i == 0) ? log : NULL, nl);
165 /* maxnri is influenced by the number of shifts (maximum is 8)
166 * and the number of energy groups.
167 * If it is not enough, nl memory will be reallocated during the run.
168 * 4 seems to be a reasonable factor, which only causes reallocation
169 * during runs with tiny and many energygroups.
171 nl->maxnri = homenr*4;
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;
295 fprintf(log, "Note: The available nonbonded kernels do not support water optimization - disabling.\n");
299 if (fr->efep != efepNO)
301 if ((fr->bEwald) && (fr->sc_alphacoul > 0)) /* need to handle long range differently if using softcore */
303 ielecf = GMX_NBKERNEL_ELEC_EWALD;
304 ielecmodf = eintmodNONE;
309 ielecmodf = ielecmod;
312 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_FREE], &nbl->nlist_lr[eNL_VDWQQ_FREE],
313 maxsr, maxlr, ivdw, ivdwmod, ielecf, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
314 init_nblist(log, &nbl->nlist_sr[eNL_VDW_FREE], &nbl->nlist_lr[eNL_VDW_FREE],
315 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
316 init_nblist(log, &nbl->nlist_sr[eNL_QQ_FREE], &nbl->nlist_lr[eNL_QQ_FREE],
317 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielecf, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
321 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
323 init_nblist(log, &fr->QMMMlist, NULL,
324 maxsr, maxlr, 0, 0, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_STANDARD);
332 fr->ns.nblist_initialized = TRUE;
335 static void reset_nblist(t_nblist *nl)
346 static void reset_neighbor_lists(t_forcerec *fr, gmx_bool bResetSR, gmx_bool bResetLR)
352 /* only reset the short-range nblist */
353 reset_nblist(&(fr->QMMMlist));
356 for (n = 0; n < fr->nnblists; n++)
358 for (i = 0; i < eNL_NR; i++)
362 reset_nblist( &(fr->nblists[n].nlist_sr[i]) );
366 reset_nblist( &(fr->nblists[n].nlist_lr[i]) );
375 static inline void new_i_nblist(t_nblist *nlist,
376 gmx_bool bLR, atom_id i_atom, int shift, int gid)
378 int i, k, nri, nshift;
382 /* Check whether we have to increase the i counter */
384 (nlist->iinr[nri] != i_atom) ||
385 (nlist->shift[nri] != shift) ||
386 (nlist->gid[nri] != gid))
388 /* This is something else. Now see if any entries have
389 * been added in the list of the previous atom.
392 ((nlist->jindex[nri+1] > nlist->jindex[nri]) &&
393 (nlist->gid[nri] != -1)))
395 /* If so increase the counter */
398 if (nlist->nri >= nlist->maxnri)
400 nlist->maxnri += over_alloc_large(nlist->nri);
401 reallocate_nblist(nlist);
404 /* Set the number of neighbours and the atom number */
405 nlist->jindex[nri+1] = nlist->jindex[nri];
406 nlist->iinr[nri] = i_atom;
407 nlist->gid[nri] = gid;
408 nlist->shift[nri] = shift;
412 /* Adding to previous list. First remove possible previous padding */
413 if(nlist->simd_padding_width>1)
415 while(nlist->nrj>0 && nlist->jjnr[nlist->nrj-1]<0)
423 static inline void close_i_nblist(t_nblist *nlist)
425 int nri = nlist->nri;
430 /* Add elements up to padding. Since we allocate memory in units
431 * of the simd_padding width, we do not have to check for possible
432 * list reallocation here.
434 while ((nlist->nrj % nlist->simd_padding_width) != 0)
436 /* Use -4 here, so we can write forces for 4 atoms before real data */
437 nlist->jjnr[nlist->nrj++] = -4;
439 nlist->jindex[nri+1] = nlist->nrj;
441 len = nlist->nrj - nlist->jindex[nri];
443 /* nlist length for water i molecules is treated statically
446 if (len > nlist->maxlen)
453 static inline void close_nblist(t_nblist *nlist)
455 /* Only close this nblist when it has been initialized.
456 * Avoid the creation of i-lists with no j-particles.
460 /* Some assembly kernels do not support empty lists,
461 * make sure here that we don't generate any empty lists.
462 * With the current ns code this branch is taken in two cases:
463 * No i-particles at all: nri=-1 here
464 * There are i-particles, but no j-particles; nri=0 here
470 /* Close list number nri by incrementing the count */
475 static inline void close_neighbor_lists(t_forcerec *fr, gmx_bool bMakeQMMMnblist)
481 close_nblist(&(fr->QMMMlist));
484 for (n = 0; n < fr->nnblists; n++)
486 for (i = 0; (i < eNL_NR); i++)
488 close_nblist(&(fr->nblists[n].nlist_sr[i]));
489 close_nblist(&(fr->nblists[n].nlist_lr[i]));
495 static inline void add_j_to_nblist(t_nblist *nlist, atom_id j_atom, gmx_bool bLR)
497 int nrj = nlist->nrj;
499 if (nlist->nrj >= nlist->maxnrj)
501 nlist->maxnrj = round_up_to_simd_width(over_alloc_small(nlist->nrj + 1), nlist->simd_padding_width);
505 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
506 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
509 srenew(nlist->jjnr, nlist->maxnrj);
512 nlist->jjnr[nrj] = j_atom;
516 static inline void add_j_to_nblist_cg(t_nblist *nlist,
517 atom_id j_start, int j_end,
518 t_excl *bexcl, gmx_bool i_is_j,
521 int nrj = nlist->nrj;
524 if (nlist->nrj >= nlist->maxnrj)
526 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
529 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
530 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
533 srenew(nlist->jjnr, nlist->maxnrj);
534 srenew(nlist->jjnr_end, nlist->maxnrj);
535 srenew(nlist->excl, nlist->maxnrj*MAX_CGCGSIZE);
538 nlist->jjnr[nrj] = j_start;
539 nlist->jjnr_end[nrj] = j_end;
541 if (j_end - j_start > MAX_CGCGSIZE)
543 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);
546 /* Set the exclusions */
547 for (j = j_start; j < j_end; j++)
549 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
553 /* Avoid double counting of intra-cg interactions */
554 for (j = 1; j < j_end-j_start; j++)
556 nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
564 put_in_list_t (gmx_bool bHaveVdW[],
581 put_in_list_at(gmx_bool bHaveVdW[],
597 /* The a[] index has been removed,
598 * to put it back in i_atom should be a[i0] and jj should be a[jj].
603 t_nblist * vdwc_free = NULL;
604 t_nblist * vdw_free = NULL;
605 t_nblist * coul_free = NULL;
606 t_nblist * vdwc_ww = NULL;
607 t_nblist * coul_ww = NULL;
609 int i, j, jcg, igid, gid, nbl_ind, ind_ij;
610 atom_id jj, jj0, jj1, i_atom;
615 real *charge, *chargeB;
616 real qi, qiB, qq, rlj;
617 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
618 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
622 /* Copy some pointers */
624 charge = md->chargeA;
625 chargeB = md->chargeB;
628 bPert = md->bPerturbed;
632 nicg = index[icg+1]-i0;
634 /* Get the i charge group info */
635 igid = GET_CGINFO_GID(cginfo[icg]);
637 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
642 /* Check if any of the particles involved are perturbed.
643 * If not we can do the cheaper normal put_in_list
644 * and use more solvent optimization.
646 for (i = 0; i < nicg; i++)
648 bFreeEnergy |= bPert[i0+i];
650 /* Loop over the j charge groups */
651 for (j = 0; (j < nj && !bFreeEnergy); j++)
656 /* Finally loop over the atoms in the j-charge group */
657 for (jj = jj0; jj < jj1; jj++)
659 bFreeEnergy |= bPert[jj];
664 /* Unpack pointers to neighbourlist structs */
665 if (fr->nnblists == 1)
671 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
675 nlist = fr->nblists[nbl_ind].nlist_lr;
679 nlist = fr->nblists[nbl_ind].nlist_sr;
682 if (iwater != esolNO)
684 vdwc = &nlist[eNL_VDWQQ_WATER];
685 vdw = &nlist[eNL_VDW];
686 coul = &nlist[eNL_QQ_WATER];
687 #ifndef DISABLE_WATERWATER_NLIST
688 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
689 coul_ww = &nlist[eNL_QQ_WATERWATER];
694 vdwc = &nlist[eNL_VDWQQ];
695 vdw = &nlist[eNL_VDW];
696 coul = &nlist[eNL_QQ];
701 if (iwater != esolNO)
703 /* Loop over the atoms in the i charge group */
705 gid = GID(igid, jgid, ngid);
706 /* Create new i_atom for each energy group */
707 if (bDoCoul && bDoVdW)
709 new_i_nblist(vdwc, bLR, i_atom, shift, gid);
710 #ifndef DISABLE_WATERWATER_NLIST
711 new_i_nblist(vdwc_ww, bLR, i_atom, shift, gid);
716 new_i_nblist(vdw, bLR, i_atom, shift, gid);
720 new_i_nblist(coul, bLR, i_atom, shift, gid);
721 #ifndef DISABLE_WATERWATER_NLIST
722 new_i_nblist(coul_ww, bLR, i_atom, shift, gid);
725 /* Loop over the j charge groups */
726 for (j = 0; (j < nj); j++)
736 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
738 if (iwater == esolSPC && jwater == esolSPC)
740 /* Interaction between two SPC molecules */
743 /* VdW only - only first atoms in each water interact */
744 add_j_to_nblist(vdw, jj0, bLR);
748 #ifdef DISABLE_WATERWATER_NLIST
749 /* Add entries for the three atoms - only do VdW if we need to */
752 add_j_to_nblist(coul, jj0, bLR);
756 add_j_to_nblist(vdwc, jj0, bLR);
758 add_j_to_nblist(coul, jj0+1, bLR);
759 add_j_to_nblist(coul, jj0+2, bLR);
761 /* One entry for the entire water-water interaction */
764 add_j_to_nblist(coul_ww, jj0, bLR);
768 add_j_to_nblist(vdwc_ww, jj0, bLR);
773 else if (iwater == esolTIP4P && jwater == esolTIP4P)
775 /* Interaction between two TIP4p molecules */
778 /* VdW only - only first atoms in each water interact */
779 add_j_to_nblist(vdw, jj0, bLR);
783 #ifdef DISABLE_WATERWATER_NLIST
784 /* Add entries for the four atoms - only do VdW if we need to */
787 add_j_to_nblist(vdw, jj0, bLR);
789 add_j_to_nblist(coul, jj0+1, bLR);
790 add_j_to_nblist(coul, jj0+2, bLR);
791 add_j_to_nblist(coul, jj0+3, bLR);
793 /* One entry for the entire water-water interaction */
796 add_j_to_nblist(coul_ww, jj0, bLR);
800 add_j_to_nblist(vdwc_ww, jj0, bLR);
807 /* j charge group is not water, but i is.
808 * Add entries to the water-other_atom lists; the geometry of the water
809 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
810 * so we don't care if it is SPC or TIP4P...
817 for (jj = jj0; (jj < jj1); jj++)
821 add_j_to_nblist(coul, jj, bLR);
827 for (jj = jj0; (jj < jj1); jj++)
829 if (bHaveVdW[type[jj]])
831 add_j_to_nblist(vdw, jj, bLR);
837 /* _charge_ _groups_ interact with both coulomb and LJ */
838 /* Check which atoms we should add to the lists! */
839 for (jj = jj0; (jj < jj1); jj++)
841 if (bHaveVdW[type[jj]])
845 add_j_to_nblist(vdwc, jj, bLR);
849 add_j_to_nblist(vdw, jj, bLR);
852 else if (charge[jj] != 0)
854 add_j_to_nblist(coul, jj, bLR);
861 close_i_nblist(coul);
862 close_i_nblist(vdwc);
863 #ifndef DISABLE_WATERWATER_NLIST
864 close_i_nblist(coul_ww);
865 close_i_nblist(vdwc_ww);
870 /* no solvent as i charge group */
871 /* Loop over the atoms in the i charge group */
872 for (i = 0; i < nicg; i++)
875 gid = GID(igid, jgid, ngid);
878 /* Create new i_atom for each energy group */
879 if (bDoVdW && bDoCoul)
881 new_i_nblist(vdwc, bLR, i_atom, shift, gid);
885 new_i_nblist(vdw, bLR, i_atom, shift, gid);
889 new_i_nblist(coul, bLR, i_atom, shift, gid);
891 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
892 bDoCoul_i = (bDoCoul && qi != 0);
894 if (bDoVdW_i || bDoCoul_i)
896 /* Loop over the j charge groups */
897 for (j = 0; (j < nj); j++)
901 /* Check for large charge groups */
912 /* Finally loop over the atoms in the j-charge group */
913 for (jj = jj0; jj < jj1; jj++)
915 bNotEx = NOTEXCL(bExcl, i, jj);
923 add_j_to_nblist(coul, jj, bLR);
928 if (bHaveVdW[type[jj]])
930 add_j_to_nblist(vdw, jj, bLR);
935 if (bHaveVdW[type[jj]])
939 add_j_to_nblist(vdwc, jj, bLR);
943 add_j_to_nblist(vdw, jj, bLR);
946 else if (charge[jj] != 0)
948 add_j_to_nblist(coul, jj, bLR);
956 close_i_nblist(coul);
957 close_i_nblist(vdwc);
963 /* we are doing free energy */
964 vdwc_free = &nlist[eNL_VDWQQ_FREE];
965 vdw_free = &nlist[eNL_VDW_FREE];
966 coul_free = &nlist[eNL_QQ_FREE];
967 /* Loop over the atoms in the i charge group */
968 for (i = 0; i < nicg; i++)
971 gid = GID(igid, jgid, ngid);
973 qiB = chargeB[i_atom];
975 /* Create new i_atom for each energy group */
976 if (bDoVdW && bDoCoul)
978 new_i_nblist(vdwc, bLR, i_atom, shift, gid);
982 new_i_nblist(vdw, bLR, i_atom, shift, gid);
986 new_i_nblist(coul, bLR, i_atom, shift, gid);
989 new_i_nblist(vdw_free, bLR, i_atom, shift, gid);
990 new_i_nblist(coul_free, bLR, i_atom, shift, gid);
991 new_i_nblist(vdwc_free, bLR, i_atom, shift, gid);
993 bDoVdW_i = (bDoVdW &&
994 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
995 bDoCoul_i = (bDoCoul && (qi != 0 || qiB != 0));
996 /* For TIP4P the first atom does not have a charge,
997 * but the last three do. So we should still put an atom
998 * without LJ but with charge in the water-atom neighborlist
999 * for a TIP4p i charge group.
1000 * For SPC type water the first atom has LJ and charge,
1001 * so there is no such problem.
1003 if (iwater == esolNO)
1005 bDoCoul_i_sol = bDoCoul_i;
1009 bDoCoul_i_sol = bDoCoul;
1012 if (bDoVdW_i || bDoCoul_i_sol)
1014 /* Loop over the j charge groups */
1015 for (j = 0; (j < nj); j++)
1019 /* Check for large charge groups */
1030 /* Finally loop over the atoms in the j-charge group */
1031 bFree = bPert[i_atom];
1032 for (jj = jj0; (jj < jj1); jj++)
1034 bFreeJ = bFree || bPert[jj];
1035 /* Complicated if, because the water H's should also
1036 * see perturbed j-particles
1038 if (iwater == esolNO || i == 0 || bFreeJ)
1040 bNotEx = NOTEXCL(bExcl, i, jj);
1048 if (charge[jj] != 0 || chargeB[jj] != 0)
1050 add_j_to_nblist(coul_free, jj, bLR);
1053 else if (!bDoCoul_i)
1055 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1057 add_j_to_nblist(vdw_free, jj, bLR);
1062 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1064 if (charge[jj] != 0 || chargeB[jj] != 0)
1066 add_j_to_nblist(vdwc_free, jj, bLR);
1070 add_j_to_nblist(vdw_free, jj, bLR);
1073 else if (charge[jj] != 0 || chargeB[jj] != 0)
1075 add_j_to_nblist(coul_free, jj, bLR);
1081 /* This is done whether or not bWater is set */
1082 if (charge[jj] != 0)
1084 add_j_to_nblist(coul, jj, bLR);
1087 else if (!bDoCoul_i_sol)
1089 if (bHaveVdW[type[jj]])
1091 add_j_to_nblist(vdw, jj, bLR);
1096 if (bHaveVdW[type[jj]])
1098 if (charge[jj] != 0)
1100 add_j_to_nblist(vdwc, jj, bLR);
1104 add_j_to_nblist(vdw, jj, bLR);
1107 else if (charge[jj] != 0)
1109 add_j_to_nblist(coul, jj, bLR);
1117 close_i_nblist(vdw);
1118 close_i_nblist(coul);
1119 close_i_nblist(vdwc);
1120 close_i_nblist(vdw_free);
1121 close_i_nblist(coul_free);
1122 close_i_nblist(vdwc_free);
1128 put_in_list_adress(gmx_bool bHaveVdW[],
1144 /* The a[] index has been removed,
1145 * to put it back in i_atom should be a[i0] and jj should be a[jj].
1150 t_nblist * vdwc_adress = NULL;
1151 t_nblist * vdw_adress = NULL;
1152 t_nblist * coul_adress = NULL;
1153 t_nblist * vdwc_ww = NULL;
1154 t_nblist * coul_ww = NULL;
1156 int i, j, jcg, igid, gid, nbl_ind, nbl_ind_adress;
1157 atom_id jj, jj0, jj1, i_atom;
1162 real *charge, *chargeB;
1164 real qi, qiB, qq, rlj;
1165 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
1166 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
1168 gmx_bool j_all_atom;
1170 t_nblist *nlist, *nlist_adress;
1171 gmx_bool bEnergyGroupCG;
1173 /* Copy some pointers */
1174 cginfo = fr->cginfo;
1175 charge = md->chargeA;
1176 chargeB = md->chargeB;
1179 bPert = md->bPerturbed;
1182 /* Get atom range */
1184 nicg = index[icg+1]-i0;
1186 /* Get the i charge group info */
1187 igid = GET_CGINFO_GID(cginfo[icg]);
1189 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
1193 gmx_fatal(FARGS, "AdResS does not support free energy pertubation\n");
1196 /* Unpack pointers to neighbourlist structs */
1197 if (fr->nnblists == 2)
1204 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
1205 nbl_ind_adress = nbl_ind+fr->nnblists/2;
1209 nlist = fr->nblists[nbl_ind].nlist_lr;
1210 nlist_adress = fr->nblists[nbl_ind_adress].nlist_lr;
1214 nlist = fr->nblists[nbl_ind].nlist_sr;
1215 nlist_adress = fr->nblists[nbl_ind_adress].nlist_sr;
1219 vdwc = &nlist[eNL_VDWQQ];
1220 vdw = &nlist[eNL_VDW];
1221 coul = &nlist[eNL_QQ];
1223 vdwc_adress = &nlist_adress[eNL_VDWQQ];
1224 vdw_adress = &nlist_adress[eNL_VDW];
1225 coul_adress = &nlist_adress[eNL_QQ];
1227 /* We do not support solvent optimization with AdResS for now.
1228 For this we would need hybrid solvent-other kernels */
1230 /* no solvent as i charge group */
1231 /* Loop over the atoms in the i charge group */
1232 for (i = 0; i < nicg; i++)
1235 gid = GID(igid, jgid, ngid);
1236 qi = charge[i_atom];
1238 /* Create new i_atom for each energy group */
1239 if (bDoVdW && bDoCoul)
1241 new_i_nblist(vdwc, bLR, i_atom, shift, gid);
1242 new_i_nblist(vdwc_adress, bLR, i_atom, shift, gid);
1247 new_i_nblist(vdw, bLR, i_atom, shift, gid);
1248 new_i_nblist(vdw_adress, bLR, i_atom, shift, gid);
1253 new_i_nblist(coul, bLR, i_atom, shift, gid);
1254 new_i_nblist(coul_adress, bLR, i_atom, shift, gid);
1256 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
1257 bDoCoul_i = (bDoCoul && qi != 0);
1259 /* Here we find out whether the energy groups interaction belong to a
1260 * coarse-grained (vsite) or atomistic interaction. Note that, beacuse
1261 * interactions between coarse-grained and other (atomistic) energygroups
1262 * are excluded automatically by grompp, it is sufficient to check for
1263 * the group id of atom i (igid) */
1264 bEnergyGroupCG = !egp_explicit(fr, igid);
1266 if (bDoVdW_i || bDoCoul_i)
1268 /* Loop over the j charge groups */
1269 for (j = 0; (j < nj); j++)
1273 /* Check for large charge groups */
1284 /* Finally loop over the atoms in the j-charge group */
1285 for (jj = jj0; jj < jj1; jj++)
1287 bNotEx = NOTEXCL(bExcl, i, jj);
1289 /* Now we have to exclude interactions which will be zero
1290 * anyway due to the AdResS weights (in previous implementations
1291 * this was done in the force kernel). This is necessary as
1292 * pure interactions (those with b_hybrid=false, i.e. w_i*w_j==1 or 0)
1293 * are put into neighbour lists which will be passed to the
1294 * standard (optimized) kernels for speed. The interactions with
1295 * b_hybrid=true are placed into the _adress neighbour lists and
1296 * processed by the generic AdResS kernel.
1298 if ( (bEnergyGroupCG &&
1299 wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS ) ||
1300 ( !bEnergyGroupCG && wf[jj] <= GMX_REAL_EPS ) )
1305 b_hybrid = !((wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS) ||
1306 (wf[i_atom] <= GMX_REAL_EPS && wf[jj] <= GMX_REAL_EPS));
1312 if (charge[jj] != 0)
1316 add_j_to_nblist(coul, jj, bLR);
1320 add_j_to_nblist(coul_adress, jj, bLR);
1324 else if (!bDoCoul_i)
1326 if (bHaveVdW[type[jj]])
1330 add_j_to_nblist(vdw, jj, bLR);
1334 add_j_to_nblist(vdw_adress, jj, bLR);
1340 if (bHaveVdW[type[jj]])
1342 if (charge[jj] != 0)
1346 add_j_to_nblist(vdwc, jj, bLR);
1350 add_j_to_nblist(vdwc_adress, jj, bLR);
1357 add_j_to_nblist(vdw, jj, bLR);
1361 add_j_to_nblist(vdw_adress, jj, bLR);
1366 else if (charge[jj] != 0)
1370 add_j_to_nblist(coul, jj, bLR);
1374 add_j_to_nblist(coul_adress, jj, bLR);
1383 close_i_nblist(vdw);
1384 close_i_nblist(coul);
1385 close_i_nblist(vdwc);
1386 close_i_nblist(vdw_adress);
1387 close_i_nblist(coul_adress);
1388 close_i_nblist(vdwc_adress);
1394 put_in_list_qmmm(gmx_bool bHaveVdW[],
1411 int i, j, jcg, igid, gid;
1412 atom_id jj, jj0, jj1, i_atom;
1416 /* Get atom range */
1418 nicg = index[icg+1]-i0;
1420 /* Get the i charge group info */
1421 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1423 coul = &fr->QMMMlist;
1425 /* Loop over atoms in the ith charge group */
1426 for (i = 0; i < nicg; i++)
1429 gid = GID(igid, jgid, ngid);
1430 /* Create new i_atom for each energy group */
1431 new_i_nblist(coul, bLR, i_atom, shift, gid);
1433 /* Loop over the j charge groups */
1434 for (j = 0; j < nj; j++)
1438 /* Charge groups cannot have QM and MM atoms simultaneously */
1443 /* Finally loop over the atoms in the j-charge group */
1444 for (jj = jj0; jj < jj1; jj++)
1446 bNotEx = NOTEXCL(bExcl, i, jj);
1449 add_j_to_nblist(coul, jj, bLR);
1454 close_i_nblist(coul);
1459 put_in_list_cg(gmx_bool bHaveVdW[],
1476 int igid, gid, nbl_ind;
1480 cginfo = fr->cginfo[icg];
1482 igid = GET_CGINFO_GID(cginfo);
1483 gid = GID(igid, jgid, ngid);
1485 /* Unpack pointers to neighbourlist structs */
1486 if (fr->nnblists == 1)
1492 nbl_ind = fr->gid2nblists[gid];
1496 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1500 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1503 /* Make a new neighbor list for charge group icg.
1504 * Currently simply one neighbor list is made with LJ and Coulomb.
1505 * If required, zero interactions could be removed here
1506 * or in the force loop.
1508 new_i_nblist(vdwc, bLR, index[icg], shift, gid);
1509 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1511 for (j = 0; (j < nj); j++)
1514 /* Skip the icg-icg pairs if all self interactions are excluded */
1515 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1517 /* Here we add the j charge group jcg to the list,
1518 * exclusions are also added to the list.
1520 add_j_to_nblist_cg(vdwc, index[jcg], index[jcg+1], bExcl, icg == jcg, bLR);
1524 close_i_nblist(vdwc);
1527 static void setexcl(atom_id start, atom_id end, t_blocka *excl, gmx_bool b,
1534 for (i = start; i < end; i++)
1536 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1538 SETEXCL(bexcl, i-start, excl->a[k]);
1544 for (i = start; i < end; i++)
1546 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1548 RMEXCL(bexcl, i-start, excl->a[k]);
1554 int calc_naaj(int icg, int cgtot)
1558 if ((cgtot % 2) == 1)
1560 /* Odd number of charge groups, easy */
1561 naaj = 1 + (cgtot/2);
1563 else if ((cgtot % 4) == 0)
1565 /* Multiple of four is hard */
1602 fprintf(log, "naaj=%d\n", naaj);
1608 /************************************************
1610 * S I M P L E C O R E S T U F F
1612 ************************************************/
1614 static real calc_image_tric(rvec xi, rvec xj, matrix box,
1615 rvec b_inv, int *shift)
1617 /* This code assumes that the cut-off is smaller than
1618 * a half times the smallest diagonal element of the box.
1620 const real h25 = 2.5;
1625 /* Compute diff vector */
1626 dz = xj[ZZ] - xi[ZZ];
1627 dy = xj[YY] - xi[YY];
1628 dx = xj[XX] - xi[XX];
1630 /* Perform NINT operation, using trunc operation, therefore
1631 * we first add 2.5 then subtract 2 again
1633 tz = dz*b_inv[ZZ] + h25;
1635 dz -= tz*box[ZZ][ZZ];
1636 dy -= tz*box[ZZ][YY];
1637 dx -= tz*box[ZZ][XX];
1639 ty = dy*b_inv[YY] + h25;
1641 dy -= ty*box[YY][YY];
1642 dx -= ty*box[YY][XX];
1644 tx = dx*b_inv[XX]+h25;
1646 dx -= tx*box[XX][XX];
1648 /* Distance squared */
1649 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1651 *shift = XYZ2IS(tx, ty, tz);
1656 static real calc_image_rect(rvec xi, rvec xj, rvec box_size,
1657 rvec b_inv, int *shift)
1659 const real h15 = 1.5;
1665 /* Compute diff vector */
1666 dx = xj[XX] - xi[XX];
1667 dy = xj[YY] - xi[YY];
1668 dz = xj[ZZ] - xi[ZZ];
1670 /* Perform NINT operation, using trunc operation, therefore
1671 * we first add 1.5 then subtract 1 again
1673 tx = dx*b_inv[XX] + h15;
1674 ty = dy*b_inv[YY] + h15;
1675 tz = dz*b_inv[ZZ] + h15;
1680 /* Correct diff vector for translation */
1681 ddx = tx*box_size[XX] - dx;
1682 ddy = ty*box_size[YY] - dy;
1683 ddz = tz*box_size[ZZ] - dz;
1685 /* Distance squared */
1686 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1688 *shift = XYZ2IS(tx, ty, tz);
1693 static void add_simple(t_ns_buf *nsbuf, int nrj, atom_id cg_j,
1694 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1695 int icg, int jgid, t_block *cgs, t_excl bexcl[],
1696 int shift, t_forcerec *fr, put_in_list_t *put_in_list)
1698 if (nsbuf->nj + nrj > MAX_CG)
1700 put_in_list(bHaveVdW, ngid, md, icg, jgid, nsbuf->ncg, nsbuf->jcg,
1701 cgs->index, bexcl, shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1702 /* Reset buffer contents */
1703 nsbuf->ncg = nsbuf->nj = 0;
1705 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1709 static void ns_inner_tric(rvec x[], int icg, int *i_egp_flags,
1710 int njcg, atom_id jcg[],
1711 matrix box, rvec b_inv, real rcut2,
1712 t_block *cgs, t_ns_buf **ns_buf,
1713 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1714 t_excl bexcl[], t_forcerec *fr,
1715 put_in_list_t *put_in_list)
1719 int *cginfo = fr->cginfo;
1720 atom_id cg_j, *cgindex;
1723 cgindex = cgs->index;
1725 for (j = 0; (j < njcg); j++)
1728 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1729 if (calc_image_tric(x[icg], x[cg_j], box, b_inv, &shift) < rcut2)
1731 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1732 if (!(i_egp_flags[jgid] & EGP_EXCL))
1734 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1735 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1742 static void ns_inner_rect(rvec x[], int icg, int *i_egp_flags,
1743 int njcg, atom_id jcg[],
1744 gmx_bool bBox, rvec box_size, rvec b_inv, real rcut2,
1745 t_block *cgs, t_ns_buf **ns_buf,
1746 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1747 t_excl bexcl[], t_forcerec *fr,
1748 put_in_list_t *put_in_list)
1752 int *cginfo = fr->cginfo;
1753 atom_id cg_j, *cgindex;
1756 cgindex = cgs->index;
1760 for (j = 0; (j < njcg); j++)
1763 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1764 if (calc_image_rect(x[icg], x[cg_j], box_size, b_inv, &shift) < rcut2)
1766 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1767 if (!(i_egp_flags[jgid] & EGP_EXCL))
1769 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1770 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1778 for (j = 0; (j < njcg); j++)
1781 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1782 if ((rcut2 == 0) || (distance2(x[icg], x[cg_j]) < rcut2))
1784 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1785 if (!(i_egp_flags[jgid] & EGP_EXCL))
1787 add_simple(&ns_buf[jgid][CENTRAL], nrj, cg_j,
1788 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, CENTRAL, fr,
1796 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1798 static int ns_simple_core(t_forcerec *fr,
1799 gmx_localtop_t *top,
1801 matrix box, rvec box_size,
1802 t_excl bexcl[], atom_id *aaj,
1803 int ngid, t_ns_buf **ns_buf,
1804 put_in_list_t *put_in_list, gmx_bool bHaveVdW[])
1808 int nsearch, icg, jcg, igid, i0, nri, nn;
1811 /* atom_id *i_atoms; */
1812 t_block *cgs = &(top->cgs);
1813 t_blocka *excl = &(top->excls);
1816 gmx_bool bBox, bTriclinic;
1819 rlist2 = sqr(fr->rlist);
1821 bBox = (fr->ePBC != epbcNONE);
1824 for (m = 0; (m < DIM); m++)
1826 b_inv[m] = divide(1.0, box_size[m]);
1828 bTriclinic = TRICLINIC(box);
1835 cginfo = fr->cginfo;
1838 for (icg = fr->cg0; (icg < fr->hcg); icg++)
1841 i0 = cgs->index[icg];
1842 nri = cgs->index[icg+1]-i0;
1843 i_atoms = &(cgs->a[i0]);
1844 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1845 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1847 igid = GET_CGINFO_GID(cginfo[icg]);
1848 i_egp_flags = fr->egp_flags + ngid*igid;
1849 setexcl(cgs->index[icg], cgs->index[icg+1], excl, TRUE, bexcl);
1851 naaj = calc_naaj(icg, cgs->nr);
1854 ns_inner_tric(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1855 box, b_inv, rlist2, cgs, ns_buf,
1856 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1860 ns_inner_rect(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1861 bBox, box_size, b_inv, rlist2, cgs, ns_buf,
1862 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1866 for (nn = 0; (nn < ngid); nn++)
1868 for (k = 0; (k < SHIFTS); k++)
1870 nsbuf = &(ns_buf[nn][k]);
1873 put_in_list(bHaveVdW, ngid, md, icg, nn, nsbuf->ncg, nsbuf->jcg,
1874 cgs->index, bexcl, k, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1875 nsbuf->ncg = nsbuf->nj = 0;
1879 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1880 setexcl(cgs->index[icg], cgs->index[icg+1], excl, FALSE, bexcl);
1882 close_neighbor_lists(fr, FALSE);
1887 /************************************************
1889 * N S 5 G R I D S T U F F
1891 ************************************************/
1893 static inline void get_dx(int Nx, real gridx, real rc2, int xgi, real x,
1894 int *dx0, int *dx1, real *dcx2)
1922 for (i = xgi0; i >= 0; i--)
1924 dcx = (i+1)*gridx-x;
1933 for (i = xgi1; i < Nx; i++)
1946 static inline void get_dx_dd(int Nx, real gridx, real rc2, int xgi, real x,
1947 int ncpddc, int shift_min, int shift_max,
1948 int *g0, int *g1, real *dcx2)
1951 int g_min, g_max, shift_home;
1984 g_min = (shift_min == shift_home ? 0 : ncpddc);
1985 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1992 else if (shift_max < 0)
2007 /* Check one grid cell down */
2008 dcx = ((*g0 - 1) + 1)*gridx - x;
2020 /* Check one grid cell up */
2021 dcx = (*g1 + 1)*gridx - x;
2033 #define sqr(x) ((x)*(x))
2034 #define calc_dx2(XI, YI, ZI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
2035 #define calc_cyl_dx2(XI, YI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
2036 /****************************************************
2038 * F A S T N E I G H B O R S E A R C H I N G
2040 * Optimized neighboursearching routine using grid
2041 * at least 1x1x1, see GROMACS manual
2043 ****************************************************/
2046 static void get_cutoff2(t_forcerec *fr, gmx_bool bDoLongRange,
2047 real *rvdw2, real *rcoul2,
2048 real *rs2, real *rm2, real *rl2)
2050 *rs2 = sqr(fr->rlist);
2052 if (bDoLongRange && fr->bTwinRange)
2054 /* The VdW and elec. LR cut-off's could be different,
2055 * so we can not simply set them to rlistlong.
2057 if (EVDW_MIGHT_BE_ZERO_AT_CUTOFF(fr->vdwtype) &&
2058 fr->rvdw > fr->rlist)
2060 *rvdw2 = sqr(fr->rlistlong);
2064 *rvdw2 = sqr(fr->rvdw);
2066 if (EEL_MIGHT_BE_ZERO_AT_CUTOFF(fr->eeltype) &&
2067 fr->rcoulomb > fr->rlist)
2069 *rcoul2 = sqr(fr->rlistlong);
2073 *rcoul2 = sqr(fr->rcoulomb);
2078 /* Workaround for a gcc -O3 or -ffast-math problem */
2082 *rm2 = min(*rvdw2, *rcoul2);
2083 *rl2 = max(*rvdw2, *rcoul2);
2086 static void init_nsgrid_lists(t_forcerec *fr, int ngid, gmx_ns_t *ns)
2088 real rvdw2, rcoul2, rs2, rm2, rl2;
2091 get_cutoff2(fr, TRUE, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2093 /* Short range buffers */
2094 snew(ns->nl_sr, ngid);
2096 snew(ns->nsr, ngid);
2097 snew(ns->nlr_ljc, ngid);
2098 snew(ns->nlr_one, ngid);
2100 /* Always allocate both list types, since rcoulomb might now change with PME load balancing */
2101 /* Long range VdW and Coul buffers */
2102 snew(ns->nl_lr_ljc, ngid);
2103 /* Long range VdW or Coul only buffers */
2104 snew(ns->nl_lr_one, ngid);
2106 for (j = 0; (j < ngid); j++)
2108 snew(ns->nl_sr[j], MAX_CG);
2109 snew(ns->nl_lr_ljc[j], MAX_CG);
2110 snew(ns->nl_lr_one[j], MAX_CG);
2115 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
2120 static int nsgrid_core(FILE *log, t_commrec *cr, t_forcerec *fr,
2121 matrix box, rvec box_size, int ngid,
2122 gmx_localtop_t *top,
2123 t_grid *grid, rvec x[],
2124 t_excl bexcl[], gmx_bool *bExcludeAlleg,
2125 t_nrnb *nrnb, t_mdatoms *md,
2126 real *lambda, real *dvdlambda,
2127 gmx_grppairener_t *grppener,
2128 put_in_list_t *put_in_list,
2129 gmx_bool bHaveVdW[],
2130 gmx_bool bDoLongRange, gmx_bool bMakeQMMMnblist)
2133 atom_id **nl_lr_ljc, **nl_lr_one, **nl_sr;
2134 int *nlr_ljc, *nlr_one, *nsr;
2135 gmx_domdec_t *dd = NULL;
2136 t_block *cgs = &(top->cgs);
2137 int *cginfo = fr->cginfo;
2138 /* atom_id *i_atoms,*cgsindex=cgs->index; */
2140 int cell_x, cell_y, cell_z;
2141 int d, tx, ty, tz, dx, dy, dz, cj;
2142 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2143 int zsh_ty, zsh_tx, ysh_tx;
2145 int dx0, dx1, dy0, dy1, dz0, dz1;
2146 int Nx, Ny, Nz, shift = -1, j, nrj, nns, nn = -1;
2147 real gridx, gridy, gridz, grid_x, grid_y, grid_z;
2148 real *dcx2, *dcy2, *dcz2;
2150 int cg0, cg1, icg = -1, cgsnr, i0, igid, nri, naaj, max_jcg;
2151 int jcg0, jcg1, jjcg, cgj0, jgid;
2152 int *grida, *gridnra, *gridind;
2153 gmx_bool rvdw_lt_rcoul, rcoul_lt_rvdw;
2154 rvec xi, *cgcm, grid_offset;
2155 real r2, rs2, rvdw2, rcoul2, rm2, rl2, XI, YI, ZI, dcx, dcy, dcz, tmp1, tmp2;
2157 gmx_bool bDomDec, bTriclinicX, bTriclinicY;
2162 bDomDec = DOMAINDECOMP(cr);
2168 bTriclinicX = ((YY < grid->npbcdim &&
2169 (!bDomDec || dd->nc[YY] == 1) && box[YY][XX] != 0) ||
2170 (ZZ < grid->npbcdim &&
2171 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][XX] != 0));
2172 bTriclinicY = (ZZ < grid->npbcdim &&
2173 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][YY] != 0);
2177 get_cutoff2(fr, bDoLongRange, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2179 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
2180 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
2182 if (bMakeQMMMnblist)
2190 nl_lr_ljc = ns->nl_lr_ljc;
2191 nl_lr_one = ns->nl_lr_one;
2192 nlr_ljc = ns->nlr_ljc;
2193 nlr_one = ns->nlr_one;
2201 gridind = grid->index;
2202 gridnra = grid->nra;
2205 gridx = grid->cell_size[XX];
2206 gridy = grid->cell_size[YY];
2207 gridz = grid->cell_size[ZZ];
2211 copy_rvec(grid->cell_offset, grid_offset);
2212 copy_ivec(grid->ncpddc, ncpddc);
2217 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2218 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
2219 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
2220 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
2221 if (zsh_tx != 0 && ysh_tx != 0)
2223 /* This could happen due to rounding, when both ratios are 0.5 */
2232 /* We only want a list for the test particle */
2241 /* Set the shift range */
2242 for (d = 0; d < DIM; d++)
2246 /* Check if we need periodicity shifts.
2247 * Without PBC or with domain decomposition we don't need them.
2249 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2256 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
2267 /* Loop over charge groups */
2268 for (icg = cg0; (icg < cg1); icg++)
2270 igid = GET_CGINFO_GID(cginfo[icg]);
2271 /* Skip this charge group if all energy groups are excluded! */
2272 if (bExcludeAlleg[igid])
2277 i0 = cgs->index[icg];
2279 if (bMakeQMMMnblist)
2281 /* Skip this charge group if it is not a QM atom while making a
2282 * QM/MM neighbourlist
2284 if (md->bQM[i0] == FALSE)
2286 continue; /* MM particle, go to next particle */
2289 /* Compute the number of charge groups that fall within the control
2292 naaj = calc_naaj(icg, cgsnr);
2299 /* make a normal neighbourlist */
2303 /* Get the j charge-group and dd cell shift ranges */
2304 dd_get_ns_ranges(cr->dd, icg, &jcg0, &jcg1, sh0, sh1);
2309 /* Compute the number of charge groups that fall within the control
2312 naaj = calc_naaj(icg, cgsnr);
2318 /* The i-particle is awlways the test particle,
2319 * so we want all j-particles
2321 max_jcg = cgsnr - 1;
2325 max_jcg = jcg1 - cgsnr;
2330 i_egp_flags = fr->egp_flags + igid*ngid;
2332 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2333 setexcl(i0, cgs->index[icg+1], &top->excls, TRUE, bexcl);
2335 ci2xyz(grid, icg, &cell_x, &cell_y, &cell_z);
2337 /* Changed iicg to icg, DvdS 990115
2338 * (but see consistency check above, DvdS 990330)
2341 fprintf(log, "icg=%5d, naaj=%5d, cell %d %d %d\n",
2342 icg, naaj, cell_x, cell_y, cell_z);
2344 /* Loop over shift vectors in three dimensions */
2345 for (tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
2347 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2348 /* Calculate range of cells in Z direction that have the shift tz */
2349 zgi = cell_z + tz*Nz;
2352 get_dx(Nz, gridz, rl2, zgi, ZI, &dz0, &dz1, dcz2);
2354 get_dx_dd(Nz, gridz, rl2, zgi, ZI-grid_offset[ZZ],
2355 ncpddc[ZZ], sh0[ZZ], sh1[ZZ], &dz0, &dz1, dcz2);
2361 for (ty = -shp[YY]; ty <= shp[YY]; ty++)
2363 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2364 /* Calculate range of cells in Y direction that have the shift ty */
2367 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2371 ygi = cell_y + ty*Ny;
2374 get_dx(Ny, gridy, rl2, ygi, YI, &dy0, &dy1, dcy2);
2376 get_dx_dd(Ny, gridy, rl2, ygi, YI-grid_offset[YY],
2377 ncpddc[YY], sh0[YY], sh1[YY], &dy0, &dy1, dcy2);
2383 for (tx = -shp[XX]; tx <= shp[XX]; tx++)
2385 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2386 /* Calculate range of cells in X direction that have the shift tx */
2389 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2393 xgi = cell_x + tx*Nx;
2396 get_dx(Nx, gridx, rl2, xgi*Nx, XI, &dx0, &dx1, dcx2);
2398 get_dx_dd(Nx, gridx, rl2, xgi, XI-grid_offset[XX],
2399 ncpddc[XX], sh0[XX], sh1[XX], &dx0, &dx1, dcx2);
2405 /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2406 * from the neigbour list as it will not interact */
2407 if (fr->adress_type != eAdressOff)
2409 if (md->wf[cgs->index[icg]] <= GMX_REAL_EPS && egp_explicit(fr, igid))
2414 /* Get shift vector */
2415 shift = XYZ2IS(tx, ty, tz);
2417 range_check(shift, 0, SHIFTS);
2419 for (nn = 0; (nn < ngid); nn++)
2426 fprintf(log, "shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2427 shift, dx0, dx1, dy0, dy1, dz0, dz1);
2428 fprintf(log, "cgcm: %8.3f %8.3f %8.3f\n", cgcm[icg][XX],
2429 cgcm[icg][YY], cgcm[icg][ZZ]);
2430 fprintf(log, "xi: %8.3f %8.3f %8.3f\n", XI, YI, ZI);
2432 for (dx = dx0; (dx <= dx1); dx++)
2434 tmp1 = rl2 - dcx2[dx];
2435 for (dy = dy0; (dy <= dy1); dy++)
2437 tmp2 = tmp1 - dcy2[dy];
2440 for (dz = dz0; (dz <= dz1); dz++)
2442 if (tmp2 > dcz2[dz])
2444 /* Find grid-cell cj in which possible neighbours are */
2445 cj = xyz2ci(Ny, Nz, dx, dy, dz);
2447 /* Check out how many cgs (nrj) there in this cell */
2450 /* Find the offset in the cg list */
2453 /* Check if all j's are out of range so we
2454 * can skip the whole cell.
2455 * Should save some time, especially with DD.
2458 (grida[cgj0] >= max_jcg &&
2459 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2465 for (j = 0; (j < nrj); j++)
2467 jjcg = grida[cgj0+j];
2469 /* check whether this guy is in range! */
2470 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2473 r2 = calc_dx2(XI, YI, ZI, cgcm[jjcg]);
2476 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2477 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2478 /* check energy group exclusions */
2479 if (!(i_egp_flags[jgid] & EGP_EXCL))
2483 if (nsr[jgid] >= MAX_CG)
2485 /* Add to short-range list */
2486 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2487 nsr[jgid], nl_sr[jgid],
2488 cgs->index, /* cgsatoms, */ bexcl,
2489 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2492 nl_sr[jgid][nsr[jgid]++] = jjcg;
2496 if (nlr_ljc[jgid] >= MAX_CG)
2498 /* Add to LJ+coulomb long-range list */
2499 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2500 nlr_ljc[jgid], nl_lr_ljc[jgid], top->cgs.index,
2501 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2504 nl_lr_ljc[jgid][nlr_ljc[jgid]++] = jjcg;
2508 if (nlr_one[jgid] >= MAX_CG)
2510 /* Add to long-range list with only coul, or only LJ */
2511 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2512 nlr_one[jgid], nl_lr_one[jgid], top->cgs.index,
2513 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2516 nl_lr_one[jgid][nlr_one[jgid]++] = jjcg;
2528 /* CHECK whether there is anything left in the buffers */
2529 for (nn = 0; (nn < ngid); nn++)
2533 put_in_list(bHaveVdW, ngid, md, icg, nn, nsr[nn], nl_sr[nn],
2534 cgs->index, /* cgsatoms, */ bexcl,
2535 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2538 if (nlr_ljc[nn] > 0)
2540 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_ljc[nn],
2541 nl_lr_ljc[nn], top->cgs.index,
2542 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2545 if (nlr_one[nn] > 0)
2547 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_one[nn],
2548 nl_lr_one[nn], top->cgs.index,
2549 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2555 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2556 setexcl(cgs->index[icg], cgs->index[icg+1], &top->excls, FALSE, bexcl);
2558 /* No need to perform any left-over force calculations anymore (as we used to do here)
2559 * since we now save the proper long-range lists for later evaluation.
2564 /* Close neighbourlists */
2565 close_neighbor_lists(fr, bMakeQMMMnblist);
2570 void ns_realloc_natoms(gmx_ns_t *ns, int natoms)
2574 if (natoms > ns->nra_alloc)
2576 ns->nra_alloc = over_alloc_dd(natoms);
2577 srenew(ns->bexcl, ns->nra_alloc);
2578 for (i = 0; i < ns->nra_alloc; i++)
2585 void init_ns(FILE *fplog, const t_commrec *cr,
2586 gmx_ns_t *ns, t_forcerec *fr,
2587 const gmx_mtop_t *mtop,
2590 int mt, icg, nr_in_cg, maxcg, i, j, jcg, ngid, ncg;
2594 /* Compute largest charge groups size (# atoms) */
2596 for (mt = 0; mt < mtop->nmoltype; mt++)
2598 cgs = &mtop->moltype[mt].cgs;
2599 for (icg = 0; (icg < cgs->nr); icg++)
2601 nr_in_cg = max(nr_in_cg, (int)(cgs->index[icg+1]-cgs->index[icg]));
2605 /* Verify whether largest charge group is <= max cg.
2606 * This is determined by the type of the local exclusion type
2607 * Exclusions are stored in bits. (If the type is not large
2608 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2610 maxcg = sizeof(t_excl)*8;
2611 if (nr_in_cg > maxcg)
2613 gmx_fatal(FARGS, "Max #atoms in a charge group: %d > %d\n",
2617 ngid = mtop->groups.grps[egcENER].nr;
2618 snew(ns->bExcludeAlleg, ngid);
2619 for (i = 0; i < ngid; i++)
2621 ns->bExcludeAlleg[i] = TRUE;
2622 for (j = 0; j < ngid; j++)
2624 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2626 ns->bExcludeAlleg[i] = FALSE;
2634 ns->grid = init_grid(fplog, fr);
2635 init_nsgrid_lists(fr, ngid, ns);
2640 snew(ns->ns_buf, ngid);
2641 for (i = 0; (i < ngid); i++)
2643 snew(ns->ns_buf[i], SHIFTS);
2645 ncg = ncg_mtop(mtop);
2646 snew(ns->simple_aaj, 2*ncg);
2647 for (jcg = 0; (jcg < ncg); jcg++)
2649 ns->simple_aaj[jcg] = jcg;
2650 ns->simple_aaj[jcg+ncg] = jcg;
2654 /* Create array that determines whether or not atoms have VdW */
2655 snew(ns->bHaveVdW, fr->ntype);
2656 for (i = 0; (i < fr->ntype); i++)
2658 for (j = 0; (j < fr->ntype); j++)
2660 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2662 ((BHAMA(fr->nbfp, fr->ntype, i, j) != 0) ||
2663 (BHAMB(fr->nbfp, fr->ntype, i, j) != 0) ||
2664 (BHAMC(fr->nbfp, fr->ntype, i, j) != 0)) :
2665 ((C6(fr->nbfp, fr->ntype, i, j) != 0) ||
2666 (C12(fr->nbfp, fr->ntype, i, j) != 0))));
2671 pr_bvec(debug, 0, "bHaveVdW", ns->bHaveVdW, fr->ntype, TRUE);
2676 if (!DOMAINDECOMP(cr))
2678 /* This could be reduced with particle decomposition */
2679 ns_realloc_natoms(ns, mtop->natoms);
2682 ns->nblist_initialized = FALSE;
2684 /* nbr list debug dump */
2686 char *ptr = getenv("GMX_DUMP_NL");
2689 ns->dump_nl = strtol(ptr, NULL, 10);
2692 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2703 int search_neighbours(FILE *log, t_forcerec *fr,
2704 rvec x[], matrix box,
2705 gmx_localtop_t *top,
2706 gmx_groups_t *groups,
2708 t_nrnb *nrnb, t_mdatoms *md,
2709 real *lambda, real *dvdlambda,
2710 gmx_grppairener_t *grppener,
2712 gmx_bool bDoLongRangeNS,
2713 gmx_bool bPadListsForKernels)
2715 t_block *cgs = &(top->cgs);
2716 rvec box_size, grid_x0, grid_x1;
2718 real min_size, grid_dens;
2722 gmx_bool *i_egp_flags;
2723 int cg_start, cg_end, start, end;
2726 gmx_domdec_zones_t *dd_zones;
2727 put_in_list_t *put_in_list;
2731 /* Set some local variables */
2733 ngid = groups->grps[egcENER].nr;
2735 for (m = 0; (m < DIM); m++)
2737 box_size[m] = box[m][m];
2740 if (fr->ePBC != epbcNONE)
2742 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC, box))
2744 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.");
2748 min_size = min(box_size[XX], min(box_size[YY], box_size[ZZ]));
2749 if (2*fr->rlistlong >= min_size)
2751 gmx_fatal(FARGS, "One of the box diagonal elements has become smaller than twice the cut-off length.");
2756 if (DOMAINDECOMP(cr))
2758 ns_realloc_natoms(ns, cgs->index[cgs->nr]);
2762 /* Reset the neighbourlists */
2763 reset_neighbor_lists(fr, TRUE, TRUE);
2765 if (bGrid && bFillGrid)
2769 if (DOMAINDECOMP(cr))
2771 dd_zones = domdec_zones(cr->dd);
2777 get_nsgrid_boundaries(grid->nboundeddim, box, NULL, NULL, NULL, NULL,
2778 cgs->nr, fr->cg_cm, grid_x0, grid_x1, &grid_dens);
2780 grid_first(log, grid, NULL, NULL, fr->ePBC, box, grid_x0, grid_x1,
2781 fr->rlistlong, grid_dens);
2785 /* Don't know why this all is... (DvdS 3/99) */
2791 end = (cgs->nr+1)/2;
2794 if (DOMAINDECOMP(cr))
2797 fill_grid(log, dd_zones, grid, end, -1, end, fr->cg_cm);
2799 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2803 fill_grid(log, NULL, grid, cgs->nr, fr->cg0, fr->hcg, fr->cg_cm);
2804 grid->icg0 = fr->cg0;
2805 grid->icg1 = fr->hcg;
2815 calc_elemnr(log, grid, start, end, cgs->nr);
2817 grid_last(log, grid, start, end, cgs->nr);
2821 check_grid(debug, grid);
2822 print_grid(debug, grid);
2827 /* Set the grid cell index for the test particle only.
2828 * The cell to cg index is not corrected, but that does not matter.
2830 fill_grid(log, NULL, ns->grid, fr->hcg, fr->hcg-1, fr->hcg, fr->cg_cm);
2834 if (fr->adress_type == eAdressOff)
2836 if (!fr->ns.bCGlist)
2838 put_in_list = put_in_list_at;
2842 put_in_list = put_in_list_cg;
2847 put_in_list = put_in_list_adress;
2854 nsearch = nsgrid_core(log, cr, fr, box, box_size, ngid, top,
2855 grid, x, ns->bexcl, ns->bExcludeAlleg,
2856 nrnb, md, lambda, dvdlambda, grppener,
2857 put_in_list, ns->bHaveVdW,
2858 bDoLongRangeNS, FALSE);
2860 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2861 * the classical calculation. The charge-charge interaction
2862 * between QM and MM atoms is handled in the QMMM core calculation
2863 * (see QMMM.c). The VDW however, we'd like to compute classically
2864 * and the QM MM atom pairs have just been put in the
2865 * corresponding neighbourlists. in case of QMMM we still need to
2866 * fill a special QMMM neighbourlist that contains all neighbours
2867 * of the QM atoms. If bQMMM is true, this list will now be made:
2869 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
2871 nsearch += nsgrid_core(log, cr, fr, box, box_size, ngid, top,
2872 grid, x, ns->bexcl, ns->bExcludeAlleg,
2873 nrnb, md, lambda, dvdlambda, grppener,
2874 put_in_list_qmmm, ns->bHaveVdW,
2875 bDoLongRangeNS, TRUE);
2880 nsearch = ns_simple_core(fr, top, md, box, box_size,
2881 ns->bexcl, ns->simple_aaj,
2882 ngid, ns->ns_buf, put_in_list, ns->bHaveVdW);
2890 inc_nrnb(nrnb, eNR_NS, nsearch);
2891 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2896 int natoms_beyond_ns_buffer(t_inputrec *ir, t_forcerec *fr, t_block *cgs,
2897 matrix scale_tot, rvec *x)
2899 int cg0, cg1, cg, a0, a1, a, i, j;
2900 real rint, hbuf2, scale;
2902 gmx_bool bIsotropic;
2907 rint = max(ir->rcoulomb, ir->rvdw);
2908 if (ir->rlist < rint)
2910 gmx_fatal(FARGS, "The neighbor search buffer has negative size: %f nm",
2918 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2920 hbuf2 = sqr(0.5*(ir->rlist - rint));
2921 for (cg = cg0; cg < cg1; cg++)
2923 a0 = cgs->index[cg];
2924 a1 = cgs->index[cg+1];
2925 for (a = a0; a < a1; a++)
2927 if (distance2(cg_cm[cg], x[a]) > hbuf2)
2937 scale = scale_tot[0][0];
2938 for (i = 1; i < DIM; i++)
2940 /* With anisotropic scaling, the original spherical ns volumes become
2941 * ellipsoids. To avoid costly transformations we use the minimum
2942 * eigenvalue of the scaling matrix for determining the buffer size.
2943 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2945 scale = min(scale, scale_tot[i][i]);
2946 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2950 for (j = 0; j < i; j++)
2952 if (scale_tot[i][j] != 0)
2958 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2961 for (cg = cg0; cg < cg1; cg++)
2963 svmul(scale, cg_cm[cg], cgsc);
2964 a0 = cgs->index[cg];
2965 a1 = cgs->index[cg+1];
2966 for (a = a0; a < a1; a++)
2968 if (distance2(cgsc, x[a]) > hbuf2)
2977 /* Anistropic scaling */
2978 for (cg = cg0; cg < cg1; cg++)
2980 /* Since scale_tot contains the transpose of the scaling matrix,
2981 * we need to multiply with the transpose.
2983 tmvmul_ur0(scale_tot, cg_cm[cg], cgsc);
2984 a0 = cgs->index[cg];
2985 a1 = cgs->index[cg+1];
2986 for (a = a0; a < a1; a++)
2988 if (distance2(cgsc, x[a]) > hbuf2)