1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
50 #include "nonbonded.h"
54 #include "gmx_fatal.h"
57 #include "mtop_util.h"
64 * E X C L U S I O N H A N D L I N G
68 static void SETEXCL_(t_excl e[], atom_id i, atom_id j)
72 static void RMEXCL_(t_excl e[], atom_id i, atom_id j)
74 e[j] = e[j] & ~(1<<i);
76 static gmx_bool ISEXCL_(t_excl e[], atom_id i, atom_id j)
78 return (gmx_bool)(e[j] & (1<<i));
80 static gmx_bool NOTEXCL_(t_excl e[], atom_id i, atom_id j)
82 return !(ISEXCL(e, i, j));
85 #define SETEXCL(e, i, j) (e)[((atom_id) (j))] |= (1<<((atom_id) (i)))
86 #define RMEXCL(e, i, j) (e)[((atom_id) (j))] &= (~(1<<((atom_id) (i))))
87 #define ISEXCL(e, i, j) (gmx_bool) ((e)[((atom_id) (j))] & (1<<((atom_id) (i))))
88 #define NOTEXCL(e, i, j) !(ISEXCL(e, i, j))
92 round_up_to_simd_width(int length, int simd_width)
94 int offset, newlength;
96 offset = (simd_width > 0) ? length % simd_width : 0;
98 return (offset == 0) ? length : length-offset+simd_width;
100 /************************************************
102 * U T I L I T I E S F O R N S
104 ************************************************/
106 static void reallocate_nblist(t_nblist *nl)
110 fprintf(debug, "reallocating neigborlist (ielec=%d, ivdw=%d, igeometry=%d, type=%d), maxnri=%d\n",
111 nl->ielec, nl->ivdw, nl->igeometry, nl->type, nl->maxnri);
113 srenew(nl->iinr, nl->maxnri);
114 if (nl->igeometry == GMX_NBLIST_GEOMETRY_CG_CG)
116 srenew(nl->iinr_end, nl->maxnri);
118 srenew(nl->gid, nl->maxnri);
119 srenew(nl->shift, nl->maxnri);
120 srenew(nl->jindex, nl->maxnri+1);
124 static void init_nblist(FILE *log, t_nblist *nl_sr, t_nblist *nl_lr,
125 int maxsr, int maxlr,
126 int ivdw, int ivdwmod,
127 int ielec, int ielecmod,
128 int igeometry, int type)
134 for (i = 0; (i < 2); i++)
136 nl = (i == 0) ? nl_sr : nl_lr;
137 homenr = (i == 0) ? maxsr : maxlr;
145 /* Set coul/vdw in neighborlist, and for the normal loops we determine
146 * an index of which one to call.
149 nl->ivdwmod = ivdwmod;
151 nl->ielecmod = ielecmod;
153 nl->igeometry = igeometry;
155 if (nl->type == GMX_NBLIST_INTERACTION_FREE_ENERGY)
157 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
160 /* This will also set the simd_padding_width field */
161 gmx_nonbonded_set_kernel_pointers( (i == 0) ? log : NULL, nl);
163 /* maxnri is influenced by the number of shifts (maximum is 8)
164 * and the number of energy groups.
165 * If it is not enough, nl memory will be reallocated during the run.
166 * 4 seems to be a reasonable factor, which only causes reallocation
167 * during runs with tiny and many energygroups.
169 nl->maxnri = homenr*4;
178 reallocate_nblist(nl);
183 fprintf(debug, "Initiating neighbourlist (ielec=%d, ivdw=%d, type=%d) for %s interactions,\nwith %d SR, %d LR atoms.\n",
184 nl->ielec, nl->ivdw, nl->type, gmx_nblist_geometry_names[nl->igeometry], maxsr, maxlr);
189 void init_neighbor_list(FILE *log, t_forcerec *fr, int homenr)
191 /* Make maxlr tunable! (does not seem to be a big difference though)
192 * This parameter determines the number of i particles in a long range
193 * neighbourlist. Too few means many function calls, too many means
196 int maxsr, maxsr_wat, maxlr, maxlr_wat;
197 int ielec, ielecf, ivdw, ielecmod, ielecmodf, ivdwmod, type;
199 int igeometry_def, igeometry_w, igeometry_ww;
203 /* maxsr = homenr-fr->nWatMol*3; */
208 gmx_fatal(FARGS, "%s, %d: Negative number of short range atoms.\n"
209 "Call your Gromacs dealer for assistance.", __FILE__, __LINE__);
211 /* This is just for initial allocation, so we do not reallocate
212 * all the nlist arrays many times in a row.
213 * The numbers seem very accurate, but they are uncritical.
215 maxsr_wat = min(fr->nWatMol, (homenr+2)/3);
219 maxlr_wat = min(maxsr_wat, maxlr);
223 maxlr = maxlr_wat = 0;
226 /* Determine the values for ielec/ivdw. */
227 ielec = fr->nbkernel_elec_interaction;
228 ivdw = fr->nbkernel_vdw_interaction;
229 ielecmod = fr->nbkernel_elec_modifier;
230 ivdwmod = fr->nbkernel_vdw_modifier;
231 type = GMX_NBLIST_INTERACTION_STANDARD;
233 fr->ns.bCGlist = (getenv("GMX_NBLISTCG") != 0);
236 igeometry_def = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
240 igeometry_def = GMX_NBLIST_GEOMETRY_CG_CG;
243 fprintf(log, "\nUsing charge-group - charge-group neighbor lists and kernels\n\n");
247 if (fr->solvent_opt == esolTIP4P)
249 igeometry_w = GMX_NBLIST_GEOMETRY_WATER4_PARTICLE;
250 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER4_WATER4;
254 igeometry_w = GMX_NBLIST_GEOMETRY_WATER3_PARTICLE;
255 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER3_WATER3;
258 for (i = 0; i < fr->nnblists; i++)
260 nbl = &(fr->nblists[i]);
262 if ((fr->adress_type != eAdressOff) && (i >= fr->nnblists/2))
264 type = GMX_NBLIST_INTERACTION_ADRESS;
266 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ], &nbl->nlist_lr[eNL_VDWQQ],
267 maxsr, maxlr, ivdw, ivdwmod, ielec, ielecmod, igeometry_def, type);
268 init_nblist(log, &nbl->nlist_sr[eNL_VDW], &nbl->nlist_lr[eNL_VDW],
269 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, igeometry_def, type);
270 init_nblist(log, &nbl->nlist_sr[eNL_QQ], &nbl->nlist_lr[eNL_QQ],
271 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_def, type);
272 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATER], &nbl->nlist_lr[eNL_VDWQQ_WATER],
273 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_w, type);
274 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATER], &nbl->nlist_lr[eNL_QQ_WATER],
275 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_w, type);
276 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATERWATER], &nbl->nlist_lr[eNL_VDWQQ_WATERWATER],
277 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_ww, type);
278 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATERWATER], &nbl->nlist_lr[eNL_QQ_WATERWATER],
279 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_ww, type);
281 /* Did we get the solvent loops so we can use optimized water kernels? */
282 if (nbl->nlist_sr[eNL_VDWQQ_WATER].kernelptr_vf == NULL
283 || nbl->nlist_sr[eNL_QQ_WATER].kernelptr_vf == NULL
284 #ifndef DISABLE_WATERWATER_NLIST
285 || nbl->nlist_sr[eNL_VDWQQ_WATERWATER].kernelptr_vf == NULL
286 || nbl->nlist_sr[eNL_QQ_WATERWATER].kernelptr_vf == NULL
290 fr->solvent_opt = esolNO;
291 fprintf(log, "Note: The available nonbonded kernels do not support water optimization - disabling.\n");
294 if (fr->efep != efepNO)
296 if ((fr->bEwald) && (fr->sc_alphacoul > 0)) /* need to handle long range differently if using softcore */
298 ielecf = GMX_NBKERNEL_ELEC_EWALD;
299 ielecmodf = eintmodNONE;
304 ielecmodf = ielecmod;
307 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_FREE], &nbl->nlist_lr[eNL_VDWQQ_FREE],
308 maxsr, maxlr, ivdw, ivdwmod, ielecf, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
309 init_nblist(log, &nbl->nlist_sr[eNL_VDW_FREE], &nbl->nlist_lr[eNL_VDW_FREE],
310 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
311 init_nblist(log, &nbl->nlist_sr[eNL_QQ_FREE], &nbl->nlist_lr[eNL_QQ_FREE],
312 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielecf, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
316 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
318 init_nblist(log, &fr->QMMMlist, NULL,
319 maxsr, maxlr, 0, 0, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_STANDARD);
327 fr->ns.nblist_initialized = TRUE;
330 static void reset_nblist(t_nblist *nl)
341 static void reset_neighbor_lists(t_forcerec *fr, gmx_bool bResetSR, gmx_bool bResetLR)
347 /* only reset the short-range nblist */
348 reset_nblist(&(fr->QMMMlist));
351 for (n = 0; n < fr->nnblists; n++)
353 for (i = 0; i < eNL_NR; i++)
357 reset_nblist( &(fr->nblists[n].nlist_sr[i]) );
361 reset_nblist( &(fr->nblists[n].nlist_lr[i]) );
370 static inline void new_i_nblist(t_nblist *nlist, atom_id i_atom, int shift, int gid)
372 int i, k, nri, nshift;
376 /* Check whether we have to increase the i counter */
378 (nlist->iinr[nri] != i_atom) ||
379 (nlist->shift[nri] != shift) ||
380 (nlist->gid[nri] != gid))
382 /* This is something else. Now see if any entries have
383 * been added in the list of the previous atom.
386 ((nlist->jindex[nri+1] > nlist->jindex[nri]) &&
387 (nlist->gid[nri] != -1)))
389 /* If so increase the counter */
392 if (nlist->nri >= nlist->maxnri)
394 nlist->maxnri += over_alloc_large(nlist->nri);
395 reallocate_nblist(nlist);
398 /* Set the number of neighbours and the atom number */
399 nlist->jindex[nri+1] = nlist->jindex[nri];
400 nlist->iinr[nri] = i_atom;
401 nlist->gid[nri] = gid;
402 nlist->shift[nri] = shift;
406 /* Adding to previous list. First remove possible previous padding */
407 if(nlist->simd_padding_width>1)
409 while(nlist->nrj>0 && nlist->jjnr[nlist->nrj-1]<0)
417 static inline void close_i_nblist(t_nblist *nlist)
419 int nri = nlist->nri;
424 /* Add elements up to padding. Since we allocate memory in units
425 * of the simd_padding width, we do not have to check for possible
426 * list reallocation here.
428 while ((nlist->nrj % nlist->simd_padding_width) != 0)
430 /* Use -4 here, so we can write forces for 4 atoms before real data */
431 nlist->jjnr[nlist->nrj++] = -4;
433 nlist->jindex[nri+1] = nlist->nrj;
435 len = nlist->nrj - nlist->jindex[nri];
437 /* nlist length for water i molecules is treated statically
440 if (len > nlist->maxlen)
447 static inline void close_nblist(t_nblist *nlist)
449 /* Only close this nblist when it has been initialized.
450 * Avoid the creation of i-lists with no j-particles.
454 /* Some assembly kernels do not support empty lists,
455 * make sure here that we don't generate any empty lists.
456 * With the current ns code this branch is taken in two cases:
457 * No i-particles at all: nri=-1 here
458 * There are i-particles, but no j-particles; nri=0 here
464 /* Close list number nri by incrementing the count */
469 static inline void close_neighbor_lists(t_forcerec *fr, gmx_bool bMakeQMMMnblist)
475 close_nblist(&(fr->QMMMlist));
478 for (n = 0; n < fr->nnblists; n++)
480 for (i = 0; (i < eNL_NR); i++)
482 close_nblist(&(fr->nblists[n].nlist_sr[i]));
483 close_nblist(&(fr->nblists[n].nlist_lr[i]));
489 static inline void add_j_to_nblist(t_nblist *nlist, atom_id j_atom, gmx_bool bLR)
491 int nrj = nlist->nrj;
493 if (nlist->nrj >= nlist->maxnrj)
495 nlist->maxnrj = round_up_to_simd_width(over_alloc_small(nlist->nrj + 1), nlist->simd_padding_width);
499 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
500 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
503 srenew(nlist->jjnr, nlist->maxnrj);
506 nlist->jjnr[nrj] = j_atom;
510 static inline void add_j_to_nblist_cg(t_nblist *nlist,
511 atom_id j_start, int j_end,
512 t_excl *bexcl, gmx_bool i_is_j,
515 int nrj = nlist->nrj;
518 if (nlist->nrj >= nlist->maxnrj)
520 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
523 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
524 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
527 srenew(nlist->jjnr, nlist->maxnrj);
528 srenew(nlist->jjnr_end, nlist->maxnrj);
529 srenew(nlist->excl, nlist->maxnrj*MAX_CGCGSIZE);
532 nlist->jjnr[nrj] = j_start;
533 nlist->jjnr_end[nrj] = j_end;
535 if (j_end - j_start > MAX_CGCGSIZE)
537 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);
540 /* Set the exclusions */
541 for (j = j_start; j < j_end; j++)
543 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
547 /* Avoid double counting of intra-cg interactions */
548 for (j = 1; j < j_end-j_start; j++)
550 nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
558 put_in_list_t (gmx_bool bHaveVdW[],
575 put_in_list_at(gmx_bool bHaveVdW[],
591 /* The a[] index has been removed,
592 * to put it back in i_atom should be a[i0] and jj should be a[jj].
597 t_nblist * vdwc_free = NULL;
598 t_nblist * vdw_free = NULL;
599 t_nblist * coul_free = NULL;
600 t_nblist * vdwc_ww = NULL;
601 t_nblist * coul_ww = NULL;
603 int i, j, jcg, igid, gid, nbl_ind, ind_ij;
604 atom_id jj, jj0, jj1, i_atom;
609 real *charge, *chargeB;
610 real qi, qiB, qq, rlj;
611 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
612 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
616 /* Copy some pointers */
618 charge = md->chargeA;
619 chargeB = md->chargeB;
622 bPert = md->bPerturbed;
626 nicg = index[icg+1]-i0;
628 /* Get the i charge group info */
629 igid = GET_CGINFO_GID(cginfo[icg]);
631 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
636 /* Check if any of the particles involved are perturbed.
637 * If not we can do the cheaper normal put_in_list
638 * and use more solvent optimization.
640 for (i = 0; i < nicg; i++)
642 bFreeEnergy |= bPert[i0+i];
644 /* Loop over the j charge groups */
645 for (j = 0; (j < nj && !bFreeEnergy); j++)
650 /* Finally loop over the atoms in the j-charge group */
651 for (jj = jj0; jj < jj1; jj++)
653 bFreeEnergy |= bPert[jj];
658 /* Unpack pointers to neighbourlist structs */
659 if (fr->nnblists == 1)
665 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
669 nlist = fr->nblists[nbl_ind].nlist_lr;
673 nlist = fr->nblists[nbl_ind].nlist_sr;
676 if (iwater != esolNO)
678 vdwc = &nlist[eNL_VDWQQ_WATER];
679 vdw = &nlist[eNL_VDW];
680 coul = &nlist[eNL_QQ_WATER];
681 #ifndef DISABLE_WATERWATER_NLIST
682 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
683 coul_ww = &nlist[eNL_QQ_WATERWATER];
688 vdwc = &nlist[eNL_VDWQQ];
689 vdw = &nlist[eNL_VDW];
690 coul = &nlist[eNL_QQ];
695 if (iwater != esolNO)
697 /* Loop over the atoms in the i charge group */
699 gid = GID(igid, jgid, ngid);
700 /* Create new i_atom for each energy group */
701 if (bDoCoul && bDoVdW)
703 new_i_nblist(vdwc, i_atom, shift, gid);
704 #ifndef DISABLE_WATERWATER_NLIST
705 new_i_nblist(vdwc_ww, i_atom, shift, gid);
710 new_i_nblist(vdw, i_atom, shift, gid);
714 new_i_nblist(coul, i_atom, shift, gid);
715 #ifndef DISABLE_WATERWATER_NLIST
716 new_i_nblist(coul_ww, i_atom, shift, gid);
719 /* Loop over the j charge groups */
720 for (j = 0; (j < nj); j++)
730 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
732 if (iwater == esolSPC && jwater == esolSPC)
734 /* Interaction between two SPC molecules */
737 /* VdW only - only first atoms in each water interact */
738 add_j_to_nblist(vdw, jj0, bLR);
742 #ifdef DISABLE_WATERWATER_NLIST
743 /* Add entries for the three atoms - only do VdW if we need to */
746 add_j_to_nblist(coul, jj0, bLR);
750 add_j_to_nblist(vdwc, jj0, bLR);
752 add_j_to_nblist(coul, jj0+1, bLR);
753 add_j_to_nblist(coul, jj0+2, bLR);
755 /* One entry for the entire water-water interaction */
758 add_j_to_nblist(coul_ww, jj0, bLR);
762 add_j_to_nblist(vdwc_ww, jj0, bLR);
767 else if (iwater == esolTIP4P && jwater == esolTIP4P)
769 /* Interaction between two TIP4p molecules */
772 /* VdW only - only first atoms in each water interact */
773 add_j_to_nblist(vdw, jj0, bLR);
777 #ifdef DISABLE_WATERWATER_NLIST
778 /* Add entries for the four atoms - only do VdW if we need to */
781 add_j_to_nblist(vdw, jj0, bLR);
783 add_j_to_nblist(coul, jj0+1, bLR);
784 add_j_to_nblist(coul, jj0+2, bLR);
785 add_j_to_nblist(coul, jj0+3, bLR);
787 /* One entry for the entire water-water interaction */
790 add_j_to_nblist(coul_ww, jj0, bLR);
794 add_j_to_nblist(vdwc_ww, jj0, bLR);
801 /* j charge group is not water, but i is.
802 * Add entries to the water-other_atom lists; the geometry of the water
803 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
804 * so we don't care if it is SPC or TIP4P...
811 for (jj = jj0; (jj < jj1); jj++)
815 add_j_to_nblist(coul, jj, bLR);
821 for (jj = jj0; (jj < jj1); jj++)
823 if (bHaveVdW[type[jj]])
825 add_j_to_nblist(vdw, jj, bLR);
831 /* _charge_ _groups_ interact with both coulomb and LJ */
832 /* Check which atoms we should add to the lists! */
833 for (jj = jj0; (jj < jj1); jj++)
835 if (bHaveVdW[type[jj]])
839 add_j_to_nblist(vdwc, jj, bLR);
843 add_j_to_nblist(vdw, jj, bLR);
846 else if (charge[jj] != 0)
848 add_j_to_nblist(coul, jj, bLR);
855 close_i_nblist(coul);
856 close_i_nblist(vdwc);
857 #ifndef DISABLE_WATERWATER_NLIST
858 close_i_nblist(coul_ww);
859 close_i_nblist(vdwc_ww);
864 /* no solvent as i charge group */
865 /* Loop over the atoms in the i charge group */
866 for (i = 0; i < nicg; i++)
869 gid = GID(igid, jgid, ngid);
872 /* Create new i_atom for each energy group */
873 if (bDoVdW && bDoCoul)
875 new_i_nblist(vdwc, i_atom, shift, gid);
879 new_i_nblist(vdw, i_atom, shift, gid);
883 new_i_nblist(coul, i_atom, shift, gid);
885 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
886 bDoCoul_i = (bDoCoul && qi != 0);
888 if (bDoVdW_i || bDoCoul_i)
890 /* Loop over the j charge groups */
891 for (j = 0; (j < nj); j++)
895 /* Check for large charge groups */
906 /* Finally loop over the atoms in the j-charge group */
907 for (jj = jj0; jj < jj1; jj++)
909 bNotEx = NOTEXCL(bExcl, i, jj);
917 add_j_to_nblist(coul, jj, bLR);
922 if (bHaveVdW[type[jj]])
924 add_j_to_nblist(vdw, jj, bLR);
929 if (bHaveVdW[type[jj]])
933 add_j_to_nblist(vdwc, jj, bLR);
937 add_j_to_nblist(vdw, jj, bLR);
940 else if (charge[jj] != 0)
942 add_j_to_nblist(coul, jj, bLR);
950 close_i_nblist(coul);
951 close_i_nblist(vdwc);
957 /* we are doing free energy */
958 vdwc_free = &nlist[eNL_VDWQQ_FREE];
959 vdw_free = &nlist[eNL_VDW_FREE];
960 coul_free = &nlist[eNL_QQ_FREE];
961 /* Loop over the atoms in the i charge group */
962 for (i = 0; i < nicg; i++)
965 gid = GID(igid, jgid, ngid);
967 qiB = chargeB[i_atom];
969 /* Create new i_atom for each energy group */
970 if (bDoVdW && bDoCoul)
972 new_i_nblist(vdwc, i_atom, shift, gid);
976 new_i_nblist(vdw, i_atom, shift, gid);
980 new_i_nblist(coul, i_atom, shift, gid);
983 new_i_nblist(vdw_free, i_atom, shift, gid);
984 new_i_nblist(coul_free, i_atom, shift, gid);
985 new_i_nblist(vdwc_free, i_atom, shift, gid);
987 bDoVdW_i = (bDoVdW &&
988 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
989 bDoCoul_i = (bDoCoul && (qi != 0 || qiB != 0));
990 /* For TIP4P the first atom does not have a charge,
991 * but the last three do. So we should still put an atom
992 * without LJ but with charge in the water-atom neighborlist
993 * for a TIP4p i charge group.
994 * For SPC type water the first atom has LJ and charge,
995 * so there is no such problem.
997 if (iwater == esolNO)
999 bDoCoul_i_sol = bDoCoul_i;
1003 bDoCoul_i_sol = bDoCoul;
1006 if (bDoVdW_i || bDoCoul_i_sol)
1008 /* Loop over the j charge groups */
1009 for (j = 0; (j < nj); j++)
1013 /* Check for large charge groups */
1024 /* Finally loop over the atoms in the j-charge group */
1025 bFree = bPert[i_atom];
1026 for (jj = jj0; (jj < jj1); jj++)
1028 bFreeJ = bFree || bPert[jj];
1029 /* Complicated if, because the water H's should also
1030 * see perturbed j-particles
1032 if (iwater == esolNO || i == 0 || bFreeJ)
1034 bNotEx = NOTEXCL(bExcl, i, jj);
1042 if (charge[jj] != 0 || chargeB[jj] != 0)
1044 add_j_to_nblist(coul_free, jj, bLR);
1047 else if (!bDoCoul_i)
1049 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1051 add_j_to_nblist(vdw_free, jj, bLR);
1056 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1058 if (charge[jj] != 0 || chargeB[jj] != 0)
1060 add_j_to_nblist(vdwc_free, jj, bLR);
1064 add_j_to_nblist(vdw_free, jj, bLR);
1067 else if (charge[jj] != 0 || chargeB[jj] != 0)
1069 add_j_to_nblist(coul_free, jj, bLR);
1075 /* This is done whether or not bWater is set */
1076 if (charge[jj] != 0)
1078 add_j_to_nblist(coul, jj, bLR);
1081 else if (!bDoCoul_i_sol)
1083 if (bHaveVdW[type[jj]])
1085 add_j_to_nblist(vdw, jj, bLR);
1090 if (bHaveVdW[type[jj]])
1092 if (charge[jj] != 0)
1094 add_j_to_nblist(vdwc, jj, bLR);
1098 add_j_to_nblist(vdw, jj, bLR);
1101 else if (charge[jj] != 0)
1103 add_j_to_nblist(coul, jj, bLR);
1111 close_i_nblist(vdw);
1112 close_i_nblist(coul);
1113 close_i_nblist(vdwc);
1114 close_i_nblist(vdw_free);
1115 close_i_nblist(coul_free);
1116 close_i_nblist(vdwc_free);
1122 put_in_list_adress(gmx_bool bHaveVdW[],
1138 /* The a[] index has been removed,
1139 * to put it back in i_atom should be a[i0] and jj should be a[jj].
1144 t_nblist * vdwc_adress = NULL;
1145 t_nblist * vdw_adress = NULL;
1146 t_nblist * coul_adress = NULL;
1147 t_nblist * vdwc_ww = NULL;
1148 t_nblist * coul_ww = NULL;
1150 int i, j, jcg, igid, gid, nbl_ind, nbl_ind_adress;
1151 atom_id jj, jj0, jj1, i_atom;
1156 real *charge, *chargeB;
1158 real qi, qiB, qq, rlj;
1159 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
1160 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
1162 gmx_bool j_all_atom;
1164 t_nblist *nlist, *nlist_adress;
1165 gmx_bool bEnergyGroupCG;
1167 /* Copy some pointers */
1168 cginfo = fr->cginfo;
1169 charge = md->chargeA;
1170 chargeB = md->chargeB;
1173 bPert = md->bPerturbed;
1176 /* Get atom range */
1178 nicg = index[icg+1]-i0;
1180 /* Get the i charge group info */
1181 igid = GET_CGINFO_GID(cginfo[icg]);
1183 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
1187 gmx_fatal(FARGS, "AdResS does not support free energy pertubation\n");
1190 /* Unpack pointers to neighbourlist structs */
1191 if (fr->nnblists == 2)
1198 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
1199 nbl_ind_adress = nbl_ind+fr->nnblists/2;
1203 nlist = fr->nblists[nbl_ind].nlist_lr;
1204 nlist_adress = fr->nblists[nbl_ind_adress].nlist_lr;
1208 nlist = fr->nblists[nbl_ind].nlist_sr;
1209 nlist_adress = fr->nblists[nbl_ind_adress].nlist_sr;
1213 vdwc = &nlist[eNL_VDWQQ];
1214 vdw = &nlist[eNL_VDW];
1215 coul = &nlist[eNL_QQ];
1217 vdwc_adress = &nlist_adress[eNL_VDWQQ];
1218 vdw_adress = &nlist_adress[eNL_VDW];
1219 coul_adress = &nlist_adress[eNL_QQ];
1221 /* We do not support solvent optimization with AdResS for now.
1222 For this we would need hybrid solvent-other kernels */
1224 /* no solvent as i charge group */
1225 /* Loop over the atoms in the i charge group */
1226 for (i = 0; i < nicg; i++)
1229 gid = GID(igid, jgid, ngid);
1230 qi = charge[i_atom];
1232 /* Create new i_atom for each energy group */
1233 if (bDoVdW && bDoCoul)
1235 new_i_nblist(vdwc, i_atom, shift, gid);
1236 new_i_nblist(vdwc_adress, i_atom, shift, gid);
1241 new_i_nblist(vdw, i_atom, shift, gid);
1242 new_i_nblist(vdw_adress, i_atom, shift, gid);
1247 new_i_nblist(coul, i_atom, shift, gid);
1248 new_i_nblist(coul_adress, i_atom, shift, gid);
1250 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
1251 bDoCoul_i = (bDoCoul && qi != 0);
1253 /* Here we find out whether the energy groups interaction belong to a
1254 * coarse-grained (vsite) or atomistic interaction. Note that, beacuse
1255 * interactions between coarse-grained and other (atomistic) energygroups
1256 * are excluded automatically by grompp, it is sufficient to check for
1257 * the group id of atom i (igid) */
1258 bEnergyGroupCG = !egp_explicit(fr, igid);
1260 if (bDoVdW_i || bDoCoul_i)
1262 /* Loop over the j charge groups */
1263 for (j = 0; (j < nj); j++)
1267 /* Check for large charge groups */
1278 /* Finally loop over the atoms in the j-charge group */
1279 for (jj = jj0; jj < jj1; jj++)
1281 bNotEx = NOTEXCL(bExcl, i, jj);
1283 /* Now we have to exclude interactions which will be zero
1284 * anyway due to the AdResS weights (in previous implementations
1285 * this was done in the force kernel). This is necessary as
1286 * pure interactions (those with b_hybrid=false, i.e. w_i*w_j==1 or 0)
1287 * are put into neighbour lists which will be passed to the
1288 * standard (optimized) kernels for speed. The interactions with
1289 * b_hybrid=true are placed into the _adress neighbour lists and
1290 * processed by the generic AdResS kernel.
1292 if ( (bEnergyGroupCG &&
1293 wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS ) ||
1294 ( !bEnergyGroupCG && wf[jj] <= GMX_REAL_EPS ) )
1299 b_hybrid = !((wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS) ||
1300 (wf[i_atom] <= GMX_REAL_EPS && wf[jj] <= GMX_REAL_EPS));
1306 if (charge[jj] != 0)
1310 add_j_to_nblist(coul, jj, bLR);
1314 add_j_to_nblist(coul_adress, jj, bLR);
1318 else if (!bDoCoul_i)
1320 if (bHaveVdW[type[jj]])
1324 add_j_to_nblist(vdw, jj, bLR);
1328 add_j_to_nblist(vdw_adress, jj, bLR);
1334 if (bHaveVdW[type[jj]])
1336 if (charge[jj] != 0)
1340 add_j_to_nblist(vdwc, jj, bLR);
1344 add_j_to_nblist(vdwc_adress, jj, bLR);
1351 add_j_to_nblist(vdw, jj, bLR);
1355 add_j_to_nblist(vdw_adress, jj, bLR);
1360 else if (charge[jj] != 0)
1364 add_j_to_nblist(coul, jj, bLR);
1368 add_j_to_nblist(coul_adress, jj, bLR);
1377 close_i_nblist(vdw);
1378 close_i_nblist(coul);
1379 close_i_nblist(vdwc);
1380 close_i_nblist(vdw_adress);
1381 close_i_nblist(coul_adress);
1382 close_i_nblist(vdwc_adress);
1388 put_in_list_qmmm(gmx_bool gmx_unused bHaveVdW[],
1390 t_mdatoms gmx_unused * md,
1400 gmx_bool gmx_unused bDoVdW,
1401 gmx_bool gmx_unused bDoCoul,
1402 int gmx_unused solvent_opt)
1405 int i, j, jcg, igid, gid;
1406 atom_id jj, jj0, jj1, i_atom;
1410 /* Get atom range */
1412 nicg = index[icg+1]-i0;
1414 /* Get the i charge group info */
1415 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1417 coul = &fr->QMMMlist;
1419 /* Loop over atoms in the ith charge group */
1420 for (i = 0; i < nicg; i++)
1423 gid = GID(igid, jgid, ngid);
1424 /* Create new i_atom for each energy group */
1425 new_i_nblist(coul, i_atom, shift, gid);
1427 /* Loop over the j charge groups */
1428 for (j = 0; j < nj; j++)
1432 /* Charge groups cannot have QM and MM atoms simultaneously */
1437 /* Finally loop over the atoms in the j-charge group */
1438 for (jj = jj0; jj < jj1; jj++)
1440 bNotEx = NOTEXCL(bExcl, i, jj);
1443 add_j_to_nblist(coul, jj, bLR);
1448 close_i_nblist(coul);
1453 put_in_list_cg(gmx_bool gmx_unused bHaveVdW[],
1455 t_mdatoms gmx_unused * md,
1465 gmx_bool gmx_unused bDoVdW,
1466 gmx_bool gmx_unused bDoCoul,
1467 int gmx_unused solvent_opt)
1470 int igid, gid, nbl_ind;
1474 cginfo = fr->cginfo[icg];
1476 igid = GET_CGINFO_GID(cginfo);
1477 gid = GID(igid, jgid, ngid);
1479 /* Unpack pointers to neighbourlist structs */
1480 if (fr->nnblists == 1)
1486 nbl_ind = fr->gid2nblists[gid];
1490 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1494 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1497 /* Make a new neighbor list for charge group icg.
1498 * Currently simply one neighbor list is made with LJ and Coulomb.
1499 * If required, zero interactions could be removed here
1500 * or in the force loop.
1502 new_i_nblist(vdwc, index[icg], shift, gid);
1503 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1505 for (j = 0; (j < nj); j++)
1508 /* Skip the icg-icg pairs if all self interactions are excluded */
1509 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1511 /* Here we add the j charge group jcg to the list,
1512 * exclusions are also added to the list.
1514 add_j_to_nblist_cg(vdwc, index[jcg], index[jcg+1], bExcl, icg == jcg, bLR);
1518 close_i_nblist(vdwc);
1521 static void setexcl(atom_id start, atom_id end, t_blocka *excl, gmx_bool b,
1528 for (i = start; i < end; i++)
1530 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1532 SETEXCL(bexcl, i-start, excl->a[k]);
1538 for (i = start; i < end; i++)
1540 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1542 RMEXCL(bexcl, i-start, excl->a[k]);
1548 int calc_naaj(int icg, int cgtot)
1552 if ((cgtot % 2) == 1)
1554 /* Odd number of charge groups, easy */
1555 naaj = 1 + (cgtot/2);
1557 else if ((cgtot % 4) == 0)
1559 /* Multiple of four is hard */
1596 fprintf(log, "naaj=%d\n", naaj);
1602 /************************************************
1604 * S I M P L E C O R E S T U F F
1606 ************************************************/
1608 static real calc_image_tric(rvec xi, rvec xj, matrix box,
1609 rvec b_inv, int *shift)
1611 /* This code assumes that the cut-off is smaller than
1612 * a half times the smallest diagonal element of the box.
1614 const real h25 = 2.5;
1619 /* Compute diff vector */
1620 dz = xj[ZZ] - xi[ZZ];
1621 dy = xj[YY] - xi[YY];
1622 dx = xj[XX] - xi[XX];
1624 /* Perform NINT operation, using trunc operation, therefore
1625 * we first add 2.5 then subtract 2 again
1627 tz = dz*b_inv[ZZ] + h25;
1629 dz -= tz*box[ZZ][ZZ];
1630 dy -= tz*box[ZZ][YY];
1631 dx -= tz*box[ZZ][XX];
1633 ty = dy*b_inv[YY] + h25;
1635 dy -= ty*box[YY][YY];
1636 dx -= ty*box[YY][XX];
1638 tx = dx*b_inv[XX]+h25;
1640 dx -= tx*box[XX][XX];
1642 /* Distance squared */
1643 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1645 *shift = XYZ2IS(tx, ty, tz);
1650 static real calc_image_rect(rvec xi, rvec xj, rvec box_size,
1651 rvec b_inv, int *shift)
1653 const real h15 = 1.5;
1659 /* Compute diff vector */
1660 dx = xj[XX] - xi[XX];
1661 dy = xj[YY] - xi[YY];
1662 dz = xj[ZZ] - xi[ZZ];
1664 /* Perform NINT operation, using trunc operation, therefore
1665 * we first add 1.5 then subtract 1 again
1667 tx = dx*b_inv[XX] + h15;
1668 ty = dy*b_inv[YY] + h15;
1669 tz = dz*b_inv[ZZ] + h15;
1674 /* Correct diff vector for translation */
1675 ddx = tx*box_size[XX] - dx;
1676 ddy = ty*box_size[YY] - dy;
1677 ddz = tz*box_size[ZZ] - dz;
1679 /* Distance squared */
1680 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1682 *shift = XYZ2IS(tx, ty, tz);
1687 static void add_simple(t_ns_buf *nsbuf, int nrj, atom_id cg_j,
1688 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1689 int icg, int jgid, t_block *cgs, t_excl bexcl[],
1690 int shift, t_forcerec *fr, put_in_list_t *put_in_list)
1692 if (nsbuf->nj + nrj > MAX_CG)
1694 put_in_list(bHaveVdW, ngid, md, icg, jgid, nsbuf->ncg, nsbuf->jcg,
1695 cgs->index, bexcl, shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1696 /* Reset buffer contents */
1697 nsbuf->ncg = nsbuf->nj = 0;
1699 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1703 static void ns_inner_tric(rvec x[], int icg, int *i_egp_flags,
1704 int njcg, atom_id jcg[],
1705 matrix box, rvec b_inv, real rcut2,
1706 t_block *cgs, t_ns_buf **ns_buf,
1707 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1708 t_excl bexcl[], t_forcerec *fr,
1709 put_in_list_t *put_in_list)
1713 int *cginfo = fr->cginfo;
1714 atom_id cg_j, *cgindex;
1717 cgindex = cgs->index;
1719 for (j = 0; (j < njcg); j++)
1722 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1723 if (calc_image_tric(x[icg], x[cg_j], box, b_inv, &shift) < rcut2)
1725 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1726 if (!(i_egp_flags[jgid] & EGP_EXCL))
1728 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1729 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1736 static void ns_inner_rect(rvec x[], int icg, int *i_egp_flags,
1737 int njcg, atom_id jcg[],
1738 gmx_bool bBox, rvec box_size, rvec b_inv, real rcut2,
1739 t_block *cgs, t_ns_buf **ns_buf,
1740 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1741 t_excl bexcl[], t_forcerec *fr,
1742 put_in_list_t *put_in_list)
1746 int *cginfo = fr->cginfo;
1747 atom_id cg_j, *cgindex;
1750 cgindex = cgs->index;
1754 for (j = 0; (j < njcg); j++)
1757 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1758 if (calc_image_rect(x[icg], x[cg_j], box_size, b_inv, &shift) < rcut2)
1760 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1761 if (!(i_egp_flags[jgid] & EGP_EXCL))
1763 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1764 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1772 for (j = 0; (j < njcg); j++)
1775 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1776 if ((rcut2 == 0) || (distance2(x[icg], x[cg_j]) < rcut2))
1778 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1779 if (!(i_egp_flags[jgid] & EGP_EXCL))
1781 add_simple(&ns_buf[jgid][CENTRAL], nrj, cg_j,
1782 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, CENTRAL, fr,
1790 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1792 static int ns_simple_core(t_forcerec *fr,
1793 gmx_localtop_t *top,
1795 matrix box, rvec box_size,
1796 t_excl bexcl[], atom_id *aaj,
1797 int ngid, t_ns_buf **ns_buf,
1798 put_in_list_t *put_in_list, gmx_bool bHaveVdW[])
1802 int nsearch, icg, jcg, igid, i0, nri, nn;
1805 /* atom_id *i_atoms; */
1806 t_block *cgs = &(top->cgs);
1807 t_blocka *excl = &(top->excls);
1810 gmx_bool bBox, bTriclinic;
1813 rlist2 = sqr(fr->rlist);
1815 bBox = (fr->ePBC != epbcNONE);
1818 for (m = 0; (m < DIM); m++)
1820 b_inv[m] = divide_err(1.0, box_size[m]);
1822 bTriclinic = TRICLINIC(box);
1829 cginfo = fr->cginfo;
1832 for (icg = fr->cg0; (icg < fr->hcg); icg++)
1835 i0 = cgs->index[icg];
1836 nri = cgs->index[icg+1]-i0;
1837 i_atoms = &(cgs->a[i0]);
1838 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1839 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1841 igid = GET_CGINFO_GID(cginfo[icg]);
1842 i_egp_flags = fr->egp_flags + ngid*igid;
1843 setexcl(cgs->index[icg], cgs->index[icg+1], excl, TRUE, bexcl);
1845 naaj = calc_naaj(icg, cgs->nr);
1848 ns_inner_tric(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1849 box, b_inv, rlist2, cgs, ns_buf,
1850 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1854 ns_inner_rect(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1855 bBox, box_size, b_inv, rlist2, cgs, ns_buf,
1856 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1860 for (nn = 0; (nn < ngid); nn++)
1862 for (k = 0; (k < SHIFTS); k++)
1864 nsbuf = &(ns_buf[nn][k]);
1867 put_in_list(bHaveVdW, ngid, md, icg, nn, nsbuf->ncg, nsbuf->jcg,
1868 cgs->index, bexcl, k, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1869 nsbuf->ncg = nsbuf->nj = 0;
1873 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1874 setexcl(cgs->index[icg], cgs->index[icg+1], excl, FALSE, bexcl);
1876 close_neighbor_lists(fr, FALSE);
1881 /************************************************
1883 * N S 5 G R I D S T U F F
1885 ************************************************/
1887 static inline void get_dx(int Nx, real gridx, real rc2, int xgi, real x,
1888 int *dx0, int *dx1, real *dcx2)
1916 for (i = xgi0; i >= 0; i--)
1918 dcx = (i+1)*gridx-x;
1927 for (i = xgi1; i < Nx; i++)
1940 static inline void get_dx_dd(int Nx, real gridx, real rc2, int xgi, real x,
1941 int ncpddc, int shift_min, int shift_max,
1942 int *g0, int *g1, real *dcx2)
1945 int g_min, g_max, shift_home;
1978 g_min = (shift_min == shift_home ? 0 : ncpddc);
1979 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1986 else if (shift_max < 0)
2001 /* Check one grid cell down */
2002 dcx = ((*g0 - 1) + 1)*gridx - x;
2014 /* Check one grid cell up */
2015 dcx = (*g1 + 1)*gridx - x;
2027 #define sqr(x) ((x)*(x))
2028 #define calc_dx2(XI, YI, ZI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
2029 #define calc_cyl_dx2(XI, YI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
2030 /****************************************************
2032 * F A S T N E I G H B O R S E A R C H I N G
2034 * Optimized neighboursearching routine using grid
2035 * at least 1x1x1, see GROMACS manual
2037 ****************************************************/
2040 static void get_cutoff2(t_forcerec *fr, gmx_bool bDoLongRange,
2041 real *rvdw2, real *rcoul2,
2042 real *rs2, real *rm2, real *rl2)
2044 *rs2 = sqr(fr->rlist);
2046 if (bDoLongRange && fr->bTwinRange)
2048 /* The VdW and elec. LR cut-off's could be different,
2049 * so we can not simply set them to rlistlong.
2051 if (EVDW_MIGHT_BE_ZERO_AT_CUTOFF(fr->vdwtype) &&
2052 fr->rvdw > fr->rlist)
2054 *rvdw2 = sqr(fr->rlistlong);
2058 *rvdw2 = sqr(fr->rvdw);
2060 if (EEL_MIGHT_BE_ZERO_AT_CUTOFF(fr->eeltype) &&
2061 fr->rcoulomb > fr->rlist)
2063 *rcoul2 = sqr(fr->rlistlong);
2067 *rcoul2 = sqr(fr->rcoulomb);
2072 /* Workaround for a gcc -O3 or -ffast-math problem */
2076 *rm2 = min(*rvdw2, *rcoul2);
2077 *rl2 = max(*rvdw2, *rcoul2);
2080 static void init_nsgrid_lists(t_forcerec *fr, int ngid, gmx_ns_t *ns)
2082 real rvdw2, rcoul2, rs2, rm2, rl2;
2085 get_cutoff2(fr, TRUE, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2087 /* Short range buffers */
2088 snew(ns->nl_sr, ngid);
2090 snew(ns->nsr, ngid);
2091 snew(ns->nlr_ljc, ngid);
2092 snew(ns->nlr_one, ngid);
2094 /* Always allocate both list types, since rcoulomb might now change with PME load balancing */
2095 /* Long range VdW and Coul buffers */
2096 snew(ns->nl_lr_ljc, ngid);
2097 /* Long range VdW or Coul only buffers */
2098 snew(ns->nl_lr_one, ngid);
2100 for (j = 0; (j < ngid); j++)
2102 snew(ns->nl_sr[j], MAX_CG);
2103 snew(ns->nl_lr_ljc[j], MAX_CG);
2104 snew(ns->nl_lr_one[j], MAX_CG);
2109 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
2114 static int nsgrid_core(t_commrec *cr, t_forcerec *fr,
2115 matrix box, int ngid,
2116 gmx_localtop_t *top,
2118 t_excl bexcl[], gmx_bool *bExcludeAlleg,
2120 put_in_list_t *put_in_list,
2121 gmx_bool bHaveVdW[],
2122 gmx_bool bDoLongRange, gmx_bool bMakeQMMMnblist)
2125 atom_id **nl_lr_ljc, **nl_lr_one, **nl_sr;
2126 int *nlr_ljc, *nlr_one, *nsr;
2127 gmx_domdec_t *dd = NULL;
2128 t_block *cgs = &(top->cgs);
2129 int *cginfo = fr->cginfo;
2130 /* atom_id *i_atoms,*cgsindex=cgs->index; */
2132 int cell_x, cell_y, cell_z;
2133 int d, tx, ty, tz, dx, dy, dz, cj;
2134 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2135 int zsh_ty, zsh_tx, ysh_tx;
2137 int dx0, dx1, dy0, dy1, dz0, dz1;
2138 int Nx, Ny, Nz, shift = -1, j, nrj, nns, nn = -1;
2139 real gridx, gridy, gridz, grid_x, grid_y, grid_z;
2140 real *dcx2, *dcy2, *dcz2;
2142 int cg0, cg1, icg = -1, cgsnr, i0, igid, nri, naaj, max_jcg;
2143 int jcg0, jcg1, jjcg, cgj0, jgid;
2144 int *grida, *gridnra, *gridind;
2145 gmx_bool rvdw_lt_rcoul, rcoul_lt_rvdw;
2146 rvec xi, *cgcm, grid_offset;
2147 real r2, rs2, rvdw2, rcoul2, rm2, rl2, XI, YI, ZI, dcx, dcy, dcz, tmp1, tmp2;
2149 gmx_bool bDomDec, bTriclinicX, bTriclinicY;
2154 bDomDec = DOMAINDECOMP(cr);
2160 bTriclinicX = ((YY < grid->npbcdim &&
2161 (!bDomDec || dd->nc[YY] == 1) && box[YY][XX] != 0) ||
2162 (ZZ < grid->npbcdim &&
2163 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][XX] != 0));
2164 bTriclinicY = (ZZ < grid->npbcdim &&
2165 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][YY] != 0);
2169 get_cutoff2(fr, bDoLongRange, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2171 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
2172 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
2174 if (bMakeQMMMnblist)
2182 nl_lr_ljc = ns->nl_lr_ljc;
2183 nl_lr_one = ns->nl_lr_one;
2184 nlr_ljc = ns->nlr_ljc;
2185 nlr_one = ns->nlr_one;
2193 gridind = grid->index;
2194 gridnra = grid->nra;
2197 gridx = grid->cell_size[XX];
2198 gridy = grid->cell_size[YY];
2199 gridz = grid->cell_size[ZZ];
2203 copy_rvec(grid->cell_offset, grid_offset);
2204 copy_ivec(grid->ncpddc, ncpddc);
2209 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2210 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
2211 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
2212 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
2213 if (zsh_tx != 0 && ysh_tx != 0)
2215 /* This could happen due to rounding, when both ratios are 0.5 */
2224 /* We only want a list for the test particle */
2233 /* Set the shift range */
2234 for (d = 0; d < DIM; d++)
2238 /* Check if we need periodicity shifts.
2239 * Without PBC or with domain decomposition we don't need them.
2241 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2248 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
2259 /* Loop over charge groups */
2260 for (icg = cg0; (icg < cg1); icg++)
2262 igid = GET_CGINFO_GID(cginfo[icg]);
2263 /* Skip this charge group if all energy groups are excluded! */
2264 if (bExcludeAlleg[igid])
2269 i0 = cgs->index[icg];
2271 if (bMakeQMMMnblist)
2273 /* Skip this charge group if it is not a QM atom while making a
2274 * QM/MM neighbourlist
2276 if (md->bQM[i0] == FALSE)
2278 continue; /* MM particle, go to next particle */
2281 /* Compute the number of charge groups that fall within the control
2284 naaj = calc_naaj(icg, cgsnr);
2291 /* make a normal neighbourlist */
2295 /* Get the j charge-group and dd cell shift ranges */
2296 dd_get_ns_ranges(cr->dd, icg, &jcg0, &jcg1, sh0, sh1);
2301 /* Compute the number of charge groups that fall within the control
2304 naaj = calc_naaj(icg, cgsnr);
2310 /* The i-particle is awlways the test particle,
2311 * so we want all j-particles
2313 max_jcg = cgsnr - 1;
2317 max_jcg = jcg1 - cgsnr;
2322 i_egp_flags = fr->egp_flags + igid*ngid;
2324 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2325 setexcl(i0, cgs->index[icg+1], &top->excls, TRUE, bexcl);
2327 ci2xyz(grid, icg, &cell_x, &cell_y, &cell_z);
2329 /* Changed iicg to icg, DvdS 990115
2330 * (but see consistency check above, DvdS 990330)
2333 fprintf(log, "icg=%5d, naaj=%5d, cell %d %d %d\n",
2334 icg, naaj, cell_x, cell_y, cell_z);
2336 /* Loop over shift vectors in three dimensions */
2337 for (tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
2339 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2340 /* Calculate range of cells in Z direction that have the shift tz */
2341 zgi = cell_z + tz*Nz;
2344 get_dx(Nz, gridz, rl2, zgi, ZI, &dz0, &dz1, dcz2);
2346 get_dx_dd(Nz, gridz, rl2, zgi, ZI-grid_offset[ZZ],
2347 ncpddc[ZZ], sh0[ZZ], sh1[ZZ], &dz0, &dz1, dcz2);
2353 for (ty = -shp[YY]; ty <= shp[YY]; ty++)
2355 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2356 /* Calculate range of cells in Y direction that have the shift ty */
2359 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2363 ygi = cell_y + ty*Ny;
2366 get_dx(Ny, gridy, rl2, ygi, YI, &dy0, &dy1, dcy2);
2368 get_dx_dd(Ny, gridy, rl2, ygi, YI-grid_offset[YY],
2369 ncpddc[YY], sh0[YY], sh1[YY], &dy0, &dy1, dcy2);
2375 for (tx = -shp[XX]; tx <= shp[XX]; tx++)
2377 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2378 /* Calculate range of cells in X direction that have the shift tx */
2381 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2385 xgi = cell_x + tx*Nx;
2388 get_dx(Nx, gridx, rl2, xgi*Nx, XI, &dx0, &dx1, dcx2);
2390 get_dx_dd(Nx, gridx, rl2, xgi, XI-grid_offset[XX],
2391 ncpddc[XX], sh0[XX], sh1[XX], &dx0, &dx1, dcx2);
2397 /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2398 * from the neigbour list as it will not interact */
2399 if (fr->adress_type != eAdressOff)
2401 if (md->wf[cgs->index[icg]] <= GMX_REAL_EPS && egp_explicit(fr, igid))
2406 /* Get shift vector */
2407 shift = XYZ2IS(tx, ty, tz);
2409 range_check(shift, 0, SHIFTS);
2411 for (nn = 0; (nn < ngid); nn++)
2418 fprintf(log, "shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2419 shift, dx0, dx1, dy0, dy1, dz0, dz1);
2420 fprintf(log, "cgcm: %8.3f %8.3f %8.3f\n", cgcm[icg][XX],
2421 cgcm[icg][YY], cgcm[icg][ZZ]);
2422 fprintf(log, "xi: %8.3f %8.3f %8.3f\n", XI, YI, ZI);
2424 for (dx = dx0; (dx <= dx1); dx++)
2426 tmp1 = rl2 - dcx2[dx];
2427 for (dy = dy0; (dy <= dy1); dy++)
2429 tmp2 = tmp1 - dcy2[dy];
2432 for (dz = dz0; (dz <= dz1); dz++)
2434 if (tmp2 > dcz2[dz])
2436 /* Find grid-cell cj in which possible neighbours are */
2437 cj = xyz2ci(Ny, Nz, dx, dy, dz);
2439 /* Check out how many cgs (nrj) there in this cell */
2442 /* Find the offset in the cg list */
2445 /* Check if all j's are out of range so we
2446 * can skip the whole cell.
2447 * Should save some time, especially with DD.
2450 (grida[cgj0] >= max_jcg &&
2451 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2457 for (j = 0; (j < nrj); j++)
2459 jjcg = grida[cgj0+j];
2461 /* check whether this guy is in range! */
2462 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2465 r2 = calc_dx2(XI, YI, ZI, cgcm[jjcg]);
2468 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2469 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2470 /* check energy group exclusions */
2471 if (!(i_egp_flags[jgid] & EGP_EXCL))
2475 if (nsr[jgid] >= MAX_CG)
2477 /* Add to short-range list */
2478 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2479 nsr[jgid], nl_sr[jgid],
2480 cgs->index, /* cgsatoms, */ bexcl,
2481 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2484 nl_sr[jgid][nsr[jgid]++] = jjcg;
2488 if (nlr_ljc[jgid] >= MAX_CG)
2490 /* Add to LJ+coulomb long-range list */
2491 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2492 nlr_ljc[jgid], nl_lr_ljc[jgid], top->cgs.index,
2493 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2496 nl_lr_ljc[jgid][nlr_ljc[jgid]++] = jjcg;
2500 if (nlr_one[jgid] >= MAX_CG)
2502 /* Add to long-range list with only coul, or only LJ */
2503 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2504 nlr_one[jgid], nl_lr_one[jgid], top->cgs.index,
2505 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2508 nl_lr_one[jgid][nlr_one[jgid]++] = jjcg;
2520 /* CHECK whether there is anything left in the buffers */
2521 for (nn = 0; (nn < ngid); nn++)
2525 put_in_list(bHaveVdW, ngid, md, icg, nn, nsr[nn], nl_sr[nn],
2526 cgs->index, /* cgsatoms, */ bexcl,
2527 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2530 if (nlr_ljc[nn] > 0)
2532 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_ljc[nn],
2533 nl_lr_ljc[nn], top->cgs.index,
2534 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2537 if (nlr_one[nn] > 0)
2539 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_one[nn],
2540 nl_lr_one[nn], top->cgs.index,
2541 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2547 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2548 setexcl(cgs->index[icg], cgs->index[icg+1], &top->excls, FALSE, bexcl);
2550 /* No need to perform any left-over force calculations anymore (as we used to do here)
2551 * since we now save the proper long-range lists for later evaluation.
2556 /* Close neighbourlists */
2557 close_neighbor_lists(fr, bMakeQMMMnblist);
2562 void ns_realloc_natoms(gmx_ns_t *ns, int natoms)
2566 if (natoms > ns->nra_alloc)
2568 ns->nra_alloc = over_alloc_dd(natoms);
2569 srenew(ns->bexcl, ns->nra_alloc);
2570 for (i = 0; i < ns->nra_alloc; i++)
2577 void init_ns(FILE *fplog, const t_commrec *cr,
2578 gmx_ns_t *ns, t_forcerec *fr,
2579 const gmx_mtop_t *mtop)
2581 int mt, icg, nr_in_cg, maxcg, i, j, jcg, ngid, ncg;
2585 /* Compute largest charge groups size (# atoms) */
2587 for (mt = 0; mt < mtop->nmoltype; mt++)
2589 cgs = &mtop->moltype[mt].cgs;
2590 for (icg = 0; (icg < cgs->nr); icg++)
2592 nr_in_cg = max(nr_in_cg, (int)(cgs->index[icg+1]-cgs->index[icg]));
2596 /* Verify whether largest charge group is <= max cg.
2597 * This is determined by the type of the local exclusion type
2598 * Exclusions are stored in bits. (If the type is not large
2599 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2601 maxcg = sizeof(t_excl)*8;
2602 if (nr_in_cg > maxcg)
2604 gmx_fatal(FARGS, "Max #atoms in a charge group: %d > %d\n",
2608 ngid = mtop->groups.grps[egcENER].nr;
2609 snew(ns->bExcludeAlleg, ngid);
2610 for (i = 0; i < ngid; i++)
2612 ns->bExcludeAlleg[i] = TRUE;
2613 for (j = 0; j < ngid; j++)
2615 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2617 ns->bExcludeAlleg[i] = FALSE;
2625 ns->grid = init_grid(fplog, fr);
2626 init_nsgrid_lists(fr, ngid, ns);
2631 snew(ns->ns_buf, ngid);
2632 for (i = 0; (i < ngid); i++)
2634 snew(ns->ns_buf[i], SHIFTS);
2636 ncg = ncg_mtop(mtop);
2637 snew(ns->simple_aaj, 2*ncg);
2638 for (jcg = 0; (jcg < ncg); jcg++)
2640 ns->simple_aaj[jcg] = jcg;
2641 ns->simple_aaj[jcg+ncg] = jcg;
2645 /* Create array that determines whether or not atoms have VdW */
2646 snew(ns->bHaveVdW, fr->ntype);
2647 for (i = 0; (i < fr->ntype); i++)
2649 for (j = 0; (j < fr->ntype); j++)
2651 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2653 ((BHAMA(fr->nbfp, fr->ntype, i, j) != 0) ||
2654 (BHAMB(fr->nbfp, fr->ntype, i, j) != 0) ||
2655 (BHAMC(fr->nbfp, fr->ntype, i, j) != 0)) :
2656 ((C6(fr->nbfp, fr->ntype, i, j) != 0) ||
2657 (C12(fr->nbfp, fr->ntype, i, j) != 0))));
2662 pr_bvec(debug, 0, "bHaveVdW", ns->bHaveVdW, fr->ntype, TRUE);
2667 if (!DOMAINDECOMP(cr))
2669 /* This could be reduced with particle decomposition */
2670 ns_realloc_natoms(ns, mtop->natoms);
2673 ns->nblist_initialized = FALSE;
2675 /* nbr list debug dump */
2677 char *ptr = getenv("GMX_DUMP_NL");
2680 ns->dump_nl = strtol(ptr, NULL, 10);
2683 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2694 int search_neighbours(FILE *log, t_forcerec *fr,
2696 gmx_localtop_t *top,
2697 gmx_groups_t *groups,
2699 t_nrnb *nrnb, t_mdatoms *md,
2701 gmx_bool bDoLongRangeNS)
2703 t_block *cgs = &(top->cgs);
2704 rvec box_size, grid_x0, grid_x1;
2706 real min_size, grid_dens;
2710 gmx_bool *i_egp_flags;
2711 int cg_start, cg_end, start, end;
2714 gmx_domdec_zones_t *dd_zones;
2715 put_in_list_t *put_in_list;
2719 /* Set some local variables */
2721 ngid = groups->grps[egcENER].nr;
2723 for (m = 0; (m < DIM); m++)
2725 box_size[m] = box[m][m];
2728 if (fr->ePBC != epbcNONE)
2730 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC, box))
2732 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.");
2736 min_size = min(box_size[XX], min(box_size[YY], box_size[ZZ]));
2737 if (2*fr->rlistlong >= min_size)
2739 gmx_fatal(FARGS, "One of the box diagonal elements has become smaller than twice the cut-off length.");
2744 if (DOMAINDECOMP(cr))
2746 ns_realloc_natoms(ns, cgs->index[cgs->nr]);
2750 /* Reset the neighbourlists */
2751 reset_neighbor_lists(fr, TRUE, TRUE);
2753 if (bGrid && bFillGrid)
2757 if (DOMAINDECOMP(cr))
2759 dd_zones = domdec_zones(cr->dd);
2765 get_nsgrid_boundaries(grid->nboundeddim, box, NULL, NULL, NULL, NULL,
2766 cgs->nr, fr->cg_cm, grid_x0, grid_x1, &grid_dens);
2768 grid_first(log, grid, NULL, NULL, box, grid_x0, grid_x1,
2769 fr->rlistlong, grid_dens);
2773 /* Don't know why this all is... (DvdS 3/99) */
2779 end = (cgs->nr+1)/2;
2782 if (DOMAINDECOMP(cr))
2785 fill_grid(dd_zones, grid, end, -1, end, fr->cg_cm);
2787 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2791 fill_grid(NULL, grid, cgs->nr, fr->cg0, fr->hcg, fr->cg_cm);
2792 grid->icg0 = fr->cg0;
2793 grid->icg1 = fr->hcg;
2803 calc_elemnr(grid, start, end, cgs->nr);
2805 grid_last(grid, start, end, cgs->nr);
2810 print_grid(debug, grid);
2815 /* Set the grid cell index for the test particle only.
2816 * The cell to cg index is not corrected, but that does not matter.
2818 fill_grid(NULL, ns->grid, fr->hcg, fr->hcg-1, fr->hcg, fr->cg_cm);
2822 if (fr->adress_type == eAdressOff)
2824 if (!fr->ns.bCGlist)
2826 put_in_list = put_in_list_at;
2830 put_in_list = put_in_list_cg;
2835 put_in_list = put_in_list_adress;
2842 nsearch = nsgrid_core(cr, fr, box, ngid, top,
2843 grid, ns->bexcl, ns->bExcludeAlleg,
2844 md, put_in_list, ns->bHaveVdW,
2845 bDoLongRangeNS, FALSE);
2847 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2848 * the classical calculation. The charge-charge interaction
2849 * between QM and MM atoms is handled in the QMMM core calculation
2850 * (see QMMM.c). The VDW however, we'd like to compute classically
2851 * and the QM MM atom pairs have just been put in the
2852 * corresponding neighbourlists. in case of QMMM we still need to
2853 * fill a special QMMM neighbourlist that contains all neighbours
2854 * of the QM atoms. If bQMMM is true, this list will now be made:
2856 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
2858 nsearch += nsgrid_core(cr, fr, box, ngid, top,
2859 grid, ns->bexcl, ns->bExcludeAlleg,
2860 md, put_in_list_qmmm, ns->bHaveVdW,
2861 bDoLongRangeNS, TRUE);
2866 nsearch = ns_simple_core(fr, top, md, box, box_size,
2867 ns->bexcl, ns->simple_aaj,
2868 ngid, ns->ns_buf, put_in_list, ns->bHaveVdW);
2876 inc_nrnb(nrnb, eNR_NS, nsearch);
2877 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2882 int natoms_beyond_ns_buffer(t_inputrec *ir, t_forcerec *fr, t_block *cgs,
2883 matrix scale_tot, rvec *x)
2885 int cg0, cg1, cg, a0, a1, a, i, j;
2886 real rint, hbuf2, scale;
2888 gmx_bool bIsotropic;
2893 rint = max(ir->rcoulomb, ir->rvdw);
2894 if (ir->rlist < rint)
2896 gmx_fatal(FARGS, "The neighbor search buffer has negative size: %f nm",
2904 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2906 hbuf2 = sqr(0.5*(ir->rlist - rint));
2907 for (cg = cg0; cg < cg1; cg++)
2909 a0 = cgs->index[cg];
2910 a1 = cgs->index[cg+1];
2911 for (a = a0; a < a1; a++)
2913 if (distance2(cg_cm[cg], x[a]) > hbuf2)
2923 scale = scale_tot[0][0];
2924 for (i = 1; i < DIM; i++)
2926 /* With anisotropic scaling, the original spherical ns volumes become
2927 * ellipsoids. To avoid costly transformations we use the minimum
2928 * eigenvalue of the scaling matrix for determining the buffer size.
2929 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2931 scale = min(scale, scale_tot[i][i]);
2932 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2936 for (j = 0; j < i; j++)
2938 if (scale_tot[i][j] != 0)
2944 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2947 for (cg = cg0; cg < cg1; cg++)
2949 svmul(scale, cg_cm[cg], cgsc);
2950 a0 = cgs->index[cg];
2951 a1 = cgs->index[cg+1];
2952 for (a = a0; a < a1; a++)
2954 if (distance2(cgsc, x[a]) > hbuf2)
2963 /* Anistropic scaling */
2964 for (cg = cg0; cg < cg1; cg++)
2966 /* Since scale_tot contains the transpose of the scaling matrix,
2967 * we need to multiply with the transpose.
2969 tmvmul_ur0(scale_tot, cg_cm[cg], cgsc);
2970 a0 = cgs->index[cg];
2971 a1 = cgs->index[cg+1];
2972 for (a = a0; a < a1; a++)
2974 if (distance2(cgsc, x[a]) > hbuf2)