2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015, 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.
39 #include "gromacs/legacyheaders/ns.h"
49 #include "gromacs/domdec/domdec.h"
50 #include "gromacs/legacyheaders/force.h"
51 #include "gromacs/legacyheaders/macros.h"
52 #include "gromacs/legacyheaders/names.h"
53 #include "gromacs/legacyheaders/network.h"
54 #include "gromacs/legacyheaders/nonbonded.h"
55 #include "gromacs/legacyheaders/nrnb.h"
56 #include "gromacs/legacyheaders/nsgrid.h"
57 #include "gromacs/legacyheaders/txtdump.h"
58 #include "gromacs/legacyheaders/types/commrec.h"
59 #include "gromacs/math/utilities.h"
60 #include "gromacs/math/vec.h"
61 #include "gromacs/pbcutil/ishift.h"
62 #include "gromacs/pbcutil/pbc.h"
63 #include "gromacs/topology/mtop_util.h"
64 #include "gromacs/utility/fatalerror.h"
65 #include "gromacs/utility/smalloc.h"
70 * E X C L U S I O N H A N D L I N G
74 static void SETEXCL_(t_excl e[], atom_id i, atom_id j)
78 static void RMEXCL_(t_excl e[], atom_id i, atom_id j)
80 e[j] = e[j] & ~(1<<i);
82 static gmx_bool ISEXCL_(t_excl e[], atom_id i, atom_id j)
84 return (gmx_bool)(e[j] & (1<<i));
86 static gmx_bool NOTEXCL_(t_excl e[], atom_id i, atom_id j)
88 return !(ISEXCL(e, i, j));
91 #define SETEXCL(e, i, j) (e)[((atom_id) (j))] |= (1<<((atom_id) (i)))
92 #define RMEXCL(e, i, j) (e)[((atom_id) (j))] &= (~(1<<((atom_id) (i))))
93 #define ISEXCL(e, i, j) (gmx_bool) ((e)[((atom_id) (j))] & (1<<((atom_id) (i))))
94 #define NOTEXCL(e, i, j) !(ISEXCL(e, i, j))
98 round_up_to_simd_width(int length, int simd_width)
102 offset = (simd_width > 0) ? length % simd_width : 0;
104 return (offset == 0) ? length : length-offset+simd_width;
106 /************************************************
108 * U T I L I T I E S F O R N S
110 ************************************************/
112 void reallocate_nblist(t_nblist *nl)
116 fprintf(debug, "reallocating neigborlist (ielec=%d, ivdw=%d, igeometry=%d, type=%d), maxnri=%d\n",
117 nl->ielec, nl->ivdw, nl->igeometry, nl->type, nl->maxnri);
119 srenew(nl->iinr, nl->maxnri);
120 if (nl->igeometry == GMX_NBLIST_GEOMETRY_CG_CG)
122 srenew(nl->iinr_end, nl->maxnri);
124 srenew(nl->gid, nl->maxnri);
125 srenew(nl->shift, nl->maxnri);
126 srenew(nl->jindex, nl->maxnri+1);
130 static void init_nblist(FILE *log, t_nblist *nl_sr, t_nblist *nl_lr,
131 int maxsr, int maxlr,
132 int ivdw, int ivdwmod,
133 int ielec, int ielecmod,
134 int igeometry, int type,
135 gmx_bool bElecAndVdwSwitchDiffers)
141 for (i = 0; (i < 2); i++)
143 nl = (i == 0) ? nl_sr : nl_lr;
144 homenr = (i == 0) ? maxsr : maxlr;
152 /* Set coul/vdw in neighborlist, and for the normal loops we determine
153 * an index of which one to call.
156 nl->ivdwmod = ivdwmod;
158 nl->ielecmod = ielecmod;
160 nl->igeometry = igeometry;
162 if (nl->type == GMX_NBLIST_INTERACTION_FREE_ENERGY)
164 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
167 /* This will also set the simd_padding_width field */
168 gmx_nonbonded_set_kernel_pointers( (i == 0) ? log : NULL, nl, bElecAndVdwSwitchDiffers);
170 /* maxnri is influenced by the number of shifts (maximum is 8)
171 * and the number of energy groups.
172 * If it is not enough, nl memory will be reallocated during the run.
173 * 4 seems to be a reasonable factor, which only causes reallocation
174 * during runs with tiny and many energygroups.
176 nl->maxnri = homenr*4;
186 reallocate_nblist(nl);
191 fprintf(debug, "Initiating neighbourlist (ielec=%d, ivdw=%d, type=%d) for %s interactions,\nwith %d SR, %d LR atoms.\n",
192 nl->ielec, nl->ivdw, nl->type, gmx_nblist_geometry_names[nl->igeometry], maxsr, maxlr);
197 void init_neighbor_list(FILE *log, t_forcerec *fr, int homenr)
199 /* Make maxlr tunable! (does not seem to be a big difference though)
200 * This parameter determines the number of i particles in a long range
201 * neighbourlist. Too few means many function calls, too many means
204 int maxsr, maxsr_wat, maxlr, maxlr_wat;
205 int ielec, ivdw, ielecmod, ivdwmod, type;
206 int igeometry_def, igeometry_w, igeometry_ww;
208 gmx_bool bElecAndVdwSwitchDiffers;
211 /* maxsr = homenr-fr->nWatMol*3; */
216 gmx_fatal(FARGS, "%s, %d: Negative number of short range atoms.\n"
217 "Call your GROMACS dealer for assistance.", __FILE__, __LINE__);
219 /* This is just for initial allocation, so we do not reallocate
220 * all the nlist arrays many times in a row.
221 * The numbers seem very accurate, but they are uncritical.
223 maxsr_wat = std::min(fr->nWatMol, (homenr+2)/3);
227 maxlr_wat = std::min(maxsr_wat, maxlr);
231 maxlr = maxlr_wat = 0;
234 /* Determine the values for ielec/ivdw. */
235 ielec = fr->nbkernel_elec_interaction;
236 ivdw = fr->nbkernel_vdw_interaction;
237 ielecmod = fr->nbkernel_elec_modifier;
238 ivdwmod = fr->nbkernel_vdw_modifier;
239 type = GMX_NBLIST_INTERACTION_STANDARD;
240 bElecAndVdwSwitchDiffers = ( (fr->rcoulomb_switch != fr->rvdw_switch) || (fr->rcoulomb != fr->rvdw));
242 fr->ns.bCGlist = (getenv("GMX_NBLISTCG") != 0);
245 igeometry_def = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
249 igeometry_def = GMX_NBLIST_GEOMETRY_CG_CG;
252 fprintf(log, "\nUsing charge-group - charge-group neighbor lists and kernels\n\n");
256 if (fr->solvent_opt == esolTIP4P)
258 igeometry_w = GMX_NBLIST_GEOMETRY_WATER4_PARTICLE;
259 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER4_WATER4;
263 igeometry_w = GMX_NBLIST_GEOMETRY_WATER3_PARTICLE;
264 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER3_WATER3;
267 for (i = 0; i < fr->nnblists; i++)
269 nbl = &(fr->nblists[i]);
271 if ((fr->adress_type != eAdressOff) && (i >= fr->nnblists/2))
273 type = GMX_NBLIST_INTERACTION_ADRESS;
275 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ], &nbl->nlist_lr[eNL_VDWQQ],
276 maxsr, maxlr, ivdw, ivdwmod, ielec, ielecmod, igeometry_def, type, bElecAndVdwSwitchDiffers);
277 init_nblist(log, &nbl->nlist_sr[eNL_VDW], &nbl->nlist_lr[eNL_VDW],
278 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, igeometry_def, type, bElecAndVdwSwitchDiffers);
279 init_nblist(log, &nbl->nlist_sr[eNL_QQ], &nbl->nlist_lr[eNL_QQ],
280 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_def, type, bElecAndVdwSwitchDiffers);
281 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATER], &nbl->nlist_lr[eNL_VDWQQ_WATER],
282 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_w, type, bElecAndVdwSwitchDiffers);
283 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATER], &nbl->nlist_lr[eNL_QQ_WATER],
284 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_w, type, bElecAndVdwSwitchDiffers);
285 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATERWATER], &nbl->nlist_lr[eNL_VDWQQ_WATERWATER],
286 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_ww, type, bElecAndVdwSwitchDiffers);
287 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATERWATER], &nbl->nlist_lr[eNL_QQ_WATERWATER],
288 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_ww, type, bElecAndVdwSwitchDiffers);
290 /* Did we get the solvent loops so we can use optimized water kernels? */
291 if (nbl->nlist_sr[eNL_VDWQQ_WATER].kernelptr_vf == NULL
292 || nbl->nlist_sr[eNL_QQ_WATER].kernelptr_vf == NULL
293 || nbl->nlist_sr[eNL_VDWQQ_WATERWATER].kernelptr_vf == NULL
294 || nbl->nlist_sr[eNL_QQ_WATERWATER].kernelptr_vf == NULL)
296 fr->solvent_opt = esolNO;
299 fprintf(log, "Note: The available nonbonded kernels do not support water optimization - disabling.\n");
303 if (fr->efep != efepNO)
305 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_FREE], &nbl->nlist_lr[eNL_VDWQQ_FREE],
306 maxsr, maxlr, ivdw, ivdwmod, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY, bElecAndVdwSwitchDiffers);
307 init_nblist(log, &nbl->nlist_sr[eNL_VDW_FREE], &nbl->nlist_lr[eNL_VDW_FREE],
308 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY, bElecAndVdwSwitchDiffers);
309 init_nblist(log, &nbl->nlist_sr[eNL_QQ_FREE], &nbl->nlist_lr[eNL_QQ_FREE],
310 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY, bElecAndVdwSwitchDiffers);
314 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
316 init_nblist(log, &fr->QMMMlist, NULL,
317 maxsr, maxlr, 0, 0, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_STANDARD, bElecAndVdwSwitchDiffers);
325 fr->ns.nblist_initialized = TRUE;
328 static void reset_nblist(t_nblist *nl)
338 static void reset_neighbor_lists(t_forcerec *fr, gmx_bool bResetSR, gmx_bool bResetLR)
344 /* only reset the short-range nblist */
345 reset_nblist(&(fr->QMMMlist));
348 for (n = 0; n < fr->nnblists; n++)
350 for (i = 0; i < eNL_NR; i++)
354 reset_nblist( &(fr->nblists[n].nlist_sr[i]) );
358 reset_nblist( &(fr->nblists[n].nlist_lr[i]) );
367 static gmx_inline void new_i_nblist(t_nblist *nlist, atom_id i_atom, int shift, int gid)
369 int nri = nlist->nri;
371 /* Check whether we have to increase the i counter */
373 (nlist->iinr[nri] != i_atom) ||
374 (nlist->shift[nri] != shift) ||
375 (nlist->gid[nri] != gid))
377 /* This is something else. Now see if any entries have
378 * been added in the list of the previous atom.
381 ((nlist->jindex[nri+1] > nlist->jindex[nri]) &&
382 (nlist->gid[nri] != -1)))
384 /* If so increase the counter */
387 if (nlist->nri >= nlist->maxnri)
389 nlist->maxnri += over_alloc_large(nlist->nri);
390 reallocate_nblist(nlist);
393 /* Set the number of neighbours and the atom number */
394 nlist->jindex[nri+1] = nlist->jindex[nri];
395 nlist->iinr[nri] = i_atom;
396 nlist->gid[nri] = gid;
397 nlist->shift[nri] = shift;
401 /* Adding to previous list. First remove possible previous padding */
402 if (nlist->simd_padding_width > 1)
404 while (nlist->nrj > 0 && nlist->jjnr[nlist->nrj-1] < 0)
412 static gmx_inline void close_i_nblist(t_nblist *nlist)
414 int nri = nlist->nri;
419 /* Add elements up to padding. Since we allocate memory in units
420 * of the simd_padding width, we do not have to check for possible
421 * list reallocation here.
423 while ((nlist->nrj % nlist->simd_padding_width) != 0)
425 /* Use -4 here, so we can write forces for 4 atoms before real data */
426 nlist->jjnr[nlist->nrj++] = -4;
428 nlist->jindex[nri+1] = nlist->nrj;
430 len = nlist->nrj - nlist->jindex[nri];
431 /* If there are no j-particles we have to reduce the
432 * number of i-particles again, to prevent errors in the
435 if ((len == 0) && (nlist->nri > 0))
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;
599 atom_id jj, jj0, jj1, i_atom;
604 real *charge, *chargeB;
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 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
677 coul_ww = &nlist[eNL_QQ_WATERWATER];
681 vdwc = &nlist[eNL_VDWQQ];
682 vdw = &nlist[eNL_VDW];
683 coul = &nlist[eNL_QQ];
688 if (iwater != esolNO)
690 /* Loop over the atoms in the i charge group */
692 gid = GID(igid, jgid, ngid);
693 /* Create new i_atom for each energy group */
694 if (bDoCoul && bDoVdW)
696 new_i_nblist(vdwc, i_atom, shift, gid);
697 new_i_nblist(vdwc_ww, i_atom, shift, gid);
701 new_i_nblist(vdw, i_atom, shift, gid);
705 new_i_nblist(coul, i_atom, shift, gid);
706 new_i_nblist(coul_ww, i_atom, shift, gid);
708 /* Loop over the j charge groups */
709 for (j = 0; (j < nj); j++)
719 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
721 if (iwater == esolSPC && jwater == esolSPC)
723 /* Interaction between two SPC molecules */
726 /* VdW only - only first atoms in each water interact */
727 add_j_to_nblist(vdw, jj0, bLR);
731 /* One entry for the entire water-water interaction */
734 add_j_to_nblist(coul_ww, jj0, bLR);
738 add_j_to_nblist(vdwc_ww, jj0, bLR);
742 else if (iwater == esolTIP4P && jwater == esolTIP4P)
744 /* Interaction between two TIP4p molecules */
747 /* VdW only - only first atoms in each water interact */
748 add_j_to_nblist(vdw, jj0, bLR);
752 /* One entry for the entire water-water interaction */
755 add_j_to_nblist(coul_ww, jj0, bLR);
759 add_j_to_nblist(vdwc_ww, jj0, bLR);
765 /* j charge group is not water, but i is.
766 * Add entries to the water-other_atom lists; the geometry of the water
767 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
768 * so we don't care if it is SPC or TIP4P...
775 for (jj = jj0; (jj < jj1); jj++)
779 add_j_to_nblist(coul, jj, bLR);
785 for (jj = jj0; (jj < jj1); jj++)
787 if (bHaveVdW[type[jj]])
789 add_j_to_nblist(vdw, jj, bLR);
795 /* _charge_ _groups_ interact with both coulomb and LJ */
796 /* Check which atoms we should add to the lists! */
797 for (jj = jj0; (jj < jj1); jj++)
799 if (bHaveVdW[type[jj]])
803 add_j_to_nblist(vdwc, jj, bLR);
807 add_j_to_nblist(vdw, jj, bLR);
810 else if (charge[jj] != 0)
812 add_j_to_nblist(coul, jj, bLR);
819 close_i_nblist(coul);
820 close_i_nblist(vdwc);
821 close_i_nblist(coul_ww);
822 close_i_nblist(vdwc_ww);
826 /* no solvent as i charge group */
827 /* Loop over the atoms in the i charge group */
828 for (i = 0; i < nicg; i++)
831 gid = GID(igid, jgid, ngid);
834 /* Create new i_atom for each energy group */
835 if (bDoVdW && bDoCoul)
837 new_i_nblist(vdwc, i_atom, shift, gid);
841 new_i_nblist(vdw, i_atom, shift, gid);
845 new_i_nblist(coul, i_atom, shift, gid);
847 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
848 bDoCoul_i = (bDoCoul && qi != 0);
850 if (bDoVdW_i || bDoCoul_i)
852 /* Loop over the j charge groups */
853 for (j = 0; (j < nj); j++)
857 /* Check for large charge groups */
868 /* Finally loop over the atoms in the j-charge group */
869 for (jj = jj0; jj < jj1; jj++)
871 bNotEx = NOTEXCL(bExcl, i, jj);
879 add_j_to_nblist(coul, jj, bLR);
884 if (bHaveVdW[type[jj]])
886 add_j_to_nblist(vdw, jj, bLR);
891 if (bHaveVdW[type[jj]])
895 add_j_to_nblist(vdwc, jj, bLR);
899 add_j_to_nblist(vdw, jj, bLR);
902 else if (charge[jj] != 0)
904 add_j_to_nblist(coul, jj, bLR);
912 close_i_nblist(coul);
913 close_i_nblist(vdwc);
919 /* we are doing free energy */
920 vdwc_free = &nlist[eNL_VDWQQ_FREE];
921 vdw_free = &nlist[eNL_VDW_FREE];
922 coul_free = &nlist[eNL_QQ_FREE];
923 /* Loop over the atoms in the i charge group */
924 for (i = 0; i < nicg; i++)
927 gid = GID(igid, jgid, ngid);
929 qiB = chargeB[i_atom];
931 /* Create new i_atom for each energy group */
932 if (bDoVdW && bDoCoul)
934 new_i_nblist(vdwc, i_atom, shift, gid);
938 new_i_nblist(vdw, i_atom, shift, gid);
942 new_i_nblist(coul, i_atom, shift, gid);
945 new_i_nblist(vdw_free, i_atom, shift, gid);
946 new_i_nblist(coul_free, i_atom, shift, gid);
947 new_i_nblist(vdwc_free, i_atom, shift, gid);
949 bDoVdW_i = (bDoVdW &&
950 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
951 bDoCoul_i = (bDoCoul && (qi != 0 || qiB != 0));
952 /* For TIP4P the first atom does not have a charge,
953 * but the last three do. So we should still put an atom
954 * without LJ but with charge in the water-atom neighborlist
955 * for a TIP4p i charge group.
956 * For SPC type water the first atom has LJ and charge,
957 * so there is no such problem.
959 if (iwater == esolNO)
961 bDoCoul_i_sol = bDoCoul_i;
965 bDoCoul_i_sol = bDoCoul;
968 if (bDoVdW_i || bDoCoul_i_sol)
970 /* Loop over the j charge groups */
971 for (j = 0; (j < nj); j++)
975 /* Check for large charge groups */
986 /* Finally loop over the atoms in the j-charge group */
987 bFree = bPert[i_atom];
988 for (jj = jj0; (jj < jj1); jj++)
990 bFreeJ = bFree || bPert[jj];
991 /* Complicated if, because the water H's should also
992 * see perturbed j-particles
994 if (iwater == esolNO || i == 0 || bFreeJ)
996 bNotEx = NOTEXCL(bExcl, i, jj);
1004 if (charge[jj] != 0 || chargeB[jj] != 0)
1006 add_j_to_nblist(coul_free, jj, bLR);
1009 else if (!bDoCoul_i)
1011 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1013 add_j_to_nblist(vdw_free, jj, bLR);
1018 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1020 if (charge[jj] != 0 || chargeB[jj] != 0)
1022 add_j_to_nblist(vdwc_free, jj, bLR);
1026 add_j_to_nblist(vdw_free, jj, bLR);
1029 else if (charge[jj] != 0 || chargeB[jj] != 0)
1031 add_j_to_nblist(coul_free, jj, bLR);
1037 /* This is done whether or not bWater is set */
1038 if (charge[jj] != 0)
1040 add_j_to_nblist(coul, jj, bLR);
1043 else if (!bDoCoul_i_sol)
1045 if (bHaveVdW[type[jj]])
1047 add_j_to_nblist(vdw, jj, bLR);
1052 if (bHaveVdW[type[jj]])
1054 if (charge[jj] != 0)
1056 add_j_to_nblist(vdwc, jj, bLR);
1060 add_j_to_nblist(vdw, jj, bLR);
1063 else if (charge[jj] != 0)
1065 add_j_to_nblist(coul, jj, bLR);
1073 close_i_nblist(vdw);
1074 close_i_nblist(coul);
1075 close_i_nblist(vdwc);
1076 close_i_nblist(vdw_free);
1077 close_i_nblist(coul_free);
1078 close_i_nblist(vdwc_free);
1084 put_in_list_adress(gmx_bool bHaveVdW[],
1098 int gmx_unused solvent_opt)
1100 /* The a[] index has been removed,
1101 * to put it back in i_atom should be a[i0] and jj should be a[jj].
1106 t_nblist * vdwc_adress = NULL;
1107 t_nblist * vdw_adress = NULL;
1108 t_nblist * coul_adress = NULL;
1110 int i, j, jcg, igid, gid, nbl_ind, nbl_ind_adress;
1111 atom_id jj, jj0, jj1, i_atom;
1120 gmx_bool bDoVdW_i, bDoCoul_i;
1122 t_nblist *nlist, *nlist_adress;
1123 gmx_bool bEnergyGroupCG;
1125 /* Copy some pointers */
1126 cginfo = fr->cginfo;
1127 charge = md->chargeA;
1131 /* Get atom range */
1133 nicg = index[icg+1]-i0;
1135 /* Get the i charge group info */
1136 igid = GET_CGINFO_GID(cginfo[icg]);
1140 gmx_fatal(FARGS, "AdResS does not support free energy pertubation\n");
1143 /* Unpack pointers to neighbourlist structs */
1144 if (fr->nnblists == 2)
1151 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
1152 nbl_ind_adress = nbl_ind+fr->nnblists/2;
1156 nlist = fr->nblists[nbl_ind].nlist_lr;
1157 nlist_adress = fr->nblists[nbl_ind_adress].nlist_lr;
1161 nlist = fr->nblists[nbl_ind].nlist_sr;
1162 nlist_adress = fr->nblists[nbl_ind_adress].nlist_sr;
1166 vdwc = &nlist[eNL_VDWQQ];
1167 vdw = &nlist[eNL_VDW];
1168 coul = &nlist[eNL_QQ];
1170 vdwc_adress = &nlist_adress[eNL_VDWQQ];
1171 vdw_adress = &nlist_adress[eNL_VDW];
1172 coul_adress = &nlist_adress[eNL_QQ];
1174 /* We do not support solvent optimization with AdResS for now.
1175 For this we would need hybrid solvent-other kernels */
1177 /* no solvent as i charge group */
1178 /* Loop over the atoms in the i charge group */
1179 for (i = 0; i < nicg; i++)
1182 gid = GID(igid, jgid, ngid);
1183 qi = charge[i_atom];
1185 /* Create new i_atom for each energy group */
1186 if (bDoVdW && bDoCoul)
1188 new_i_nblist(vdwc, i_atom, shift, gid);
1189 new_i_nblist(vdwc_adress, i_atom, shift, gid);
1194 new_i_nblist(vdw, i_atom, shift, gid);
1195 new_i_nblist(vdw_adress, i_atom, shift, gid);
1200 new_i_nblist(coul, i_atom, shift, gid);
1201 new_i_nblist(coul_adress, i_atom, shift, gid);
1203 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
1204 bDoCoul_i = (bDoCoul && qi != 0);
1206 /* Here we find out whether the energy groups interaction belong to a
1207 * coarse-grained (vsite) or atomistic interaction. Note that, beacuse
1208 * interactions between coarse-grained and other (atomistic) energygroups
1209 * are excluded automatically by grompp, it is sufficient to check for
1210 * the group id of atom i (igid) */
1211 bEnergyGroupCG = !egp_explicit(fr, igid);
1213 if (bDoVdW_i || bDoCoul_i)
1215 /* Loop over the j charge groups */
1216 for (j = 0; (j < nj); j++)
1220 /* Check for large charge groups */
1231 /* Finally loop over the atoms in the j-charge group */
1232 for (jj = jj0; jj < jj1; jj++)
1234 bNotEx = NOTEXCL(bExcl, i, jj);
1236 /* Now we have to exclude interactions which will be zero
1237 * anyway due to the AdResS weights (in previous implementations
1238 * this was done in the force kernel). This is necessary as
1239 * pure interactions (those with b_hybrid=false, i.e. w_i*w_j==1 or 0)
1240 * are put into neighbour lists which will be passed to the
1241 * standard (optimized) kernels for speed. The interactions with
1242 * b_hybrid=true are placed into the _adress neighbour lists and
1243 * processed by the generic AdResS kernel.
1245 if ( (bEnergyGroupCG &&
1246 wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS ) ||
1247 ( !bEnergyGroupCG && wf[jj] <= GMX_REAL_EPS ) )
1252 b_hybrid = !((wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS) ||
1253 (wf[i_atom] <= GMX_REAL_EPS && wf[jj] <= GMX_REAL_EPS));
1259 if (charge[jj] != 0)
1263 add_j_to_nblist(coul, jj, bLR);
1267 add_j_to_nblist(coul_adress, jj, bLR);
1271 else if (!bDoCoul_i)
1273 if (bHaveVdW[type[jj]])
1277 add_j_to_nblist(vdw, jj, bLR);
1281 add_j_to_nblist(vdw_adress, jj, bLR);
1287 if (bHaveVdW[type[jj]])
1289 if (charge[jj] != 0)
1293 add_j_to_nblist(vdwc, jj, bLR);
1297 add_j_to_nblist(vdwc_adress, jj, bLR);
1304 add_j_to_nblist(vdw, jj, bLR);
1308 add_j_to_nblist(vdw_adress, jj, bLR);
1313 else if (charge[jj] != 0)
1317 add_j_to_nblist(coul, jj, bLR);
1321 add_j_to_nblist(coul_adress, jj, bLR);
1330 close_i_nblist(vdw);
1331 close_i_nblist(coul);
1332 close_i_nblist(vdwc);
1333 close_i_nblist(vdw_adress);
1334 close_i_nblist(coul_adress);
1335 close_i_nblist(vdwc_adress);
1341 put_in_list_qmmm(gmx_bool gmx_unused bHaveVdW[],
1343 t_mdatoms gmx_unused * md,
1353 gmx_bool gmx_unused bDoVdW,
1354 gmx_bool gmx_unused bDoCoul,
1355 int gmx_unused solvent_opt)
1358 int i, j, jcg, igid, gid;
1359 atom_id jj, jj0, jj1, i_atom;
1363 /* Get atom range */
1365 nicg = index[icg+1]-i0;
1367 /* Get the i charge group info */
1368 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1370 coul = &fr->QMMMlist;
1372 /* Loop over atoms in the ith charge group */
1373 for (i = 0; i < nicg; i++)
1376 gid = GID(igid, jgid, ngid);
1377 /* Create new i_atom for each energy group */
1378 new_i_nblist(coul, i_atom, shift, gid);
1380 /* Loop over the j charge groups */
1381 for (j = 0; j < nj; j++)
1385 /* Charge groups cannot have QM and MM atoms simultaneously */
1390 /* Finally loop over the atoms in the j-charge group */
1391 for (jj = jj0; jj < jj1; jj++)
1393 bNotEx = NOTEXCL(bExcl, i, jj);
1396 add_j_to_nblist(coul, jj, bLR);
1401 close_i_nblist(coul);
1406 put_in_list_cg(gmx_bool gmx_unused bHaveVdW[],
1408 t_mdatoms gmx_unused * md,
1418 gmx_bool gmx_unused bDoVdW,
1419 gmx_bool gmx_unused bDoCoul,
1420 int gmx_unused solvent_opt)
1423 int igid, gid, nbl_ind;
1427 cginfo = fr->cginfo[icg];
1429 igid = GET_CGINFO_GID(cginfo);
1430 gid = GID(igid, jgid, ngid);
1432 /* Unpack pointers to neighbourlist structs */
1433 if (fr->nnblists == 1)
1439 nbl_ind = fr->gid2nblists[gid];
1443 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1447 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1450 /* Make a new neighbor list for charge group icg.
1451 * Currently simply one neighbor list is made with LJ and Coulomb.
1452 * If required, zero interactions could be removed here
1453 * or in the force loop.
1455 new_i_nblist(vdwc, index[icg], shift, gid);
1456 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1458 for (j = 0; (j < nj); j++)
1461 /* Skip the icg-icg pairs if all self interactions are excluded */
1462 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1464 /* Here we add the j charge group jcg to the list,
1465 * exclusions are also added to the list.
1467 add_j_to_nblist_cg(vdwc, index[jcg], index[jcg+1], bExcl, icg == jcg, bLR);
1471 close_i_nblist(vdwc);
1474 static void setexcl(atom_id start, atom_id end, t_blocka *excl, gmx_bool b,
1481 for (i = start; i < end; i++)
1483 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1485 SETEXCL(bexcl, i-start, excl->a[k]);
1491 for (i = start; i < end; i++)
1493 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1495 RMEXCL(bexcl, i-start, excl->a[k]);
1501 int calc_naaj(int icg, int cgtot)
1505 if ((cgtot % 2) == 1)
1507 /* Odd number of charge groups, easy */
1508 naaj = 1 + (cgtot/2);
1510 else if ((cgtot % 4) == 0)
1512 /* Multiple of four is hard */
1549 fprintf(log, "naaj=%d\n", naaj);
1555 /************************************************
1557 * S I M P L E C O R E S T U F F
1559 ************************************************/
1561 static real calc_image_tric(rvec xi, rvec xj, matrix box,
1562 rvec b_inv, int *shift)
1564 /* This code assumes that the cut-off is smaller than
1565 * a half times the smallest diagonal element of the box.
1567 const real h25 = 2.5;
1572 /* Compute diff vector */
1573 dz = xj[ZZ] - xi[ZZ];
1574 dy = xj[YY] - xi[YY];
1575 dx = xj[XX] - xi[XX];
1577 /* Perform NINT operation, using trunc operation, therefore
1578 * we first add 2.5 then subtract 2 again
1580 tz = static_cast<int>(dz*b_inv[ZZ] + h25);
1582 dz -= tz*box[ZZ][ZZ];
1583 dy -= tz*box[ZZ][YY];
1584 dx -= tz*box[ZZ][XX];
1586 ty = static_cast<int>(dy*b_inv[YY] + h25);
1588 dy -= ty*box[YY][YY];
1589 dx -= ty*box[YY][XX];
1591 tx = static_cast<int>(dx*b_inv[XX]+h25);
1593 dx -= tx*box[XX][XX];
1595 /* Distance squared */
1596 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1598 *shift = XYZ2IS(tx, ty, tz);
1603 static real calc_image_rect(rvec xi, rvec xj, rvec box_size,
1604 rvec b_inv, int *shift)
1606 const real h15 = 1.5;
1612 /* Compute diff vector */
1613 dx = xj[XX] - xi[XX];
1614 dy = xj[YY] - xi[YY];
1615 dz = xj[ZZ] - xi[ZZ];
1617 /* Perform NINT operation, using trunc operation, therefore
1618 * we first add 1.5 then subtract 1 again
1620 tx = static_cast<int>(dx*b_inv[XX] + h15);
1621 ty = static_cast<int>(dy*b_inv[YY] + h15);
1622 tz = static_cast<int>(dz*b_inv[ZZ] + h15);
1627 /* Correct diff vector for translation */
1628 ddx = tx*box_size[XX] - dx;
1629 ddy = ty*box_size[YY] - dy;
1630 ddz = tz*box_size[ZZ] - dz;
1632 /* Distance squared */
1633 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1635 *shift = XYZ2IS(tx, ty, tz);
1640 static void add_simple(t_ns_buf * nsbuf, int nrj, atom_id cg_j,
1641 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1642 int icg, int jgid, t_block *cgs, t_excl bexcl[],
1643 int shift, t_forcerec *fr, put_in_list_t *put_in_list)
1645 if (nsbuf->nj + nrj > MAX_CG)
1647 put_in_list(bHaveVdW, ngid, md, icg, jgid, nsbuf->ncg, nsbuf->jcg,
1648 cgs->index, bexcl, shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1649 /* Reset buffer contents */
1650 nsbuf->ncg = nsbuf->nj = 0;
1652 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1656 static void ns_inner_tric(rvec x[], int icg, int *i_egp_flags,
1657 int njcg, atom_id jcg[],
1658 matrix box, rvec b_inv, real rcut2,
1659 t_block *cgs, t_ns_buf **ns_buf,
1660 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1661 t_excl bexcl[], t_forcerec *fr,
1662 put_in_list_t *put_in_list)
1666 int *cginfo = fr->cginfo;
1667 atom_id cg_j, *cgindex;
1669 cgindex = cgs->index;
1671 for (j = 0; (j < njcg); j++)
1674 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1675 if (calc_image_tric(x[icg], x[cg_j], box, b_inv, &shift) < rcut2)
1677 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1678 if (!(i_egp_flags[jgid] & EGP_EXCL))
1680 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1681 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1688 static void ns_inner_rect(rvec x[], int icg, int *i_egp_flags,
1689 int njcg, atom_id jcg[],
1690 gmx_bool bBox, rvec box_size, rvec b_inv, real rcut2,
1691 t_block *cgs, t_ns_buf **ns_buf,
1692 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1693 t_excl bexcl[], t_forcerec *fr,
1694 put_in_list_t *put_in_list)
1698 int *cginfo = fr->cginfo;
1699 atom_id cg_j, *cgindex;
1701 cgindex = cgs->index;
1705 for (j = 0; (j < njcg); j++)
1708 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1709 if (calc_image_rect(x[icg], x[cg_j], box_size, b_inv, &shift) < rcut2)
1711 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1712 if (!(i_egp_flags[jgid] & EGP_EXCL))
1714 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1715 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1723 for (j = 0; (j < njcg); j++)
1726 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1727 if ((rcut2 == 0) || (distance2(x[icg], x[cg_j]) < rcut2))
1729 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1730 if (!(i_egp_flags[jgid] & EGP_EXCL))
1732 add_simple(&ns_buf[jgid][CENTRAL], nrj, cg_j,
1733 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, CENTRAL, fr,
1741 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1743 static int ns_simple_core(t_forcerec *fr,
1744 gmx_localtop_t *top,
1746 matrix box, rvec box_size,
1747 t_excl bexcl[], atom_id *aaj,
1748 int ngid, t_ns_buf **ns_buf,
1749 put_in_list_t *put_in_list, gmx_bool bHaveVdW[])
1753 int nsearch, icg, igid, nn;
1756 /* atom_id *i_atoms; */
1757 t_block *cgs = &(top->cgs);
1758 t_blocka *excl = &(top->excls);
1761 gmx_bool bBox, bTriclinic;
1764 rlist2 = sqr(fr->rlist);
1766 bBox = (fr->ePBC != epbcNONE);
1769 for (m = 0; (m < DIM); m++)
1771 b_inv[m] = divide_err(1.0, box_size[m]);
1773 bTriclinic = TRICLINIC(box);
1780 cginfo = fr->cginfo;
1783 for (icg = fr->cg0; (icg < fr->hcg); icg++)
1786 i0 = cgs->index[icg];
1787 nri = cgs->index[icg+1]-i0;
1788 i_atoms = &(cgs->a[i0]);
1789 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1790 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1792 igid = GET_CGINFO_GID(cginfo[icg]);
1793 i_egp_flags = fr->egp_flags + ngid*igid;
1794 setexcl(cgs->index[icg], cgs->index[icg+1], excl, TRUE, bexcl);
1796 naaj = calc_naaj(icg, cgs->nr);
1799 ns_inner_tric(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1800 box, b_inv, rlist2, cgs, ns_buf,
1801 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1805 ns_inner_rect(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1806 bBox, box_size, b_inv, rlist2, cgs, ns_buf,
1807 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1811 for (nn = 0; (nn < ngid); nn++)
1813 for (k = 0; (k < SHIFTS); k++)
1815 nsbuf = &(ns_buf[nn][k]);
1818 put_in_list(bHaveVdW, ngid, md, icg, nn, nsbuf->ncg, nsbuf->jcg,
1819 cgs->index, bexcl, k, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1820 nsbuf->ncg = nsbuf->nj = 0;
1824 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1825 setexcl(cgs->index[icg], cgs->index[icg+1], excl, FALSE, bexcl);
1827 close_neighbor_lists(fr, FALSE);
1832 /************************************************
1834 * N S 5 G R I D S T U F F
1836 ************************************************/
1838 static gmx_inline void get_dx(int Nx, real gridx, real rc2, int xgi, real x,
1839 int *dx0, int *dx1, real *dcx2)
1867 for (i = xgi0; i >= 0; i--)
1869 dcx = (i+1)*gridx-x;
1878 for (i = xgi1; i < Nx; i++)
1891 static gmx_inline void get_dx_dd(int Nx, real gridx, real rc2, int xgi, real x,
1892 int ncpddc, int shift_min, int shift_max,
1893 int *g0, int *g1, real *dcx2)
1896 int g_min, g_max, shift_home;
1929 g_min = (shift_min == shift_home ? 0 : ncpddc);
1930 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1937 else if (shift_max < 0)
1952 /* Check one grid cell down */
1953 dcx = ((*g0 - 1) + 1)*gridx - x;
1965 /* Check one grid cell up */
1966 dcx = (*g1 + 1)*gridx - x;
1978 #define sqr(x) ((x)*(x))
1979 #define calc_dx2(XI, YI, ZI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
1980 #define calc_cyl_dx2(XI, YI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
1981 /****************************************************
1983 * F A S T N E I G H B O R S E A R C H I N G
1985 * Optimized neighboursearching routine using grid
1986 * at least 1x1x1, see GROMACS manual
1988 ****************************************************/
1991 static void get_cutoff2(t_forcerec *fr, gmx_bool bDoLongRange,
1992 real *rvdw2, real *rcoul2,
1993 real *rs2, real *rm2, real *rl2)
1995 *rs2 = sqr(fr->rlist);
1997 if (bDoLongRange && fr->bTwinRange)
1999 /* With plain cut-off or RF we need to make the list exactly
2000 * up to the cut-off and the cut-off's can be different,
2001 * so we can not simply set them to rlistlong.
2002 * To keep this code compatible with (exotic) old cases,
2003 * we also create lists up to rvdw/rcoulomb for PME and Ewald.
2004 * The interaction check should correspond to:
2005 * !ir_vdw/coulomb_might_be_zero_at_cutoff from inputrec.c.
2007 if (((fr->vdwtype == evdwCUT || fr->vdwtype == evdwPME) &&
2008 fr->vdw_modifier == eintmodNONE) ||
2009 fr->rvdw <= fr->rlist)
2011 *rvdw2 = sqr(fr->rvdw);
2015 *rvdw2 = sqr(fr->rlistlong);
2017 if (((fr->eeltype == eelCUT ||
2018 (EEL_RF(fr->eeltype) && fr->eeltype != eelRF_ZERO) ||
2019 fr->eeltype == eelPME ||
2020 fr->eeltype == eelEWALD) &&
2021 fr->coulomb_modifier == eintmodNONE) ||
2022 fr->rcoulomb <= fr->rlist)
2024 *rcoul2 = sqr(fr->rcoulomb);
2028 *rcoul2 = sqr(fr->rlistlong);
2033 /* Workaround for a gcc -O3 or -ffast-math problem */
2037 *rm2 = std::min(*rvdw2, *rcoul2);
2038 *rl2 = std::max(*rvdw2, *rcoul2);
2041 static void init_nsgrid_lists(t_forcerec *fr, int ngid, gmx_ns_t *ns)
2043 real rvdw2, rcoul2, rs2, rm2, rl2;
2046 get_cutoff2(fr, TRUE, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2048 /* Short range buffers */
2049 snew(ns->nl_sr, ngid);
2051 snew(ns->nsr, ngid);
2052 snew(ns->nlr_ljc, ngid);
2053 snew(ns->nlr_one, ngid);
2055 /* Always allocate both list types, since rcoulomb might now change with PME load balancing */
2056 /* Long range VdW and Coul buffers */
2057 snew(ns->nl_lr_ljc, ngid);
2058 /* Long range VdW or Coul only buffers */
2059 snew(ns->nl_lr_one, ngid);
2061 for (j = 0; (j < ngid); j++)
2063 snew(ns->nl_sr[j], MAX_CG);
2064 snew(ns->nl_lr_ljc[j], MAX_CG);
2065 snew(ns->nl_lr_one[j], MAX_CG);
2070 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
2075 static int nsgrid_core(t_commrec *cr, t_forcerec *fr,
2076 matrix box, int ngid,
2077 gmx_localtop_t *top,
2079 t_excl bexcl[], gmx_bool *bExcludeAlleg,
2081 put_in_list_t *put_in_list,
2082 gmx_bool bHaveVdW[],
2083 gmx_bool bDoLongRange, gmx_bool bMakeQMMMnblist)
2086 atom_id **nl_lr_ljc, **nl_lr_one, **nl_sr;
2087 int *nlr_ljc, *nlr_one, *nsr;
2089 t_block *cgs = &(top->cgs);
2090 int *cginfo = fr->cginfo;
2091 /* atom_id *i_atoms,*cgsindex=cgs->index; */
2093 int cell_x, cell_y, cell_z;
2094 int d, tx, ty, tz, dx, dy, dz, cj;
2095 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2096 int zsh_ty, zsh_tx, ysh_tx;
2098 int dx0, dx1, dy0, dy1, dz0, dz1;
2099 int Nx, Ny, Nz, shift = -1, j, nrj, nns, nn = -1;
2100 real gridx, gridy, gridz, grid_x, grid_y;
2101 real *dcx2, *dcy2, *dcz2;
2103 int cg0, cg1, icg = -1, cgsnr, i0, igid, naaj, max_jcg;
2104 int jcg0, jcg1, jjcg, cgj0, jgid;
2105 int *grida, *gridnra, *gridind;
2106 gmx_bool rvdw_lt_rcoul, rcoul_lt_rvdw;
2107 rvec *cgcm, grid_offset;
2108 real r2, rs2, rvdw2, rcoul2, rm2, rl2, XI, YI, ZI, tmp1, tmp2;
2110 gmx_bool bDomDec, bTriclinicX, bTriclinicY;
2115 bDomDec = DOMAINDECOMP(cr);
2118 bTriclinicX = ((YY < grid->npbcdim &&
2119 (!bDomDec || dd->nc[YY] == 1) && box[YY][XX] != 0) ||
2120 (ZZ < grid->npbcdim &&
2121 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][XX] != 0));
2122 bTriclinicY = (ZZ < grid->npbcdim &&
2123 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][YY] != 0);
2127 get_cutoff2(fr, bDoLongRange, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2129 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
2130 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
2132 if (bMakeQMMMnblist)
2140 nl_lr_ljc = ns->nl_lr_ljc;
2141 nl_lr_one = ns->nl_lr_one;
2142 nlr_ljc = ns->nlr_ljc;
2143 nlr_one = ns->nlr_one;
2151 gridind = grid->index;
2152 gridnra = grid->nra;
2155 gridx = grid->cell_size[XX];
2156 gridy = grid->cell_size[YY];
2157 gridz = grid->cell_size[ZZ];
2160 copy_rvec(grid->cell_offset, grid_offset);
2161 copy_ivec(grid->ncpddc, ncpddc);
2166 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2167 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
2168 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
2169 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
2170 if (zsh_tx != 0 && ysh_tx != 0)
2172 /* This could happen due to rounding, when both ratios are 0.5 */
2181 /* We only want a list for the test particle */
2190 /* Set the shift range */
2191 for (d = 0; d < DIM; d++)
2195 /* Check if we need periodicity shifts.
2196 * Without PBC or with domain decomposition we don't need them.
2198 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2205 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < std::sqrt(rl2))
2216 /* Loop over charge groups */
2217 for (icg = cg0; (icg < cg1); icg++)
2219 igid = GET_CGINFO_GID(cginfo[icg]);
2220 /* Skip this charge group if all energy groups are excluded! */
2221 if (bExcludeAlleg[igid])
2226 i0 = cgs->index[icg];
2228 if (bMakeQMMMnblist)
2230 /* Skip this charge group if it is not a QM atom while making a
2231 * QM/MM neighbourlist
2233 if (md->bQM[i0] == FALSE)
2235 continue; /* MM particle, go to next particle */
2238 /* Compute the number of charge groups that fall within the control
2241 naaj = calc_naaj(icg, cgsnr);
2248 /* make a normal neighbourlist */
2252 /* Get the j charge-group and dd cell shift ranges */
2253 dd_get_ns_ranges(cr->dd, icg, &jcg0, &jcg1, sh0, sh1);
2258 /* Compute the number of charge groups that fall within the control
2261 naaj = calc_naaj(icg, cgsnr);
2267 /* The i-particle is awlways the test particle,
2268 * so we want all j-particles
2270 max_jcg = cgsnr - 1;
2274 max_jcg = jcg1 - cgsnr;
2279 i_egp_flags = fr->egp_flags + igid*ngid;
2281 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2282 setexcl(i0, cgs->index[icg+1], &top->excls, TRUE, bexcl);
2284 ci2xyz(grid, icg, &cell_x, &cell_y, &cell_z);
2286 /* Changed iicg to icg, DvdS 990115
2287 * (but see consistency check above, DvdS 990330)
2290 fprintf(log, "icg=%5d, naaj=%5d, cell %d %d %d\n",
2291 icg, naaj, cell_x, cell_y, cell_z);
2293 /* Loop over shift vectors in three dimensions */
2294 for (tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
2296 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2297 /* Calculate range of cells in Z direction that have the shift tz */
2298 zgi = cell_z + tz*Nz;
2301 get_dx(Nz, gridz, rl2, zgi, ZI, &dz0, &dz1, dcz2);
2303 get_dx_dd(Nz, gridz, rl2, zgi, ZI-grid_offset[ZZ],
2304 ncpddc[ZZ], sh0[ZZ], sh1[ZZ], &dz0, &dz1, dcz2);
2310 for (ty = -shp[YY]; ty <= shp[YY]; ty++)
2312 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2313 /* Calculate range of cells in Y direction that have the shift ty */
2316 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2320 ygi = cell_y + ty*Ny;
2323 get_dx(Ny, gridy, rl2, ygi, YI, &dy0, &dy1, dcy2);
2325 get_dx_dd(Ny, gridy, rl2, ygi, YI-grid_offset[YY],
2326 ncpddc[YY], sh0[YY], sh1[YY], &dy0, &dy1, dcy2);
2332 for (tx = -shp[XX]; tx <= shp[XX]; tx++)
2334 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2335 /* Calculate range of cells in X direction that have the shift tx */
2338 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2342 xgi = cell_x + tx*Nx;
2345 get_dx(Nx, gridx, rl2, xgi*Nx, XI, &dx0, &dx1, dcx2);
2347 get_dx_dd(Nx, gridx, rl2, xgi, XI-grid_offset[XX],
2348 ncpddc[XX], sh0[XX], sh1[XX], &dx0, &dx1, dcx2);
2354 /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2355 * from the neigbour list as it will not interact */
2356 if (fr->adress_type != eAdressOff)
2358 if (md->wf[cgs->index[icg]] <= GMX_REAL_EPS && egp_explicit(fr, igid))
2363 /* Get shift vector */
2364 shift = XYZ2IS(tx, ty, tz);
2366 range_check(shift, 0, SHIFTS);
2368 for (nn = 0; (nn < ngid); nn++)
2375 fprintf(log, "shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2376 shift, dx0, dx1, dy0, dy1, dz0, dz1);
2377 fprintf(log, "cgcm: %8.3f %8.3f %8.3f\n", cgcm[icg][XX],
2378 cgcm[icg][YY], cgcm[icg][ZZ]);
2379 fprintf(log, "xi: %8.3f %8.3f %8.3f\n", XI, YI, ZI);
2381 for (dx = dx0; (dx <= dx1); dx++)
2383 tmp1 = rl2 - dcx2[dx];
2384 for (dy = dy0; (dy <= dy1); dy++)
2386 tmp2 = tmp1 - dcy2[dy];
2389 for (dz = dz0; (dz <= dz1); dz++)
2391 if (tmp2 > dcz2[dz])
2393 /* Find grid-cell cj in which possible neighbours are */
2394 cj = xyz2ci(Ny, Nz, dx, dy, dz);
2396 /* Check out how many cgs (nrj) there in this cell */
2399 /* Find the offset in the cg list */
2402 /* Check if all j's are out of range so we
2403 * can skip the whole cell.
2404 * Should save some time, especially with DD.
2407 (grida[cgj0] >= max_jcg &&
2408 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2414 for (j = 0; (j < nrj); j++)
2416 jjcg = grida[cgj0+j];
2418 /* check whether this guy is in range! */
2419 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2422 r2 = calc_dx2(XI, YI, ZI, cgcm[jjcg]);
2425 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2426 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2427 /* check energy group exclusions */
2428 if (!(i_egp_flags[jgid] & EGP_EXCL))
2432 if (nsr[jgid] >= MAX_CG)
2434 /* Add to short-range list */
2435 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2436 nsr[jgid], nl_sr[jgid],
2437 cgs->index, /* cgsatoms, */ bexcl,
2438 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2441 nl_sr[jgid][nsr[jgid]++] = jjcg;
2445 if (nlr_ljc[jgid] >= MAX_CG)
2447 /* Add to LJ+coulomb long-range list */
2448 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2449 nlr_ljc[jgid], nl_lr_ljc[jgid], top->cgs.index,
2450 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2453 nl_lr_ljc[jgid][nlr_ljc[jgid]++] = jjcg;
2457 if (nlr_one[jgid] >= MAX_CG)
2459 /* Add to long-range list with only coul, or only LJ */
2460 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2461 nlr_one[jgid], nl_lr_one[jgid], top->cgs.index,
2462 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2465 nl_lr_one[jgid][nlr_one[jgid]++] = jjcg;
2477 /* CHECK whether there is anything left in the buffers */
2478 for (nn = 0; (nn < ngid); nn++)
2482 put_in_list(bHaveVdW, ngid, md, icg, nn, nsr[nn], nl_sr[nn],
2483 cgs->index, /* cgsatoms, */ bexcl,
2484 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2487 if (nlr_ljc[nn] > 0)
2489 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_ljc[nn],
2490 nl_lr_ljc[nn], top->cgs.index,
2491 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2494 if (nlr_one[nn] > 0)
2496 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_one[nn],
2497 nl_lr_one[nn], top->cgs.index,
2498 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2504 setexcl(cgs->index[icg], cgs->index[icg+1], &top->excls, FALSE, bexcl);
2506 /* No need to perform any left-over force calculations anymore (as we used to do here)
2507 * since we now save the proper long-range lists for later evaluation.
2512 /* Close neighbourlists */
2513 close_neighbor_lists(fr, bMakeQMMMnblist);
2518 void ns_realloc_natoms(gmx_ns_t *ns, int natoms)
2522 if (natoms > ns->nra_alloc)
2524 ns->nra_alloc = over_alloc_dd(natoms);
2525 srenew(ns->bexcl, ns->nra_alloc);
2526 for (i = 0; i < ns->nra_alloc; i++)
2533 void init_ns(FILE *fplog, const t_commrec *cr,
2534 gmx_ns_t *ns, t_forcerec *fr,
2535 const gmx_mtop_t *mtop)
2537 int mt, icg, nr_in_cg, maxcg, i, j, jcg, ngid, ncg;
2540 /* Compute largest charge groups size (# atoms) */
2542 for (mt = 0; mt < mtop->nmoltype; mt++)
2544 cgs = &mtop->moltype[mt].cgs;
2545 for (icg = 0; (icg < cgs->nr); icg++)
2547 nr_in_cg = std::max(nr_in_cg, (int)(cgs->index[icg+1]-cgs->index[icg]));
2551 /* Verify whether largest charge group is <= max cg.
2552 * This is determined by the type of the local exclusion type
2553 * Exclusions are stored in bits. (If the type is not large
2554 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2556 maxcg = sizeof(t_excl)*8;
2557 if (nr_in_cg > maxcg)
2559 gmx_fatal(FARGS, "Max #atoms in a charge group: %d > %d\n",
2563 ngid = mtop->groups.grps[egcENER].nr;
2564 snew(ns->bExcludeAlleg, ngid);
2565 for (i = 0; i < ngid; i++)
2567 ns->bExcludeAlleg[i] = TRUE;
2568 for (j = 0; j < ngid; j++)
2570 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2572 ns->bExcludeAlleg[i] = FALSE;
2580 ns->grid = init_grid(fplog, fr);
2581 init_nsgrid_lists(fr, ngid, ns);
2586 snew(ns->ns_buf, ngid);
2587 for (i = 0; (i < ngid); i++)
2589 snew(ns->ns_buf[i], SHIFTS);
2591 ncg = ncg_mtop(mtop);
2592 snew(ns->simple_aaj, 2*ncg);
2593 for (jcg = 0; (jcg < ncg); jcg++)
2595 ns->simple_aaj[jcg] = jcg;
2596 ns->simple_aaj[jcg+ncg] = jcg;
2600 /* Create array that determines whether or not atoms have VdW */
2601 snew(ns->bHaveVdW, fr->ntype);
2602 for (i = 0; (i < fr->ntype); i++)
2604 for (j = 0; (j < fr->ntype); j++)
2606 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2608 ((BHAMA(fr->nbfp, fr->ntype, i, j) != 0) ||
2609 (BHAMB(fr->nbfp, fr->ntype, i, j) != 0) ||
2610 (BHAMC(fr->nbfp, fr->ntype, i, j) != 0)) :
2611 ((C6(fr->nbfp, fr->ntype, i, j) != 0) ||
2612 (C12(fr->nbfp, fr->ntype, i, j) != 0))));
2617 pr_bvec(debug, 0, "bHaveVdW", ns->bHaveVdW, fr->ntype, TRUE);
2622 if (!DOMAINDECOMP(cr))
2624 ns_realloc_natoms(ns, mtop->natoms);
2627 ns->nblist_initialized = FALSE;
2629 /* nbr list debug dump */
2631 char *ptr = getenv("GMX_DUMP_NL");
2634 ns->dump_nl = strtol(ptr, NULL, 10);
2637 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2648 int search_neighbours(FILE *log, t_forcerec *fr,
2650 gmx_localtop_t *top,
2651 gmx_groups_t *groups,
2653 t_nrnb *nrnb, t_mdatoms *md,
2655 gmx_bool bDoLongRangeNS)
2657 t_block *cgs = &(top->cgs);
2658 rvec box_size, grid_x0, grid_x1;
2660 real min_size, grid_dens;
2666 gmx_domdec_zones_t *dd_zones;
2667 put_in_list_t *put_in_list;
2671 /* Set some local variables */
2673 ngid = groups->grps[egcENER].nr;
2675 for (m = 0; (m < DIM); m++)
2677 box_size[m] = box[m][m];
2680 if (fr->ePBC != epbcNONE)
2682 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC, box))
2684 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.");
2688 min_size = std::min(box_size[XX], std::min(box_size[YY], box_size[ZZ]));
2689 if (2*fr->rlistlong >= min_size)
2691 gmx_fatal(FARGS, "One of the box diagonal elements has become smaller than twice the cut-off length.");
2696 if (DOMAINDECOMP(cr))
2698 ns_realloc_natoms(ns, cgs->index[cgs->nr]);
2702 /* Reset the neighbourlists */
2703 reset_neighbor_lists(fr, TRUE, TRUE);
2705 if (bGrid && bFillGrid)
2709 if (DOMAINDECOMP(cr))
2711 dd_zones = domdec_zones(cr->dd);
2717 get_nsgrid_boundaries(grid->nboundeddim, box, NULL, NULL, NULL, NULL,
2718 cgs->nr, fr->cg_cm, grid_x0, grid_x1, &grid_dens);
2720 grid_first(log, grid, NULL, NULL, box, grid_x0, grid_x1,
2721 fr->rlistlong, grid_dens);
2728 if (DOMAINDECOMP(cr))
2731 fill_grid(dd_zones, grid, end, -1, end, fr->cg_cm);
2733 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2737 fill_grid(NULL, grid, cgs->nr, fr->cg0, fr->hcg, fr->cg_cm);
2738 grid->icg0 = fr->cg0;
2739 grid->icg1 = fr->hcg;
2743 calc_elemnr(grid, start, end, cgs->nr);
2745 grid_last(grid, start, end, cgs->nr);
2750 print_grid(debug, grid);
2755 /* Set the grid cell index for the test particle only.
2756 * The cell to cg index is not corrected, but that does not matter.
2758 fill_grid(NULL, ns->grid, fr->hcg, fr->hcg-1, fr->hcg, fr->cg_cm);
2762 if (fr->adress_type == eAdressOff)
2764 if (!fr->ns.bCGlist)
2766 put_in_list = put_in_list_at;
2770 put_in_list = put_in_list_cg;
2775 put_in_list = put_in_list_adress;
2782 nsearch = nsgrid_core(cr, fr, box, ngid, top,
2783 grid, ns->bexcl, ns->bExcludeAlleg,
2784 md, put_in_list, ns->bHaveVdW,
2785 bDoLongRangeNS, FALSE);
2787 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2788 * the classical calculation. The charge-charge interaction
2789 * between QM and MM atoms is handled in the QMMM core calculation
2790 * (see QMMM.c). The VDW however, we'd like to compute classically
2791 * and the QM MM atom pairs have just been put in the
2792 * corresponding neighbourlists. in case of QMMM we still need to
2793 * fill a special QMMM neighbourlist that contains all neighbours
2794 * of the QM atoms. If bQMMM is true, this list will now be made:
2796 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
2798 nsearch += nsgrid_core(cr, fr, box, ngid, top,
2799 grid, ns->bexcl, ns->bExcludeAlleg,
2800 md, put_in_list_qmmm, ns->bHaveVdW,
2801 bDoLongRangeNS, TRUE);
2806 nsearch = ns_simple_core(fr, top, md, box, box_size,
2807 ns->bexcl, ns->simple_aaj,
2808 ngid, ns->ns_buf, put_in_list, ns->bHaveVdW);
2816 inc_nrnb(nrnb, eNR_NS, nsearch);
2817 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */