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"
67 * E X C L U S I O N H A N D L I N G
71 static void SETEXCL_(t_excl e[],atom_id i,atom_id j)
72 { e[j] = e[j] | (1<<i); }
73 static void RMEXCL_(t_excl e[],atom_id i,atom_id j)
74 { e[j]=e[j] & ~(1<<i); }
75 static gmx_bool ISEXCL_(t_excl e[],atom_id i,atom_id j)
76 { return (gmx_bool)(e[j] & (1<<i)); }
77 static gmx_bool NOTEXCL_(t_excl e[],atom_id i,atom_id j)
78 { return !(ISEXCL(e,i,j)); }
80 #define SETEXCL(e,i,j) (e)[((atom_id) (j))] |= (1<<((atom_id) (i)))
81 #define RMEXCL(e,i,j) (e)[((atom_id) (j))] &= (~(1<<((atom_id) (i))))
82 #define ISEXCL(e,i,j) (gmx_bool) ((e)[((atom_id) (j))] & (1<<((atom_id) (i))))
83 #define NOTEXCL(e,i,j) !(ISEXCL(e,i,j))
86 /************************************************
88 * U T I L I T I E S F O R N S
90 ************************************************/
92 static void reallocate_nblist(t_nblist *nl)
96 fprintf(debug,"reallocating neigborlist il_code=%d, maxnri=%d\n",
97 nl->il_code,nl->maxnri);
99 srenew(nl->iinr, nl->maxnri);
100 if (nl->enlist == enlistCG_CG)
102 srenew(nl->iinr_end,nl->maxnri);
104 srenew(nl->gid, nl->maxnri);
105 srenew(nl->shift, nl->maxnri);
106 srenew(nl->jindex, nl->maxnri+1);
109 /* ivdw/icoul are used to determine the type of interaction, so we
110 * can set an innerloop index here. The obvious choice for this would have
111 * been the vdwtype/coultype values in the forcerecord, but unfortunately
112 * those types are braindead - for instance both Buckingham and normal
113 * Lennard-Jones use the same value (evdwCUT), and a separate gmx_boolean variable
114 * to determine which interaction is used. There is further no special value
115 * for 'no interaction'. For backward compatibility with old TPR files we won't
116 * change this in the 3.x series, so when calling this routine you should use:
118 * icoul=0 no coulomb interaction
119 * icoul=1 cutoff standard coulomb
120 * icoul=2 reaction-field coulomb
121 * icoul=3 tabulated coulomb
123 * ivdw=0 no vdw interaction
124 * ivdw=1 standard L-J interaction
126 * ivdw=3 tabulated vdw.
128 * Kind of ugly, but it works.
130 static void init_nblist(t_nblist *nl_sr,t_nblist *nl_lr,
133 gmx_bool bfree, int enlist)
165 nl = (i == 0) ? nl_sr : nl_lr;
166 homenr = (i == 0) ? maxsr : maxlr;
173 /* Set coul/vdw in neighborlist, and for the normal loops we determine
174 * an index of which one to call.
178 nl->free_energy = bfree;
182 nl->enlist = enlistATOM_ATOM;
183 nl->il_code = eNR_NBKERNEL_FREE_ENERGY;
189 nn = inloop[4*icoul + ivdw];
191 /* solvent loops follow directly after the corresponding
192 * ordinary loops, in the order:
194 * SPC, SPC-SPC, TIP4p, TIP4p-TIP4p
198 case enlistATOM_ATOM:
201 case enlistSPC_ATOM: nn += 1; break;
202 case enlistSPC_SPC: nn += 2; break;
203 case enlistTIP4P_ATOM: nn += 3; break;
204 case enlistTIP4P_TIP4P: nn += 4; break;
211 fprintf(debug,"Initiating neighbourlist type %d for %s interactions,\nwith %d SR, %d LR atoms.\n",
212 nl->il_code,ENLISTTYPE(enlist),maxsr,maxlr);
214 /* maxnri is influenced by the number of shifts (maximum is 8)
215 * and the number of energy groups.
216 * If it is not enough, nl memory will be reallocated during the run.
217 * 4 seems to be a reasonable factor, which only causes reallocation
218 * during runs with tiny and many energygroups.
220 nl->maxnri = homenr*4;
229 reallocate_nblist(nl);
231 #ifdef GMX_THREAD_SHM_FDECOMP
234 pthread_mutex_init(nl->mtx,NULL);
239 void init_neighbor_list(FILE *log,t_forcerec *fr,int homenr)
241 /* Make maxlr tunable! (does not seem to be a big difference though)
242 * This parameter determines the number of i particles in a long range
243 * neighbourlist. Too few means many function calls, too many means
246 int maxsr,maxsr_wat,maxlr,maxlr_wat;
247 int icoul,icoulf,ivdw;
249 int enlist_def,enlist_w,enlist_ww;
253 /* maxsr = homenr-fr->nWatMol*3; */
258 gmx_fatal(FARGS,"%s, %d: Negative number of short range atoms.\n"
259 "Call your Gromacs dealer for assistance.",__FILE__,__LINE__);
261 /* This is just for initial allocation, so we do not reallocate
262 * all the nlist arrays many times in a row.
263 * The numbers seem very accurate, but they are uncritical.
265 maxsr_wat = min(fr->nWatMol,(homenr+2)/3);
269 maxlr_wat = min(maxsr_wat,maxlr);
273 maxlr = maxlr_wat = 0;
276 /* Determine the values for icoul/ivdw. */
282 else if (fr->bcoultab)
286 else if (EEL_RF(fr->eeltype))
308 fr->ns.bCGlist = (getenv("GMX_NBLISTCG") != 0);
311 enlist_def = enlistATOM_ATOM;
315 enlist_def = enlistCG_CG;
318 fprintf(log,"\nUsing charge-group - charge-group neighbor lists and kernels\n\n");
322 if (fr->solvent_opt == esolTIP4P) {
323 enlist_w = enlistTIP4P_ATOM;
324 enlist_ww = enlistTIP4P_TIP4P;
326 enlist_w = enlistSPC_ATOM;
327 enlist_ww = enlistSPC_SPC;
330 for(i=0; i<fr->nnblists; i++)
332 nbl = &(fr->nblists[i]);
333 init_nblist(&nbl->nlist_sr[eNL_VDWQQ],&nbl->nlist_lr[eNL_VDWQQ],
334 maxsr,maxlr,ivdw,icoul,FALSE,enlist_def);
335 init_nblist(&nbl->nlist_sr[eNL_VDW],&nbl->nlist_lr[eNL_VDW],
336 maxsr,maxlr,ivdw,0,FALSE,enlist_def);
337 init_nblist(&nbl->nlist_sr[eNL_QQ],&nbl->nlist_lr[eNL_QQ],
338 maxsr,maxlr,0,icoul,FALSE,enlist_def);
339 init_nblist(&nbl->nlist_sr[eNL_VDWQQ_WATER],&nbl->nlist_lr[eNL_VDWQQ_WATER],
340 maxsr_wat,maxlr_wat,ivdw,icoul, FALSE,enlist_w);
341 init_nblist(&nbl->nlist_sr[eNL_QQ_WATER],&nbl->nlist_lr[eNL_QQ_WATER],
342 maxsr_wat,maxlr_wat,0,icoul, FALSE,enlist_w);
343 init_nblist(&nbl->nlist_sr[eNL_VDWQQ_WATERWATER],&nbl->nlist_lr[eNL_VDWQQ_WATERWATER],
344 maxsr_wat,maxlr_wat,ivdw,icoul, FALSE,enlist_ww);
345 init_nblist(&nbl->nlist_sr[eNL_QQ_WATERWATER],&nbl->nlist_lr[eNL_QQ_WATERWATER],
346 maxsr_wat,maxlr_wat,0,icoul, FALSE,enlist_ww);
348 if (fr->efep != efepNO)
359 init_nblist(&nbl->nlist_sr[eNL_VDWQQ_FREE],&nbl->nlist_lr[eNL_VDWQQ_FREE],
360 maxsr,maxlr,ivdw,icoulf,TRUE,enlistATOM_ATOM);
361 init_nblist(&nbl->nlist_sr[eNL_VDW_FREE],&nbl->nlist_lr[eNL_VDW_FREE],
362 maxsr,maxlr,ivdw,0,TRUE,enlistATOM_ATOM);
363 init_nblist(&nbl->nlist_sr[eNL_QQ_FREE],&nbl->nlist_lr[eNL_QQ_FREE],
364 maxsr,maxlr,0,icoulf,TRUE,enlistATOM_ATOM);
368 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
370 init_nblist(&fr->QMMMlist,NULL,
371 maxsr,maxlr,0,icoul,FALSE,enlistATOM_ATOM);
374 fr->ns.nblist_initialized=TRUE;
377 static void reset_nblist(t_nblist *nl)
388 static void reset_neighbor_list(t_forcerec *fr,gmx_bool bLR,int nls,int eNL)
394 reset_nblist(&(fr->nblists[nls].nlist_lr[eNL]));
398 for(n=0; n<fr->nnblists; n++)
400 for(i=0; i<eNL_NR; i++)
402 reset_nblist(&(fr->nblists[n].nlist_sr[i]));
407 /* only reset the short-range nblist */
408 reset_nblist(&(fr->QMMMlist));
416 static inline void new_i_nblist(t_nblist *nlist,
417 gmx_bool bLR,atom_id i_atom,int shift,int gid)
423 /* Check whether we have to increase the i counter */
425 (nlist->iinr[nri] != i_atom) ||
426 (nlist->shift[nri] != shift) ||
427 (nlist->gid[nri] != gid))
429 /* This is something else. Now see if any entries have
430 * been added in the list of the previous atom.
433 ((nlist->jindex[nri+1] > nlist->jindex[nri]) &&
434 (nlist->gid[nri] != -1)))
436 /* If so increase the counter */
439 if (nlist->nri >= nlist->maxnri)
441 nlist->maxnri += over_alloc_large(nlist->nri);
442 reallocate_nblist(nlist);
445 /* Set the number of neighbours and the atom number */
446 nlist->jindex[nri+1] = nlist->jindex[nri];
447 nlist->iinr[nri] = i_atom;
448 nlist->gid[nri] = gid;
449 nlist->shift[nri] = shift;
453 static inline void close_i_nblist(t_nblist *nlist)
455 int nri = nlist->nri;
460 nlist->jindex[nri+1] = nlist->nrj;
462 len=nlist->nrj - nlist->jindex[nri];
464 /* nlist length for water i molecules is treated statically
467 if (len > nlist->maxlen)
474 static inline void close_nblist(t_nblist *nlist)
476 /* Only close this nblist when it has been initialized.
477 * Avoid the creation of i-lists with no j-particles.
481 /* Some assembly kernels do not support empty lists,
482 * make sure here that we don't generate any empty lists.
483 * With the current ns code this branch is taken in two cases:
484 * No i-particles at all: nri=-1 here
485 * There are i-particles, but no j-particles; nri=0 here
491 /* Close list number nri by incrementing the count */
496 static inline void close_neighbor_list(t_forcerec *fr,gmx_bool bLR,int nls,int eNL,
497 gmx_bool bMakeQMMMnblist)
501 if (bMakeQMMMnblist) {
504 close_nblist(&(fr->QMMMlist));
511 close_nblist(&(fr->nblists[nls].nlist_lr[eNL]));
515 for(n=0; n<fr->nnblists; n++)
517 for(i=0; (i<eNL_NR); i++)
519 close_nblist(&(fr->nblists[n].nlist_sr[i]));
526 static inline void add_j_to_nblist(t_nblist *nlist,atom_id j_atom,gmx_bool bLR)
530 if (nlist->nrj >= nlist->maxnrj)
532 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
534 fprintf(debug,"Increasing %s nblist %s j size to %d\n",
535 bLR ? "LR" : "SR",nrnb_str(nlist->il_code),nlist->maxnrj);
537 srenew(nlist->jjnr,nlist->maxnrj);
540 nlist->jjnr[nrj] = j_atom;
544 static inline void add_j_to_nblist_cg(t_nblist *nlist,
545 atom_id j_start,int j_end,
546 t_excl *bexcl,gmx_bool i_is_j,
552 if (nlist->nrj >= nlist->maxnrj)
554 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
556 fprintf(debug,"Increasing %s nblist %s j size to %d\n",
557 bLR ? "LR" : "SR",nrnb_str(nlist->il_code),nlist->maxnrj);
559 srenew(nlist->jjnr ,nlist->maxnrj);
560 srenew(nlist->jjnr_end,nlist->maxnrj);
561 srenew(nlist->excl ,nlist->maxnrj*MAX_CGCGSIZE);
564 nlist->jjnr[nrj] = j_start;
565 nlist->jjnr_end[nrj] = j_end;
567 if (j_end - j_start > MAX_CGCGSIZE)
569 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);
572 /* Set the exclusions */
573 for(j=j_start; j<j_end; j++)
575 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
579 /* Avoid double counting of intra-cg interactions */
580 for(j=1; j<j_end-j_start; j++)
582 nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
590 put_in_list_t(gmx_bool bHaveVdW[],
606 put_in_list_at(gmx_bool bHaveVdW[],
621 /* The a[] index has been removed,
622 * to put it back in i_atom should be a[i0] and jj should be a[jj].
627 t_nblist * vdwc_free = NULL;
628 t_nblist * vdw_free = NULL;
629 t_nblist * coul_free = NULL;
630 t_nblist * vdwc_ww = NULL;
631 t_nblist * coul_ww = NULL;
633 int i,j,jcg,igid,gid,nbl_ind,ind_ij;
634 atom_id jj,jj0,jj1,i_atom;
639 real *charge,*chargeB;
641 gmx_bool bFreeEnergy,bFree,bFreeJ,bNotEx,*bPert;
642 gmx_bool bDoVdW_i,bDoCoul_i,bDoCoul_i_sol;
646 /* Copy some pointers */
648 charge = md->chargeA;
649 chargeB = md->chargeB;
652 bPert = md->bPerturbed;
656 nicg = index[icg+1]-i0;
658 /* Get the i charge group info */
659 igid = GET_CGINFO_GID(cginfo[icg]);
660 iwater = GET_CGINFO_SOLOPT(cginfo[icg]);
665 /* Check if any of the particles involved are perturbed.
666 * If not we can do the cheaper normal put_in_list
667 * and use more solvent optimization.
669 for(i=0; i<nicg; i++)
671 bFreeEnergy |= bPert[i0+i];
673 /* Loop over the j charge groups */
674 for(j=0; (j<nj && !bFreeEnergy); j++)
679 /* Finally loop over the atoms in the j-charge group */
680 for(jj=jj0; jj<jj1; jj++)
682 bFreeEnergy |= bPert[jj];
687 /* Unpack pointers to neighbourlist structs */
688 if (fr->nnblists == 1)
694 nbl_ind = fr->gid2nblists[GID(igid,jgid,ngid)];
698 nlist = fr->nblists[nbl_ind].nlist_lr;
702 nlist = fr->nblists[nbl_ind].nlist_sr;
705 if (iwater != esolNO)
707 vdwc = &nlist[eNL_VDWQQ_WATER];
708 vdw = &nlist[eNL_VDW];
709 coul = &nlist[eNL_QQ_WATER];
710 #ifndef DISABLE_WATERWATER_NLIST
711 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
712 coul_ww = &nlist[eNL_QQ_WATERWATER];
717 vdwc = &nlist[eNL_VDWQQ];
718 vdw = &nlist[eNL_VDW];
719 coul = &nlist[eNL_QQ];
724 if (iwater != esolNO)
726 /* Loop over the atoms in the i charge group */
728 gid = GID(igid,jgid,ngid);
729 /* Create new i_atom for each energy group */
730 if (bDoCoul && bDoVdW)
732 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
733 #ifndef DISABLE_WATERWATER_NLIST
734 new_i_nblist(vdwc_ww,bLR,i_atom,shift,gid);
739 new_i_nblist(vdw,bLR,i_atom,shift,gid);
743 new_i_nblist(coul,bLR,i_atom,shift,gid);
744 #ifndef DISABLE_WATERWATER_NLIST
745 new_i_nblist(coul_ww,bLR,i_atom,shift,gid);
748 /* Loop over the j charge groups */
749 for(j=0; (j<nj); j++)
759 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
761 if (iwater == esolSPC && jwater == esolSPC)
763 /* Interaction between two SPC molecules */
766 /* VdW only - only first atoms in each water interact */
767 add_j_to_nblist(vdw,jj0,bLR);
771 #ifdef DISABLE_WATERWATER_NLIST
772 /* Add entries for the three atoms - only do VdW if we need to */
775 add_j_to_nblist(coul,jj0,bLR);
779 add_j_to_nblist(vdwc,jj0,bLR);
781 add_j_to_nblist(coul,jj0+1,bLR);
782 add_j_to_nblist(coul,jj0+2,bLR);
784 /* One entry for the entire water-water interaction */
787 add_j_to_nblist(coul_ww,jj0,bLR);
791 add_j_to_nblist(vdwc_ww,jj0,bLR);
796 else if (iwater == esolTIP4P && jwater == esolTIP4P)
798 /* Interaction between two TIP4p molecules */
801 /* VdW only - only first atoms in each water interact */
802 add_j_to_nblist(vdw,jj0,bLR);
806 #ifdef DISABLE_WATERWATER_NLIST
807 /* Add entries for the four atoms - only do VdW if we need to */
810 add_j_to_nblist(vdw,jj0,bLR);
812 add_j_to_nblist(coul,jj0+1,bLR);
813 add_j_to_nblist(coul,jj0+2,bLR);
814 add_j_to_nblist(coul,jj0+3,bLR);
816 /* One entry for the entire water-water interaction */
819 add_j_to_nblist(coul_ww,jj0,bLR);
823 add_j_to_nblist(vdwc_ww,jj0,bLR);
830 /* j charge group is not water, but i is.
831 * Add entries to the water-other_atom lists; the geometry of the water
832 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
833 * so we don't care if it is SPC or TIP4P...
840 for(jj=jj0; (jj<jj1); jj++)
844 add_j_to_nblist(coul,jj,bLR);
850 for(jj=jj0; (jj<jj1); jj++)
852 if (bHaveVdW[type[jj]])
854 add_j_to_nblist(vdw,jj,bLR);
860 /* _charge_ _groups_ interact with both coulomb and LJ */
861 /* Check which atoms we should add to the lists! */
862 for(jj=jj0; (jj<jj1); jj++)
864 if (bHaveVdW[type[jj]])
868 add_j_to_nblist(vdwc,jj,bLR);
872 add_j_to_nblist(vdw,jj,bLR);
875 else if (charge[jj] != 0)
877 add_j_to_nblist(coul,jj,bLR);
884 close_i_nblist(coul);
885 close_i_nblist(vdwc);
886 #ifndef DISABLE_WATERWATER_NLIST
887 close_i_nblist(coul_ww);
888 close_i_nblist(vdwc_ww);
893 /* no solvent as i charge group */
894 /* Loop over the atoms in the i charge group */
895 for(i=0; i<nicg; i++)
898 gid = GID(igid,jgid,ngid);
901 /* Create new i_atom for each energy group */
902 if (bDoVdW && bDoCoul)
904 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
908 new_i_nblist(vdw,bLR,i_atom,shift,gid);
912 new_i_nblist(coul,bLR,i_atom,shift,gid);
914 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
915 bDoCoul_i = (bDoCoul && qi!=0);
917 if (bDoVdW_i || bDoCoul_i)
919 /* Loop over the j charge groups */
920 for(j=0; (j<nj); j++)
924 /* Check for large charge groups */
935 /* Finally loop over the atoms in the j-charge group */
936 for(jj=jj0; jj<jj1; jj++)
938 bNotEx = NOTEXCL(bExcl,i,jj);
946 add_j_to_nblist(coul,jj,bLR);
951 if (bHaveVdW[type[jj]])
953 add_j_to_nblist(vdw,jj,bLR);
958 if (bHaveVdW[type[jj]])
962 add_j_to_nblist(vdwc,jj,bLR);
966 add_j_to_nblist(vdw,jj,bLR);
969 else if (charge[jj] != 0)
971 add_j_to_nblist(coul,jj,bLR);
979 close_i_nblist(coul);
980 close_i_nblist(vdwc);
986 /* we are doing free energy */
987 vdwc_free = &nlist[eNL_VDWQQ_FREE];
988 vdw_free = &nlist[eNL_VDW_FREE];
989 coul_free = &nlist[eNL_QQ_FREE];
990 /* Loop over the atoms in the i charge group */
991 for(i=0; i<nicg; i++)
994 gid = GID(igid,jgid,ngid);
996 qiB = chargeB[i_atom];
998 /* Create new i_atom for each energy group */
999 if (bDoVdW && bDoCoul)
1000 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
1002 new_i_nblist(vdw,bLR,i_atom,shift,gid);
1004 new_i_nblist(coul,bLR,i_atom,shift,gid);
1006 new_i_nblist(vdw_free,bLR,i_atom,shift,gid);
1007 new_i_nblist(coul_free,bLR,i_atom,shift,gid);
1008 new_i_nblist(vdwc_free,bLR,i_atom,shift,gid);
1010 bDoVdW_i = (bDoVdW &&
1011 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
1012 bDoCoul_i = (bDoCoul && (qi!=0 || qiB!=0));
1013 /* For TIP4P the first atom does not have a charge,
1014 * but the last three do. So we should still put an atom
1015 * without LJ but with charge in the water-atom neighborlist
1016 * for a TIP4p i charge group.
1017 * For SPC type water the first atom has LJ and charge,
1018 * so there is no such problem.
1020 if (iwater == esolNO)
1022 bDoCoul_i_sol = bDoCoul_i;
1026 bDoCoul_i_sol = bDoCoul;
1029 if (bDoVdW_i || bDoCoul_i_sol)
1031 /* Loop over the j charge groups */
1032 for(j=0; (j<nj); j++)
1036 /* Check for large charge groups */
1047 /* Finally loop over the atoms in the j-charge group */
1048 bFree = bPert[i_atom];
1049 for(jj=jj0; (jj<jj1); jj++)
1051 bFreeJ = bFree || bPert[jj];
1052 /* Complicated if, because the water H's should also
1053 * see perturbed j-particles
1055 if (iwater==esolNO || i==0 || bFreeJ)
1057 bNotEx = NOTEXCL(bExcl,i,jj);
1065 if (charge[jj]!=0 || chargeB[jj]!=0)
1067 add_j_to_nblist(coul_free,jj,bLR);
1070 else if (!bDoCoul_i)
1072 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1074 add_j_to_nblist(vdw_free,jj,bLR);
1079 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1081 if (charge[jj]!=0 || chargeB[jj]!=0)
1083 add_j_to_nblist(vdwc_free,jj,bLR);
1087 add_j_to_nblist(vdw_free,jj,bLR);
1090 else if (charge[jj]!=0 || chargeB[jj]!=0)
1091 add_j_to_nblist(coul_free,jj,bLR);
1096 /* This is done whether or not bWater is set */
1097 if (charge[jj] != 0)
1099 add_j_to_nblist(coul,jj,bLR);
1102 else if (!bDoCoul_i_sol)
1104 if (bHaveVdW[type[jj]])
1106 add_j_to_nblist(vdw,jj,bLR);
1111 if (bHaveVdW[type[jj]])
1113 if (charge[jj] != 0)
1115 add_j_to_nblist(vdwc,jj,bLR);
1119 add_j_to_nblist(vdw,jj,bLR);
1122 else if (charge[jj] != 0)
1124 add_j_to_nblist(coul,jj,bLR);
1132 close_i_nblist(vdw);
1133 close_i_nblist(coul);
1134 close_i_nblist(vdwc);
1135 close_i_nblist(vdw_free);
1136 close_i_nblist(coul_free);
1137 close_i_nblist(vdwc_free);
1143 put_in_list_qmmm(gmx_bool bHaveVdW[],
1159 int i,j,jcg,igid,gid;
1160 atom_id jj,jj0,jj1,i_atom;
1164 /* Get atom range */
1166 nicg = index[icg+1]-i0;
1168 /* Get the i charge group info */
1169 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1171 coul = &fr->QMMMlist;
1173 /* Loop over atoms in the ith charge group */
1174 for (i=0;i<nicg;i++)
1177 gid = GID(igid,jgid,ngid);
1178 /* Create new i_atom for each energy group */
1179 new_i_nblist(coul,bLR,i_atom,shift,gid);
1181 /* Loop over the j charge groups */
1186 /* Charge groups cannot have QM and MM atoms simultaneously */
1191 /* Finally loop over the atoms in the j-charge group */
1192 for(jj=jj0; jj<jj1; jj++)
1194 bNotEx = NOTEXCL(bExcl,i,jj);
1196 add_j_to_nblist(coul,jj,bLR);
1200 close_i_nblist(coul);
1205 put_in_list_cg(gmx_bool bHaveVdW[],
1221 int igid,gid,nbl_ind;
1225 cginfo = fr->cginfo[icg];
1227 igid = GET_CGINFO_GID(cginfo);
1228 gid = GID(igid,jgid,ngid);
1230 /* Unpack pointers to neighbourlist structs */
1231 if (fr->nnblists == 1)
1237 nbl_ind = fr->gid2nblists[gid];
1241 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1245 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1248 /* Make a new neighbor list for charge group icg.
1249 * Currently simply one neighbor list is made with LJ and Coulomb.
1250 * If required, zero interactions could be removed here
1251 * or in the force loop.
1253 new_i_nblist(vdwc,bLR,index[icg],shift,gid);
1254 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1256 for(j=0; (j<nj); j++)
1259 /* Skip the icg-icg pairs if all self interactions are excluded */
1260 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1262 /* Here we add the j charge group jcg to the list,
1263 * exclusions are also added to the list.
1265 add_j_to_nblist_cg(vdwc,index[jcg],index[jcg+1],bExcl,icg==jcg,bLR);
1269 close_i_nblist(vdwc);
1272 static void setexcl(atom_id start,atom_id end,t_blocka *excl,gmx_bool b,
1279 for(i=start; i<end; i++)
1281 for(k=excl->index[i]; k<excl->index[i+1]; k++)
1283 SETEXCL(bexcl,i-start,excl->a[k]);
1289 for(i=start; i<end; i++)
1291 for(k=excl->index[i]; k<excl->index[i+1]; k++)
1293 RMEXCL(bexcl,i-start,excl->a[k]);
1299 int calc_naaj(int icg,int cgtot)
1303 if ((cgtot % 2) == 1)
1305 /* Odd number of charge groups, easy */
1306 naaj = 1 + (cgtot/2);
1308 else if ((cgtot % 4) == 0)
1310 /* Multiple of four is hard */
1347 fprintf(log,"naaj=%d\n",naaj);
1353 /************************************************
1355 * S I M P L E C O R E S T U F F
1357 ************************************************/
1359 static real calc_image_tric(rvec xi,rvec xj,matrix box,
1360 rvec b_inv,int *shift)
1362 /* This code assumes that the cut-off is smaller than
1363 * a half times the smallest diagonal element of the box.
1370 /* Compute diff vector */
1371 dz = xj[ZZ] - xi[ZZ];
1372 dy = xj[YY] - xi[YY];
1373 dx = xj[XX] - xi[XX];
1375 /* Perform NINT operation, using trunc operation, therefore
1376 * we first add 2.5 then subtract 2 again
1378 tz = dz*b_inv[ZZ] + h25;
1380 dz -= tz*box[ZZ][ZZ];
1381 dy -= tz*box[ZZ][YY];
1382 dx -= tz*box[ZZ][XX];
1384 ty = dy*b_inv[YY] + h25;
1386 dy -= ty*box[YY][YY];
1387 dx -= ty*box[YY][XX];
1389 tx = dx*b_inv[XX]+h25;
1391 dx -= tx*box[XX][XX];
1393 /* Distance squared */
1394 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1396 *shift = XYZ2IS(tx,ty,tz);
1401 static real calc_image_rect(rvec xi,rvec xj,rvec box_size,
1402 rvec b_inv,int *shift)
1410 /* Compute diff vector */
1411 dx = xj[XX] - xi[XX];
1412 dy = xj[YY] - xi[YY];
1413 dz = xj[ZZ] - xi[ZZ];
1415 /* Perform NINT operation, using trunc operation, therefore
1416 * we first add 1.5 then subtract 1 again
1418 tx = dx*b_inv[XX] + h15;
1419 ty = dy*b_inv[YY] + h15;
1420 tz = dz*b_inv[ZZ] + h15;
1425 /* Correct diff vector for translation */
1426 ddx = tx*box_size[XX] - dx;
1427 ddy = ty*box_size[YY] - dy;
1428 ddz = tz*box_size[ZZ] - dz;
1430 /* Distance squared */
1431 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1433 *shift = XYZ2IS(tx,ty,tz);
1438 static void add_simple(t_ns_buf *nsbuf,int nrj,atom_id cg_j,
1439 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1440 int icg,int jgid,t_block *cgs,t_excl bexcl[],
1441 int shift,t_forcerec *fr,put_in_list_t *put_in_list)
1443 if (nsbuf->nj + nrj > MAX_CG)
1445 put_in_list(bHaveVdW,ngid,md,icg,jgid,nsbuf->ncg,nsbuf->jcg,
1446 cgs->index,bexcl,shift,fr,FALSE,TRUE,TRUE);
1447 /* Reset buffer contents */
1448 nsbuf->ncg = nsbuf->nj = 0;
1450 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1454 static void ns_inner_tric(rvec x[],int icg,int *i_egp_flags,
1455 int njcg,atom_id jcg[],
1456 matrix box,rvec b_inv,real rcut2,
1457 t_block *cgs,t_ns_buf **ns_buf,
1458 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1459 t_excl bexcl[],t_forcerec *fr,
1460 put_in_list_t *put_in_list)
1464 int *cginfo=fr->cginfo;
1465 atom_id cg_j,*cgindex;
1468 cgindex = cgs->index;
1470 for(j=0; (j<njcg); j++)
1473 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1474 if (calc_image_tric(x[icg],x[cg_j],box,b_inv,&shift) < rcut2)
1476 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1477 if (!(i_egp_flags[jgid] & EGP_EXCL))
1479 add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1480 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1487 static void ns_inner_rect(rvec x[],int icg,int *i_egp_flags,
1488 int njcg,atom_id jcg[],
1489 gmx_bool bBox,rvec box_size,rvec b_inv,real rcut2,
1490 t_block *cgs,t_ns_buf **ns_buf,
1491 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1492 t_excl bexcl[],t_forcerec *fr,
1493 put_in_list_t *put_in_list)
1497 int *cginfo=fr->cginfo;
1498 atom_id cg_j,*cgindex;
1501 cgindex = cgs->index;
1505 for(j=0; (j<njcg); j++)
1508 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1509 if (calc_image_rect(x[icg],x[cg_j],box_size,b_inv,&shift) < rcut2)
1511 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1512 if (!(i_egp_flags[jgid] & EGP_EXCL))
1514 add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1515 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1523 for(j=0; (j<njcg); j++)
1526 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1527 if ((rcut2 == 0) || (distance2(x[icg],x[cg_j]) < rcut2)) {
1528 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1529 if (!(i_egp_flags[jgid] & EGP_EXCL))
1531 add_simple(&ns_buf[jgid][CENTRAL],nrj,cg_j,
1532 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,CENTRAL,fr,
1540 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1542 static int ns_simple_core(t_forcerec *fr,
1543 gmx_localtop_t *top,
1545 matrix box,rvec box_size,
1546 t_excl bexcl[],atom_id *aaj,
1547 int ngid,t_ns_buf **ns_buf,
1548 put_in_list_t *put_in_list,gmx_bool bHaveVdW[])
1552 int nsearch,icg,jcg,igid,i0,nri,nn;
1555 /* atom_id *i_atoms; */
1556 t_block *cgs=&(top->cgs);
1557 t_blocka *excl=&(top->excls);
1560 gmx_bool bBox,bTriclinic;
1563 rlist2 = sqr(fr->rlist);
1565 bBox = (fr->ePBC != epbcNONE);
1568 for(m=0; (m<DIM); m++)
1570 b_inv[m] = divide_err(1.0,box_size[m]);
1572 bTriclinic = TRICLINIC(box);
1579 cginfo = fr->cginfo;
1582 for (icg=fr->cg0; (icg<fr->hcg); icg++)
1585 i0 = cgs->index[icg];
1586 nri = cgs->index[icg+1]-i0;
1587 i_atoms = &(cgs->a[i0]);
1588 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1589 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1591 igid = GET_CGINFO_GID(cginfo[icg]);
1592 i_egp_flags = fr->egp_flags + ngid*igid;
1593 setexcl(cgs->index[icg],cgs->index[icg+1],excl,TRUE,bexcl);
1595 naaj=calc_naaj(icg,cgs->nr);
1598 ns_inner_tric(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1599 box,b_inv,rlist2,cgs,ns_buf,
1600 bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1604 ns_inner_rect(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1605 bBox,box_size,b_inv,rlist2,cgs,ns_buf,
1606 bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1610 for(nn=0; (nn<ngid); nn++)
1612 for(k=0; (k<SHIFTS); k++)
1614 nsbuf = &(ns_buf[nn][k]);
1617 put_in_list(bHaveVdW,ngid,md,icg,nn,nsbuf->ncg,nsbuf->jcg,
1618 cgs->index,bexcl,k,fr,FALSE,TRUE,TRUE);
1619 nsbuf->ncg=nsbuf->nj=0;
1623 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1624 setexcl(cgs->index[icg],cgs->index[icg+1],excl,FALSE,bexcl);
1626 close_neighbor_list(fr,FALSE,-1,-1,FALSE);
1631 /************************************************
1633 * N S 5 G R I D S T U F F
1635 ************************************************/
1637 static inline void get_dx(int Nx,real gridx,real rc2,int xgi,real x,
1638 int *dx0,int *dx1,real *dcx2)
1666 for(i=xgi0; i>=0; i--)
1668 dcx = (i+1)*gridx-x;
1675 for(i=xgi1; i<Nx; i++)
1688 static inline void get_dx_dd(int Nx,real gridx,real rc2,int xgi,real x,
1689 int ncpddc,int shift_min,int shift_max,
1690 int *g0,int *g1,real *dcx2)
1693 int g_min,g_max,shift_home;
1726 g_min = (shift_min == shift_home ? 0 : ncpddc);
1727 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1734 else if (shift_max < 0)
1749 /* Check one grid cell down */
1750 dcx = ((*g0 - 1) + 1)*gridx - x;
1762 /* Check one grid cell up */
1763 dcx = (*g1 + 1)*gridx - x;
1775 #define sqr(x) ((x)*(x))
1776 #define calc_dx2(XI,YI,ZI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
1777 #define calc_cyl_dx2(XI,YI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
1778 /****************************************************
1780 * F A S T N E I G H B O R S E A R C H I N G
1782 * Optimized neighboursearching routine using grid
1783 * at least 1x1x1, see GROMACS manual
1785 ****************************************************/
1787 static void do_longrange(t_commrec *cr,gmx_localtop_t *top,t_forcerec *fr,
1788 int ngid,t_mdatoms *md,int icg,
1790 atom_id lr[],t_excl bexcl[],int shift,
1791 rvec x[],rvec box_size,t_nrnb *nrnb,
1792 real lambda,real *dvdlambda,
1793 gmx_grppairener_t *grppener,
1794 gmx_bool bDoVdW,gmx_bool bDoCoul,
1795 gmx_bool bEvaluateNow,put_in_list_t *put_in_list,
1796 gmx_bool bHaveVdW[],
1797 gmx_bool bDoForces,rvec *f)
1802 for(n=0; n<fr->nnblists; n++)
1804 for(i=0; (i<eNL_NR); i++)
1806 nl = &fr->nblists[n].nlist_lr[i];
1807 if ((nl->nri > nl->maxnri-32) || bEvaluateNow)
1809 close_neighbor_list(fr,TRUE,n,i,FALSE);
1810 /* Evaluate the energies and forces */
1811 do_nonbonded(cr,fr,x,f,md,NULL,
1812 grppener->ener[fr->bBHAM ? egBHAMLR : egLJLR],
1813 grppener->ener[egCOULLR],
1814 grppener->ener[egGB],box_size,
1815 nrnb,lambda,dvdlambda,n,i,
1816 GMX_DONB_LR | GMX_DONB_FORCES);
1818 reset_neighbor_list(fr,TRUE,n,i);
1825 /* Put the long range particles in a list */
1826 /* do_longrange is never called for QMMM */
1827 put_in_list(bHaveVdW,ngid,md,icg,jgid,nlr,lr,top->cgs.index,
1828 bexcl,shift,fr,TRUE,bDoVdW,bDoCoul);
1832 static void get_cutoff2(t_forcerec *fr,gmx_bool bDoLongRange,
1833 real *rvdw2,real *rcoul2,
1834 real *rs2,real *rm2,real *rl2)
1836 *rs2 = sqr(fr->rlist);
1837 if (bDoLongRange && fr->bTwinRange)
1839 /* The VdW and elec. LR cut-off's could be different,
1840 * so we can not simply set them to rlistlong.
1842 if (EVDW_MIGHT_BE_ZERO_AT_CUTOFF(fr->vdwtype) &&
1843 fr->rvdw > fr->rlist)
1845 *rvdw2 = sqr(fr->rlistlong);
1849 *rvdw2 = sqr(fr->rvdw);
1851 if (EEL_MIGHT_BE_ZERO_AT_CUTOFF(fr->eeltype) &&
1852 fr->rcoulomb > fr->rlist)
1854 *rcoul2 = sqr(fr->rlistlong);
1858 *rcoul2 = sqr(fr->rcoulomb);
1863 /* Workaround for a gcc -O3 or -ffast-math problem */
1867 *rm2 = min(*rvdw2,*rcoul2);
1868 *rl2 = max(*rvdw2,*rcoul2);
1871 static void init_nsgrid_lists(t_forcerec *fr,int ngid,gmx_ns_t *ns)
1873 real rvdw2,rcoul2,rs2,rm2,rl2;
1876 get_cutoff2(fr,TRUE,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
1878 /* Short range buffers */
1879 snew(ns->nl_sr,ngid);
1882 snew(ns->nlr_ljc,ngid);
1883 snew(ns->nlr_one,ngid);
1887 /* Long range VdW and Coul buffers */
1888 snew(ns->nl_lr_ljc,ngid);
1892 /* Long range VdW or Coul only buffers */
1893 snew(ns->nl_lr_one,ngid);
1895 for(j=0; (j<ngid); j++) {
1896 snew(ns->nl_sr[j],MAX_CG);
1899 snew(ns->nl_lr_ljc[j],MAX_CG);
1903 snew(ns->nl_lr_one[j],MAX_CG);
1909 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
1914 static int nsgrid_core(FILE *log,t_commrec *cr,t_forcerec *fr,
1915 matrix box,rvec box_size,int ngid,
1916 gmx_localtop_t *top,
1917 t_grid *grid,rvec x[],
1918 t_excl bexcl[],gmx_bool *bExcludeAlleg,
1919 t_nrnb *nrnb,t_mdatoms *md,
1920 real lambda,real *dvdlambda,
1921 gmx_grppairener_t *grppener,
1922 put_in_list_t *put_in_list,
1923 gmx_bool bHaveVdW[],
1924 gmx_bool bDoLongRange,gmx_bool bDoForces,rvec *f,
1925 gmx_bool bMakeQMMMnblist)
1928 atom_id **nl_lr_ljc,**nl_lr_one,**nl_sr;
1929 int *nlr_ljc,*nlr_one,*nsr;
1930 gmx_domdec_t *dd=NULL;
1931 t_block *cgs=&(top->cgs);
1932 int *cginfo=fr->cginfo;
1933 /* atom_id *i_atoms,*cgsindex=cgs->index; */
1935 int cell_x,cell_y,cell_z;
1936 int d,tx,ty,tz,dx,dy,dz,cj;
1937 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
1938 int zsh_ty,zsh_tx,ysh_tx;
1940 int dx0,dx1,dy0,dy1,dz0,dz1;
1941 int Nx,Ny,Nz,shift=-1,j,nrj,nns,nn=-1;
1942 real gridx,gridy,gridz,grid_x,grid_y,grid_z;
1943 real *dcx2,*dcy2,*dcz2;
1945 int cg0,cg1,icg=-1,cgsnr,i0,igid,nri,naaj,max_jcg;
1946 int jcg0,jcg1,jjcg,cgj0,jgid;
1947 int *grida,*gridnra,*gridind;
1948 gmx_bool rvdw_lt_rcoul,rcoul_lt_rvdw;
1949 rvec xi,*cgcm,grid_offset;
1950 real r2,rs2,rvdw2,rcoul2,rm2,rl2,XI,YI,ZI,dcx,dcy,dcz,tmp1,tmp2;
1952 gmx_bool bDomDec,bTriclinicX,bTriclinicY;
1957 bDomDec = DOMAINDECOMP(cr);
1963 bTriclinicX = ((YY < grid->npbcdim &&
1964 (!bDomDec || dd->nc[YY]==1) && box[YY][XX] != 0) ||
1965 (ZZ < grid->npbcdim &&
1966 (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][XX] != 0));
1967 bTriclinicY = (ZZ < grid->npbcdim &&
1968 (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][YY] != 0);
1972 get_cutoff2(fr,bDoLongRange,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
1974 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
1975 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
1977 if (bMakeQMMMnblist)
1985 nl_lr_ljc = ns->nl_lr_ljc;
1986 nl_lr_one = ns->nl_lr_one;
1987 nlr_ljc = ns->nlr_ljc;
1988 nlr_one = ns->nlr_one;
1996 gridind = grid->index;
1997 gridnra = grid->nra;
2000 gridx = grid->cell_size[XX];
2001 gridy = grid->cell_size[YY];
2002 gridz = grid->cell_size[ZZ];
2006 copy_rvec(grid->cell_offset,grid_offset);
2007 copy_ivec(grid->ncpddc,ncpddc);
2012 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2013 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
2014 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
2015 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
2016 if (zsh_tx!=0 && ysh_tx!=0)
2018 /* This could happen due to rounding, when both ratios are 0.5 */
2027 /* We only want a list for the test particle */
2036 /* Set the shift range */
2037 for(d=0; d<DIM; d++)
2041 /* Check if we need periodicity shifts.
2042 * Without PBC or with domain decomposition we don't need them.
2044 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2051 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
2062 /* Loop over charge groups */
2063 for(icg=cg0; (icg < cg1); icg++)
2065 igid = GET_CGINFO_GID(cginfo[icg]);
2066 /* Skip this charge group if all energy groups are excluded! */
2067 if (bExcludeAlleg[igid])
2072 i0 = cgs->index[icg];
2074 if (bMakeQMMMnblist)
2076 /* Skip this charge group if it is not a QM atom while making a
2077 * QM/MM neighbourlist
2079 if (md->bQM[i0]==FALSE)
2081 continue; /* MM particle, go to next particle */
2084 /* Compute the number of charge groups that fall within the control
2087 naaj = calc_naaj(icg,cgsnr);
2094 /* make a normal neighbourlist */
2098 /* Get the j charge-group and dd cell shift ranges */
2099 dd_get_ns_ranges(cr->dd,icg,&jcg0,&jcg1,sh0,sh1);
2104 /* Compute the number of charge groups that fall within the control
2107 naaj = calc_naaj(icg,cgsnr);
2113 /* The i-particle is awlways the test particle,
2114 * so we want all j-particles
2116 max_jcg = cgsnr - 1;
2120 max_jcg = jcg1 - cgsnr;
2125 i_egp_flags = fr->egp_flags + igid*ngid;
2127 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2128 setexcl(i0,cgs->index[icg+1],&top->excls,TRUE,bexcl);
2130 ci2xyz(grid,icg,&cell_x,&cell_y,&cell_z);
2132 /* Changed iicg to icg, DvdS 990115
2133 * (but see consistency check above, DvdS 990330)
2136 fprintf(log,"icg=%5d, naaj=%5d, cell %d %d %d\n",
2137 icg,naaj,cell_x,cell_y,cell_z);
2139 /* Loop over shift vectors in three dimensions */
2140 for (tz=-shp[ZZ]; tz<=shp[ZZ]; tz++)
2142 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2143 /* Calculate range of cells in Z direction that have the shift tz */
2144 zgi = cell_z + tz*Nz;
2147 get_dx(Nz,gridz,rl2,zgi,ZI,&dz0,&dz1,dcz2);
2149 get_dx_dd(Nz,gridz,rl2,zgi,ZI-grid_offset[ZZ],
2150 ncpddc[ZZ],sh0[ZZ],sh1[ZZ],&dz0,&dz1,dcz2);
2156 for (ty=-shp[YY]; ty<=shp[YY]; ty++)
2158 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2159 /* Calculate range of cells in Y direction that have the shift ty */
2162 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2166 ygi = cell_y + ty*Ny;
2169 get_dx(Ny,gridy,rl2,ygi,YI,&dy0,&dy1,dcy2);
2171 get_dx_dd(Ny,gridy,rl2,ygi,YI-grid_offset[YY],
2172 ncpddc[YY],sh0[YY],sh1[YY],&dy0,&dy1,dcy2);
2178 for (tx=-shp[XX]; tx<=shp[XX]; tx++)
2180 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2181 /* Calculate range of cells in X direction that have the shift tx */
2184 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2188 xgi = cell_x + tx*Nx;
2191 get_dx(Nx,gridx,rl2,xgi*Nx,XI,&dx0,&dx1,dcx2);
2193 get_dx_dd(Nx,gridx,rl2,xgi,XI-grid_offset[XX],
2194 ncpddc[XX],sh0[XX],sh1[XX],&dx0,&dx1,dcx2);
2200 /* Get shift vector */
2201 shift=XYZ2IS(tx,ty,tz);
2203 range_check(shift,0,SHIFTS);
2205 for(nn=0; (nn<ngid); nn++)
2212 fprintf(log,"shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2213 shift,dx0,dx1,dy0,dy1,dz0,dz1);
2214 fprintf(log,"cgcm: %8.3f %8.3f %8.3f\n",cgcm[icg][XX],
2215 cgcm[icg][YY],cgcm[icg][ZZ]);
2216 fprintf(log,"xi: %8.3f %8.3f %8.3f\n",XI,YI,ZI);
2218 for (dx=dx0; (dx<=dx1); dx++)
2220 tmp1 = rl2 - dcx2[dx];
2221 for (dy=dy0; (dy<=dy1); dy++)
2223 tmp2 = tmp1 - dcy2[dy];
2226 for (dz=dz0; (dz<=dz1); dz++) {
2227 if (tmp2 > dcz2[dz]) {
2228 /* Find grid-cell cj in which possible neighbours are */
2229 cj = xyz2ci(Ny,Nz,dx,dy,dz);
2231 /* Check out how many cgs (nrj) there in this cell */
2234 /* Find the offset in the cg list */
2237 /* Check if all j's are out of range so we
2238 * can skip the whole cell.
2239 * Should save some time, especially with DD.
2242 (grida[cgj0] >= max_jcg &&
2243 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2249 for (j=0; (j<nrj); j++)
2251 jjcg = grida[cgj0+j];
2253 /* check whether this guy is in range! */
2254 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2257 r2=calc_dx2(XI,YI,ZI,cgcm[jjcg]);
2259 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2260 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2261 /* check energy group exclusions */
2262 if (!(i_egp_flags[jgid] & EGP_EXCL))
2266 if (nsr[jgid] >= MAX_CG)
2268 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2269 nsr[jgid],nl_sr[jgid],
2270 cgs->index,/* cgsatoms, */ bexcl,
2271 shift,fr,FALSE,TRUE,TRUE);
2274 nl_sr[jgid][nsr[jgid]++]=jjcg;
2278 if (nlr_ljc[jgid] >= MAX_CG)
2280 do_longrange(cr,top,fr,ngid,md,icg,jgid,
2282 nl_lr_ljc[jgid],bexcl,shift,x,
2292 nl_lr_ljc[jgid][nlr_ljc[jgid]++]=jjcg;
2296 if (nlr_one[jgid] >= MAX_CG) {
2297 do_longrange(cr,top,fr,ngid,md,icg,jgid,
2299 nl_lr_one[jgid],bexcl,shift,x,
2303 rvdw_lt_rcoul,rcoul_lt_rvdw,FALSE,
2309 nl_lr_one[jgid][nlr_one[jgid]++]=jjcg;
2321 /* CHECK whether there is anything left in the buffers */
2322 for(nn=0; (nn<ngid); nn++)
2326 put_in_list(bHaveVdW,ngid,md,icg,nn,nsr[nn],nl_sr[nn],
2327 cgs->index, /* cgsatoms, */ bexcl,
2328 shift,fr,FALSE,TRUE,TRUE);
2331 if (nlr_ljc[nn] > 0)
2333 do_longrange(cr,top,fr,ngid,md,icg,nn,nlr_ljc[nn],
2334 nl_lr_ljc[nn],bexcl,shift,x,box_size,nrnb,
2335 lambda,dvdlambda,grppener,TRUE,TRUE,FALSE,
2336 put_in_list,bHaveVdW,bDoForces,f);
2339 if (nlr_one[nn] > 0)
2341 do_longrange(cr,top,fr,ngid,md,icg,nn,nlr_one[nn],
2342 nl_lr_one[nn],bexcl,shift,x,box_size,nrnb,
2343 lambda,dvdlambda,grppener,
2344 rvdw_lt_rcoul,rcoul_lt_rvdw,FALSE,
2345 put_in_list,bHaveVdW,bDoForces,f);
2351 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2352 setexcl(cgs->index[icg],cgs->index[icg+1],&top->excls,FALSE,bexcl);
2354 /* Perform any left over force calculations */
2355 for (nn=0; (nn<ngid); nn++)
2359 do_longrange(cr,top,fr,0,md,icg,nn,nlr_ljc[nn],
2360 nl_lr_ljc[nn],bexcl,shift,x,box_size,nrnb,
2361 lambda,dvdlambda,grppener,
2362 TRUE,TRUE,TRUE,put_in_list,bHaveVdW,bDoForces,f);
2365 do_longrange(cr,top,fr,0,md,icg,nn,nlr_one[nn],
2366 nl_lr_one[nn],bexcl,shift,x,box_size,nrnb,
2367 lambda,dvdlambda,grppener,
2368 rvdw_lt_rcoul,rcoul_lt_rvdw,
2369 TRUE,put_in_list,bHaveVdW,bDoForces,f);
2374 /* Close off short range neighbourlists */
2375 close_neighbor_list(fr,FALSE,-1,-1,bMakeQMMMnblist);
2380 void ns_realloc_natoms(gmx_ns_t *ns,int natoms)
2384 if (natoms > ns->nra_alloc)
2386 ns->nra_alloc = over_alloc_dd(natoms);
2387 srenew(ns->bexcl,ns->nra_alloc);
2388 for(i=0; i<ns->nra_alloc; i++)
2395 void init_ns(FILE *fplog,const t_commrec *cr,
2396 gmx_ns_t *ns,t_forcerec *fr,
2397 const gmx_mtop_t *mtop,
2400 int mt,icg,nr_in_cg,maxcg,i,j,jcg,ngid,ncg;
2404 /* Compute largest charge groups size (# atoms) */
2406 for(mt=0; mt<mtop->nmoltype; mt++) {
2407 cgs = &mtop->moltype[mt].cgs;
2408 for (icg=0; (icg < cgs->nr); icg++)
2410 nr_in_cg=max(nr_in_cg,(int)(cgs->index[icg+1]-cgs->index[icg]));
2414 /* Verify whether largest charge group is <= max cg.
2415 * This is determined by the type of the local exclusion type
2416 * Exclusions are stored in bits. (If the type is not large
2417 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2419 maxcg = sizeof(t_excl)*8;
2420 if (nr_in_cg > maxcg)
2422 gmx_fatal(FARGS,"Max #atoms in a charge group: %d > %d\n",
2426 ngid = mtop->groups.grps[egcENER].nr;
2427 snew(ns->bExcludeAlleg,ngid);
2428 for(i=0; i<ngid; i++) {
2429 ns->bExcludeAlleg[i] = TRUE;
2430 for(j=0; j<ngid; j++)
2432 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2434 ns->bExcludeAlleg[i] = FALSE;
2441 ns->grid = init_grid(fplog,fr);
2442 init_nsgrid_lists(fr,ngid,ns);
2447 snew(ns->ns_buf,ngid);
2448 for(i=0; (i<ngid); i++)
2450 snew(ns->ns_buf[i],SHIFTS);
2452 ncg = ncg_mtop(mtop);
2453 snew(ns->simple_aaj,2*ncg);
2454 for(jcg=0; (jcg<ncg); jcg++)
2456 ns->simple_aaj[jcg] = jcg;
2457 ns->simple_aaj[jcg+ncg] = jcg;
2461 /* Create array that determines whether or not atoms have VdW */
2462 snew(ns->bHaveVdW,fr->ntype);
2463 for(i=0; (i<fr->ntype); i++)
2465 for(j=0; (j<fr->ntype); j++)
2467 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2469 ((BHAMA(fr->nbfp,fr->ntype,i,j) != 0) ||
2470 (BHAMB(fr->nbfp,fr->ntype,i,j) != 0) ||
2471 (BHAMC(fr->nbfp,fr->ntype,i,j) != 0)) :
2472 ((C6(fr->nbfp,fr->ntype,i,j) != 0) ||
2473 (C12(fr->nbfp,fr->ntype,i,j) != 0))));
2477 pr_bvec(debug,0,"bHaveVdW",ns->bHaveVdW,fr->ntype,TRUE);
2481 if (!DOMAINDECOMP(cr))
2483 /* This could be reduced with particle decomposition */
2484 ns_realloc_natoms(ns,mtop->natoms);
2487 ns->nblist_initialized=FALSE;
2489 /* nbr list debug dump */
2491 char *ptr=getenv("GMX_DUMP_NL");
2494 ns->dump_nl=strtol(ptr,NULL,10);
2497 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2508 int search_neighbours(FILE *log,t_forcerec *fr,
2509 rvec x[],matrix box,
2510 gmx_localtop_t *top,
2511 gmx_groups_t *groups,
2513 t_nrnb *nrnb,t_mdatoms *md,
2514 real lambda,real *dvdlambda,
2515 gmx_grppairener_t *grppener,
2517 gmx_bool bDoLongRange,
2518 gmx_bool bDoForces,rvec *f)
2520 t_block *cgs=&(top->cgs);
2521 rvec box_size,grid_x0,grid_x1;
2523 real min_size,grid_dens;
2527 gmx_bool *i_egp_flags;
2528 int cg_start,cg_end,start,end;
2531 gmx_domdec_zones_t *dd_zones;
2532 put_in_list_t *put_in_list;
2536 /* Set some local variables */
2538 ngid = groups->grps[egcENER].nr;
2540 for(m=0; (m<DIM); m++)
2542 box_size[m] = box[m][m];
2545 if (fr->ePBC != epbcNONE)
2547 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC,box))
2549 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.");
2553 min_size = min(box_size[XX],min(box_size[YY],box_size[ZZ]));
2554 if (2*fr->rlistlong >= min_size)
2555 gmx_fatal(FARGS,"One of the box diagonal elements has become smaller than twice the cut-off length.");
2559 if (DOMAINDECOMP(cr))
2561 ns_realloc_natoms(ns,cgs->index[cgs->nr]);
2565 /* Reset the neighbourlists */
2566 reset_neighbor_list(fr,FALSE,-1,-1);
2568 if (bGrid && bFillGrid)
2572 if (DOMAINDECOMP(cr))
2574 dd_zones = domdec_zones(cr->dd);
2580 get_nsgrid_boundaries(grid,NULL,box,NULL,NULL,NULL,
2581 cgs->nr,fr->cg_cm,grid_x0,grid_x1,&grid_dens);
2583 grid_first(log,grid,NULL,NULL,fr->ePBC,box,grid_x0,grid_x1,
2584 fr->rlistlong,grid_dens);
2588 /* Don't know why this all is... (DvdS 3/99) */
2594 end = (cgs->nr+1)/2;
2597 if (DOMAINDECOMP(cr))
2600 fill_grid(log,dd_zones,grid,end,-1,end,fr->cg_cm);
2602 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2606 fill_grid(log,NULL,grid,cgs->nr,fr->cg0,fr->hcg,fr->cg_cm);
2607 grid->icg0 = fr->cg0;
2608 grid->icg1 = fr->hcg;
2616 calc_elemnr(log,grid,start,end,cgs->nr);
2618 grid_last(log,grid,start,end,cgs->nr);
2622 check_grid(debug,grid);
2623 print_grid(debug,grid);
2628 /* Set the grid cell index for the test particle only.
2629 * The cell to cg index is not corrected, but that does not matter.
2631 fill_grid(log,NULL,ns->grid,fr->hcg,fr->hcg-1,fr->hcg,fr->cg_cm);
2635 if (!fr->ns.bCGlist)
2637 put_in_list = put_in_list_at;
2641 put_in_list = put_in_list_cg;
2648 nsearch = nsgrid_core(log,cr,fr,box,box_size,ngid,top,
2649 grid,x,ns->bexcl,ns->bExcludeAlleg,
2650 nrnb,md,lambda,dvdlambda,grppener,
2651 put_in_list,ns->bHaveVdW,
2652 bDoLongRange,bDoForces,f,
2655 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2656 * the classical calculation. The charge-charge interaction
2657 * between QM and MM atoms is handled in the QMMM core calculation
2658 * (see QMMM.c). The VDW however, we'd like to compute classically
2659 * and the QM MM atom pairs have just been put in the
2660 * corresponding neighbourlists. in case of QMMM we still need to
2661 * fill a special QMMM neighbourlist that contains all neighbours
2662 * of the QM atoms. If bQMMM is true, this list will now be made:
2664 if (fr->bQMMM && fr->qr->QMMMscheme!=eQMMMschemeoniom)
2666 nsearch += nsgrid_core(log,cr,fr,box,box_size,ngid,top,
2667 grid,x,ns->bexcl,ns->bExcludeAlleg,
2668 nrnb,md,lambda,dvdlambda,grppener,
2669 put_in_list_qmmm,ns->bHaveVdW,
2670 bDoLongRange,bDoForces,f,
2676 nsearch = ns_simple_core(fr,top,md,box,box_size,
2677 ns->bexcl,ns->simple_aaj,
2678 ngid,ns->ns_buf,put_in_list,ns->bHaveVdW);
2686 inc_nrnb(nrnb,eNR_NS,nsearch);
2687 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2692 int natoms_beyond_ns_buffer(t_inputrec *ir,t_forcerec *fr,t_block *cgs,
2693 matrix scale_tot,rvec *x)
2695 int cg0,cg1,cg,a0,a1,a,i,j;
2696 real rint,hbuf2,scale;
2698 gmx_bool bIsotropic;
2703 rint = max(ir->rcoulomb,ir->rvdw);
2704 if (ir->rlist < rint)
2706 gmx_fatal(FARGS,"The neighbor search buffer has negative size: %f nm",
2714 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2716 hbuf2 = sqr(0.5*(ir->rlist - rint));
2717 for(cg=cg0; cg<cg1; cg++)
2719 a0 = cgs->index[cg];
2720 a1 = cgs->index[cg+1];
2721 for(a=a0; a<a1; a++)
2723 if (distance2(cg_cm[cg],x[a]) > hbuf2)
2733 scale = scale_tot[0][0];
2734 for(i=1; i<DIM; i++)
2736 /* With anisotropic scaling, the original spherical ns volumes become
2737 * ellipsoids. To avoid costly transformations we use the minimum
2738 * eigenvalue of the scaling matrix for determining the buffer size.
2739 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2741 scale = min(scale,scale_tot[i][i]);
2742 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2748 if (scale_tot[i][j] != 0)
2754 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2757 for(cg=cg0; cg<cg1; cg++)
2759 svmul(scale,cg_cm[cg],cgsc);
2760 a0 = cgs->index[cg];
2761 a1 = cgs->index[cg+1];
2762 for(a=a0; a<a1; a++)
2764 if (distance2(cgsc,x[a]) > hbuf2)
2773 /* Anistropic scaling */
2774 for(cg=cg0; cg<cg1; cg++)
2776 /* Since scale_tot contains the transpose of the scaling matrix,
2777 * we need to multiply with the transpose.
2779 tmvmul_ur0(scale_tot,cg_cm[cg],cgsc);
2780 a0 = cgs->index[cg];
2781 a1 = cgs->index[cg+1];
2782 for(a=a0; a<a1; a++)
2784 if (distance2(cgsc,x[a]) > hbuf2)