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,2015,2017,2018, 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.
47 #include "gromacs/domdec/domdec.h"
48 #include "gromacs/domdec/domdec_struct.h"
49 #include "gromacs/gmxlib/network.h"
50 #include "gromacs/gmxlib/nrnb.h"
51 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/utilities.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/math/vecdump.h"
56 #include "gromacs/mdlib/nsgrid.h"
57 #include "gromacs/mdlib/qmmm.h"
58 #include "gromacs/mdtypes/commrec.h"
59 #include "gromacs/mdtypes/forcerec.h"
60 #include "gromacs/mdtypes/group.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/mdtypes/md_enums.h"
63 #include "gromacs/pbcutil/ishift.h"
64 #include "gromacs/pbcutil/pbc.h"
65 #include "gromacs/topology/mtop_util.h"
66 #include "gromacs/utility/fatalerror.h"
67 #include "gromacs/utility/smalloc.h"
70 * E X C L U S I O N H A N D L I N G
74 static void SETEXCL_(t_excl e[], int i, int j)
78 static void RMEXCL_(t_excl e[], int i, int j)
80 e[j] = e[j] & ~(1<<i);
82 static gmx_bool ISEXCL_(t_excl e[], int i, int j)
84 return (gmx_bool)(e[j] & (1<<i));
86 static gmx_bool NOTEXCL_(t_excl e[], int i, int j)
88 return !(ISEXCL(e, i, j));
91 #define SETEXCL(e, i, j) (e)[int(j)] |= (1<<(int(i)))
92 #define RMEXCL(e, i, j) (e)[int(j)] &= (~(1<<(int(i))))
93 #define ISEXCL(e, i, j) static_cast<gmx_bool>((e)[(int(j))] & (1<<(int(i))))
94 #define NOTEXCL(e, i, j) !(ISEXCL(e, i, j))
98 round_up_to_simd_width(int length, int simd_width)
102 offset = (simd_width > 0) ? length % simd_width : 0;
104 return (offset == 0) ? length : length-offset+simd_width;
106 /************************************************
108 * U T I L I T I E S F O R N S
110 ************************************************/
112 void reallocate_nblist(t_nblist *nl)
116 fprintf(debug, "reallocating neigborlist (ielec=%d, ivdw=%d, igeometry=%d, type=%d), maxnri=%d\n",
117 nl->ielec, nl->ivdw, nl->igeometry, nl->type, nl->maxnri);
119 srenew(nl->iinr, nl->maxnri);
120 if (nl->igeometry == GMX_NBLIST_GEOMETRY_CG_CG)
122 srenew(nl->iinr_end, nl->maxnri);
124 srenew(nl->gid, nl->maxnri);
125 srenew(nl->shift, nl->maxnri);
126 srenew(nl->jindex, nl->maxnri+1);
130 static void init_nblist(FILE *log, t_nblist *nl_sr,
132 int ivdw, int ivdwmod,
133 int ielec, int ielecmod,
134 int igeometry, int type,
135 gmx_bool bElecAndVdwSwitchDiffers)
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(log, 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;
181 nl->jindex = nullptr;
183 nl->excl_fep = nullptr;
184 reallocate_nblist(nl);
189 fprintf(debug, "Initiating neighbourlist (ielec=%d, ivdw=%d, type=%d) for %s interactions,\nwith %d SR atoms.\n",
190 nl->ielec, nl->ivdw, nl->type, gmx_nblist_geometry_names[nl->igeometry], maxsr);
195 void init_neighbor_list(FILE *log, t_forcerec *fr, int homenr)
197 int maxsr, maxsr_wat;
198 int ielec, ivdw, ielecmod, ivdwmod, type;
199 int igeometry_def, igeometry_w, igeometry_ww;
201 gmx_bool bElecAndVdwSwitchDiffers;
204 /* maxsr = homenr-fr->nWatMol*3; */
209 gmx_fatal(FARGS, "%s, %d: Negative number of short range atoms.\n"
210 "Call your GROMACS dealer for assistance.", __FILE__, __LINE__);
212 /* This is just for initial allocation, so we do not reallocate
213 * all the nlist arrays many times in a row.
214 * The numbers seem very accurate, but they are uncritical.
216 maxsr_wat = std::min(fr->nWatMol, (homenr+2)/3);
218 /* Determine the values for ielec/ivdw. */
219 ielec = fr->nbkernel_elec_interaction;
220 ivdw = fr->nbkernel_vdw_interaction;
221 ielecmod = fr->nbkernel_elec_modifier;
222 ivdwmod = fr->nbkernel_vdw_modifier;
223 type = GMX_NBLIST_INTERACTION_STANDARD;
224 bElecAndVdwSwitchDiffers = ( (fr->ic->rcoulomb_switch != fr->ic->rvdw_switch) || (fr->ic->rcoulomb != fr->ic->rvdw));
226 fr->ns->bCGlist = (getenv("GMX_NBLISTCG") != nullptr);
227 if (!fr->ns->bCGlist)
229 igeometry_def = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
233 igeometry_def = GMX_NBLIST_GEOMETRY_CG_CG;
236 fprintf(log, "\nUsing charge-group - charge-group neighbor lists and kernels\n\n");
240 if (fr->solvent_opt == esolTIP4P)
242 igeometry_w = GMX_NBLIST_GEOMETRY_WATER4_PARTICLE;
243 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER4_WATER4;
247 igeometry_w = GMX_NBLIST_GEOMETRY_WATER3_PARTICLE;
248 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER3_WATER3;
251 for (i = 0; i < fr->nnblists; i++)
253 nbl = &(fr->nblists[i]);
255 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ],
256 maxsr, ivdw, ivdwmod, ielec, ielecmod, igeometry_def, type, bElecAndVdwSwitchDiffers);
257 init_nblist(log, &nbl->nlist_sr[eNL_VDW],
258 maxsr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, igeometry_def, type, bElecAndVdwSwitchDiffers);
259 init_nblist(log, &nbl->nlist_sr[eNL_QQ],
260 maxsr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_def, type, bElecAndVdwSwitchDiffers);
261 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATER],
262 maxsr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_w, type, bElecAndVdwSwitchDiffers);
263 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATER],
264 maxsr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_w, type, bElecAndVdwSwitchDiffers);
265 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATERWATER],
266 maxsr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_ww, type, bElecAndVdwSwitchDiffers);
267 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATERWATER],
268 maxsr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_ww, type, bElecAndVdwSwitchDiffers);
270 /* Did we get the solvent loops so we can use optimized water kernels? */
271 if (nbl->nlist_sr[eNL_VDWQQ_WATER].kernelptr_vf == nullptr
272 || nbl->nlist_sr[eNL_QQ_WATER].kernelptr_vf == nullptr
273 || nbl->nlist_sr[eNL_VDWQQ_WATERWATER].kernelptr_vf == nullptr
274 || nbl->nlist_sr[eNL_QQ_WATERWATER].kernelptr_vf == nullptr)
276 fr->solvent_opt = esolNO;
279 fprintf(log, "Note: The available nonbonded kernels do not support water optimization - disabling.\n");
283 if (fr->efep != efepNO)
285 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_FREE],
286 maxsr, ivdw, ivdwmod, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY, bElecAndVdwSwitchDiffers);
287 init_nblist(log, &nbl->nlist_sr[eNL_VDW_FREE],
288 maxsr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY, bElecAndVdwSwitchDiffers);
289 init_nblist(log, &nbl->nlist_sr[eNL_QQ_FREE],
290 maxsr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY, bElecAndVdwSwitchDiffers);
294 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
296 if (nullptr == fr->QMMMlist)
298 snew(fr->QMMMlist, 1);
300 init_nblist(log, fr->QMMMlist,
301 maxsr, 0, 0, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_STANDARD, bElecAndVdwSwitchDiffers);
309 fr->ns->nblist_initialized = TRUE;
312 static void reset_nblist(t_nblist *nl)
322 static void reset_neighbor_lists(t_forcerec *fr)
328 /* only reset the short-range nblist */
329 reset_nblist(fr->QMMMlist);
332 for (n = 0; n < fr->nnblists; n++)
334 for (i = 0; i < eNL_NR; i++)
336 reset_nblist( &(fr->nblists[n].nlist_sr[i]) );
344 static inline void new_i_nblist(t_nblist *nlist, int i_atom, int shift, int gid)
346 int nri = nlist->nri;
348 /* Check whether we have to increase the i counter */
350 (nlist->iinr[nri] != i_atom) ||
351 (nlist->shift[nri] != shift) ||
352 (nlist->gid[nri] != gid))
354 /* This is something else. Now see if any entries have
355 * been added in the list of the previous atom.
358 ((nlist->jindex[nri+1] > nlist->jindex[nri]) &&
359 (nlist->gid[nri] != -1)))
361 /* If so increase the counter */
364 if (nlist->nri >= nlist->maxnri)
366 nlist->maxnri += over_alloc_large(nlist->nri);
367 reallocate_nblist(nlist);
370 /* Set the number of neighbours and the atom number */
371 nlist->jindex[nri+1] = nlist->jindex[nri];
372 nlist->iinr[nri] = i_atom;
373 nlist->gid[nri] = gid;
374 nlist->shift[nri] = shift;
378 /* Adding to previous list. First remove possible previous padding */
379 if (nlist->simd_padding_width > 1)
381 while (nlist->nrj > 0 && nlist->jjnr[nlist->nrj-1] < 0)
389 static inline void close_i_nblist(t_nblist *nlist)
391 int nri = nlist->nri;
396 /* Add elements up to padding. Since we allocate memory in units
397 * of the simd_padding width, we do not have to check for possible
398 * list reallocation here.
400 while ((nlist->nrj % nlist->simd_padding_width) != 0)
402 /* Use -4 here, so we can write forces for 4 atoms before real data */
403 nlist->jjnr[nlist->nrj++] = -4;
405 nlist->jindex[nri+1] = nlist->nrj;
407 len = nlist->nrj - nlist->jindex[nri];
408 /* If there are no j-particles we have to reduce the
409 * number of i-particles again, to prevent errors in the
412 if ((len == 0) && (nlist->nri > 0))
419 static inline void close_nblist(t_nblist *nlist)
421 /* Only close this nblist when it has been initialized.
422 * Avoid the creation of i-lists with no j-particles.
426 /* Some assembly kernels do not support empty lists,
427 * make sure here that we don't generate any empty lists.
428 * With the current ns code this branch is taken in two cases:
429 * No i-particles at all: nri=-1 here
430 * There are i-particles, but no j-particles; nri=0 here
436 /* Close list number nri by incrementing the count */
441 static inline void close_neighbor_lists(t_forcerec *fr, gmx_bool bMakeQMMMnblist)
447 close_nblist(fr->QMMMlist);
450 for (n = 0; n < fr->nnblists; n++)
452 for (i = 0; (i < eNL_NR); i++)
454 close_nblist(&(fr->nblists[n].nlist_sr[i]));
460 static inline void add_j_to_nblist(t_nblist *nlist, int j_atom)
462 int nrj = nlist->nrj;
464 if (nlist->nrj >= nlist->maxnrj)
466 nlist->maxnrj = round_up_to_simd_width(over_alloc_small(nlist->nrj + 1), nlist->simd_padding_width);
470 fprintf(debug, "Increasing SR nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
471 nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
474 srenew(nlist->jjnr, nlist->maxnrj);
477 nlist->jjnr[nrj] = j_atom;
481 static inline void add_j_to_nblist_cg(t_nblist *nlist,
482 int j_start, int j_end,
483 const t_excl *bexcl, gmx_bool i_is_j)
485 int nrj = nlist->nrj;
488 if (nlist->nrj >= nlist->maxnrj)
490 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
493 fprintf(debug, "Increasing SR nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
494 nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
497 srenew(nlist->jjnr, nlist->maxnrj);
498 srenew(nlist->jjnr_end, nlist->maxnrj);
499 srenew(nlist->excl, nlist->maxnrj*MAX_CGCGSIZE);
502 nlist->jjnr[nrj] = j_start;
503 nlist->jjnr_end[nrj] = j_end;
505 if (j_end - j_start > MAX_CGCGSIZE)
507 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);
510 /* Set the exclusions */
511 for (j = j_start; j < j_end; j++)
513 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
517 /* Avoid double counting of intra-cg interactions */
518 for (j = 1; j < j_end-j_start; j++)
520 nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
528 put_in_list_t (const gmx_bool bHaveVdW[],
536 const t_excl bExcl[],
544 put_in_list_at(const gmx_bool bHaveVdW[],
552 const t_excl bExcl[],
559 /* The a[] index has been removed,
560 * to put it back in i_atom should be a[i0] and jj should be a[jj].
565 t_nblist * vdwc_free = nullptr;
566 t_nblist * vdw_free = nullptr;
567 t_nblist * coul_free = nullptr;
568 t_nblist * vdwc_ww = nullptr;
569 t_nblist * coul_ww = nullptr;
571 int i, j, jcg, igid, gid, nbl_ind;
572 int jj, jj0, jj1, i_atom;
577 real *charge, *chargeB;
579 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
580 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
584 /* Copy some pointers */
586 charge = md->chargeA;
587 chargeB = md->chargeB;
590 bPert = md->bPerturbed;
594 nicg = index[icg+1]-i0;
596 /* Get the i charge group info */
597 igid = GET_CGINFO_GID(cginfo[icg]);
599 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
604 /* Check if any of the particles involved are perturbed.
605 * If not we can do the cheaper normal put_in_list
606 * and use more solvent optimization.
608 for (i = 0; i < nicg; i++)
610 bFreeEnergy |= bPert[i0+i];
612 /* Loop over the j charge groups */
613 for (j = 0; (j < nj && !bFreeEnergy); j++)
618 /* Finally loop over the atoms in the j-charge group */
619 for (jj = jj0; jj < jj1; jj++)
621 bFreeEnergy |= bPert[jj];
626 /* Unpack pointers to neighbourlist structs */
627 if (fr->nnblists == 1)
633 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
635 nlist = fr->nblists[nbl_ind].nlist_sr;
637 if (iwater != esolNO)
639 vdwc = &nlist[eNL_VDWQQ_WATER];
640 vdw = &nlist[eNL_VDW];
641 coul = &nlist[eNL_QQ_WATER];
642 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
643 coul_ww = &nlist[eNL_QQ_WATERWATER];
647 vdwc = &nlist[eNL_VDWQQ];
648 vdw = &nlist[eNL_VDW];
649 coul = &nlist[eNL_QQ];
654 if (iwater != esolNO)
656 /* Loop over the atoms in the i charge group */
658 gid = GID(igid, jgid, ngid);
659 /* Create new i_atom for each energy group */
660 if (bDoCoul && bDoVdW)
662 new_i_nblist(vdwc, i_atom, shift, gid);
663 new_i_nblist(vdwc_ww, i_atom, shift, gid);
667 new_i_nblist(vdw, i_atom, shift, gid);
671 new_i_nblist(coul, i_atom, shift, gid);
672 new_i_nblist(coul_ww, i_atom, shift, gid);
674 /* Loop over the j charge groups */
675 for (j = 0; (j < nj); j++)
685 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
687 if (iwater == esolSPC && jwater == esolSPC)
689 /* Interaction between two SPC molecules */
692 /* VdW only - only first atoms in each water interact */
693 add_j_to_nblist(vdw, jj0);
697 /* One entry for the entire water-water interaction */
700 add_j_to_nblist(coul_ww, jj0);
704 add_j_to_nblist(vdwc_ww, jj0);
708 else if (iwater == esolTIP4P && jwater == esolTIP4P)
710 /* Interaction between two TIP4p molecules */
713 /* VdW only - only first atoms in each water interact */
714 add_j_to_nblist(vdw, jj0);
718 /* One entry for the entire water-water interaction */
721 add_j_to_nblist(coul_ww, jj0);
725 add_j_to_nblist(vdwc_ww, jj0);
731 /* j charge group is not water, but i is.
732 * Add entries to the water-other_atom lists; the geometry of the water
733 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
734 * so we don't care if it is SPC or TIP4P...
741 for (jj = jj0; (jj < jj1); jj++)
745 add_j_to_nblist(coul, jj);
751 for (jj = jj0; (jj < jj1); jj++)
753 if (bHaveVdW[type[jj]])
755 add_j_to_nblist(vdw, jj);
761 /* _charge_ _groups_ interact with both coulomb and LJ */
762 /* Check which atoms we should add to the lists! */
763 for (jj = jj0; (jj < jj1); jj++)
765 if (bHaveVdW[type[jj]])
769 add_j_to_nblist(vdwc, jj);
773 add_j_to_nblist(vdw, jj);
776 else if (charge[jj] != 0)
778 add_j_to_nblist(coul, jj);
785 close_i_nblist(coul);
786 close_i_nblist(vdwc);
787 close_i_nblist(coul_ww);
788 close_i_nblist(vdwc_ww);
792 /* no solvent as i charge group */
793 /* Loop over the atoms in the i charge group */
794 for (i = 0; i < nicg; i++)
797 gid = GID(igid, jgid, ngid);
800 /* Create new i_atom for each energy group */
801 if (bDoVdW && bDoCoul)
803 new_i_nblist(vdwc, i_atom, shift, gid);
807 new_i_nblist(vdw, i_atom, shift, gid);
811 new_i_nblist(coul, i_atom, shift, gid);
813 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
814 bDoCoul_i = (bDoCoul && qi != 0);
816 if (bDoVdW_i || bDoCoul_i)
818 /* Loop over the j charge groups */
819 for (j = 0; (j < nj); j++)
823 /* Check for large charge groups */
834 /* Finally loop over the atoms in the j-charge group */
835 for (jj = jj0; jj < jj1; jj++)
837 bNotEx = NOTEXCL(bExcl, i, jj);
845 add_j_to_nblist(coul, jj);
850 if (bHaveVdW[type[jj]])
852 add_j_to_nblist(vdw, jj);
857 if (bHaveVdW[type[jj]])
861 add_j_to_nblist(vdwc, jj);
865 add_j_to_nblist(vdw, jj);
868 else if (charge[jj] != 0)
870 add_j_to_nblist(coul, jj);
878 close_i_nblist(coul);
879 close_i_nblist(vdwc);
885 /* we are doing free energy */
886 vdwc_free = &nlist[eNL_VDWQQ_FREE];
887 vdw_free = &nlist[eNL_VDW_FREE];
888 coul_free = &nlist[eNL_QQ_FREE];
889 /* Loop over the atoms in the i charge group */
890 for (i = 0; i < nicg; i++)
893 gid = GID(igid, jgid, ngid);
895 qiB = chargeB[i_atom];
897 /* Create new i_atom for each energy group */
898 if (bDoVdW && bDoCoul)
900 new_i_nblist(vdwc, i_atom, shift, gid);
904 new_i_nblist(vdw, i_atom, shift, gid);
908 new_i_nblist(coul, i_atom, shift, gid);
911 new_i_nblist(vdw_free, i_atom, shift, gid);
912 new_i_nblist(coul_free, i_atom, shift, gid);
913 new_i_nblist(vdwc_free, i_atom, shift, gid);
915 bDoVdW_i = (bDoVdW &&
916 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
917 bDoCoul_i = (bDoCoul && (qi != 0 || qiB != 0));
918 /* For TIP4P the first atom does not have a charge,
919 * but the last three do. So we should still put an atom
920 * without LJ but with charge in the water-atom neighborlist
921 * for a TIP4p i charge group.
922 * For SPC type water the first atom has LJ and charge,
923 * so there is no such problem.
925 if (iwater == esolNO)
927 bDoCoul_i_sol = bDoCoul_i;
931 bDoCoul_i_sol = bDoCoul;
934 if (bDoVdW_i || bDoCoul_i_sol)
936 /* Loop over the j charge groups */
937 for (j = 0; (j < nj); j++)
941 /* Check for large charge groups */
952 /* Finally loop over the atoms in the j-charge group */
953 bFree = bPert[i_atom];
954 for (jj = jj0; (jj < jj1); jj++)
956 bFreeJ = bFree || bPert[jj];
957 /* Complicated if, because the water H's should also
958 * see perturbed j-particles
960 if (iwater == esolNO || i == 0 || bFreeJ)
962 bNotEx = NOTEXCL(bExcl, i, jj);
970 if (charge[jj] != 0 || chargeB[jj] != 0)
972 add_j_to_nblist(coul_free, jj);
977 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
979 add_j_to_nblist(vdw_free, jj);
984 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
986 if (charge[jj] != 0 || chargeB[jj] != 0)
988 add_j_to_nblist(vdwc_free, jj);
992 add_j_to_nblist(vdw_free, jj);
995 else if (charge[jj] != 0 || chargeB[jj] != 0)
997 add_j_to_nblist(coul_free, jj);
1003 /* This is done whether or not bWater is set */
1004 if (charge[jj] != 0)
1006 add_j_to_nblist(coul, jj);
1009 else if (!bDoCoul_i_sol)
1011 if (bHaveVdW[type[jj]])
1013 add_j_to_nblist(vdw, jj);
1018 if (bHaveVdW[type[jj]])
1020 if (charge[jj] != 0)
1022 add_j_to_nblist(vdwc, jj);
1026 add_j_to_nblist(vdw, jj);
1029 else if (charge[jj] != 0)
1031 add_j_to_nblist(coul, jj);
1039 close_i_nblist(vdw);
1040 close_i_nblist(coul);
1041 close_i_nblist(vdwc);
1042 close_i_nblist(vdw_free);
1043 close_i_nblist(coul_free);
1044 close_i_nblist(vdwc_free);
1050 put_in_list_qmmm(const gmx_bool gmx_unused bHaveVdW[],
1052 const t_mdatoms * /* md */,
1058 const t_excl bExcl[],
1061 gmx_bool gmx_unused bDoVdW,
1062 gmx_bool gmx_unused bDoCoul,
1063 int gmx_unused solvent_opt)
1066 int i, j, jcg, igid, gid;
1067 int jj, jj0, jj1, i_atom;
1071 /* Get atom range */
1073 nicg = index[icg+1]-i0;
1075 /* Get the i charge group info */
1076 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1078 coul = fr->QMMMlist;
1080 /* Loop over atoms in the ith charge group */
1081 for (i = 0; i < nicg; i++)
1084 gid = GID(igid, jgid, ngid);
1085 /* Create new i_atom for each energy group */
1086 new_i_nblist(coul, i_atom, shift, gid);
1088 /* Loop over the j charge groups */
1089 for (j = 0; j < nj; j++)
1093 /* Charge groups cannot have QM and MM atoms simultaneously */
1098 /* Finally loop over the atoms in the j-charge group */
1099 for (jj = jj0; jj < jj1; jj++)
1101 bNotEx = NOTEXCL(bExcl, i, jj);
1104 add_j_to_nblist(coul, jj);
1109 close_i_nblist(coul);
1114 put_in_list_cg(const gmx_bool gmx_unused bHaveVdW[],
1116 const t_mdatoms * /* md */,
1122 const t_excl bExcl[],
1125 gmx_bool gmx_unused bDoVdW,
1126 gmx_bool gmx_unused bDoCoul,
1127 int gmx_unused solvent_opt)
1130 int igid, gid, nbl_ind;
1134 cginfo = fr->cginfo[icg];
1136 igid = GET_CGINFO_GID(cginfo);
1137 gid = GID(igid, jgid, ngid);
1139 /* Unpack pointers to neighbourlist structs */
1140 if (fr->nnblists == 1)
1146 nbl_ind = fr->gid2nblists[gid];
1148 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1150 /* Make a new neighbor list for charge group icg.
1151 * Currently simply one neighbor list is made with LJ and Coulomb.
1152 * If required, zero interactions could be removed here
1153 * or in the force loop.
1155 new_i_nblist(vdwc, index[icg], shift, gid);
1156 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1158 for (j = 0; (j < nj); j++)
1161 /* Skip the icg-icg pairs if all self interactions are excluded */
1162 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1164 /* Here we add the j charge group jcg to the list,
1165 * exclusions are also added to the list.
1167 add_j_to_nblist_cg(vdwc, index[jcg], index[jcg+1], bExcl, icg == jcg);
1171 close_i_nblist(vdwc);
1174 static void setexcl(int start, int end, t_blocka *excl, gmx_bool b,
1181 for (i = start; i < end; i++)
1183 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1185 SETEXCL(bexcl, i-start, excl->a[k]);
1191 for (i = start; i < end; i++)
1193 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1195 RMEXCL(bexcl, i-start, excl->a[k]);
1201 int calc_naaj(int icg, int cgtot)
1205 if ((cgtot % 2) == 1)
1207 /* Odd number of charge groups, easy */
1208 naaj = 1 + (cgtot/2);
1210 else if ((cgtot % 4) == 0)
1212 /* Multiple of four is hard */
1249 fprintf(log, "naaj=%d\n", naaj);
1255 /************************************************
1257 * S I M P L E C O R E S T U F F
1259 ************************************************/
1261 static real calc_image_tric(const rvec xi, const rvec xj, matrix box,
1262 const rvec b_inv, int *shift)
1264 /* This code assumes that the cut-off is smaller than
1265 * a half times the smallest diagonal element of the box.
1267 const real h25 = 2.5;
1272 /* Compute diff vector */
1273 dz = xj[ZZ] - xi[ZZ];
1274 dy = xj[YY] - xi[YY];
1275 dx = xj[XX] - xi[XX];
1277 /* Perform NINT operation, using trunc operation, therefore
1278 * we first add 2.5 then subtract 2 again
1280 tz = static_cast<int>(dz*b_inv[ZZ] + h25);
1282 dz -= tz*box[ZZ][ZZ];
1283 dy -= tz*box[ZZ][YY];
1284 dx -= tz*box[ZZ][XX];
1286 ty = static_cast<int>(dy*b_inv[YY] + h25);
1288 dy -= ty*box[YY][YY];
1289 dx -= ty*box[YY][XX];
1291 tx = static_cast<int>(dx*b_inv[XX]+h25);
1293 dx -= tx*box[XX][XX];
1295 /* Distance squared */
1296 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1298 *shift = XYZ2IS(tx, ty, tz);
1303 static real calc_image_rect(const rvec xi, const rvec xj, const rvec box_size,
1304 const rvec b_inv, int *shift)
1306 const real h15 = 1.5;
1312 /* Compute diff vector */
1313 dx = xj[XX] - xi[XX];
1314 dy = xj[YY] - xi[YY];
1315 dz = xj[ZZ] - xi[ZZ];
1317 /* Perform NINT operation, using trunc operation, therefore
1318 * we first add 1.5 then subtract 1 again
1320 tx = static_cast<int>(dx*b_inv[XX] + h15);
1321 ty = static_cast<int>(dy*b_inv[YY] + h15);
1322 tz = static_cast<int>(dz*b_inv[ZZ] + h15);
1327 /* Correct diff vector for translation */
1328 ddx = tx*box_size[XX] - dx;
1329 ddy = ty*box_size[YY] - dy;
1330 ddz = tz*box_size[ZZ] - dz;
1332 /* Distance squared */
1333 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1335 *shift = XYZ2IS(tx, ty, tz);
1340 static void add_simple(t_ns_buf * nsbuf,
1343 gmx_bool bHaveVdW[],
1345 const t_mdatoms *md,
1352 put_in_list_t *put_in_list)
1354 if (nsbuf->nj + nrj > MAX_CG)
1356 put_in_list(bHaveVdW, ngid, md, icg, jgid, nsbuf->ncg, nsbuf->jcg,
1357 cgs->index, bexcl, shift, fr, TRUE, TRUE, fr->solvent_opt);
1358 /* Reset buffer contents */
1359 nsbuf->ncg = nsbuf->nj = 0;
1361 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1365 static void ns_inner_tric(rvec x[],
1367 const int *i_egp_flags,
1375 gmx_bool bHaveVdW[],
1377 const t_mdatoms *md,
1380 put_in_list_t *put_in_list)
1384 int *cginfo = fr->cginfo;
1387 cgindex = cgs->index;
1389 for (j = 0; (j < njcg); j++)
1392 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1393 if (calc_image_tric(x[icg], x[cg_j], box, b_inv, &shift) < rcut2)
1395 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1396 if (!(i_egp_flags[jgid] & EGP_EXCL))
1398 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1399 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1406 static void ns_inner_rect(rvec x[],
1408 const int *i_egp_flags,
1417 gmx_bool bHaveVdW[],
1419 const t_mdatoms *md,
1422 put_in_list_t *put_in_list)
1426 int *cginfo = fr->cginfo;
1429 cgindex = cgs->index;
1433 for (j = 0; (j < njcg); j++)
1436 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1437 if (calc_image_rect(x[icg], x[cg_j], box_size, b_inv, &shift) < rcut2)
1439 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1440 if (!(i_egp_flags[jgid] & EGP_EXCL))
1442 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1443 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1451 for (j = 0; (j < njcg); j++)
1454 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1455 if ((rcut2 == 0) || (distance2(x[icg], x[cg_j]) < rcut2))
1457 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1458 if (!(i_egp_flags[jgid] & EGP_EXCL))
1460 add_simple(&ns_buf[jgid][CENTRAL], nrj, cg_j,
1461 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, CENTRAL, fr,
1469 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1471 static int ns_simple_core(t_forcerec *fr,
1472 gmx_localtop_t *top,
1473 const t_mdatoms *md,
1480 put_in_list_t *put_in_list,
1481 gmx_bool bHaveVdW[])
1485 int nsearch, icg, igid, nn;
1489 t_block *cgs = &(top->cgs);
1490 t_blocka *excl = &(top->excls);
1493 gmx_bool bBox, bTriclinic;
1496 rlist2 = gmx::square(fr->rlist);
1498 bBox = (fr->ePBC != epbcNONE);
1501 for (m = 0; (m < DIM); m++)
1503 if (gmx_numzero(box_size[m]))
1505 gmx_fatal(FARGS, "Dividing by zero box size!");
1507 b_inv[m] = 1.0/box_size[m];
1509 bTriclinic = TRICLINIC(box);
1516 cginfo = fr->cginfo;
1519 for (icg = fr->cg0; (icg < fr->hcg); icg++)
1522 i0 = cgs->index[icg];
1523 nri = cgs->index[icg+1]-i0;
1524 i_atoms = &(cgs->a[i0]);
1525 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1526 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1528 igid = GET_CGINFO_GID(cginfo[icg]);
1529 i_egp_flags = fr->egp_flags + ngid*igid;
1530 setexcl(cgs->index[icg], cgs->index[icg+1], excl, TRUE, bexcl);
1532 naaj = calc_naaj(icg, cgs->nr);
1535 ns_inner_tric(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1536 box, b_inv, rlist2, cgs, ns_buf,
1537 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1541 ns_inner_rect(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1542 bBox, box_size, b_inv, rlist2, cgs, ns_buf,
1543 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1547 for (nn = 0; (nn < ngid); nn++)
1549 for (k = 0; (k < SHIFTS); k++)
1551 nsbuf = &(ns_buf[nn][k]);
1554 put_in_list(bHaveVdW, ngid, md, icg, nn, nsbuf->ncg, nsbuf->jcg,
1555 cgs->index, bexcl, k, fr, TRUE, TRUE, fr->solvent_opt);
1556 nsbuf->ncg = nsbuf->nj = 0;
1560 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1561 setexcl(cgs->index[icg], cgs->index[icg+1], excl, FALSE, bexcl);
1563 close_neighbor_lists(fr, FALSE);
1568 /************************************************
1570 * N S 5 G R I D S T U F F
1572 ************************************************/
1574 static inline void get_dx_dd(int Nx, real gridx, real rc2, int xgi, real x,
1575 int ncpddc, int shift_min, int shift_max,
1576 int *g0, int *g1, real *dcx2)
1579 int g_min, g_max, shift_home;
1612 g_min = (shift_min == shift_home ? 0 : ncpddc);
1613 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1620 else if (shift_max < 0)
1635 /* Check one grid cell down */
1636 dcx = ((*g0 - 1) + 1)*gridx - x;
1648 /* Check one grid cell up */
1649 dcx = (*g1 + 1)*gridx - x;
1661 #define calc_dx2(XI, YI, ZI, y) (gmx::square((XI)-(y)[XX]) + gmx::square((YI)-(y)[YY]) + gmx::square((ZI)-(y)[ZZ]))
1662 #define calc_cyl_dx2(XI, YI, y) (gmx::square((XI)-(y)[XX]) + gmx::square((YI)-(y)[YY]))
1663 /****************************************************
1665 * F A S T N E I G H B O R S E A R C H I N G
1667 * Optimized neighboursearching routine using grid
1668 * at least 1x1x1, see GROMACS manual
1670 ****************************************************/
1673 static void get_cutoff2(t_forcerec *fr, real *rs2)
1675 *rs2 = gmx::square(fr->rlist);
1678 static void init_nsgrid_lists(t_forcerec *fr, int ngid, gmx_ns_t *ns)
1683 get_cutoff2(fr, &rs2);
1685 /* Short range buffers */
1686 snew(ns->nl_sr, ngid);
1688 snew(ns->nsr, ngid);
1690 for (j = 0; (j < ngid); j++)
1692 snew(ns->nl_sr[j], MAX_CG);
1697 "ns5_core: rs2 = %g (nm^2)\n",
1702 static int nsgrid_core(const t_commrec *cr,
1706 gmx_localtop_t *top,
1709 const gmx_bool *bExcludeAlleg,
1710 const t_mdatoms *md,
1711 put_in_list_t *put_in_list,
1712 gmx_bool bHaveVdW[],
1713 gmx_bool bMakeQMMMnblist)
1719 const t_block *cgs = &(top->cgs);
1720 int *cginfo = fr->cginfo;
1721 /* int *i_atoms,*cgsindex=cgs->index; */
1723 int cell_x, cell_y, cell_z;
1724 int d, tx, ty, tz, dx, dy, dz, cj;
1725 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
1726 int zsh_ty, zsh_tx, ysh_tx;
1728 int dx0, dx1, dy0, dy1, dz0, dz1;
1729 int Nx, Ny, Nz, shift = -1, j, nrj, nns, nn = -1;
1730 real gridx, gridy, gridz, grid_x, grid_y;
1731 real *dcx2, *dcy2, *dcz2;
1733 int cg0, cg1, icg = -1, cgsnr, i0, igid, naaj, max_jcg;
1734 int jcg0, jcg1, jjcg, cgj0, jgid;
1735 int *grida, *gridnra, *gridind;
1736 rvec *cgcm, grid_offset;
1737 real r2, rs2, XI, YI, ZI, tmp1, tmp2;
1739 gmx_bool bDomDec, bTriclinicX, bTriclinicY;
1744 bDomDec = DOMAINDECOMP(cr);
1747 bTriclinicX = ((YY < grid->npbcdim &&
1748 (!bDomDec || dd->nc[YY] == 1) && box[YY][XX] != 0) ||
1749 (ZZ < grid->npbcdim &&
1750 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][XX] != 0));
1751 bTriclinicY = (ZZ < grid->npbcdim &&
1752 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][YY] != 0);
1756 get_cutoff2(fr, &rs2);
1767 gridind = grid->index;
1768 gridnra = grid->nra;
1771 gridx = grid->cell_size[XX];
1772 gridy = grid->cell_size[YY];
1773 gridz = grid->cell_size[ZZ];
1776 copy_rvec(grid->cell_offset, grid_offset);
1777 copy_ivec(grid->ncpddc, ncpddc);
1782 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
1783 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
1784 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
1785 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
1786 if (zsh_tx != 0 && ysh_tx != 0)
1788 /* This could happen due to rounding, when both ratios are 0.5 */
1795 /* We only want a list for the test particle */
1804 /* Set the shift range */
1805 for (d = 0; d < DIM; d++)
1809 /* Check if we need periodicity shifts.
1810 * Without PBC or with domain decomposition we don't need them.
1812 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
1819 box[XX][XX] - std::abs(box[YY][XX]) - std::abs(box[ZZ][XX]) < std::sqrt(rs2))
1830 /* Loop over charge groups */
1831 for (icg = cg0; (icg < cg1); icg++)
1833 igid = GET_CGINFO_GID(cginfo[icg]);
1834 /* Skip this charge group if all energy groups are excluded! */
1835 if (bExcludeAlleg[igid])
1840 i0 = cgs->index[icg];
1842 if (bMakeQMMMnblist)
1844 /* Skip this charge group if it is not a QM atom while making a
1845 * QM/MM neighbourlist
1847 if (md->bQM[i0] == FALSE)
1849 continue; /* MM particle, go to next particle */
1852 /* Compute the number of charge groups that fall within the control
1855 naaj = calc_naaj(icg, cgsnr);
1862 /* make a normal neighbourlist */
1866 /* Get the j charge-group and dd cell shift ranges */
1867 dd_get_ns_ranges(cr->dd, icg, &jcg0, &jcg1, sh0, sh1);
1872 /* Compute the number of charge groups that fall within the control
1875 naaj = calc_naaj(icg, cgsnr);
1881 /* The i-particle is awlways the test particle,
1882 * so we want all j-particles
1884 max_jcg = cgsnr - 1;
1888 max_jcg = jcg1 - cgsnr;
1893 i_egp_flags = fr->egp_flags + igid*ngid;
1895 /* Set the exclusions for the atoms in charge group icg using a bitmask */
1896 setexcl(i0, cgs->index[icg+1], &top->excls, TRUE, bexcl);
1898 ci2xyz(grid, icg, &cell_x, &cell_y, &cell_z);
1900 /* Changed iicg to icg, DvdS 990115
1901 * (but see consistency check above, DvdS 990330)
1904 fprintf(log, "icg=%5d, naaj=%5d, cell %d %d %d\n",
1905 icg, naaj, cell_x, cell_y, cell_z);
1907 /* Loop over shift vectors in three dimensions */
1908 for (tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
1910 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
1911 /* Calculate range of cells in Z direction that have the shift tz */
1912 zgi = cell_z + tz*Nz;
1913 get_dx_dd(Nz, gridz, rs2, zgi, ZI-grid_offset[ZZ],
1914 ncpddc[ZZ], sh0[ZZ], sh1[ZZ], &dz0, &dz1, dcz2);
1919 for (ty = -shp[YY]; ty <= shp[YY]; ty++)
1921 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
1922 /* Calculate range of cells in Y direction that have the shift ty */
1925 ygi = static_cast<int>(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
1929 ygi = cell_y + ty*Ny;
1931 get_dx_dd(Ny, gridy, rs2, ygi, YI-grid_offset[YY],
1932 ncpddc[YY], sh0[YY], sh1[YY], &dy0, &dy1, dcy2);
1937 for (tx = -shp[XX]; tx <= shp[XX]; tx++)
1939 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
1940 /* Calculate range of cells in X direction that have the shift tx */
1943 xgi = static_cast<int>(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
1947 xgi = cell_x + tx*Nx;
1949 get_dx_dd(Nx, gridx, rs2, xgi, XI-grid_offset[XX],
1950 ncpddc[XX], sh0[XX], sh1[XX], &dx0, &dx1, dcx2);
1955 /* Get shift vector */
1956 shift = XYZ2IS(tx, ty, tz);
1958 range_check(shift, 0, SHIFTS);
1960 for (nn = 0; (nn < ngid); nn++)
1965 fprintf(log, "shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
1966 shift, dx0, dx1, dy0, dy1, dz0, dz1);
1967 fprintf(log, "cgcm: %8.3f %8.3f %8.3f\n", cgcm[icg][XX],
1968 cgcm[icg][YY], cgcm[icg][ZZ]);
1969 fprintf(log, "xi: %8.3f %8.3f %8.3f\n", XI, YI, ZI);
1971 for (dx = dx0; (dx <= dx1); dx++)
1973 tmp1 = rs2 - dcx2[dx];
1974 for (dy = dy0; (dy <= dy1); dy++)
1976 tmp2 = tmp1 - dcy2[dy];
1979 for (dz = dz0; (dz <= dz1); dz++)
1981 if (tmp2 > dcz2[dz])
1983 /* Find grid-cell cj in which possible neighbours are */
1984 cj = xyz2ci(Ny, Nz, dx, dy, dz);
1986 /* Check out how many cgs (nrj) there in this cell */
1989 /* Find the offset in the cg list */
1992 /* Check if all j's are out of range so we
1993 * can skip the whole cell.
1994 * Should save some time, especially with DD.
1997 (grida[cgj0] >= max_jcg &&
1998 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2004 for (j = 0; (j < nrj); j++)
2006 jjcg = grida[cgj0+j];
2008 /* check whether this guy is in range! */
2009 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2012 r2 = calc_dx2(XI, YI, ZI, cgcm[jjcg]);
2015 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2016 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2017 /* check energy group exclusions */
2018 if (!(i_egp_flags[jgid] & EGP_EXCL))
2020 if (nsr[jgid] >= MAX_CG)
2022 /* Add to short-range list */
2023 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2024 nsr[jgid], nl_sr[jgid],
2025 cgs->index, /* cgsatoms, */ bexcl,
2026 shift, fr, TRUE, TRUE, fr->solvent_opt);
2029 nl_sr[jgid][nsr[jgid]++] = jjcg;
2040 /* CHECK whether there is anything left in the buffers */
2041 for (nn = 0; (nn < ngid); nn++)
2045 put_in_list(bHaveVdW, ngid, md, icg, nn, nsr[nn], nl_sr[nn],
2046 cgs->index, /* cgsatoms, */ bexcl,
2047 shift, fr, TRUE, TRUE, fr->solvent_opt);
2053 setexcl(cgs->index[icg], cgs->index[icg+1], &top->excls, FALSE, bexcl);
2055 /* No need to perform any left-over force calculations anymore (as we used to do here)
2056 * since we now save the proper long-range lists for later evaluation.
2059 /* Close neighbourlists */
2060 close_neighbor_lists(fr, bMakeQMMMnblist);
2065 static void ns_realloc_natoms(gmx_ns_t *ns, int natoms)
2069 if (natoms > ns->nra_alloc)
2071 ns->nra_alloc = over_alloc_dd(natoms);
2072 srenew(ns->bexcl, ns->nra_alloc);
2073 for (i = 0; i < ns->nra_alloc; i++)
2080 void init_ns(FILE *fplog, const t_commrec *cr,
2081 gmx_ns_t *ns, t_forcerec *fr,
2082 const gmx_mtop_t *mtop)
2084 int icg, nr_in_cg, maxcg, i, j, jcg, ngid, ncg;
2086 /* Compute largest charge groups size (# atoms) */
2088 for (const gmx_moltype_t &molt : mtop->moltype)
2090 const t_block *cgs = &molt.cgs;
2091 for (icg = 0; (icg < cgs->nr); icg++)
2093 nr_in_cg = std::max(nr_in_cg, (cgs->index[icg+1]-cgs->index[icg]));
2097 /* Verify whether largest charge group is <= max cg.
2098 * This is determined by the type of the local exclusion type
2099 * Exclusions are stored in bits. (If the type is not large
2100 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2102 maxcg = sizeof(t_excl)*8;
2103 if (nr_in_cg > maxcg)
2105 gmx_fatal(FARGS, "Max #atoms in a charge group: %d > %d\n",
2109 ngid = mtop->groups.grps[egcENER].nr;
2110 snew(ns->bExcludeAlleg, ngid);
2111 for (i = 0; i < ngid; i++)
2113 ns->bExcludeAlleg[i] = TRUE;
2114 for (j = 0; j < ngid; j++)
2116 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2118 ns->bExcludeAlleg[i] = FALSE;
2126 ns->grid = init_grid(fplog, fr);
2127 init_nsgrid_lists(fr, ngid, ns);
2132 snew(ns->ns_buf, ngid);
2133 for (i = 0; (i < ngid); i++)
2135 snew(ns->ns_buf[i], SHIFTS);
2137 ncg = ncg_mtop(mtop);
2138 snew(ns->simple_aaj, 2*ncg);
2139 for (jcg = 0; (jcg < ncg); jcg++)
2141 ns->simple_aaj[jcg] = jcg;
2142 ns->simple_aaj[jcg+ncg] = jcg;
2146 /* Create array that determines whether or not atoms have VdW */
2147 snew(ns->bHaveVdW, fr->ntype);
2148 for (i = 0; (i < fr->ntype); i++)
2150 for (j = 0; (j < fr->ntype); j++)
2152 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2154 ((BHAMA(fr->nbfp, fr->ntype, i, j) != 0) ||
2155 (BHAMB(fr->nbfp, fr->ntype, i, j) != 0) ||
2156 (BHAMC(fr->nbfp, fr->ntype, i, j) != 0)) :
2157 ((C6(fr->nbfp, fr->ntype, i, j) != 0) ||
2158 (C12(fr->nbfp, fr->ntype, i, j) != 0))));
2163 pr_bvec(debug, 0, "bHaveVdW", ns->bHaveVdW, fr->ntype, TRUE);
2167 ns->bexcl = nullptr;
2168 if (!DOMAINDECOMP(cr))
2170 ns_realloc_natoms(ns, mtop->natoms);
2173 ns->nblist_initialized = FALSE;
2175 /* nbr list debug dump */
2177 char *ptr = getenv("GMX_DUMP_NL");
2180 ns->dump_nl = strtol(ptr, nullptr, 10);
2183 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2193 void done_ns(gmx_ns_t *ns, int numEnergyGroups)
2195 sfree(ns->bExcludeAlleg);
2198 for (int i = 0; i < numEnergyGroups; i++)
2200 sfree(ns->ns_buf[i]);
2204 sfree(ns->simple_aaj);
2205 sfree(ns->bHaveVdW);
2206 done_grid(ns->grid);
2210 int search_neighbours(FILE *log,
2213 gmx_localtop_t *top,
2214 const gmx_groups_t *groups,
2215 const t_commrec *cr,
2217 const t_mdatoms *md,
2220 const t_block *cgs = &(top->cgs);
2221 rvec box_size, grid_x0, grid_x1;
2223 real min_size, grid_dens;
2229 gmx_domdec_zones_t *dd_zones;
2230 put_in_list_t *put_in_list;
2234 /* Set some local variables */
2236 ngid = groups->grps[egcENER].nr;
2238 for (m = 0; (m < DIM); m++)
2240 box_size[m] = box[m][m];
2243 if (fr->ePBC != epbcNONE)
2245 if (gmx::square(fr->rlist) >= max_cutoff2(fr->ePBC, box))
2247 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.");
2251 min_size = std::min(box_size[XX], std::min(box_size[YY], box_size[ZZ]));
2252 if (2*fr->rlist >= min_size)
2254 gmx_fatal(FARGS, "One of the box diagonal elements has become smaller than twice the cut-off length.");
2259 if (DOMAINDECOMP(cr))
2261 ns_realloc_natoms(ns, cgs->index[cgs->nr]);
2264 /* Reset the neighbourlists */
2265 reset_neighbor_lists(fr);
2267 if (bGrid && bFillGrid)
2271 if (DOMAINDECOMP(cr))
2273 dd_zones = domdec_zones(cr->dd);
2279 get_nsgrid_boundaries(grid->nboundeddim, box, nullptr, nullptr, nullptr, nullptr,
2280 cgs->nr, fr->cg_cm, grid_x0, grid_x1, &grid_dens);
2282 grid_first(log, grid, nullptr, nullptr, box, grid_x0, grid_x1,
2283 fr->rlist, grid_dens);
2289 if (DOMAINDECOMP(cr))
2292 fill_grid(dd_zones, grid, end, -1, end, fr->cg_cm);
2294 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2298 fill_grid(nullptr, grid, cgs->nr, fr->cg0, fr->hcg, fr->cg_cm);
2299 grid->icg0 = fr->cg0;
2300 grid->icg1 = fr->hcg;
2303 calc_elemnr(grid, start, end, cgs->nr);
2305 grid_last(grid, start, end, cgs->nr);
2310 print_grid(debug, grid);
2315 /* Set the grid cell index for the test particle only.
2316 * The cell to cg index is not corrected, but that does not matter.
2318 fill_grid(nullptr, ns->grid, fr->hcg, fr->hcg-1, fr->hcg, fr->cg_cm);
2321 if (!fr->ns->bCGlist)
2323 put_in_list = put_in_list_at;
2327 put_in_list = put_in_list_cg;
2334 nsearch = nsgrid_core(cr, fr, box, ngid, top,
2335 grid, ns->bexcl, ns->bExcludeAlleg,
2336 md, put_in_list, ns->bHaveVdW,
2339 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2340 * the classical calculation. The charge-charge interaction
2341 * between QM and MM atoms is handled in the QMMM core calculation
2342 * (see QMMM.c). The VDW however, we'd like to compute classically
2343 * and the QM MM atom pairs have just been put in the
2344 * corresponding neighbourlists. in case of QMMM we still need to
2345 * fill a special QMMM neighbourlist that contains all neighbours
2346 * of the QM atoms. If bQMMM is true, this list will now be made:
2348 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
2350 nsearch += nsgrid_core(cr, fr, box, ngid, top,
2351 grid, ns->bexcl, ns->bExcludeAlleg,
2352 md, put_in_list_qmmm, ns->bHaveVdW,
2358 nsearch = ns_simple_core(fr, top, md, box, box_size,
2359 ns->bexcl, ns->simple_aaj,
2360 ngid, ns->ns_buf, put_in_list, ns->bHaveVdW);
2363 inc_nrnb(nrnb, eNR_NS, nsearch);