2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "gromacs/legacyheaders/macros.h"
46 #include "gromacs/math/utilities.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/legacyheaders/types/commrec.h"
49 #include "gromacs/legacyheaders/network.h"
50 #include "gromacs/legacyheaders/nsgrid.h"
51 #include "gromacs/legacyheaders/force.h"
52 #include "gromacs/legacyheaders/nonbonded.h"
53 #include "gromacs/legacyheaders/ns.h"
54 #include "gromacs/legacyheaders/names.h"
55 #include "gromacs/legacyheaders/nrnb.h"
56 #include "gromacs/legacyheaders/txtdump.h"
57 #include "gromacs/topology/mtop_util.h"
59 #include "gromacs/legacyheaders/domdec.h"
62 #include "gromacs/pbcutil/ishift.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/utility/fatalerror.h"
65 #include "gromacs/utility/smalloc.h"
68 * E X C L U S I O N H A N D L I N G
72 static void SETEXCL_(t_excl e[], atom_id i, atom_id j)
76 static void RMEXCL_(t_excl e[], atom_id i, atom_id j)
78 e[j] = e[j] & ~(1<<i);
80 static gmx_bool ISEXCL_(t_excl e[], atom_id i, atom_id j)
82 return (gmx_bool)(e[j] & (1<<i));
84 static gmx_bool NOTEXCL_(t_excl e[], atom_id i, atom_id j)
86 return !(ISEXCL(e, i, j));
89 #define SETEXCL(e, i, j) (e)[((atom_id) (j))] |= (1<<((atom_id) (i)))
90 #define RMEXCL(e, i, j) (e)[((atom_id) (j))] &= (~(1<<((atom_id) (i))))
91 #define ISEXCL(e, i, j) (gmx_bool) ((e)[((atom_id) (j))] & (1<<((atom_id) (i))))
92 #define NOTEXCL(e, i, j) !(ISEXCL(e, i, j))
96 round_up_to_simd_width(int length, int simd_width)
98 int offset, newlength;
100 offset = (simd_width > 0) ? length % simd_width : 0;
102 return (offset == 0) ? length : length-offset+simd_width;
104 /************************************************
106 * U T I L I T I E S F O R N S
108 ************************************************/
110 void reallocate_nblist(t_nblist *nl)
114 fprintf(debug, "reallocating neigborlist (ielec=%d, ivdw=%d, igeometry=%d, type=%d), maxnri=%d\n",
115 nl->ielec, nl->ivdw, nl->igeometry, nl->type, nl->maxnri);
117 srenew(nl->iinr, nl->maxnri);
118 if (nl->igeometry == GMX_NBLIST_GEOMETRY_CG_CG)
120 srenew(nl->iinr_end, nl->maxnri);
122 srenew(nl->gid, nl->maxnri);
123 srenew(nl->shift, nl->maxnri);
124 srenew(nl->jindex, nl->maxnri+1);
128 static void init_nblist(FILE *log, t_nblist *nl_sr, t_nblist *nl_lr,
129 int maxsr, int maxlr,
130 int ivdw, int ivdwmod,
131 int ielec, int ielecmod,
132 int igeometry, int type,
133 gmx_bool bElecAndVdwSwitchDiffers)
139 for (i = 0; (i < 2); i++)
141 nl = (i == 0) ? nl_sr : nl_lr;
142 homenr = (i == 0) ? maxsr : maxlr;
150 /* Set coul/vdw in neighborlist, and for the normal loops we determine
151 * an index of which one to call.
154 nl->ivdwmod = ivdwmod;
156 nl->ielecmod = ielecmod;
158 nl->igeometry = igeometry;
160 if (nl->type == GMX_NBLIST_INTERACTION_FREE_ENERGY)
162 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
165 /* This will also set the simd_padding_width field */
166 gmx_nonbonded_set_kernel_pointers( (i == 0) ? log : NULL, nl, bElecAndVdwSwitchDiffers);
168 /* maxnri is influenced by the number of shifts (maximum is 8)
169 * and the number of energy groups.
170 * If it is not enough, nl memory will be reallocated during the run.
171 * 4 seems to be a reasonable factor, which only causes reallocation
172 * during runs with tiny and many energygroups.
174 nl->maxnri = homenr*4;
184 reallocate_nblist(nl);
189 fprintf(debug, "Initiating neighbourlist (ielec=%d, ivdw=%d, type=%d) for %s interactions,\nwith %d SR, %d LR atoms.\n",
190 nl->ielec, nl->ivdw, nl->type, gmx_nblist_geometry_names[nl->igeometry], maxsr, maxlr);
195 void init_neighbor_list(FILE *log, t_forcerec *fr, int homenr)
197 /* Make maxlr tunable! (does not seem to be a big difference though)
198 * This parameter determines the number of i particles in a long range
199 * neighbourlist. Too few means many function calls, too many means
202 int maxsr, maxsr_wat, maxlr, maxlr_wat;
203 int ielec, ivdw, ielecmod, ivdwmod, type;
205 int igeometry_def, igeometry_w, igeometry_ww;
207 gmx_bool bElecAndVdwSwitchDiffers;
210 /* maxsr = homenr-fr->nWatMol*3; */
215 gmx_fatal(FARGS, "%s, %d: Negative number of short range atoms.\n"
216 "Call your Gromacs dealer for assistance.", __FILE__, __LINE__);
218 /* This is just for initial allocation, so we do not reallocate
219 * all the nlist arrays many times in a row.
220 * The numbers seem very accurate, but they are uncritical.
222 maxsr_wat = min(fr->nWatMol, (homenr+2)/3);
226 maxlr_wat = min(maxsr_wat, maxlr);
230 maxlr = maxlr_wat = 0;
233 /* Determine the values for ielec/ivdw. */
234 ielec = fr->nbkernel_elec_interaction;
235 ivdw = fr->nbkernel_vdw_interaction;
236 ielecmod = fr->nbkernel_elec_modifier;
237 ivdwmod = fr->nbkernel_vdw_modifier;
238 type = GMX_NBLIST_INTERACTION_STANDARD;
239 bElecAndVdwSwitchDiffers = ( (fr->rcoulomb_switch != fr->rvdw_switch) || (fr->rcoulomb != fr->rvdw));
241 fr->ns.bCGlist = (getenv("GMX_NBLISTCG") != 0);
244 igeometry_def = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
248 igeometry_def = GMX_NBLIST_GEOMETRY_CG_CG;
251 fprintf(log, "\nUsing charge-group - charge-group neighbor lists and kernels\n\n");
255 if (fr->solvent_opt == esolTIP4P)
257 igeometry_w = GMX_NBLIST_GEOMETRY_WATER4_PARTICLE;
258 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER4_WATER4;
262 igeometry_w = GMX_NBLIST_GEOMETRY_WATER3_PARTICLE;
263 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER3_WATER3;
266 for (i = 0; i < fr->nnblists; i++)
268 nbl = &(fr->nblists[i]);
270 if ((fr->adress_type != eAdressOff) && (i >= fr->nnblists/2))
272 type = GMX_NBLIST_INTERACTION_ADRESS;
274 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ], &nbl->nlist_lr[eNL_VDWQQ],
275 maxsr, maxlr, ivdw, ivdwmod, ielec, ielecmod, igeometry_def, type, bElecAndVdwSwitchDiffers);
276 init_nblist(log, &nbl->nlist_sr[eNL_VDW], &nbl->nlist_lr[eNL_VDW],
277 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, igeometry_def, type, bElecAndVdwSwitchDiffers);
278 init_nblist(log, &nbl->nlist_sr[eNL_QQ], &nbl->nlist_lr[eNL_QQ],
279 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_def, type, bElecAndVdwSwitchDiffers);
280 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATER], &nbl->nlist_lr[eNL_VDWQQ_WATER],
281 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_w, type, bElecAndVdwSwitchDiffers);
282 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATER], &nbl->nlist_lr[eNL_QQ_WATER],
283 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_w, type, bElecAndVdwSwitchDiffers);
284 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATERWATER], &nbl->nlist_lr[eNL_VDWQQ_WATERWATER],
285 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_ww, type, bElecAndVdwSwitchDiffers);
286 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATERWATER], &nbl->nlist_lr[eNL_QQ_WATERWATER],
287 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_ww, type, bElecAndVdwSwitchDiffers);
289 /* Did we get the solvent loops so we can use optimized water kernels? */
290 if (nbl->nlist_sr[eNL_VDWQQ_WATER].kernelptr_vf == NULL
291 || nbl->nlist_sr[eNL_QQ_WATER].kernelptr_vf == NULL
292 || nbl->nlist_sr[eNL_VDWQQ_WATERWATER].kernelptr_vf == NULL
293 || nbl->nlist_sr[eNL_QQ_WATERWATER].kernelptr_vf == NULL)
295 fr->solvent_opt = esolNO;
298 fprintf(log, "Note: The available nonbonded kernels do not support water optimization - disabling.\n");
302 if (fr->efep != efepNO)
304 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_FREE], &nbl->nlist_lr[eNL_VDWQQ_FREE],
305 maxsr, maxlr, ivdw, ivdwmod, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY, bElecAndVdwSwitchDiffers);
306 init_nblist(log, &nbl->nlist_sr[eNL_VDW_FREE], &nbl->nlist_lr[eNL_VDW_FREE],
307 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY, bElecAndVdwSwitchDiffers);
308 init_nblist(log, &nbl->nlist_sr[eNL_QQ_FREE], &nbl->nlist_lr[eNL_QQ_FREE],
309 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY, bElecAndVdwSwitchDiffers);
313 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
315 init_nblist(log, &fr->QMMMlist, NULL,
316 maxsr, maxlr, 0, 0, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_STANDARD, bElecAndVdwSwitchDiffers);
324 fr->ns.nblist_initialized = TRUE;
327 static void reset_nblist(t_nblist *nl)
337 static void reset_neighbor_lists(t_forcerec *fr, gmx_bool bResetSR, gmx_bool bResetLR)
343 /* only reset the short-range nblist */
344 reset_nblist(&(fr->QMMMlist));
347 for (n = 0; n < fr->nnblists; n++)
349 for (i = 0; i < eNL_NR; i++)
353 reset_nblist( &(fr->nblists[n].nlist_sr[i]) );
357 reset_nblist( &(fr->nblists[n].nlist_lr[i]) );
366 static gmx_inline void new_i_nblist(t_nblist *nlist, atom_id i_atom, int shift, int gid)
368 int i, k, nri, nshift;
372 /* Check whether we have to increase the i counter */
374 (nlist->iinr[nri] != i_atom) ||
375 (nlist->shift[nri] != shift) ||
376 (nlist->gid[nri] != gid))
378 /* This is something else. Now see if any entries have
379 * been added in the list of the previous atom.
382 ((nlist->jindex[nri+1] > nlist->jindex[nri]) &&
383 (nlist->gid[nri] != -1)))
385 /* If so increase the counter */
388 if (nlist->nri >= nlist->maxnri)
390 nlist->maxnri += over_alloc_large(nlist->nri);
391 reallocate_nblist(nlist);
394 /* Set the number of neighbours and the atom number */
395 nlist->jindex[nri+1] = nlist->jindex[nri];
396 nlist->iinr[nri] = i_atom;
397 nlist->gid[nri] = gid;
398 nlist->shift[nri] = shift;
402 /* Adding to previous list. First remove possible previous padding */
403 if (nlist->simd_padding_width > 1)
405 while (nlist->nrj > 0 && nlist->jjnr[nlist->nrj-1] < 0)
413 static gmx_inline void close_i_nblist(t_nblist *nlist)
415 int nri = nlist->nri;
420 /* Add elements up to padding. Since we allocate memory in units
421 * of the simd_padding width, we do not have to check for possible
422 * list reallocation here.
424 while ((nlist->nrj % nlist->simd_padding_width) != 0)
426 /* Use -4 here, so we can write forces for 4 atoms before real data */
427 nlist->jjnr[nlist->nrj++] = -4;
429 nlist->jindex[nri+1] = nlist->nrj;
431 len = nlist->nrj - nlist->jindex[nri];
435 static gmx_inline void close_nblist(t_nblist *nlist)
437 /* Only close this nblist when it has been initialized.
438 * Avoid the creation of i-lists with no j-particles.
442 /* Some assembly kernels do not support empty lists,
443 * make sure here that we don't generate any empty lists.
444 * With the current ns code this branch is taken in two cases:
445 * No i-particles at all: nri=-1 here
446 * There are i-particles, but no j-particles; nri=0 here
452 /* Close list number nri by incrementing the count */
457 static gmx_inline void close_neighbor_lists(t_forcerec *fr, gmx_bool bMakeQMMMnblist)
463 close_nblist(&(fr->QMMMlist));
466 for (n = 0; n < fr->nnblists; n++)
468 for (i = 0; (i < eNL_NR); i++)
470 close_nblist(&(fr->nblists[n].nlist_sr[i]));
471 close_nblist(&(fr->nblists[n].nlist_lr[i]));
477 static gmx_inline void add_j_to_nblist(t_nblist *nlist, atom_id j_atom, gmx_bool bLR)
479 int nrj = nlist->nrj;
481 if (nlist->nrj >= nlist->maxnrj)
483 nlist->maxnrj = round_up_to_simd_width(over_alloc_small(nlist->nrj + 1), nlist->simd_padding_width);
487 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
488 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
491 srenew(nlist->jjnr, nlist->maxnrj);
494 nlist->jjnr[nrj] = j_atom;
498 static gmx_inline void add_j_to_nblist_cg(t_nblist *nlist,
499 atom_id j_start, int j_end,
500 t_excl *bexcl, gmx_bool i_is_j,
503 int nrj = nlist->nrj;
506 if (nlist->nrj >= nlist->maxnrj)
508 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
511 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
512 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
515 srenew(nlist->jjnr, nlist->maxnrj);
516 srenew(nlist->jjnr_end, nlist->maxnrj);
517 srenew(nlist->excl, nlist->maxnrj*MAX_CGCGSIZE);
520 nlist->jjnr[nrj] = j_start;
521 nlist->jjnr_end[nrj] = j_end;
523 if (j_end - j_start > MAX_CGCGSIZE)
525 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);
528 /* Set the exclusions */
529 for (j = j_start; j < j_end; j++)
531 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
535 /* Avoid double counting of intra-cg interactions */
536 for (j = 1; j < j_end-j_start; j++)
538 nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
546 put_in_list_t (gmx_bool bHaveVdW[],
563 put_in_list_at(gmx_bool bHaveVdW[],
579 /* The a[] index has been removed,
580 * to put it back in i_atom should be a[i0] and jj should be a[jj].
585 t_nblist * vdwc_free = NULL;
586 t_nblist * vdw_free = NULL;
587 t_nblist * coul_free = NULL;
588 t_nblist * vdwc_ww = NULL;
589 t_nblist * coul_ww = NULL;
591 int i, j, jcg, igid, gid, nbl_ind, ind_ij;
592 atom_id jj, jj0, jj1, i_atom;
597 real *charge, *chargeB;
598 real qi, qiB, qq, rlj;
599 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
600 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
604 /* Copy some pointers */
606 charge = md->chargeA;
607 chargeB = md->chargeB;
610 bPert = md->bPerturbed;
614 nicg = index[icg+1]-i0;
616 /* Get the i charge group info */
617 igid = GET_CGINFO_GID(cginfo[icg]);
619 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
624 /* Check if any of the particles involved are perturbed.
625 * If not we can do the cheaper normal put_in_list
626 * and use more solvent optimization.
628 for (i = 0; i < nicg; i++)
630 bFreeEnergy |= bPert[i0+i];
632 /* Loop over the j charge groups */
633 for (j = 0; (j < nj && !bFreeEnergy); j++)
638 /* Finally loop over the atoms in the j-charge group */
639 for (jj = jj0; jj < jj1; jj++)
641 bFreeEnergy |= bPert[jj];
646 /* Unpack pointers to neighbourlist structs */
647 if (fr->nnblists == 1)
653 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
657 nlist = fr->nblists[nbl_ind].nlist_lr;
661 nlist = fr->nblists[nbl_ind].nlist_sr;
664 if (iwater != esolNO)
666 vdwc = &nlist[eNL_VDWQQ_WATER];
667 vdw = &nlist[eNL_VDW];
668 coul = &nlist[eNL_QQ_WATER];
669 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
670 coul_ww = &nlist[eNL_QQ_WATERWATER];
674 vdwc = &nlist[eNL_VDWQQ];
675 vdw = &nlist[eNL_VDW];
676 coul = &nlist[eNL_QQ];
681 if (iwater != esolNO)
683 /* Loop over the atoms in the i charge group */
685 gid = GID(igid, jgid, ngid);
686 /* Create new i_atom for each energy group */
687 if (bDoCoul && bDoVdW)
689 new_i_nblist(vdwc, i_atom, shift, gid);
690 new_i_nblist(vdwc_ww, i_atom, shift, gid);
694 new_i_nblist(vdw, i_atom, shift, gid);
698 new_i_nblist(coul, i_atom, shift, gid);
699 new_i_nblist(coul_ww, i_atom, shift, gid);
701 /* Loop over the j charge groups */
702 for (j = 0; (j < nj); j++)
712 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
714 if (iwater == esolSPC && jwater == esolSPC)
716 /* Interaction between two SPC molecules */
719 /* VdW only - only first atoms in each water interact */
720 add_j_to_nblist(vdw, jj0, bLR);
724 /* One entry for the entire water-water interaction */
727 add_j_to_nblist(coul_ww, jj0, bLR);
731 add_j_to_nblist(vdwc_ww, jj0, bLR);
735 else if (iwater == esolTIP4P && jwater == esolTIP4P)
737 /* Interaction between two TIP4p molecules */
740 /* VdW only - only first atoms in each water interact */
741 add_j_to_nblist(vdw, jj0, bLR);
745 /* One entry for the entire water-water interaction */
748 add_j_to_nblist(coul_ww, jj0, bLR);
752 add_j_to_nblist(vdwc_ww, jj0, bLR);
758 /* j charge group is not water, but i is.
759 * Add entries to the water-other_atom lists; the geometry of the water
760 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
761 * so we don't care if it is SPC or TIP4P...
768 for (jj = jj0; (jj < jj1); jj++)
772 add_j_to_nblist(coul, jj, bLR);
778 for (jj = jj0; (jj < jj1); jj++)
780 if (bHaveVdW[type[jj]])
782 add_j_to_nblist(vdw, jj, bLR);
788 /* _charge_ _groups_ interact with both coulomb and LJ */
789 /* Check which atoms we should add to the lists! */
790 for (jj = jj0; (jj < jj1); jj++)
792 if (bHaveVdW[type[jj]])
796 add_j_to_nblist(vdwc, jj, bLR);
800 add_j_to_nblist(vdw, jj, bLR);
803 else if (charge[jj] != 0)
805 add_j_to_nblist(coul, jj, bLR);
812 close_i_nblist(coul);
813 close_i_nblist(vdwc);
814 close_i_nblist(coul_ww);
815 close_i_nblist(vdwc_ww);
819 /* no solvent as i charge group */
820 /* Loop over the atoms in the i charge group */
821 for (i = 0; i < nicg; i++)
824 gid = GID(igid, jgid, ngid);
827 /* Create new i_atom for each energy group */
828 if (bDoVdW && bDoCoul)
830 new_i_nblist(vdwc, i_atom, shift, gid);
834 new_i_nblist(vdw, i_atom, shift, gid);
838 new_i_nblist(coul, i_atom, shift, gid);
840 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
841 bDoCoul_i = (bDoCoul && qi != 0);
843 if (bDoVdW_i || bDoCoul_i)
845 /* Loop over the j charge groups */
846 for (j = 0; (j < nj); j++)
850 /* Check for large charge groups */
861 /* Finally loop over the atoms in the j-charge group */
862 for (jj = jj0; jj < jj1; jj++)
864 bNotEx = NOTEXCL(bExcl, i, jj);
872 add_j_to_nblist(coul, jj, bLR);
877 if (bHaveVdW[type[jj]])
879 add_j_to_nblist(vdw, jj, bLR);
884 if (bHaveVdW[type[jj]])
888 add_j_to_nblist(vdwc, jj, bLR);
892 add_j_to_nblist(vdw, jj, bLR);
895 else if (charge[jj] != 0)
897 add_j_to_nblist(coul, jj, bLR);
905 close_i_nblist(coul);
906 close_i_nblist(vdwc);
912 /* we are doing free energy */
913 vdwc_free = &nlist[eNL_VDWQQ_FREE];
914 vdw_free = &nlist[eNL_VDW_FREE];
915 coul_free = &nlist[eNL_QQ_FREE];
916 /* Loop over the atoms in the i charge group */
917 for (i = 0; i < nicg; i++)
920 gid = GID(igid, jgid, ngid);
922 qiB = chargeB[i_atom];
924 /* Create new i_atom for each energy group */
925 if (bDoVdW && bDoCoul)
927 new_i_nblist(vdwc, i_atom, shift, gid);
931 new_i_nblist(vdw, i_atom, shift, gid);
935 new_i_nblist(coul, i_atom, shift, gid);
938 new_i_nblist(vdw_free, i_atom, shift, gid);
939 new_i_nblist(coul_free, i_atom, shift, gid);
940 new_i_nblist(vdwc_free, i_atom, shift, gid);
942 bDoVdW_i = (bDoVdW &&
943 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
944 bDoCoul_i = (bDoCoul && (qi != 0 || qiB != 0));
945 /* For TIP4P the first atom does not have a charge,
946 * but the last three do. So we should still put an atom
947 * without LJ but with charge in the water-atom neighborlist
948 * for a TIP4p i charge group.
949 * For SPC type water the first atom has LJ and charge,
950 * so there is no such problem.
952 if (iwater == esolNO)
954 bDoCoul_i_sol = bDoCoul_i;
958 bDoCoul_i_sol = bDoCoul;
961 if (bDoVdW_i || bDoCoul_i_sol)
963 /* Loop over the j charge groups */
964 for (j = 0; (j < nj); j++)
968 /* Check for large charge groups */
979 /* Finally loop over the atoms in the j-charge group */
980 bFree = bPert[i_atom];
981 for (jj = jj0; (jj < jj1); jj++)
983 bFreeJ = bFree || bPert[jj];
984 /* Complicated if, because the water H's should also
985 * see perturbed j-particles
987 if (iwater == esolNO || i == 0 || bFreeJ)
989 bNotEx = NOTEXCL(bExcl, i, jj);
997 if (charge[jj] != 0 || chargeB[jj] != 0)
999 add_j_to_nblist(coul_free, jj, bLR);
1002 else if (!bDoCoul_i)
1004 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1006 add_j_to_nblist(vdw_free, jj, bLR);
1011 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1013 if (charge[jj] != 0 || chargeB[jj] != 0)
1015 add_j_to_nblist(vdwc_free, jj, bLR);
1019 add_j_to_nblist(vdw_free, jj, bLR);
1022 else if (charge[jj] != 0 || chargeB[jj] != 0)
1024 add_j_to_nblist(coul_free, jj, bLR);
1030 /* This is done whether or not bWater is set */
1031 if (charge[jj] != 0)
1033 add_j_to_nblist(coul, jj, bLR);
1036 else if (!bDoCoul_i_sol)
1038 if (bHaveVdW[type[jj]])
1040 add_j_to_nblist(vdw, jj, bLR);
1045 if (bHaveVdW[type[jj]])
1047 if (charge[jj] != 0)
1049 add_j_to_nblist(vdwc, jj, bLR);
1053 add_j_to_nblist(vdw, jj, bLR);
1056 else if (charge[jj] != 0)
1058 add_j_to_nblist(coul, jj, bLR);
1066 close_i_nblist(vdw);
1067 close_i_nblist(coul);
1068 close_i_nblist(vdwc);
1069 close_i_nblist(vdw_free);
1070 close_i_nblist(coul_free);
1071 close_i_nblist(vdwc_free);
1077 put_in_list_adress(gmx_bool bHaveVdW[],
1093 /* The a[] index has been removed,
1094 * to put it back in i_atom should be a[i0] and jj should be a[jj].
1099 t_nblist * vdwc_adress = NULL;
1100 t_nblist * vdw_adress = NULL;
1101 t_nblist * coul_adress = NULL;
1102 t_nblist * vdwc_ww = NULL;
1103 t_nblist * coul_ww = NULL;
1105 int i, j, jcg, igid, gid, nbl_ind, nbl_ind_adress;
1106 atom_id jj, jj0, jj1, i_atom;
1111 real *charge, *chargeB;
1113 real qi, qiB, qq, rlj;
1114 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
1115 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
1117 gmx_bool j_all_atom;
1119 t_nblist *nlist, *nlist_adress;
1120 gmx_bool bEnergyGroupCG;
1122 /* Copy some pointers */
1123 cginfo = fr->cginfo;
1124 charge = md->chargeA;
1125 chargeB = md->chargeB;
1128 bPert = md->bPerturbed;
1131 /* Get atom range */
1133 nicg = index[icg+1]-i0;
1135 /* Get the i charge group info */
1136 igid = GET_CGINFO_GID(cginfo[icg]);
1138 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
1142 gmx_fatal(FARGS, "AdResS does not support free energy pertubation\n");
1145 /* Unpack pointers to neighbourlist structs */
1146 if (fr->nnblists == 2)
1153 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
1154 nbl_ind_adress = nbl_ind+fr->nnblists/2;
1158 nlist = fr->nblists[nbl_ind].nlist_lr;
1159 nlist_adress = fr->nblists[nbl_ind_adress].nlist_lr;
1163 nlist = fr->nblists[nbl_ind].nlist_sr;
1164 nlist_adress = fr->nblists[nbl_ind_adress].nlist_sr;
1168 vdwc = &nlist[eNL_VDWQQ];
1169 vdw = &nlist[eNL_VDW];
1170 coul = &nlist[eNL_QQ];
1172 vdwc_adress = &nlist_adress[eNL_VDWQQ];
1173 vdw_adress = &nlist_adress[eNL_VDW];
1174 coul_adress = &nlist_adress[eNL_QQ];
1176 /* We do not support solvent optimization with AdResS for now.
1177 For this we would need hybrid solvent-other kernels */
1179 /* no solvent as i charge group */
1180 /* Loop over the atoms in the i charge group */
1181 for (i = 0; i < nicg; i++)
1184 gid = GID(igid, jgid, ngid);
1185 qi = charge[i_atom];
1187 /* Create new i_atom for each energy group */
1188 if (bDoVdW && bDoCoul)
1190 new_i_nblist(vdwc, i_atom, shift, gid);
1191 new_i_nblist(vdwc_adress, i_atom, shift, gid);
1196 new_i_nblist(vdw, i_atom, shift, gid);
1197 new_i_nblist(vdw_adress, i_atom, shift, gid);
1202 new_i_nblist(coul, i_atom, shift, gid);
1203 new_i_nblist(coul_adress, i_atom, shift, gid);
1205 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
1206 bDoCoul_i = (bDoCoul && qi != 0);
1208 /* Here we find out whether the energy groups interaction belong to a
1209 * coarse-grained (vsite) or atomistic interaction. Note that, beacuse
1210 * interactions between coarse-grained and other (atomistic) energygroups
1211 * are excluded automatically by grompp, it is sufficient to check for
1212 * the group id of atom i (igid) */
1213 bEnergyGroupCG = !egp_explicit(fr, igid);
1215 if (bDoVdW_i || bDoCoul_i)
1217 /* Loop over the j charge groups */
1218 for (j = 0; (j < nj); j++)
1222 /* Check for large charge groups */
1233 /* Finally loop over the atoms in the j-charge group */
1234 for (jj = jj0; jj < jj1; jj++)
1236 bNotEx = NOTEXCL(bExcl, i, jj);
1238 /* Now we have to exclude interactions which will be zero
1239 * anyway due to the AdResS weights (in previous implementations
1240 * this was done in the force kernel). This is necessary as
1241 * pure interactions (those with b_hybrid=false, i.e. w_i*w_j==1 or 0)
1242 * are put into neighbour lists which will be passed to the
1243 * standard (optimized) kernels for speed. The interactions with
1244 * b_hybrid=true are placed into the _adress neighbour lists and
1245 * processed by the generic AdResS kernel.
1247 if ( (bEnergyGroupCG &&
1248 wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS ) ||
1249 ( !bEnergyGroupCG && wf[jj] <= GMX_REAL_EPS ) )
1254 b_hybrid = !((wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS) ||
1255 (wf[i_atom] <= GMX_REAL_EPS && wf[jj] <= GMX_REAL_EPS));
1261 if (charge[jj] != 0)
1265 add_j_to_nblist(coul, jj, bLR);
1269 add_j_to_nblist(coul_adress, jj, bLR);
1273 else if (!bDoCoul_i)
1275 if (bHaveVdW[type[jj]])
1279 add_j_to_nblist(vdw, jj, bLR);
1283 add_j_to_nblist(vdw_adress, jj, bLR);
1289 if (bHaveVdW[type[jj]])
1291 if (charge[jj] != 0)
1295 add_j_to_nblist(vdwc, jj, bLR);
1299 add_j_to_nblist(vdwc_adress, jj, bLR);
1306 add_j_to_nblist(vdw, jj, bLR);
1310 add_j_to_nblist(vdw_adress, jj, bLR);
1315 else if (charge[jj] != 0)
1319 add_j_to_nblist(coul, jj, bLR);
1323 add_j_to_nblist(coul_adress, jj, bLR);
1332 close_i_nblist(vdw);
1333 close_i_nblist(coul);
1334 close_i_nblist(vdwc);
1335 close_i_nblist(vdw_adress);
1336 close_i_nblist(coul_adress);
1337 close_i_nblist(vdwc_adress);
1343 put_in_list_qmmm(gmx_bool gmx_unused bHaveVdW[],
1345 t_mdatoms gmx_unused * md,
1355 gmx_bool gmx_unused bDoVdW,
1356 gmx_bool gmx_unused bDoCoul,
1357 int gmx_unused solvent_opt)
1360 int i, j, jcg, igid, gid;
1361 atom_id jj, jj0, jj1, i_atom;
1365 /* Get atom range */
1367 nicg = index[icg+1]-i0;
1369 /* Get the i charge group info */
1370 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1372 coul = &fr->QMMMlist;
1374 /* Loop over atoms in the ith charge group */
1375 for (i = 0; i < nicg; i++)
1378 gid = GID(igid, jgid, ngid);
1379 /* Create new i_atom for each energy group */
1380 new_i_nblist(coul, i_atom, shift, gid);
1382 /* Loop over the j charge groups */
1383 for (j = 0; j < nj; j++)
1387 /* Charge groups cannot have QM and MM atoms simultaneously */
1392 /* Finally loop over the atoms in the j-charge group */
1393 for (jj = jj0; jj < jj1; jj++)
1395 bNotEx = NOTEXCL(bExcl, i, jj);
1398 add_j_to_nblist(coul, jj, bLR);
1403 close_i_nblist(coul);
1408 put_in_list_cg(gmx_bool gmx_unused bHaveVdW[],
1410 t_mdatoms gmx_unused * md,
1420 gmx_bool gmx_unused bDoVdW,
1421 gmx_bool gmx_unused bDoCoul,
1422 int gmx_unused solvent_opt)
1425 int igid, gid, nbl_ind;
1429 cginfo = fr->cginfo[icg];
1431 igid = GET_CGINFO_GID(cginfo);
1432 gid = GID(igid, jgid, ngid);
1434 /* Unpack pointers to neighbourlist structs */
1435 if (fr->nnblists == 1)
1441 nbl_ind = fr->gid2nblists[gid];
1445 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1449 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1452 /* Make a new neighbor list for charge group icg.
1453 * Currently simply one neighbor list is made with LJ and Coulomb.
1454 * If required, zero interactions could be removed here
1455 * or in the force loop.
1457 new_i_nblist(vdwc, index[icg], shift, gid);
1458 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1460 for (j = 0; (j < nj); j++)
1463 /* Skip the icg-icg pairs if all self interactions are excluded */
1464 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1466 /* Here we add the j charge group jcg to the list,
1467 * exclusions are also added to the list.
1469 add_j_to_nblist_cg(vdwc, index[jcg], index[jcg+1], bExcl, icg == jcg, bLR);
1473 close_i_nblist(vdwc);
1476 static void setexcl(atom_id start, atom_id end, t_blocka *excl, gmx_bool b,
1483 for (i = start; i < end; i++)
1485 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1487 SETEXCL(bexcl, i-start, excl->a[k]);
1493 for (i = start; i < end; i++)
1495 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1497 RMEXCL(bexcl, i-start, excl->a[k]);
1503 int calc_naaj(int icg, int cgtot)
1507 if ((cgtot % 2) == 1)
1509 /* Odd number of charge groups, easy */
1510 naaj = 1 + (cgtot/2);
1512 else if ((cgtot % 4) == 0)
1514 /* Multiple of four is hard */
1551 fprintf(log, "naaj=%d\n", naaj);
1557 /************************************************
1559 * S I M P L E C O R E S T U F F
1561 ************************************************/
1563 static real calc_image_tric(rvec xi, rvec xj, matrix box,
1564 rvec b_inv, int *shift)
1566 /* This code assumes that the cut-off is smaller than
1567 * a half times the smallest diagonal element of the box.
1569 const real h25 = 2.5;
1574 /* Compute diff vector */
1575 dz = xj[ZZ] - xi[ZZ];
1576 dy = xj[YY] - xi[YY];
1577 dx = xj[XX] - xi[XX];
1579 /* Perform NINT operation, using trunc operation, therefore
1580 * we first add 2.5 then subtract 2 again
1582 tz = dz*b_inv[ZZ] + h25;
1584 dz -= tz*box[ZZ][ZZ];
1585 dy -= tz*box[ZZ][YY];
1586 dx -= tz*box[ZZ][XX];
1588 ty = dy*b_inv[YY] + h25;
1590 dy -= ty*box[YY][YY];
1591 dx -= ty*box[YY][XX];
1593 tx = dx*b_inv[XX]+h25;
1595 dx -= tx*box[XX][XX];
1597 /* Distance squared */
1598 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1600 *shift = XYZ2IS(tx, ty, tz);
1605 static real calc_image_rect(rvec xi, rvec xj, rvec box_size,
1606 rvec b_inv, int *shift)
1608 const real h15 = 1.5;
1614 /* Compute diff vector */
1615 dx = xj[XX] - xi[XX];
1616 dy = xj[YY] - xi[YY];
1617 dz = xj[ZZ] - xi[ZZ];
1619 /* Perform NINT operation, using trunc operation, therefore
1620 * we first add 1.5 then subtract 1 again
1622 tx = dx*b_inv[XX] + h15;
1623 ty = dy*b_inv[YY] + h15;
1624 tz = dz*b_inv[ZZ] + h15;
1629 /* Correct diff vector for translation */
1630 ddx = tx*box_size[XX] - dx;
1631 ddy = ty*box_size[YY] - dy;
1632 ddz = tz*box_size[ZZ] - dz;
1634 /* Distance squared */
1635 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1637 *shift = XYZ2IS(tx, ty, tz);
1642 static void add_simple(t_ns_buf *nsbuf, int nrj, atom_id cg_j,
1643 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1644 int icg, int jgid, t_block *cgs, t_excl bexcl[],
1645 int shift, t_forcerec *fr, put_in_list_t *put_in_list)
1647 if (nsbuf->nj + nrj > MAX_CG)
1649 put_in_list(bHaveVdW, ngid, md, icg, jgid, nsbuf->ncg, nsbuf->jcg,
1650 cgs->index, bexcl, shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1651 /* Reset buffer contents */
1652 nsbuf->ncg = nsbuf->nj = 0;
1654 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1658 static void ns_inner_tric(rvec x[], int icg, int *i_egp_flags,
1659 int njcg, atom_id jcg[],
1660 matrix box, rvec b_inv, real rcut2,
1661 t_block *cgs, t_ns_buf **ns_buf,
1662 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1663 t_excl bexcl[], t_forcerec *fr,
1664 put_in_list_t *put_in_list)
1668 int *cginfo = fr->cginfo;
1669 atom_id cg_j, *cgindex;
1672 cgindex = cgs->index;
1674 for (j = 0; (j < njcg); j++)
1677 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1678 if (calc_image_tric(x[icg], x[cg_j], box, b_inv, &shift) < rcut2)
1680 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1681 if (!(i_egp_flags[jgid] & EGP_EXCL))
1683 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1684 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1691 static void ns_inner_rect(rvec x[], int icg, int *i_egp_flags,
1692 int njcg, atom_id jcg[],
1693 gmx_bool bBox, rvec box_size, rvec b_inv, real rcut2,
1694 t_block *cgs, t_ns_buf **ns_buf,
1695 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1696 t_excl bexcl[], t_forcerec *fr,
1697 put_in_list_t *put_in_list)
1701 int *cginfo = fr->cginfo;
1702 atom_id cg_j, *cgindex;
1705 cgindex = cgs->index;
1709 for (j = 0; (j < njcg); j++)
1712 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1713 if (calc_image_rect(x[icg], x[cg_j], box_size, b_inv, &shift) < rcut2)
1715 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1716 if (!(i_egp_flags[jgid] & EGP_EXCL))
1718 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1719 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1727 for (j = 0; (j < njcg); j++)
1730 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1731 if ((rcut2 == 0) || (distance2(x[icg], x[cg_j]) < rcut2))
1733 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1734 if (!(i_egp_flags[jgid] & EGP_EXCL))
1736 add_simple(&ns_buf[jgid][CENTRAL], nrj, cg_j,
1737 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, CENTRAL, fr,
1745 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1747 static int ns_simple_core(t_forcerec *fr,
1748 gmx_localtop_t *top,
1750 matrix box, rvec box_size,
1751 t_excl bexcl[], atom_id *aaj,
1752 int ngid, t_ns_buf **ns_buf,
1753 put_in_list_t *put_in_list, gmx_bool bHaveVdW[])
1757 int nsearch, icg, jcg, igid, i0, nri, nn;
1760 /* atom_id *i_atoms; */
1761 t_block *cgs = &(top->cgs);
1762 t_blocka *excl = &(top->excls);
1765 gmx_bool bBox, bTriclinic;
1768 rlist2 = sqr(fr->rlist);
1770 bBox = (fr->ePBC != epbcNONE);
1773 for (m = 0; (m < DIM); m++)
1775 b_inv[m] = divide_err(1.0, box_size[m]);
1777 bTriclinic = TRICLINIC(box);
1784 cginfo = fr->cginfo;
1787 for (icg = fr->cg0; (icg < fr->hcg); icg++)
1790 i0 = cgs->index[icg];
1791 nri = cgs->index[icg+1]-i0;
1792 i_atoms = &(cgs->a[i0]);
1793 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1794 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1796 igid = GET_CGINFO_GID(cginfo[icg]);
1797 i_egp_flags = fr->egp_flags + ngid*igid;
1798 setexcl(cgs->index[icg], cgs->index[icg+1], excl, TRUE, bexcl);
1800 naaj = calc_naaj(icg, cgs->nr);
1803 ns_inner_tric(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1804 box, b_inv, rlist2, cgs, ns_buf,
1805 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1809 ns_inner_rect(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1810 bBox, box_size, b_inv, rlist2, cgs, ns_buf,
1811 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1815 for (nn = 0; (nn < ngid); nn++)
1817 for (k = 0; (k < SHIFTS); k++)
1819 nsbuf = &(ns_buf[nn][k]);
1822 put_in_list(bHaveVdW, ngid, md, icg, nn, nsbuf->ncg, nsbuf->jcg,
1823 cgs->index, bexcl, k, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1824 nsbuf->ncg = nsbuf->nj = 0;
1828 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1829 setexcl(cgs->index[icg], cgs->index[icg+1], excl, FALSE, bexcl);
1831 close_neighbor_lists(fr, FALSE);
1836 /************************************************
1838 * N S 5 G R I D S T U F F
1840 ************************************************/
1842 static gmx_inline void get_dx(int Nx, real gridx, real rc2, int xgi, real x,
1843 int *dx0, int *dx1, real *dcx2)
1871 for (i = xgi0; i >= 0; i--)
1873 dcx = (i+1)*gridx-x;
1882 for (i = xgi1; i < Nx; i++)
1895 static gmx_inline void get_dx_dd(int Nx, real gridx, real rc2, int xgi, real x,
1896 int ncpddc, int shift_min, int shift_max,
1897 int *g0, int *g1, real *dcx2)
1900 int g_min, g_max, shift_home;
1933 g_min = (shift_min == shift_home ? 0 : ncpddc);
1934 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1941 else if (shift_max < 0)
1956 /* Check one grid cell down */
1957 dcx = ((*g0 - 1) + 1)*gridx - x;
1969 /* Check one grid cell up */
1970 dcx = (*g1 + 1)*gridx - x;
1982 #define sqr(x) ((x)*(x))
1983 #define calc_dx2(XI, YI, ZI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
1984 #define calc_cyl_dx2(XI, YI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
1985 /****************************************************
1987 * F A S T N E I G H B O R S E A R C H I N G
1989 * Optimized neighboursearching routine using grid
1990 * at least 1x1x1, see GROMACS manual
1992 ****************************************************/
1995 static void get_cutoff2(t_forcerec *fr, gmx_bool bDoLongRange,
1996 real *rvdw2, real *rcoul2,
1997 real *rs2, real *rm2, real *rl2)
1999 *rs2 = sqr(fr->rlist);
2001 if (bDoLongRange && fr->bTwinRange)
2003 /* With plain cut-off or RF we need to make the list exactly
2004 * up to the cut-off and the cut-off's can be different,
2005 * so we can not simply set them to rlistlong.
2006 * To keep this code compatible with (exotic) old cases,
2007 * we also create lists up to rvdw/rcoulomb for PME and Ewald.
2008 * The interaction check should correspond to:
2009 * !ir_vdw/coulomb_might_be_zero_at_cutoff from inputrec.c.
2011 if (((fr->vdwtype == evdwCUT || fr->vdwtype == evdwPME) &&
2012 fr->vdw_modifier == eintmodNONE) ||
2013 fr->rvdw <= fr->rlist)
2015 *rvdw2 = sqr(fr->rvdw);
2019 *rvdw2 = sqr(fr->rlistlong);
2021 if (((fr->eeltype == eelCUT ||
2022 (EEL_RF(fr->eeltype) && fr->eeltype != eelRF_ZERO) ||
2023 fr->eeltype == eelPME ||
2024 fr->eeltype == eelEWALD) &&
2025 fr->coulomb_modifier == eintmodNONE) ||
2026 fr->rcoulomb <= fr->rlist)
2028 *rcoul2 = sqr(fr->rcoulomb);
2032 *rcoul2 = sqr(fr->rlistlong);
2037 /* Workaround for a gcc -O3 or -ffast-math problem */
2041 *rm2 = min(*rvdw2, *rcoul2);
2042 *rl2 = max(*rvdw2, *rcoul2);
2045 static void init_nsgrid_lists(t_forcerec *fr, int ngid, gmx_ns_t *ns)
2047 real rvdw2, rcoul2, rs2, rm2, rl2;
2050 get_cutoff2(fr, TRUE, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2052 /* Short range buffers */
2053 snew(ns->nl_sr, ngid);
2055 snew(ns->nsr, ngid);
2056 snew(ns->nlr_ljc, ngid);
2057 snew(ns->nlr_one, ngid);
2059 /* Always allocate both list types, since rcoulomb might now change with PME load balancing */
2060 /* Long range VdW and Coul buffers */
2061 snew(ns->nl_lr_ljc, ngid);
2062 /* Long range VdW or Coul only buffers */
2063 snew(ns->nl_lr_one, ngid);
2065 for (j = 0; (j < ngid); j++)
2067 snew(ns->nl_sr[j], MAX_CG);
2068 snew(ns->nl_lr_ljc[j], MAX_CG);
2069 snew(ns->nl_lr_one[j], MAX_CG);
2074 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
2079 static int nsgrid_core(t_commrec *cr, t_forcerec *fr,
2080 matrix box, int ngid,
2081 gmx_localtop_t *top,
2083 t_excl bexcl[], gmx_bool *bExcludeAlleg,
2085 put_in_list_t *put_in_list,
2086 gmx_bool bHaveVdW[],
2087 gmx_bool bDoLongRange, gmx_bool bMakeQMMMnblist)
2090 atom_id **nl_lr_ljc, **nl_lr_one, **nl_sr;
2091 int *nlr_ljc, *nlr_one, *nsr;
2093 t_block *cgs = &(top->cgs);
2094 int *cginfo = fr->cginfo;
2095 /* atom_id *i_atoms,*cgsindex=cgs->index; */
2097 int cell_x, cell_y, cell_z;
2098 int d, tx, ty, tz, dx, dy, dz, cj;
2099 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2100 int zsh_ty, zsh_tx, ysh_tx;
2102 int dx0, dx1, dy0, dy1, dz0, dz1;
2103 int Nx, Ny, Nz, shift = -1, j, nrj, nns, nn = -1;
2104 real gridx, gridy, gridz, grid_x, grid_y, grid_z;
2105 real *dcx2, *dcy2, *dcz2;
2107 int cg0, cg1, icg = -1, cgsnr, i0, igid, nri, naaj, max_jcg;
2108 int jcg0, jcg1, jjcg, cgj0, jgid;
2109 int *grida, *gridnra, *gridind;
2110 gmx_bool rvdw_lt_rcoul, rcoul_lt_rvdw;
2111 rvec xi, *cgcm, grid_offset;
2112 real r2, rs2, rvdw2, rcoul2, rm2, rl2, XI, YI, ZI, dcx, dcy, dcz, tmp1, tmp2;
2114 gmx_bool bDomDec, bTriclinicX, bTriclinicY;
2119 bDomDec = DOMAINDECOMP(cr);
2122 bTriclinicX = ((YY < grid->npbcdim &&
2123 (!bDomDec || dd->nc[YY] == 1) && box[YY][XX] != 0) ||
2124 (ZZ < grid->npbcdim &&
2125 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][XX] != 0));
2126 bTriclinicY = (ZZ < grid->npbcdim &&
2127 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][YY] != 0);
2131 get_cutoff2(fr, bDoLongRange, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2133 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
2134 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
2136 if (bMakeQMMMnblist)
2144 nl_lr_ljc = ns->nl_lr_ljc;
2145 nl_lr_one = ns->nl_lr_one;
2146 nlr_ljc = ns->nlr_ljc;
2147 nlr_one = ns->nlr_one;
2155 gridind = grid->index;
2156 gridnra = grid->nra;
2159 gridx = grid->cell_size[XX];
2160 gridy = grid->cell_size[YY];
2161 gridz = grid->cell_size[ZZ];
2165 copy_rvec(grid->cell_offset, grid_offset);
2166 copy_ivec(grid->ncpddc, ncpddc);
2171 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2172 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
2173 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
2174 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
2175 if (zsh_tx != 0 && ysh_tx != 0)
2177 /* This could happen due to rounding, when both ratios are 0.5 */
2186 /* We only want a list for the test particle */
2195 /* Set the shift range */
2196 for (d = 0; d < DIM; d++)
2200 /* Check if we need periodicity shifts.
2201 * Without PBC or with domain decomposition we don't need them.
2203 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2210 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
2221 /* Loop over charge groups */
2222 for (icg = cg0; (icg < cg1); icg++)
2224 igid = GET_CGINFO_GID(cginfo[icg]);
2225 /* Skip this charge group if all energy groups are excluded! */
2226 if (bExcludeAlleg[igid])
2231 i0 = cgs->index[icg];
2233 if (bMakeQMMMnblist)
2235 /* Skip this charge group if it is not a QM atom while making a
2236 * QM/MM neighbourlist
2238 if (md->bQM[i0] == FALSE)
2240 continue; /* MM particle, go to next particle */
2243 /* Compute the number of charge groups that fall within the control
2246 naaj = calc_naaj(icg, cgsnr);
2253 /* make a normal neighbourlist */
2257 /* Get the j charge-group and dd cell shift ranges */
2258 dd_get_ns_ranges(cr->dd, icg, &jcg0, &jcg1, sh0, sh1);
2263 /* Compute the number of charge groups that fall within the control
2266 naaj = calc_naaj(icg, cgsnr);
2272 /* The i-particle is awlways the test particle,
2273 * so we want all j-particles
2275 max_jcg = cgsnr - 1;
2279 max_jcg = jcg1 - cgsnr;
2284 i_egp_flags = fr->egp_flags + igid*ngid;
2286 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2287 setexcl(i0, cgs->index[icg+1], &top->excls, TRUE, bexcl);
2289 ci2xyz(grid, icg, &cell_x, &cell_y, &cell_z);
2291 /* Changed iicg to icg, DvdS 990115
2292 * (but see consistency check above, DvdS 990330)
2295 fprintf(log, "icg=%5d, naaj=%5d, cell %d %d %d\n",
2296 icg, naaj, cell_x, cell_y, cell_z);
2298 /* Loop over shift vectors in three dimensions */
2299 for (tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
2301 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2302 /* Calculate range of cells in Z direction that have the shift tz */
2303 zgi = cell_z + tz*Nz;
2306 get_dx(Nz, gridz, rl2, zgi, ZI, &dz0, &dz1, dcz2);
2308 get_dx_dd(Nz, gridz, rl2, zgi, ZI-grid_offset[ZZ],
2309 ncpddc[ZZ], sh0[ZZ], sh1[ZZ], &dz0, &dz1, dcz2);
2315 for (ty = -shp[YY]; ty <= shp[YY]; ty++)
2317 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2318 /* Calculate range of cells in Y direction that have the shift ty */
2321 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2325 ygi = cell_y + ty*Ny;
2328 get_dx(Ny, gridy, rl2, ygi, YI, &dy0, &dy1, dcy2);
2330 get_dx_dd(Ny, gridy, rl2, ygi, YI-grid_offset[YY],
2331 ncpddc[YY], sh0[YY], sh1[YY], &dy0, &dy1, dcy2);
2337 for (tx = -shp[XX]; tx <= shp[XX]; tx++)
2339 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2340 /* Calculate range of cells in X direction that have the shift tx */
2343 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2347 xgi = cell_x + tx*Nx;
2350 get_dx(Nx, gridx, rl2, xgi*Nx, XI, &dx0, &dx1, dcx2);
2352 get_dx_dd(Nx, gridx, rl2, xgi, XI-grid_offset[XX],
2353 ncpddc[XX], sh0[XX], sh1[XX], &dx0, &dx1, dcx2);
2359 /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2360 * from the neigbour list as it will not interact */
2361 if (fr->adress_type != eAdressOff)
2363 if (md->wf[cgs->index[icg]] <= GMX_REAL_EPS && egp_explicit(fr, igid))
2368 /* Get shift vector */
2369 shift = XYZ2IS(tx, ty, tz);
2371 range_check(shift, 0, SHIFTS);
2373 for (nn = 0; (nn < ngid); nn++)
2380 fprintf(log, "shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2381 shift, dx0, dx1, dy0, dy1, dz0, dz1);
2382 fprintf(log, "cgcm: %8.3f %8.3f %8.3f\n", cgcm[icg][XX],
2383 cgcm[icg][YY], cgcm[icg][ZZ]);
2384 fprintf(log, "xi: %8.3f %8.3f %8.3f\n", XI, YI, ZI);
2386 for (dx = dx0; (dx <= dx1); dx++)
2388 tmp1 = rl2 - dcx2[dx];
2389 for (dy = dy0; (dy <= dy1); dy++)
2391 tmp2 = tmp1 - dcy2[dy];
2394 for (dz = dz0; (dz <= dz1); dz++)
2396 if (tmp2 > dcz2[dz])
2398 /* Find grid-cell cj in which possible neighbours are */
2399 cj = xyz2ci(Ny, Nz, dx, dy, dz);
2401 /* Check out how many cgs (nrj) there in this cell */
2404 /* Find the offset in the cg list */
2407 /* Check if all j's are out of range so we
2408 * can skip the whole cell.
2409 * Should save some time, especially with DD.
2412 (grida[cgj0] >= max_jcg &&
2413 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2419 for (j = 0; (j < nrj); j++)
2421 jjcg = grida[cgj0+j];
2423 /* check whether this guy is in range! */
2424 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2427 r2 = calc_dx2(XI, YI, ZI, cgcm[jjcg]);
2430 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2431 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2432 /* check energy group exclusions */
2433 if (!(i_egp_flags[jgid] & EGP_EXCL))
2437 if (nsr[jgid] >= MAX_CG)
2439 /* Add to short-range list */
2440 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2441 nsr[jgid], nl_sr[jgid],
2442 cgs->index, /* cgsatoms, */ bexcl,
2443 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2446 nl_sr[jgid][nsr[jgid]++] = jjcg;
2450 if (nlr_ljc[jgid] >= MAX_CG)
2452 /* Add to LJ+coulomb long-range list */
2453 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2454 nlr_ljc[jgid], nl_lr_ljc[jgid], top->cgs.index,
2455 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2458 nl_lr_ljc[jgid][nlr_ljc[jgid]++] = jjcg;
2462 if (nlr_one[jgid] >= MAX_CG)
2464 /* Add to long-range list with only coul, or only LJ */
2465 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2466 nlr_one[jgid], nl_lr_one[jgid], top->cgs.index,
2467 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2470 nl_lr_one[jgid][nlr_one[jgid]++] = jjcg;
2482 /* CHECK whether there is anything left in the buffers */
2483 for (nn = 0; (nn < ngid); nn++)
2487 put_in_list(bHaveVdW, ngid, md, icg, nn, nsr[nn], nl_sr[nn],
2488 cgs->index, /* cgsatoms, */ bexcl,
2489 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2492 if (nlr_ljc[nn] > 0)
2494 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_ljc[nn],
2495 nl_lr_ljc[nn], top->cgs.index,
2496 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2499 if (nlr_one[nn] > 0)
2501 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_one[nn],
2502 nl_lr_one[nn], top->cgs.index,
2503 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2509 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2510 setexcl(cgs->index[icg], cgs->index[icg+1], &top->excls, FALSE, bexcl);
2512 /* No need to perform any left-over force calculations anymore (as we used to do here)
2513 * since we now save the proper long-range lists for later evaluation.
2518 /* Close neighbourlists */
2519 close_neighbor_lists(fr, bMakeQMMMnblist);
2524 void ns_realloc_natoms(gmx_ns_t *ns, int natoms)
2528 if (natoms > ns->nra_alloc)
2530 ns->nra_alloc = over_alloc_dd(natoms);
2531 srenew(ns->bexcl, ns->nra_alloc);
2532 for (i = 0; i < ns->nra_alloc; i++)
2539 void init_ns(FILE *fplog, const t_commrec *cr,
2540 gmx_ns_t *ns, t_forcerec *fr,
2541 const gmx_mtop_t *mtop)
2543 int mt, icg, nr_in_cg, maxcg, i, j, jcg, ngid, ncg;
2547 /* Compute largest charge groups size (# atoms) */
2549 for (mt = 0; mt < mtop->nmoltype; mt++)
2551 cgs = &mtop->moltype[mt].cgs;
2552 for (icg = 0; (icg < cgs->nr); icg++)
2554 nr_in_cg = max(nr_in_cg, (int)(cgs->index[icg+1]-cgs->index[icg]));
2558 /* Verify whether largest charge group is <= max cg.
2559 * This is determined by the type of the local exclusion type
2560 * Exclusions are stored in bits. (If the type is not large
2561 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2563 maxcg = sizeof(t_excl)*8;
2564 if (nr_in_cg > maxcg)
2566 gmx_fatal(FARGS, "Max #atoms in a charge group: %d > %d\n",
2570 ngid = mtop->groups.grps[egcENER].nr;
2571 snew(ns->bExcludeAlleg, ngid);
2572 for (i = 0; i < ngid; i++)
2574 ns->bExcludeAlleg[i] = TRUE;
2575 for (j = 0; j < ngid; j++)
2577 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2579 ns->bExcludeAlleg[i] = FALSE;
2587 ns->grid = init_grid(fplog, fr);
2588 init_nsgrid_lists(fr, ngid, ns);
2593 snew(ns->ns_buf, ngid);
2594 for (i = 0; (i < ngid); i++)
2596 snew(ns->ns_buf[i], SHIFTS);
2598 ncg = ncg_mtop(mtop);
2599 snew(ns->simple_aaj, 2*ncg);
2600 for (jcg = 0; (jcg < ncg); jcg++)
2602 ns->simple_aaj[jcg] = jcg;
2603 ns->simple_aaj[jcg+ncg] = jcg;
2607 /* Create array that determines whether or not atoms have VdW */
2608 snew(ns->bHaveVdW, fr->ntype);
2609 for (i = 0; (i < fr->ntype); i++)
2611 for (j = 0; (j < fr->ntype); j++)
2613 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2615 ((BHAMA(fr->nbfp, fr->ntype, i, j) != 0) ||
2616 (BHAMB(fr->nbfp, fr->ntype, i, j) != 0) ||
2617 (BHAMC(fr->nbfp, fr->ntype, i, j) != 0)) :
2618 ((C6(fr->nbfp, fr->ntype, i, j) != 0) ||
2619 (C12(fr->nbfp, fr->ntype, i, j) != 0))));
2624 pr_bvec(debug, 0, "bHaveVdW", ns->bHaveVdW, fr->ntype, TRUE);
2629 if (!DOMAINDECOMP(cr))
2631 ns_realloc_natoms(ns, mtop->natoms);
2634 ns->nblist_initialized = FALSE;
2636 /* nbr list debug dump */
2638 char *ptr = getenv("GMX_DUMP_NL");
2641 ns->dump_nl = strtol(ptr, NULL, 10);
2644 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2655 int search_neighbours(FILE *log, t_forcerec *fr,
2657 gmx_localtop_t *top,
2658 gmx_groups_t *groups,
2660 t_nrnb *nrnb, t_mdatoms *md,
2662 gmx_bool bDoLongRangeNS)
2664 t_block *cgs = &(top->cgs);
2665 rvec box_size, grid_x0, grid_x1;
2667 real min_size, grid_dens;
2671 gmx_bool *i_egp_flags;
2672 int cg_start, cg_end, start, end;
2675 gmx_domdec_zones_t *dd_zones;
2676 put_in_list_t *put_in_list;
2680 /* Set some local variables */
2682 ngid = groups->grps[egcENER].nr;
2684 for (m = 0; (m < DIM); m++)
2686 box_size[m] = box[m][m];
2689 if (fr->ePBC != epbcNONE)
2691 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC, box))
2693 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.");
2697 min_size = min(box_size[XX], min(box_size[YY], box_size[ZZ]));
2698 if (2*fr->rlistlong >= min_size)
2700 gmx_fatal(FARGS, "One of the box diagonal elements has become smaller than twice the cut-off length.");
2705 if (DOMAINDECOMP(cr))
2707 ns_realloc_natoms(ns, cgs->index[cgs->nr]);
2711 /* Reset the neighbourlists */
2712 reset_neighbor_lists(fr, TRUE, TRUE);
2714 if (bGrid && bFillGrid)
2718 if (DOMAINDECOMP(cr))
2720 dd_zones = domdec_zones(cr->dd);
2726 get_nsgrid_boundaries(grid->nboundeddim, box, NULL, NULL, NULL, NULL,
2727 cgs->nr, fr->cg_cm, grid_x0, grid_x1, &grid_dens);
2729 grid_first(log, grid, NULL, NULL, box, grid_x0, grid_x1,
2730 fr->rlistlong, grid_dens);
2737 if (DOMAINDECOMP(cr))
2740 fill_grid(dd_zones, grid, end, -1, end, fr->cg_cm);
2742 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2746 fill_grid(NULL, grid, cgs->nr, fr->cg0, fr->hcg, fr->cg_cm);
2747 grid->icg0 = fr->cg0;
2748 grid->icg1 = fr->hcg;
2752 calc_elemnr(grid, start, end, cgs->nr);
2754 grid_last(grid, start, end, cgs->nr);
2759 print_grid(debug, grid);
2764 /* Set the grid cell index for the test particle only.
2765 * The cell to cg index is not corrected, but that does not matter.
2767 fill_grid(NULL, ns->grid, fr->hcg, fr->hcg-1, fr->hcg, fr->cg_cm);
2771 if (fr->adress_type == eAdressOff)
2773 if (!fr->ns.bCGlist)
2775 put_in_list = put_in_list_at;
2779 put_in_list = put_in_list_cg;
2784 put_in_list = put_in_list_adress;
2791 nsearch = nsgrid_core(cr, fr, box, ngid, top,
2792 grid, ns->bexcl, ns->bExcludeAlleg,
2793 md, put_in_list, ns->bHaveVdW,
2794 bDoLongRangeNS, FALSE);
2796 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2797 * the classical calculation. The charge-charge interaction
2798 * between QM and MM atoms is handled in the QMMM core calculation
2799 * (see QMMM.c). The VDW however, we'd like to compute classically
2800 * and the QM MM atom pairs have just been put in the
2801 * corresponding neighbourlists. in case of QMMM we still need to
2802 * fill a special QMMM neighbourlist that contains all neighbours
2803 * of the QM atoms. If bQMMM is true, this list will now be made:
2805 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
2807 nsearch += nsgrid_core(cr, fr, box, ngid, top,
2808 grid, ns->bexcl, ns->bExcludeAlleg,
2809 md, put_in_list_qmmm, ns->bHaveVdW,
2810 bDoLongRangeNS, TRUE);
2815 nsearch = ns_simple_core(fr, top, md, box, box_size,
2816 ns->bexcl, ns->simple_aaj,
2817 ngid, ns->ns_buf, put_in_list, ns->bHaveVdW);
2825 inc_nrnb(nrnb, eNR_NS, nsearch);
2826 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2831 int natoms_beyond_ns_buffer(t_inputrec *ir, t_forcerec *fr, t_block *cgs,
2832 matrix scale_tot, rvec *x)
2834 int cg0, cg1, cg, a0, a1, a, i, j;
2835 real rint, hbuf2, scale;
2837 gmx_bool bIsotropic;
2842 rint = max(ir->rcoulomb, ir->rvdw);
2843 if (ir->rlist < rint)
2845 gmx_fatal(FARGS, "The neighbor search buffer has negative size: %f nm",
2853 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2855 hbuf2 = sqr(0.5*(ir->rlist - rint));
2856 for (cg = cg0; cg < cg1; cg++)
2858 a0 = cgs->index[cg];
2859 a1 = cgs->index[cg+1];
2860 for (a = a0; a < a1; a++)
2862 if (distance2(cg_cm[cg], x[a]) > hbuf2)
2872 scale = scale_tot[0][0];
2873 for (i = 1; i < DIM; i++)
2875 /* With anisotropic scaling, the original spherical ns volumes become
2876 * ellipsoids. To avoid costly transformations we use the minimum
2877 * eigenvalue of the scaling matrix for determining the buffer size.
2878 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2880 scale = min(scale, scale_tot[i][i]);
2881 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2885 for (j = 0; j < i; j++)
2887 if (scale_tot[i][j] != 0)
2893 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2896 for (cg = cg0; cg < cg1; cg++)
2898 svmul(scale, cg_cm[cg], cgsc);
2899 a0 = cgs->index[cg];
2900 a1 = cgs->index[cg+1];
2901 for (a = a0; a < a1; a++)
2903 if (distance2(cgsc, x[a]) > hbuf2)
2912 /* Anistropic scaling */
2913 for (cg = cg0; cg < cg1; cg++)
2915 /* Since scale_tot contains the transpose of the scaling matrix,
2916 * we need to multiply with the transpose.
2918 tmvmul_ur0(scale_tot, cg_cm[cg], cgsc);
2919 a0 = cgs->index[cg];
2920 a1 = cgs->index[cg+1];
2921 for (a = a0; a < a1; a++)
2923 if (distance2(cgsc, x[a]) > hbuf2)