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/utility/smalloc.h"
47 #include "gromacs/math/utilities.h"
49 #include "types/commrec.h"
53 #include "nonbonded.h"
57 #include "gromacs/utility/fatalerror.h"
60 #include "mtop_util.h"
67 * E X C L U S I O N H A N D L I N G
71 static void SETEXCL_(t_excl e[], atom_id i, atom_id j)
75 static void RMEXCL_(t_excl e[], atom_id i, atom_id j)
77 e[j] = e[j] & ~(1<<i);
79 static gmx_bool ISEXCL_(t_excl e[], atom_id i, atom_id j)
81 return (gmx_bool)(e[j] & (1<<i));
83 static gmx_bool NOTEXCL_(t_excl e[], atom_id i, atom_id j)
85 return !(ISEXCL(e, i, j));
88 #define SETEXCL(e, i, j) (e)[((atom_id) (j))] |= (1<<((atom_id) (i)))
89 #define RMEXCL(e, i, j) (e)[((atom_id) (j))] &= (~(1<<((atom_id) (i))))
90 #define ISEXCL(e, i, j) (gmx_bool) ((e)[((atom_id) (j))] & (1<<((atom_id) (i))))
91 #define NOTEXCL(e, i, j) !(ISEXCL(e, i, j))
95 round_up_to_simd_width(int length, int simd_width)
97 int offset, newlength;
99 offset = (simd_width > 0) ? length % simd_width : 0;
101 return (offset == 0) ? length : length-offset+simd_width;
103 /************************************************
105 * U T I L I T I E S F O R N S
107 ************************************************/
109 void reallocate_nblist(t_nblist *nl)
113 fprintf(debug, "reallocating neigborlist (ielec=%d, ivdw=%d, igeometry=%d, type=%d), maxnri=%d\n",
114 nl->ielec, nl->ivdw, nl->igeometry, nl->type, nl->maxnri);
116 srenew(nl->iinr, nl->maxnri);
117 if (nl->igeometry == GMX_NBLIST_GEOMETRY_CG_CG)
119 srenew(nl->iinr_end, nl->maxnri);
121 srenew(nl->gid, nl->maxnri);
122 srenew(nl->shift, nl->maxnri);
123 srenew(nl->jindex, nl->maxnri+1);
127 static void init_nblist(FILE *log, t_nblist *nl_sr, t_nblist *nl_lr,
128 int maxsr, int maxlr,
129 int ivdw, int ivdwmod,
130 int ielec, int ielecmod,
131 int igeometry, int type)
137 for (i = 0; (i < 2); i++)
139 nl = (i == 0) ? nl_sr : nl_lr;
140 homenr = (i == 0) ? maxsr : maxlr;
148 /* Set coul/vdw in neighborlist, and for the normal loops we determine
149 * an index of which one to call.
152 nl->ivdwmod = ivdwmod;
154 nl->ielecmod = ielecmod;
156 nl->igeometry = igeometry;
158 if (nl->type == GMX_NBLIST_INTERACTION_FREE_ENERGY)
160 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
163 /* This will also set the simd_padding_width field */
164 gmx_nonbonded_set_kernel_pointers( (i == 0) ? log : NULL, nl);
166 /* maxnri is influenced by the number of shifts (maximum is 8)
167 * and the number of energy groups.
168 * If it is not enough, nl memory will be reallocated during the run.
169 * 4 seems to be a reasonable factor, which only causes reallocation
170 * during runs with tiny and many energygroups.
172 nl->maxnri = homenr*4;
182 reallocate_nblist(nl);
187 fprintf(debug, "Initiating neighbourlist (ielec=%d, ivdw=%d, type=%d) for %s interactions,\nwith %d SR, %d LR atoms.\n",
188 nl->ielec, nl->ivdw, nl->type, gmx_nblist_geometry_names[nl->igeometry], maxsr, maxlr);
193 void init_neighbor_list(FILE *log, t_forcerec *fr, int homenr)
195 /* Make maxlr tunable! (does not seem to be a big difference though)
196 * This parameter determines the number of i particles in a long range
197 * neighbourlist. Too few means many function calls, too many means
200 int maxsr, maxsr_wat, maxlr, maxlr_wat;
201 int ielec, ielecf, ivdw, ielecmod, ielecmodf, ivdwmod, type;
203 int igeometry_def, igeometry_w, igeometry_ww;
207 /* maxsr = homenr-fr->nWatMol*3; */
212 gmx_fatal(FARGS, "%s, %d: Negative number of short range atoms.\n"
213 "Call your Gromacs dealer for assistance.", __FILE__, __LINE__);
215 /* This is just for initial allocation, so we do not reallocate
216 * all the nlist arrays many times in a row.
217 * The numbers seem very accurate, but they are uncritical.
219 maxsr_wat = min(fr->nWatMol, (homenr+2)/3);
223 maxlr_wat = min(maxsr_wat, maxlr);
227 maxlr = maxlr_wat = 0;
230 /* Determine the values for ielec/ivdw. */
231 ielec = fr->nbkernel_elec_interaction;
232 ivdw = fr->nbkernel_vdw_interaction;
233 ielecmod = fr->nbkernel_elec_modifier;
234 ivdwmod = fr->nbkernel_vdw_modifier;
235 type = GMX_NBLIST_INTERACTION_STANDARD;
237 fr->ns.bCGlist = (getenv("GMX_NBLISTCG") != 0);
240 igeometry_def = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
244 igeometry_def = GMX_NBLIST_GEOMETRY_CG_CG;
247 fprintf(log, "\nUsing charge-group - charge-group neighbor lists and kernels\n\n");
251 if (fr->solvent_opt == esolTIP4P)
253 igeometry_w = GMX_NBLIST_GEOMETRY_WATER4_PARTICLE;
254 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER4_WATER4;
258 igeometry_w = GMX_NBLIST_GEOMETRY_WATER3_PARTICLE;
259 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER3_WATER3;
262 for (i = 0; i < fr->nnblists; i++)
264 nbl = &(fr->nblists[i]);
266 if ((fr->adress_type != eAdressOff) && (i >= fr->nnblists/2))
268 type = GMX_NBLIST_INTERACTION_ADRESS;
270 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ], &nbl->nlist_lr[eNL_VDWQQ],
271 maxsr, maxlr, ivdw, ivdwmod, ielec, ielecmod, igeometry_def, type);
272 init_nblist(log, &nbl->nlist_sr[eNL_VDW], &nbl->nlist_lr[eNL_VDW],
273 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, igeometry_def, type);
274 init_nblist(log, &nbl->nlist_sr[eNL_QQ], &nbl->nlist_lr[eNL_QQ],
275 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_def, type);
276 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATER], &nbl->nlist_lr[eNL_VDWQQ_WATER],
277 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_w, type);
278 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATER], &nbl->nlist_lr[eNL_QQ_WATER],
279 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_w, type);
280 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATERWATER], &nbl->nlist_lr[eNL_VDWQQ_WATERWATER],
281 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_ww, type);
282 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATERWATER], &nbl->nlist_lr[eNL_QQ_WATERWATER],
283 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_ww, type);
285 /* Did we get the solvent loops so we can use optimized water kernels? */
286 if (nbl->nlist_sr[eNL_VDWQQ_WATER].kernelptr_vf == NULL
287 || nbl->nlist_sr[eNL_QQ_WATER].kernelptr_vf == NULL
288 #ifndef DISABLE_WATERWATER_NLIST
289 || nbl->nlist_sr[eNL_VDWQQ_WATERWATER].kernelptr_vf == NULL
290 || nbl->nlist_sr[eNL_QQ_WATERWATER].kernelptr_vf == NULL
294 fr->solvent_opt = esolNO;
295 fprintf(log, "Note: The available nonbonded kernels do not support water optimization - disabling.\n");
298 if (fr->efep != efepNO)
300 if ((fr->bEwald) && (fr->sc_alphacoul > 0)) /* need to handle long range differently if using softcore */
302 ielecf = GMX_NBKERNEL_ELEC_EWALD;
303 ielecmodf = eintmodNONE;
308 ielecmodf = ielecmod;
311 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_FREE], &nbl->nlist_lr[eNL_VDWQQ_FREE],
312 maxsr, maxlr, ivdw, ivdwmod, ielecf, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
313 init_nblist(log, &nbl->nlist_sr[eNL_VDW_FREE], &nbl->nlist_lr[eNL_VDW_FREE],
314 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
315 init_nblist(log, &nbl->nlist_sr[eNL_QQ_FREE], &nbl->nlist_lr[eNL_QQ_FREE],
316 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielecf, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
320 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
322 init_nblist(log, &fr->QMMMlist, NULL,
323 maxsr, maxlr, 0, 0, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_STANDARD);
331 fr->ns.nblist_initialized = TRUE;
334 static void reset_nblist(t_nblist *nl)
344 static void reset_neighbor_lists(t_forcerec *fr, gmx_bool bResetSR, gmx_bool bResetLR)
350 /* only reset the short-range nblist */
351 reset_nblist(&(fr->QMMMlist));
354 for (n = 0; n < fr->nnblists; n++)
356 for (i = 0; i < eNL_NR; i++)
360 reset_nblist( &(fr->nblists[n].nlist_sr[i]) );
364 reset_nblist( &(fr->nblists[n].nlist_lr[i]) );
373 static gmx_inline void new_i_nblist(t_nblist *nlist, atom_id i_atom, int shift, int gid)
375 int i, k, nri, nshift;
379 /* Check whether we have to increase the i counter */
381 (nlist->iinr[nri] != i_atom) ||
382 (nlist->shift[nri] != shift) ||
383 (nlist->gid[nri] != gid))
385 /* This is something else. Now see if any entries have
386 * been added in the list of the previous atom.
389 ((nlist->jindex[nri+1] > nlist->jindex[nri]) &&
390 (nlist->gid[nri] != -1)))
392 /* If so increase the counter */
395 if (nlist->nri >= nlist->maxnri)
397 nlist->maxnri += over_alloc_large(nlist->nri);
398 reallocate_nblist(nlist);
401 /* Set the number of neighbours and the atom number */
402 nlist->jindex[nri+1] = nlist->jindex[nri];
403 nlist->iinr[nri] = i_atom;
404 nlist->gid[nri] = gid;
405 nlist->shift[nri] = shift;
409 /* Adding to previous list. First remove possible previous padding */
410 if (nlist->simd_padding_width > 1)
412 while (nlist->nrj > 0 && nlist->jjnr[nlist->nrj-1] < 0)
420 static gmx_inline void close_i_nblist(t_nblist *nlist)
422 int nri = nlist->nri;
427 /* Add elements up to padding. Since we allocate memory in units
428 * of the simd_padding width, we do not have to check for possible
429 * list reallocation here.
431 while ((nlist->nrj % nlist->simd_padding_width) != 0)
433 /* Use -4 here, so we can write forces for 4 atoms before real data */
434 nlist->jjnr[nlist->nrj++] = -4;
436 nlist->jindex[nri+1] = nlist->nrj;
438 len = nlist->nrj - nlist->jindex[nri];
442 static gmx_inline void close_nblist(t_nblist *nlist)
444 /* Only close this nblist when it has been initialized.
445 * Avoid the creation of i-lists with no j-particles.
449 /* Some assembly kernels do not support empty lists,
450 * make sure here that we don't generate any empty lists.
451 * With the current ns code this branch is taken in two cases:
452 * No i-particles at all: nri=-1 here
453 * There are i-particles, but no j-particles; nri=0 here
459 /* Close list number nri by incrementing the count */
464 static gmx_inline void close_neighbor_lists(t_forcerec *fr, gmx_bool bMakeQMMMnblist)
470 close_nblist(&(fr->QMMMlist));
473 for (n = 0; n < fr->nnblists; n++)
475 for (i = 0; (i < eNL_NR); i++)
477 close_nblist(&(fr->nblists[n].nlist_sr[i]));
478 close_nblist(&(fr->nblists[n].nlist_lr[i]));
484 static gmx_inline void add_j_to_nblist(t_nblist *nlist, atom_id j_atom, gmx_bool bLR)
486 int nrj = nlist->nrj;
488 if (nlist->nrj >= nlist->maxnrj)
490 nlist->maxnrj = round_up_to_simd_width(over_alloc_small(nlist->nrj + 1), nlist->simd_padding_width);
494 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
495 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
498 srenew(nlist->jjnr, nlist->maxnrj);
501 nlist->jjnr[nrj] = j_atom;
505 static gmx_inline void add_j_to_nblist_cg(t_nblist *nlist,
506 atom_id j_start, int j_end,
507 t_excl *bexcl, gmx_bool i_is_j,
510 int nrj = nlist->nrj;
513 if (nlist->nrj >= nlist->maxnrj)
515 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
518 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
519 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
522 srenew(nlist->jjnr, nlist->maxnrj);
523 srenew(nlist->jjnr_end, nlist->maxnrj);
524 srenew(nlist->excl, nlist->maxnrj*MAX_CGCGSIZE);
527 nlist->jjnr[nrj] = j_start;
528 nlist->jjnr_end[nrj] = j_end;
530 if (j_end - j_start > MAX_CGCGSIZE)
532 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);
535 /* Set the exclusions */
536 for (j = j_start; j < j_end; j++)
538 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
542 /* Avoid double counting of intra-cg interactions */
543 for (j = 1; j < j_end-j_start; j++)
545 nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
553 put_in_list_t (gmx_bool bHaveVdW[],
570 put_in_list_at(gmx_bool bHaveVdW[],
586 /* The a[] index has been removed,
587 * to put it back in i_atom should be a[i0] and jj should be a[jj].
592 t_nblist * vdwc_free = NULL;
593 t_nblist * vdw_free = NULL;
594 t_nblist * coul_free = NULL;
595 t_nblist * vdwc_ww = NULL;
596 t_nblist * coul_ww = NULL;
598 int i, j, jcg, igid, gid, nbl_ind, ind_ij;
599 atom_id jj, jj0, jj1, i_atom;
604 real *charge, *chargeB;
605 real qi, qiB, qq, rlj;
606 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
607 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
611 /* Copy some pointers */
613 charge = md->chargeA;
614 chargeB = md->chargeB;
617 bPert = md->bPerturbed;
621 nicg = index[icg+1]-i0;
623 /* Get the i charge group info */
624 igid = GET_CGINFO_GID(cginfo[icg]);
626 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
631 /* Check if any of the particles involved are perturbed.
632 * If not we can do the cheaper normal put_in_list
633 * and use more solvent optimization.
635 for (i = 0; i < nicg; i++)
637 bFreeEnergy |= bPert[i0+i];
639 /* Loop over the j charge groups */
640 for (j = 0; (j < nj && !bFreeEnergy); j++)
645 /* Finally loop over the atoms in the j-charge group */
646 for (jj = jj0; jj < jj1; jj++)
648 bFreeEnergy |= bPert[jj];
653 /* Unpack pointers to neighbourlist structs */
654 if (fr->nnblists == 1)
660 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
664 nlist = fr->nblists[nbl_ind].nlist_lr;
668 nlist = fr->nblists[nbl_ind].nlist_sr;
671 if (iwater != esolNO)
673 vdwc = &nlist[eNL_VDWQQ_WATER];
674 vdw = &nlist[eNL_VDW];
675 coul = &nlist[eNL_QQ_WATER];
676 #ifndef DISABLE_WATERWATER_NLIST
677 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
678 coul_ww = &nlist[eNL_QQ_WATERWATER];
683 vdwc = &nlist[eNL_VDWQQ];
684 vdw = &nlist[eNL_VDW];
685 coul = &nlist[eNL_QQ];
690 if (iwater != esolNO)
692 /* Loop over the atoms in the i charge group */
694 gid = GID(igid, jgid, ngid);
695 /* Create new i_atom for each energy group */
696 if (bDoCoul && bDoVdW)
698 new_i_nblist(vdwc, i_atom, shift, gid);
699 #ifndef DISABLE_WATERWATER_NLIST
700 new_i_nblist(vdwc_ww, i_atom, shift, gid);
705 new_i_nblist(vdw, i_atom, shift, gid);
709 new_i_nblist(coul, i_atom, shift, gid);
710 #ifndef DISABLE_WATERWATER_NLIST
711 new_i_nblist(coul_ww, i_atom, shift, gid);
714 /* Loop over the j charge groups */
715 for (j = 0; (j < nj); j++)
725 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
727 if (iwater == esolSPC && jwater == esolSPC)
729 /* Interaction between two SPC molecules */
732 /* VdW only - only first atoms in each water interact */
733 add_j_to_nblist(vdw, jj0, bLR);
737 #ifdef DISABLE_WATERWATER_NLIST
738 /* Add entries for the three atoms - only do VdW if we need to */
741 add_j_to_nblist(coul, jj0, bLR);
745 add_j_to_nblist(vdwc, jj0, bLR);
747 add_j_to_nblist(coul, jj0+1, bLR);
748 add_j_to_nblist(coul, jj0+2, bLR);
750 /* One entry for the entire water-water interaction */
753 add_j_to_nblist(coul_ww, jj0, bLR);
757 add_j_to_nblist(vdwc_ww, jj0, bLR);
762 else if (iwater == esolTIP4P && jwater == esolTIP4P)
764 /* Interaction between two TIP4p molecules */
767 /* VdW only - only first atoms in each water interact */
768 add_j_to_nblist(vdw, jj0, bLR);
772 #ifdef DISABLE_WATERWATER_NLIST
773 /* Add entries for the four atoms - only do VdW if we need to */
776 add_j_to_nblist(vdw, jj0, bLR);
778 add_j_to_nblist(coul, jj0+1, bLR);
779 add_j_to_nblist(coul, jj0+2, bLR);
780 add_j_to_nblist(coul, jj0+3, bLR);
782 /* One entry for the entire water-water interaction */
785 add_j_to_nblist(coul_ww, jj0, bLR);
789 add_j_to_nblist(vdwc_ww, jj0, bLR);
796 /* j charge group is not water, but i is.
797 * Add entries to the water-other_atom lists; the geometry of the water
798 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
799 * so we don't care if it is SPC or TIP4P...
806 for (jj = jj0; (jj < jj1); jj++)
810 add_j_to_nblist(coul, jj, bLR);
816 for (jj = jj0; (jj < jj1); jj++)
818 if (bHaveVdW[type[jj]])
820 add_j_to_nblist(vdw, jj, bLR);
826 /* _charge_ _groups_ interact with both coulomb and LJ */
827 /* Check which atoms we should add to the lists! */
828 for (jj = jj0; (jj < jj1); jj++)
830 if (bHaveVdW[type[jj]])
834 add_j_to_nblist(vdwc, jj, bLR);
838 add_j_to_nblist(vdw, jj, bLR);
841 else if (charge[jj] != 0)
843 add_j_to_nblist(coul, jj, bLR);
850 close_i_nblist(coul);
851 close_i_nblist(vdwc);
852 #ifndef DISABLE_WATERWATER_NLIST
853 close_i_nblist(coul_ww);
854 close_i_nblist(vdwc_ww);
859 /* no solvent as i charge group */
860 /* Loop over the atoms in the i charge group */
861 for (i = 0; i < nicg; i++)
864 gid = GID(igid, jgid, ngid);
867 /* Create new i_atom for each energy group */
868 if (bDoVdW && bDoCoul)
870 new_i_nblist(vdwc, i_atom, shift, gid);
874 new_i_nblist(vdw, i_atom, shift, gid);
878 new_i_nblist(coul, i_atom, shift, gid);
880 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
881 bDoCoul_i = (bDoCoul && qi != 0);
883 if (bDoVdW_i || bDoCoul_i)
885 /* Loop over the j charge groups */
886 for (j = 0; (j < nj); j++)
890 /* Check for large charge groups */
901 /* Finally loop over the atoms in the j-charge group */
902 for (jj = jj0; jj < jj1; jj++)
904 bNotEx = NOTEXCL(bExcl, i, jj);
912 add_j_to_nblist(coul, jj, bLR);
917 if (bHaveVdW[type[jj]])
919 add_j_to_nblist(vdw, jj, bLR);
924 if (bHaveVdW[type[jj]])
928 add_j_to_nblist(vdwc, jj, bLR);
932 add_j_to_nblist(vdw, jj, bLR);
935 else if (charge[jj] != 0)
937 add_j_to_nblist(coul, jj, bLR);
945 close_i_nblist(coul);
946 close_i_nblist(vdwc);
952 /* we are doing free energy */
953 vdwc_free = &nlist[eNL_VDWQQ_FREE];
954 vdw_free = &nlist[eNL_VDW_FREE];
955 coul_free = &nlist[eNL_QQ_FREE];
956 /* Loop over the atoms in the i charge group */
957 for (i = 0; i < nicg; i++)
960 gid = GID(igid, jgid, ngid);
962 qiB = chargeB[i_atom];
964 /* Create new i_atom for each energy group */
965 if (bDoVdW && bDoCoul)
967 new_i_nblist(vdwc, i_atom, shift, gid);
971 new_i_nblist(vdw, i_atom, shift, gid);
975 new_i_nblist(coul, i_atom, shift, gid);
978 new_i_nblist(vdw_free, i_atom, shift, gid);
979 new_i_nblist(coul_free, i_atom, shift, gid);
980 new_i_nblist(vdwc_free, i_atom, shift, gid);
982 bDoVdW_i = (bDoVdW &&
983 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
984 bDoCoul_i = (bDoCoul && (qi != 0 || qiB != 0));
985 /* For TIP4P the first atom does not have a charge,
986 * but the last three do. So we should still put an atom
987 * without LJ but with charge in the water-atom neighborlist
988 * for a TIP4p i charge group.
989 * For SPC type water the first atom has LJ and charge,
990 * so there is no such problem.
992 if (iwater == esolNO)
994 bDoCoul_i_sol = bDoCoul_i;
998 bDoCoul_i_sol = bDoCoul;
1001 if (bDoVdW_i || bDoCoul_i_sol)
1003 /* Loop over the j charge groups */
1004 for (j = 0; (j < nj); j++)
1008 /* Check for large charge groups */
1019 /* Finally loop over the atoms in the j-charge group */
1020 bFree = bPert[i_atom];
1021 for (jj = jj0; (jj < jj1); jj++)
1023 bFreeJ = bFree || bPert[jj];
1024 /* Complicated if, because the water H's should also
1025 * see perturbed j-particles
1027 if (iwater == esolNO || i == 0 || bFreeJ)
1029 bNotEx = NOTEXCL(bExcl, i, jj);
1037 if (charge[jj] != 0 || chargeB[jj] != 0)
1039 add_j_to_nblist(coul_free, jj, bLR);
1042 else if (!bDoCoul_i)
1044 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1046 add_j_to_nblist(vdw_free, jj, bLR);
1051 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1053 if (charge[jj] != 0 || chargeB[jj] != 0)
1055 add_j_to_nblist(vdwc_free, jj, bLR);
1059 add_j_to_nblist(vdw_free, jj, bLR);
1062 else if (charge[jj] != 0 || chargeB[jj] != 0)
1064 add_j_to_nblist(coul_free, jj, bLR);
1070 /* This is done whether or not bWater is set */
1071 if (charge[jj] != 0)
1073 add_j_to_nblist(coul, jj, bLR);
1076 else if (!bDoCoul_i_sol)
1078 if (bHaveVdW[type[jj]])
1080 add_j_to_nblist(vdw, jj, bLR);
1085 if (bHaveVdW[type[jj]])
1087 if (charge[jj] != 0)
1089 add_j_to_nblist(vdwc, jj, bLR);
1093 add_j_to_nblist(vdw, jj, bLR);
1096 else if (charge[jj] != 0)
1098 add_j_to_nblist(coul, jj, bLR);
1106 close_i_nblist(vdw);
1107 close_i_nblist(coul);
1108 close_i_nblist(vdwc);
1109 close_i_nblist(vdw_free);
1110 close_i_nblist(coul_free);
1111 close_i_nblist(vdwc_free);
1117 put_in_list_adress(gmx_bool bHaveVdW[],
1133 /* The a[] index has been removed,
1134 * to put it back in i_atom should be a[i0] and jj should be a[jj].
1139 t_nblist * vdwc_adress = NULL;
1140 t_nblist * vdw_adress = NULL;
1141 t_nblist * coul_adress = NULL;
1142 t_nblist * vdwc_ww = NULL;
1143 t_nblist * coul_ww = NULL;
1145 int i, j, jcg, igid, gid, nbl_ind, nbl_ind_adress;
1146 atom_id jj, jj0, jj1, i_atom;
1151 real *charge, *chargeB;
1153 real qi, qiB, qq, rlj;
1154 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
1155 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
1157 gmx_bool j_all_atom;
1159 t_nblist *nlist, *nlist_adress;
1160 gmx_bool bEnergyGroupCG;
1162 /* Copy some pointers */
1163 cginfo = fr->cginfo;
1164 charge = md->chargeA;
1165 chargeB = md->chargeB;
1168 bPert = md->bPerturbed;
1171 /* Get atom range */
1173 nicg = index[icg+1]-i0;
1175 /* Get the i charge group info */
1176 igid = GET_CGINFO_GID(cginfo[icg]);
1178 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
1182 gmx_fatal(FARGS, "AdResS does not support free energy pertubation\n");
1185 /* Unpack pointers to neighbourlist structs */
1186 if (fr->nnblists == 2)
1193 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
1194 nbl_ind_adress = nbl_ind+fr->nnblists/2;
1198 nlist = fr->nblists[nbl_ind].nlist_lr;
1199 nlist_adress = fr->nblists[nbl_ind_adress].nlist_lr;
1203 nlist = fr->nblists[nbl_ind].nlist_sr;
1204 nlist_adress = fr->nblists[nbl_ind_adress].nlist_sr;
1208 vdwc = &nlist[eNL_VDWQQ];
1209 vdw = &nlist[eNL_VDW];
1210 coul = &nlist[eNL_QQ];
1212 vdwc_adress = &nlist_adress[eNL_VDWQQ];
1213 vdw_adress = &nlist_adress[eNL_VDW];
1214 coul_adress = &nlist_adress[eNL_QQ];
1216 /* We do not support solvent optimization with AdResS for now.
1217 For this we would need hybrid solvent-other kernels */
1219 /* no solvent as i charge group */
1220 /* Loop over the atoms in the i charge group */
1221 for (i = 0; i < nicg; i++)
1224 gid = GID(igid, jgid, ngid);
1225 qi = charge[i_atom];
1227 /* Create new i_atom for each energy group */
1228 if (bDoVdW && bDoCoul)
1230 new_i_nblist(vdwc, i_atom, shift, gid);
1231 new_i_nblist(vdwc_adress, i_atom, shift, gid);
1236 new_i_nblist(vdw, i_atom, shift, gid);
1237 new_i_nblist(vdw_adress, i_atom, shift, gid);
1242 new_i_nblist(coul, i_atom, shift, gid);
1243 new_i_nblist(coul_adress, i_atom, shift, gid);
1245 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
1246 bDoCoul_i = (bDoCoul && qi != 0);
1248 /* Here we find out whether the energy groups interaction belong to a
1249 * coarse-grained (vsite) or atomistic interaction. Note that, beacuse
1250 * interactions between coarse-grained and other (atomistic) energygroups
1251 * are excluded automatically by grompp, it is sufficient to check for
1252 * the group id of atom i (igid) */
1253 bEnergyGroupCG = !egp_explicit(fr, igid);
1255 if (bDoVdW_i || bDoCoul_i)
1257 /* Loop over the j charge groups */
1258 for (j = 0; (j < nj); j++)
1262 /* Check for large charge groups */
1273 /* Finally loop over the atoms in the j-charge group */
1274 for (jj = jj0; jj < jj1; jj++)
1276 bNotEx = NOTEXCL(bExcl, i, jj);
1278 /* Now we have to exclude interactions which will be zero
1279 * anyway due to the AdResS weights (in previous implementations
1280 * this was done in the force kernel). This is necessary as
1281 * pure interactions (those with b_hybrid=false, i.e. w_i*w_j==1 or 0)
1282 * are put into neighbour lists which will be passed to the
1283 * standard (optimized) kernels for speed. The interactions with
1284 * b_hybrid=true are placed into the _adress neighbour lists and
1285 * processed by the generic AdResS kernel.
1287 if ( (bEnergyGroupCG &&
1288 wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS ) ||
1289 ( !bEnergyGroupCG && wf[jj] <= GMX_REAL_EPS ) )
1294 b_hybrid = !((wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS) ||
1295 (wf[i_atom] <= GMX_REAL_EPS && wf[jj] <= GMX_REAL_EPS));
1301 if (charge[jj] != 0)
1305 add_j_to_nblist(coul, jj, bLR);
1309 add_j_to_nblist(coul_adress, jj, bLR);
1313 else if (!bDoCoul_i)
1315 if (bHaveVdW[type[jj]])
1319 add_j_to_nblist(vdw, jj, bLR);
1323 add_j_to_nblist(vdw_adress, jj, bLR);
1329 if (bHaveVdW[type[jj]])
1331 if (charge[jj] != 0)
1335 add_j_to_nblist(vdwc, jj, bLR);
1339 add_j_to_nblist(vdwc_adress, jj, bLR);
1346 add_j_to_nblist(vdw, jj, bLR);
1350 add_j_to_nblist(vdw_adress, jj, bLR);
1355 else if (charge[jj] != 0)
1359 add_j_to_nblist(coul, jj, bLR);
1363 add_j_to_nblist(coul_adress, jj, bLR);
1372 close_i_nblist(vdw);
1373 close_i_nblist(coul);
1374 close_i_nblist(vdwc);
1375 close_i_nblist(vdw_adress);
1376 close_i_nblist(coul_adress);
1377 close_i_nblist(vdwc_adress);
1383 put_in_list_qmmm(gmx_bool gmx_unused bHaveVdW[],
1385 t_mdatoms gmx_unused * md,
1395 gmx_bool gmx_unused bDoVdW,
1396 gmx_bool gmx_unused bDoCoul,
1397 int gmx_unused solvent_opt)
1400 int i, j, jcg, igid, gid;
1401 atom_id jj, jj0, jj1, i_atom;
1405 /* Get atom range */
1407 nicg = index[icg+1]-i0;
1409 /* Get the i charge group info */
1410 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1412 coul = &fr->QMMMlist;
1414 /* Loop over atoms in the ith charge group */
1415 for (i = 0; i < nicg; i++)
1418 gid = GID(igid, jgid, ngid);
1419 /* Create new i_atom for each energy group */
1420 new_i_nblist(coul, i_atom, shift, gid);
1422 /* Loop over the j charge groups */
1423 for (j = 0; j < nj; j++)
1427 /* Charge groups cannot have QM and MM atoms simultaneously */
1432 /* Finally loop over the atoms in the j-charge group */
1433 for (jj = jj0; jj < jj1; jj++)
1435 bNotEx = NOTEXCL(bExcl, i, jj);
1438 add_j_to_nblist(coul, jj, bLR);
1443 close_i_nblist(coul);
1448 put_in_list_cg(gmx_bool gmx_unused bHaveVdW[],
1450 t_mdatoms gmx_unused * md,
1460 gmx_bool gmx_unused bDoVdW,
1461 gmx_bool gmx_unused bDoCoul,
1462 int gmx_unused solvent_opt)
1465 int igid, gid, nbl_ind;
1469 cginfo = fr->cginfo[icg];
1471 igid = GET_CGINFO_GID(cginfo);
1472 gid = GID(igid, jgid, ngid);
1474 /* Unpack pointers to neighbourlist structs */
1475 if (fr->nnblists == 1)
1481 nbl_ind = fr->gid2nblists[gid];
1485 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1489 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1492 /* Make a new neighbor list for charge group icg.
1493 * Currently simply one neighbor list is made with LJ and Coulomb.
1494 * If required, zero interactions could be removed here
1495 * or in the force loop.
1497 new_i_nblist(vdwc, index[icg], shift, gid);
1498 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1500 for (j = 0; (j < nj); j++)
1503 /* Skip the icg-icg pairs if all self interactions are excluded */
1504 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1506 /* Here we add the j charge group jcg to the list,
1507 * exclusions are also added to the list.
1509 add_j_to_nblist_cg(vdwc, index[jcg], index[jcg+1], bExcl, icg == jcg, bLR);
1513 close_i_nblist(vdwc);
1516 static void setexcl(atom_id start, atom_id end, t_blocka *excl, gmx_bool b,
1523 for (i = start; i < end; i++)
1525 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1527 SETEXCL(bexcl, i-start, excl->a[k]);
1533 for (i = start; i < end; i++)
1535 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1537 RMEXCL(bexcl, i-start, excl->a[k]);
1543 int calc_naaj(int icg, int cgtot)
1547 if ((cgtot % 2) == 1)
1549 /* Odd number of charge groups, easy */
1550 naaj = 1 + (cgtot/2);
1552 else if ((cgtot % 4) == 0)
1554 /* Multiple of four is hard */
1591 fprintf(log, "naaj=%d\n", naaj);
1597 /************************************************
1599 * S I M P L E C O R E S T U F F
1601 ************************************************/
1603 static real calc_image_tric(rvec xi, rvec xj, matrix box,
1604 rvec b_inv, int *shift)
1606 /* This code assumes that the cut-off is smaller than
1607 * a half times the smallest diagonal element of the box.
1609 const real h25 = 2.5;
1614 /* Compute diff vector */
1615 dz = xj[ZZ] - xi[ZZ];
1616 dy = xj[YY] - xi[YY];
1617 dx = xj[XX] - xi[XX];
1619 /* Perform NINT operation, using trunc operation, therefore
1620 * we first add 2.5 then subtract 2 again
1622 tz = dz*b_inv[ZZ] + h25;
1624 dz -= tz*box[ZZ][ZZ];
1625 dy -= tz*box[ZZ][YY];
1626 dx -= tz*box[ZZ][XX];
1628 ty = dy*b_inv[YY] + h25;
1630 dy -= ty*box[YY][YY];
1631 dx -= ty*box[YY][XX];
1633 tx = dx*b_inv[XX]+h25;
1635 dx -= tx*box[XX][XX];
1637 /* Distance squared */
1638 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1640 *shift = XYZ2IS(tx, ty, tz);
1645 static real calc_image_rect(rvec xi, rvec xj, rvec box_size,
1646 rvec b_inv, int *shift)
1648 const real h15 = 1.5;
1654 /* Compute diff vector */
1655 dx = xj[XX] - xi[XX];
1656 dy = xj[YY] - xi[YY];
1657 dz = xj[ZZ] - xi[ZZ];
1659 /* Perform NINT operation, using trunc operation, therefore
1660 * we first add 1.5 then subtract 1 again
1662 tx = dx*b_inv[XX] + h15;
1663 ty = dy*b_inv[YY] + h15;
1664 tz = dz*b_inv[ZZ] + h15;
1669 /* Correct diff vector for translation */
1670 ddx = tx*box_size[XX] - dx;
1671 ddy = ty*box_size[YY] - dy;
1672 ddz = tz*box_size[ZZ] - dz;
1674 /* Distance squared */
1675 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1677 *shift = XYZ2IS(tx, ty, tz);
1682 static void add_simple(t_ns_buf *nsbuf, int nrj, atom_id cg_j,
1683 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1684 int icg, int jgid, t_block *cgs, t_excl bexcl[],
1685 int shift, t_forcerec *fr, put_in_list_t *put_in_list)
1687 if (nsbuf->nj + nrj > MAX_CG)
1689 put_in_list(bHaveVdW, ngid, md, icg, jgid, nsbuf->ncg, nsbuf->jcg,
1690 cgs->index, bexcl, shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1691 /* Reset buffer contents */
1692 nsbuf->ncg = nsbuf->nj = 0;
1694 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1698 static void ns_inner_tric(rvec x[], int icg, int *i_egp_flags,
1699 int njcg, atom_id jcg[],
1700 matrix box, rvec b_inv, real rcut2,
1701 t_block *cgs, t_ns_buf **ns_buf,
1702 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1703 t_excl bexcl[], t_forcerec *fr,
1704 put_in_list_t *put_in_list)
1708 int *cginfo = fr->cginfo;
1709 atom_id cg_j, *cgindex;
1712 cgindex = cgs->index;
1714 for (j = 0; (j < njcg); j++)
1717 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1718 if (calc_image_tric(x[icg], x[cg_j], box, b_inv, &shift) < rcut2)
1720 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1721 if (!(i_egp_flags[jgid] & EGP_EXCL))
1723 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1724 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1731 static void ns_inner_rect(rvec x[], int icg, int *i_egp_flags,
1732 int njcg, atom_id jcg[],
1733 gmx_bool bBox, rvec box_size, rvec b_inv, real rcut2,
1734 t_block *cgs, t_ns_buf **ns_buf,
1735 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1736 t_excl bexcl[], t_forcerec *fr,
1737 put_in_list_t *put_in_list)
1741 int *cginfo = fr->cginfo;
1742 atom_id cg_j, *cgindex;
1745 cgindex = cgs->index;
1749 for (j = 0; (j < njcg); j++)
1752 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1753 if (calc_image_rect(x[icg], x[cg_j], box_size, b_inv, &shift) < rcut2)
1755 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1756 if (!(i_egp_flags[jgid] & EGP_EXCL))
1758 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1759 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1767 for (j = 0; (j < njcg); j++)
1770 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1771 if ((rcut2 == 0) || (distance2(x[icg], x[cg_j]) < rcut2))
1773 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1774 if (!(i_egp_flags[jgid] & EGP_EXCL))
1776 add_simple(&ns_buf[jgid][CENTRAL], nrj, cg_j,
1777 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, CENTRAL, fr,
1785 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1787 static int ns_simple_core(t_forcerec *fr,
1788 gmx_localtop_t *top,
1790 matrix box, rvec box_size,
1791 t_excl bexcl[], atom_id *aaj,
1792 int ngid, t_ns_buf **ns_buf,
1793 put_in_list_t *put_in_list, gmx_bool bHaveVdW[])
1797 int nsearch, icg, jcg, igid, i0, nri, nn;
1800 /* atom_id *i_atoms; */
1801 t_block *cgs = &(top->cgs);
1802 t_blocka *excl = &(top->excls);
1805 gmx_bool bBox, bTriclinic;
1808 rlist2 = sqr(fr->rlist);
1810 bBox = (fr->ePBC != epbcNONE);
1813 for (m = 0; (m < DIM); m++)
1815 b_inv[m] = divide_err(1.0, box_size[m]);
1817 bTriclinic = TRICLINIC(box);
1824 cginfo = fr->cginfo;
1827 for (icg = fr->cg0; (icg < fr->hcg); icg++)
1830 i0 = cgs->index[icg];
1831 nri = cgs->index[icg+1]-i0;
1832 i_atoms = &(cgs->a[i0]);
1833 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1834 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1836 igid = GET_CGINFO_GID(cginfo[icg]);
1837 i_egp_flags = fr->egp_flags + ngid*igid;
1838 setexcl(cgs->index[icg], cgs->index[icg+1], excl, TRUE, bexcl);
1840 naaj = calc_naaj(icg, cgs->nr);
1843 ns_inner_tric(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1844 box, b_inv, rlist2, cgs, ns_buf,
1845 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1849 ns_inner_rect(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1850 bBox, box_size, b_inv, rlist2, cgs, ns_buf,
1851 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1855 for (nn = 0; (nn < ngid); nn++)
1857 for (k = 0; (k < SHIFTS); k++)
1859 nsbuf = &(ns_buf[nn][k]);
1862 put_in_list(bHaveVdW, ngid, md, icg, nn, nsbuf->ncg, nsbuf->jcg,
1863 cgs->index, bexcl, k, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1864 nsbuf->ncg = nsbuf->nj = 0;
1868 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1869 setexcl(cgs->index[icg], cgs->index[icg+1], excl, FALSE, bexcl);
1871 close_neighbor_lists(fr, FALSE);
1876 /************************************************
1878 * N S 5 G R I D S T U F F
1880 ************************************************/
1882 static gmx_inline void get_dx(int Nx, real gridx, real rc2, int xgi, real x,
1883 int *dx0, int *dx1, real *dcx2)
1911 for (i = xgi0; i >= 0; i--)
1913 dcx = (i+1)*gridx-x;
1922 for (i = xgi1; i < Nx; i++)
1935 static gmx_inline void get_dx_dd(int Nx, real gridx, real rc2, int xgi, real x,
1936 int ncpddc, int shift_min, int shift_max,
1937 int *g0, int *g1, real *dcx2)
1940 int g_min, g_max, shift_home;
1973 g_min = (shift_min == shift_home ? 0 : ncpddc);
1974 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1981 else if (shift_max < 0)
1996 /* Check one grid cell down */
1997 dcx = ((*g0 - 1) + 1)*gridx - x;
2009 /* Check one grid cell up */
2010 dcx = (*g1 + 1)*gridx - x;
2022 #define sqr(x) ((x)*(x))
2023 #define calc_dx2(XI, YI, ZI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
2024 #define calc_cyl_dx2(XI, YI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
2025 /****************************************************
2027 * F A S T N E I G H B O R S E A R C H I N G
2029 * Optimized neighboursearching routine using grid
2030 * at least 1x1x1, see GROMACS manual
2032 ****************************************************/
2035 static void get_cutoff2(t_forcerec *fr, gmx_bool bDoLongRange,
2036 real *rvdw2, real *rcoul2,
2037 real *rs2, real *rm2, real *rl2)
2039 *rs2 = sqr(fr->rlist);
2041 if (bDoLongRange && fr->bTwinRange)
2043 /* With plain cut-off or RF we need to make the list exactly
2044 * up to the cut-off and the cut-off's can be different,
2045 * so we can not simply set them to rlistlong.
2046 * To keep this code compatible with (exotic) old cases,
2047 * we also create lists up to rvdw/rcoulomb for PME and Ewald.
2048 * The interaction check should correspond to:
2049 * !ir_vdw/coulomb_might_be_zero_at_cutoff from inputrec.c.
2051 if (((fr->vdwtype == evdwCUT || fr->vdwtype == evdwPME) &&
2052 fr->vdw_modifier == eintmodNONE) ||
2053 fr->rvdw <= fr->rlist)
2055 *rvdw2 = sqr(fr->rvdw);
2059 *rvdw2 = sqr(fr->rlistlong);
2061 if (((fr->eeltype == eelCUT ||
2062 (EEL_RF(fr->eeltype) && fr->eeltype != eelRF_ZERO) ||
2063 fr->eeltype == eelPME ||
2064 fr->eeltype == eelEWALD) &&
2065 fr->coulomb_modifier == eintmodNONE) ||
2066 fr->rcoulomb <= fr->rlist)
2068 *rcoul2 = sqr(fr->rcoulomb);
2072 *rcoul2 = sqr(fr->rlistlong);
2077 /* Workaround for a gcc -O3 or -ffast-math problem */
2081 *rm2 = min(*rvdw2, *rcoul2);
2082 *rl2 = max(*rvdw2, *rcoul2);
2085 static void init_nsgrid_lists(t_forcerec *fr, int ngid, gmx_ns_t *ns)
2087 real rvdw2, rcoul2, rs2, rm2, rl2;
2090 get_cutoff2(fr, TRUE, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2092 /* Short range buffers */
2093 snew(ns->nl_sr, ngid);
2095 snew(ns->nsr, ngid);
2096 snew(ns->nlr_ljc, ngid);
2097 snew(ns->nlr_one, ngid);
2099 /* Always allocate both list types, since rcoulomb might now change with PME load balancing */
2100 /* Long range VdW and Coul buffers */
2101 snew(ns->nl_lr_ljc, ngid);
2102 /* Long range VdW or Coul only buffers */
2103 snew(ns->nl_lr_one, ngid);
2105 for (j = 0; (j < ngid); j++)
2107 snew(ns->nl_sr[j], MAX_CG);
2108 snew(ns->nl_lr_ljc[j], MAX_CG);
2109 snew(ns->nl_lr_one[j], MAX_CG);
2114 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
2119 static int nsgrid_core(t_commrec *cr, t_forcerec *fr,
2120 matrix box, int ngid,
2121 gmx_localtop_t *top,
2123 t_excl bexcl[], gmx_bool *bExcludeAlleg,
2125 put_in_list_t *put_in_list,
2126 gmx_bool bHaveVdW[],
2127 gmx_bool bDoLongRange, gmx_bool bMakeQMMMnblist)
2130 atom_id **nl_lr_ljc, **nl_lr_one, **nl_sr;
2131 int *nlr_ljc, *nlr_one, *nsr;
2132 gmx_domdec_t *dd = NULL;
2133 t_block *cgs = &(top->cgs);
2134 int *cginfo = fr->cginfo;
2135 /* atom_id *i_atoms,*cgsindex=cgs->index; */
2137 int cell_x, cell_y, cell_z;
2138 int d, tx, ty, tz, dx, dy, dz, cj;
2139 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2140 int zsh_ty, zsh_tx, ysh_tx;
2142 int dx0, dx1, dy0, dy1, dz0, dz1;
2143 int Nx, Ny, Nz, shift = -1, j, nrj, nns, nn = -1;
2144 real gridx, gridy, gridz, grid_x, grid_y, grid_z;
2145 real *dcx2, *dcy2, *dcz2;
2147 int cg0, cg1, icg = -1, cgsnr, i0, igid, nri, naaj, max_jcg;
2148 int jcg0, jcg1, jjcg, cgj0, jgid;
2149 int *grida, *gridnra, *gridind;
2150 gmx_bool rvdw_lt_rcoul, rcoul_lt_rvdw;
2151 rvec xi, *cgcm, grid_offset;
2152 real r2, rs2, rvdw2, rcoul2, rm2, rl2, XI, YI, ZI, dcx, dcy, dcz, tmp1, tmp2;
2154 gmx_bool bDomDec, bTriclinicX, bTriclinicY;
2159 bDomDec = DOMAINDECOMP(cr);
2165 bTriclinicX = ((YY < grid->npbcdim &&
2166 (!bDomDec || dd->nc[YY] == 1) && box[YY][XX] != 0) ||
2167 (ZZ < grid->npbcdim &&
2168 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][XX] != 0));
2169 bTriclinicY = (ZZ < grid->npbcdim &&
2170 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][YY] != 0);
2174 get_cutoff2(fr, bDoLongRange, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2176 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
2177 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
2179 if (bMakeQMMMnblist)
2187 nl_lr_ljc = ns->nl_lr_ljc;
2188 nl_lr_one = ns->nl_lr_one;
2189 nlr_ljc = ns->nlr_ljc;
2190 nlr_one = ns->nlr_one;
2198 gridind = grid->index;
2199 gridnra = grid->nra;
2202 gridx = grid->cell_size[XX];
2203 gridy = grid->cell_size[YY];
2204 gridz = grid->cell_size[ZZ];
2208 copy_rvec(grid->cell_offset, grid_offset);
2209 copy_ivec(grid->ncpddc, ncpddc);
2214 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2215 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
2216 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
2217 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
2218 if (zsh_tx != 0 && ysh_tx != 0)
2220 /* This could happen due to rounding, when both ratios are 0.5 */
2229 /* We only want a list for the test particle */
2238 /* Set the shift range */
2239 for (d = 0; d < DIM; d++)
2243 /* Check if we need periodicity shifts.
2244 * Without PBC or with domain decomposition we don't need them.
2246 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2253 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
2264 /* Loop over charge groups */
2265 for (icg = cg0; (icg < cg1); icg++)
2267 igid = GET_CGINFO_GID(cginfo[icg]);
2268 /* Skip this charge group if all energy groups are excluded! */
2269 if (bExcludeAlleg[igid])
2274 i0 = cgs->index[icg];
2276 if (bMakeQMMMnblist)
2278 /* Skip this charge group if it is not a QM atom while making a
2279 * QM/MM neighbourlist
2281 if (md->bQM[i0] == FALSE)
2283 continue; /* MM particle, go to next particle */
2286 /* Compute the number of charge groups that fall within the control
2289 naaj = calc_naaj(icg, cgsnr);
2296 /* make a normal neighbourlist */
2300 /* Get the j charge-group and dd cell shift ranges */
2301 dd_get_ns_ranges(cr->dd, icg, &jcg0, &jcg1, sh0, sh1);
2306 /* Compute the number of charge groups that fall within the control
2309 naaj = calc_naaj(icg, cgsnr);
2315 /* The i-particle is awlways the test particle,
2316 * so we want all j-particles
2318 max_jcg = cgsnr - 1;
2322 max_jcg = jcg1 - cgsnr;
2327 i_egp_flags = fr->egp_flags + igid*ngid;
2329 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2330 setexcl(i0, cgs->index[icg+1], &top->excls, TRUE, bexcl);
2332 ci2xyz(grid, icg, &cell_x, &cell_y, &cell_z);
2334 /* Changed iicg to icg, DvdS 990115
2335 * (but see consistency check above, DvdS 990330)
2338 fprintf(log, "icg=%5d, naaj=%5d, cell %d %d %d\n",
2339 icg, naaj, cell_x, cell_y, cell_z);
2341 /* Loop over shift vectors in three dimensions */
2342 for (tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
2344 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2345 /* Calculate range of cells in Z direction that have the shift tz */
2346 zgi = cell_z + tz*Nz;
2349 get_dx(Nz, gridz, rl2, zgi, ZI, &dz0, &dz1, dcz2);
2351 get_dx_dd(Nz, gridz, rl2, zgi, ZI-grid_offset[ZZ],
2352 ncpddc[ZZ], sh0[ZZ], sh1[ZZ], &dz0, &dz1, dcz2);
2358 for (ty = -shp[YY]; ty <= shp[YY]; ty++)
2360 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2361 /* Calculate range of cells in Y direction that have the shift ty */
2364 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2368 ygi = cell_y + ty*Ny;
2371 get_dx(Ny, gridy, rl2, ygi, YI, &dy0, &dy1, dcy2);
2373 get_dx_dd(Ny, gridy, rl2, ygi, YI-grid_offset[YY],
2374 ncpddc[YY], sh0[YY], sh1[YY], &dy0, &dy1, dcy2);
2380 for (tx = -shp[XX]; tx <= shp[XX]; tx++)
2382 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2383 /* Calculate range of cells in X direction that have the shift tx */
2386 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2390 xgi = cell_x + tx*Nx;
2393 get_dx(Nx, gridx, rl2, xgi*Nx, XI, &dx0, &dx1, dcx2);
2395 get_dx_dd(Nx, gridx, rl2, xgi, XI-grid_offset[XX],
2396 ncpddc[XX], sh0[XX], sh1[XX], &dx0, &dx1, dcx2);
2402 /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2403 * from the neigbour list as it will not interact */
2404 if (fr->adress_type != eAdressOff)
2406 if (md->wf[cgs->index[icg]] <= GMX_REAL_EPS && egp_explicit(fr, igid))
2411 /* Get shift vector */
2412 shift = XYZ2IS(tx, ty, tz);
2414 range_check(shift, 0, SHIFTS);
2416 for (nn = 0; (nn < ngid); nn++)
2423 fprintf(log, "shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2424 shift, dx0, dx1, dy0, dy1, dz0, dz1);
2425 fprintf(log, "cgcm: %8.3f %8.3f %8.3f\n", cgcm[icg][XX],
2426 cgcm[icg][YY], cgcm[icg][ZZ]);
2427 fprintf(log, "xi: %8.3f %8.3f %8.3f\n", XI, YI, ZI);
2429 for (dx = dx0; (dx <= dx1); dx++)
2431 tmp1 = rl2 - dcx2[dx];
2432 for (dy = dy0; (dy <= dy1); dy++)
2434 tmp2 = tmp1 - dcy2[dy];
2437 for (dz = dz0; (dz <= dz1); dz++)
2439 if (tmp2 > dcz2[dz])
2441 /* Find grid-cell cj in which possible neighbours are */
2442 cj = xyz2ci(Ny, Nz, dx, dy, dz);
2444 /* Check out how many cgs (nrj) there in this cell */
2447 /* Find the offset in the cg list */
2450 /* Check if all j's are out of range so we
2451 * can skip the whole cell.
2452 * Should save some time, especially with DD.
2455 (grida[cgj0] >= max_jcg &&
2456 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2462 for (j = 0; (j < nrj); j++)
2464 jjcg = grida[cgj0+j];
2466 /* check whether this guy is in range! */
2467 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2470 r2 = calc_dx2(XI, YI, ZI, cgcm[jjcg]);
2473 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2474 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2475 /* check energy group exclusions */
2476 if (!(i_egp_flags[jgid] & EGP_EXCL))
2480 if (nsr[jgid] >= MAX_CG)
2482 /* Add to short-range list */
2483 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2484 nsr[jgid], nl_sr[jgid],
2485 cgs->index, /* cgsatoms, */ bexcl,
2486 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2489 nl_sr[jgid][nsr[jgid]++] = jjcg;
2493 if (nlr_ljc[jgid] >= MAX_CG)
2495 /* Add to LJ+coulomb long-range list */
2496 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2497 nlr_ljc[jgid], nl_lr_ljc[jgid], top->cgs.index,
2498 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2501 nl_lr_ljc[jgid][nlr_ljc[jgid]++] = jjcg;
2505 if (nlr_one[jgid] >= MAX_CG)
2507 /* Add to long-range list with only coul, or only LJ */
2508 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2509 nlr_one[jgid], nl_lr_one[jgid], top->cgs.index,
2510 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2513 nl_lr_one[jgid][nlr_one[jgid]++] = jjcg;
2525 /* CHECK whether there is anything left in the buffers */
2526 for (nn = 0; (nn < ngid); nn++)
2530 put_in_list(bHaveVdW, ngid, md, icg, nn, nsr[nn], nl_sr[nn],
2531 cgs->index, /* cgsatoms, */ bexcl,
2532 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2535 if (nlr_ljc[nn] > 0)
2537 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_ljc[nn],
2538 nl_lr_ljc[nn], top->cgs.index,
2539 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2542 if (nlr_one[nn] > 0)
2544 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_one[nn],
2545 nl_lr_one[nn], top->cgs.index,
2546 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2552 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2553 setexcl(cgs->index[icg], cgs->index[icg+1], &top->excls, FALSE, bexcl);
2555 /* No need to perform any left-over force calculations anymore (as we used to do here)
2556 * since we now save the proper long-range lists for later evaluation.
2561 /* Close neighbourlists */
2562 close_neighbor_lists(fr, bMakeQMMMnblist);
2567 void ns_realloc_natoms(gmx_ns_t *ns, int natoms)
2571 if (natoms > ns->nra_alloc)
2573 ns->nra_alloc = over_alloc_dd(natoms);
2574 srenew(ns->bexcl, ns->nra_alloc);
2575 for (i = 0; i < ns->nra_alloc; i++)
2582 void init_ns(FILE *fplog, const t_commrec *cr,
2583 gmx_ns_t *ns, t_forcerec *fr,
2584 const gmx_mtop_t *mtop)
2586 int mt, icg, nr_in_cg, maxcg, i, j, jcg, ngid, ncg;
2590 /* Compute largest charge groups size (# atoms) */
2592 for (mt = 0; mt < mtop->nmoltype; mt++)
2594 cgs = &mtop->moltype[mt].cgs;
2595 for (icg = 0; (icg < cgs->nr); icg++)
2597 nr_in_cg = max(nr_in_cg, (int)(cgs->index[icg+1]-cgs->index[icg]));
2601 /* Verify whether largest charge group is <= max cg.
2602 * This is determined by the type of the local exclusion type
2603 * Exclusions are stored in bits. (If the type is not large
2604 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2606 maxcg = sizeof(t_excl)*8;
2607 if (nr_in_cg > maxcg)
2609 gmx_fatal(FARGS, "Max #atoms in a charge group: %d > %d\n",
2613 ngid = mtop->groups.grps[egcENER].nr;
2614 snew(ns->bExcludeAlleg, ngid);
2615 for (i = 0; i < ngid; i++)
2617 ns->bExcludeAlleg[i] = TRUE;
2618 for (j = 0; j < ngid; j++)
2620 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2622 ns->bExcludeAlleg[i] = FALSE;
2630 ns->grid = init_grid(fplog, fr);
2631 init_nsgrid_lists(fr, ngid, ns);
2636 snew(ns->ns_buf, ngid);
2637 for (i = 0; (i < ngid); i++)
2639 snew(ns->ns_buf[i], SHIFTS);
2641 ncg = ncg_mtop(mtop);
2642 snew(ns->simple_aaj, 2*ncg);
2643 for (jcg = 0; (jcg < ncg); jcg++)
2645 ns->simple_aaj[jcg] = jcg;
2646 ns->simple_aaj[jcg+ncg] = jcg;
2650 /* Create array that determines whether or not atoms have VdW */
2651 snew(ns->bHaveVdW, fr->ntype);
2652 for (i = 0; (i < fr->ntype); i++)
2654 for (j = 0; (j < fr->ntype); j++)
2656 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2658 ((BHAMA(fr->nbfp, fr->ntype, i, j) != 0) ||
2659 (BHAMB(fr->nbfp, fr->ntype, i, j) != 0) ||
2660 (BHAMC(fr->nbfp, fr->ntype, i, j) != 0)) :
2661 ((C6(fr->nbfp, fr->ntype, i, j) != 0) ||
2662 (C12(fr->nbfp, fr->ntype, i, j) != 0))));
2667 pr_bvec(debug, 0, "bHaveVdW", ns->bHaveVdW, fr->ntype, TRUE);
2672 if (!DOMAINDECOMP(cr))
2674 ns_realloc_natoms(ns, mtop->natoms);
2677 ns->nblist_initialized = FALSE;
2679 /* nbr list debug dump */
2681 char *ptr = getenv("GMX_DUMP_NL");
2684 ns->dump_nl = strtol(ptr, NULL, 10);
2687 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2698 int search_neighbours(FILE *log, t_forcerec *fr,
2700 gmx_localtop_t *top,
2701 gmx_groups_t *groups,
2703 t_nrnb *nrnb, t_mdatoms *md,
2705 gmx_bool bDoLongRangeNS)
2707 t_block *cgs = &(top->cgs);
2708 rvec box_size, grid_x0, grid_x1;
2710 real min_size, grid_dens;
2714 gmx_bool *i_egp_flags;
2715 int cg_start, cg_end, start, end;
2718 gmx_domdec_zones_t *dd_zones;
2719 put_in_list_t *put_in_list;
2723 /* Set some local variables */
2725 ngid = groups->grps[egcENER].nr;
2727 for (m = 0; (m < DIM); m++)
2729 box_size[m] = box[m][m];
2732 if (fr->ePBC != epbcNONE)
2734 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC, box))
2736 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.");
2740 min_size = min(box_size[XX], min(box_size[YY], box_size[ZZ]));
2741 if (2*fr->rlistlong >= min_size)
2743 gmx_fatal(FARGS, "One of the box diagonal elements has become smaller than twice the cut-off length.");
2748 if (DOMAINDECOMP(cr))
2750 ns_realloc_natoms(ns, cgs->index[cgs->nr]);
2754 /* Reset the neighbourlists */
2755 reset_neighbor_lists(fr, TRUE, TRUE);
2757 if (bGrid && bFillGrid)
2761 if (DOMAINDECOMP(cr))
2763 dd_zones = domdec_zones(cr->dd);
2769 get_nsgrid_boundaries(grid->nboundeddim, box, NULL, NULL, NULL, NULL,
2770 cgs->nr, fr->cg_cm, grid_x0, grid_x1, &grid_dens);
2772 grid_first(log, grid, NULL, NULL, box, grid_x0, grid_x1,
2773 fr->rlistlong, grid_dens);
2780 if (DOMAINDECOMP(cr))
2783 fill_grid(dd_zones, grid, end, -1, end, fr->cg_cm);
2785 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2789 fill_grid(NULL, grid, cgs->nr, fr->cg0, fr->hcg, fr->cg_cm);
2790 grid->icg0 = fr->cg0;
2791 grid->icg1 = fr->hcg;
2795 calc_elemnr(grid, start, end, cgs->nr);
2797 grid_last(grid, start, end, cgs->nr);
2802 print_grid(debug, grid);
2807 /* Set the grid cell index for the test particle only.
2808 * The cell to cg index is not corrected, but that does not matter.
2810 fill_grid(NULL, ns->grid, fr->hcg, fr->hcg-1, fr->hcg, fr->cg_cm);
2814 if (fr->adress_type == eAdressOff)
2816 if (!fr->ns.bCGlist)
2818 put_in_list = put_in_list_at;
2822 put_in_list = put_in_list_cg;
2827 put_in_list = put_in_list_adress;
2834 nsearch = nsgrid_core(cr, fr, box, ngid, top,
2835 grid, ns->bexcl, ns->bExcludeAlleg,
2836 md, put_in_list, ns->bHaveVdW,
2837 bDoLongRangeNS, FALSE);
2839 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2840 * the classical calculation. The charge-charge interaction
2841 * between QM and MM atoms is handled in the QMMM core calculation
2842 * (see QMMM.c). The VDW however, we'd like to compute classically
2843 * and the QM MM atom pairs have just been put in the
2844 * corresponding neighbourlists. in case of QMMM we still need to
2845 * fill a special QMMM neighbourlist that contains all neighbours
2846 * of the QM atoms. If bQMMM is true, this list will now be made:
2848 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
2850 nsearch += nsgrid_core(cr, fr, box, ngid, top,
2851 grid, ns->bexcl, ns->bExcludeAlleg,
2852 md, put_in_list_qmmm, ns->bHaveVdW,
2853 bDoLongRangeNS, TRUE);
2858 nsearch = ns_simple_core(fr, top, md, box, box_size,
2859 ns->bexcl, ns->simple_aaj,
2860 ngid, ns->ns_buf, put_in_list, ns->bHaveVdW);
2868 inc_nrnb(nrnb, eNR_NS, nsearch);
2869 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2874 int natoms_beyond_ns_buffer(t_inputrec *ir, t_forcerec *fr, t_block *cgs,
2875 matrix scale_tot, rvec *x)
2877 int cg0, cg1, cg, a0, a1, a, i, j;
2878 real rint, hbuf2, scale;
2880 gmx_bool bIsotropic;
2885 rint = max(ir->rcoulomb, ir->rvdw);
2886 if (ir->rlist < rint)
2888 gmx_fatal(FARGS, "The neighbor search buffer has negative size: %f nm",
2896 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2898 hbuf2 = sqr(0.5*(ir->rlist - rint));
2899 for (cg = cg0; cg < cg1; cg++)
2901 a0 = cgs->index[cg];
2902 a1 = cgs->index[cg+1];
2903 for (a = a0; a < a1; a++)
2905 if (distance2(cg_cm[cg], x[a]) > hbuf2)
2915 scale = scale_tot[0][0];
2916 for (i = 1; i < DIM; i++)
2918 /* With anisotropic scaling, the original spherical ns volumes become
2919 * ellipsoids. To avoid costly transformations we use the minimum
2920 * eigenvalue of the scaling matrix for determining the buffer size.
2921 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2923 scale = min(scale, scale_tot[i][i]);
2924 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2928 for (j = 0; j < i; j++)
2930 if (scale_tot[i][j] != 0)
2936 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2939 for (cg = cg0; cg < cg1; cg++)
2941 svmul(scale, cg_cm[cg], cgsc);
2942 a0 = cgs->index[cg];
2943 a1 = cgs->index[cg+1];
2944 for (a = a0; a < a1; a++)
2946 if (distance2(cgsc, x[a]) > hbuf2)
2955 /* Anistropic scaling */
2956 for (cg = cg0; cg < cg1; cg++)
2958 /* Since scale_tot contains the transpose of the scaling matrix,
2959 * we need to multiply with the transpose.
2961 tmvmul_ur0(scale_tot, cg_cm[cg], cgsc);
2962 a0 = cgs->index[cg];
2963 a1 = cgs->index[cg+1];
2964 for (a = a0; a < a1; a++)
2966 if (distance2(cgsc, x[a]) > hbuf2)