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)
69 { e[j] = e[j] | (1<<i); }
70 static void RMEXCL_(t_excl e[],atom_id i,atom_id j)
71 { e[j]=e[j] & ~(1<<i); }
72 static gmx_bool ISEXCL_(t_excl e[],atom_id i,atom_id j)
73 { return (gmx_bool)(e[j] & (1<<i)); }
74 static gmx_bool NOTEXCL_(t_excl e[],atom_id i,atom_id j)
75 { return !(ISEXCL(e,i,j)); }
77 #define SETEXCL(e,i,j) (e)[((atom_id) (j))] |= (1<<((atom_id) (i)))
78 #define RMEXCL(e,i,j) (e)[((atom_id) (j))] &= (~(1<<((atom_id) (i))))
79 #define ISEXCL(e,i,j) (gmx_bool) ((e)[((atom_id) (j))] & (1<<((atom_id) (i))))
80 #define NOTEXCL(e,i,j) !(ISEXCL(e,i,j))
84 round_up_to_simd_width(int length, int simd_width)
88 offset = (simd_width>0) ? length % simd_width : 0;
90 return (offset==0) ? length : length-offset+simd_width;
92 /************************************************
94 * U T I L I T I E S F O R N S
96 ************************************************/
98 static void reallocate_nblist(t_nblist *nl)
102 fprintf(debug,"reallocating neigborlist (ielec=%d, ivdw=%d, igeometry=%d, free_energy=%d), maxnri=%d\n",
103 nl->ielec,nl->ivdw,nl->igeometry,nl->free_energy,nl->maxnri);
105 srenew(nl->iinr, nl->maxnri);
106 if (nl->igeometry == GMX_NBLIST_GEOMETRY_CG_CG)
108 srenew(nl->iinr_end,nl->maxnri);
110 srenew(nl->gid, nl->maxnri);
111 srenew(nl->shift, nl->maxnri);
112 srenew(nl->jindex, nl->maxnri+1);
116 static void init_nblist(FILE *log, t_nblist *nl_sr,t_nblist *nl_lr,
118 int ivdw, int ivdwmod,
119 int ielec, int ielecmod,
120 gmx_bool bfree, int igeometry)
128 nl = (i == 0) ? nl_sr : nl_lr;
129 homenr = (i == 0) ? maxsr : maxlr;
137 /* Set coul/vdw in neighborlist, and for the normal loops we determine
138 * an index of which one to call.
141 nl->ivdwmod = ivdwmod;
143 nl->ielecmod = ielecmod;
144 nl->free_energy = bfree;
145 nl->igeometry = igeometry;
149 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
152 /* This will also set the simd_padding_width field */
153 gmx_nonbonded_set_kernel_pointers( (i==0) ? log : NULL,nl);
155 /* maxnri is influenced by the number of shifts (maximum is 8)
156 * and the number of energy groups.
157 * If it is not enough, nl memory will be reallocated during the run.
158 * 4 seems to be a reasonable factor, which only causes reallocation
159 * during runs with tiny and many energygroups.
161 nl->maxnri = homenr*4;
170 reallocate_nblist(nl);
175 fprintf(debug,"Initiating neighbourlist (ielec=%d, ivdw=%d, free=%d) for %s interactions,\nwith %d SR, %d LR atoms.\n",
176 nl->ielec,nl->ivdw,nl->free_energy,gmx_nblist_geometry_names[nl->igeometry],maxsr,maxlr);
181 void init_neighbor_list(FILE *log,t_forcerec *fr,int homenr)
183 /* Make maxlr tunable! (does not seem to be a big difference though)
184 * This parameter determines the number of i particles in a long range
185 * neighbourlist. Too few means many function calls, too many means
188 int maxsr,maxsr_wat,maxlr,maxlr_wat;
189 int ielec,ielecf,ivdw,ielecmod,ielecmodf,ivdwmod;
191 int igeometry_def,igeometry_w,igeometry_ww;
195 /* maxsr = homenr-fr->nWatMol*3; */
200 gmx_fatal(FARGS,"%s, %d: Negative number of short range atoms.\n"
201 "Call your Gromacs dealer for assistance.",__FILE__,__LINE__);
203 /* This is just for initial allocation, so we do not reallocate
204 * all the nlist arrays many times in a row.
205 * The numbers seem very accurate, but they are uncritical.
207 maxsr_wat = min(fr->nWatMol,(homenr+2)/3);
211 maxlr_wat = min(maxsr_wat,maxlr);
215 maxlr = maxlr_wat = 0;
218 /* Determine the values for ielec/ivdw. */
219 ielec = fr->nbkernel_elec_interaction;
220 ivdw = fr->nbkernel_vdw_interaction;
221 ielecmod = fr->nbkernel_elec_modifier;
222 ivdwmod = fr->nbkernel_vdw_modifier;
224 fr->ns.bCGlist = (getenv("GMX_NBLISTCG") != 0);
227 igeometry_def = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
231 igeometry_def = GMX_NBLIST_GEOMETRY_CG_CG;
234 fprintf(log,"\nUsing charge-group - charge-group neighbor lists and kernels\n\n");
238 if (fr->solvent_opt == esolTIP4P) {
239 igeometry_w = GMX_NBLIST_GEOMETRY_WATER4_PARTICLE;
240 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER4_WATER4;
242 igeometry_w = GMX_NBLIST_GEOMETRY_WATER3_PARTICLE;
243 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER3_WATER3;
246 for(i=0; i<fr->nnblists; i++)
248 nbl = &(fr->nblists[i]);
249 init_nblist(log,&nbl->nlist_sr[eNL_VDWQQ],&nbl->nlist_lr[eNL_VDWQQ],
250 maxsr,maxlr,ivdw,ivdwmod,ielec,ielecmod,FALSE,igeometry_def);
251 init_nblist(log,&nbl->nlist_sr[eNL_VDW],&nbl->nlist_lr[eNL_VDW],
252 maxsr,maxlr,ivdw,ivdwmod,GMX_NBKERNEL_ELEC_NONE,eintmodNONE,FALSE,igeometry_def);
253 init_nblist(log,&nbl->nlist_sr[eNL_QQ],&nbl->nlist_lr[eNL_QQ],
254 maxsr,maxlr,GMX_NBKERNEL_VDW_NONE,eintmodNONE,ielec,ielecmod,FALSE,igeometry_def);
255 init_nblist(log,&nbl->nlist_sr[eNL_VDWQQ_WATER],&nbl->nlist_lr[eNL_VDWQQ_WATER],
256 maxsr_wat,maxlr_wat,ivdw,ivdwmod,ielec,ielecmod, FALSE,igeometry_w);
257 init_nblist(log,&nbl->nlist_sr[eNL_QQ_WATER],&nbl->nlist_lr[eNL_QQ_WATER],
258 maxsr_wat,maxlr_wat,GMX_NBKERNEL_VDW_NONE,eintmodNONE,ielec,ielecmod, FALSE,igeometry_w);
259 init_nblist(log,&nbl->nlist_sr[eNL_VDWQQ_WATERWATER],&nbl->nlist_lr[eNL_VDWQQ_WATERWATER],
260 maxsr_wat,maxlr_wat,ivdw,ivdwmod,ielec,ielecmod, FALSE,igeometry_ww);
261 init_nblist(log,&nbl->nlist_sr[eNL_QQ_WATERWATER],&nbl->nlist_lr[eNL_QQ_WATERWATER],
262 maxsr_wat,maxlr_wat,GMX_NBKERNEL_VDW_NONE,eintmodNONE,ielec,ielecmod, FALSE,igeometry_ww);
264 /* Did we get the solvent loops so we can use optimized water kernels? */
265 if(nbl->nlist_sr[eNL_VDWQQ_WATER].kernelptr_vf==NULL
266 || nbl->nlist_sr[eNL_QQ_WATER].kernelptr_vf==NULL
267 #ifndef DISABLE_WATERWATER_NLIST
268 || nbl->nlist_sr[eNL_VDWQQ_WATERWATER].kernelptr_vf==NULL
269 || nbl->nlist_sr[eNL_QQ_WATERWATER].kernelptr_vf==NULL
273 fr->solvent_opt = esolNO;
274 fprintf(log,"Note: The available nonbonded kernels do not support water optimization - disabling.\n");
277 if (fr->efep != efepNO)
279 if ((fr->bEwald) && (fr->sc_alphacoul > 0)) /* need to handle long range differently if using softcore */
281 ielecf = GMX_NBKERNEL_ELEC_EWALD;
282 ielecmodf = eintmodNONE;
287 ielecmodf = ielecmod;
290 init_nblist(log,&nbl->nlist_sr[eNL_VDWQQ_FREE],&nbl->nlist_lr[eNL_VDWQQ_FREE],
291 maxsr,maxlr,ivdw,ivdwmod,ielecf,ielecmod,TRUE,GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE);
292 init_nblist(log,&nbl->nlist_sr[eNL_VDW_FREE],&nbl->nlist_lr[eNL_VDW_FREE],
293 maxsr,maxlr,ivdw,ivdwmod,GMX_NBKERNEL_ELEC_NONE,eintmodNONE,TRUE,GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE);
294 init_nblist(log,&nbl->nlist_sr[eNL_QQ_FREE],&nbl->nlist_lr[eNL_QQ_FREE],
295 maxsr,maxlr,GMX_NBKERNEL_VDW_NONE,eintmodNONE,ielecf,ielecmod,TRUE,GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE);
299 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
301 init_nblist(log,&fr->QMMMlist,NULL,
302 maxsr,maxlr,0,0,ielec,ielecmod,FALSE,GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE);
310 fr->ns.nblist_initialized=TRUE;
313 static void reset_nblist(t_nblist *nl)
324 static void reset_neighbor_lists(t_forcerec *fr,gmx_bool bResetSR, gmx_bool bResetLR)
330 /* only reset the short-range nblist */
331 reset_nblist(&(fr->QMMMlist));
334 for(n=0; n<fr->nnblists; n++)
336 for(i=0; i<eNL_NR; i++)
340 reset_nblist( &(fr->nblists[n].nlist_sr[i]) );
344 reset_nblist( &(fr->nblists[n].nlist_lr[i]) );
353 static inline void new_i_nblist(t_nblist *nlist,
354 gmx_bool bLR,atom_id i_atom,int shift,int gid)
360 /* Check whether we have to increase the i counter */
362 (nlist->iinr[nri] != i_atom) ||
363 (nlist->shift[nri] != shift) ||
364 (nlist->gid[nri] != gid))
366 /* This is something else. Now see if any entries have
367 * been added in the list of the previous atom.
370 ((nlist->jindex[nri+1] > nlist->jindex[nri]) &&
371 (nlist->gid[nri] != -1)))
373 /* If so increase the counter */
376 if (nlist->nri >= nlist->maxnri)
378 nlist->maxnri += over_alloc_large(nlist->nri);
379 reallocate_nblist(nlist);
382 /* Set the number of neighbours and the atom number */
383 nlist->jindex[nri+1] = nlist->jindex[nri];
384 nlist->iinr[nri] = i_atom;
385 nlist->gid[nri] = gid;
386 nlist->shift[nri] = shift;
390 static inline void close_i_nblist(t_nblist *nlist)
392 int nri = nlist->nri;
397 /* Add elements up to padding. Since we allocate memory in units
398 * of the simd_padding width, we do not have to check for possible
399 * list reallocation here.
401 while((nlist->nrj % nlist->simd_padding_width)!=0)
403 /* Use -4 here, so we can write forces for 4 atoms before real data */
404 nlist->jjnr[nlist->nrj++]=-4;
406 nlist->jindex[nri+1] = nlist->nrj;
408 len=nlist->nrj - nlist->jindex[nri];
410 /* nlist length for water i molecules is treated statically
413 if (len > nlist->maxlen)
420 static inline void close_nblist(t_nblist *nlist)
422 /* Only close this nblist when it has been initialized.
423 * Avoid the creation of i-lists with no j-particles.
427 /* Some assembly kernels do not support empty lists,
428 * make sure here that we don't generate any empty lists.
429 * With the current ns code this branch is taken in two cases:
430 * No i-particles at all: nri=-1 here
431 * There are i-particles, but no j-particles; nri=0 here
437 /* Close list number nri by incrementing the count */
442 static inline void close_neighbor_lists(t_forcerec *fr,gmx_bool bMakeQMMMnblist)
448 close_nblist(&(fr->QMMMlist));
451 for(n=0; n<fr->nnblists; n++)
453 for(i=0; (i<eNL_NR); i++)
455 close_nblist(&(fr->nblists[n].nlist_sr[i]));
456 close_nblist(&(fr->nblists[n].nlist_lr[i]));
462 static inline void add_j_to_nblist(t_nblist *nlist,atom_id j_atom,gmx_bool bLR)
466 if (nlist->nrj >= nlist->maxnrj)
468 nlist->maxnrj = round_up_to_simd_width(over_alloc_small(nlist->nrj + 1),nlist->simd_padding_width);
471 fprintf(debug,"Increasing %s nblist (ielec=%d,ivdw=%d,free=%d,igeometry=%d) j size to %d\n",
472 bLR ? "LR" : "SR",nlist->ielec,nlist->ivdw,nlist->free_energy,nlist->igeometry,nlist->maxnrj);
474 srenew(nlist->jjnr,nlist->maxnrj);
477 nlist->jjnr[nrj] = j_atom;
481 static inline void add_j_to_nblist_cg(t_nblist *nlist,
482 atom_id j_start,int j_end,
483 t_excl *bexcl,gmx_bool i_is_j,
489 if (nlist->nrj >= nlist->maxnrj)
491 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
493 fprintf(debug,"Increasing %s nblist (ielec=%d,ivdw=%d,free=%d,igeometry=%d) j size to %d\n",
494 bLR ? "LR" : "SR",nlist->ielec,nlist->ivdw,nlist->free_energy,nlist->igeometry,nlist->maxnrj);
496 srenew(nlist->jjnr ,nlist->maxnrj);
497 srenew(nlist->jjnr_end,nlist->maxnrj);
498 srenew(nlist->excl ,nlist->maxnrj*MAX_CGCGSIZE);
501 nlist->jjnr[nrj] = j_start;
502 nlist->jjnr_end[nrj] = j_end;
504 if (j_end - j_start > MAX_CGCGSIZE)
506 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);
509 /* Set the exclusions */
510 for(j=j_start; j<j_end; j++)
512 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
516 /* Avoid double counting of intra-cg interactions */
517 for(j=1; j<j_end-j_start; j++)
519 nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
527 put_in_list_t(gmx_bool bHaveVdW[],
544 put_in_list_at(gmx_bool bHaveVdW[],
560 /* The a[] index has been removed,
561 * to put it back in i_atom should be a[i0] and jj should be a[jj].
566 t_nblist * vdwc_free = NULL;
567 t_nblist * vdw_free = NULL;
568 t_nblist * coul_free = NULL;
569 t_nblist * vdwc_ww = NULL;
570 t_nblist * coul_ww = NULL;
572 int i,j,jcg,igid,gid,nbl_ind,ind_ij;
573 atom_id jj,jj0,jj1,i_atom;
578 real *charge,*chargeB;
580 gmx_bool bFreeEnergy,bFree,bFreeJ,bNotEx,*bPert;
581 gmx_bool bDoVdW_i,bDoCoul_i,bDoCoul_i_sol;
585 /* Copy some pointers */
587 charge = md->chargeA;
588 chargeB = md->chargeB;
591 bPert = md->bPerturbed;
595 nicg = index[icg+1]-i0;
597 /* Get the i charge group info */
598 igid = GET_CGINFO_GID(cginfo[icg]);
600 iwater = (solvent_opt!=esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
605 /* Check if any of the particles involved are perturbed.
606 * If not we can do the cheaper normal put_in_list
607 * and use more solvent optimization.
609 for(i=0; i<nicg; i++)
611 bFreeEnergy |= bPert[i0+i];
613 /* Loop over the j charge groups */
614 for(j=0; (j<nj && !bFreeEnergy); j++)
619 /* Finally loop over the atoms in the j-charge group */
620 for(jj=jj0; jj<jj1; jj++)
622 bFreeEnergy |= bPert[jj];
627 /* Unpack pointers to neighbourlist structs */
628 if (fr->nnblists == 1)
634 nbl_ind = fr->gid2nblists[GID(igid,jgid,ngid)];
638 nlist = fr->nblists[nbl_ind].nlist_lr;
642 nlist = fr->nblists[nbl_ind].nlist_sr;
645 if (iwater != esolNO)
647 vdwc = &nlist[eNL_VDWQQ_WATER];
648 vdw = &nlist[eNL_VDW];
649 coul = &nlist[eNL_QQ_WATER];
650 #ifndef DISABLE_WATERWATER_NLIST
651 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
652 coul_ww = &nlist[eNL_QQ_WATERWATER];
657 vdwc = &nlist[eNL_VDWQQ];
658 vdw = &nlist[eNL_VDW];
659 coul = &nlist[eNL_QQ];
664 if (iwater != esolNO)
666 /* Loop over the atoms in the i charge group */
668 gid = GID(igid,jgid,ngid);
669 /* Create new i_atom for each energy group */
670 if (bDoCoul && bDoVdW)
672 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
673 #ifndef DISABLE_WATERWATER_NLIST
674 new_i_nblist(vdwc_ww,bLR,i_atom,shift,gid);
679 new_i_nblist(vdw,bLR,i_atom,shift,gid);
683 new_i_nblist(coul,bLR,i_atom,shift,gid);
684 #ifndef DISABLE_WATERWATER_NLIST
685 new_i_nblist(coul_ww,bLR,i_atom,shift,gid);
688 /* Loop over the j charge groups */
689 for(j=0; (j<nj); j++)
699 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
701 if (iwater == esolSPC && jwater == esolSPC)
703 /* Interaction between two SPC molecules */
706 /* VdW only - only first atoms in each water interact */
707 add_j_to_nblist(vdw,jj0,bLR);
711 #ifdef DISABLE_WATERWATER_NLIST
712 /* Add entries for the three atoms - only do VdW if we need to */
715 add_j_to_nblist(coul,jj0,bLR);
719 add_j_to_nblist(vdwc,jj0,bLR);
721 add_j_to_nblist(coul,jj0+1,bLR);
722 add_j_to_nblist(coul,jj0+2,bLR);
724 /* One entry for the entire water-water interaction */
727 add_j_to_nblist(coul_ww,jj0,bLR);
731 add_j_to_nblist(vdwc_ww,jj0,bLR);
736 else if (iwater == esolTIP4P && jwater == esolTIP4P)
738 /* Interaction between two TIP4p molecules */
741 /* VdW only - only first atoms in each water interact */
742 add_j_to_nblist(vdw,jj0,bLR);
746 #ifdef DISABLE_WATERWATER_NLIST
747 /* Add entries for the four atoms - only do VdW if we need to */
750 add_j_to_nblist(vdw,jj0,bLR);
752 add_j_to_nblist(coul,jj0+1,bLR);
753 add_j_to_nblist(coul,jj0+2,bLR);
754 add_j_to_nblist(coul,jj0+3,bLR);
756 /* One entry for the entire water-water interaction */
759 add_j_to_nblist(coul_ww,jj0,bLR);
763 add_j_to_nblist(vdwc_ww,jj0,bLR);
770 /* j charge group is not water, but i is.
771 * Add entries to the water-other_atom lists; the geometry of the water
772 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
773 * so we don't care if it is SPC or TIP4P...
780 for(jj=jj0; (jj<jj1); jj++)
784 add_j_to_nblist(coul,jj,bLR);
790 for(jj=jj0; (jj<jj1); jj++)
792 if (bHaveVdW[type[jj]])
794 add_j_to_nblist(vdw,jj,bLR);
800 /* _charge_ _groups_ interact with both coulomb and LJ */
801 /* Check which atoms we should add to the lists! */
802 for(jj=jj0; (jj<jj1); jj++)
804 if (bHaveVdW[type[jj]])
808 add_j_to_nblist(vdwc,jj,bLR);
812 add_j_to_nblist(vdw,jj,bLR);
815 else if (charge[jj] != 0)
817 add_j_to_nblist(coul,jj,bLR);
824 close_i_nblist(coul);
825 close_i_nblist(vdwc);
826 #ifndef DISABLE_WATERWATER_NLIST
827 close_i_nblist(coul_ww);
828 close_i_nblist(vdwc_ww);
833 /* no solvent as i charge group */
834 /* Loop over the atoms in the i charge group */
835 for(i=0; i<nicg; i++)
838 gid = GID(igid,jgid,ngid);
841 /* Create new i_atom for each energy group */
842 if (bDoVdW && bDoCoul)
844 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
848 new_i_nblist(vdw,bLR,i_atom,shift,gid);
852 new_i_nblist(coul,bLR,i_atom,shift,gid);
854 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
855 bDoCoul_i = (bDoCoul && qi!=0);
857 if (bDoVdW_i || bDoCoul_i)
859 /* Loop over the j charge groups */
860 for(j=0; (j<nj); j++)
864 /* Check for large charge groups */
875 /* Finally loop over the atoms in the j-charge group */
876 for(jj=jj0; jj<jj1; jj++)
878 bNotEx = NOTEXCL(bExcl,i,jj);
886 add_j_to_nblist(coul,jj,bLR);
891 if (bHaveVdW[type[jj]])
893 add_j_to_nblist(vdw,jj,bLR);
898 if (bHaveVdW[type[jj]])
902 add_j_to_nblist(vdwc,jj,bLR);
906 add_j_to_nblist(vdw,jj,bLR);
909 else if (charge[jj] != 0)
911 add_j_to_nblist(coul,jj,bLR);
919 close_i_nblist(coul);
920 close_i_nblist(vdwc);
926 /* we are doing free energy */
927 vdwc_free = &nlist[eNL_VDWQQ_FREE];
928 vdw_free = &nlist[eNL_VDW_FREE];
929 coul_free = &nlist[eNL_QQ_FREE];
930 /* Loop over the atoms in the i charge group */
931 for(i=0; i<nicg; i++)
934 gid = GID(igid,jgid,ngid);
936 qiB = chargeB[i_atom];
938 /* Create new i_atom for each energy group */
939 if (bDoVdW && bDoCoul)
940 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
942 new_i_nblist(vdw,bLR,i_atom,shift,gid);
944 new_i_nblist(coul,bLR,i_atom,shift,gid);
946 new_i_nblist(vdw_free,bLR,i_atom,shift,gid);
947 new_i_nblist(coul_free,bLR,i_atom,shift,gid);
948 new_i_nblist(vdwc_free,bLR,i_atom,shift,gid);
950 bDoVdW_i = (bDoVdW &&
951 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
952 bDoCoul_i = (bDoCoul && (qi!=0 || qiB!=0));
953 /* For TIP4P the first atom does not have a charge,
954 * but the last three do. So we should still put an atom
955 * without LJ but with charge in the water-atom neighborlist
956 * for a TIP4p i charge group.
957 * For SPC type water the first atom has LJ and charge,
958 * so there is no such problem.
960 if (iwater == esolNO)
962 bDoCoul_i_sol = bDoCoul_i;
966 bDoCoul_i_sol = bDoCoul;
969 if (bDoVdW_i || bDoCoul_i_sol)
971 /* Loop over the j charge groups */
972 for(j=0; (j<nj); j++)
976 /* Check for large charge groups */
987 /* Finally loop over the atoms in the j-charge group */
988 bFree = bPert[i_atom];
989 for(jj=jj0; (jj<jj1); jj++)
991 bFreeJ = bFree || bPert[jj];
992 /* Complicated if, because the water H's should also
993 * see perturbed j-particles
995 if (iwater==esolNO || i==0 || bFreeJ)
997 bNotEx = NOTEXCL(bExcl,i,jj);
1005 if (charge[jj]!=0 || chargeB[jj]!=0)
1007 add_j_to_nblist(coul_free,jj,bLR);
1010 else if (!bDoCoul_i)
1012 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1014 add_j_to_nblist(vdw_free,jj,bLR);
1019 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1021 if (charge[jj]!=0 || chargeB[jj]!=0)
1023 add_j_to_nblist(vdwc_free,jj,bLR);
1027 add_j_to_nblist(vdw_free,jj,bLR);
1030 else if (charge[jj]!=0 || chargeB[jj]!=0)
1031 add_j_to_nblist(coul_free,jj,bLR);
1036 /* This is done whether or not bWater is set */
1037 if (charge[jj] != 0)
1039 add_j_to_nblist(coul,jj,bLR);
1042 else if (!bDoCoul_i_sol)
1044 if (bHaveVdW[type[jj]])
1046 add_j_to_nblist(vdw,jj,bLR);
1051 if (bHaveVdW[type[jj]])
1053 if (charge[jj] != 0)
1055 add_j_to_nblist(vdwc,jj,bLR);
1059 add_j_to_nblist(vdw,jj,bLR);
1062 else if (charge[jj] != 0)
1064 add_j_to_nblist(coul,jj,bLR);
1072 close_i_nblist(vdw);
1073 close_i_nblist(coul);
1074 close_i_nblist(vdwc);
1075 close_i_nblist(vdw_free);
1076 close_i_nblist(coul_free);
1077 close_i_nblist(vdwc_free);
1083 put_in_list_qmmm(gmx_bool bHaveVdW[],
1100 int i,j,jcg,igid,gid;
1101 atom_id jj,jj0,jj1,i_atom;
1105 /* Get atom range */
1107 nicg = index[icg+1]-i0;
1109 /* Get the i charge group info */
1110 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1112 coul = &fr->QMMMlist;
1114 /* Loop over atoms in the ith charge group */
1115 for (i=0;i<nicg;i++)
1118 gid = GID(igid,jgid,ngid);
1119 /* Create new i_atom for each energy group */
1120 new_i_nblist(coul,bLR,i_atom,shift,gid);
1122 /* Loop over the j charge groups */
1127 /* Charge groups cannot have QM and MM atoms simultaneously */
1132 /* Finally loop over the atoms in the j-charge group */
1133 for(jj=jj0; jj<jj1; jj++)
1135 bNotEx = NOTEXCL(bExcl,i,jj);
1137 add_j_to_nblist(coul,jj,bLR);
1141 close_i_nblist(coul);
1146 put_in_list_cg(gmx_bool bHaveVdW[],
1163 int igid,gid,nbl_ind;
1167 cginfo = fr->cginfo[icg];
1169 igid = GET_CGINFO_GID(cginfo);
1170 gid = GID(igid,jgid,ngid);
1172 /* Unpack pointers to neighbourlist structs */
1173 if (fr->nnblists == 1)
1179 nbl_ind = fr->gid2nblists[gid];
1183 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1187 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1190 /* Make a new neighbor list for charge group icg.
1191 * Currently simply one neighbor list is made with LJ and Coulomb.
1192 * If required, zero interactions could be removed here
1193 * or in the force loop.
1195 new_i_nblist(vdwc,bLR,index[icg],shift,gid);
1196 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1198 for(j=0; (j<nj); j++)
1201 /* Skip the icg-icg pairs if all self interactions are excluded */
1202 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1204 /* Here we add the j charge group jcg to the list,
1205 * exclusions are also added to the list.
1207 add_j_to_nblist_cg(vdwc,index[jcg],index[jcg+1],bExcl,icg==jcg,bLR);
1211 close_i_nblist(vdwc);
1214 static void setexcl(atom_id start,atom_id end,t_blocka *excl,gmx_bool b,
1221 for(i=start; i<end; i++)
1223 for(k=excl->index[i]; k<excl->index[i+1]; k++)
1225 SETEXCL(bexcl,i-start,excl->a[k]);
1231 for(i=start; i<end; i++)
1233 for(k=excl->index[i]; k<excl->index[i+1]; k++)
1235 RMEXCL(bexcl,i-start,excl->a[k]);
1241 int calc_naaj(int icg,int cgtot)
1245 if ((cgtot % 2) == 1)
1247 /* Odd number of charge groups, easy */
1248 naaj = 1 + (cgtot/2);
1250 else if ((cgtot % 4) == 0)
1252 /* Multiple of four is hard */
1289 fprintf(log,"naaj=%d\n",naaj);
1295 /************************************************
1297 * S I M P L E C O R E S T U F F
1299 ************************************************/
1301 static real calc_image_tric(rvec xi,rvec xj,matrix box,
1302 rvec b_inv,int *shift)
1304 /* This code assumes that the cut-off is smaller than
1305 * a half times the smallest diagonal element of the box.
1312 /* Compute diff vector */
1313 dz = xj[ZZ] - xi[ZZ];
1314 dy = xj[YY] - xi[YY];
1315 dx = xj[XX] - xi[XX];
1317 /* Perform NINT operation, using trunc operation, therefore
1318 * we first add 2.5 then subtract 2 again
1320 tz = dz*b_inv[ZZ] + h25;
1322 dz -= tz*box[ZZ][ZZ];
1323 dy -= tz*box[ZZ][YY];
1324 dx -= tz*box[ZZ][XX];
1326 ty = dy*b_inv[YY] + h25;
1328 dy -= ty*box[YY][YY];
1329 dx -= ty*box[YY][XX];
1331 tx = dx*b_inv[XX]+h25;
1333 dx -= tx*box[XX][XX];
1335 /* Distance squared */
1336 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1338 *shift = XYZ2IS(tx,ty,tz);
1343 static real calc_image_rect(rvec xi,rvec xj,rvec box_size,
1344 rvec b_inv,int *shift)
1352 /* Compute diff vector */
1353 dx = xj[XX] - xi[XX];
1354 dy = xj[YY] - xi[YY];
1355 dz = xj[ZZ] - xi[ZZ];
1357 /* Perform NINT operation, using trunc operation, therefore
1358 * we first add 1.5 then subtract 1 again
1360 tx = dx*b_inv[XX] + h15;
1361 ty = dy*b_inv[YY] + h15;
1362 tz = dz*b_inv[ZZ] + h15;
1367 /* Correct diff vector for translation */
1368 ddx = tx*box_size[XX] - dx;
1369 ddy = ty*box_size[YY] - dy;
1370 ddz = tz*box_size[ZZ] - dz;
1372 /* Distance squared */
1373 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1375 *shift = XYZ2IS(tx,ty,tz);
1380 static void add_simple(t_ns_buf *nsbuf,int nrj,atom_id cg_j,
1381 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1382 int icg,int jgid,t_block *cgs,t_excl bexcl[],
1383 int shift,t_forcerec *fr,put_in_list_t *put_in_list)
1385 if (nsbuf->nj + nrj > MAX_CG)
1387 put_in_list(bHaveVdW,ngid,md,icg,jgid,nsbuf->ncg,nsbuf->jcg,
1388 cgs->index,bexcl,shift,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
1389 /* Reset buffer contents */
1390 nsbuf->ncg = nsbuf->nj = 0;
1392 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1396 static void ns_inner_tric(rvec x[],int icg,int *i_egp_flags,
1397 int njcg,atom_id jcg[],
1398 matrix box,rvec b_inv,real rcut2,
1399 t_block *cgs,t_ns_buf **ns_buf,
1400 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1401 t_excl bexcl[],t_forcerec *fr,
1402 put_in_list_t *put_in_list)
1406 int *cginfo=fr->cginfo;
1407 atom_id cg_j,*cgindex;
1410 cgindex = cgs->index;
1412 for(j=0; (j<njcg); j++)
1415 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1416 if (calc_image_tric(x[icg],x[cg_j],box,b_inv,&shift) < rcut2)
1418 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1419 if (!(i_egp_flags[jgid] & EGP_EXCL))
1421 add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1422 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1429 static void ns_inner_rect(rvec x[],int icg,int *i_egp_flags,
1430 int njcg,atom_id jcg[],
1431 gmx_bool bBox,rvec box_size,rvec b_inv,real rcut2,
1432 t_block *cgs,t_ns_buf **ns_buf,
1433 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1434 t_excl bexcl[],t_forcerec *fr,
1435 put_in_list_t *put_in_list)
1439 int *cginfo=fr->cginfo;
1440 atom_id cg_j,*cgindex;
1443 cgindex = cgs->index;
1447 for(j=0; (j<njcg); j++)
1450 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1451 if (calc_image_rect(x[icg],x[cg_j],box_size,b_inv,&shift) < rcut2)
1453 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1454 if (!(i_egp_flags[jgid] & EGP_EXCL))
1456 add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1457 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1465 for(j=0; (j<njcg); j++)
1468 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1469 if ((rcut2 == 0) || (distance2(x[icg],x[cg_j]) < rcut2)) {
1470 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1471 if (!(i_egp_flags[jgid] & EGP_EXCL))
1473 add_simple(&ns_buf[jgid][CENTRAL],nrj,cg_j,
1474 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,CENTRAL,fr,
1482 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1484 static int ns_simple_core(t_forcerec *fr,
1485 gmx_localtop_t *top,
1487 matrix box,rvec box_size,
1488 t_excl bexcl[],atom_id *aaj,
1489 int ngid,t_ns_buf **ns_buf,
1490 put_in_list_t *put_in_list,gmx_bool bHaveVdW[])
1494 int nsearch,icg,jcg,igid,i0,nri,nn;
1497 /* atom_id *i_atoms; */
1498 t_block *cgs=&(top->cgs);
1499 t_blocka *excl=&(top->excls);
1502 gmx_bool bBox,bTriclinic;
1505 rlist2 = sqr(fr->rlist);
1507 bBox = (fr->ePBC != epbcNONE);
1510 for(m=0; (m<DIM); m++)
1512 b_inv[m] = divide_err(1.0,box_size[m]);
1514 bTriclinic = TRICLINIC(box);
1521 cginfo = fr->cginfo;
1524 for (icg=fr->cg0; (icg<fr->hcg); icg++)
1527 i0 = cgs->index[icg];
1528 nri = cgs->index[icg+1]-i0;
1529 i_atoms = &(cgs->a[i0]);
1530 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1531 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1533 igid = GET_CGINFO_GID(cginfo[icg]);
1534 i_egp_flags = fr->egp_flags + ngid*igid;
1535 setexcl(cgs->index[icg],cgs->index[icg+1],excl,TRUE,bexcl);
1537 naaj=calc_naaj(icg,cgs->nr);
1540 ns_inner_tric(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1541 box,b_inv,rlist2,cgs,ns_buf,
1542 bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1546 ns_inner_rect(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1547 bBox,box_size,b_inv,rlist2,cgs,ns_buf,
1548 bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1552 for(nn=0; (nn<ngid); nn++)
1554 for(k=0; (k<SHIFTS); k++)
1556 nsbuf = &(ns_buf[nn][k]);
1559 put_in_list(bHaveVdW,ngid,md,icg,nn,nsbuf->ncg,nsbuf->jcg,
1560 cgs->index,bexcl,k,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
1561 nsbuf->ncg=nsbuf->nj=0;
1565 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1566 setexcl(cgs->index[icg],cgs->index[icg+1],excl,FALSE,bexcl);
1568 close_neighbor_lists(fr,FALSE);
1573 /************************************************
1575 * N S 5 G R I D S T U F F
1577 ************************************************/
1579 static inline void get_dx(int Nx,real gridx,real rc2,int xgi,real x,
1580 int *dx0,int *dx1,real *dcx2)
1608 for(i=xgi0; i>=0; i--)
1610 dcx = (i+1)*gridx-x;
1617 for(i=xgi1; i<Nx; i++)
1630 static inline void get_dx_dd(int Nx,real gridx,real rc2,int xgi,real x,
1631 int ncpddc,int shift_min,int shift_max,
1632 int *g0,int *g1,real *dcx2)
1635 int g_min,g_max,shift_home;
1668 g_min = (shift_min == shift_home ? 0 : ncpddc);
1669 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1676 else if (shift_max < 0)
1691 /* Check one grid cell down */
1692 dcx = ((*g0 - 1) + 1)*gridx - x;
1704 /* Check one grid cell up */
1705 dcx = (*g1 + 1)*gridx - x;
1717 #define sqr(x) ((x)*(x))
1718 #define calc_dx2(XI,YI,ZI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
1719 #define calc_cyl_dx2(XI,YI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
1720 /****************************************************
1722 * F A S T N E I G H B O R S E A R C H I N G
1724 * Optimized neighboursearching routine using grid
1725 * at least 1x1x1, see GROMACS manual
1727 ****************************************************/
1730 static void get_cutoff2(t_forcerec *fr,gmx_bool bDoLongRange,
1731 real *rvdw2,real *rcoul2,
1732 real *rs2,real *rm2,real *rl2)
1734 *rs2 = sqr(fr->rlist);
1736 if (bDoLongRange && fr->bTwinRange)
1738 /* The VdW and elec. LR cut-off's could be different,
1739 * so we can not simply set them to rlistlong.
1741 if (EVDW_MIGHT_BE_ZERO_AT_CUTOFF(fr->vdwtype) &&
1742 fr->rvdw > fr->rlist)
1744 *rvdw2 = sqr(fr->rlistlong);
1748 *rvdw2 = sqr(fr->rvdw);
1750 if (EEL_MIGHT_BE_ZERO_AT_CUTOFF(fr->eeltype) &&
1751 fr->rcoulomb > fr->rlist)
1753 *rcoul2 = sqr(fr->rlistlong);
1757 *rcoul2 = sqr(fr->rcoulomb);
1762 /* Workaround for a gcc -O3 or -ffast-math problem */
1766 *rm2 = min(*rvdw2,*rcoul2);
1767 *rl2 = max(*rvdw2,*rcoul2);
1770 static void init_nsgrid_lists(t_forcerec *fr,int ngid,gmx_ns_t *ns)
1772 real rvdw2,rcoul2,rs2,rm2,rl2;
1775 get_cutoff2(fr,TRUE,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
1777 /* Short range buffers */
1778 snew(ns->nl_sr,ngid);
1781 snew(ns->nlr_ljc,ngid);
1782 snew(ns->nlr_one,ngid);
1784 /* Always allocate both list types, since rcoulomb might now change with PME load balancing */
1785 /* Long range VdW and Coul buffers */
1786 snew(ns->nl_lr_ljc,ngid);
1787 /* Long range VdW or Coul only buffers */
1788 snew(ns->nl_lr_one,ngid);
1790 for(j=0; (j<ngid); j++) {
1791 snew(ns->nl_sr[j],MAX_CG);
1792 snew(ns->nl_lr_ljc[j],MAX_CG);
1793 snew(ns->nl_lr_one[j],MAX_CG);
1798 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
1803 static int nsgrid_core(FILE *log,t_commrec *cr,t_forcerec *fr,
1804 matrix box,rvec box_size,int ngid,
1805 gmx_localtop_t *top,
1806 t_grid *grid,rvec x[],
1807 t_excl bexcl[],gmx_bool *bExcludeAlleg,
1808 t_nrnb *nrnb,t_mdatoms *md,
1809 real *lambda,real *dvdlambda,
1810 gmx_grppairener_t *grppener,
1811 put_in_list_t *put_in_list,
1812 gmx_bool bHaveVdW[],
1813 gmx_bool bDoLongRange,gmx_bool bMakeQMMMnblist)
1816 atom_id **nl_lr_ljc,**nl_lr_one,**nl_sr;
1817 int *nlr_ljc,*nlr_one,*nsr;
1818 gmx_domdec_t *dd=NULL;
1819 t_block *cgs=&(top->cgs);
1820 int *cginfo=fr->cginfo;
1821 /* atom_id *i_atoms,*cgsindex=cgs->index; */
1823 int cell_x,cell_y,cell_z;
1824 int d,tx,ty,tz,dx,dy,dz,cj;
1825 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
1826 int zsh_ty,zsh_tx,ysh_tx;
1828 int dx0,dx1,dy0,dy1,dz0,dz1;
1829 int Nx,Ny,Nz,shift=-1,j,nrj,nns,nn=-1;
1830 real gridx,gridy,gridz,grid_x,grid_y,grid_z;
1831 real *dcx2,*dcy2,*dcz2;
1833 int cg0,cg1,icg=-1,cgsnr,i0,igid,nri,naaj,max_jcg;
1834 int jcg0,jcg1,jjcg,cgj0,jgid;
1835 int *grida,*gridnra,*gridind;
1836 gmx_bool rvdw_lt_rcoul,rcoul_lt_rvdw;
1837 rvec xi,*cgcm,grid_offset;
1838 real r2,rs2,rvdw2,rcoul2,rm2,rl2,XI,YI,ZI,dcx,dcy,dcz,tmp1,tmp2;
1840 gmx_bool bDomDec,bTriclinicX,bTriclinicY;
1845 bDomDec = DOMAINDECOMP(cr);
1851 bTriclinicX = ((YY < grid->npbcdim &&
1852 (!bDomDec || dd->nc[YY]==1) && box[YY][XX] != 0) ||
1853 (ZZ < grid->npbcdim &&
1854 (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][XX] != 0));
1855 bTriclinicY = (ZZ < grid->npbcdim &&
1856 (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][YY] != 0);
1860 get_cutoff2(fr,bDoLongRange,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
1862 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
1863 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
1865 if (bMakeQMMMnblist)
1873 nl_lr_ljc = ns->nl_lr_ljc;
1874 nl_lr_one = ns->nl_lr_one;
1875 nlr_ljc = ns->nlr_ljc;
1876 nlr_one = ns->nlr_one;
1884 gridind = grid->index;
1885 gridnra = grid->nra;
1888 gridx = grid->cell_size[XX];
1889 gridy = grid->cell_size[YY];
1890 gridz = grid->cell_size[ZZ];
1894 copy_rvec(grid->cell_offset,grid_offset);
1895 copy_ivec(grid->ncpddc,ncpddc);
1900 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
1901 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
1902 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
1903 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
1904 if (zsh_tx!=0 && ysh_tx!=0)
1906 /* This could happen due to rounding, when both ratios are 0.5 */
1915 /* We only want a list for the test particle */
1924 /* Set the shift range */
1925 for(d=0; d<DIM; d++)
1929 /* Check if we need periodicity shifts.
1930 * Without PBC or with domain decomposition we don't need them.
1932 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
1939 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
1950 /* Loop over charge groups */
1951 for(icg=cg0; (icg < cg1); icg++)
1953 igid = GET_CGINFO_GID(cginfo[icg]);
1954 /* Skip this charge group if all energy groups are excluded! */
1955 if (bExcludeAlleg[igid])
1960 i0 = cgs->index[icg];
1962 if (bMakeQMMMnblist)
1964 /* Skip this charge group if it is not a QM atom while making a
1965 * QM/MM neighbourlist
1967 if (md->bQM[i0]==FALSE)
1969 continue; /* MM particle, go to next particle */
1972 /* Compute the number of charge groups that fall within the control
1975 naaj = calc_naaj(icg,cgsnr);
1982 /* make a normal neighbourlist */
1986 /* Get the j charge-group and dd cell shift ranges */
1987 dd_get_ns_ranges(cr->dd,icg,&jcg0,&jcg1,sh0,sh1);
1992 /* Compute the number of charge groups that fall within the control
1995 naaj = calc_naaj(icg,cgsnr);
2001 /* The i-particle is awlways the test particle,
2002 * so we want all j-particles
2004 max_jcg = cgsnr - 1;
2008 max_jcg = jcg1 - cgsnr;
2013 i_egp_flags = fr->egp_flags + igid*ngid;
2015 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2016 setexcl(i0,cgs->index[icg+1],&top->excls,TRUE,bexcl);
2018 ci2xyz(grid,icg,&cell_x,&cell_y,&cell_z);
2020 /* Changed iicg to icg, DvdS 990115
2021 * (but see consistency check above, DvdS 990330)
2024 fprintf(log,"icg=%5d, naaj=%5d, cell %d %d %d\n",
2025 icg,naaj,cell_x,cell_y,cell_z);
2027 /* Loop over shift vectors in three dimensions */
2028 for (tz=-shp[ZZ]; tz<=shp[ZZ]; tz++)
2030 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2031 /* Calculate range of cells in Z direction that have the shift tz */
2032 zgi = cell_z + tz*Nz;
2035 get_dx(Nz,gridz,rl2,zgi,ZI,&dz0,&dz1,dcz2);
2037 get_dx_dd(Nz,gridz,rl2,zgi,ZI-grid_offset[ZZ],
2038 ncpddc[ZZ],sh0[ZZ],sh1[ZZ],&dz0,&dz1,dcz2);
2044 for (ty=-shp[YY]; ty<=shp[YY]; ty++)
2046 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2047 /* Calculate range of cells in Y direction that have the shift ty */
2050 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2054 ygi = cell_y + ty*Ny;
2057 get_dx(Ny,gridy,rl2,ygi,YI,&dy0,&dy1,dcy2);
2059 get_dx_dd(Ny,gridy,rl2,ygi,YI-grid_offset[YY],
2060 ncpddc[YY],sh0[YY],sh1[YY],&dy0,&dy1,dcy2);
2066 for (tx=-shp[XX]; tx<=shp[XX]; tx++)
2068 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2069 /* Calculate range of cells in X direction that have the shift tx */
2072 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2076 xgi = cell_x + tx*Nx;
2079 get_dx(Nx,gridx,rl2,xgi*Nx,XI,&dx0,&dx1,dcx2);
2081 get_dx_dd(Nx,gridx,rl2,xgi,XI-grid_offset[XX],
2082 ncpddc[XX],sh0[XX],sh1[XX],&dx0,&dx1,dcx2);
2088 /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2089 * from the neigbour list as it will not interact */
2090 if (fr->adress_type != eAdressOff){
2091 if (md->wf[cgs->index[icg]]==0 && egp_explicit(fr, igid)){
2095 /* Get shift vector */
2096 shift=XYZ2IS(tx,ty,tz);
2098 range_check(shift,0,SHIFTS);
2100 for(nn=0; (nn<ngid); nn++)
2107 fprintf(log,"shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2108 shift,dx0,dx1,dy0,dy1,dz0,dz1);
2109 fprintf(log,"cgcm: %8.3f %8.3f %8.3f\n",cgcm[icg][XX],
2110 cgcm[icg][YY],cgcm[icg][ZZ]);
2111 fprintf(log,"xi: %8.3f %8.3f %8.3f\n",XI,YI,ZI);
2113 for (dx=dx0; (dx<=dx1); dx++)
2115 tmp1 = rl2 - dcx2[dx];
2116 for (dy=dy0; (dy<=dy1); dy++)
2118 tmp2 = tmp1 - dcy2[dy];
2121 for (dz=dz0; (dz<=dz1); dz++) {
2122 if (tmp2 > dcz2[dz]) {
2123 /* Find grid-cell cj in which possible neighbours are */
2124 cj = xyz2ci(Ny,Nz,dx,dy,dz);
2126 /* Check out how many cgs (nrj) there in this cell */
2129 /* Find the offset in the cg list */
2132 /* Check if all j's are out of range so we
2133 * can skip the whole cell.
2134 * Should save some time, especially with DD.
2137 (grida[cgj0] >= max_jcg &&
2138 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2144 for (j=0; (j<nrj); j++)
2146 jjcg = grida[cgj0+j];
2148 /* check whether this guy is in range! */
2149 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2152 r2=calc_dx2(XI,YI,ZI,cgcm[jjcg]);
2154 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2155 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2156 /* check energy group exclusions */
2157 if (!(i_egp_flags[jgid] & EGP_EXCL))
2161 if (nsr[jgid] >= MAX_CG)
2163 /* Add to short-range list */
2164 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2165 nsr[jgid],nl_sr[jgid],
2166 cgs->index,/* cgsatoms, */ bexcl,
2167 shift,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
2170 nl_sr[jgid][nsr[jgid]++]=jjcg;
2174 if (nlr_ljc[jgid] >= MAX_CG)
2176 /* Add to LJ+coulomb long-range list */
2177 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2178 nlr_ljc[jgid],nl_lr_ljc[jgid],top->cgs.index,
2179 bexcl,shift,fr,TRUE,TRUE,TRUE,fr->solvent_opt);
2182 nl_lr_ljc[jgid][nlr_ljc[jgid]++]=jjcg;
2186 if (nlr_one[jgid] >= MAX_CG)
2188 /* Add to long-range list with only coul, or only LJ */
2189 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2190 nlr_one[jgid],nl_lr_one[jgid],top->cgs.index,
2191 bexcl,shift,fr,TRUE,rvdw_lt_rcoul,rcoul_lt_rvdw,fr->solvent_opt);
2194 nl_lr_one[jgid][nlr_one[jgid]++]=jjcg;
2206 /* CHECK whether there is anything left in the buffers */
2207 for(nn=0; (nn<ngid); nn++)
2211 put_in_list(bHaveVdW,ngid,md,icg,nn,nsr[nn],nl_sr[nn],
2212 cgs->index, /* cgsatoms, */ bexcl,
2213 shift,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
2216 if (nlr_ljc[nn] > 0)
2218 put_in_list(bHaveVdW,ngid,md,icg,nn,nlr_ljc[nn],
2219 nl_lr_ljc[nn],top->cgs.index,
2220 bexcl,shift,fr,TRUE,TRUE,TRUE,fr->solvent_opt);
2223 if (nlr_one[nn] > 0)
2225 put_in_list(bHaveVdW,ngid,md,icg,nn,nlr_one[nn],
2226 nl_lr_one[nn],top->cgs.index,
2227 bexcl,shift,fr,TRUE,rvdw_lt_rcoul,rcoul_lt_rvdw,fr->solvent_opt);
2233 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2234 setexcl(cgs->index[icg],cgs->index[icg+1],&top->excls,FALSE,bexcl);
2236 /* No need to perform any left-over force calculations anymore (as we used to do here)
2237 * since we now save the proper long-range lists for later evaluation.
2242 /* Close neighbourlists */
2243 close_neighbor_lists(fr,bMakeQMMMnblist);
2248 void ns_realloc_natoms(gmx_ns_t *ns,int natoms)
2252 if (natoms > ns->nra_alloc)
2254 ns->nra_alloc = over_alloc_dd(natoms);
2255 srenew(ns->bexcl,ns->nra_alloc);
2256 for(i=0; i<ns->nra_alloc; i++)
2263 void init_ns(FILE *fplog,const t_commrec *cr,
2264 gmx_ns_t *ns,t_forcerec *fr,
2265 const gmx_mtop_t *mtop,
2268 int mt,icg,nr_in_cg,maxcg,i,j,jcg,ngid,ncg;
2272 /* Compute largest charge groups size (# atoms) */
2274 for(mt=0; mt<mtop->nmoltype; mt++) {
2275 cgs = &mtop->moltype[mt].cgs;
2276 for (icg=0; (icg < cgs->nr); icg++)
2278 nr_in_cg=max(nr_in_cg,(int)(cgs->index[icg+1]-cgs->index[icg]));
2282 /* Verify whether largest charge group is <= max cg.
2283 * This is determined by the type of the local exclusion type
2284 * Exclusions are stored in bits. (If the type is not large
2285 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2287 maxcg = sizeof(t_excl)*8;
2288 if (nr_in_cg > maxcg)
2290 gmx_fatal(FARGS,"Max #atoms in a charge group: %d > %d\n",
2294 ngid = mtop->groups.grps[egcENER].nr;
2295 snew(ns->bExcludeAlleg,ngid);
2296 for(i=0; i<ngid; i++) {
2297 ns->bExcludeAlleg[i] = TRUE;
2298 for(j=0; j<ngid; j++)
2300 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2302 ns->bExcludeAlleg[i] = FALSE;
2309 ns->grid = init_grid(fplog,fr);
2310 init_nsgrid_lists(fr,ngid,ns);
2315 snew(ns->ns_buf,ngid);
2316 for(i=0; (i<ngid); i++)
2318 snew(ns->ns_buf[i],SHIFTS);
2320 ncg = ncg_mtop(mtop);
2321 snew(ns->simple_aaj,2*ncg);
2322 for(jcg=0; (jcg<ncg); jcg++)
2324 ns->simple_aaj[jcg] = jcg;
2325 ns->simple_aaj[jcg+ncg] = jcg;
2329 /* Create array that determines whether or not atoms have VdW */
2330 snew(ns->bHaveVdW,fr->ntype);
2331 for(i=0; (i<fr->ntype); i++)
2333 for(j=0; (j<fr->ntype); j++)
2335 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2337 ((BHAMA(fr->nbfp,fr->ntype,i,j) != 0) ||
2338 (BHAMB(fr->nbfp,fr->ntype,i,j) != 0) ||
2339 (BHAMC(fr->nbfp,fr->ntype,i,j) != 0)) :
2340 ((C6(fr->nbfp,fr->ntype,i,j) != 0) ||
2341 (C12(fr->nbfp,fr->ntype,i,j) != 0))));
2345 pr_bvec(debug,0,"bHaveVdW",ns->bHaveVdW,fr->ntype,TRUE);
2349 if (!DOMAINDECOMP(cr))
2351 /* This could be reduced with particle decomposition */
2352 ns_realloc_natoms(ns,mtop->natoms);
2355 ns->nblist_initialized=FALSE;
2357 /* nbr list debug dump */
2359 char *ptr=getenv("GMX_DUMP_NL");
2362 ns->dump_nl=strtol(ptr,NULL,10);
2365 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2376 int search_neighbours(FILE *log,t_forcerec *fr,
2377 rvec x[],matrix box,
2378 gmx_localtop_t *top,
2379 gmx_groups_t *groups,
2381 t_nrnb *nrnb,t_mdatoms *md,
2382 real *lambda,real *dvdlambda,
2383 gmx_grppairener_t *grppener,
2385 gmx_bool bDoLongRangeNS,
2386 gmx_bool bPadListsForKernels)
2388 t_block *cgs=&(top->cgs);
2389 rvec box_size,grid_x0,grid_x1;
2391 real min_size,grid_dens;
2395 gmx_bool *i_egp_flags;
2396 int cg_start,cg_end,start,end;
2399 gmx_domdec_zones_t *dd_zones;
2400 put_in_list_t *put_in_list;
2404 /* Set some local variables */
2406 ngid = groups->grps[egcENER].nr;
2408 for(m=0; (m<DIM); m++)
2410 box_size[m] = box[m][m];
2413 if (fr->ePBC != epbcNONE)
2415 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC,box))
2417 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.");
2421 min_size = min(box_size[XX],min(box_size[YY],box_size[ZZ]));
2422 if (2*fr->rlistlong >= min_size)
2423 gmx_fatal(FARGS,"One of the box diagonal elements has become smaller than twice the cut-off length.");
2427 if (DOMAINDECOMP(cr))
2429 ns_realloc_natoms(ns,cgs->index[cgs->nr]);
2433 /* Reset the neighbourlists */
2434 reset_neighbor_lists(fr,TRUE,TRUE);
2436 if (bGrid && bFillGrid)
2440 if (DOMAINDECOMP(cr))
2442 dd_zones = domdec_zones(cr->dd);
2448 get_nsgrid_boundaries(grid->nboundeddim,box,NULL,NULL,NULL,NULL,
2449 cgs->nr,fr->cg_cm,grid_x0,grid_x1,&grid_dens);
2451 grid_first(log,grid,NULL,NULL,fr->ePBC,box,grid_x0,grid_x1,
2452 fr->rlistlong,grid_dens);
2456 /* Don't know why this all is... (DvdS 3/99) */
2462 end = (cgs->nr+1)/2;
2465 if (DOMAINDECOMP(cr))
2468 fill_grid(log,dd_zones,grid,end,-1,end,fr->cg_cm);
2470 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2474 fill_grid(log,NULL,grid,cgs->nr,fr->cg0,fr->hcg,fr->cg_cm);
2475 grid->icg0 = fr->cg0;
2476 grid->icg1 = fr->hcg;
2484 calc_elemnr(log,grid,start,end,cgs->nr);
2486 grid_last(log,grid,start,end,cgs->nr);
2490 check_grid(debug,grid);
2491 print_grid(debug,grid);
2496 /* Set the grid cell index for the test particle only.
2497 * The cell to cg index is not corrected, but that does not matter.
2499 fill_grid(log,NULL,ns->grid,fr->hcg,fr->hcg-1,fr->hcg,fr->cg_cm);
2503 if (!fr->ns.bCGlist)
2505 put_in_list = put_in_list_at;
2509 put_in_list = put_in_list_cg;
2516 nsearch = nsgrid_core(log,cr,fr,box,box_size,ngid,top,
2517 grid,x,ns->bexcl,ns->bExcludeAlleg,
2518 nrnb,md,lambda,dvdlambda,grppener,
2519 put_in_list,ns->bHaveVdW,
2520 bDoLongRangeNS,FALSE);
2522 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2523 * the classical calculation. The charge-charge interaction
2524 * between QM and MM atoms is handled in the QMMM core calculation
2525 * (see QMMM.c). The VDW however, we'd like to compute classically
2526 * and the QM MM atom pairs have just been put in the
2527 * corresponding neighbourlists. in case of QMMM we still need to
2528 * fill a special QMMM neighbourlist that contains all neighbours
2529 * of the QM atoms. If bQMMM is true, this list will now be made:
2531 if (fr->bQMMM && fr->qr->QMMMscheme!=eQMMMschemeoniom)
2533 nsearch += nsgrid_core(log,cr,fr,box,box_size,ngid,top,
2534 grid,x,ns->bexcl,ns->bExcludeAlleg,
2535 nrnb,md,lambda,dvdlambda,grppener,
2536 put_in_list_qmmm,ns->bHaveVdW,
2537 bDoLongRangeNS,TRUE);
2542 nsearch = ns_simple_core(fr,top,md,box,box_size,
2543 ns->bexcl,ns->simple_aaj,
2544 ngid,ns->ns_buf,put_in_list,ns->bHaveVdW);
2552 inc_nrnb(nrnb,eNR_NS,nsearch);
2553 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2558 int natoms_beyond_ns_buffer(t_inputrec *ir,t_forcerec *fr,t_block *cgs,
2559 matrix scale_tot,rvec *x)
2561 int cg0,cg1,cg,a0,a1,a,i,j;
2562 real rint,hbuf2,scale;
2564 gmx_bool bIsotropic;
2569 rint = max(ir->rcoulomb,ir->rvdw);
2570 if (ir->rlist < rint)
2572 gmx_fatal(FARGS,"The neighbor search buffer has negative size: %f nm",
2580 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2582 hbuf2 = sqr(0.5*(ir->rlist - rint));
2583 for(cg=cg0; cg<cg1; cg++)
2585 a0 = cgs->index[cg];
2586 a1 = cgs->index[cg+1];
2587 for(a=a0; a<a1; a++)
2589 if (distance2(cg_cm[cg],x[a]) > hbuf2)
2599 scale = scale_tot[0][0];
2600 for(i=1; i<DIM; i++)
2602 /* With anisotropic scaling, the original spherical ns volumes become
2603 * ellipsoids. To avoid costly transformations we use the minimum
2604 * eigenvalue of the scaling matrix for determining the buffer size.
2605 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2607 scale = min(scale,scale_tot[i][i]);
2608 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2614 if (scale_tot[i][j] != 0)
2620 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2623 for(cg=cg0; cg<cg1; cg++)
2625 svmul(scale,cg_cm[cg],cgsc);
2626 a0 = cgs->index[cg];
2627 a1 = cgs->index[cg+1];
2628 for(a=a0; a<a1; a++)
2630 if (distance2(cgsc,x[a]) > hbuf2)
2639 /* Anistropic scaling */
2640 for(cg=cg0; cg<cg1; cg++)
2642 /* Since scale_tot contains the transpose of the scaling matrix,
2643 * we need to multiply with the transpose.
2645 tmvmul_ur0(scale_tot,cg_cm[cg],cgsc);
2646 a0 = cgs->index[cg];
2647 a1 = cgs->index[cg+1];
2648 for(a=a0; a<a1; a++)
2650 if (distance2(cgsc,x[a]) > hbuf2)