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))
83 /************************************************
85 * U T I L I T I E S F O R N S
87 ************************************************/
89 static void reallocate_nblist(t_nblist *nl)
93 fprintf(debug,"reallocating neigborlist (ielec=%d, ivdw=%d, igeometry=%d, free_energy=%d), maxnri=%d\n",
94 nl->ielec,nl->ivdw,nl->igeometry,nl->free_energy,nl->maxnri);
96 srenew(nl->iinr, nl->maxnri);
97 if (nl->igeometry == GMX_NBLIST_GEOMETRY_CG_CG)
99 srenew(nl->iinr_end,nl->maxnri);
101 srenew(nl->gid, nl->maxnri);
102 srenew(nl->shift, nl->maxnri);
103 srenew(nl->jindex, nl->maxnri+1);
107 static void init_nblist(FILE *log, t_nblist *nl_sr,t_nblist *nl_lr,
109 int ivdw, int ivdwmod,
110 int ielec, int ielecmod,
111 gmx_bool bfree, int igeometry)
119 nl = (i == 0) ? nl_sr : nl_lr;
120 homenr = (i == 0) ? maxsr : maxlr;
128 /* Set coul/vdw in neighborlist, and for the normal loops we determine
129 * an index of which one to call.
132 nl->ivdwmod = ivdwmod;
134 nl->ielecmod = ielecmod;
135 nl->free_energy = bfree;
136 nl->igeometry = igeometry;
140 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
143 gmx_nonbonded_set_kernel_pointers( (i==0) ? log : NULL,nl);
145 /* maxnri is influenced by the number of shifts (maximum is 8)
146 * and the number of energy groups.
147 * If it is not enough, nl memory will be reallocated during the run.
148 * 4 seems to be a reasonable factor, which only causes reallocation
149 * during runs with tiny and many energygroups.
151 nl->maxnri = homenr*4;
160 reallocate_nblist(nl);
165 fprintf(debug,"Initiating neighbourlist (ielec=%d, ivdw=%d, free=%d) for %s interactions,\nwith %d SR, %d LR atoms.\n",
166 nl->ielec,nl->ivdw,nl->free_energy,gmx_nblist_geometry_names[nl->igeometry],maxsr,maxlr);
171 void init_neighbor_list(FILE *log,t_forcerec *fr,int homenr)
173 /* Make maxlr tunable! (does not seem to be a big difference though)
174 * This parameter determines the number of i particles in a long range
175 * neighbourlist. Too few means many function calls, too many means
178 int maxsr,maxsr_wat,maxlr,maxlr_wat;
179 int ielec,ielecf,ivdw,ielecmod,ielecmodf,ivdwmod;
181 int igeometry_def,igeometry_w,igeometry_ww;
185 /* maxsr = homenr-fr->nWatMol*3; */
190 gmx_fatal(FARGS,"%s, %d: Negative number of short range atoms.\n"
191 "Call your Gromacs dealer for assistance.",__FILE__,__LINE__);
193 /* This is just for initial allocation, so we do not reallocate
194 * all the nlist arrays many times in a row.
195 * The numbers seem very accurate, but they are uncritical.
197 maxsr_wat = min(fr->nWatMol,(homenr+2)/3);
201 maxlr_wat = min(maxsr_wat,maxlr);
205 maxlr = maxlr_wat = 0;
208 /* Determine the values for ielec/ivdw. */
209 ielec = fr->nbkernel_elec_interaction;
210 ivdw = fr->nbkernel_vdw_interaction;
211 ielecmod = fr->nbkernel_elec_modifier;
212 ivdwmod = fr->nbkernel_vdw_modifier;
214 fr->ns.bCGlist = (getenv("GMX_NBLISTCG") != 0);
217 igeometry_def = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
221 igeometry_def = GMX_NBLIST_GEOMETRY_CG_CG;
224 fprintf(log,"\nUsing charge-group - charge-group neighbor lists and kernels\n\n");
228 if (fr->solvent_opt == esolTIP4P) {
229 igeometry_w = GMX_NBLIST_GEOMETRY_WATER4_PARTICLE;
230 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER4_WATER4;
232 igeometry_w = GMX_NBLIST_GEOMETRY_WATER3_PARTICLE;
233 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER3_WATER3;
236 for(i=0; i<fr->nnblists; i++)
238 nbl = &(fr->nblists[i]);
239 init_nblist(log,&nbl->nlist_sr[eNL_VDWQQ],&nbl->nlist_lr[eNL_VDWQQ],
240 maxsr,maxlr,ivdw,ivdwmod,ielec,ielecmod,FALSE,igeometry_def);
241 init_nblist(log,&nbl->nlist_sr[eNL_VDW],&nbl->nlist_lr[eNL_VDW],
242 maxsr,maxlr,ivdw,ivdwmod,GMX_NBKERNEL_ELEC_NONE,eintmodNONE,FALSE,igeometry_def);
243 init_nblist(log,&nbl->nlist_sr[eNL_QQ],&nbl->nlist_lr[eNL_QQ],
244 maxsr,maxlr,GMX_NBKERNEL_VDW_NONE,eintmodNONE,ielec,ielecmod,FALSE,igeometry_def);
245 init_nblist(log,&nbl->nlist_sr[eNL_VDWQQ_WATER],&nbl->nlist_lr[eNL_VDWQQ_WATER],
246 maxsr_wat,maxlr_wat,ivdw,ivdwmod,ielec,ielecmod, FALSE,igeometry_w);
247 init_nblist(log,&nbl->nlist_sr[eNL_QQ_WATER],&nbl->nlist_lr[eNL_QQ_WATER],
248 maxsr_wat,maxlr_wat,GMX_NBKERNEL_VDW_NONE,eintmodNONE,ielec,ielecmod, FALSE,igeometry_w);
249 init_nblist(log,&nbl->nlist_sr[eNL_VDWQQ_WATERWATER],&nbl->nlist_lr[eNL_VDWQQ_WATERWATER],
250 maxsr_wat,maxlr_wat,ivdw,ivdwmod,ielec,ielecmod, FALSE,igeometry_ww);
251 init_nblist(log,&nbl->nlist_sr[eNL_QQ_WATERWATER],&nbl->nlist_lr[eNL_QQ_WATERWATER],
252 maxsr_wat,maxlr_wat,GMX_NBKERNEL_VDW_NONE,eintmodNONE,ielec,ielecmod, FALSE,igeometry_ww);
254 /* Did we get the solvent loops so we can use optimized water kernels? */
255 if(nbl->nlist_sr[eNL_VDWQQ_WATER].kernelptr_vf==NULL
256 || nbl->nlist_sr[eNL_QQ_WATER].kernelptr_vf==NULL
257 #ifndef DISABLE_WATERWATER_NLIST
258 || nbl->nlist_sr[eNL_VDWQQ_WATERWATER].kernelptr_vf==NULL
259 || nbl->nlist_sr[eNL_QQ_WATERWATER].kernelptr_vf==NULL
263 fr->solvent_opt = esolNO;
264 fprintf(log,"Note: The available nonbonded kernels do not support water optimization - disabling.\n");
267 if (fr->efep != efepNO)
269 if ((fr->bEwald) && (fr->sc_alphacoul > 0)) /* need to handle long range differently if using softcore */
271 ielecf = GMX_NBKERNEL_ELEC_EWALD;
272 ielecmodf = eintmodNONE;
277 ielecmodf = ielecmod;
280 init_nblist(log,&nbl->nlist_sr[eNL_VDWQQ_FREE],&nbl->nlist_lr[eNL_VDWQQ_FREE],
281 maxsr,maxlr,ivdw,ivdwmod,ielecf,ielecmod,TRUE,GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE);
282 init_nblist(log,&nbl->nlist_sr[eNL_VDW_FREE],&nbl->nlist_lr[eNL_VDW_FREE],
283 maxsr,maxlr,ivdw,ivdwmod,GMX_NBKERNEL_ELEC_NONE,eintmodNONE,TRUE,GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE);
284 init_nblist(log,&nbl->nlist_sr[eNL_QQ_FREE],&nbl->nlist_lr[eNL_QQ_FREE],
285 maxsr,maxlr,GMX_NBKERNEL_VDW_NONE,eintmodNONE,ielecf,ielecmod,TRUE,GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE);
289 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
291 init_nblist(log,&fr->QMMMlist,NULL,
292 maxsr,maxlr,0,0,ielec,ielecmod,FALSE,GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE);
300 fr->ns.nblist_initialized=TRUE;
303 static void reset_nblist(t_nblist *nl)
314 static void reset_neighbor_lists(t_forcerec *fr,gmx_bool bResetSR, gmx_bool bResetLR)
320 /* only reset the short-range nblist */
321 reset_nblist(&(fr->QMMMlist));
324 for(n=0; n<fr->nnblists; n++)
326 for(i=0; i<eNL_NR; i++)
330 reset_nblist( &(fr->nblists[n].nlist_sr[i]) );
334 reset_nblist( &(fr->nblists[n].nlist_lr[i]) );
343 static inline void new_i_nblist(t_nblist *nlist,
344 gmx_bool bLR,atom_id i_atom,int shift,int gid)
350 /* Check whether we have to increase the i counter */
352 (nlist->iinr[nri] != i_atom) ||
353 (nlist->shift[nri] != shift) ||
354 (nlist->gid[nri] != gid))
356 /* This is something else. Now see if any entries have
357 * been added in the list of the previous atom.
360 ((nlist->jindex[nri+1] > nlist->jindex[nri]) &&
361 (nlist->gid[nri] != -1)))
363 /* If so increase the counter */
366 if (nlist->nri >= nlist->maxnri)
368 nlist->maxnri += over_alloc_large(nlist->nri);
369 reallocate_nblist(nlist);
372 /* Set the number of neighbours and the atom number */
373 nlist->jindex[nri+1] = nlist->jindex[nri];
374 nlist->iinr[nri] = i_atom;
375 nlist->gid[nri] = gid;
376 nlist->shift[nri] = shift;
380 static inline void close_i_nblist(t_nblist *nlist)
382 int nri = nlist->nri;
387 nlist->jindex[nri+1] = nlist->nrj;
389 len=nlist->nrj - nlist->jindex[nri];
391 /* nlist length for water i molecules is treated statically
394 if (len > nlist->maxlen)
401 static inline void close_nblist(t_nblist *nlist)
403 /* Only close this nblist when it has been initialized.
404 * Avoid the creation of i-lists with no j-particles.
408 /* Some assembly kernels do not support empty lists,
409 * make sure here that we don't generate any empty lists.
410 * With the current ns code this branch is taken in two cases:
411 * No i-particles at all: nri=-1 here
412 * There are i-particles, but no j-particles; nri=0 here
418 /* Close list number nri by incrementing the count */
423 static inline void close_neighbor_lists(t_forcerec *fr,gmx_bool bMakeQMMMnblist)
429 close_nblist(&(fr->QMMMlist));
432 for(n=0; n<fr->nnblists; n++)
434 for(i=0; (i<eNL_NR); i++)
436 close_nblist(&(fr->nblists[n].nlist_sr[i]));
437 close_nblist(&(fr->nblists[n].nlist_lr[i]));
443 static inline void add_j_to_nblist(t_nblist *nlist,atom_id j_atom,gmx_bool bLR)
447 if (nlist->nrj >= nlist->maxnrj)
449 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
451 fprintf(debug,"Increasing %s nblist (ielec=%d,ivdw=%d,free=%d,igeometry=%d) j size to %d\n",
452 bLR ? "LR" : "SR",nlist->ielec,nlist->ivdw,nlist->free_energy,nlist->igeometry,nlist->maxnrj);
454 srenew(nlist->jjnr,nlist->maxnrj);
457 nlist->jjnr[nrj] = j_atom;
461 static inline void add_j_to_nblist_cg(t_nblist *nlist,
462 atom_id j_start,int j_end,
463 t_excl *bexcl,gmx_bool i_is_j,
469 if (nlist->nrj >= nlist->maxnrj)
471 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
473 fprintf(debug,"Increasing %s nblist (ielec=%d,ivdw=%d,free=%d,igeometry=%d) j size to %d\n",
474 bLR ? "LR" : "SR",nlist->ielec,nlist->ivdw,nlist->free_energy,nlist->igeometry,nlist->maxnrj);
476 srenew(nlist->jjnr ,nlist->maxnrj);
477 srenew(nlist->jjnr_end,nlist->maxnrj);
478 srenew(nlist->excl ,nlist->maxnrj*MAX_CGCGSIZE);
481 nlist->jjnr[nrj] = j_start;
482 nlist->jjnr_end[nrj] = j_end;
484 if (j_end - j_start > MAX_CGCGSIZE)
486 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);
489 /* Set the exclusions */
490 for(j=j_start; j<j_end; j++)
492 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
496 /* Avoid double counting of intra-cg interactions */
497 for(j=1; j<j_end-j_start; j++)
499 nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
507 put_in_list_t(gmx_bool bHaveVdW[],
524 put_in_list_at(gmx_bool bHaveVdW[],
540 /* The a[] index has been removed,
541 * to put it back in i_atom should be a[i0] and jj should be a[jj].
546 t_nblist * vdwc_free = NULL;
547 t_nblist * vdw_free = NULL;
548 t_nblist * coul_free = NULL;
549 t_nblist * vdwc_ww = NULL;
550 t_nblist * coul_ww = NULL;
552 int i,j,jcg,igid,gid,nbl_ind,ind_ij;
553 atom_id jj,jj0,jj1,i_atom;
558 real *charge,*chargeB;
560 gmx_bool bFreeEnergy,bFree,bFreeJ,bNotEx,*bPert;
561 gmx_bool bDoVdW_i,bDoCoul_i,bDoCoul_i_sol;
565 /* Copy some pointers */
567 charge = md->chargeA;
568 chargeB = md->chargeB;
571 bPert = md->bPerturbed;
575 nicg = index[icg+1]-i0;
577 /* Get the i charge group info */
578 igid = GET_CGINFO_GID(cginfo[icg]);
580 iwater = (solvent_opt!=esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
585 /* Check if any of the particles involved are perturbed.
586 * If not we can do the cheaper normal put_in_list
587 * and use more solvent optimization.
589 for(i=0; i<nicg; i++)
591 bFreeEnergy |= bPert[i0+i];
593 /* Loop over the j charge groups */
594 for(j=0; (j<nj && !bFreeEnergy); j++)
599 /* Finally loop over the atoms in the j-charge group */
600 for(jj=jj0; jj<jj1; jj++)
602 bFreeEnergy |= bPert[jj];
607 /* Unpack pointers to neighbourlist structs */
608 if (fr->nnblists == 1)
614 nbl_ind = fr->gid2nblists[GID(igid,jgid,ngid)];
618 nlist = fr->nblists[nbl_ind].nlist_lr;
622 nlist = fr->nblists[nbl_ind].nlist_sr;
625 if (iwater != esolNO)
627 vdwc = &nlist[eNL_VDWQQ_WATER];
628 vdw = &nlist[eNL_VDW];
629 coul = &nlist[eNL_QQ_WATER];
630 #ifndef DISABLE_WATERWATER_NLIST
631 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
632 coul_ww = &nlist[eNL_QQ_WATERWATER];
637 vdwc = &nlist[eNL_VDWQQ];
638 vdw = &nlist[eNL_VDW];
639 coul = &nlist[eNL_QQ];
644 if (iwater != esolNO)
646 /* Loop over the atoms in the i charge group */
648 gid = GID(igid,jgid,ngid);
649 /* Create new i_atom for each energy group */
650 if (bDoCoul && bDoVdW)
652 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
653 #ifndef DISABLE_WATERWATER_NLIST
654 new_i_nblist(vdwc_ww,bLR,i_atom,shift,gid);
659 new_i_nblist(vdw,bLR,i_atom,shift,gid);
663 new_i_nblist(coul,bLR,i_atom,shift,gid);
664 #ifndef DISABLE_WATERWATER_NLIST
665 new_i_nblist(coul_ww,bLR,i_atom,shift,gid);
668 /* Loop over the j charge groups */
669 for(j=0; (j<nj); j++)
679 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
681 if (iwater == esolSPC && jwater == esolSPC)
683 /* Interaction between two SPC molecules */
686 /* VdW only - only first atoms in each water interact */
687 add_j_to_nblist(vdw,jj0,bLR);
691 #ifdef DISABLE_WATERWATER_NLIST
692 /* Add entries for the three atoms - only do VdW if we need to */
695 add_j_to_nblist(coul,jj0,bLR);
699 add_j_to_nblist(vdwc,jj0,bLR);
701 add_j_to_nblist(coul,jj0+1,bLR);
702 add_j_to_nblist(coul,jj0+2,bLR);
704 /* One entry for the entire water-water interaction */
707 add_j_to_nblist(coul_ww,jj0,bLR);
711 add_j_to_nblist(vdwc_ww,jj0,bLR);
716 else if (iwater == esolTIP4P && jwater == esolTIP4P)
718 /* Interaction between two TIP4p molecules */
721 /* VdW only - only first atoms in each water interact */
722 add_j_to_nblist(vdw,jj0,bLR);
726 #ifdef DISABLE_WATERWATER_NLIST
727 /* Add entries for the four atoms - only do VdW if we need to */
730 add_j_to_nblist(vdw,jj0,bLR);
732 add_j_to_nblist(coul,jj0+1,bLR);
733 add_j_to_nblist(coul,jj0+2,bLR);
734 add_j_to_nblist(coul,jj0+3,bLR);
736 /* One entry for the entire water-water interaction */
739 add_j_to_nblist(coul_ww,jj0,bLR);
743 add_j_to_nblist(vdwc_ww,jj0,bLR);
750 /* j charge group is not water, but i is.
751 * Add entries to the water-other_atom lists; the geometry of the water
752 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
753 * so we don't care if it is SPC or TIP4P...
760 for(jj=jj0; (jj<jj1); jj++)
764 add_j_to_nblist(coul,jj,bLR);
770 for(jj=jj0; (jj<jj1); jj++)
772 if (bHaveVdW[type[jj]])
774 add_j_to_nblist(vdw,jj,bLR);
780 /* _charge_ _groups_ interact with both coulomb and LJ */
781 /* Check which atoms we should add to the lists! */
782 for(jj=jj0; (jj<jj1); jj++)
784 if (bHaveVdW[type[jj]])
788 add_j_to_nblist(vdwc,jj,bLR);
792 add_j_to_nblist(vdw,jj,bLR);
795 else if (charge[jj] != 0)
797 add_j_to_nblist(coul,jj,bLR);
804 close_i_nblist(coul);
805 close_i_nblist(vdwc);
806 #ifndef DISABLE_WATERWATER_NLIST
807 close_i_nblist(coul_ww);
808 close_i_nblist(vdwc_ww);
813 /* no solvent as i charge group */
814 /* Loop over the atoms in the i charge group */
815 for(i=0; i<nicg; i++)
818 gid = GID(igid,jgid,ngid);
821 /* Create new i_atom for each energy group */
822 if (bDoVdW && bDoCoul)
824 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
828 new_i_nblist(vdw,bLR,i_atom,shift,gid);
832 new_i_nblist(coul,bLR,i_atom,shift,gid);
834 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
835 bDoCoul_i = (bDoCoul && qi!=0);
837 if (bDoVdW_i || bDoCoul_i)
839 /* Loop over the j charge groups */
840 for(j=0; (j<nj); j++)
844 /* Check for large charge groups */
855 /* Finally loop over the atoms in the j-charge group */
856 for(jj=jj0; jj<jj1; jj++)
858 bNotEx = NOTEXCL(bExcl,i,jj);
866 add_j_to_nblist(coul,jj,bLR);
871 if (bHaveVdW[type[jj]])
873 add_j_to_nblist(vdw,jj,bLR);
878 if (bHaveVdW[type[jj]])
882 add_j_to_nblist(vdwc,jj,bLR);
886 add_j_to_nblist(vdw,jj,bLR);
889 else if (charge[jj] != 0)
891 add_j_to_nblist(coul,jj,bLR);
899 close_i_nblist(coul);
900 close_i_nblist(vdwc);
906 /* we are doing free energy */
907 vdwc_free = &nlist[eNL_VDWQQ_FREE];
908 vdw_free = &nlist[eNL_VDW_FREE];
909 coul_free = &nlist[eNL_QQ_FREE];
910 /* Loop over the atoms in the i charge group */
911 for(i=0; i<nicg; i++)
914 gid = GID(igid,jgid,ngid);
916 qiB = chargeB[i_atom];
918 /* Create new i_atom for each energy group */
919 if (bDoVdW && bDoCoul)
920 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
922 new_i_nblist(vdw,bLR,i_atom,shift,gid);
924 new_i_nblist(coul,bLR,i_atom,shift,gid);
926 new_i_nblist(vdw_free,bLR,i_atom,shift,gid);
927 new_i_nblist(coul_free,bLR,i_atom,shift,gid);
928 new_i_nblist(vdwc_free,bLR,i_atom,shift,gid);
930 bDoVdW_i = (bDoVdW &&
931 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
932 bDoCoul_i = (bDoCoul && (qi!=0 || qiB!=0));
933 /* For TIP4P the first atom does not have a charge,
934 * but the last three do. So we should still put an atom
935 * without LJ but with charge in the water-atom neighborlist
936 * for a TIP4p i charge group.
937 * For SPC type water the first atom has LJ and charge,
938 * so there is no such problem.
940 if (iwater == esolNO)
942 bDoCoul_i_sol = bDoCoul_i;
946 bDoCoul_i_sol = bDoCoul;
949 if (bDoVdW_i || bDoCoul_i_sol)
951 /* Loop over the j charge groups */
952 for(j=0; (j<nj); j++)
956 /* Check for large charge groups */
967 /* Finally loop over the atoms in the j-charge group */
968 bFree = bPert[i_atom];
969 for(jj=jj0; (jj<jj1); jj++)
971 bFreeJ = bFree || bPert[jj];
972 /* Complicated if, because the water H's should also
973 * see perturbed j-particles
975 if (iwater==esolNO || i==0 || bFreeJ)
977 bNotEx = NOTEXCL(bExcl,i,jj);
985 if (charge[jj]!=0 || chargeB[jj]!=0)
987 add_j_to_nblist(coul_free,jj,bLR);
992 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
994 add_j_to_nblist(vdw_free,jj,bLR);
999 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1001 if (charge[jj]!=0 || chargeB[jj]!=0)
1003 add_j_to_nblist(vdwc_free,jj,bLR);
1007 add_j_to_nblist(vdw_free,jj,bLR);
1010 else if (charge[jj]!=0 || chargeB[jj]!=0)
1011 add_j_to_nblist(coul_free,jj,bLR);
1016 /* This is done whether or not bWater is set */
1017 if (charge[jj] != 0)
1019 add_j_to_nblist(coul,jj,bLR);
1022 else if (!bDoCoul_i_sol)
1024 if (bHaveVdW[type[jj]])
1026 add_j_to_nblist(vdw,jj,bLR);
1031 if (bHaveVdW[type[jj]])
1033 if (charge[jj] != 0)
1035 add_j_to_nblist(vdwc,jj,bLR);
1039 add_j_to_nblist(vdw,jj,bLR);
1042 else if (charge[jj] != 0)
1044 add_j_to_nblist(coul,jj,bLR);
1052 close_i_nblist(vdw);
1053 close_i_nblist(coul);
1054 close_i_nblist(vdwc);
1055 close_i_nblist(vdw_free);
1056 close_i_nblist(coul_free);
1057 close_i_nblist(vdwc_free);
1063 put_in_list_qmmm(gmx_bool bHaveVdW[],
1080 int i,j,jcg,igid,gid;
1081 atom_id jj,jj0,jj1,i_atom;
1085 /* Get atom range */
1087 nicg = index[icg+1]-i0;
1089 /* Get the i charge group info */
1090 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1092 coul = &fr->QMMMlist;
1094 /* Loop over atoms in the ith charge group */
1095 for (i=0;i<nicg;i++)
1098 gid = GID(igid,jgid,ngid);
1099 /* Create new i_atom for each energy group */
1100 new_i_nblist(coul,bLR,i_atom,shift,gid);
1102 /* Loop over the j charge groups */
1107 /* Charge groups cannot have QM and MM atoms simultaneously */
1112 /* Finally loop over the atoms in the j-charge group */
1113 for(jj=jj0; jj<jj1; jj++)
1115 bNotEx = NOTEXCL(bExcl,i,jj);
1117 add_j_to_nblist(coul,jj,bLR);
1121 close_i_nblist(coul);
1126 put_in_list_cg(gmx_bool bHaveVdW[],
1143 int igid,gid,nbl_ind;
1147 cginfo = fr->cginfo[icg];
1149 igid = GET_CGINFO_GID(cginfo);
1150 gid = GID(igid,jgid,ngid);
1152 /* Unpack pointers to neighbourlist structs */
1153 if (fr->nnblists == 1)
1159 nbl_ind = fr->gid2nblists[gid];
1163 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1167 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1170 /* Make a new neighbor list for charge group icg.
1171 * Currently simply one neighbor list is made with LJ and Coulomb.
1172 * If required, zero interactions could be removed here
1173 * or in the force loop.
1175 new_i_nblist(vdwc,bLR,index[icg],shift,gid);
1176 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1178 for(j=0; (j<nj); j++)
1181 /* Skip the icg-icg pairs if all self interactions are excluded */
1182 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1184 /* Here we add the j charge group jcg to the list,
1185 * exclusions are also added to the list.
1187 add_j_to_nblist_cg(vdwc,index[jcg],index[jcg+1],bExcl,icg==jcg,bLR);
1191 close_i_nblist(vdwc);
1194 static void setexcl(atom_id start,atom_id end,t_blocka *excl,gmx_bool b,
1201 for(i=start; i<end; i++)
1203 for(k=excl->index[i]; k<excl->index[i+1]; k++)
1205 SETEXCL(bexcl,i-start,excl->a[k]);
1211 for(i=start; i<end; i++)
1213 for(k=excl->index[i]; k<excl->index[i+1]; k++)
1215 RMEXCL(bexcl,i-start,excl->a[k]);
1221 int calc_naaj(int icg,int cgtot)
1225 if ((cgtot % 2) == 1)
1227 /* Odd number of charge groups, easy */
1228 naaj = 1 + (cgtot/2);
1230 else if ((cgtot % 4) == 0)
1232 /* Multiple of four is hard */
1269 fprintf(log,"naaj=%d\n",naaj);
1275 /************************************************
1277 * S I M P L E C O R E S T U F F
1279 ************************************************/
1281 static real calc_image_tric(rvec xi,rvec xj,matrix box,
1282 rvec b_inv,int *shift)
1284 /* This code assumes that the cut-off is smaller than
1285 * a half times the smallest diagonal element of the box.
1292 /* Compute diff vector */
1293 dz = xj[ZZ] - xi[ZZ];
1294 dy = xj[YY] - xi[YY];
1295 dx = xj[XX] - xi[XX];
1297 /* Perform NINT operation, using trunc operation, therefore
1298 * we first add 2.5 then subtract 2 again
1300 tz = dz*b_inv[ZZ] + h25;
1302 dz -= tz*box[ZZ][ZZ];
1303 dy -= tz*box[ZZ][YY];
1304 dx -= tz*box[ZZ][XX];
1306 ty = dy*b_inv[YY] + h25;
1308 dy -= ty*box[YY][YY];
1309 dx -= ty*box[YY][XX];
1311 tx = dx*b_inv[XX]+h25;
1313 dx -= tx*box[XX][XX];
1315 /* Distance squared */
1316 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1318 *shift = XYZ2IS(tx,ty,tz);
1323 static real calc_image_rect(rvec xi,rvec xj,rvec box_size,
1324 rvec b_inv,int *shift)
1332 /* Compute diff vector */
1333 dx = xj[XX] - xi[XX];
1334 dy = xj[YY] - xi[YY];
1335 dz = xj[ZZ] - xi[ZZ];
1337 /* Perform NINT operation, using trunc operation, therefore
1338 * we first add 1.5 then subtract 1 again
1340 tx = dx*b_inv[XX] + h15;
1341 ty = dy*b_inv[YY] + h15;
1342 tz = dz*b_inv[ZZ] + h15;
1347 /* Correct diff vector for translation */
1348 ddx = tx*box_size[XX] - dx;
1349 ddy = ty*box_size[YY] - dy;
1350 ddz = tz*box_size[ZZ] - dz;
1352 /* Distance squared */
1353 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1355 *shift = XYZ2IS(tx,ty,tz);
1360 static void add_simple(t_ns_buf *nsbuf,int nrj,atom_id cg_j,
1361 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1362 int icg,int jgid,t_block *cgs,t_excl bexcl[],
1363 int shift,t_forcerec *fr,put_in_list_t *put_in_list)
1365 if (nsbuf->nj + nrj > MAX_CG)
1367 put_in_list(bHaveVdW,ngid,md,icg,jgid,nsbuf->ncg,nsbuf->jcg,
1368 cgs->index,bexcl,shift,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
1369 /* Reset buffer contents */
1370 nsbuf->ncg = nsbuf->nj = 0;
1372 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1376 static void ns_inner_tric(rvec x[],int icg,int *i_egp_flags,
1377 int njcg,atom_id jcg[],
1378 matrix box,rvec b_inv,real rcut2,
1379 t_block *cgs,t_ns_buf **ns_buf,
1380 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1381 t_excl bexcl[],t_forcerec *fr,
1382 put_in_list_t *put_in_list)
1386 int *cginfo=fr->cginfo;
1387 atom_id cg_j,*cgindex;
1390 cgindex = cgs->index;
1392 for(j=0; (j<njcg); j++)
1395 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1396 if (calc_image_tric(x[icg],x[cg_j],box,b_inv,&shift) < rcut2)
1398 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1399 if (!(i_egp_flags[jgid] & EGP_EXCL))
1401 add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1402 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1409 static void ns_inner_rect(rvec x[],int icg,int *i_egp_flags,
1410 int njcg,atom_id jcg[],
1411 gmx_bool bBox,rvec box_size,rvec b_inv,real rcut2,
1412 t_block *cgs,t_ns_buf **ns_buf,
1413 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1414 t_excl bexcl[],t_forcerec *fr,
1415 put_in_list_t *put_in_list)
1419 int *cginfo=fr->cginfo;
1420 atom_id cg_j,*cgindex;
1423 cgindex = cgs->index;
1427 for(j=0; (j<njcg); j++)
1430 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1431 if (calc_image_rect(x[icg],x[cg_j],box_size,b_inv,&shift) < rcut2)
1433 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1434 if (!(i_egp_flags[jgid] & EGP_EXCL))
1436 add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1437 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1445 for(j=0; (j<njcg); j++)
1448 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1449 if ((rcut2 == 0) || (distance2(x[icg],x[cg_j]) < rcut2)) {
1450 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1451 if (!(i_egp_flags[jgid] & EGP_EXCL))
1453 add_simple(&ns_buf[jgid][CENTRAL],nrj,cg_j,
1454 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,CENTRAL,fr,
1462 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1464 static int ns_simple_core(t_forcerec *fr,
1465 gmx_localtop_t *top,
1467 matrix box,rvec box_size,
1468 t_excl bexcl[],atom_id *aaj,
1469 int ngid,t_ns_buf **ns_buf,
1470 put_in_list_t *put_in_list,gmx_bool bHaveVdW[])
1474 int nsearch,icg,jcg,igid,i0,nri,nn;
1477 /* atom_id *i_atoms; */
1478 t_block *cgs=&(top->cgs);
1479 t_blocka *excl=&(top->excls);
1482 gmx_bool bBox,bTriclinic;
1485 rlist2 = sqr(fr->rlist);
1487 bBox = (fr->ePBC != epbcNONE);
1490 for(m=0; (m<DIM); m++)
1492 b_inv[m] = divide_err(1.0,box_size[m]);
1494 bTriclinic = TRICLINIC(box);
1501 cginfo = fr->cginfo;
1504 for (icg=fr->cg0; (icg<fr->hcg); icg++)
1507 i0 = cgs->index[icg];
1508 nri = cgs->index[icg+1]-i0;
1509 i_atoms = &(cgs->a[i0]);
1510 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1511 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1513 igid = GET_CGINFO_GID(cginfo[icg]);
1514 i_egp_flags = fr->egp_flags + ngid*igid;
1515 setexcl(cgs->index[icg],cgs->index[icg+1],excl,TRUE,bexcl);
1517 naaj=calc_naaj(icg,cgs->nr);
1520 ns_inner_tric(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1521 box,b_inv,rlist2,cgs,ns_buf,
1522 bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1526 ns_inner_rect(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1527 bBox,box_size,b_inv,rlist2,cgs,ns_buf,
1528 bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1532 for(nn=0; (nn<ngid); nn++)
1534 for(k=0; (k<SHIFTS); k++)
1536 nsbuf = &(ns_buf[nn][k]);
1539 put_in_list(bHaveVdW,ngid,md,icg,nn,nsbuf->ncg,nsbuf->jcg,
1540 cgs->index,bexcl,k,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
1541 nsbuf->ncg=nsbuf->nj=0;
1545 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1546 setexcl(cgs->index[icg],cgs->index[icg+1],excl,FALSE,bexcl);
1548 close_neighbor_lists(fr,FALSE);
1553 /************************************************
1555 * N S 5 G R I D S T U F F
1557 ************************************************/
1559 static inline void get_dx(int Nx,real gridx,real rc2,int xgi,real x,
1560 int *dx0,int *dx1,real *dcx2)
1588 for(i=xgi0; i>=0; i--)
1590 dcx = (i+1)*gridx-x;
1597 for(i=xgi1; i<Nx; i++)
1610 static inline void get_dx_dd(int Nx,real gridx,real rc2,int xgi,real x,
1611 int ncpddc,int shift_min,int shift_max,
1612 int *g0,int *g1,real *dcx2)
1615 int g_min,g_max,shift_home;
1648 g_min = (shift_min == shift_home ? 0 : ncpddc);
1649 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1656 else if (shift_max < 0)
1671 /* Check one grid cell down */
1672 dcx = ((*g0 - 1) + 1)*gridx - x;
1684 /* Check one grid cell up */
1685 dcx = (*g1 + 1)*gridx - x;
1697 #define sqr(x) ((x)*(x))
1698 #define calc_dx2(XI,YI,ZI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
1699 #define calc_cyl_dx2(XI,YI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
1700 /****************************************************
1702 * F A S T N E I G H B O R S E A R C H I N G
1704 * Optimized neighboursearching routine using grid
1705 * at least 1x1x1, see GROMACS manual
1707 ****************************************************/
1710 static void get_cutoff2(t_forcerec *fr,gmx_bool bDoLongRange,
1711 real *rvdw2,real *rcoul2,
1712 real *rs2,real *rm2,real *rl2)
1714 *rs2 = sqr(fr->rlist);
1716 if (bDoLongRange && fr->bTwinRange)
1718 /* The VdW and elec. LR cut-off's could be different,
1719 * so we can not simply set them to rlistlong.
1721 if (EVDW_MIGHT_BE_ZERO_AT_CUTOFF(fr->vdwtype) &&
1722 fr->rvdw > fr->rlist)
1724 *rvdw2 = sqr(fr->rlistlong);
1728 *rvdw2 = sqr(fr->rvdw);
1730 if (EEL_MIGHT_BE_ZERO_AT_CUTOFF(fr->eeltype) &&
1731 fr->rcoulomb > fr->rlist)
1733 *rcoul2 = sqr(fr->rlistlong);
1737 *rcoul2 = sqr(fr->rcoulomb);
1742 /* Workaround for a gcc -O3 or -ffast-math problem */
1746 *rm2 = min(*rvdw2,*rcoul2);
1747 *rl2 = max(*rvdw2,*rcoul2);
1750 static void init_nsgrid_lists(t_forcerec *fr,int ngid,gmx_ns_t *ns)
1752 real rvdw2,rcoul2,rs2,rm2,rl2;
1755 get_cutoff2(fr,TRUE,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
1757 /* Short range buffers */
1758 snew(ns->nl_sr,ngid);
1761 snew(ns->nlr_ljc,ngid);
1762 snew(ns->nlr_one,ngid);
1764 /* Always allocate both list types, since rcoulomb might now change with PME load balancing */
1765 /* Long range VdW and Coul buffers */
1766 snew(ns->nl_lr_ljc,ngid);
1767 /* Long range VdW or Coul only buffers */
1768 snew(ns->nl_lr_one,ngid);
1770 for(j=0; (j<ngid); j++) {
1771 snew(ns->nl_sr[j],MAX_CG);
1772 snew(ns->nl_lr_ljc[j],MAX_CG);
1773 snew(ns->nl_lr_one[j],MAX_CG);
1778 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
1783 static int nsgrid_core(FILE *log,t_commrec *cr,t_forcerec *fr,
1784 matrix box,rvec box_size,int ngid,
1785 gmx_localtop_t *top,
1786 t_grid *grid,rvec x[],
1787 t_excl bexcl[],gmx_bool *bExcludeAlleg,
1788 t_nrnb *nrnb,t_mdatoms *md,
1789 real *lambda,real *dvdlambda,
1790 gmx_grppairener_t *grppener,
1791 put_in_list_t *put_in_list,
1792 gmx_bool bHaveVdW[],
1793 gmx_bool bDoLongRange,gmx_bool bMakeQMMMnblist)
1796 atom_id **nl_lr_ljc,**nl_lr_one,**nl_sr;
1797 int *nlr_ljc,*nlr_one,*nsr;
1798 gmx_domdec_t *dd=NULL;
1799 t_block *cgs=&(top->cgs);
1800 int *cginfo=fr->cginfo;
1801 /* atom_id *i_atoms,*cgsindex=cgs->index; */
1803 int cell_x,cell_y,cell_z;
1804 int d,tx,ty,tz,dx,dy,dz,cj;
1805 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
1806 int zsh_ty,zsh_tx,ysh_tx;
1808 int dx0,dx1,dy0,dy1,dz0,dz1;
1809 int Nx,Ny,Nz,shift=-1,j,nrj,nns,nn=-1;
1810 real gridx,gridy,gridz,grid_x,grid_y,grid_z;
1811 real *dcx2,*dcy2,*dcz2;
1813 int cg0,cg1,icg=-1,cgsnr,i0,igid,nri,naaj,max_jcg;
1814 int jcg0,jcg1,jjcg,cgj0,jgid;
1815 int *grida,*gridnra,*gridind;
1816 gmx_bool rvdw_lt_rcoul,rcoul_lt_rvdw;
1817 rvec xi,*cgcm,grid_offset;
1818 real r2,rs2,rvdw2,rcoul2,rm2,rl2,XI,YI,ZI,dcx,dcy,dcz,tmp1,tmp2;
1820 gmx_bool bDomDec,bTriclinicX,bTriclinicY;
1825 bDomDec = DOMAINDECOMP(cr);
1831 bTriclinicX = ((YY < grid->npbcdim &&
1832 (!bDomDec || dd->nc[YY]==1) && box[YY][XX] != 0) ||
1833 (ZZ < grid->npbcdim &&
1834 (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][XX] != 0));
1835 bTriclinicY = (ZZ < grid->npbcdim &&
1836 (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][YY] != 0);
1840 get_cutoff2(fr,bDoLongRange,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
1842 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
1843 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
1845 if (bMakeQMMMnblist)
1853 nl_lr_ljc = ns->nl_lr_ljc;
1854 nl_lr_one = ns->nl_lr_one;
1855 nlr_ljc = ns->nlr_ljc;
1856 nlr_one = ns->nlr_one;
1864 gridind = grid->index;
1865 gridnra = grid->nra;
1868 gridx = grid->cell_size[XX];
1869 gridy = grid->cell_size[YY];
1870 gridz = grid->cell_size[ZZ];
1874 copy_rvec(grid->cell_offset,grid_offset);
1875 copy_ivec(grid->ncpddc,ncpddc);
1880 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
1881 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
1882 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
1883 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
1884 if (zsh_tx!=0 && ysh_tx!=0)
1886 /* This could happen due to rounding, when both ratios are 0.5 */
1895 /* We only want a list for the test particle */
1904 /* Set the shift range */
1905 for(d=0; d<DIM; d++)
1909 /* Check if we need periodicity shifts.
1910 * Without PBC or with domain decomposition we don't need them.
1912 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
1919 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
1930 /* Loop over charge groups */
1931 for(icg=cg0; (icg < cg1); icg++)
1933 igid = GET_CGINFO_GID(cginfo[icg]);
1934 /* Skip this charge group if all energy groups are excluded! */
1935 if (bExcludeAlleg[igid])
1940 i0 = cgs->index[icg];
1942 if (bMakeQMMMnblist)
1944 /* Skip this charge group if it is not a QM atom while making a
1945 * QM/MM neighbourlist
1947 if (md->bQM[i0]==FALSE)
1949 continue; /* MM particle, go to next particle */
1952 /* Compute the number of charge groups that fall within the control
1955 naaj = calc_naaj(icg,cgsnr);
1962 /* make a normal neighbourlist */
1966 /* Get the j charge-group and dd cell shift ranges */
1967 dd_get_ns_ranges(cr->dd,icg,&jcg0,&jcg1,sh0,sh1);
1972 /* Compute the number of charge groups that fall within the control
1975 naaj = calc_naaj(icg,cgsnr);
1981 /* The i-particle is awlways the test particle,
1982 * so we want all j-particles
1984 max_jcg = cgsnr - 1;
1988 max_jcg = jcg1 - cgsnr;
1993 i_egp_flags = fr->egp_flags + igid*ngid;
1995 /* Set the exclusions for the atoms in charge group icg using a bitmask */
1996 setexcl(i0,cgs->index[icg+1],&top->excls,TRUE,bexcl);
1998 ci2xyz(grid,icg,&cell_x,&cell_y,&cell_z);
2000 /* Changed iicg to icg, DvdS 990115
2001 * (but see consistency check above, DvdS 990330)
2004 fprintf(log,"icg=%5d, naaj=%5d, cell %d %d %d\n",
2005 icg,naaj,cell_x,cell_y,cell_z);
2007 /* Loop over shift vectors in three dimensions */
2008 for (tz=-shp[ZZ]; tz<=shp[ZZ]; tz++)
2010 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2011 /* Calculate range of cells in Z direction that have the shift tz */
2012 zgi = cell_z + tz*Nz;
2015 get_dx(Nz,gridz,rl2,zgi,ZI,&dz0,&dz1,dcz2);
2017 get_dx_dd(Nz,gridz,rl2,zgi,ZI-grid_offset[ZZ],
2018 ncpddc[ZZ],sh0[ZZ],sh1[ZZ],&dz0,&dz1,dcz2);
2024 for (ty=-shp[YY]; ty<=shp[YY]; ty++)
2026 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2027 /* Calculate range of cells in Y direction that have the shift ty */
2030 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2034 ygi = cell_y + ty*Ny;
2037 get_dx(Ny,gridy,rl2,ygi,YI,&dy0,&dy1,dcy2);
2039 get_dx_dd(Ny,gridy,rl2,ygi,YI-grid_offset[YY],
2040 ncpddc[YY],sh0[YY],sh1[YY],&dy0,&dy1,dcy2);
2046 for (tx=-shp[XX]; tx<=shp[XX]; tx++)
2048 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2049 /* Calculate range of cells in X direction that have the shift tx */
2052 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2056 xgi = cell_x + tx*Nx;
2059 get_dx(Nx,gridx,rl2,xgi*Nx,XI,&dx0,&dx1,dcx2);
2061 get_dx_dd(Nx,gridx,rl2,xgi,XI-grid_offset[XX],
2062 ncpddc[XX],sh0[XX],sh1[XX],&dx0,&dx1,dcx2);
2068 /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2069 * from the neigbour list as it will not interact */
2070 if (fr->adress_type != eAdressOff){
2071 if (md->wf[cgs->index[icg]]==0 && egp_explicit(fr, igid)){
2075 /* Get shift vector */
2076 shift=XYZ2IS(tx,ty,tz);
2078 range_check(shift,0,SHIFTS);
2080 for(nn=0; (nn<ngid); nn++)
2087 fprintf(log,"shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2088 shift,dx0,dx1,dy0,dy1,dz0,dz1);
2089 fprintf(log,"cgcm: %8.3f %8.3f %8.3f\n",cgcm[icg][XX],
2090 cgcm[icg][YY],cgcm[icg][ZZ]);
2091 fprintf(log,"xi: %8.3f %8.3f %8.3f\n",XI,YI,ZI);
2093 for (dx=dx0; (dx<=dx1); dx++)
2095 tmp1 = rl2 - dcx2[dx];
2096 for (dy=dy0; (dy<=dy1); dy++)
2098 tmp2 = tmp1 - dcy2[dy];
2101 for (dz=dz0; (dz<=dz1); dz++) {
2102 if (tmp2 > dcz2[dz]) {
2103 /* Find grid-cell cj in which possible neighbours are */
2104 cj = xyz2ci(Ny,Nz,dx,dy,dz);
2106 /* Check out how many cgs (nrj) there in this cell */
2109 /* Find the offset in the cg list */
2112 /* Check if all j's are out of range so we
2113 * can skip the whole cell.
2114 * Should save some time, especially with DD.
2117 (grida[cgj0] >= max_jcg &&
2118 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2124 for (j=0; (j<nrj); j++)
2126 jjcg = grida[cgj0+j];
2128 /* check whether this guy is in range! */
2129 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2132 r2=calc_dx2(XI,YI,ZI,cgcm[jjcg]);
2134 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2135 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2136 /* check energy group exclusions */
2137 if (!(i_egp_flags[jgid] & EGP_EXCL))
2141 if (nsr[jgid] >= MAX_CG)
2143 /* Add to short-range list */
2144 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2145 nsr[jgid],nl_sr[jgid],
2146 cgs->index,/* cgsatoms, */ bexcl,
2147 shift,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
2150 nl_sr[jgid][nsr[jgid]++]=jjcg;
2154 if (nlr_ljc[jgid] >= MAX_CG)
2156 /* Add to LJ+coulomb long-range list */
2157 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2158 nlr_ljc[jgid],nl_lr_ljc[jgid],top->cgs.index,
2159 bexcl,shift,fr,TRUE,TRUE,TRUE,fr->solvent_opt);
2162 nl_lr_ljc[jgid][nlr_ljc[jgid]++]=jjcg;
2166 if (nlr_one[jgid] >= MAX_CG)
2168 /* Add to long-range list with only coul, or only LJ */
2169 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2170 nlr_one[jgid],nl_lr_one[jgid],top->cgs.index,
2171 bexcl,shift,fr,TRUE,rvdw_lt_rcoul,rcoul_lt_rvdw,fr->solvent_opt);
2174 nl_lr_one[jgid][nlr_one[jgid]++]=jjcg;
2186 /* CHECK whether there is anything left in the buffers */
2187 for(nn=0; (nn<ngid); nn++)
2191 put_in_list(bHaveVdW,ngid,md,icg,nn,nsr[nn],nl_sr[nn],
2192 cgs->index, /* cgsatoms, */ bexcl,
2193 shift,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
2196 if (nlr_ljc[nn] > 0)
2198 put_in_list(bHaveVdW,ngid,md,icg,nn,nlr_ljc[nn],
2199 nl_lr_ljc[nn],top->cgs.index,
2200 bexcl,shift,fr,TRUE,TRUE,TRUE,fr->solvent_opt);
2203 if (nlr_one[nn] > 0)
2205 put_in_list(bHaveVdW,ngid,md,icg,nn,nlr_one[nn],
2206 nl_lr_one[nn],top->cgs.index,
2207 bexcl,shift,fr,TRUE,rvdw_lt_rcoul,rcoul_lt_rvdw,fr->solvent_opt);
2213 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2214 setexcl(cgs->index[icg],cgs->index[icg+1],&top->excls,FALSE,bexcl);
2216 /* No need to perform any left-over force calculations anymore (as we used to do here)
2217 * since we now save the proper long-range lists for later evaluation.
2222 /* Close neighbourlists */
2223 close_neighbor_lists(fr,bMakeQMMMnblist);
2228 void ns_realloc_natoms(gmx_ns_t *ns,int natoms)
2232 if (natoms > ns->nra_alloc)
2234 ns->nra_alloc = over_alloc_dd(natoms);
2235 srenew(ns->bexcl,ns->nra_alloc);
2236 for(i=0; i<ns->nra_alloc; i++)
2243 void init_ns(FILE *fplog,const t_commrec *cr,
2244 gmx_ns_t *ns,t_forcerec *fr,
2245 const gmx_mtop_t *mtop,
2248 int mt,icg,nr_in_cg,maxcg,i,j,jcg,ngid,ncg;
2252 /* Compute largest charge groups size (# atoms) */
2254 for(mt=0; mt<mtop->nmoltype; mt++) {
2255 cgs = &mtop->moltype[mt].cgs;
2256 for (icg=0; (icg < cgs->nr); icg++)
2258 nr_in_cg=max(nr_in_cg,(int)(cgs->index[icg+1]-cgs->index[icg]));
2262 /* Verify whether largest charge group is <= max cg.
2263 * This is determined by the type of the local exclusion type
2264 * Exclusions are stored in bits. (If the type is not large
2265 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2267 maxcg = sizeof(t_excl)*8;
2268 if (nr_in_cg > maxcg)
2270 gmx_fatal(FARGS,"Max #atoms in a charge group: %d > %d\n",
2274 ngid = mtop->groups.grps[egcENER].nr;
2275 snew(ns->bExcludeAlleg,ngid);
2276 for(i=0; i<ngid; i++) {
2277 ns->bExcludeAlleg[i] = TRUE;
2278 for(j=0; j<ngid; j++)
2280 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2282 ns->bExcludeAlleg[i] = FALSE;
2289 ns->grid = init_grid(fplog,fr);
2290 init_nsgrid_lists(fr,ngid,ns);
2295 snew(ns->ns_buf,ngid);
2296 for(i=0; (i<ngid); i++)
2298 snew(ns->ns_buf[i],SHIFTS);
2300 ncg = ncg_mtop(mtop);
2301 snew(ns->simple_aaj,2*ncg);
2302 for(jcg=0; (jcg<ncg); jcg++)
2304 ns->simple_aaj[jcg] = jcg;
2305 ns->simple_aaj[jcg+ncg] = jcg;
2309 /* Create array that determines whether or not atoms have VdW */
2310 snew(ns->bHaveVdW,fr->ntype);
2311 for(i=0; (i<fr->ntype); i++)
2313 for(j=0; (j<fr->ntype); j++)
2315 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2317 ((BHAMA(fr->nbfp,fr->ntype,i,j) != 0) ||
2318 (BHAMB(fr->nbfp,fr->ntype,i,j) != 0) ||
2319 (BHAMC(fr->nbfp,fr->ntype,i,j) != 0)) :
2320 ((C6(fr->nbfp,fr->ntype,i,j) != 0) ||
2321 (C12(fr->nbfp,fr->ntype,i,j) != 0))));
2325 pr_bvec(debug,0,"bHaveVdW",ns->bHaveVdW,fr->ntype,TRUE);
2329 if (!DOMAINDECOMP(cr))
2331 /* This could be reduced with particle decomposition */
2332 ns_realloc_natoms(ns,mtop->natoms);
2335 ns->nblist_initialized=FALSE;
2337 /* nbr list debug dump */
2339 char *ptr=getenv("GMX_DUMP_NL");
2342 ns->dump_nl=strtol(ptr,NULL,10);
2345 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2356 int search_neighbours(FILE *log,t_forcerec *fr,
2357 rvec x[],matrix box,
2358 gmx_localtop_t *top,
2359 gmx_groups_t *groups,
2361 t_nrnb *nrnb,t_mdatoms *md,
2362 real *lambda,real *dvdlambda,
2363 gmx_grppairener_t *grppener,
2365 gmx_bool bDoLongRangeNS)
2367 t_block *cgs=&(top->cgs);
2368 rvec box_size,grid_x0,grid_x1;
2370 real min_size,grid_dens;
2374 gmx_bool *i_egp_flags;
2375 int cg_start,cg_end,start,end;
2378 gmx_domdec_zones_t *dd_zones;
2379 put_in_list_t *put_in_list;
2383 /* Set some local variables */
2385 ngid = groups->grps[egcENER].nr;
2387 for(m=0; (m<DIM); m++)
2389 box_size[m] = box[m][m];
2392 if (fr->ePBC != epbcNONE)
2394 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC,box))
2396 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.");
2400 min_size = min(box_size[XX],min(box_size[YY],box_size[ZZ]));
2401 if (2*fr->rlistlong >= min_size)
2402 gmx_fatal(FARGS,"One of the box diagonal elements has become smaller than twice the cut-off length.");
2406 if (DOMAINDECOMP(cr))
2408 ns_realloc_natoms(ns,cgs->index[cgs->nr]);
2412 /* Reset the neighbourlists */
2413 reset_neighbor_lists(fr,TRUE,TRUE);
2415 if (bGrid && bFillGrid)
2419 if (DOMAINDECOMP(cr))
2421 dd_zones = domdec_zones(cr->dd);
2427 get_nsgrid_boundaries(grid->nboundeddim,box,NULL,NULL,NULL,NULL,
2428 cgs->nr,fr->cg_cm,grid_x0,grid_x1,&grid_dens);
2430 grid_first(log,grid,NULL,NULL,fr->ePBC,box,grid_x0,grid_x1,
2431 fr->rlistlong,grid_dens);
2435 /* Don't know why this all is... (DvdS 3/99) */
2441 end = (cgs->nr+1)/2;
2444 if (DOMAINDECOMP(cr))
2447 fill_grid(log,dd_zones,grid,end,-1,end,fr->cg_cm);
2449 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2453 fill_grid(log,NULL,grid,cgs->nr,fr->cg0,fr->hcg,fr->cg_cm);
2454 grid->icg0 = fr->cg0;
2455 grid->icg1 = fr->hcg;
2463 calc_elemnr(log,grid,start,end,cgs->nr);
2465 grid_last(log,grid,start,end,cgs->nr);
2469 check_grid(debug,grid);
2470 print_grid(debug,grid);
2475 /* Set the grid cell index for the test particle only.
2476 * The cell to cg index is not corrected, but that does not matter.
2478 fill_grid(log,NULL,ns->grid,fr->hcg,fr->hcg-1,fr->hcg,fr->cg_cm);
2482 if (!fr->ns.bCGlist)
2484 put_in_list = put_in_list_at;
2488 put_in_list = put_in_list_cg;
2495 nsearch = nsgrid_core(log,cr,fr,box,box_size,ngid,top,
2496 grid,x,ns->bexcl,ns->bExcludeAlleg,
2497 nrnb,md,lambda,dvdlambda,grppener,
2498 put_in_list,ns->bHaveVdW,
2499 bDoLongRangeNS,FALSE);
2501 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2502 * the classical calculation. The charge-charge interaction
2503 * between QM and MM atoms is handled in the QMMM core calculation
2504 * (see QMMM.c). The VDW however, we'd like to compute classically
2505 * and the QM MM atom pairs have just been put in the
2506 * corresponding neighbourlists. in case of QMMM we still need to
2507 * fill a special QMMM neighbourlist that contains all neighbours
2508 * of the QM atoms. If bQMMM is true, this list will now be made:
2510 if (fr->bQMMM && fr->qr->QMMMscheme!=eQMMMschemeoniom)
2512 nsearch += nsgrid_core(log,cr,fr,box,box_size,ngid,top,
2513 grid,x,ns->bexcl,ns->bExcludeAlleg,
2514 nrnb,md,lambda,dvdlambda,grppener,
2515 put_in_list_qmmm,ns->bHaveVdW,
2516 bDoLongRangeNS,TRUE);
2521 nsearch = ns_simple_core(fr,top,md,box,box_size,
2522 ns->bexcl,ns->simple_aaj,
2523 ngid,ns->ns_buf,put_in_list,ns->bHaveVdW);
2531 inc_nrnb(nrnb,eNR_NS,nsearch);
2532 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2537 int natoms_beyond_ns_buffer(t_inputrec *ir,t_forcerec *fr,t_block *cgs,
2538 matrix scale_tot,rvec *x)
2540 int cg0,cg1,cg,a0,a1,a,i,j;
2541 real rint,hbuf2,scale;
2543 gmx_bool bIsotropic;
2548 rint = max(ir->rcoulomb,ir->rvdw);
2549 if (ir->rlist < rint)
2551 gmx_fatal(FARGS,"The neighbor search buffer has negative size: %f nm",
2559 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2561 hbuf2 = sqr(0.5*(ir->rlist - rint));
2562 for(cg=cg0; cg<cg1; cg++)
2564 a0 = cgs->index[cg];
2565 a1 = cgs->index[cg+1];
2566 for(a=a0; a<a1; a++)
2568 if (distance2(cg_cm[cg],x[a]) > hbuf2)
2578 scale = scale_tot[0][0];
2579 for(i=1; i<DIM; i++)
2581 /* With anisotropic scaling, the original spherical ns volumes become
2582 * ellipsoids. To avoid costly transformations we use the minimum
2583 * eigenvalue of the scaling matrix for determining the buffer size.
2584 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2586 scale = min(scale,scale_tot[i][i]);
2587 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2593 if (scale_tot[i][j] != 0)
2599 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2602 for(cg=cg0; cg<cg1; cg++)
2604 svmul(scale,cg_cm[cg],cgsc);
2605 a0 = cgs->index[cg];
2606 a1 = cgs->index[cg+1];
2607 for(a=a0; a<a1; a++)
2609 if (distance2(cgsc,x[a]) > hbuf2)
2618 /* Anistropic scaling */
2619 for(cg=cg0; cg<cg1; cg++)
2621 /* Since scale_tot contains the transpose of the scaling matrix,
2622 * we need to multiply with the transpose.
2624 tmvmul_ur0(scale_tot,cg_cm[cg],cgsc);
2625 a0 = cgs->index[cg];
2626 a1 = cgs->index[cg+1];
2627 for(a=a0; a<a1; a++)
2629 if (distance2(cgsc,x[a]) > hbuf2)