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, type=%d), maxnri=%d\n",
103 nl->ielec,nl->ivdw,nl->igeometry,nl->type,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 int igeometry, int type)
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;
145 nl->igeometry = igeometry;
147 if (nl->type==GMX_NBLIST_INTERACTION_FREE_ENERGY)
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, type=%d) for %s interactions,\nwith %d SR, %d LR atoms.\n",
176 nl->ielec,nl->ivdw,nl->type,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,type;
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;
223 type = GMX_NBLIST_INTERACTION_STANDARD;
225 fr->ns.bCGlist = (getenv("GMX_NBLISTCG") != 0);
228 igeometry_def = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
232 igeometry_def = GMX_NBLIST_GEOMETRY_CG_CG;
235 fprintf(log,"\nUsing charge-group - charge-group neighbor lists and kernels\n\n");
239 if (fr->solvent_opt == esolTIP4P) {
240 igeometry_w = GMX_NBLIST_GEOMETRY_WATER4_PARTICLE;
241 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER4_WATER4;
243 igeometry_w = GMX_NBLIST_GEOMETRY_WATER3_PARTICLE;
244 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER3_WATER3;
247 for(i=0; i<fr->nnblists; i++)
249 nbl = &(fr->nblists[i]);
251 if ((fr->adress_type!=eAdressOff) && (i>=fr->nnblists/2)){
252 type=GMX_NBLIST_INTERACTION_ADRESS;
254 init_nblist(log,&nbl->nlist_sr[eNL_VDWQQ],&nbl->nlist_lr[eNL_VDWQQ],
255 maxsr,maxlr,ivdw,ivdwmod,ielec,ielecmod,igeometry_def, type);
256 init_nblist(log,&nbl->nlist_sr[eNL_VDW],&nbl->nlist_lr[eNL_VDW],
257 maxsr,maxlr,ivdw,ivdwmod,GMX_NBKERNEL_ELEC_NONE,eintmodNONE,igeometry_def, type);
258 init_nblist(log,&nbl->nlist_sr[eNL_QQ],&nbl->nlist_lr[eNL_QQ],
259 maxsr,maxlr,GMX_NBKERNEL_VDW_NONE,eintmodNONE,ielec,ielecmod,igeometry_def, type);
260 init_nblist(log,&nbl->nlist_sr[eNL_VDWQQ_WATER],&nbl->nlist_lr[eNL_VDWQQ_WATER],
261 maxsr_wat,maxlr_wat,ivdw,ivdwmod,ielec,ielecmod, igeometry_w, type);
262 init_nblist(log,&nbl->nlist_sr[eNL_QQ_WATER],&nbl->nlist_lr[eNL_QQ_WATER],
263 maxsr_wat,maxlr_wat,GMX_NBKERNEL_VDW_NONE,eintmodNONE,ielec,ielecmod, igeometry_w, type);
264 init_nblist(log,&nbl->nlist_sr[eNL_VDWQQ_WATERWATER],&nbl->nlist_lr[eNL_VDWQQ_WATERWATER],
265 maxsr_wat,maxlr_wat,ivdw,ivdwmod,ielec,ielecmod, igeometry_ww, type);
266 init_nblist(log,&nbl->nlist_sr[eNL_QQ_WATERWATER],&nbl->nlist_lr[eNL_QQ_WATERWATER],
267 maxsr_wat,maxlr_wat,GMX_NBKERNEL_VDW_NONE,eintmodNONE,ielec,ielecmod, igeometry_ww, type);
269 /* Did we get the solvent loops so we can use optimized water kernels? */
270 if(nbl->nlist_sr[eNL_VDWQQ_WATER].kernelptr_vf==NULL
271 || nbl->nlist_sr[eNL_QQ_WATER].kernelptr_vf==NULL
272 #ifndef DISABLE_WATERWATER_NLIST
273 || nbl->nlist_sr[eNL_VDWQQ_WATERWATER].kernelptr_vf==NULL
274 || nbl->nlist_sr[eNL_QQ_WATERWATER].kernelptr_vf==NULL
278 fr->solvent_opt = esolNO;
279 fprintf(log,"Note: The available nonbonded kernels do not support water optimization - disabling.\n");
282 if (fr->efep != efepNO)
284 if ((fr->bEwald) && (fr->sc_alphacoul > 0)) /* need to handle long range differently if using softcore */
286 ielecf = GMX_NBKERNEL_ELEC_EWALD;
287 ielecmodf = eintmodNONE;
292 ielecmodf = ielecmod;
295 init_nblist(log,&nbl->nlist_sr[eNL_VDWQQ_FREE],&nbl->nlist_lr[eNL_VDWQQ_FREE],
296 maxsr,maxlr,ivdw,ivdwmod,ielecf,ielecmod,GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
297 init_nblist(log,&nbl->nlist_sr[eNL_VDW_FREE],&nbl->nlist_lr[eNL_VDW_FREE],
298 maxsr,maxlr,ivdw,ivdwmod,GMX_NBKERNEL_ELEC_NONE,eintmodNONE,GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
299 init_nblist(log,&nbl->nlist_sr[eNL_QQ_FREE],&nbl->nlist_lr[eNL_QQ_FREE],
300 maxsr,maxlr,GMX_NBKERNEL_VDW_NONE,eintmodNONE,ielecf,ielecmod,GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
304 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
306 init_nblist(log,&fr->QMMMlist,NULL,
307 maxsr,maxlr,0,0,ielec,ielecmod,GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_STANDARD);
315 fr->ns.nblist_initialized=TRUE;
318 static void reset_nblist(t_nblist *nl)
329 static void reset_neighbor_lists(t_forcerec *fr,gmx_bool bResetSR, gmx_bool bResetLR)
335 /* only reset the short-range nblist */
336 reset_nblist(&(fr->QMMMlist));
339 for(n=0; n<fr->nnblists; n++)
341 for(i=0; i<eNL_NR; i++)
345 reset_nblist( &(fr->nblists[n].nlist_sr[i]) );
349 reset_nblist( &(fr->nblists[n].nlist_lr[i]) );
358 static inline void new_i_nblist(t_nblist *nlist,
359 gmx_bool bLR,atom_id i_atom,int shift,int gid)
365 /* Check whether we have to increase the i counter */
367 (nlist->iinr[nri] != i_atom) ||
368 (nlist->shift[nri] != shift) ||
369 (nlist->gid[nri] != gid))
371 /* This is something else. Now see if any entries have
372 * been added in the list of the previous atom.
375 ((nlist->jindex[nri+1] > nlist->jindex[nri]) &&
376 (nlist->gid[nri] != -1)))
378 /* If so increase the counter */
381 if (nlist->nri >= nlist->maxnri)
383 nlist->maxnri += over_alloc_large(nlist->nri);
384 reallocate_nblist(nlist);
387 /* Set the number of neighbours and the atom number */
388 nlist->jindex[nri+1] = nlist->jindex[nri];
389 nlist->iinr[nri] = i_atom;
390 nlist->gid[nri] = gid;
391 nlist->shift[nri] = shift;
395 static inline void close_i_nblist(t_nblist *nlist)
397 int nri = nlist->nri;
402 /* Add elements up to padding. Since we allocate memory in units
403 * of the simd_padding width, we do not have to check for possible
404 * list reallocation here.
406 while((nlist->nrj % nlist->simd_padding_width)!=0)
408 /* Use -4 here, so we can write forces for 4 atoms before real data */
409 nlist->jjnr[nlist->nrj++]=-4;
411 nlist->jindex[nri+1] = nlist->nrj;
413 len=nlist->nrj - nlist->jindex[nri];
415 /* nlist length for water i molecules is treated statically
418 if (len > nlist->maxlen)
425 static inline void close_nblist(t_nblist *nlist)
427 /* Only close this nblist when it has been initialized.
428 * Avoid the creation of i-lists with no j-particles.
432 /* Some assembly kernels do not support empty lists,
433 * make sure here that we don't generate any empty lists.
434 * With the current ns code this branch is taken in two cases:
435 * No i-particles at all: nri=-1 here
436 * There are i-particles, but no j-particles; nri=0 here
442 /* Close list number nri by incrementing the count */
447 static inline void close_neighbor_lists(t_forcerec *fr,gmx_bool bMakeQMMMnblist)
453 close_nblist(&(fr->QMMMlist));
456 for(n=0; n<fr->nnblists; n++)
458 for(i=0; (i<eNL_NR); i++)
460 close_nblist(&(fr->nblists[n].nlist_sr[i]));
461 close_nblist(&(fr->nblists[n].nlist_lr[i]));
467 static inline void add_j_to_nblist(t_nblist *nlist,atom_id j_atom,gmx_bool bLR)
471 if (nlist->nrj >= nlist->maxnrj)
473 nlist->maxnrj = round_up_to_simd_width(over_alloc_small(nlist->nrj + 1),nlist->simd_padding_width);
476 fprintf(debug,"Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
477 bLR ? "LR" : "SR",nlist->ielec,nlist->ivdw,nlist->type,nlist->igeometry,nlist->maxnrj);
479 srenew(nlist->jjnr,nlist->maxnrj);
482 nlist->jjnr[nrj] = j_atom;
486 static inline void add_j_to_nblist_cg(t_nblist *nlist,
487 atom_id j_start,int j_end,
488 t_excl *bexcl,gmx_bool i_is_j,
494 if (nlist->nrj >= nlist->maxnrj)
496 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
498 fprintf(debug,"Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
499 bLR ? "LR" : "SR",nlist->ielec,nlist->ivdw,nlist->type,nlist->igeometry,nlist->maxnrj);
501 srenew(nlist->jjnr ,nlist->maxnrj);
502 srenew(nlist->jjnr_end,nlist->maxnrj);
503 srenew(nlist->excl ,nlist->maxnrj*MAX_CGCGSIZE);
506 nlist->jjnr[nrj] = j_start;
507 nlist->jjnr_end[nrj] = j_end;
509 if (j_end - j_start > MAX_CGCGSIZE)
511 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);
514 /* Set the exclusions */
515 for(j=j_start; j<j_end; j++)
517 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
521 /* Avoid double counting of intra-cg interactions */
522 for(j=1; j<j_end-j_start; j++)
524 nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
532 put_in_list_t(gmx_bool bHaveVdW[],
549 put_in_list_at(gmx_bool bHaveVdW[],
565 /* The a[] index has been removed,
566 * to put it back in i_atom should be a[i0] and jj should be a[jj].
571 t_nblist * vdwc_free = NULL;
572 t_nblist * vdw_free = NULL;
573 t_nblist * coul_free = NULL;
574 t_nblist * vdwc_ww = NULL;
575 t_nblist * coul_ww = NULL;
577 int i,j,jcg,igid,gid,nbl_ind,ind_ij;
578 atom_id jj,jj0,jj1,i_atom;
583 real *charge,*chargeB;
585 gmx_bool bFreeEnergy,bFree,bFreeJ,bNotEx,*bPert;
586 gmx_bool bDoVdW_i,bDoCoul_i,bDoCoul_i_sol;
590 /* Copy some pointers */
592 charge = md->chargeA;
593 chargeB = md->chargeB;
596 bPert = md->bPerturbed;
600 nicg = index[icg+1]-i0;
602 /* Get the i charge group info */
603 igid = GET_CGINFO_GID(cginfo[icg]);
605 iwater = (solvent_opt!=esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
610 /* Check if any of the particles involved are perturbed.
611 * If not we can do the cheaper normal put_in_list
612 * and use more solvent optimization.
614 for(i=0; i<nicg; i++)
616 bFreeEnergy |= bPert[i0+i];
618 /* Loop over the j charge groups */
619 for(j=0; (j<nj && !bFreeEnergy); j++)
624 /* Finally loop over the atoms in the j-charge group */
625 for(jj=jj0; jj<jj1; jj++)
627 bFreeEnergy |= bPert[jj];
632 /* Unpack pointers to neighbourlist structs */
633 if (fr->nnblists == 1)
639 nbl_ind = fr->gid2nblists[GID(igid,jgid,ngid)];
643 nlist = fr->nblists[nbl_ind].nlist_lr;
647 nlist = fr->nblists[nbl_ind].nlist_sr;
650 if (iwater != esolNO)
652 vdwc = &nlist[eNL_VDWQQ_WATER];
653 vdw = &nlist[eNL_VDW];
654 coul = &nlist[eNL_QQ_WATER];
655 #ifndef DISABLE_WATERWATER_NLIST
656 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
657 coul_ww = &nlist[eNL_QQ_WATERWATER];
662 vdwc = &nlist[eNL_VDWQQ];
663 vdw = &nlist[eNL_VDW];
664 coul = &nlist[eNL_QQ];
669 if (iwater != esolNO)
671 /* Loop over the atoms in the i charge group */
673 gid = GID(igid,jgid,ngid);
674 /* Create new i_atom for each energy group */
675 if (bDoCoul && bDoVdW)
677 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
678 #ifndef DISABLE_WATERWATER_NLIST
679 new_i_nblist(vdwc_ww,bLR,i_atom,shift,gid);
684 new_i_nblist(vdw,bLR,i_atom,shift,gid);
688 new_i_nblist(coul,bLR,i_atom,shift,gid);
689 #ifndef DISABLE_WATERWATER_NLIST
690 new_i_nblist(coul_ww,bLR,i_atom,shift,gid);
693 /* Loop over the j charge groups */
694 for(j=0; (j<nj); j++)
704 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
706 if (iwater == esolSPC && jwater == esolSPC)
708 /* Interaction between two SPC molecules */
711 /* VdW only - only first atoms in each water interact */
712 add_j_to_nblist(vdw,jj0,bLR);
716 #ifdef DISABLE_WATERWATER_NLIST
717 /* Add entries for the three atoms - only do VdW if we need to */
720 add_j_to_nblist(coul,jj0,bLR);
724 add_j_to_nblist(vdwc,jj0,bLR);
726 add_j_to_nblist(coul,jj0+1,bLR);
727 add_j_to_nblist(coul,jj0+2,bLR);
729 /* One entry for the entire water-water interaction */
732 add_j_to_nblist(coul_ww,jj0,bLR);
736 add_j_to_nblist(vdwc_ww,jj0,bLR);
741 else if (iwater == esolTIP4P && jwater == esolTIP4P)
743 /* Interaction between two TIP4p molecules */
746 /* VdW only - only first atoms in each water interact */
747 add_j_to_nblist(vdw,jj0,bLR);
751 #ifdef DISABLE_WATERWATER_NLIST
752 /* Add entries for the four atoms - only do VdW if we need to */
755 add_j_to_nblist(vdw,jj0,bLR);
757 add_j_to_nblist(coul,jj0+1,bLR);
758 add_j_to_nblist(coul,jj0+2,bLR);
759 add_j_to_nblist(coul,jj0+3,bLR);
761 /* One entry for the entire water-water interaction */
764 add_j_to_nblist(coul_ww,jj0,bLR);
768 add_j_to_nblist(vdwc_ww,jj0,bLR);
775 /* j charge group is not water, but i is.
776 * Add entries to the water-other_atom lists; the geometry of the water
777 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
778 * so we don't care if it is SPC or TIP4P...
785 for(jj=jj0; (jj<jj1); jj++)
789 add_j_to_nblist(coul,jj,bLR);
795 for(jj=jj0; (jj<jj1); jj++)
797 if (bHaveVdW[type[jj]])
799 add_j_to_nblist(vdw,jj,bLR);
805 /* _charge_ _groups_ interact with both coulomb and LJ */
806 /* Check which atoms we should add to the lists! */
807 for(jj=jj0; (jj<jj1); jj++)
809 if (bHaveVdW[type[jj]])
813 add_j_to_nblist(vdwc,jj,bLR);
817 add_j_to_nblist(vdw,jj,bLR);
820 else if (charge[jj] != 0)
822 add_j_to_nblist(coul,jj,bLR);
829 close_i_nblist(coul);
830 close_i_nblist(vdwc);
831 #ifndef DISABLE_WATERWATER_NLIST
832 close_i_nblist(coul_ww);
833 close_i_nblist(vdwc_ww);
838 /* no solvent as i charge group */
839 /* Loop over the atoms in the i charge group */
840 for(i=0; i<nicg; i++)
843 gid = GID(igid,jgid,ngid);
846 /* Create new i_atom for each energy group */
847 if (bDoVdW && bDoCoul)
849 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
853 new_i_nblist(vdw,bLR,i_atom,shift,gid);
857 new_i_nblist(coul,bLR,i_atom,shift,gid);
859 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
860 bDoCoul_i = (bDoCoul && qi!=0);
862 if (bDoVdW_i || bDoCoul_i)
864 /* Loop over the j charge groups */
865 for(j=0; (j<nj); j++)
869 /* Check for large charge groups */
880 /* Finally loop over the atoms in the j-charge group */
881 for(jj=jj0; jj<jj1; jj++)
883 bNotEx = NOTEXCL(bExcl,i,jj);
891 add_j_to_nblist(coul,jj,bLR);
896 if (bHaveVdW[type[jj]])
898 add_j_to_nblist(vdw,jj,bLR);
903 if (bHaveVdW[type[jj]])
907 add_j_to_nblist(vdwc,jj,bLR);
911 add_j_to_nblist(vdw,jj,bLR);
914 else if (charge[jj] != 0)
916 add_j_to_nblist(coul,jj,bLR);
924 close_i_nblist(coul);
925 close_i_nblist(vdwc);
931 /* we are doing free energy */
932 vdwc_free = &nlist[eNL_VDWQQ_FREE];
933 vdw_free = &nlist[eNL_VDW_FREE];
934 coul_free = &nlist[eNL_QQ_FREE];
935 /* Loop over the atoms in the i charge group */
936 for(i=0; i<nicg; i++)
939 gid = GID(igid,jgid,ngid);
941 qiB = chargeB[i_atom];
943 /* Create new i_atom for each energy group */
944 if (bDoVdW && bDoCoul)
945 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
947 new_i_nblist(vdw,bLR,i_atom,shift,gid);
949 new_i_nblist(coul,bLR,i_atom,shift,gid);
951 new_i_nblist(vdw_free,bLR,i_atom,shift,gid);
952 new_i_nblist(coul_free,bLR,i_atom,shift,gid);
953 new_i_nblist(vdwc_free,bLR,i_atom,shift,gid);
955 bDoVdW_i = (bDoVdW &&
956 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
957 bDoCoul_i = (bDoCoul && (qi!=0 || qiB!=0));
958 /* For TIP4P the first atom does not have a charge,
959 * but the last three do. So we should still put an atom
960 * without LJ but with charge in the water-atom neighborlist
961 * for a TIP4p i charge group.
962 * For SPC type water the first atom has LJ and charge,
963 * so there is no such problem.
965 if (iwater == esolNO)
967 bDoCoul_i_sol = bDoCoul_i;
971 bDoCoul_i_sol = bDoCoul;
974 if (bDoVdW_i || bDoCoul_i_sol)
976 /* Loop over the j charge groups */
977 for(j=0; (j<nj); j++)
981 /* Check for large charge groups */
992 /* Finally loop over the atoms in the j-charge group */
993 bFree = bPert[i_atom];
994 for(jj=jj0; (jj<jj1); jj++)
996 bFreeJ = bFree || bPert[jj];
997 /* Complicated if, because the water H's should also
998 * see perturbed j-particles
1000 if (iwater==esolNO || i==0 || bFreeJ)
1002 bNotEx = NOTEXCL(bExcl,i,jj);
1010 if (charge[jj]!=0 || chargeB[jj]!=0)
1012 add_j_to_nblist(coul_free,jj,bLR);
1015 else if (!bDoCoul_i)
1017 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1019 add_j_to_nblist(vdw_free,jj,bLR);
1024 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1026 if (charge[jj]!=0 || chargeB[jj]!=0)
1028 add_j_to_nblist(vdwc_free,jj,bLR);
1032 add_j_to_nblist(vdw_free,jj,bLR);
1035 else if (charge[jj]!=0 || chargeB[jj]!=0)
1036 add_j_to_nblist(coul_free,jj,bLR);
1041 /* This is done whether or not bWater is set */
1042 if (charge[jj] != 0)
1044 add_j_to_nblist(coul,jj,bLR);
1047 else if (!bDoCoul_i_sol)
1049 if (bHaveVdW[type[jj]])
1051 add_j_to_nblist(vdw,jj,bLR);
1056 if (bHaveVdW[type[jj]])
1058 if (charge[jj] != 0)
1060 add_j_to_nblist(vdwc,jj,bLR);
1064 add_j_to_nblist(vdw,jj,bLR);
1067 else if (charge[jj] != 0)
1069 add_j_to_nblist(coul,jj,bLR);
1077 close_i_nblist(vdw);
1078 close_i_nblist(coul);
1079 close_i_nblist(vdwc);
1080 close_i_nblist(vdw_free);
1081 close_i_nblist(coul_free);
1082 close_i_nblist(vdwc_free);
1088 put_in_list_adress(gmx_bool bHaveVdW[],
1104 /* The a[] index has been removed,
1105 * to put it back in i_atom should be a[i0] and jj should be a[jj].
1110 t_nblist * vdwc_adress = NULL;
1111 t_nblist * vdw_adress = NULL;
1112 t_nblist * coul_adress = NULL;
1113 t_nblist * vdwc_ww = NULL;
1114 t_nblist * coul_ww = NULL;
1116 int i,j,jcg,igid,gid,nbl_ind,nbl_ind_adress;
1117 atom_id jj,jj0,jj1,i_atom;
1122 real *charge,*chargeB;
1125 gmx_bool bFreeEnergy,bFree,bFreeJ,bNotEx,*bPert;
1126 gmx_bool bDoVdW_i,bDoCoul_i,bDoCoul_i_sol;
1128 gmx_bool j_all_atom;
1130 t_nblist *nlist, *nlist_adress;
1132 /* Copy some pointers */
1133 cginfo = fr->cginfo;
1134 charge = md->chargeA;
1135 chargeB = md->chargeB;
1138 bPert = md->bPerturbed;
1141 /* Get atom range */
1143 nicg = index[icg+1]-i0;
1145 /* Get the i charge group info */
1146 igid = GET_CGINFO_GID(cginfo[icg]);
1148 iwater = (solvent_opt!=esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
1152 gmx_fatal(FARGS,"AdResS does not support free energy pertubation\n");
1155 /* Unpack pointers to neighbourlist structs */
1156 if (fr->nnblists == 2)
1163 nbl_ind = fr->gid2nblists[GID(igid,jgid,ngid)];
1164 nbl_ind_adress = nbl_ind+fr->nnblists/2;
1168 nlist = fr->nblists[nbl_ind].nlist_lr;
1169 nlist_adress= fr->nblists[nbl_ind_adress].nlist_lr;
1173 nlist = fr->nblists[nbl_ind].nlist_sr;
1174 nlist_adress = fr->nblists[nbl_ind_adress].nlist_sr;
1178 vdwc = &nlist[eNL_VDWQQ];
1179 vdw = &nlist[eNL_VDW];
1180 coul = &nlist[eNL_QQ];
1182 vdwc_adress = &nlist_adress[eNL_VDWQQ];
1183 vdw_adress = &nlist_adress[eNL_VDW];
1184 coul_adress = &nlist_adress[eNL_QQ];
1186 /* We do not support solvent optimization with AdResS for now.
1187 For this we would need hybrid solvent-other kernels */
1189 /* no solvent as i charge group */
1190 /* Loop over the atoms in the i charge group */
1191 for(i=0; i<nicg; i++)
1194 gid = GID(igid,jgid,ngid);
1195 qi = charge[i_atom];
1197 /* Create new i_atom for each energy group */
1198 if (bDoVdW && bDoCoul)
1200 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
1201 new_i_nblist(vdwc_adress,bLR,i_atom,shift,gid);
1206 new_i_nblist(vdw,bLR,i_atom,shift,gid);
1207 new_i_nblist(vdw_adress,bLR,i_atom,shift,gid);
1212 new_i_nblist(coul,bLR,i_atom,shift,gid);
1213 new_i_nblist(coul_adress,bLR,i_atom,shift,gid);
1215 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
1216 bDoCoul_i = (bDoCoul && qi!=0);
1218 if (bDoVdW_i || bDoCoul_i)
1220 /* Loop over the j charge groups */
1221 for(j=0; (j<nj); j++)
1225 /* Check for large charge groups */
1236 /* Finally loop over the atoms in the j-charge group */
1237 for(jj=jj0; jj<jj1; jj++)
1239 bNotEx = NOTEXCL(bExcl,i,jj);
1241 b_hybrid=!((wf[i_atom]==1&&wf[jj]==1)||(wf[i_atom] ==0 && wf[jj]==0));
1247 if (charge[jj] != 0)
1250 add_j_to_nblist(coul,jj,bLR);
1252 add_j_to_nblist(coul_adress,jj,bLR);
1256 else if (!bDoCoul_i)
1258 if (bHaveVdW[type[jj]])
1261 add_j_to_nblist(vdw,jj,bLR);
1263 add_j_to_nblist(vdw_adress,jj,bLR);
1269 if (bHaveVdW[type[jj]])
1271 if (charge[jj] != 0)
1274 add_j_to_nblist(vdwc,jj,bLR);
1276 add_j_to_nblist(vdwc_adress,jj,bLR);
1282 add_j_to_nblist(vdw,jj,bLR);
1284 add_j_to_nblist(vdw_adress,jj,bLR);
1289 else if (charge[jj] != 0)
1292 add_j_to_nblist(coul,jj,bLR);
1294 add_j_to_nblist(coul_adress,jj,bLR);
1303 close_i_nblist(vdw);
1304 close_i_nblist(coul);
1305 close_i_nblist(vdwc);
1306 close_i_nblist(vdw_adress);
1307 close_i_nblist(coul_adress);
1308 close_i_nblist(vdwc_adress);
1314 put_in_list_qmmm(gmx_bool bHaveVdW[],
1331 int i,j,jcg,igid,gid;
1332 atom_id jj,jj0,jj1,i_atom;
1336 /* Get atom range */
1338 nicg = index[icg+1]-i0;
1340 /* Get the i charge group info */
1341 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1343 coul = &fr->QMMMlist;
1345 /* Loop over atoms in the ith charge group */
1346 for (i=0;i<nicg;i++)
1349 gid = GID(igid,jgid,ngid);
1350 /* Create new i_atom for each energy group */
1351 new_i_nblist(coul,bLR,i_atom,shift,gid);
1353 /* Loop over the j charge groups */
1358 /* Charge groups cannot have QM and MM atoms simultaneously */
1363 /* Finally loop over the atoms in the j-charge group */
1364 for(jj=jj0; jj<jj1; jj++)
1366 bNotEx = NOTEXCL(bExcl,i,jj);
1368 add_j_to_nblist(coul,jj,bLR);
1372 close_i_nblist(coul);
1377 put_in_list_cg(gmx_bool bHaveVdW[],
1394 int igid,gid,nbl_ind;
1398 cginfo = fr->cginfo[icg];
1400 igid = GET_CGINFO_GID(cginfo);
1401 gid = GID(igid,jgid,ngid);
1403 /* Unpack pointers to neighbourlist structs */
1404 if (fr->nnblists == 1)
1410 nbl_ind = fr->gid2nblists[gid];
1414 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1418 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1421 /* Make a new neighbor list for charge group icg.
1422 * Currently simply one neighbor list is made with LJ and Coulomb.
1423 * If required, zero interactions could be removed here
1424 * or in the force loop.
1426 new_i_nblist(vdwc,bLR,index[icg],shift,gid);
1427 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1429 for(j=0; (j<nj); j++)
1432 /* Skip the icg-icg pairs if all self interactions are excluded */
1433 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1435 /* Here we add the j charge group jcg to the list,
1436 * exclusions are also added to the list.
1438 add_j_to_nblist_cg(vdwc,index[jcg],index[jcg+1],bExcl,icg==jcg,bLR);
1442 close_i_nblist(vdwc);
1445 static void setexcl(atom_id start,atom_id end,t_blocka *excl,gmx_bool b,
1452 for(i=start; i<end; i++)
1454 for(k=excl->index[i]; k<excl->index[i+1]; k++)
1456 SETEXCL(bexcl,i-start,excl->a[k]);
1462 for(i=start; i<end; i++)
1464 for(k=excl->index[i]; k<excl->index[i+1]; k++)
1466 RMEXCL(bexcl,i-start,excl->a[k]);
1472 int calc_naaj(int icg,int cgtot)
1476 if ((cgtot % 2) == 1)
1478 /* Odd number of charge groups, easy */
1479 naaj = 1 + (cgtot/2);
1481 else if ((cgtot % 4) == 0)
1483 /* Multiple of four is hard */
1520 fprintf(log,"naaj=%d\n",naaj);
1526 /************************************************
1528 * S I M P L E C O R E S T U F F
1530 ************************************************/
1532 static real calc_image_tric(rvec xi,rvec xj,matrix box,
1533 rvec b_inv,int *shift)
1535 /* This code assumes that the cut-off is smaller than
1536 * a half times the smallest diagonal element of the box.
1543 /* Compute diff vector */
1544 dz = xj[ZZ] - xi[ZZ];
1545 dy = xj[YY] - xi[YY];
1546 dx = xj[XX] - xi[XX];
1548 /* Perform NINT operation, using trunc operation, therefore
1549 * we first add 2.5 then subtract 2 again
1551 tz = dz*b_inv[ZZ] + h25;
1553 dz -= tz*box[ZZ][ZZ];
1554 dy -= tz*box[ZZ][YY];
1555 dx -= tz*box[ZZ][XX];
1557 ty = dy*b_inv[YY] + h25;
1559 dy -= ty*box[YY][YY];
1560 dx -= ty*box[YY][XX];
1562 tx = dx*b_inv[XX]+h25;
1564 dx -= tx*box[XX][XX];
1566 /* Distance squared */
1567 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1569 *shift = XYZ2IS(tx,ty,tz);
1574 static real calc_image_rect(rvec xi,rvec xj,rvec box_size,
1575 rvec b_inv,int *shift)
1583 /* Compute diff vector */
1584 dx = xj[XX] - xi[XX];
1585 dy = xj[YY] - xi[YY];
1586 dz = xj[ZZ] - xi[ZZ];
1588 /* Perform NINT operation, using trunc operation, therefore
1589 * we first add 1.5 then subtract 1 again
1591 tx = dx*b_inv[XX] + h15;
1592 ty = dy*b_inv[YY] + h15;
1593 tz = dz*b_inv[ZZ] + h15;
1598 /* Correct diff vector for translation */
1599 ddx = tx*box_size[XX] - dx;
1600 ddy = ty*box_size[YY] - dy;
1601 ddz = tz*box_size[ZZ] - dz;
1603 /* Distance squared */
1604 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1606 *shift = XYZ2IS(tx,ty,tz);
1611 static void add_simple(t_ns_buf *nsbuf,int nrj,atom_id cg_j,
1612 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1613 int icg,int jgid,t_block *cgs,t_excl bexcl[],
1614 int shift,t_forcerec *fr,put_in_list_t *put_in_list)
1616 if (nsbuf->nj + nrj > MAX_CG)
1618 put_in_list(bHaveVdW,ngid,md,icg,jgid,nsbuf->ncg,nsbuf->jcg,
1619 cgs->index,bexcl,shift,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
1620 /* Reset buffer contents */
1621 nsbuf->ncg = nsbuf->nj = 0;
1623 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1627 static void ns_inner_tric(rvec x[],int icg,int *i_egp_flags,
1628 int njcg,atom_id jcg[],
1629 matrix box,rvec b_inv,real rcut2,
1630 t_block *cgs,t_ns_buf **ns_buf,
1631 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1632 t_excl bexcl[],t_forcerec *fr,
1633 put_in_list_t *put_in_list)
1637 int *cginfo=fr->cginfo;
1638 atom_id cg_j,*cgindex;
1641 cgindex = cgs->index;
1643 for(j=0; (j<njcg); j++)
1646 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1647 if (calc_image_tric(x[icg],x[cg_j],box,b_inv,&shift) < rcut2)
1649 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1650 if (!(i_egp_flags[jgid] & EGP_EXCL))
1652 add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1653 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1660 static void ns_inner_rect(rvec x[],int icg,int *i_egp_flags,
1661 int njcg,atom_id jcg[],
1662 gmx_bool bBox,rvec box_size,rvec b_inv,real rcut2,
1663 t_block *cgs,t_ns_buf **ns_buf,
1664 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1665 t_excl bexcl[],t_forcerec *fr,
1666 put_in_list_t *put_in_list)
1670 int *cginfo=fr->cginfo;
1671 atom_id cg_j,*cgindex;
1674 cgindex = cgs->index;
1678 for(j=0; (j<njcg); j++)
1681 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1682 if (calc_image_rect(x[icg],x[cg_j],box_size,b_inv,&shift) < rcut2)
1684 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1685 if (!(i_egp_flags[jgid] & EGP_EXCL))
1687 add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1688 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1696 for(j=0; (j<njcg); j++)
1699 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1700 if ((rcut2 == 0) || (distance2(x[icg],x[cg_j]) < rcut2)) {
1701 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1702 if (!(i_egp_flags[jgid] & EGP_EXCL))
1704 add_simple(&ns_buf[jgid][CENTRAL],nrj,cg_j,
1705 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,CENTRAL,fr,
1713 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1715 static int ns_simple_core(t_forcerec *fr,
1716 gmx_localtop_t *top,
1718 matrix box,rvec box_size,
1719 t_excl bexcl[],atom_id *aaj,
1720 int ngid,t_ns_buf **ns_buf,
1721 put_in_list_t *put_in_list,gmx_bool bHaveVdW[])
1725 int nsearch,icg,jcg,igid,i0,nri,nn;
1728 /* atom_id *i_atoms; */
1729 t_block *cgs=&(top->cgs);
1730 t_blocka *excl=&(top->excls);
1733 gmx_bool bBox,bTriclinic;
1736 rlist2 = sqr(fr->rlist);
1738 bBox = (fr->ePBC != epbcNONE);
1741 for(m=0; (m<DIM); m++)
1743 b_inv[m] = divide_err(1.0,box_size[m]);
1745 bTriclinic = TRICLINIC(box);
1752 cginfo = fr->cginfo;
1755 for (icg=fr->cg0; (icg<fr->hcg); icg++)
1758 i0 = cgs->index[icg];
1759 nri = cgs->index[icg+1]-i0;
1760 i_atoms = &(cgs->a[i0]);
1761 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1762 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1764 igid = GET_CGINFO_GID(cginfo[icg]);
1765 i_egp_flags = fr->egp_flags + ngid*igid;
1766 setexcl(cgs->index[icg],cgs->index[icg+1],excl,TRUE,bexcl);
1768 naaj=calc_naaj(icg,cgs->nr);
1771 ns_inner_tric(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1772 box,b_inv,rlist2,cgs,ns_buf,
1773 bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1777 ns_inner_rect(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1778 bBox,box_size,b_inv,rlist2,cgs,ns_buf,
1779 bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1783 for(nn=0; (nn<ngid); nn++)
1785 for(k=0; (k<SHIFTS); k++)
1787 nsbuf = &(ns_buf[nn][k]);
1790 put_in_list(bHaveVdW,ngid,md,icg,nn,nsbuf->ncg,nsbuf->jcg,
1791 cgs->index,bexcl,k,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
1792 nsbuf->ncg=nsbuf->nj=0;
1796 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1797 setexcl(cgs->index[icg],cgs->index[icg+1],excl,FALSE,bexcl);
1799 close_neighbor_lists(fr,FALSE);
1804 /************************************************
1806 * N S 5 G R I D S T U F F
1808 ************************************************/
1810 static inline void get_dx(int Nx,real gridx,real rc2,int xgi,real x,
1811 int *dx0,int *dx1,real *dcx2)
1839 for(i=xgi0; i>=0; i--)
1841 dcx = (i+1)*gridx-x;
1848 for(i=xgi1; i<Nx; i++)
1861 static inline void get_dx_dd(int Nx,real gridx,real rc2,int xgi,real x,
1862 int ncpddc,int shift_min,int shift_max,
1863 int *g0,int *g1,real *dcx2)
1866 int g_min,g_max,shift_home;
1899 g_min = (shift_min == shift_home ? 0 : ncpddc);
1900 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1907 else if (shift_max < 0)
1922 /* Check one grid cell down */
1923 dcx = ((*g0 - 1) + 1)*gridx - x;
1935 /* Check one grid cell up */
1936 dcx = (*g1 + 1)*gridx - x;
1948 #define sqr(x) ((x)*(x))
1949 #define calc_dx2(XI,YI,ZI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
1950 #define calc_cyl_dx2(XI,YI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
1951 /****************************************************
1953 * F A S T N E I G H B O R S E A R C H I N G
1955 * Optimized neighboursearching routine using grid
1956 * at least 1x1x1, see GROMACS manual
1958 ****************************************************/
1961 static void get_cutoff2(t_forcerec *fr,gmx_bool bDoLongRange,
1962 real *rvdw2,real *rcoul2,
1963 real *rs2,real *rm2,real *rl2)
1965 *rs2 = sqr(fr->rlist);
1967 if (bDoLongRange && fr->bTwinRange)
1969 /* The VdW and elec. LR cut-off's could be different,
1970 * so we can not simply set them to rlistlong.
1972 if (EVDW_MIGHT_BE_ZERO_AT_CUTOFF(fr->vdwtype) &&
1973 fr->rvdw > fr->rlist)
1975 *rvdw2 = sqr(fr->rlistlong);
1979 *rvdw2 = sqr(fr->rvdw);
1981 if (EEL_MIGHT_BE_ZERO_AT_CUTOFF(fr->eeltype) &&
1982 fr->rcoulomb > fr->rlist)
1984 *rcoul2 = sqr(fr->rlistlong);
1988 *rcoul2 = sqr(fr->rcoulomb);
1993 /* Workaround for a gcc -O3 or -ffast-math problem */
1997 *rm2 = min(*rvdw2,*rcoul2);
1998 *rl2 = max(*rvdw2,*rcoul2);
2001 static void init_nsgrid_lists(t_forcerec *fr,int ngid,gmx_ns_t *ns)
2003 real rvdw2,rcoul2,rs2,rm2,rl2;
2006 get_cutoff2(fr,TRUE,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
2008 /* Short range buffers */
2009 snew(ns->nl_sr,ngid);
2012 snew(ns->nlr_ljc,ngid);
2013 snew(ns->nlr_one,ngid);
2015 /* Always allocate both list types, since rcoulomb might now change with PME load balancing */
2016 /* Long range VdW and Coul buffers */
2017 snew(ns->nl_lr_ljc,ngid);
2018 /* Long range VdW or Coul only buffers */
2019 snew(ns->nl_lr_one,ngid);
2021 for(j=0; (j<ngid); j++) {
2022 snew(ns->nl_sr[j],MAX_CG);
2023 snew(ns->nl_lr_ljc[j],MAX_CG);
2024 snew(ns->nl_lr_one[j],MAX_CG);
2029 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
2034 static int nsgrid_core(FILE *log,t_commrec *cr,t_forcerec *fr,
2035 matrix box,rvec box_size,int ngid,
2036 gmx_localtop_t *top,
2037 t_grid *grid,rvec x[],
2038 t_excl bexcl[],gmx_bool *bExcludeAlleg,
2039 t_nrnb *nrnb,t_mdatoms *md,
2040 real *lambda,real *dvdlambda,
2041 gmx_grppairener_t *grppener,
2042 put_in_list_t *put_in_list,
2043 gmx_bool bHaveVdW[],
2044 gmx_bool bDoLongRange,gmx_bool bMakeQMMMnblist)
2047 atom_id **nl_lr_ljc,**nl_lr_one,**nl_sr;
2048 int *nlr_ljc,*nlr_one,*nsr;
2049 gmx_domdec_t *dd=NULL;
2050 t_block *cgs=&(top->cgs);
2051 int *cginfo=fr->cginfo;
2052 /* atom_id *i_atoms,*cgsindex=cgs->index; */
2054 int cell_x,cell_y,cell_z;
2055 int d,tx,ty,tz,dx,dy,dz,cj;
2056 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2057 int zsh_ty,zsh_tx,ysh_tx;
2059 int dx0,dx1,dy0,dy1,dz0,dz1;
2060 int Nx,Ny,Nz,shift=-1,j,nrj,nns,nn=-1;
2061 real gridx,gridy,gridz,grid_x,grid_y,grid_z;
2062 real *dcx2,*dcy2,*dcz2;
2064 int cg0,cg1,icg=-1,cgsnr,i0,igid,nri,naaj,max_jcg;
2065 int jcg0,jcg1,jjcg,cgj0,jgid;
2066 int *grida,*gridnra,*gridind;
2067 gmx_bool rvdw_lt_rcoul,rcoul_lt_rvdw;
2068 rvec xi,*cgcm,grid_offset;
2069 real r2,rs2,rvdw2,rcoul2,rm2,rl2,XI,YI,ZI,dcx,dcy,dcz,tmp1,tmp2;
2071 gmx_bool bDomDec,bTriclinicX,bTriclinicY;
2076 bDomDec = DOMAINDECOMP(cr);
2082 bTriclinicX = ((YY < grid->npbcdim &&
2083 (!bDomDec || dd->nc[YY]==1) && box[YY][XX] != 0) ||
2084 (ZZ < grid->npbcdim &&
2085 (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][XX] != 0));
2086 bTriclinicY = (ZZ < grid->npbcdim &&
2087 (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][YY] != 0);
2091 get_cutoff2(fr,bDoLongRange,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
2093 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
2094 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
2096 if (bMakeQMMMnblist)
2104 nl_lr_ljc = ns->nl_lr_ljc;
2105 nl_lr_one = ns->nl_lr_one;
2106 nlr_ljc = ns->nlr_ljc;
2107 nlr_one = ns->nlr_one;
2115 gridind = grid->index;
2116 gridnra = grid->nra;
2119 gridx = grid->cell_size[XX];
2120 gridy = grid->cell_size[YY];
2121 gridz = grid->cell_size[ZZ];
2125 copy_rvec(grid->cell_offset,grid_offset);
2126 copy_ivec(grid->ncpddc,ncpddc);
2131 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2132 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
2133 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
2134 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
2135 if (zsh_tx!=0 && ysh_tx!=0)
2137 /* This could happen due to rounding, when both ratios are 0.5 */
2146 /* We only want a list for the test particle */
2155 /* Set the shift range */
2156 for(d=0; d<DIM; d++)
2160 /* Check if we need periodicity shifts.
2161 * Without PBC or with domain decomposition we don't need them.
2163 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2170 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
2181 /* Loop over charge groups */
2182 for(icg=cg0; (icg < cg1); icg++)
2184 igid = GET_CGINFO_GID(cginfo[icg]);
2185 /* Skip this charge group if all energy groups are excluded! */
2186 if (bExcludeAlleg[igid])
2191 i0 = cgs->index[icg];
2193 if (bMakeQMMMnblist)
2195 /* Skip this charge group if it is not a QM atom while making a
2196 * QM/MM neighbourlist
2198 if (md->bQM[i0]==FALSE)
2200 continue; /* MM particle, go to next particle */
2203 /* Compute the number of charge groups that fall within the control
2206 naaj = calc_naaj(icg,cgsnr);
2213 /* make a normal neighbourlist */
2217 /* Get the j charge-group and dd cell shift ranges */
2218 dd_get_ns_ranges(cr->dd,icg,&jcg0,&jcg1,sh0,sh1);
2223 /* Compute the number of charge groups that fall within the control
2226 naaj = calc_naaj(icg,cgsnr);
2232 /* The i-particle is awlways the test particle,
2233 * so we want all j-particles
2235 max_jcg = cgsnr - 1;
2239 max_jcg = jcg1 - cgsnr;
2244 i_egp_flags = fr->egp_flags + igid*ngid;
2246 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2247 setexcl(i0,cgs->index[icg+1],&top->excls,TRUE,bexcl);
2249 ci2xyz(grid,icg,&cell_x,&cell_y,&cell_z);
2251 /* Changed iicg to icg, DvdS 990115
2252 * (but see consistency check above, DvdS 990330)
2255 fprintf(log,"icg=%5d, naaj=%5d, cell %d %d %d\n",
2256 icg,naaj,cell_x,cell_y,cell_z);
2258 /* Loop over shift vectors in three dimensions */
2259 for (tz=-shp[ZZ]; tz<=shp[ZZ]; tz++)
2261 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2262 /* Calculate range of cells in Z direction that have the shift tz */
2263 zgi = cell_z + tz*Nz;
2266 get_dx(Nz,gridz,rl2,zgi,ZI,&dz0,&dz1,dcz2);
2268 get_dx_dd(Nz,gridz,rl2,zgi,ZI-grid_offset[ZZ],
2269 ncpddc[ZZ],sh0[ZZ],sh1[ZZ],&dz0,&dz1,dcz2);
2275 for (ty=-shp[YY]; ty<=shp[YY]; ty++)
2277 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2278 /* Calculate range of cells in Y direction that have the shift ty */
2281 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2285 ygi = cell_y + ty*Ny;
2288 get_dx(Ny,gridy,rl2,ygi,YI,&dy0,&dy1,dcy2);
2290 get_dx_dd(Ny,gridy,rl2,ygi,YI-grid_offset[YY],
2291 ncpddc[YY],sh0[YY],sh1[YY],&dy0,&dy1,dcy2);
2297 for (tx=-shp[XX]; tx<=shp[XX]; tx++)
2299 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2300 /* Calculate range of cells in X direction that have the shift tx */
2303 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2307 xgi = cell_x + tx*Nx;
2310 get_dx(Nx,gridx,rl2,xgi*Nx,XI,&dx0,&dx1,dcx2);
2312 get_dx_dd(Nx,gridx,rl2,xgi,XI-grid_offset[XX],
2313 ncpddc[XX],sh0[XX],sh1[XX],&dx0,&dx1,dcx2);
2319 /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2320 * from the neigbour list as it will not interact */
2321 if (fr->adress_type != eAdressOff){
2322 if (md->wf[cgs->index[icg]]==0 && egp_explicit(fr, igid)){
2326 /* Get shift vector */
2327 shift=XYZ2IS(tx,ty,tz);
2329 range_check(shift,0,SHIFTS);
2331 for(nn=0; (nn<ngid); nn++)
2338 fprintf(log,"shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2339 shift,dx0,dx1,dy0,dy1,dz0,dz1);
2340 fprintf(log,"cgcm: %8.3f %8.3f %8.3f\n",cgcm[icg][XX],
2341 cgcm[icg][YY],cgcm[icg][ZZ]);
2342 fprintf(log,"xi: %8.3f %8.3f %8.3f\n",XI,YI,ZI);
2344 for (dx=dx0; (dx<=dx1); dx++)
2346 tmp1 = rl2 - dcx2[dx];
2347 for (dy=dy0; (dy<=dy1); dy++)
2349 tmp2 = tmp1 - dcy2[dy];
2352 for (dz=dz0; (dz<=dz1); dz++) {
2353 if (tmp2 > dcz2[dz]) {
2354 /* Find grid-cell cj in which possible neighbours are */
2355 cj = xyz2ci(Ny,Nz,dx,dy,dz);
2357 /* Check out how many cgs (nrj) there in this cell */
2360 /* Find the offset in the cg list */
2363 /* Check if all j's are out of range so we
2364 * can skip the whole cell.
2365 * Should save some time, especially with DD.
2368 (grida[cgj0] >= max_jcg &&
2369 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2375 for (j=0; (j<nrj); j++)
2377 jjcg = grida[cgj0+j];
2379 /* check whether this guy is in range! */
2380 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2383 r2=calc_dx2(XI,YI,ZI,cgcm[jjcg]);
2385 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2386 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2387 /* check energy group exclusions */
2388 if (!(i_egp_flags[jgid] & EGP_EXCL))
2392 if (nsr[jgid] >= MAX_CG)
2394 /* Add to short-range list */
2395 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2396 nsr[jgid],nl_sr[jgid],
2397 cgs->index,/* cgsatoms, */ bexcl,
2398 shift,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
2401 nl_sr[jgid][nsr[jgid]++]=jjcg;
2405 if (nlr_ljc[jgid] >= MAX_CG)
2407 /* Add to LJ+coulomb long-range list */
2408 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2409 nlr_ljc[jgid],nl_lr_ljc[jgid],top->cgs.index,
2410 bexcl,shift,fr,TRUE,TRUE,TRUE,fr->solvent_opt);
2413 nl_lr_ljc[jgid][nlr_ljc[jgid]++]=jjcg;
2417 if (nlr_one[jgid] >= MAX_CG)
2419 /* Add to long-range list with only coul, or only LJ */
2420 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2421 nlr_one[jgid],nl_lr_one[jgid],top->cgs.index,
2422 bexcl,shift,fr,TRUE,rvdw_lt_rcoul,rcoul_lt_rvdw,fr->solvent_opt);
2425 nl_lr_one[jgid][nlr_one[jgid]++]=jjcg;
2437 /* CHECK whether there is anything left in the buffers */
2438 for(nn=0; (nn<ngid); nn++)
2442 put_in_list(bHaveVdW,ngid,md,icg,nn,nsr[nn],nl_sr[nn],
2443 cgs->index, /* cgsatoms, */ bexcl,
2444 shift,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
2447 if (nlr_ljc[nn] > 0)
2449 put_in_list(bHaveVdW,ngid,md,icg,nn,nlr_ljc[nn],
2450 nl_lr_ljc[nn],top->cgs.index,
2451 bexcl,shift,fr,TRUE,TRUE,TRUE,fr->solvent_opt);
2454 if (nlr_one[nn] > 0)
2456 put_in_list(bHaveVdW,ngid,md,icg,nn,nlr_one[nn],
2457 nl_lr_one[nn],top->cgs.index,
2458 bexcl,shift,fr,TRUE,rvdw_lt_rcoul,rcoul_lt_rvdw,fr->solvent_opt);
2464 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2465 setexcl(cgs->index[icg],cgs->index[icg+1],&top->excls,FALSE,bexcl);
2467 /* No need to perform any left-over force calculations anymore (as we used to do here)
2468 * since we now save the proper long-range lists for later evaluation.
2473 /* Close neighbourlists */
2474 close_neighbor_lists(fr,bMakeQMMMnblist);
2479 void ns_realloc_natoms(gmx_ns_t *ns,int natoms)
2483 if (natoms > ns->nra_alloc)
2485 ns->nra_alloc = over_alloc_dd(natoms);
2486 srenew(ns->bexcl,ns->nra_alloc);
2487 for(i=0; i<ns->nra_alloc; i++)
2494 void init_ns(FILE *fplog,const t_commrec *cr,
2495 gmx_ns_t *ns,t_forcerec *fr,
2496 const gmx_mtop_t *mtop,
2499 int mt,icg,nr_in_cg,maxcg,i,j,jcg,ngid,ncg;
2503 /* Compute largest charge groups size (# atoms) */
2505 for(mt=0; mt<mtop->nmoltype; mt++) {
2506 cgs = &mtop->moltype[mt].cgs;
2507 for (icg=0; (icg < cgs->nr); icg++)
2509 nr_in_cg=max(nr_in_cg,(int)(cgs->index[icg+1]-cgs->index[icg]));
2513 /* Verify whether largest charge group is <= max cg.
2514 * This is determined by the type of the local exclusion type
2515 * Exclusions are stored in bits. (If the type is not large
2516 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2518 maxcg = sizeof(t_excl)*8;
2519 if (nr_in_cg > maxcg)
2521 gmx_fatal(FARGS,"Max #atoms in a charge group: %d > %d\n",
2525 ngid = mtop->groups.grps[egcENER].nr;
2526 snew(ns->bExcludeAlleg,ngid);
2527 for(i=0; i<ngid; i++) {
2528 ns->bExcludeAlleg[i] = TRUE;
2529 for(j=0; j<ngid; j++)
2531 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2533 ns->bExcludeAlleg[i] = FALSE;
2540 ns->grid = init_grid(fplog,fr);
2541 init_nsgrid_lists(fr,ngid,ns);
2546 snew(ns->ns_buf,ngid);
2547 for(i=0; (i<ngid); i++)
2549 snew(ns->ns_buf[i],SHIFTS);
2551 ncg = ncg_mtop(mtop);
2552 snew(ns->simple_aaj,2*ncg);
2553 for(jcg=0; (jcg<ncg); jcg++)
2555 ns->simple_aaj[jcg] = jcg;
2556 ns->simple_aaj[jcg+ncg] = jcg;
2560 /* Create array that determines whether or not atoms have VdW */
2561 snew(ns->bHaveVdW,fr->ntype);
2562 for(i=0; (i<fr->ntype); i++)
2564 for(j=0; (j<fr->ntype); j++)
2566 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2568 ((BHAMA(fr->nbfp,fr->ntype,i,j) != 0) ||
2569 (BHAMB(fr->nbfp,fr->ntype,i,j) != 0) ||
2570 (BHAMC(fr->nbfp,fr->ntype,i,j) != 0)) :
2571 ((C6(fr->nbfp,fr->ntype,i,j) != 0) ||
2572 (C12(fr->nbfp,fr->ntype,i,j) != 0))));
2576 pr_bvec(debug,0,"bHaveVdW",ns->bHaveVdW,fr->ntype,TRUE);
2580 if (!DOMAINDECOMP(cr))
2582 /* This could be reduced with particle decomposition */
2583 ns_realloc_natoms(ns,mtop->natoms);
2586 ns->nblist_initialized=FALSE;
2588 /* nbr list debug dump */
2590 char *ptr=getenv("GMX_DUMP_NL");
2593 ns->dump_nl=strtol(ptr,NULL,10);
2596 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2607 int search_neighbours(FILE *log,t_forcerec *fr,
2608 rvec x[],matrix box,
2609 gmx_localtop_t *top,
2610 gmx_groups_t *groups,
2612 t_nrnb *nrnb,t_mdatoms *md,
2613 real *lambda,real *dvdlambda,
2614 gmx_grppairener_t *grppener,
2616 gmx_bool bDoLongRangeNS,
2617 gmx_bool bPadListsForKernels)
2619 t_block *cgs=&(top->cgs);
2620 rvec box_size,grid_x0,grid_x1;
2622 real min_size,grid_dens;
2626 gmx_bool *i_egp_flags;
2627 int cg_start,cg_end,start,end;
2630 gmx_domdec_zones_t *dd_zones;
2631 put_in_list_t *put_in_list;
2635 /* Set some local variables */
2637 ngid = groups->grps[egcENER].nr;
2639 for(m=0; (m<DIM); m++)
2641 box_size[m] = box[m][m];
2644 if (fr->ePBC != epbcNONE)
2646 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC,box))
2648 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.");
2652 min_size = min(box_size[XX],min(box_size[YY],box_size[ZZ]));
2653 if (2*fr->rlistlong >= min_size)
2654 gmx_fatal(FARGS,"One of the box diagonal elements has become smaller than twice the cut-off length.");
2658 if (DOMAINDECOMP(cr))
2660 ns_realloc_natoms(ns,cgs->index[cgs->nr]);
2664 /* Reset the neighbourlists */
2665 reset_neighbor_lists(fr,TRUE,TRUE);
2667 if (bGrid && bFillGrid)
2671 if (DOMAINDECOMP(cr))
2673 dd_zones = domdec_zones(cr->dd);
2679 get_nsgrid_boundaries(grid->nboundeddim,box,NULL,NULL,NULL,NULL,
2680 cgs->nr,fr->cg_cm,grid_x0,grid_x1,&grid_dens);
2682 grid_first(log,grid,NULL,NULL,fr->ePBC,box,grid_x0,grid_x1,
2683 fr->rlistlong,grid_dens);
2687 /* Don't know why this all is... (DvdS 3/99) */
2693 end = (cgs->nr+1)/2;
2696 if (DOMAINDECOMP(cr))
2699 fill_grid(log,dd_zones,grid,end,-1,end,fr->cg_cm);
2701 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2705 fill_grid(log,NULL,grid,cgs->nr,fr->cg0,fr->hcg,fr->cg_cm);
2706 grid->icg0 = fr->cg0;
2707 grid->icg1 = fr->hcg;
2715 calc_elemnr(log,grid,start,end,cgs->nr);
2717 grid_last(log,grid,start,end,cgs->nr);
2721 check_grid(debug,grid);
2722 print_grid(debug,grid);
2727 /* Set the grid cell index for the test particle only.
2728 * The cell to cg index is not corrected, but that does not matter.
2730 fill_grid(log,NULL,ns->grid,fr->hcg,fr->hcg-1,fr->hcg,fr->cg_cm);
2734 if (fr->adress_type == eAdressOff){
2735 if (!fr->ns.bCGlist)
2737 put_in_list = put_in_list_at;
2741 put_in_list = put_in_list_cg;
2744 put_in_list = put_in_list_adress;
2751 nsearch = nsgrid_core(log,cr,fr,box,box_size,ngid,top,
2752 grid,x,ns->bexcl,ns->bExcludeAlleg,
2753 nrnb,md,lambda,dvdlambda,grppener,
2754 put_in_list,ns->bHaveVdW,
2755 bDoLongRangeNS,FALSE);
2757 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2758 * the classical calculation. The charge-charge interaction
2759 * between QM and MM atoms is handled in the QMMM core calculation
2760 * (see QMMM.c). The VDW however, we'd like to compute classically
2761 * and the QM MM atom pairs have just been put in the
2762 * corresponding neighbourlists. in case of QMMM we still need to
2763 * fill a special QMMM neighbourlist that contains all neighbours
2764 * of the QM atoms. If bQMMM is true, this list will now be made:
2766 if (fr->bQMMM && fr->qr->QMMMscheme!=eQMMMschemeoniom)
2768 nsearch += nsgrid_core(log,cr,fr,box,box_size,ngid,top,
2769 grid,x,ns->bexcl,ns->bExcludeAlleg,
2770 nrnb,md,lambda,dvdlambda,grppener,
2771 put_in_list_qmmm,ns->bHaveVdW,
2772 bDoLongRangeNS,TRUE);
2777 nsearch = ns_simple_core(fr,top,md,box,box_size,
2778 ns->bexcl,ns->simple_aaj,
2779 ngid,ns->ns_buf,put_in_list,ns->bHaveVdW);
2787 inc_nrnb(nrnb,eNR_NS,nsearch);
2788 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2793 int natoms_beyond_ns_buffer(t_inputrec *ir,t_forcerec *fr,t_block *cgs,
2794 matrix scale_tot,rvec *x)
2796 int cg0,cg1,cg,a0,a1,a,i,j;
2797 real rint,hbuf2,scale;
2799 gmx_bool bIsotropic;
2804 rint = max(ir->rcoulomb,ir->rvdw);
2805 if (ir->rlist < rint)
2807 gmx_fatal(FARGS,"The neighbor search buffer has negative size: %f nm",
2815 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2817 hbuf2 = sqr(0.5*(ir->rlist - rint));
2818 for(cg=cg0; cg<cg1; cg++)
2820 a0 = cgs->index[cg];
2821 a1 = cgs->index[cg+1];
2822 for(a=a0; a<a1; a++)
2824 if (distance2(cg_cm[cg],x[a]) > hbuf2)
2834 scale = scale_tot[0][0];
2835 for(i=1; i<DIM; i++)
2837 /* With anisotropic scaling, the original spherical ns volumes become
2838 * ellipsoids. To avoid costly transformations we use the minimum
2839 * eigenvalue of the scaling matrix for determining the buffer size.
2840 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2842 scale = min(scale,scale_tot[i][i]);
2843 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2849 if (scale_tot[i][j] != 0)
2855 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2858 for(cg=cg0; cg<cg1; cg++)
2860 svmul(scale,cg_cm[cg],cgsc);
2861 a0 = cgs->index[cg];
2862 a1 = cgs->index[cg+1];
2863 for(a=a0; a<a1; a++)
2865 if (distance2(cgsc,x[a]) > hbuf2)
2874 /* Anistropic scaling */
2875 for(cg=cg0; cg<cg1; cg++)
2877 /* Since scale_tot contains the transpose of the scaling matrix,
2878 * we need to multiply with the transpose.
2880 tmvmul_ur0(scale_tot,cg_cm[cg],cgsc);
2881 a0 = cgs->index[cg];
2882 a1 = cgs->index[cg+1];
2883 for(a=a0; a<a1; a++)
2885 if (distance2(cgsc,x[a]) > hbuf2)