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;
293 fprintf(log, "Note: The available nonbonded kernels do not support water optimization - disabling.\n");
296 if (fr->efep != efepNO)
298 if ((fr->bEwald) && (fr->sc_alphacoul > 0)) /* need to handle long range differently if using softcore */
300 ielecf = GMX_NBKERNEL_ELEC_EWALD;
301 ielecmodf = eintmodNONE;
306 ielecmodf = ielecmod;
309 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_FREE], &nbl->nlist_lr[eNL_VDWQQ_FREE],
310 maxsr, maxlr, ivdw, ivdwmod, ielecf, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
311 init_nblist(log, &nbl->nlist_sr[eNL_VDW_FREE], &nbl->nlist_lr[eNL_VDW_FREE],
312 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
313 init_nblist(log, &nbl->nlist_sr[eNL_QQ_FREE], &nbl->nlist_lr[eNL_QQ_FREE],
314 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielecf, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
318 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
320 init_nblist(log, &fr->QMMMlist, NULL,
321 maxsr, maxlr, 0, 0, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_STANDARD);
329 fr->ns.nblist_initialized = TRUE;
332 static void reset_nblist(t_nblist *nl)
343 static void reset_neighbor_lists(t_forcerec *fr, gmx_bool bResetSR, gmx_bool bResetLR)
349 /* only reset the short-range nblist */
350 reset_nblist(&(fr->QMMMlist));
353 for (n = 0; n < fr->nnblists; n++)
355 for (i = 0; i < eNL_NR; i++)
359 reset_nblist( &(fr->nblists[n].nlist_sr[i]) );
363 reset_nblist( &(fr->nblists[n].nlist_lr[i]) );
372 static inline void new_i_nblist(t_nblist *nlist,
373 gmx_bool bLR, atom_id i_atom, int shift, int gid)
375 int i, k, nri, nshift;
379 /* Check whether we have to increase the i counter */
381 (nlist->iinr[nri] != i_atom) ||
382 (nlist->shift[nri] != shift) ||
383 (nlist->gid[nri] != gid))
385 /* This is something else. Now see if any entries have
386 * been added in the list of the previous atom.
389 ((nlist->jindex[nri+1] > nlist->jindex[nri]) &&
390 (nlist->gid[nri] != -1)))
392 /* If so increase the counter */
395 if (nlist->nri >= nlist->maxnri)
397 nlist->maxnri += over_alloc_large(nlist->nri);
398 reallocate_nblist(nlist);
401 /* Set the number of neighbours and the atom number */
402 nlist->jindex[nri+1] = nlist->jindex[nri];
403 nlist->iinr[nri] = i_atom;
404 nlist->gid[nri] = gid;
405 nlist->shift[nri] = shift;
409 /* Adding to previous list. First remove possible previous padding */
410 if(nlist->simd_padding_width>1)
412 while(nlist->nrj>0 && nlist->jjnr[nlist->nrj-1]<0)
420 static inline void close_i_nblist(t_nblist *nlist)
422 int nri = nlist->nri;
427 /* Add elements up to padding. Since we allocate memory in units
428 * of the simd_padding width, we do not have to check for possible
429 * list reallocation here.
431 while ((nlist->nrj % nlist->simd_padding_width) != 0)
433 /* Use -4 here, so we can write forces for 4 atoms before real data */
434 nlist->jjnr[nlist->nrj++] = -4;
436 nlist->jindex[nri+1] = nlist->nrj;
438 len = nlist->nrj - nlist->jindex[nri];
440 /* nlist length for water i molecules is treated statically
443 if (len > nlist->maxlen)
450 static inline void close_nblist(t_nblist *nlist)
452 /* Only close this nblist when it has been initialized.
453 * Avoid the creation of i-lists with no j-particles.
457 /* Some assembly kernels do not support empty lists,
458 * make sure here that we don't generate any empty lists.
459 * With the current ns code this branch is taken in two cases:
460 * No i-particles at all: nri=-1 here
461 * There are i-particles, but no j-particles; nri=0 here
467 /* Close list number nri by incrementing the count */
472 static inline void close_neighbor_lists(t_forcerec *fr, gmx_bool bMakeQMMMnblist)
478 close_nblist(&(fr->QMMMlist));
481 for (n = 0; n < fr->nnblists; n++)
483 for (i = 0; (i < eNL_NR); i++)
485 close_nblist(&(fr->nblists[n].nlist_sr[i]));
486 close_nblist(&(fr->nblists[n].nlist_lr[i]));
492 static inline void add_j_to_nblist(t_nblist *nlist, atom_id j_atom, gmx_bool bLR)
494 int nrj = nlist->nrj;
496 if (nlist->nrj >= nlist->maxnrj)
498 nlist->maxnrj = round_up_to_simd_width(over_alloc_small(nlist->nrj + 1), nlist->simd_padding_width);
502 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
503 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
506 srenew(nlist->jjnr, nlist->maxnrj);
509 nlist->jjnr[nrj] = j_atom;
513 static inline void add_j_to_nblist_cg(t_nblist *nlist,
514 atom_id j_start, int j_end,
515 t_excl *bexcl, gmx_bool i_is_j,
518 int nrj = nlist->nrj;
521 if (nlist->nrj >= nlist->maxnrj)
523 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
526 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
527 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
530 srenew(nlist->jjnr, nlist->maxnrj);
531 srenew(nlist->jjnr_end, nlist->maxnrj);
532 srenew(nlist->excl, nlist->maxnrj*MAX_CGCGSIZE);
535 nlist->jjnr[nrj] = j_start;
536 nlist->jjnr_end[nrj] = j_end;
538 if (j_end - j_start > MAX_CGCGSIZE)
540 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);
543 /* Set the exclusions */
544 for (j = j_start; j < j_end; j++)
546 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
550 /* Avoid double counting of intra-cg interactions */
551 for (j = 1; j < j_end-j_start; j++)
553 nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
561 put_in_list_t (gmx_bool bHaveVdW[],
578 put_in_list_at(gmx_bool bHaveVdW[],
594 /* The a[] index has been removed,
595 * to put it back in i_atom should be a[i0] and jj should be a[jj].
600 t_nblist * vdwc_free = NULL;
601 t_nblist * vdw_free = NULL;
602 t_nblist * coul_free = NULL;
603 t_nblist * vdwc_ww = NULL;
604 t_nblist * coul_ww = NULL;
606 int i, j, jcg, igid, gid, nbl_ind, ind_ij;
607 atom_id jj, jj0, jj1, i_atom;
612 real *charge, *chargeB;
613 real qi, qiB, qq, rlj;
614 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
615 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
619 /* Copy some pointers */
621 charge = md->chargeA;
622 chargeB = md->chargeB;
625 bPert = md->bPerturbed;
629 nicg = index[icg+1]-i0;
631 /* Get the i charge group info */
632 igid = GET_CGINFO_GID(cginfo[icg]);
634 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
639 /* Check if any of the particles involved are perturbed.
640 * If not we can do the cheaper normal put_in_list
641 * and use more solvent optimization.
643 for (i = 0; i < nicg; i++)
645 bFreeEnergy |= bPert[i0+i];
647 /* Loop over the j charge groups */
648 for (j = 0; (j < nj && !bFreeEnergy); j++)
653 /* Finally loop over the atoms in the j-charge group */
654 for (jj = jj0; jj < jj1; jj++)
656 bFreeEnergy |= bPert[jj];
661 /* Unpack pointers to neighbourlist structs */
662 if (fr->nnblists == 1)
668 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
672 nlist = fr->nblists[nbl_ind].nlist_lr;
676 nlist = fr->nblists[nbl_ind].nlist_sr;
679 if (iwater != esolNO)
681 vdwc = &nlist[eNL_VDWQQ_WATER];
682 vdw = &nlist[eNL_VDW];
683 coul = &nlist[eNL_QQ_WATER];
684 #ifndef DISABLE_WATERWATER_NLIST
685 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
686 coul_ww = &nlist[eNL_QQ_WATERWATER];
691 vdwc = &nlist[eNL_VDWQQ];
692 vdw = &nlist[eNL_VDW];
693 coul = &nlist[eNL_QQ];
698 if (iwater != esolNO)
700 /* Loop over the atoms in the i charge group */
702 gid = GID(igid, jgid, ngid);
703 /* Create new i_atom for each energy group */
704 if (bDoCoul && bDoVdW)
706 new_i_nblist(vdwc, bLR, i_atom, shift, gid);
707 #ifndef DISABLE_WATERWATER_NLIST
708 new_i_nblist(vdwc_ww, bLR, i_atom, shift, gid);
713 new_i_nblist(vdw, bLR, i_atom, shift, gid);
717 new_i_nblist(coul, bLR, i_atom, shift, gid);
718 #ifndef DISABLE_WATERWATER_NLIST
719 new_i_nblist(coul_ww, bLR, i_atom, shift, gid);
722 /* Loop over the j charge groups */
723 for (j = 0; (j < nj); j++)
733 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
735 if (iwater == esolSPC && jwater == esolSPC)
737 /* Interaction between two SPC molecules */
740 /* VdW only - only first atoms in each water interact */
741 add_j_to_nblist(vdw, jj0, bLR);
745 #ifdef DISABLE_WATERWATER_NLIST
746 /* Add entries for the three atoms - only do VdW if we need to */
749 add_j_to_nblist(coul, jj0, bLR);
753 add_j_to_nblist(vdwc, jj0, bLR);
755 add_j_to_nblist(coul, jj0+1, bLR);
756 add_j_to_nblist(coul, jj0+2, bLR);
758 /* One entry for the entire water-water interaction */
761 add_j_to_nblist(coul_ww, jj0, bLR);
765 add_j_to_nblist(vdwc_ww, jj0, bLR);
770 else if (iwater == esolTIP4P && jwater == esolTIP4P)
772 /* Interaction between two TIP4p molecules */
775 /* VdW only - only first atoms in each water interact */
776 add_j_to_nblist(vdw, jj0, bLR);
780 #ifdef DISABLE_WATERWATER_NLIST
781 /* Add entries for the four atoms - only do VdW if we need to */
784 add_j_to_nblist(vdw, jj0, bLR);
786 add_j_to_nblist(coul, jj0+1, bLR);
787 add_j_to_nblist(coul, jj0+2, bLR);
788 add_j_to_nblist(coul, jj0+3, bLR);
790 /* One entry for the entire water-water interaction */
793 add_j_to_nblist(coul_ww, jj0, bLR);
797 add_j_to_nblist(vdwc_ww, jj0, bLR);
804 /* j charge group is not water, but i is.
805 * Add entries to the water-other_atom lists; the geometry of the water
806 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
807 * so we don't care if it is SPC or TIP4P...
814 for (jj = jj0; (jj < jj1); jj++)
818 add_j_to_nblist(coul, jj, bLR);
824 for (jj = jj0; (jj < jj1); jj++)
826 if (bHaveVdW[type[jj]])
828 add_j_to_nblist(vdw, jj, bLR);
834 /* _charge_ _groups_ interact with both coulomb and LJ */
835 /* Check which atoms we should add to the lists! */
836 for (jj = jj0; (jj < jj1); jj++)
838 if (bHaveVdW[type[jj]])
842 add_j_to_nblist(vdwc, jj, bLR);
846 add_j_to_nblist(vdw, jj, bLR);
849 else if (charge[jj] != 0)
851 add_j_to_nblist(coul, jj, bLR);
858 close_i_nblist(coul);
859 close_i_nblist(vdwc);
860 #ifndef DISABLE_WATERWATER_NLIST
861 close_i_nblist(coul_ww);
862 close_i_nblist(vdwc_ww);
867 /* no solvent as i charge group */
868 /* Loop over the atoms in the i charge group */
869 for (i = 0; i < nicg; i++)
872 gid = GID(igid, jgid, ngid);
875 /* Create new i_atom for each energy group */
876 if (bDoVdW && bDoCoul)
878 new_i_nblist(vdwc, bLR, i_atom, shift, gid);
882 new_i_nblist(vdw, bLR, i_atom, shift, gid);
886 new_i_nblist(coul, bLR, i_atom, shift, gid);
888 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
889 bDoCoul_i = (bDoCoul && qi != 0);
891 if (bDoVdW_i || bDoCoul_i)
893 /* Loop over the j charge groups */
894 for (j = 0; (j < nj); j++)
898 /* Check for large charge groups */
909 /* Finally loop over the atoms in the j-charge group */
910 for (jj = jj0; jj < jj1; jj++)
912 bNotEx = NOTEXCL(bExcl, i, jj);
920 add_j_to_nblist(coul, jj, bLR);
925 if (bHaveVdW[type[jj]])
927 add_j_to_nblist(vdw, jj, bLR);
932 if (bHaveVdW[type[jj]])
936 add_j_to_nblist(vdwc, jj, bLR);
940 add_j_to_nblist(vdw, jj, bLR);
943 else if (charge[jj] != 0)
945 add_j_to_nblist(coul, jj, bLR);
953 close_i_nblist(coul);
954 close_i_nblist(vdwc);
960 /* we are doing free energy */
961 vdwc_free = &nlist[eNL_VDWQQ_FREE];
962 vdw_free = &nlist[eNL_VDW_FREE];
963 coul_free = &nlist[eNL_QQ_FREE];
964 /* Loop over the atoms in the i charge group */
965 for (i = 0; i < nicg; i++)
968 gid = GID(igid, jgid, ngid);
970 qiB = chargeB[i_atom];
972 /* Create new i_atom for each energy group */
973 if (bDoVdW && bDoCoul)
975 new_i_nblist(vdwc, bLR, i_atom, shift, gid);
979 new_i_nblist(vdw, bLR, i_atom, shift, gid);
983 new_i_nblist(coul, bLR, i_atom, shift, gid);
986 new_i_nblist(vdw_free, bLR, i_atom, shift, gid);
987 new_i_nblist(coul_free, bLR, i_atom, shift, gid);
988 new_i_nblist(vdwc_free, bLR, i_atom, shift, gid);
990 bDoVdW_i = (bDoVdW &&
991 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
992 bDoCoul_i = (bDoCoul && (qi != 0 || qiB != 0));
993 /* For TIP4P the first atom does not have a charge,
994 * but the last three do. So we should still put an atom
995 * without LJ but with charge in the water-atom neighborlist
996 * for a TIP4p i charge group.
997 * For SPC type water the first atom has LJ and charge,
998 * so there is no such problem.
1000 if (iwater == esolNO)
1002 bDoCoul_i_sol = bDoCoul_i;
1006 bDoCoul_i_sol = bDoCoul;
1009 if (bDoVdW_i || bDoCoul_i_sol)
1011 /* Loop over the j charge groups */
1012 for (j = 0; (j < nj); j++)
1016 /* Check for large charge groups */
1027 /* Finally loop over the atoms in the j-charge group */
1028 bFree = bPert[i_atom];
1029 for (jj = jj0; (jj < jj1); jj++)
1031 bFreeJ = bFree || bPert[jj];
1032 /* Complicated if, because the water H's should also
1033 * see perturbed j-particles
1035 if (iwater == esolNO || i == 0 || bFreeJ)
1037 bNotEx = NOTEXCL(bExcl, i, jj);
1045 if (charge[jj] != 0 || chargeB[jj] != 0)
1047 add_j_to_nblist(coul_free, jj, bLR);
1050 else if (!bDoCoul_i)
1052 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1054 add_j_to_nblist(vdw_free, jj, bLR);
1059 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1061 if (charge[jj] != 0 || chargeB[jj] != 0)
1063 add_j_to_nblist(vdwc_free, jj, bLR);
1067 add_j_to_nblist(vdw_free, jj, bLR);
1070 else if (charge[jj] != 0 || chargeB[jj] != 0)
1072 add_j_to_nblist(coul_free, jj, bLR);
1078 /* This is done whether or not bWater is set */
1079 if (charge[jj] != 0)
1081 add_j_to_nblist(coul, jj, bLR);
1084 else if (!bDoCoul_i_sol)
1086 if (bHaveVdW[type[jj]])
1088 add_j_to_nblist(vdw, jj, bLR);
1093 if (bHaveVdW[type[jj]])
1095 if (charge[jj] != 0)
1097 add_j_to_nblist(vdwc, jj, bLR);
1101 add_j_to_nblist(vdw, jj, bLR);
1104 else if (charge[jj] != 0)
1106 add_j_to_nblist(coul, jj, bLR);
1114 close_i_nblist(vdw);
1115 close_i_nblist(coul);
1116 close_i_nblist(vdwc);
1117 close_i_nblist(vdw_free);
1118 close_i_nblist(coul_free);
1119 close_i_nblist(vdwc_free);
1125 put_in_list_adress(gmx_bool bHaveVdW[],
1141 /* The a[] index has been removed,
1142 * to put it back in i_atom should be a[i0] and jj should be a[jj].
1147 t_nblist * vdwc_adress = NULL;
1148 t_nblist * vdw_adress = NULL;
1149 t_nblist * coul_adress = NULL;
1150 t_nblist * vdwc_ww = NULL;
1151 t_nblist * coul_ww = NULL;
1153 int i, j, jcg, igid, gid, nbl_ind, nbl_ind_adress;
1154 atom_id jj, jj0, jj1, i_atom;
1159 real *charge, *chargeB;
1161 real qi, qiB, qq, rlj;
1162 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
1163 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
1165 gmx_bool j_all_atom;
1167 t_nblist *nlist, *nlist_adress;
1168 gmx_bool bEnergyGroupCG;
1170 /* Copy some pointers */
1171 cginfo = fr->cginfo;
1172 charge = md->chargeA;
1173 chargeB = md->chargeB;
1176 bPert = md->bPerturbed;
1179 /* Get atom range */
1181 nicg = index[icg+1]-i0;
1183 /* Get the i charge group info */
1184 igid = GET_CGINFO_GID(cginfo[icg]);
1186 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
1190 gmx_fatal(FARGS, "AdResS does not support free energy pertubation\n");
1193 /* Unpack pointers to neighbourlist structs */
1194 if (fr->nnblists == 2)
1201 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
1202 nbl_ind_adress = nbl_ind+fr->nnblists/2;
1206 nlist = fr->nblists[nbl_ind].nlist_lr;
1207 nlist_adress = fr->nblists[nbl_ind_adress].nlist_lr;
1211 nlist = fr->nblists[nbl_ind].nlist_sr;
1212 nlist_adress = fr->nblists[nbl_ind_adress].nlist_sr;
1216 vdwc = &nlist[eNL_VDWQQ];
1217 vdw = &nlist[eNL_VDW];
1218 coul = &nlist[eNL_QQ];
1220 vdwc_adress = &nlist_adress[eNL_VDWQQ];
1221 vdw_adress = &nlist_adress[eNL_VDW];
1222 coul_adress = &nlist_adress[eNL_QQ];
1224 /* We do not support solvent optimization with AdResS for now.
1225 For this we would need hybrid solvent-other kernels */
1227 /* no solvent as i charge group */
1228 /* Loop over the atoms in the i charge group */
1229 for (i = 0; i < nicg; i++)
1232 gid = GID(igid, jgid, ngid);
1233 qi = charge[i_atom];
1235 /* Create new i_atom for each energy group */
1236 if (bDoVdW && bDoCoul)
1238 new_i_nblist(vdwc, bLR, i_atom, shift, gid);
1239 new_i_nblist(vdwc_adress, bLR, i_atom, shift, gid);
1244 new_i_nblist(vdw, bLR, i_atom, shift, gid);
1245 new_i_nblist(vdw_adress, bLR, i_atom, shift, gid);
1250 new_i_nblist(coul, bLR, i_atom, shift, gid);
1251 new_i_nblist(coul_adress, bLR, i_atom, shift, gid);
1253 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
1254 bDoCoul_i = (bDoCoul && qi != 0);
1256 /* Here we find out whether the energy groups interaction belong to a
1257 * coarse-grained (vsite) or atomistic interaction. Note that, beacuse
1258 * interactions between coarse-grained and other (atomistic) energygroups
1259 * are excluded automatically by grompp, it is sufficient to check for
1260 * the group id of atom i (igid) */
1261 bEnergyGroupCG = !egp_explicit(fr, igid);
1263 if (bDoVdW_i || bDoCoul_i)
1265 /* Loop over the j charge groups */
1266 for (j = 0; (j < nj); j++)
1270 /* Check for large charge groups */
1281 /* Finally loop over the atoms in the j-charge group */
1282 for (jj = jj0; jj < jj1; jj++)
1284 bNotEx = NOTEXCL(bExcl, i, jj);
1286 /* Now we have to exclude interactions which will be zero
1287 * anyway due to the AdResS weights (in previous implementations
1288 * this was done in the force kernel). This is necessary as
1289 * pure interactions (those with b_hybrid=false, i.e. w_i*w_j==1 or 0)
1290 * are put into neighbour lists which will be passed to the
1291 * standard (optimized) kernels for speed. The interactions with
1292 * b_hybrid=true are placed into the _adress neighbour lists and
1293 * processed by the generic AdResS kernel.
1295 if ( (bEnergyGroupCG &&
1296 wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS ) ||
1297 ( !bEnergyGroupCG && wf[jj] <= GMX_REAL_EPS ) )
1302 b_hybrid = !((wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS) ||
1303 (wf[i_atom] <= GMX_REAL_EPS && wf[jj] <= GMX_REAL_EPS));
1309 if (charge[jj] != 0)
1313 add_j_to_nblist(coul, jj, bLR);
1317 add_j_to_nblist(coul_adress, jj, bLR);
1321 else if (!bDoCoul_i)
1323 if (bHaveVdW[type[jj]])
1327 add_j_to_nblist(vdw, jj, bLR);
1331 add_j_to_nblist(vdw_adress, jj, bLR);
1337 if (bHaveVdW[type[jj]])
1339 if (charge[jj] != 0)
1343 add_j_to_nblist(vdwc, jj, bLR);
1347 add_j_to_nblist(vdwc_adress, jj, bLR);
1354 add_j_to_nblist(vdw, jj, bLR);
1358 add_j_to_nblist(vdw_adress, jj, bLR);
1363 else if (charge[jj] != 0)
1367 add_j_to_nblist(coul, jj, bLR);
1371 add_j_to_nblist(coul_adress, jj, bLR);
1380 close_i_nblist(vdw);
1381 close_i_nblist(coul);
1382 close_i_nblist(vdwc);
1383 close_i_nblist(vdw_adress);
1384 close_i_nblist(coul_adress);
1385 close_i_nblist(vdwc_adress);
1391 put_in_list_qmmm(gmx_bool bHaveVdW[],
1408 int i, j, jcg, igid, gid;
1409 atom_id jj, jj0, jj1, i_atom;
1413 /* Get atom range */
1415 nicg = index[icg+1]-i0;
1417 /* Get the i charge group info */
1418 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1420 coul = &fr->QMMMlist;
1422 /* Loop over atoms in the ith charge group */
1423 for (i = 0; i < nicg; i++)
1426 gid = GID(igid, jgid, ngid);
1427 /* Create new i_atom for each energy group */
1428 new_i_nblist(coul, bLR, i_atom, shift, gid);
1430 /* Loop over the j charge groups */
1431 for (j = 0; j < nj; j++)
1435 /* Charge groups cannot have QM and MM atoms simultaneously */
1440 /* Finally loop over the atoms in the j-charge group */
1441 for (jj = jj0; jj < jj1; jj++)
1443 bNotEx = NOTEXCL(bExcl, i, jj);
1446 add_j_to_nblist(coul, jj, bLR);
1451 close_i_nblist(coul);
1456 put_in_list_cg(gmx_bool bHaveVdW[],
1473 int igid, gid, nbl_ind;
1477 cginfo = fr->cginfo[icg];
1479 igid = GET_CGINFO_GID(cginfo);
1480 gid = GID(igid, jgid, ngid);
1482 /* Unpack pointers to neighbourlist structs */
1483 if (fr->nnblists == 1)
1489 nbl_ind = fr->gid2nblists[gid];
1493 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1497 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1500 /* Make a new neighbor list for charge group icg.
1501 * Currently simply one neighbor list is made with LJ and Coulomb.
1502 * If required, zero interactions could be removed here
1503 * or in the force loop.
1505 new_i_nblist(vdwc, bLR, index[icg], shift, gid);
1506 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1508 for (j = 0; (j < nj); j++)
1511 /* Skip the icg-icg pairs if all self interactions are excluded */
1512 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1514 /* Here we add the j charge group jcg to the list,
1515 * exclusions are also added to the list.
1517 add_j_to_nblist_cg(vdwc, index[jcg], index[jcg+1], bExcl, icg == jcg, bLR);
1521 close_i_nblist(vdwc);
1524 static void setexcl(atom_id start, atom_id end, t_blocka *excl, gmx_bool b,
1531 for (i = start; i < end; i++)
1533 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1535 SETEXCL(bexcl, i-start, excl->a[k]);
1541 for (i = start; i < end; i++)
1543 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1545 RMEXCL(bexcl, i-start, excl->a[k]);
1551 int calc_naaj(int icg, int cgtot)
1555 if ((cgtot % 2) == 1)
1557 /* Odd number of charge groups, easy */
1558 naaj = 1 + (cgtot/2);
1560 else if ((cgtot % 4) == 0)
1562 /* Multiple of four is hard */
1599 fprintf(log, "naaj=%d\n", naaj);
1605 /************************************************
1607 * S I M P L E C O R E S T U F F
1609 ************************************************/
1611 static real calc_image_tric(rvec xi, rvec xj, matrix box,
1612 rvec b_inv, int *shift)
1614 /* This code assumes that the cut-off is smaller than
1615 * a half times the smallest diagonal element of the box.
1617 const real h25 = 2.5;
1622 /* Compute diff vector */
1623 dz = xj[ZZ] - xi[ZZ];
1624 dy = xj[YY] - xi[YY];
1625 dx = xj[XX] - xi[XX];
1627 /* Perform NINT operation, using trunc operation, therefore
1628 * we first add 2.5 then subtract 2 again
1630 tz = dz*b_inv[ZZ] + h25;
1632 dz -= tz*box[ZZ][ZZ];
1633 dy -= tz*box[ZZ][YY];
1634 dx -= tz*box[ZZ][XX];
1636 ty = dy*b_inv[YY] + h25;
1638 dy -= ty*box[YY][YY];
1639 dx -= ty*box[YY][XX];
1641 tx = dx*b_inv[XX]+h25;
1643 dx -= tx*box[XX][XX];
1645 /* Distance squared */
1646 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1648 *shift = XYZ2IS(tx, ty, tz);
1653 static real calc_image_rect(rvec xi, rvec xj, rvec box_size,
1654 rvec b_inv, int *shift)
1656 const real h15 = 1.5;
1662 /* Compute diff vector */
1663 dx = xj[XX] - xi[XX];
1664 dy = xj[YY] - xi[YY];
1665 dz = xj[ZZ] - xi[ZZ];
1667 /* Perform NINT operation, using trunc operation, therefore
1668 * we first add 1.5 then subtract 1 again
1670 tx = dx*b_inv[XX] + h15;
1671 ty = dy*b_inv[YY] + h15;
1672 tz = dz*b_inv[ZZ] + h15;
1677 /* Correct diff vector for translation */
1678 ddx = tx*box_size[XX] - dx;
1679 ddy = ty*box_size[YY] - dy;
1680 ddz = tz*box_size[ZZ] - dz;
1682 /* Distance squared */
1683 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1685 *shift = XYZ2IS(tx, ty, tz);
1690 static void add_simple(t_ns_buf *nsbuf, int nrj, atom_id cg_j,
1691 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1692 int icg, int jgid, t_block *cgs, t_excl bexcl[],
1693 int shift, t_forcerec *fr, put_in_list_t *put_in_list)
1695 if (nsbuf->nj + nrj > MAX_CG)
1697 put_in_list(bHaveVdW, ngid, md, icg, jgid, nsbuf->ncg, nsbuf->jcg,
1698 cgs->index, bexcl, shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1699 /* Reset buffer contents */
1700 nsbuf->ncg = nsbuf->nj = 0;
1702 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1706 static void ns_inner_tric(rvec x[], int icg, int *i_egp_flags,
1707 int njcg, atom_id jcg[],
1708 matrix box, rvec b_inv, real rcut2,
1709 t_block *cgs, t_ns_buf **ns_buf,
1710 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1711 t_excl bexcl[], t_forcerec *fr,
1712 put_in_list_t *put_in_list)
1716 int *cginfo = fr->cginfo;
1717 atom_id cg_j, *cgindex;
1720 cgindex = cgs->index;
1722 for (j = 0; (j < njcg); j++)
1725 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1726 if (calc_image_tric(x[icg], x[cg_j], box, b_inv, &shift) < rcut2)
1728 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1729 if (!(i_egp_flags[jgid] & EGP_EXCL))
1731 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1732 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1739 static void ns_inner_rect(rvec x[], int icg, int *i_egp_flags,
1740 int njcg, atom_id jcg[],
1741 gmx_bool bBox, rvec box_size, rvec b_inv, real rcut2,
1742 t_block *cgs, t_ns_buf **ns_buf,
1743 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1744 t_excl bexcl[], t_forcerec *fr,
1745 put_in_list_t *put_in_list)
1749 int *cginfo = fr->cginfo;
1750 atom_id cg_j, *cgindex;
1753 cgindex = cgs->index;
1757 for (j = 0; (j < njcg); j++)
1760 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1761 if (calc_image_rect(x[icg], x[cg_j], box_size, b_inv, &shift) < rcut2)
1763 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1764 if (!(i_egp_flags[jgid] & EGP_EXCL))
1766 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1767 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1775 for (j = 0; (j < njcg); j++)
1778 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1779 if ((rcut2 == 0) || (distance2(x[icg], x[cg_j]) < rcut2))
1781 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1782 if (!(i_egp_flags[jgid] & EGP_EXCL))
1784 add_simple(&ns_buf[jgid][CENTRAL], nrj, cg_j,
1785 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, CENTRAL, fr,
1793 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1795 static int ns_simple_core(t_forcerec *fr,
1796 gmx_localtop_t *top,
1798 matrix box, rvec box_size,
1799 t_excl bexcl[], atom_id *aaj,
1800 int ngid, t_ns_buf **ns_buf,
1801 put_in_list_t *put_in_list, gmx_bool bHaveVdW[])
1805 int nsearch, icg, jcg, igid, i0, nri, nn;
1808 /* atom_id *i_atoms; */
1809 t_block *cgs = &(top->cgs);
1810 t_blocka *excl = &(top->excls);
1813 gmx_bool bBox, bTriclinic;
1816 rlist2 = sqr(fr->rlist);
1818 bBox = (fr->ePBC != epbcNONE);
1821 for (m = 0; (m < DIM); m++)
1823 b_inv[m] = divide(1.0, box_size[m]);
1825 bTriclinic = TRICLINIC(box);
1832 cginfo = fr->cginfo;
1835 for (icg = fr->cg0; (icg < fr->hcg); icg++)
1838 i0 = cgs->index[icg];
1839 nri = cgs->index[icg+1]-i0;
1840 i_atoms = &(cgs->a[i0]);
1841 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1842 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1844 igid = GET_CGINFO_GID(cginfo[icg]);
1845 i_egp_flags = fr->egp_flags + ngid*igid;
1846 setexcl(cgs->index[icg], cgs->index[icg+1], excl, TRUE, bexcl);
1848 naaj = calc_naaj(icg, cgs->nr);
1851 ns_inner_tric(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1852 box, b_inv, rlist2, cgs, ns_buf,
1853 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1857 ns_inner_rect(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1858 bBox, box_size, b_inv, rlist2, cgs, ns_buf,
1859 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1863 for (nn = 0; (nn < ngid); nn++)
1865 for (k = 0; (k < SHIFTS); k++)
1867 nsbuf = &(ns_buf[nn][k]);
1870 put_in_list(bHaveVdW, ngid, md, icg, nn, nsbuf->ncg, nsbuf->jcg,
1871 cgs->index, bexcl, k, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1872 nsbuf->ncg = nsbuf->nj = 0;
1876 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1877 setexcl(cgs->index[icg], cgs->index[icg+1], excl, FALSE, bexcl);
1879 close_neighbor_lists(fr, FALSE);
1884 /************************************************
1886 * N S 5 G R I D S T U F F
1888 ************************************************/
1890 static inline void get_dx(int Nx, real gridx, real rc2, int xgi, real x,
1891 int *dx0, int *dx1, real *dcx2)
1919 for (i = xgi0; i >= 0; i--)
1921 dcx = (i+1)*gridx-x;
1930 for (i = xgi1; i < Nx; i++)
1943 static inline void get_dx_dd(int Nx, real gridx, real rc2, int xgi, real x,
1944 int ncpddc, int shift_min, int shift_max,
1945 int *g0, int *g1, real *dcx2)
1948 int g_min, g_max, shift_home;
1981 g_min = (shift_min == shift_home ? 0 : ncpddc);
1982 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1989 else if (shift_max < 0)
2004 /* Check one grid cell down */
2005 dcx = ((*g0 - 1) + 1)*gridx - x;
2017 /* Check one grid cell up */
2018 dcx = (*g1 + 1)*gridx - x;
2030 #define sqr(x) ((x)*(x))
2031 #define calc_dx2(XI, YI, ZI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
2032 #define calc_cyl_dx2(XI, YI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
2033 /****************************************************
2035 * F A S T N E I G H B O R S E A R C H I N G
2037 * Optimized neighboursearching routine using grid
2038 * at least 1x1x1, see GROMACS manual
2040 ****************************************************/
2043 static void get_cutoff2(t_forcerec *fr, gmx_bool bDoLongRange,
2044 real *rvdw2, real *rcoul2,
2045 real *rs2, real *rm2, real *rl2)
2047 *rs2 = sqr(fr->rlist);
2049 if (bDoLongRange && fr->bTwinRange)
2051 /* The VdW and elec. LR cut-off's could be different,
2052 * so we can not simply set them to rlistlong.
2054 if (EVDW_MIGHT_BE_ZERO_AT_CUTOFF(fr->vdwtype) &&
2055 fr->rvdw > fr->rlist)
2057 *rvdw2 = sqr(fr->rlistlong);
2061 *rvdw2 = sqr(fr->rvdw);
2063 if (EEL_MIGHT_BE_ZERO_AT_CUTOFF(fr->eeltype) &&
2064 fr->rcoulomb > fr->rlist)
2066 *rcoul2 = sqr(fr->rlistlong);
2070 *rcoul2 = sqr(fr->rcoulomb);
2075 /* Workaround for a gcc -O3 or -ffast-math problem */
2079 *rm2 = min(*rvdw2, *rcoul2);
2080 *rl2 = max(*rvdw2, *rcoul2);
2083 static void init_nsgrid_lists(t_forcerec *fr, int ngid, gmx_ns_t *ns)
2085 real rvdw2, rcoul2, rs2, rm2, rl2;
2088 get_cutoff2(fr, TRUE, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2090 /* Short range buffers */
2091 snew(ns->nl_sr, ngid);
2093 snew(ns->nsr, ngid);
2094 snew(ns->nlr_ljc, ngid);
2095 snew(ns->nlr_one, ngid);
2097 /* Always allocate both list types, since rcoulomb might now change with PME load balancing */
2098 /* Long range VdW and Coul buffers */
2099 snew(ns->nl_lr_ljc, ngid);
2100 /* Long range VdW or Coul only buffers */
2101 snew(ns->nl_lr_one, ngid);
2103 for (j = 0; (j < ngid); j++)
2105 snew(ns->nl_sr[j], MAX_CG);
2106 snew(ns->nl_lr_ljc[j], MAX_CG);
2107 snew(ns->nl_lr_one[j], MAX_CG);
2112 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
2117 static int nsgrid_core(FILE *log, t_commrec *cr, t_forcerec *fr,
2118 matrix box, rvec box_size, int ngid,
2119 gmx_localtop_t *top,
2120 t_grid *grid, rvec x[],
2121 t_excl bexcl[], gmx_bool *bExcludeAlleg,
2122 t_nrnb *nrnb, t_mdatoms *md,
2123 real *lambda, real *dvdlambda,
2124 gmx_grppairener_t *grppener,
2125 put_in_list_t *put_in_list,
2126 gmx_bool bHaveVdW[],
2127 gmx_bool bDoLongRange, gmx_bool bMakeQMMMnblist)
2130 atom_id **nl_lr_ljc, **nl_lr_one, **nl_sr;
2131 int *nlr_ljc, *nlr_one, *nsr;
2132 gmx_domdec_t *dd = NULL;
2133 t_block *cgs = &(top->cgs);
2134 int *cginfo = fr->cginfo;
2135 /* atom_id *i_atoms,*cgsindex=cgs->index; */
2137 int cell_x, cell_y, cell_z;
2138 int d, tx, ty, tz, dx, dy, dz, cj;
2139 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2140 int zsh_ty, zsh_tx, ysh_tx;
2142 int dx0, dx1, dy0, dy1, dz0, dz1;
2143 int Nx, Ny, Nz, shift = -1, j, nrj, nns, nn = -1;
2144 real gridx, gridy, gridz, grid_x, grid_y, grid_z;
2145 real *dcx2, *dcy2, *dcz2;
2147 int cg0, cg1, icg = -1, cgsnr, i0, igid, nri, naaj, max_jcg;
2148 int jcg0, jcg1, jjcg, cgj0, jgid;
2149 int *grida, *gridnra, *gridind;
2150 gmx_bool rvdw_lt_rcoul, rcoul_lt_rvdw;
2151 rvec xi, *cgcm, grid_offset;
2152 real r2, rs2, rvdw2, rcoul2, rm2, rl2, XI, YI, ZI, dcx, dcy, dcz, tmp1, tmp2;
2154 gmx_bool bDomDec, bTriclinicX, bTriclinicY;
2159 bDomDec = DOMAINDECOMP(cr);
2165 bTriclinicX = ((YY < grid->npbcdim &&
2166 (!bDomDec || dd->nc[YY] == 1) && box[YY][XX] != 0) ||
2167 (ZZ < grid->npbcdim &&
2168 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][XX] != 0));
2169 bTriclinicY = (ZZ < grid->npbcdim &&
2170 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][YY] != 0);
2174 get_cutoff2(fr, bDoLongRange, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2176 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
2177 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
2179 if (bMakeQMMMnblist)
2187 nl_lr_ljc = ns->nl_lr_ljc;
2188 nl_lr_one = ns->nl_lr_one;
2189 nlr_ljc = ns->nlr_ljc;
2190 nlr_one = ns->nlr_one;
2198 gridind = grid->index;
2199 gridnra = grid->nra;
2202 gridx = grid->cell_size[XX];
2203 gridy = grid->cell_size[YY];
2204 gridz = grid->cell_size[ZZ];
2208 copy_rvec(grid->cell_offset, grid_offset);
2209 copy_ivec(grid->ncpddc, ncpddc);
2214 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2215 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
2216 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
2217 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
2218 if (zsh_tx != 0 && ysh_tx != 0)
2220 /* This could happen due to rounding, when both ratios are 0.5 */
2229 /* We only want a list for the test particle */
2238 /* Set the shift range */
2239 for (d = 0; d < DIM; d++)
2243 /* Check if we need periodicity shifts.
2244 * Without PBC or with domain decomposition we don't need them.
2246 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2253 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
2264 /* Loop over charge groups */
2265 for (icg = cg0; (icg < cg1); icg++)
2267 igid = GET_CGINFO_GID(cginfo[icg]);
2268 /* Skip this charge group if all energy groups are excluded! */
2269 if (bExcludeAlleg[igid])
2274 i0 = cgs->index[icg];
2276 if (bMakeQMMMnblist)
2278 /* Skip this charge group if it is not a QM atom while making a
2279 * QM/MM neighbourlist
2281 if (md->bQM[i0] == FALSE)
2283 continue; /* MM particle, go to next particle */
2286 /* Compute the number of charge groups that fall within the control
2289 naaj = calc_naaj(icg, cgsnr);
2296 /* make a normal neighbourlist */
2300 /* Get the j charge-group and dd cell shift ranges */
2301 dd_get_ns_ranges(cr->dd, icg, &jcg0, &jcg1, sh0, sh1);
2306 /* Compute the number of charge groups that fall within the control
2309 naaj = calc_naaj(icg, cgsnr);
2315 /* The i-particle is awlways the test particle,
2316 * so we want all j-particles
2318 max_jcg = cgsnr - 1;
2322 max_jcg = jcg1 - cgsnr;
2327 i_egp_flags = fr->egp_flags + igid*ngid;
2329 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2330 setexcl(i0, cgs->index[icg+1], &top->excls, TRUE, bexcl);
2332 ci2xyz(grid, icg, &cell_x, &cell_y, &cell_z);
2334 /* Changed iicg to icg, DvdS 990115
2335 * (but see consistency check above, DvdS 990330)
2338 fprintf(log, "icg=%5d, naaj=%5d, cell %d %d %d\n",
2339 icg, naaj, cell_x, cell_y, cell_z);
2341 /* Loop over shift vectors in three dimensions */
2342 for (tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
2344 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2345 /* Calculate range of cells in Z direction that have the shift tz */
2346 zgi = cell_z + tz*Nz;
2349 get_dx(Nz, gridz, rl2, zgi, ZI, &dz0, &dz1, dcz2);
2351 get_dx_dd(Nz, gridz, rl2, zgi, ZI-grid_offset[ZZ],
2352 ncpddc[ZZ], sh0[ZZ], sh1[ZZ], &dz0, &dz1, dcz2);
2358 for (ty = -shp[YY]; ty <= shp[YY]; ty++)
2360 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2361 /* Calculate range of cells in Y direction that have the shift ty */
2364 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2368 ygi = cell_y + ty*Ny;
2371 get_dx(Ny, gridy, rl2, ygi, YI, &dy0, &dy1, dcy2);
2373 get_dx_dd(Ny, gridy, rl2, ygi, YI-grid_offset[YY],
2374 ncpddc[YY], sh0[YY], sh1[YY], &dy0, &dy1, dcy2);
2380 for (tx = -shp[XX]; tx <= shp[XX]; tx++)
2382 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2383 /* Calculate range of cells in X direction that have the shift tx */
2386 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2390 xgi = cell_x + tx*Nx;
2393 get_dx(Nx, gridx, rl2, xgi*Nx, XI, &dx0, &dx1, dcx2);
2395 get_dx_dd(Nx, gridx, rl2, xgi, XI-grid_offset[XX],
2396 ncpddc[XX], sh0[XX], sh1[XX], &dx0, &dx1, dcx2);
2402 /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2403 * from the neigbour list as it will not interact */
2404 if (fr->adress_type != eAdressOff)
2406 if (md->wf[cgs->index[icg]] <= GMX_REAL_EPS && egp_explicit(fr, igid))
2411 /* Get shift vector */
2412 shift = XYZ2IS(tx, ty, tz);
2414 range_check(shift, 0, SHIFTS);
2416 for (nn = 0; (nn < ngid); nn++)
2423 fprintf(log, "shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2424 shift, dx0, dx1, dy0, dy1, dz0, dz1);
2425 fprintf(log, "cgcm: %8.3f %8.3f %8.3f\n", cgcm[icg][XX],
2426 cgcm[icg][YY], cgcm[icg][ZZ]);
2427 fprintf(log, "xi: %8.3f %8.3f %8.3f\n", XI, YI, ZI);
2429 for (dx = dx0; (dx <= dx1); dx++)
2431 tmp1 = rl2 - dcx2[dx];
2432 for (dy = dy0; (dy <= dy1); dy++)
2434 tmp2 = tmp1 - dcy2[dy];
2437 for (dz = dz0; (dz <= dz1); dz++)
2439 if (tmp2 > dcz2[dz])
2441 /* Find grid-cell cj in which possible neighbours are */
2442 cj = xyz2ci(Ny, Nz, dx, dy, dz);
2444 /* Check out how many cgs (nrj) there in this cell */
2447 /* Find the offset in the cg list */
2450 /* Check if all j's are out of range so we
2451 * can skip the whole cell.
2452 * Should save some time, especially with DD.
2455 (grida[cgj0] >= max_jcg &&
2456 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2462 for (j = 0; (j < nrj); j++)
2464 jjcg = grida[cgj0+j];
2466 /* check whether this guy is in range! */
2467 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2470 r2 = calc_dx2(XI, YI, ZI, cgcm[jjcg]);
2473 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2474 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2475 /* check energy group exclusions */
2476 if (!(i_egp_flags[jgid] & EGP_EXCL))
2480 if (nsr[jgid] >= MAX_CG)
2482 /* Add to short-range list */
2483 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2484 nsr[jgid], nl_sr[jgid],
2485 cgs->index, /* cgsatoms, */ bexcl,
2486 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2489 nl_sr[jgid][nsr[jgid]++] = jjcg;
2493 if (nlr_ljc[jgid] >= MAX_CG)
2495 /* Add to LJ+coulomb long-range list */
2496 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2497 nlr_ljc[jgid], nl_lr_ljc[jgid], top->cgs.index,
2498 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2501 nl_lr_ljc[jgid][nlr_ljc[jgid]++] = jjcg;
2505 if (nlr_one[jgid] >= MAX_CG)
2507 /* Add to long-range list with only coul, or only LJ */
2508 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2509 nlr_one[jgid], nl_lr_one[jgid], top->cgs.index,
2510 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2513 nl_lr_one[jgid][nlr_one[jgid]++] = jjcg;
2525 /* CHECK whether there is anything left in the buffers */
2526 for (nn = 0; (nn < ngid); nn++)
2530 put_in_list(bHaveVdW, ngid, md, icg, nn, nsr[nn], nl_sr[nn],
2531 cgs->index, /* cgsatoms, */ bexcl,
2532 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2535 if (nlr_ljc[nn] > 0)
2537 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_ljc[nn],
2538 nl_lr_ljc[nn], top->cgs.index,
2539 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2542 if (nlr_one[nn] > 0)
2544 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_one[nn],
2545 nl_lr_one[nn], top->cgs.index,
2546 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2552 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2553 setexcl(cgs->index[icg], cgs->index[icg+1], &top->excls, FALSE, bexcl);
2555 /* No need to perform any left-over force calculations anymore (as we used to do here)
2556 * since we now save the proper long-range lists for later evaluation.
2561 /* Close neighbourlists */
2562 close_neighbor_lists(fr, bMakeQMMMnblist);
2567 void ns_realloc_natoms(gmx_ns_t *ns, int natoms)
2571 if (natoms > ns->nra_alloc)
2573 ns->nra_alloc = over_alloc_dd(natoms);
2574 srenew(ns->bexcl, ns->nra_alloc);
2575 for (i = 0; i < ns->nra_alloc; i++)
2582 void init_ns(FILE *fplog, const t_commrec *cr,
2583 gmx_ns_t *ns, t_forcerec *fr,
2584 const gmx_mtop_t *mtop,
2587 int mt, icg, nr_in_cg, maxcg, i, j, jcg, ngid, ncg;
2591 /* Compute largest charge groups size (# atoms) */
2593 for (mt = 0; mt < mtop->nmoltype; mt++)
2595 cgs = &mtop->moltype[mt].cgs;
2596 for (icg = 0; (icg < cgs->nr); icg++)
2598 nr_in_cg = max(nr_in_cg, (int)(cgs->index[icg+1]-cgs->index[icg]));
2602 /* Verify whether largest charge group is <= max cg.
2603 * This is determined by the type of the local exclusion type
2604 * Exclusions are stored in bits. (If the type is not large
2605 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2607 maxcg = sizeof(t_excl)*8;
2608 if (nr_in_cg > maxcg)
2610 gmx_fatal(FARGS, "Max #atoms in a charge group: %d > %d\n",
2614 ngid = mtop->groups.grps[egcENER].nr;
2615 snew(ns->bExcludeAlleg, ngid);
2616 for (i = 0; i < ngid; i++)
2618 ns->bExcludeAlleg[i] = TRUE;
2619 for (j = 0; j < ngid; j++)
2621 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2623 ns->bExcludeAlleg[i] = FALSE;
2631 ns->grid = init_grid(fplog, fr);
2632 init_nsgrid_lists(fr, ngid, ns);
2637 snew(ns->ns_buf, ngid);
2638 for (i = 0; (i < ngid); i++)
2640 snew(ns->ns_buf[i], SHIFTS);
2642 ncg = ncg_mtop(mtop);
2643 snew(ns->simple_aaj, 2*ncg);
2644 for (jcg = 0; (jcg < ncg); jcg++)
2646 ns->simple_aaj[jcg] = jcg;
2647 ns->simple_aaj[jcg+ncg] = jcg;
2651 /* Create array that determines whether or not atoms have VdW */
2652 snew(ns->bHaveVdW, fr->ntype);
2653 for (i = 0; (i < fr->ntype); i++)
2655 for (j = 0; (j < fr->ntype); j++)
2657 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2659 ((BHAMA(fr->nbfp, fr->ntype, i, j) != 0) ||
2660 (BHAMB(fr->nbfp, fr->ntype, i, j) != 0) ||
2661 (BHAMC(fr->nbfp, fr->ntype, i, j) != 0)) :
2662 ((C6(fr->nbfp, fr->ntype, i, j) != 0) ||
2663 (C12(fr->nbfp, fr->ntype, i, j) != 0))));
2668 pr_bvec(debug, 0, "bHaveVdW", ns->bHaveVdW, fr->ntype, TRUE);
2673 if (!DOMAINDECOMP(cr))
2675 /* This could be reduced with particle decomposition */
2676 ns_realloc_natoms(ns, mtop->natoms);
2679 ns->nblist_initialized = FALSE;
2681 /* nbr list debug dump */
2683 char *ptr = getenv("GMX_DUMP_NL");
2686 ns->dump_nl = strtol(ptr, NULL, 10);
2689 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2700 int search_neighbours(FILE *log, t_forcerec *fr,
2701 rvec x[], matrix box,
2702 gmx_localtop_t *top,
2703 gmx_groups_t *groups,
2705 t_nrnb *nrnb, t_mdatoms *md,
2706 real *lambda, real *dvdlambda,
2707 gmx_grppairener_t *grppener,
2709 gmx_bool bDoLongRangeNS,
2710 gmx_bool bPadListsForKernels)
2712 t_block *cgs = &(top->cgs);
2713 rvec box_size, grid_x0, grid_x1;
2715 real min_size, grid_dens;
2719 gmx_bool *i_egp_flags;
2720 int cg_start, cg_end, start, end;
2723 gmx_domdec_zones_t *dd_zones;
2724 put_in_list_t *put_in_list;
2728 /* Set some local variables */
2730 ngid = groups->grps[egcENER].nr;
2732 for (m = 0; (m < DIM); m++)
2734 box_size[m] = box[m][m];
2737 if (fr->ePBC != epbcNONE)
2739 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC, box))
2741 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.");
2745 min_size = min(box_size[XX], min(box_size[YY], box_size[ZZ]));
2746 if (2*fr->rlistlong >= min_size)
2748 gmx_fatal(FARGS, "One of the box diagonal elements has become smaller than twice the cut-off length.");
2753 if (DOMAINDECOMP(cr))
2755 ns_realloc_natoms(ns, cgs->index[cgs->nr]);
2759 /* Reset the neighbourlists */
2760 reset_neighbor_lists(fr, TRUE, TRUE);
2762 if (bGrid && bFillGrid)
2766 if (DOMAINDECOMP(cr))
2768 dd_zones = domdec_zones(cr->dd);
2774 get_nsgrid_boundaries(grid->nboundeddim, box, NULL, NULL, NULL, NULL,
2775 cgs->nr, fr->cg_cm, grid_x0, grid_x1, &grid_dens);
2777 grid_first(log, grid, NULL, NULL, fr->ePBC, box, grid_x0, grid_x1,
2778 fr->rlistlong, grid_dens);
2782 /* Don't know why this all is... (DvdS 3/99) */
2788 end = (cgs->nr+1)/2;
2791 if (DOMAINDECOMP(cr))
2794 fill_grid(log, dd_zones, grid, end, -1, end, fr->cg_cm);
2796 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2800 fill_grid(log, NULL, grid, cgs->nr, fr->cg0, fr->hcg, fr->cg_cm);
2801 grid->icg0 = fr->cg0;
2802 grid->icg1 = fr->hcg;
2812 calc_elemnr(log, grid, start, end, cgs->nr);
2814 grid_last(log, grid, start, end, cgs->nr);
2818 check_grid(debug, grid);
2819 print_grid(debug, grid);
2824 /* Set the grid cell index for the test particle only.
2825 * The cell to cg index is not corrected, but that does not matter.
2827 fill_grid(log, NULL, ns->grid, fr->hcg, fr->hcg-1, fr->hcg, fr->cg_cm);
2831 if (fr->adress_type == eAdressOff)
2833 if (!fr->ns.bCGlist)
2835 put_in_list = put_in_list_at;
2839 put_in_list = put_in_list_cg;
2844 put_in_list = put_in_list_adress;
2851 nsearch = nsgrid_core(log, cr, fr, box, box_size, ngid, top,
2852 grid, x, ns->bexcl, ns->bExcludeAlleg,
2853 nrnb, md, lambda, dvdlambda, grppener,
2854 put_in_list, ns->bHaveVdW,
2855 bDoLongRangeNS, FALSE);
2857 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2858 * the classical calculation. The charge-charge interaction
2859 * between QM and MM atoms is handled in the QMMM core calculation
2860 * (see QMMM.c). The VDW however, we'd like to compute classically
2861 * and the QM MM atom pairs have just been put in the
2862 * corresponding neighbourlists. in case of QMMM we still need to
2863 * fill a special QMMM neighbourlist that contains all neighbours
2864 * of the QM atoms. If bQMMM is true, this list will now be made:
2866 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
2868 nsearch += nsgrid_core(log, cr, fr, box, box_size, ngid, top,
2869 grid, x, ns->bexcl, ns->bExcludeAlleg,
2870 nrnb, md, lambda, dvdlambda, grppener,
2871 put_in_list_qmmm, ns->bHaveVdW,
2872 bDoLongRangeNS, TRUE);
2877 nsearch = ns_simple_core(fr, top, md, box, box_size,
2878 ns->bexcl, ns->simple_aaj,
2879 ngid, ns->ns_buf, put_in_list, ns->bHaveVdW);
2887 inc_nrnb(nrnb, eNR_NS, nsearch);
2888 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2893 int natoms_beyond_ns_buffer(t_inputrec *ir, t_forcerec *fr, t_block *cgs,
2894 matrix scale_tot, rvec *x)
2896 int cg0, cg1, cg, a0, a1, a, i, j;
2897 real rint, hbuf2, scale;
2899 gmx_bool bIsotropic;
2904 rint = max(ir->rcoulomb, ir->rvdw);
2905 if (ir->rlist < rint)
2907 gmx_fatal(FARGS, "The neighbor search buffer has negative size: %f nm",
2915 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2917 hbuf2 = sqr(0.5*(ir->rlist - rint));
2918 for (cg = cg0; cg < cg1; cg++)
2920 a0 = cgs->index[cg];
2921 a1 = cgs->index[cg+1];
2922 for (a = a0; a < a1; a++)
2924 if (distance2(cg_cm[cg], x[a]) > hbuf2)
2934 scale = scale_tot[0][0];
2935 for (i = 1; i < DIM; i++)
2937 /* With anisotropic scaling, the original spherical ns volumes become
2938 * ellipsoids. To avoid costly transformations we use the minimum
2939 * eigenvalue of the scaling matrix for determining the buffer size.
2940 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2942 scale = min(scale, scale_tot[i][i]);
2943 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2947 for (j = 0; j < i; j++)
2949 if (scale_tot[i][j] != 0)
2955 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2958 for (cg = cg0; cg < cg1; cg++)
2960 svmul(scale, cg_cm[cg], cgsc);
2961 a0 = cgs->index[cg];
2962 a1 = cgs->index[cg+1];
2963 for (a = a0; a < a1; a++)
2965 if (distance2(cgsc, x[a]) > hbuf2)
2974 /* Anistropic scaling */
2975 for (cg = cg0; cg < cg1; cg++)
2977 /* Since scale_tot contains the transpose of the scaling matrix,
2978 * we need to multiply with the transpose.
2980 tmvmul_ur0(scale_tot, cg_cm[cg], cgsc);
2981 a0 = cgs->index[cg];
2982 a1 = cgs->index[cg+1];
2983 for (a = a0; a < a1; a++)
2985 if (distance2(cgsc, x[a]) > hbuf2)