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
40 #ifdef GMX_THREAD_SHM_FDECOMP
54 #include "nonbonded.h"
58 #include "gmx_fatal.h"
61 #include "mtop_util.h"
68 * E X C L U S I O N H A N D L I N G
72 static void SETEXCL_(t_excl e[],atom_id i,atom_id j)
73 { e[j] = e[j] | (1<<i); }
74 static void RMEXCL_(t_excl e[],atom_id i,atom_id j)
75 { e[j]=e[j] & ~(1<<i); }
76 static gmx_bool ISEXCL_(t_excl e[],atom_id i,atom_id j)
77 { return (gmx_bool)(e[j] & (1<<i)); }
78 static gmx_bool NOTEXCL_(t_excl e[],atom_id i,atom_id j)
79 { return !(ISEXCL(e,i,j)); }
81 #define SETEXCL(e,i,j) (e)[((atom_id) (j))] |= (1<<((atom_id) (i)))
82 #define RMEXCL(e,i,j) (e)[((atom_id) (j))] &= (~(1<<((atom_id) (i))))
83 #define ISEXCL(e,i,j) (gmx_bool) ((e)[((atom_id) (j))] & (1<<((atom_id) (i))))
84 #define NOTEXCL(e,i,j) !(ISEXCL(e,i,j))
87 /************************************************
89 * U T I L I T I E S F O R N S
91 ************************************************/
93 static void reallocate_nblist(t_nblist *nl)
97 fprintf(debug,"reallocating neigborlist (ielec=%d, ivdw=%d, igeometry=%d, free_energy=%d), maxnri=%d\n",
98 nl->ielec,nl->ivdw,nl->igeometry,nl->free_energy,nl->maxnri);
100 srenew(nl->iinr, nl->maxnri);
101 if (nl->igeometry == GMX_NBLIST_GEOMETRY_CG_CG)
103 srenew(nl->iinr_end,nl->maxnri);
105 srenew(nl->gid, nl->maxnri);
106 srenew(nl->shift, nl->maxnri);
107 srenew(nl->jindex, nl->maxnri+1);
111 static void init_nblist(FILE *log, t_nblist *nl_sr,t_nblist *nl_lr,
113 int ivdw, int ivdwmod,
114 int ielec, int ielecmod,
115 gmx_bool bfree, int igeometry)
123 nl = (i == 0) ? nl_sr : nl_lr;
124 homenr = (i == 0) ? maxsr : maxlr;
132 /* Set coul/vdw in neighborlist, and for the normal loops we determine
133 * an index of which one to call.
136 nl->ivdwmod = ivdwmod;
138 nl->ielecmod = ielecmod;
139 nl->free_energy = bfree;
140 nl->igeometry = igeometry;
144 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
147 gmx_nonbonded_set_kernel_pointers( (i==0) ? log : NULL,nl);
149 /* maxnri is influenced by the number of shifts (maximum is 8)
150 * and the number of energy groups.
151 * If it is not enough, nl memory will be reallocated during the run.
152 * 4 seems to be a reasonable factor, which only causes reallocation
153 * during runs with tiny and many energygroups.
155 nl->maxnri = homenr*4;
164 reallocate_nblist(nl);
169 fprintf(debug,"Initiating neighbourlist (ielec=%d, ivdw=%d, free=%d) for %s interactions,\nwith %d SR, %d LR atoms.\n",
170 nl->ielec,nl->ivdw,nl->free_energy,gmx_nblist_geometry_names[nl->igeometry],maxsr,maxlr);
173 #ifdef GMX_THREAD_SHM_FDECOMP
176 pthread_mutex_init(nl->mtx,NULL);
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 nlist->jindex[nri+1] = nlist->nrj;
399 len=nlist->nrj - nlist->jindex[nri];
401 /* nlist length for water i molecules is treated statically
404 if (len > nlist->maxlen)
411 static inline void close_nblist(t_nblist *nlist)
413 /* Only close this nblist when it has been initialized.
414 * Avoid the creation of i-lists with no j-particles.
418 /* Some assembly kernels do not support empty lists,
419 * make sure here that we don't generate any empty lists.
420 * With the current ns code this branch is taken in two cases:
421 * No i-particles at all: nri=-1 here
422 * There are i-particles, but no j-particles; nri=0 here
428 /* Close list number nri by incrementing the count */
433 static inline void close_neighbor_lists(t_forcerec *fr,gmx_bool bMakeQMMMnblist)
439 close_nblist(&(fr->QMMMlist));
442 for(n=0; n<fr->nnblists; n++)
444 for(i=0; (i<eNL_NR); i++)
446 close_nblist(&(fr->nblists[n].nlist_sr[i]));
447 close_nblist(&(fr->nblists[n].nlist_lr[i]));
453 static inline void add_j_to_nblist(t_nblist *nlist,atom_id j_atom,gmx_bool bLR)
457 if (nlist->nrj >= nlist->maxnrj)
459 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
461 fprintf(debug,"Increasing %s nblist (ielec=%d,ivdw=%d,free=%d,igeometry=%d) j size to %d\n",
462 bLR ? "LR" : "SR",nlist->ielec,nlist->ivdw,nlist->free_energy,nlist->igeometry,nlist->maxnrj);
464 srenew(nlist->jjnr,nlist->maxnrj);
467 nlist->jjnr[nrj] = j_atom;
471 static inline void add_j_to_nblist_cg(t_nblist *nlist,
472 atom_id j_start,int j_end,
473 t_excl *bexcl,gmx_bool i_is_j,
479 if (nlist->nrj >= nlist->maxnrj)
481 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
483 fprintf(debug,"Increasing %s nblist (ielec=%d,ivdw=%d,free=%d,igeometry=%d) j size to %d\n",
484 bLR ? "LR" : "SR",nlist->ielec,nlist->ivdw,nlist->free_energy,nlist->igeometry,nlist->maxnrj);
486 srenew(nlist->jjnr ,nlist->maxnrj);
487 srenew(nlist->jjnr_end,nlist->maxnrj);
488 srenew(nlist->excl ,nlist->maxnrj*MAX_CGCGSIZE);
491 nlist->jjnr[nrj] = j_start;
492 nlist->jjnr_end[nrj] = j_end;
494 if (j_end - j_start > MAX_CGCGSIZE)
496 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);
499 /* Set the exclusions */
500 for(j=j_start; j<j_end; j++)
502 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
506 /* Avoid double counting of intra-cg interactions */
507 for(j=1; j<j_end-j_start; j++)
509 nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
517 put_in_list_t(gmx_bool bHaveVdW[],
534 put_in_list_at(gmx_bool bHaveVdW[],
550 /* The a[] index has been removed,
551 * to put it back in i_atom should be a[i0] and jj should be a[jj].
556 t_nblist * vdwc_free = NULL;
557 t_nblist * vdw_free = NULL;
558 t_nblist * coul_free = NULL;
559 t_nblist * vdwc_ww = NULL;
560 t_nblist * coul_ww = NULL;
562 int i,j,jcg,igid,gid,nbl_ind,ind_ij;
563 atom_id jj,jj0,jj1,i_atom;
568 real *charge,*chargeB;
570 gmx_bool bFreeEnergy,bFree,bFreeJ,bNotEx,*bPert;
571 gmx_bool bDoVdW_i,bDoCoul_i,bDoCoul_i_sol;
575 /* Copy some pointers */
577 charge = md->chargeA;
578 chargeB = md->chargeB;
581 bPert = md->bPerturbed;
585 nicg = index[icg+1]-i0;
587 /* Get the i charge group info */
588 igid = GET_CGINFO_GID(cginfo[icg]);
590 iwater = (solvent_opt!=esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
595 /* Check if any of the particles involved are perturbed.
596 * If not we can do the cheaper normal put_in_list
597 * and use more solvent optimization.
599 for(i=0; i<nicg; i++)
601 bFreeEnergy |= bPert[i0+i];
603 /* Loop over the j charge groups */
604 for(j=0; (j<nj && !bFreeEnergy); j++)
609 /* Finally loop over the atoms in the j-charge group */
610 for(jj=jj0; jj<jj1; jj++)
612 bFreeEnergy |= bPert[jj];
617 /* Unpack pointers to neighbourlist structs */
618 if (fr->nnblists == 1)
624 nbl_ind = fr->gid2nblists[GID(igid,jgid,ngid)];
628 nlist = fr->nblists[nbl_ind].nlist_lr;
632 nlist = fr->nblists[nbl_ind].nlist_sr;
635 if (iwater != esolNO)
637 vdwc = &nlist[eNL_VDWQQ_WATER];
638 vdw = &nlist[eNL_VDW];
639 coul = &nlist[eNL_QQ_WATER];
640 #ifndef DISABLE_WATERWATER_NLIST
641 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
642 coul_ww = &nlist[eNL_QQ_WATERWATER];
647 vdwc = &nlist[eNL_VDWQQ];
648 vdw = &nlist[eNL_VDW];
649 coul = &nlist[eNL_QQ];
654 if (iwater != esolNO)
656 /* Loop over the atoms in the i charge group */
658 gid = GID(igid,jgid,ngid);
659 /* Create new i_atom for each energy group */
660 if (bDoCoul && bDoVdW)
662 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
663 #ifndef DISABLE_WATERWATER_NLIST
664 new_i_nblist(vdwc_ww,bLR,i_atom,shift,gid);
669 new_i_nblist(vdw,bLR,i_atom,shift,gid);
673 new_i_nblist(coul,bLR,i_atom,shift,gid);
674 #ifndef DISABLE_WATERWATER_NLIST
675 new_i_nblist(coul_ww,bLR,i_atom,shift,gid);
678 /* Loop over the j charge groups */
679 for(j=0; (j<nj); j++)
689 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
691 if (iwater == esolSPC && jwater == esolSPC)
693 /* Interaction between two SPC molecules */
696 /* VdW only - only first atoms in each water interact */
697 add_j_to_nblist(vdw,jj0,bLR);
701 #ifdef DISABLE_WATERWATER_NLIST
702 /* Add entries for the three atoms - only do VdW if we need to */
705 add_j_to_nblist(coul,jj0,bLR);
709 add_j_to_nblist(vdwc,jj0,bLR);
711 add_j_to_nblist(coul,jj0+1,bLR);
712 add_j_to_nblist(coul,jj0+2,bLR);
714 /* One entry for the entire water-water interaction */
717 add_j_to_nblist(coul_ww,jj0,bLR);
721 add_j_to_nblist(vdwc_ww,jj0,bLR);
726 else if (iwater == esolTIP4P && jwater == esolTIP4P)
728 /* Interaction between two TIP4p molecules */
731 /* VdW only - only first atoms in each water interact */
732 add_j_to_nblist(vdw,jj0,bLR);
736 #ifdef DISABLE_WATERWATER_NLIST
737 /* Add entries for the four atoms - only do VdW if we need to */
740 add_j_to_nblist(vdw,jj0,bLR);
742 add_j_to_nblist(coul,jj0+1,bLR);
743 add_j_to_nblist(coul,jj0+2,bLR);
744 add_j_to_nblist(coul,jj0+3,bLR);
746 /* One entry for the entire water-water interaction */
749 add_j_to_nblist(coul_ww,jj0,bLR);
753 add_j_to_nblist(vdwc_ww,jj0,bLR);
760 /* j charge group is not water, but i is.
761 * Add entries to the water-other_atom lists; the geometry of the water
762 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
763 * so we don't care if it is SPC or TIP4P...
770 for(jj=jj0; (jj<jj1); jj++)
774 add_j_to_nblist(coul,jj,bLR);
780 for(jj=jj0; (jj<jj1); jj++)
782 if (bHaveVdW[type[jj]])
784 add_j_to_nblist(vdw,jj,bLR);
790 /* _charge_ _groups_ interact with both coulomb and LJ */
791 /* Check which atoms we should add to the lists! */
792 for(jj=jj0; (jj<jj1); jj++)
794 if (bHaveVdW[type[jj]])
798 add_j_to_nblist(vdwc,jj,bLR);
802 add_j_to_nblist(vdw,jj,bLR);
805 else if (charge[jj] != 0)
807 add_j_to_nblist(coul,jj,bLR);
814 close_i_nblist(coul);
815 close_i_nblist(vdwc);
816 #ifndef DISABLE_WATERWATER_NLIST
817 close_i_nblist(coul_ww);
818 close_i_nblist(vdwc_ww);
823 /* no solvent as i charge group */
824 /* Loop over the atoms in the i charge group */
825 for(i=0; i<nicg; i++)
828 gid = GID(igid,jgid,ngid);
831 /* Create new i_atom for each energy group */
832 if (bDoVdW && bDoCoul)
834 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
838 new_i_nblist(vdw,bLR,i_atom,shift,gid);
842 new_i_nblist(coul,bLR,i_atom,shift,gid);
844 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
845 bDoCoul_i = (bDoCoul && qi!=0);
847 if (bDoVdW_i || bDoCoul_i)
849 /* Loop over the j charge groups */
850 for(j=0; (j<nj); j++)
854 /* Check for large charge groups */
865 /* Finally loop over the atoms in the j-charge group */
866 for(jj=jj0; jj<jj1; jj++)
868 bNotEx = NOTEXCL(bExcl,i,jj);
876 add_j_to_nblist(coul,jj,bLR);
881 if (bHaveVdW[type[jj]])
883 add_j_to_nblist(vdw,jj,bLR);
888 if (bHaveVdW[type[jj]])
892 add_j_to_nblist(vdwc,jj,bLR);
896 add_j_to_nblist(vdw,jj,bLR);
899 else if (charge[jj] != 0)
901 add_j_to_nblist(coul,jj,bLR);
909 close_i_nblist(coul);
910 close_i_nblist(vdwc);
916 /* we are doing free energy */
917 vdwc_free = &nlist[eNL_VDWQQ_FREE];
918 vdw_free = &nlist[eNL_VDW_FREE];
919 coul_free = &nlist[eNL_QQ_FREE];
920 /* Loop over the atoms in the i charge group */
921 for(i=0; i<nicg; i++)
924 gid = GID(igid,jgid,ngid);
926 qiB = chargeB[i_atom];
928 /* Create new i_atom for each energy group */
929 if (bDoVdW && bDoCoul)
930 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
932 new_i_nblist(vdw,bLR,i_atom,shift,gid);
934 new_i_nblist(coul,bLR,i_atom,shift,gid);
936 new_i_nblist(vdw_free,bLR,i_atom,shift,gid);
937 new_i_nblist(coul_free,bLR,i_atom,shift,gid);
938 new_i_nblist(vdwc_free,bLR,i_atom,shift,gid);
940 bDoVdW_i = (bDoVdW &&
941 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
942 bDoCoul_i = (bDoCoul && (qi!=0 || qiB!=0));
943 /* For TIP4P the first atom does not have a charge,
944 * but the last three do. So we should still put an atom
945 * without LJ but with charge in the water-atom neighborlist
946 * for a TIP4p i charge group.
947 * For SPC type water the first atom has LJ and charge,
948 * so there is no such problem.
950 if (iwater == esolNO)
952 bDoCoul_i_sol = bDoCoul_i;
956 bDoCoul_i_sol = bDoCoul;
959 if (bDoVdW_i || bDoCoul_i_sol)
961 /* Loop over the j charge groups */
962 for(j=0; (j<nj); j++)
966 /* Check for large charge groups */
977 /* Finally loop over the atoms in the j-charge group */
978 bFree = bPert[i_atom];
979 for(jj=jj0; (jj<jj1); jj++)
981 bFreeJ = bFree || bPert[jj];
982 /* Complicated if, because the water H's should also
983 * see perturbed j-particles
985 if (iwater==esolNO || i==0 || bFreeJ)
987 bNotEx = NOTEXCL(bExcl,i,jj);
995 if (charge[jj]!=0 || chargeB[jj]!=0)
997 add_j_to_nblist(coul_free,jj,bLR);
1000 else if (!bDoCoul_i)
1002 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1004 add_j_to_nblist(vdw_free,jj,bLR);
1009 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1011 if (charge[jj]!=0 || chargeB[jj]!=0)
1013 add_j_to_nblist(vdwc_free,jj,bLR);
1017 add_j_to_nblist(vdw_free,jj,bLR);
1020 else if (charge[jj]!=0 || chargeB[jj]!=0)
1021 add_j_to_nblist(coul_free,jj,bLR);
1026 /* This is done whether or not bWater is set */
1027 if (charge[jj] != 0)
1029 add_j_to_nblist(coul,jj,bLR);
1032 else if (!bDoCoul_i_sol)
1034 if (bHaveVdW[type[jj]])
1036 add_j_to_nblist(vdw,jj,bLR);
1041 if (bHaveVdW[type[jj]])
1043 if (charge[jj] != 0)
1045 add_j_to_nblist(vdwc,jj,bLR);
1049 add_j_to_nblist(vdw,jj,bLR);
1052 else if (charge[jj] != 0)
1054 add_j_to_nblist(coul,jj,bLR);
1062 close_i_nblist(vdw);
1063 close_i_nblist(coul);
1064 close_i_nblist(vdwc);
1065 close_i_nblist(vdw_free);
1066 close_i_nblist(coul_free);
1067 close_i_nblist(vdwc_free);
1073 put_in_list_qmmm(gmx_bool bHaveVdW[],
1090 int i,j,jcg,igid,gid;
1091 atom_id jj,jj0,jj1,i_atom;
1095 /* Get atom range */
1097 nicg = index[icg+1]-i0;
1099 /* Get the i charge group info */
1100 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1102 coul = &fr->QMMMlist;
1104 /* Loop over atoms in the ith charge group */
1105 for (i=0;i<nicg;i++)
1108 gid = GID(igid,jgid,ngid);
1109 /* Create new i_atom for each energy group */
1110 new_i_nblist(coul,bLR,i_atom,shift,gid);
1112 /* Loop over the j charge groups */
1117 /* Charge groups cannot have QM and MM atoms simultaneously */
1122 /* Finally loop over the atoms in the j-charge group */
1123 for(jj=jj0; jj<jj1; jj++)
1125 bNotEx = NOTEXCL(bExcl,i,jj);
1127 add_j_to_nblist(coul,jj,bLR);
1131 close_i_nblist(coul);
1136 put_in_list_cg(gmx_bool bHaveVdW[],
1153 int igid,gid,nbl_ind;
1157 cginfo = fr->cginfo[icg];
1159 igid = GET_CGINFO_GID(cginfo);
1160 gid = GID(igid,jgid,ngid);
1162 /* Unpack pointers to neighbourlist structs */
1163 if (fr->nnblists == 1)
1169 nbl_ind = fr->gid2nblists[gid];
1173 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1177 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1180 /* Make a new neighbor list for charge group icg.
1181 * Currently simply one neighbor list is made with LJ and Coulomb.
1182 * If required, zero interactions could be removed here
1183 * or in the force loop.
1185 new_i_nblist(vdwc,bLR,index[icg],shift,gid);
1186 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1188 for(j=0; (j<nj); j++)
1191 /* Skip the icg-icg pairs if all self interactions are excluded */
1192 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1194 /* Here we add the j charge group jcg to the list,
1195 * exclusions are also added to the list.
1197 add_j_to_nblist_cg(vdwc,index[jcg],index[jcg+1],bExcl,icg==jcg,bLR);
1201 close_i_nblist(vdwc);
1204 static void setexcl(atom_id start,atom_id end,t_blocka *excl,gmx_bool b,
1211 for(i=start; i<end; i++)
1213 for(k=excl->index[i]; k<excl->index[i+1]; k++)
1215 SETEXCL(bexcl,i-start,excl->a[k]);
1221 for(i=start; i<end; i++)
1223 for(k=excl->index[i]; k<excl->index[i+1]; k++)
1225 RMEXCL(bexcl,i-start,excl->a[k]);
1231 int calc_naaj(int icg,int cgtot)
1235 if ((cgtot % 2) == 1)
1237 /* Odd number of charge groups, easy */
1238 naaj = 1 + (cgtot/2);
1240 else if ((cgtot % 4) == 0)
1242 /* Multiple of four is hard */
1279 fprintf(log,"naaj=%d\n",naaj);
1285 /************************************************
1287 * S I M P L E C O R E S T U F F
1289 ************************************************/
1291 static real calc_image_tric(rvec xi,rvec xj,matrix box,
1292 rvec b_inv,int *shift)
1294 /* This code assumes that the cut-off is smaller than
1295 * a half times the smallest diagonal element of the box.
1302 /* Compute diff vector */
1303 dz = xj[ZZ] - xi[ZZ];
1304 dy = xj[YY] - xi[YY];
1305 dx = xj[XX] - xi[XX];
1307 /* Perform NINT operation, using trunc operation, therefore
1308 * we first add 2.5 then subtract 2 again
1310 tz = dz*b_inv[ZZ] + h25;
1312 dz -= tz*box[ZZ][ZZ];
1313 dy -= tz*box[ZZ][YY];
1314 dx -= tz*box[ZZ][XX];
1316 ty = dy*b_inv[YY] + h25;
1318 dy -= ty*box[YY][YY];
1319 dx -= ty*box[YY][XX];
1321 tx = dx*b_inv[XX]+h25;
1323 dx -= tx*box[XX][XX];
1325 /* Distance squared */
1326 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1328 *shift = XYZ2IS(tx,ty,tz);
1333 static real calc_image_rect(rvec xi,rvec xj,rvec box_size,
1334 rvec b_inv,int *shift)
1342 /* Compute diff vector */
1343 dx = xj[XX] - xi[XX];
1344 dy = xj[YY] - xi[YY];
1345 dz = xj[ZZ] - xi[ZZ];
1347 /* Perform NINT operation, using trunc operation, therefore
1348 * we first add 1.5 then subtract 1 again
1350 tx = dx*b_inv[XX] + h15;
1351 ty = dy*b_inv[YY] + h15;
1352 tz = dz*b_inv[ZZ] + h15;
1357 /* Correct diff vector for translation */
1358 ddx = tx*box_size[XX] - dx;
1359 ddy = ty*box_size[YY] - dy;
1360 ddz = tz*box_size[ZZ] - dz;
1362 /* Distance squared */
1363 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1365 *shift = XYZ2IS(tx,ty,tz);
1370 static void add_simple(t_ns_buf *nsbuf,int nrj,atom_id cg_j,
1371 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1372 int icg,int jgid,t_block *cgs,t_excl bexcl[],
1373 int shift,t_forcerec *fr,put_in_list_t *put_in_list)
1375 if (nsbuf->nj + nrj > MAX_CG)
1377 put_in_list(bHaveVdW,ngid,md,icg,jgid,nsbuf->ncg,nsbuf->jcg,
1378 cgs->index,bexcl,shift,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
1379 /* Reset buffer contents */
1380 nsbuf->ncg = nsbuf->nj = 0;
1382 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1386 static void ns_inner_tric(rvec x[],int icg,int *i_egp_flags,
1387 int njcg,atom_id jcg[],
1388 matrix box,rvec b_inv,real rcut2,
1389 t_block *cgs,t_ns_buf **ns_buf,
1390 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1391 t_excl bexcl[],t_forcerec *fr,
1392 put_in_list_t *put_in_list)
1396 int *cginfo=fr->cginfo;
1397 atom_id cg_j,*cgindex;
1400 cgindex = cgs->index;
1402 for(j=0; (j<njcg); j++)
1405 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1406 if (calc_image_tric(x[icg],x[cg_j],box,b_inv,&shift) < rcut2)
1408 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1409 if (!(i_egp_flags[jgid] & EGP_EXCL))
1411 add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1412 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1419 static void ns_inner_rect(rvec x[],int icg,int *i_egp_flags,
1420 int njcg,atom_id jcg[],
1421 gmx_bool bBox,rvec box_size,rvec b_inv,real rcut2,
1422 t_block *cgs,t_ns_buf **ns_buf,
1423 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1424 t_excl bexcl[],t_forcerec *fr,
1425 put_in_list_t *put_in_list)
1429 int *cginfo=fr->cginfo;
1430 atom_id cg_j,*cgindex;
1433 cgindex = cgs->index;
1437 for(j=0; (j<njcg); j++)
1440 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1441 if (calc_image_rect(x[icg],x[cg_j],box_size,b_inv,&shift) < rcut2)
1443 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1444 if (!(i_egp_flags[jgid] & EGP_EXCL))
1446 add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1447 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1455 for(j=0; (j<njcg); j++)
1458 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1459 if ((rcut2 == 0) || (distance2(x[icg],x[cg_j]) < rcut2)) {
1460 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1461 if (!(i_egp_flags[jgid] & EGP_EXCL))
1463 add_simple(&ns_buf[jgid][CENTRAL],nrj,cg_j,
1464 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,CENTRAL,fr,
1472 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1474 static int ns_simple_core(t_forcerec *fr,
1475 gmx_localtop_t *top,
1477 matrix box,rvec box_size,
1478 t_excl bexcl[],atom_id *aaj,
1479 int ngid,t_ns_buf **ns_buf,
1480 put_in_list_t *put_in_list,gmx_bool bHaveVdW[])
1484 int nsearch,icg,jcg,igid,i0,nri,nn;
1487 /* atom_id *i_atoms; */
1488 t_block *cgs=&(top->cgs);
1489 t_blocka *excl=&(top->excls);
1492 gmx_bool bBox,bTriclinic;
1495 rlist2 = sqr(fr->rlist);
1497 bBox = (fr->ePBC != epbcNONE);
1500 for(m=0; (m<DIM); m++)
1502 b_inv[m] = divide_err(1.0,box_size[m]);
1504 bTriclinic = TRICLINIC(box);
1511 cginfo = fr->cginfo;
1514 for (icg=fr->cg0; (icg<fr->hcg); icg++)
1517 i0 = cgs->index[icg];
1518 nri = cgs->index[icg+1]-i0;
1519 i_atoms = &(cgs->a[i0]);
1520 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1521 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1523 igid = GET_CGINFO_GID(cginfo[icg]);
1524 i_egp_flags = fr->egp_flags + ngid*igid;
1525 setexcl(cgs->index[icg],cgs->index[icg+1],excl,TRUE,bexcl);
1527 naaj=calc_naaj(icg,cgs->nr);
1530 ns_inner_tric(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1531 box,b_inv,rlist2,cgs,ns_buf,
1532 bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1536 ns_inner_rect(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1537 bBox,box_size,b_inv,rlist2,cgs,ns_buf,
1538 bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1542 for(nn=0; (nn<ngid); nn++)
1544 for(k=0; (k<SHIFTS); k++)
1546 nsbuf = &(ns_buf[nn][k]);
1549 put_in_list(bHaveVdW,ngid,md,icg,nn,nsbuf->ncg,nsbuf->jcg,
1550 cgs->index,bexcl,k,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
1551 nsbuf->ncg=nsbuf->nj=0;
1555 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1556 setexcl(cgs->index[icg],cgs->index[icg+1],excl,FALSE,bexcl);
1558 close_neighbor_lists(fr,FALSE);
1563 /************************************************
1565 * N S 5 G R I D S T U F F
1567 ************************************************/
1569 static inline void get_dx(int Nx,real gridx,real rc2,int xgi,real x,
1570 int *dx0,int *dx1,real *dcx2)
1598 for(i=xgi0; i>=0; i--)
1600 dcx = (i+1)*gridx-x;
1607 for(i=xgi1; i<Nx; i++)
1620 static inline void get_dx_dd(int Nx,real gridx,real rc2,int xgi,real x,
1621 int ncpddc,int shift_min,int shift_max,
1622 int *g0,int *g1,real *dcx2)
1625 int g_min,g_max,shift_home;
1658 g_min = (shift_min == shift_home ? 0 : ncpddc);
1659 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1666 else if (shift_max < 0)
1681 /* Check one grid cell down */
1682 dcx = ((*g0 - 1) + 1)*gridx - x;
1694 /* Check one grid cell up */
1695 dcx = (*g1 + 1)*gridx - x;
1707 #define sqr(x) ((x)*(x))
1708 #define calc_dx2(XI,YI,ZI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
1709 #define calc_cyl_dx2(XI,YI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
1710 /****************************************************
1712 * F A S T N E I G H B O R S E A R C H I N G
1714 * Optimized neighboursearching routine using grid
1715 * at least 1x1x1, see GROMACS manual
1717 ****************************************************/
1720 static void get_cutoff2(t_forcerec *fr,gmx_bool bDoLongRange,
1721 real *rvdw2,real *rcoul2,
1722 real *rs2,real *rm2,real *rl2)
1724 *rs2 = sqr(fr->rlist);
1726 if (bDoLongRange && fr->bTwinRange)
1728 /* The VdW and elec. LR cut-off's could be different,
1729 * so we can not simply set them to rlistlong.
1731 if (EVDW_MIGHT_BE_ZERO_AT_CUTOFF(fr->vdwtype) &&
1732 fr->rvdw > fr->rlist)
1734 *rvdw2 = sqr(fr->rlistlong);
1738 *rvdw2 = sqr(fr->rvdw);
1740 if (EEL_MIGHT_BE_ZERO_AT_CUTOFF(fr->eeltype) &&
1741 fr->rcoulomb > fr->rlist)
1743 *rcoul2 = sqr(fr->rlistlong);
1747 *rcoul2 = sqr(fr->rcoulomb);
1752 /* Workaround for a gcc -O3 or -ffast-math problem */
1756 *rm2 = min(*rvdw2,*rcoul2);
1757 *rl2 = max(*rvdw2,*rcoul2);
1760 static void init_nsgrid_lists(t_forcerec *fr,int ngid,gmx_ns_t *ns)
1762 real rvdw2,rcoul2,rs2,rm2,rl2;
1765 get_cutoff2(fr,TRUE,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
1767 /* Short range buffers */
1768 snew(ns->nl_sr,ngid);
1771 snew(ns->nlr_ljc,ngid);
1772 snew(ns->nlr_one,ngid);
1774 /* Always allocate both list types, since rcoulomb might now change with PME load balancing */
1775 /* Long range VdW and Coul buffers */
1776 snew(ns->nl_lr_ljc,ngid);
1777 /* Long range VdW or Coul only buffers */
1778 snew(ns->nl_lr_one,ngid);
1780 for(j=0; (j<ngid); j++) {
1781 snew(ns->nl_sr[j],MAX_CG);
1782 snew(ns->nl_lr_ljc[j],MAX_CG);
1783 snew(ns->nl_lr_one[j],MAX_CG);
1788 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
1793 static int nsgrid_core(FILE *log,t_commrec *cr,t_forcerec *fr,
1794 matrix box,rvec box_size,int ngid,
1795 gmx_localtop_t *top,
1796 t_grid *grid,rvec x[],
1797 t_excl bexcl[],gmx_bool *bExcludeAlleg,
1798 t_nrnb *nrnb,t_mdatoms *md,
1799 real *lambda,real *dvdlambda,
1800 gmx_grppairener_t *grppener,
1801 put_in_list_t *put_in_list,
1802 gmx_bool bHaveVdW[],
1803 gmx_bool bDoLongRange,gmx_bool bMakeQMMMnblist)
1806 atom_id **nl_lr_ljc,**nl_lr_one,**nl_sr;
1807 int *nlr_ljc,*nlr_one,*nsr;
1808 gmx_domdec_t *dd=NULL;
1809 t_block *cgs=&(top->cgs);
1810 int *cginfo=fr->cginfo;
1811 /* atom_id *i_atoms,*cgsindex=cgs->index; */
1813 int cell_x,cell_y,cell_z;
1814 int d,tx,ty,tz,dx,dy,dz,cj;
1815 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
1816 int zsh_ty,zsh_tx,ysh_tx;
1818 int dx0,dx1,dy0,dy1,dz0,dz1;
1819 int Nx,Ny,Nz,shift=-1,j,nrj,nns,nn=-1;
1820 real gridx,gridy,gridz,grid_x,grid_y,grid_z;
1821 real *dcx2,*dcy2,*dcz2;
1823 int cg0,cg1,icg=-1,cgsnr,i0,igid,nri,naaj,max_jcg;
1824 int jcg0,jcg1,jjcg,cgj0,jgid;
1825 int *grida,*gridnra,*gridind;
1826 gmx_bool rvdw_lt_rcoul,rcoul_lt_rvdw;
1827 rvec xi,*cgcm,grid_offset;
1828 real r2,rs2,rvdw2,rcoul2,rm2,rl2,XI,YI,ZI,dcx,dcy,dcz,tmp1,tmp2;
1830 gmx_bool bDomDec,bTriclinicX,bTriclinicY;
1835 bDomDec = DOMAINDECOMP(cr);
1841 bTriclinicX = ((YY < grid->npbcdim &&
1842 (!bDomDec || dd->nc[YY]==1) && box[YY][XX] != 0) ||
1843 (ZZ < grid->npbcdim &&
1844 (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][XX] != 0));
1845 bTriclinicY = (ZZ < grid->npbcdim &&
1846 (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][YY] != 0);
1850 get_cutoff2(fr,bDoLongRange,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
1852 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
1853 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
1855 if (bMakeQMMMnblist)
1863 nl_lr_ljc = ns->nl_lr_ljc;
1864 nl_lr_one = ns->nl_lr_one;
1865 nlr_ljc = ns->nlr_ljc;
1866 nlr_one = ns->nlr_one;
1874 gridind = grid->index;
1875 gridnra = grid->nra;
1878 gridx = grid->cell_size[XX];
1879 gridy = grid->cell_size[YY];
1880 gridz = grid->cell_size[ZZ];
1884 copy_rvec(grid->cell_offset,grid_offset);
1885 copy_ivec(grid->ncpddc,ncpddc);
1890 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
1891 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
1892 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
1893 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
1894 if (zsh_tx!=0 && ysh_tx!=0)
1896 /* This could happen due to rounding, when both ratios are 0.5 */
1905 /* We only want a list for the test particle */
1914 /* Set the shift range */
1915 for(d=0; d<DIM; d++)
1919 /* Check if we need periodicity shifts.
1920 * Without PBC or with domain decomposition we don't need them.
1922 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
1929 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
1940 /* Loop over charge groups */
1941 for(icg=cg0; (icg < cg1); icg++)
1943 igid = GET_CGINFO_GID(cginfo[icg]);
1944 /* Skip this charge group if all energy groups are excluded! */
1945 if (bExcludeAlleg[igid])
1950 i0 = cgs->index[icg];
1952 if (bMakeQMMMnblist)
1954 /* Skip this charge group if it is not a QM atom while making a
1955 * QM/MM neighbourlist
1957 if (md->bQM[i0]==FALSE)
1959 continue; /* MM particle, go to next particle */
1962 /* Compute the number of charge groups that fall within the control
1965 naaj = calc_naaj(icg,cgsnr);
1972 /* make a normal neighbourlist */
1976 /* Get the j charge-group and dd cell shift ranges */
1977 dd_get_ns_ranges(cr->dd,icg,&jcg0,&jcg1,sh0,sh1);
1982 /* Compute the number of charge groups that fall within the control
1985 naaj = calc_naaj(icg,cgsnr);
1991 /* The i-particle is awlways the test particle,
1992 * so we want all j-particles
1994 max_jcg = cgsnr - 1;
1998 max_jcg = jcg1 - cgsnr;
2003 i_egp_flags = fr->egp_flags + igid*ngid;
2005 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2006 setexcl(i0,cgs->index[icg+1],&top->excls,TRUE,bexcl);
2008 ci2xyz(grid,icg,&cell_x,&cell_y,&cell_z);
2010 /* Changed iicg to icg, DvdS 990115
2011 * (but see consistency check above, DvdS 990330)
2014 fprintf(log,"icg=%5d, naaj=%5d, cell %d %d %d\n",
2015 icg,naaj,cell_x,cell_y,cell_z);
2017 /* Loop over shift vectors in three dimensions */
2018 for (tz=-shp[ZZ]; tz<=shp[ZZ]; tz++)
2020 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2021 /* Calculate range of cells in Z direction that have the shift tz */
2022 zgi = cell_z + tz*Nz;
2025 get_dx(Nz,gridz,rl2,zgi,ZI,&dz0,&dz1,dcz2);
2027 get_dx_dd(Nz,gridz,rl2,zgi,ZI-grid_offset[ZZ],
2028 ncpddc[ZZ],sh0[ZZ],sh1[ZZ],&dz0,&dz1,dcz2);
2034 for (ty=-shp[YY]; ty<=shp[YY]; ty++)
2036 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2037 /* Calculate range of cells in Y direction that have the shift ty */
2040 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2044 ygi = cell_y + ty*Ny;
2047 get_dx(Ny,gridy,rl2,ygi,YI,&dy0,&dy1,dcy2);
2049 get_dx_dd(Ny,gridy,rl2,ygi,YI-grid_offset[YY],
2050 ncpddc[YY],sh0[YY],sh1[YY],&dy0,&dy1,dcy2);
2056 for (tx=-shp[XX]; tx<=shp[XX]; tx++)
2058 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2059 /* Calculate range of cells in X direction that have the shift tx */
2062 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2066 xgi = cell_x + tx*Nx;
2069 get_dx(Nx,gridx,rl2,xgi*Nx,XI,&dx0,&dx1,dcx2);
2071 get_dx_dd(Nx,gridx,rl2,xgi,XI-grid_offset[XX],
2072 ncpddc[XX],sh0[XX],sh1[XX],&dx0,&dx1,dcx2);
2078 /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2079 * from the neigbour list as it will not interact */
2080 if (fr->adress_type != eAdressOff){
2081 if (md->wf[cgs->index[icg]]==0 && egp_explicit(fr, igid)){
2085 /* Get shift vector */
2086 shift=XYZ2IS(tx,ty,tz);
2088 range_check(shift,0,SHIFTS);
2090 for(nn=0; (nn<ngid); nn++)
2097 fprintf(log,"shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2098 shift,dx0,dx1,dy0,dy1,dz0,dz1);
2099 fprintf(log,"cgcm: %8.3f %8.3f %8.3f\n",cgcm[icg][XX],
2100 cgcm[icg][YY],cgcm[icg][ZZ]);
2101 fprintf(log,"xi: %8.3f %8.3f %8.3f\n",XI,YI,ZI);
2103 for (dx=dx0; (dx<=dx1); dx++)
2105 tmp1 = rl2 - dcx2[dx];
2106 for (dy=dy0; (dy<=dy1); dy++)
2108 tmp2 = tmp1 - dcy2[dy];
2111 for (dz=dz0; (dz<=dz1); dz++) {
2112 if (tmp2 > dcz2[dz]) {
2113 /* Find grid-cell cj in which possible neighbours are */
2114 cj = xyz2ci(Ny,Nz,dx,dy,dz);
2116 /* Check out how many cgs (nrj) there in this cell */
2119 /* Find the offset in the cg list */
2122 /* Check if all j's are out of range so we
2123 * can skip the whole cell.
2124 * Should save some time, especially with DD.
2127 (grida[cgj0] >= max_jcg &&
2128 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2134 for (j=0; (j<nrj); j++)
2136 jjcg = grida[cgj0+j];
2138 /* check whether this guy is in range! */
2139 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2142 r2=calc_dx2(XI,YI,ZI,cgcm[jjcg]);
2144 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2145 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2146 /* check energy group exclusions */
2147 if (!(i_egp_flags[jgid] & EGP_EXCL))
2151 if (nsr[jgid] >= MAX_CG)
2153 /* Add to short-range list */
2154 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2155 nsr[jgid],nl_sr[jgid],
2156 cgs->index,/* cgsatoms, */ bexcl,
2157 shift,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
2160 nl_sr[jgid][nsr[jgid]++]=jjcg;
2164 if (nlr_ljc[jgid] >= MAX_CG)
2166 /* Add to LJ+coulomb long-range list */
2167 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2168 nlr_ljc[jgid],nl_lr_ljc[jgid],top->cgs.index,
2169 bexcl,shift,fr,TRUE,TRUE,TRUE,fr->solvent_opt);
2172 nl_lr_ljc[jgid][nlr_ljc[jgid]++]=jjcg;
2176 if (nlr_one[jgid] >= MAX_CG)
2178 /* Add to long-range list with only coul, or only LJ */
2179 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2180 nlr_one[jgid],nl_lr_one[jgid],top->cgs.index,
2181 bexcl,shift,fr,TRUE,rvdw_lt_rcoul,rcoul_lt_rvdw,fr->solvent_opt);
2184 nl_lr_one[jgid][nlr_one[jgid]++]=jjcg;
2196 /* CHECK whether there is anything left in the buffers */
2197 for(nn=0; (nn<ngid); nn++)
2201 put_in_list(bHaveVdW,ngid,md,icg,nn,nsr[nn],nl_sr[nn],
2202 cgs->index, /* cgsatoms, */ bexcl,
2203 shift,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
2206 if (nlr_ljc[nn] > 0)
2208 put_in_list(bHaveVdW,ngid,md,icg,nn,nlr_ljc[nn],
2209 nl_lr_ljc[nn],top->cgs.index,
2210 bexcl,shift,fr,TRUE,TRUE,TRUE,fr->solvent_opt);
2213 if (nlr_one[nn] > 0)
2215 put_in_list(bHaveVdW,ngid,md,icg,nn,nlr_one[nn],
2216 nl_lr_one[nn],top->cgs.index,
2217 bexcl,shift,fr,TRUE,rvdw_lt_rcoul,rcoul_lt_rvdw,fr->solvent_opt);
2223 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2224 setexcl(cgs->index[icg],cgs->index[icg+1],&top->excls,FALSE,bexcl);
2226 /* No need to perform any left-over force calculations anymore (as we used to do here)
2227 * since we now save the proper long-range lists for later evaluation.
2232 /* Close neighbourlists */
2233 close_neighbor_lists(fr,bMakeQMMMnblist);
2238 void ns_realloc_natoms(gmx_ns_t *ns,int natoms)
2242 if (natoms > ns->nra_alloc)
2244 ns->nra_alloc = over_alloc_dd(natoms);
2245 srenew(ns->bexcl,ns->nra_alloc);
2246 for(i=0; i<ns->nra_alloc; i++)
2253 void init_ns(FILE *fplog,const t_commrec *cr,
2254 gmx_ns_t *ns,t_forcerec *fr,
2255 const gmx_mtop_t *mtop,
2258 int mt,icg,nr_in_cg,maxcg,i,j,jcg,ngid,ncg;
2262 /* Compute largest charge groups size (# atoms) */
2264 for(mt=0; mt<mtop->nmoltype; mt++) {
2265 cgs = &mtop->moltype[mt].cgs;
2266 for (icg=0; (icg < cgs->nr); icg++)
2268 nr_in_cg=max(nr_in_cg,(int)(cgs->index[icg+1]-cgs->index[icg]));
2272 /* Verify whether largest charge group is <= max cg.
2273 * This is determined by the type of the local exclusion type
2274 * Exclusions are stored in bits. (If the type is not large
2275 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2277 maxcg = sizeof(t_excl)*8;
2278 if (nr_in_cg > maxcg)
2280 gmx_fatal(FARGS,"Max #atoms in a charge group: %d > %d\n",
2284 ngid = mtop->groups.grps[egcENER].nr;
2285 snew(ns->bExcludeAlleg,ngid);
2286 for(i=0; i<ngid; i++) {
2287 ns->bExcludeAlleg[i] = TRUE;
2288 for(j=0; j<ngid; j++)
2290 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2292 ns->bExcludeAlleg[i] = FALSE;
2299 ns->grid = init_grid(fplog,fr);
2300 init_nsgrid_lists(fr,ngid,ns);
2305 snew(ns->ns_buf,ngid);
2306 for(i=0; (i<ngid); i++)
2308 snew(ns->ns_buf[i],SHIFTS);
2310 ncg = ncg_mtop(mtop);
2311 snew(ns->simple_aaj,2*ncg);
2312 for(jcg=0; (jcg<ncg); jcg++)
2314 ns->simple_aaj[jcg] = jcg;
2315 ns->simple_aaj[jcg+ncg] = jcg;
2319 /* Create array that determines whether or not atoms have VdW */
2320 snew(ns->bHaveVdW,fr->ntype);
2321 for(i=0; (i<fr->ntype); i++)
2323 for(j=0; (j<fr->ntype); j++)
2325 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2327 ((BHAMA(fr->nbfp,fr->ntype,i,j) != 0) ||
2328 (BHAMB(fr->nbfp,fr->ntype,i,j) != 0) ||
2329 (BHAMC(fr->nbfp,fr->ntype,i,j) != 0)) :
2330 ((C6(fr->nbfp,fr->ntype,i,j) != 0) ||
2331 (C12(fr->nbfp,fr->ntype,i,j) != 0))));
2335 pr_bvec(debug,0,"bHaveVdW",ns->bHaveVdW,fr->ntype,TRUE);
2339 if (!DOMAINDECOMP(cr))
2341 /* This could be reduced with particle decomposition */
2342 ns_realloc_natoms(ns,mtop->natoms);
2345 ns->nblist_initialized=FALSE;
2347 /* nbr list debug dump */
2349 char *ptr=getenv("GMX_DUMP_NL");
2352 ns->dump_nl=strtol(ptr,NULL,10);
2355 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2366 int search_neighbours(FILE *log,t_forcerec *fr,
2367 rvec x[],matrix box,
2368 gmx_localtop_t *top,
2369 gmx_groups_t *groups,
2371 t_nrnb *nrnb,t_mdatoms *md,
2372 real *lambda,real *dvdlambda,
2373 gmx_grppairener_t *grppener,
2375 gmx_bool bDoLongRangeNS)
2377 t_block *cgs=&(top->cgs);
2378 rvec box_size,grid_x0,grid_x1;
2380 real min_size,grid_dens;
2384 gmx_bool *i_egp_flags;
2385 int cg_start,cg_end,start,end;
2388 gmx_domdec_zones_t *dd_zones;
2389 put_in_list_t *put_in_list;
2393 /* Set some local variables */
2395 ngid = groups->grps[egcENER].nr;
2397 for(m=0; (m<DIM); m++)
2399 box_size[m] = box[m][m];
2402 if (fr->ePBC != epbcNONE)
2404 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC,box))
2406 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.");
2410 min_size = min(box_size[XX],min(box_size[YY],box_size[ZZ]));
2411 if (2*fr->rlistlong >= min_size)
2412 gmx_fatal(FARGS,"One of the box diagonal elements has become smaller than twice the cut-off length.");
2416 if (DOMAINDECOMP(cr))
2418 ns_realloc_natoms(ns,cgs->index[cgs->nr]);
2422 /* Reset the neighbourlists */
2423 reset_neighbor_lists(fr,TRUE,TRUE);
2425 if (bGrid && bFillGrid)
2429 if (DOMAINDECOMP(cr))
2431 dd_zones = domdec_zones(cr->dd);
2437 get_nsgrid_boundaries(grid->nboundeddim,box,NULL,NULL,NULL,NULL,
2438 cgs->nr,fr->cg_cm,grid_x0,grid_x1,&grid_dens);
2440 grid_first(log,grid,NULL,NULL,fr->ePBC,box,grid_x0,grid_x1,
2441 fr->rlistlong,grid_dens);
2445 /* Don't know why this all is... (DvdS 3/99) */
2451 end = (cgs->nr+1)/2;
2454 if (DOMAINDECOMP(cr))
2457 fill_grid(log,dd_zones,grid,end,-1,end,fr->cg_cm);
2459 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2463 fill_grid(log,NULL,grid,cgs->nr,fr->cg0,fr->hcg,fr->cg_cm);
2464 grid->icg0 = fr->cg0;
2465 grid->icg1 = fr->hcg;
2473 calc_elemnr(log,grid,start,end,cgs->nr);
2475 grid_last(log,grid,start,end,cgs->nr);
2479 check_grid(debug,grid);
2480 print_grid(debug,grid);
2485 /* Set the grid cell index for the test particle only.
2486 * The cell to cg index is not corrected, but that does not matter.
2488 fill_grid(log,NULL,ns->grid,fr->hcg,fr->hcg-1,fr->hcg,fr->cg_cm);
2492 if (!fr->ns.bCGlist)
2494 put_in_list = put_in_list_at;
2498 put_in_list = put_in_list_cg;
2505 nsearch = nsgrid_core(log,cr,fr,box,box_size,ngid,top,
2506 grid,x,ns->bexcl,ns->bExcludeAlleg,
2507 nrnb,md,lambda,dvdlambda,grppener,
2508 put_in_list,ns->bHaveVdW,
2509 bDoLongRangeNS,FALSE);
2511 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2512 * the classical calculation. The charge-charge interaction
2513 * between QM and MM atoms is handled in the QMMM core calculation
2514 * (see QMMM.c). The VDW however, we'd like to compute classically
2515 * and the QM MM atom pairs have just been put in the
2516 * corresponding neighbourlists. in case of QMMM we still need to
2517 * fill a special QMMM neighbourlist that contains all neighbours
2518 * of the QM atoms. If bQMMM is true, this list will now be made:
2520 if (fr->bQMMM && fr->qr->QMMMscheme!=eQMMMschemeoniom)
2522 nsearch += nsgrid_core(log,cr,fr,box,box_size,ngid,top,
2523 grid,x,ns->bexcl,ns->bExcludeAlleg,
2524 nrnb,md,lambda,dvdlambda,grppener,
2525 put_in_list_qmmm,ns->bHaveVdW,
2526 bDoLongRangeNS,TRUE);
2531 nsearch = ns_simple_core(fr,top,md,box,box_size,
2532 ns->bexcl,ns->simple_aaj,
2533 ngid,ns->ns_buf,put_in_list,ns->bHaveVdW);
2541 inc_nrnb(nrnb,eNR_NS,nsearch);
2542 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2547 int natoms_beyond_ns_buffer(t_inputrec *ir,t_forcerec *fr,t_block *cgs,
2548 matrix scale_tot,rvec *x)
2550 int cg0,cg1,cg,a0,a1,a,i,j;
2551 real rint,hbuf2,scale;
2553 gmx_bool bIsotropic;
2558 rint = max(ir->rcoulomb,ir->rvdw);
2559 if (ir->rlist < rint)
2561 gmx_fatal(FARGS,"The neighbor search buffer has negative size: %f nm",
2569 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2571 hbuf2 = sqr(0.5*(ir->rlist - rint));
2572 for(cg=cg0; cg<cg1; cg++)
2574 a0 = cgs->index[cg];
2575 a1 = cgs->index[cg+1];
2576 for(a=a0; a<a1; a++)
2578 if (distance2(cg_cm[cg],x[a]) > hbuf2)
2588 scale = scale_tot[0][0];
2589 for(i=1; i<DIM; i++)
2591 /* With anisotropic scaling, the original spherical ns volumes become
2592 * ellipsoids. To avoid costly transformations we use the minimum
2593 * eigenvalue of the scaling matrix for determining the buffer size.
2594 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2596 scale = min(scale,scale_tot[i][i]);
2597 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2603 if (scale_tot[i][j] != 0)
2609 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2612 for(cg=cg0; cg<cg1; cg++)
2614 svmul(scale,cg_cm[cg],cgsc);
2615 a0 = cgs->index[cg];
2616 a1 = cgs->index[cg+1];
2617 for(a=a0; a<a1; a++)
2619 if (distance2(cgsc,x[a]) > hbuf2)
2628 /* Anistropic scaling */
2629 for(cg=cg0; cg<cg1; cg++)
2631 /* Since scale_tot contains the transpose of the scaling matrix,
2632 * we need to multiply with the transpose.
2634 tmvmul_ur0(scale_tot,cg_cm[cg],cgsc);
2635 a0 = cgs->index[cg];
2636 a1 = cgs->index[cg+1];
2637 for(a=a0; a<a1; a++)
2639 if (distance2(cgsc,x[a]) > hbuf2)