1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
40 #ifdef GMX_THREAD_SHM_FDECOMP
54 #include "nonbonded.h"
58 #include "gmx_fatal.h"
61 #include "mtop_util.h"
68 * E X C L U S I O N H A N D L I N G
72 static void SETEXCL_(t_excl e[],atom_id i,atom_id j)
73 { e[j] = e[j] | (1<<i); }
74 static void RMEXCL_(t_excl e[],atom_id i,atom_id j)
75 { e[j]=e[j] & ~(1<<i); }
76 static gmx_bool ISEXCL_(t_excl e[],atom_id i,atom_id j)
77 { return (gmx_bool)(e[j] & (1<<i)); }
78 static gmx_bool NOTEXCL_(t_excl e[],atom_id i,atom_id j)
79 { return !(ISEXCL(e,i,j)); }
81 #define SETEXCL(e,i,j) (e)[((atom_id) (j))] |= (1<<((atom_id) (i)))
82 #define RMEXCL(e,i,j) (e)[((atom_id) (j))] &= (~(1<<((atom_id) (i))))
83 #define ISEXCL(e,i,j) (gmx_bool) ((e)[((atom_id) (j))] & (1<<((atom_id) (i))))
84 #define NOTEXCL(e,i,j) !(ISEXCL(e,i,j))
87 /************************************************
89 * U T I L I T I E S F O R N S
91 ************************************************/
93 static void reallocate_nblist(t_nblist *nl)
97 fprintf(debug,"reallocating neigborlist il_code=%d, maxnri=%d\n",
98 nl->il_code,nl->maxnri);
100 srenew(nl->iinr, nl->maxnri);
101 if (nl->enlist == enlistCG_CG)
103 srenew(nl->iinr_end,nl->maxnri);
105 srenew(nl->gid, nl->maxnri);
106 srenew(nl->shift, nl->maxnri);
107 srenew(nl->jindex, nl->maxnri+1);
110 /* ivdw/icoul are used to determine the type of interaction, so we
111 * can set an innerloop index here. The obvious choice for this would have
112 * been the vdwtype/coultype values in the forcerecord, but unfortunately
113 * those types are braindead - for instance both Buckingham and normal
114 * Lennard-Jones use the same value (evdwCUT), and a separate gmx_boolean variable
115 * to determine which interaction is used. There is further no special value
116 * for 'no interaction'. For backward compatibility with old TPR files we won't
117 * change this in the 3.x series, so when calling this routine you should use:
119 * icoul=0 no coulomb interaction
120 * icoul=1 cutoff standard coulomb
121 * icoul=2 reaction-field coulomb
122 * icoul=3 tabulated coulomb
124 * ivdw=0 no vdw interaction
125 * ivdw=1 standard L-J interaction
127 * ivdw=3 tabulated vdw.
129 * Kind of ugly, but it works.
131 static void init_nblist(t_nblist *nl_sr,t_nblist *nl_lr,
134 gmx_bool bfree, int enlist)
166 nl = (i == 0) ? nl_sr : nl_lr;
167 homenr = (i == 0) ? maxsr : maxlr;
174 /* Set coul/vdw in neighborlist, and for the normal loops we determine
175 * an index of which one to call.
179 nl->free_energy = bfree;
183 nl->enlist = enlistATOM_ATOM;
184 nl->il_code = eNR_NBKERNEL_FREE_ENERGY;
190 nn = inloop[4*icoul + ivdw];
192 /* solvent loops follow directly after the corresponding
193 * ordinary loops, in the order:
195 * SPC, SPC-SPC, TIP4p, TIP4p-TIP4p
199 case enlistATOM_ATOM:
202 case enlistSPC_ATOM: nn += 1; break;
203 case enlistSPC_SPC: nn += 2; break;
204 case enlistTIP4P_ATOM: nn += 3; break;
205 case enlistTIP4P_TIP4P: nn += 4; break;
212 fprintf(debug,"Initiating neighbourlist type %d for %s interactions,\nwith %d SR, %d LR atoms.\n",
213 nl->il_code,ENLISTTYPE(enlist),maxsr,maxlr);
215 /* maxnri is influenced by the number of shifts (maximum is 8)
216 * and the number of energy groups.
217 * If it is not enough, nl memory will be reallocated during the run.
218 * 4 seems to be a reasonable factor, which only causes reallocation
219 * during runs with tiny and many energygroups.
221 nl->maxnri = homenr*4;
230 reallocate_nblist(nl);
232 #ifdef GMX_THREAD_SHM_FDECOMP
235 pthread_mutex_init(nl->mtx,NULL);
240 void init_neighbor_list(FILE *log,t_forcerec *fr,int homenr)
242 /* Make maxlr tunable! (does not seem to be a big difference though)
243 * This parameter determines the number of i particles in a long range
244 * neighbourlist. Too few means many function calls, too many means
247 int maxsr,maxsr_wat,maxlr,maxlr_wat;
248 int icoul,icoulf,ivdw;
250 int enlist_def,enlist_w,enlist_ww;
254 /* maxsr = homenr-fr->nWatMol*3; */
259 gmx_fatal(FARGS,"%s, %d: Negative number of short range atoms.\n"
260 "Call your Gromacs dealer for assistance.",__FILE__,__LINE__);
262 /* This is just for initial allocation, so we do not reallocate
263 * all the nlist arrays many times in a row.
264 * The numbers seem very accurate, but they are uncritical.
266 maxsr_wat = min(fr->nWatMol,(homenr+2)/3);
270 maxlr_wat = min(maxsr_wat,maxlr);
274 maxlr = maxlr_wat = 0;
277 /* Determine the values for icoul/ivdw. */
283 else if (fr->bcoultab)
287 else if (EEL_RF(fr->eeltype))
309 fr->ns.bCGlist = (getenv("GMX_NBLISTCG") != 0);
312 enlist_def = enlistATOM_ATOM;
316 enlist_def = enlistCG_CG;
319 fprintf(log,"\nUsing charge-group - charge-group neighbor lists and kernels\n\n");
323 if (fr->solvent_opt == esolTIP4P) {
324 enlist_w = enlistTIP4P_ATOM;
325 enlist_ww = enlistTIP4P_TIP4P;
327 enlist_w = enlistSPC_ATOM;
328 enlist_ww = enlistSPC_SPC;
331 for(i=0; i<fr->nnblists; i++)
333 nbl = &(fr->nblists[i]);
334 init_nblist(&nbl->nlist_sr[eNL_VDWQQ],&nbl->nlist_lr[eNL_VDWQQ],
335 maxsr,maxlr,ivdw,icoul,FALSE,enlist_def);
336 init_nblist(&nbl->nlist_sr[eNL_VDW],&nbl->nlist_lr[eNL_VDW],
337 maxsr,maxlr,ivdw,0,FALSE,enlist_def);
338 init_nblist(&nbl->nlist_sr[eNL_QQ],&nbl->nlist_lr[eNL_QQ],
339 maxsr,maxlr,0,icoul,FALSE,enlist_def);
340 init_nblist(&nbl->nlist_sr[eNL_VDWQQ_WATER],&nbl->nlist_lr[eNL_VDWQQ_WATER],
341 maxsr_wat,maxlr_wat,ivdw,icoul, FALSE,enlist_w);
342 init_nblist(&nbl->nlist_sr[eNL_QQ_WATER],&nbl->nlist_lr[eNL_QQ_WATER],
343 maxsr_wat,maxlr_wat,0,icoul, FALSE,enlist_w);
344 init_nblist(&nbl->nlist_sr[eNL_VDWQQ_WATERWATER],&nbl->nlist_lr[eNL_VDWQQ_WATERWATER],
345 maxsr_wat,maxlr_wat,ivdw,icoul, FALSE,enlist_ww);
346 init_nblist(&nbl->nlist_sr[eNL_QQ_WATERWATER],&nbl->nlist_lr[eNL_QQ_WATERWATER],
347 maxsr_wat,maxlr_wat,0,icoul, FALSE,enlist_ww);
349 if (fr->efep != efepNO)
360 init_nblist(&nbl->nlist_sr[eNL_VDWQQ_FREE],&nbl->nlist_lr[eNL_VDWQQ_FREE],
361 maxsr,maxlr,ivdw,icoulf,TRUE,enlistATOM_ATOM);
362 init_nblist(&nbl->nlist_sr[eNL_VDW_FREE],&nbl->nlist_lr[eNL_VDW_FREE],
363 maxsr,maxlr,ivdw,0,TRUE,enlistATOM_ATOM);
364 init_nblist(&nbl->nlist_sr[eNL_QQ_FREE],&nbl->nlist_lr[eNL_QQ_FREE],
365 maxsr,maxlr,0,icoulf,TRUE,enlistATOM_ATOM);
369 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
371 init_nblist(&fr->QMMMlist,NULL,
372 maxsr,maxlr,0,icoul,FALSE,enlistATOM_ATOM);
375 fr->ns.nblist_initialized=TRUE;
378 static void reset_nblist(t_nblist *nl)
389 static void reset_neighbor_list(t_forcerec *fr,gmx_bool bLR,int nls,int eNL)
395 reset_nblist(&(fr->nblists[nls].nlist_lr[eNL]));
399 for(n=0; n<fr->nnblists; n++)
401 for(i=0; i<eNL_NR; i++)
403 reset_nblist(&(fr->nblists[n].nlist_sr[i]));
408 /* only reset the short-range nblist */
409 reset_nblist(&(fr->QMMMlist));
417 static inline void new_i_nblist(t_nblist *nlist,
418 gmx_bool bLR,atom_id i_atom,int shift,int gid)
424 /* Check whether we have to increase the i counter */
426 (nlist->iinr[nri] != i_atom) ||
427 (nlist->shift[nri] != shift) ||
428 (nlist->gid[nri] != gid))
430 /* This is something else. Now see if any entries have
431 * been added in the list of the previous atom.
434 ((nlist->jindex[nri+1] > nlist->jindex[nri]) &&
435 (nlist->gid[nri] != -1)))
437 /* If so increase the counter */
440 if (nlist->nri >= nlist->maxnri)
442 nlist->maxnri += over_alloc_large(nlist->nri);
443 reallocate_nblist(nlist);
446 /* Set the number of neighbours and the atom number */
447 nlist->jindex[nri+1] = nlist->jindex[nri];
448 nlist->iinr[nri] = i_atom;
449 nlist->gid[nri] = gid;
450 nlist->shift[nri] = shift;
454 static inline void close_i_nblist(t_nblist *nlist)
456 int nri = nlist->nri;
461 nlist->jindex[nri+1] = nlist->nrj;
463 len=nlist->nrj - nlist->jindex[nri];
465 /* nlist length for water i molecules is treated statically
468 if (len > nlist->maxlen)
475 static inline void close_nblist(t_nblist *nlist)
477 /* Only close this nblist when it has been initialized.
478 * Avoid the creation of i-lists with no j-particles.
482 /* Some assembly kernels do not support empty lists,
483 * make sure here that we don't generate any empty lists.
484 * With the current ns code this branch is taken in two cases:
485 * No i-particles at all: nri=-1 here
486 * There are i-particles, but no j-particles; nri=0 here
492 /* Close list number nri by incrementing the count */
497 static inline void close_neighbor_list(t_forcerec *fr,gmx_bool bLR,int nls,int eNL,
498 gmx_bool bMakeQMMMnblist)
502 if (bMakeQMMMnblist) {
505 close_nblist(&(fr->QMMMlist));
512 close_nblist(&(fr->nblists[nls].nlist_lr[eNL]));
516 for(n=0; n<fr->nnblists; n++)
518 for(i=0; (i<eNL_NR); i++)
520 close_nblist(&(fr->nblists[n].nlist_sr[i]));
527 static inline void add_j_to_nblist(t_nblist *nlist,atom_id j_atom,gmx_bool bLR)
531 if (nlist->nrj >= nlist->maxnrj)
533 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
535 fprintf(debug,"Increasing %s nblist %s j size to %d\n",
536 bLR ? "LR" : "SR",nrnb_str(nlist->il_code),nlist->maxnrj);
538 srenew(nlist->jjnr,nlist->maxnrj);
541 nlist->jjnr[nrj] = j_atom;
545 static inline void add_j_to_nblist_cg(t_nblist *nlist,
546 atom_id j_start,int j_end,
547 t_excl *bexcl,gmx_bool i_is_j,
553 if (nlist->nrj >= nlist->maxnrj)
555 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
557 fprintf(debug,"Increasing %s nblist %s j size to %d\n",
558 bLR ? "LR" : "SR",nrnb_str(nlist->il_code),nlist->maxnrj);
560 srenew(nlist->jjnr ,nlist->maxnrj);
561 srenew(nlist->jjnr_end,nlist->maxnrj);
562 srenew(nlist->excl ,nlist->maxnrj*MAX_CGCGSIZE);
565 nlist->jjnr[nrj] = j_start;
566 nlist->jjnr_end[nrj] = j_end;
568 if (j_end - j_start > MAX_CGCGSIZE)
570 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);
573 /* Set the exclusions */
574 for(j=j_start; j<j_end; j++)
576 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
580 /* Avoid double counting of intra-cg interactions */
581 for(j=1; j<j_end-j_start; j++)
583 nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
591 put_in_list_t(gmx_bool bHaveVdW[],
607 put_in_list_at(gmx_bool bHaveVdW[],
622 /* The a[] index has been removed,
623 * to put it back in i_atom should be a[i0] and jj should be a[jj].
628 t_nblist * vdwc_free = NULL;
629 t_nblist * vdw_free = NULL;
630 t_nblist * coul_free = NULL;
631 t_nblist * vdwc_ww = NULL;
632 t_nblist * coul_ww = NULL;
634 int i,j,jcg,igid,gid,nbl_ind,ind_ij;
635 atom_id jj,jj0,jj1,i_atom;
640 real *charge,*chargeB;
642 gmx_bool bFreeEnergy,bFree,bFreeJ,bNotEx,*bPert;
643 gmx_bool bDoVdW_i,bDoCoul_i,bDoCoul_i_sol;
647 /* Copy some pointers */
649 charge = md->chargeA;
650 chargeB = md->chargeB;
653 bPert = md->bPerturbed;
657 nicg = index[icg+1]-i0;
659 /* Get the i charge group info */
660 igid = GET_CGINFO_GID(cginfo[icg]);
661 iwater = GET_CGINFO_SOLOPT(cginfo[icg]);
666 /* Check if any of the particles involved are perturbed.
667 * If not we can do the cheaper normal put_in_list
668 * and use more solvent optimization.
670 for(i=0; i<nicg; i++)
672 bFreeEnergy |= bPert[i0+i];
674 /* Loop over the j charge groups */
675 for(j=0; (j<nj && !bFreeEnergy); j++)
680 /* Finally loop over the atoms in the j-charge group */
681 for(jj=jj0; jj<jj1; jj++)
683 bFreeEnergy |= bPert[jj];
688 /* Unpack pointers to neighbourlist structs */
689 if (fr->nnblists == 1)
695 nbl_ind = fr->gid2nblists[GID(igid,jgid,ngid)];
699 nlist = fr->nblists[nbl_ind].nlist_lr;
703 nlist = fr->nblists[nbl_ind].nlist_sr;
706 if (iwater != esolNO)
708 vdwc = &nlist[eNL_VDWQQ_WATER];
709 vdw = &nlist[eNL_VDW];
710 coul = &nlist[eNL_QQ_WATER];
711 #ifndef DISABLE_WATERWATER_NLIST
712 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
713 coul_ww = &nlist[eNL_QQ_WATERWATER];
718 vdwc = &nlist[eNL_VDWQQ];
719 vdw = &nlist[eNL_VDW];
720 coul = &nlist[eNL_QQ];
725 if (iwater != esolNO)
727 /* Loop over the atoms in the i charge group */
729 gid = GID(igid,jgid,ngid);
730 /* Create new i_atom for each energy group */
731 if (bDoCoul && bDoVdW)
733 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
734 #ifndef DISABLE_WATERWATER_NLIST
735 new_i_nblist(vdwc_ww,bLR,i_atom,shift,gid);
740 new_i_nblist(vdw,bLR,i_atom,shift,gid);
744 new_i_nblist(coul,bLR,i_atom,shift,gid);
745 #ifndef DISABLE_WATERWATER_NLIST
746 new_i_nblist(coul_ww,bLR,i_atom,shift,gid);
749 /* Loop over the j charge groups */
750 for(j=0; (j<nj); j++)
760 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
762 if (iwater == esolSPC && jwater == esolSPC)
764 /* Interaction between two SPC molecules */
767 /* VdW only - only first atoms in each water interact */
768 add_j_to_nblist(vdw,jj0,bLR);
772 #ifdef DISABLE_WATERWATER_NLIST
773 /* Add entries for the three atoms - only do VdW if we need to */
776 add_j_to_nblist(coul,jj0,bLR);
780 add_j_to_nblist(vdwc,jj0,bLR);
782 add_j_to_nblist(coul,jj0+1,bLR);
783 add_j_to_nblist(coul,jj0+2,bLR);
785 /* One entry for the entire water-water interaction */
788 add_j_to_nblist(coul_ww,jj0,bLR);
792 add_j_to_nblist(vdwc_ww,jj0,bLR);
797 else if (iwater == esolTIP4P && jwater == esolTIP4P)
799 /* Interaction between two TIP4p molecules */
802 /* VdW only - only first atoms in each water interact */
803 add_j_to_nblist(vdw,jj0,bLR);
807 #ifdef DISABLE_WATERWATER_NLIST
808 /* Add entries for the four atoms - only do VdW if we need to */
811 add_j_to_nblist(vdw,jj0,bLR);
813 add_j_to_nblist(coul,jj0+1,bLR);
814 add_j_to_nblist(coul,jj0+2,bLR);
815 add_j_to_nblist(coul,jj0+3,bLR);
817 /* One entry for the entire water-water interaction */
820 add_j_to_nblist(coul_ww,jj0,bLR);
824 add_j_to_nblist(vdwc_ww,jj0,bLR);
831 /* j charge group is not water, but i is.
832 * Add entries to the water-other_atom lists; the geometry of the water
833 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
834 * so we don't care if it is SPC or TIP4P...
841 for(jj=jj0; (jj<jj1); jj++)
845 add_j_to_nblist(coul,jj,bLR);
851 for(jj=jj0; (jj<jj1); jj++)
853 if (bHaveVdW[type[jj]])
855 add_j_to_nblist(vdw,jj,bLR);
861 /* _charge_ _groups_ interact with both coulomb and LJ */
862 /* Check which atoms we should add to the lists! */
863 for(jj=jj0; (jj<jj1); jj++)
865 if (bHaveVdW[type[jj]])
869 add_j_to_nblist(vdwc,jj,bLR);
873 add_j_to_nblist(vdw,jj,bLR);
876 else if (charge[jj] != 0)
878 add_j_to_nblist(coul,jj,bLR);
885 close_i_nblist(coul);
886 close_i_nblist(vdwc);
887 #ifndef DISABLE_WATERWATER_NLIST
888 close_i_nblist(coul_ww);
889 close_i_nblist(vdwc_ww);
894 /* no solvent as i charge group */
895 /* Loop over the atoms in the i charge group */
896 for(i=0; i<nicg; i++)
899 gid = GID(igid,jgid,ngid);
902 /* Create new i_atom for each energy group */
903 if (bDoVdW && bDoCoul)
905 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
909 new_i_nblist(vdw,bLR,i_atom,shift,gid);
913 new_i_nblist(coul,bLR,i_atom,shift,gid);
915 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
916 bDoCoul_i = (bDoCoul && qi!=0);
918 if (bDoVdW_i || bDoCoul_i)
920 /* Loop over the j charge groups */
921 for(j=0; (j<nj); j++)
925 /* Check for large charge groups */
936 /* Finally loop over the atoms in the j-charge group */
937 for(jj=jj0; jj<jj1; jj++)
939 bNotEx = NOTEXCL(bExcl,i,jj);
947 add_j_to_nblist(coul,jj,bLR);
952 if (bHaveVdW[type[jj]])
954 add_j_to_nblist(vdw,jj,bLR);
959 if (bHaveVdW[type[jj]])
963 add_j_to_nblist(vdwc,jj,bLR);
967 add_j_to_nblist(vdw,jj,bLR);
970 else if (charge[jj] != 0)
972 add_j_to_nblist(coul,jj,bLR);
980 close_i_nblist(coul);
981 close_i_nblist(vdwc);
987 /* we are doing free energy */
988 vdwc_free = &nlist[eNL_VDWQQ_FREE];
989 vdw_free = &nlist[eNL_VDW_FREE];
990 coul_free = &nlist[eNL_QQ_FREE];
991 /* Loop over the atoms in the i charge group */
992 for(i=0; i<nicg; i++)
995 gid = GID(igid,jgid,ngid);
997 qiB = chargeB[i_atom];
999 /* Create new i_atom for each energy group */
1000 if (bDoVdW && bDoCoul)
1001 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
1003 new_i_nblist(vdw,bLR,i_atom,shift,gid);
1005 new_i_nblist(coul,bLR,i_atom,shift,gid);
1007 new_i_nblist(vdw_free,bLR,i_atom,shift,gid);
1008 new_i_nblist(coul_free,bLR,i_atom,shift,gid);
1009 new_i_nblist(vdwc_free,bLR,i_atom,shift,gid);
1011 bDoVdW_i = (bDoVdW &&
1012 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
1013 bDoCoul_i = (bDoCoul && (qi!=0 || qiB!=0));
1014 /* For TIP4P the first atom does not have a charge,
1015 * but the last three do. So we should still put an atom
1016 * without LJ but with charge in the water-atom neighborlist
1017 * for a TIP4p i charge group.
1018 * For SPC type water the first atom has LJ and charge,
1019 * so there is no such problem.
1021 if (iwater == esolNO)
1023 bDoCoul_i_sol = bDoCoul_i;
1027 bDoCoul_i_sol = bDoCoul;
1030 if (bDoVdW_i || bDoCoul_i_sol)
1032 /* Loop over the j charge groups */
1033 for(j=0; (j<nj); j++)
1037 /* Check for large charge groups */
1048 /* Finally loop over the atoms in the j-charge group */
1049 bFree = bPert[i_atom];
1050 for(jj=jj0; (jj<jj1); jj++)
1052 bFreeJ = bFree || bPert[jj];
1053 /* Complicated if, because the water H's should also
1054 * see perturbed j-particles
1056 if (iwater==esolNO || i==0 || bFreeJ)
1058 bNotEx = NOTEXCL(bExcl,i,jj);
1066 if (charge[jj]!=0 || chargeB[jj]!=0)
1068 add_j_to_nblist(coul_free,jj,bLR);
1071 else if (!bDoCoul_i)
1073 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1075 add_j_to_nblist(vdw_free,jj,bLR);
1080 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1082 if (charge[jj]!=0 || chargeB[jj]!=0)
1084 add_j_to_nblist(vdwc_free,jj,bLR);
1088 add_j_to_nblist(vdw_free,jj,bLR);
1091 else if (charge[jj]!=0 || chargeB[jj]!=0)
1092 add_j_to_nblist(coul_free,jj,bLR);
1097 /* This is done whether or not bWater is set */
1098 if (charge[jj] != 0)
1100 add_j_to_nblist(coul,jj,bLR);
1103 else if (!bDoCoul_i_sol)
1105 if (bHaveVdW[type[jj]])
1107 add_j_to_nblist(vdw,jj,bLR);
1112 if (bHaveVdW[type[jj]])
1114 if (charge[jj] != 0)
1116 add_j_to_nblist(vdwc,jj,bLR);
1120 add_j_to_nblist(vdw,jj,bLR);
1123 else if (charge[jj] != 0)
1125 add_j_to_nblist(coul,jj,bLR);
1133 close_i_nblist(vdw);
1134 close_i_nblist(coul);
1135 close_i_nblist(vdwc);
1136 close_i_nblist(vdw_free);
1137 close_i_nblist(coul_free);
1138 close_i_nblist(vdwc_free);
1144 put_in_list_qmmm(gmx_bool bHaveVdW[],
1160 int i,j,jcg,igid,gid;
1161 atom_id jj,jj0,jj1,i_atom;
1165 /* Get atom range */
1167 nicg = index[icg+1]-i0;
1169 /* Get the i charge group info */
1170 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1172 coul = &fr->QMMMlist;
1174 /* Loop over atoms in the ith charge group */
1175 for (i=0;i<nicg;i++)
1178 gid = GID(igid,jgid,ngid);
1179 /* Create new i_atom for each energy group */
1180 new_i_nblist(coul,bLR,i_atom,shift,gid);
1182 /* Loop over the j charge groups */
1187 /* Charge groups cannot have QM and MM atoms simultaneously */
1192 /* Finally loop over the atoms in the j-charge group */
1193 for(jj=jj0; jj<jj1; jj++)
1195 bNotEx = NOTEXCL(bExcl,i,jj);
1197 add_j_to_nblist(coul,jj,bLR);
1201 close_i_nblist(coul);
1206 put_in_list_cg(gmx_bool bHaveVdW[],
1222 int igid,gid,nbl_ind;
1226 cginfo = fr->cginfo[icg];
1228 igid = GET_CGINFO_GID(cginfo);
1229 gid = GID(igid,jgid,ngid);
1231 /* Unpack pointers to neighbourlist structs */
1232 if (fr->nnblists == 1)
1238 nbl_ind = fr->gid2nblists[gid];
1242 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1246 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1249 /* Make a new neighbor list for charge group icg.
1250 * Currently simply one neighbor list is made with LJ and Coulomb.
1251 * If required, zero interactions could be removed here
1252 * or in the force loop.
1254 new_i_nblist(vdwc,bLR,index[icg],shift,gid);
1255 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1257 for(j=0; (j<nj); j++)
1260 /* Skip the icg-icg pairs if all self interactions are excluded */
1261 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1263 /* Here we add the j charge group jcg to the list,
1264 * exclusions are also added to the list.
1266 add_j_to_nblist_cg(vdwc,index[jcg],index[jcg+1],bExcl,icg==jcg,bLR);
1270 close_i_nblist(vdwc);
1273 static void setexcl(atom_id start,atom_id end,t_blocka *excl,gmx_bool b,
1280 for(i=start; i<end; i++)
1282 for(k=excl->index[i]; k<excl->index[i+1]; k++)
1284 SETEXCL(bexcl,i-start,excl->a[k]);
1290 for(i=start; i<end; i++)
1292 for(k=excl->index[i]; k<excl->index[i+1]; k++)
1294 RMEXCL(bexcl,i-start,excl->a[k]);
1300 int calc_naaj(int icg,int cgtot)
1304 if ((cgtot % 2) == 1)
1306 /* Odd number of charge groups, easy */
1307 naaj = 1 + (cgtot/2);
1309 else if ((cgtot % 4) == 0)
1311 /* Multiple of four is hard */
1348 fprintf(log,"naaj=%d\n",naaj);
1354 /************************************************
1356 * S I M P L E C O R E S T U F F
1358 ************************************************/
1360 static real calc_image_tric(rvec xi,rvec xj,matrix box,
1361 rvec b_inv,int *shift)
1363 /* This code assumes that the cut-off is smaller than
1364 * a half times the smallest diagonal element of the box.
1371 /* Compute diff vector */
1372 dz = xj[ZZ] - xi[ZZ];
1373 dy = xj[YY] - xi[YY];
1374 dx = xj[XX] - xi[XX];
1376 /* Perform NINT operation, using trunc operation, therefore
1377 * we first add 2.5 then subtract 2 again
1379 tz = dz*b_inv[ZZ] + h25;
1381 dz -= tz*box[ZZ][ZZ];
1382 dy -= tz*box[ZZ][YY];
1383 dx -= tz*box[ZZ][XX];
1385 ty = dy*b_inv[YY] + h25;
1387 dy -= ty*box[YY][YY];
1388 dx -= ty*box[YY][XX];
1390 tx = dx*b_inv[XX]+h25;
1392 dx -= tx*box[XX][XX];
1394 /* Distance squared */
1395 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1397 *shift = XYZ2IS(tx,ty,tz);
1402 static real calc_image_rect(rvec xi,rvec xj,rvec box_size,
1403 rvec b_inv,int *shift)
1411 /* Compute diff vector */
1412 dx = xj[XX] - xi[XX];
1413 dy = xj[YY] - xi[YY];
1414 dz = xj[ZZ] - xi[ZZ];
1416 /* Perform NINT operation, using trunc operation, therefore
1417 * we first add 1.5 then subtract 1 again
1419 tx = dx*b_inv[XX] + h15;
1420 ty = dy*b_inv[YY] + h15;
1421 tz = dz*b_inv[ZZ] + h15;
1426 /* Correct diff vector for translation */
1427 ddx = tx*box_size[XX] - dx;
1428 ddy = ty*box_size[YY] - dy;
1429 ddz = tz*box_size[ZZ] - dz;
1431 /* Distance squared */
1432 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1434 *shift = XYZ2IS(tx,ty,tz);
1439 static void add_simple(t_ns_buf *nsbuf,int nrj,atom_id cg_j,
1440 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1441 int icg,int jgid,t_block *cgs,t_excl bexcl[],
1442 int shift,t_forcerec *fr,put_in_list_t *put_in_list)
1444 if (nsbuf->nj + nrj > MAX_CG)
1446 put_in_list(bHaveVdW,ngid,md,icg,jgid,nsbuf->ncg,nsbuf->jcg,
1447 cgs->index,bexcl,shift,fr,FALSE,TRUE,TRUE);
1448 /* Reset buffer contents */
1449 nsbuf->ncg = nsbuf->nj = 0;
1451 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1455 static void ns_inner_tric(rvec x[],int icg,int *i_egp_flags,
1456 int njcg,atom_id jcg[],
1457 matrix box,rvec b_inv,real rcut2,
1458 t_block *cgs,t_ns_buf **ns_buf,
1459 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1460 t_excl bexcl[],t_forcerec *fr,
1461 put_in_list_t *put_in_list)
1465 int *cginfo=fr->cginfo;
1466 atom_id cg_j,*cgindex;
1469 cgindex = cgs->index;
1471 for(j=0; (j<njcg); j++)
1474 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1475 if (calc_image_tric(x[icg],x[cg_j],box,b_inv,&shift) < rcut2)
1477 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1478 if (!(i_egp_flags[jgid] & EGP_EXCL))
1480 add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1481 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1488 static void ns_inner_rect(rvec x[],int icg,int *i_egp_flags,
1489 int njcg,atom_id jcg[],
1490 gmx_bool bBox,rvec box_size,rvec b_inv,real rcut2,
1491 t_block *cgs,t_ns_buf **ns_buf,
1492 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1493 t_excl bexcl[],t_forcerec *fr,
1494 put_in_list_t *put_in_list)
1498 int *cginfo=fr->cginfo;
1499 atom_id cg_j,*cgindex;
1502 cgindex = cgs->index;
1506 for(j=0; (j<njcg); j++)
1509 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1510 if (calc_image_rect(x[icg],x[cg_j],box_size,b_inv,&shift) < rcut2)
1512 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1513 if (!(i_egp_flags[jgid] & EGP_EXCL))
1515 add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1516 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1524 for(j=0; (j<njcg); j++)
1527 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1528 if ((rcut2 == 0) || (distance2(x[icg],x[cg_j]) < rcut2)) {
1529 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1530 if (!(i_egp_flags[jgid] & EGP_EXCL))
1532 add_simple(&ns_buf[jgid][CENTRAL],nrj,cg_j,
1533 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,CENTRAL,fr,
1541 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1543 static int ns_simple_core(t_forcerec *fr,
1544 gmx_localtop_t *top,
1546 matrix box,rvec box_size,
1547 t_excl bexcl[],atom_id *aaj,
1548 int ngid,t_ns_buf **ns_buf,
1549 put_in_list_t *put_in_list,gmx_bool bHaveVdW[])
1553 int nsearch,icg,jcg,igid,i0,nri,nn;
1556 /* atom_id *i_atoms; */
1557 t_block *cgs=&(top->cgs);
1558 t_blocka *excl=&(top->excls);
1561 gmx_bool bBox,bTriclinic;
1564 rlist2 = sqr(fr->rlist);
1566 bBox = (fr->ePBC != epbcNONE);
1569 for(m=0; (m<DIM); m++)
1571 b_inv[m] = divide_err(1.0,box_size[m]);
1573 bTriclinic = TRICLINIC(box);
1580 cginfo = fr->cginfo;
1583 for (icg=fr->cg0; (icg<fr->hcg); icg++)
1586 i0 = cgs->index[icg];
1587 nri = cgs->index[icg+1]-i0;
1588 i_atoms = &(cgs->a[i0]);
1589 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1590 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1592 igid = GET_CGINFO_GID(cginfo[icg]);
1593 i_egp_flags = fr->egp_flags + ngid*igid;
1594 setexcl(cgs->index[icg],cgs->index[icg+1],excl,TRUE,bexcl);
1596 naaj=calc_naaj(icg,cgs->nr);
1599 ns_inner_tric(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1600 box,b_inv,rlist2,cgs,ns_buf,
1601 bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1605 ns_inner_rect(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1606 bBox,box_size,b_inv,rlist2,cgs,ns_buf,
1607 bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1611 for(nn=0; (nn<ngid); nn++)
1613 for(k=0; (k<SHIFTS); k++)
1615 nsbuf = &(ns_buf[nn][k]);
1618 put_in_list(bHaveVdW,ngid,md,icg,nn,nsbuf->ncg,nsbuf->jcg,
1619 cgs->index,bexcl,k,fr,FALSE,TRUE,TRUE);
1620 nsbuf->ncg=nsbuf->nj=0;
1624 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1625 setexcl(cgs->index[icg],cgs->index[icg+1],excl,FALSE,bexcl);
1627 close_neighbor_list(fr,FALSE,-1,-1,FALSE);
1632 /************************************************
1634 * N S 5 G R I D S T U F F
1636 ************************************************/
1638 static inline void get_dx(int Nx,real gridx,real rc2,int xgi,real x,
1639 int *dx0,int *dx1,real *dcx2)
1667 for(i=xgi0; i>=0; i--)
1669 dcx = (i+1)*gridx-x;
1676 for(i=xgi1; i<Nx; i++)
1689 static inline void get_dx_dd(int Nx,real gridx,real rc2,int xgi,real x,
1690 int ncpddc,int shift_min,int shift_max,
1691 int *g0,int *g1,real *dcx2)
1694 int g_min,g_max,shift_home;
1727 g_min = (shift_min == shift_home ? 0 : ncpddc);
1728 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1735 else if (shift_max < 0)
1750 /* Check one grid cell down */
1751 dcx = ((*g0 - 1) + 1)*gridx - x;
1763 /* Check one grid cell up */
1764 dcx = (*g1 + 1)*gridx - x;
1776 #define sqr(x) ((x)*(x))
1777 #define calc_dx2(XI,YI,ZI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
1778 #define calc_cyl_dx2(XI,YI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
1779 /****************************************************
1781 * F A S T N E I G H B O R S E A R C H I N G
1783 * Optimized neighboursearching routine using grid
1784 * at least 1x1x1, see GROMACS manual
1786 ****************************************************/
1788 static void do_longrange(t_commrec *cr,gmx_localtop_t *top,t_forcerec *fr,
1789 int ngid,t_mdatoms *md,int icg,
1791 atom_id lr[],t_excl bexcl[],int shift,
1792 rvec x[],rvec box_size,t_nrnb *nrnb,
1793 real lambda,real *dvdlambda,
1794 gmx_grppairener_t *grppener,
1795 gmx_bool bDoVdW,gmx_bool bDoCoul,
1796 gmx_bool bEvaluateNow,put_in_list_t *put_in_list,
1797 gmx_bool bHaveVdW[],
1798 gmx_bool bDoForces,rvec *f)
1803 for(n=0; n<fr->nnblists; n++)
1805 for(i=0; (i<eNL_NR); i++)
1807 nl = &fr->nblists[n].nlist_lr[i];
1808 if ((nl->nri > nl->maxnri-32) || bEvaluateNow)
1810 close_neighbor_list(fr,TRUE,n,i,FALSE);
1811 /* Evaluate the energies and forces */
1812 do_nonbonded(cr,fr,x,f,md,NULL,
1813 grppener->ener[fr->bBHAM ? egBHAMLR : egLJLR],
1814 grppener->ener[egCOULLR],
1815 grppener->ener[egGB],box_size,
1816 nrnb,lambda,dvdlambda,n,i,
1817 GMX_DONB_LR | GMX_DONB_FORCES);
1819 reset_neighbor_list(fr,TRUE,n,i);
1826 /* Put the long range particles in a list */
1827 /* do_longrange is never called for QMMM */
1828 put_in_list(bHaveVdW,ngid,md,icg,jgid,nlr,lr,top->cgs.index,
1829 bexcl,shift,fr,TRUE,bDoVdW,bDoCoul);
1833 static void get_cutoff2(t_forcerec *fr,gmx_bool bDoLongRange,
1834 real *rvdw2,real *rcoul2,
1835 real *rs2,real *rm2,real *rl2)
1837 *rs2 = sqr(fr->rlist);
1838 if (bDoLongRange && fr->bTwinRange)
1840 /* The VdW and elec. LR cut-off's could be different,
1841 * so we can not simply set them to rlistlong.
1843 if (EVDW_MIGHT_BE_ZERO_AT_CUTOFF(fr->vdwtype) &&
1844 fr->rvdw > fr->rlist)
1846 *rvdw2 = sqr(fr->rlistlong);
1850 *rvdw2 = sqr(fr->rvdw);
1852 if (EEL_MIGHT_BE_ZERO_AT_CUTOFF(fr->eeltype) &&
1853 fr->rcoulomb > fr->rlist)
1855 *rcoul2 = sqr(fr->rlistlong);
1859 *rcoul2 = sqr(fr->rcoulomb);
1864 /* Workaround for a gcc -O3 or -ffast-math problem */
1868 *rm2 = min(*rvdw2,*rcoul2);
1869 *rl2 = max(*rvdw2,*rcoul2);
1872 static void init_nsgrid_lists(t_forcerec *fr,int ngid,gmx_ns_t *ns)
1874 real rvdw2,rcoul2,rs2,rm2,rl2;
1877 get_cutoff2(fr,TRUE,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
1879 /* Short range buffers */
1880 snew(ns->nl_sr,ngid);
1883 snew(ns->nlr_ljc,ngid);
1884 snew(ns->nlr_one,ngid);
1888 /* Long range VdW and Coul buffers */
1889 snew(ns->nl_lr_ljc,ngid);
1893 /* Long range VdW or Coul only buffers */
1894 snew(ns->nl_lr_one,ngid);
1896 for(j=0; (j<ngid); j++) {
1897 snew(ns->nl_sr[j],MAX_CG);
1900 snew(ns->nl_lr_ljc[j],MAX_CG);
1904 snew(ns->nl_lr_one[j],MAX_CG);
1910 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
1915 static int nsgrid_core(FILE *log,t_commrec *cr,t_forcerec *fr,
1916 matrix box,rvec box_size,int ngid,
1917 gmx_localtop_t *top,
1918 t_grid *grid,rvec x[],
1919 t_excl bexcl[],gmx_bool *bExcludeAlleg,
1920 t_nrnb *nrnb,t_mdatoms *md,
1921 real lambda,real *dvdlambda,
1922 gmx_grppairener_t *grppener,
1923 put_in_list_t *put_in_list,
1924 gmx_bool bHaveVdW[],
1925 gmx_bool bDoLongRange,gmx_bool bDoForces,rvec *f,
1926 gmx_bool bMakeQMMMnblist)
1929 atom_id **nl_lr_ljc,**nl_lr_one,**nl_sr;
1930 int *nlr_ljc,*nlr_one,*nsr;
1931 gmx_domdec_t *dd=NULL;
1932 t_block *cgs=&(top->cgs);
1933 int *cginfo=fr->cginfo;
1934 /* atom_id *i_atoms,*cgsindex=cgs->index; */
1936 int cell_x,cell_y,cell_z;
1937 int d,tx,ty,tz,dx,dy,dz,cj;
1938 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
1939 int zsh_ty,zsh_tx,ysh_tx;
1941 int dx0,dx1,dy0,dy1,dz0,dz1;
1942 int Nx,Ny,Nz,shift=-1,j,nrj,nns,nn=-1;
1943 real gridx,gridy,gridz,grid_x,grid_y,grid_z;
1944 real *dcx2,*dcy2,*dcz2;
1946 int cg0,cg1,icg=-1,cgsnr,i0,igid,nri,naaj,max_jcg;
1947 int jcg0,jcg1,jjcg,cgj0,jgid;
1948 int *grida,*gridnra,*gridind;
1949 gmx_bool rvdw_lt_rcoul,rcoul_lt_rvdw;
1950 rvec xi,*cgcm,grid_offset;
1951 real r2,rs2,rvdw2,rcoul2,rm2,rl2,XI,YI,ZI,dcx,dcy,dcz,tmp1,tmp2;
1953 gmx_bool bDomDec,bTriclinicX,bTriclinicY;
1958 bDomDec = DOMAINDECOMP(cr);
1964 bTriclinicX = ((YY < grid->npbcdim &&
1965 (!bDomDec || dd->nc[YY]==1) && box[YY][XX] != 0) ||
1966 (ZZ < grid->npbcdim &&
1967 (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][XX] != 0));
1968 bTriclinicY = (ZZ < grid->npbcdim &&
1969 (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][YY] != 0);
1973 get_cutoff2(fr,bDoLongRange,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
1975 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
1976 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
1978 if (bMakeQMMMnblist)
1986 nl_lr_ljc = ns->nl_lr_ljc;
1987 nl_lr_one = ns->nl_lr_one;
1988 nlr_ljc = ns->nlr_ljc;
1989 nlr_one = ns->nlr_one;
1997 gridind = grid->index;
1998 gridnra = grid->nra;
2001 gridx = grid->cell_size[XX];
2002 gridy = grid->cell_size[YY];
2003 gridz = grid->cell_size[ZZ];
2007 copy_rvec(grid->cell_offset,grid_offset);
2008 copy_ivec(grid->ncpddc,ncpddc);
2013 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2014 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
2015 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
2016 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
2017 if (zsh_tx!=0 && ysh_tx!=0)
2019 /* This could happen due to rounding, when both ratios are 0.5 */
2028 /* We only want a list for the test particle */
2037 /* Set the shift range */
2038 for(d=0; d<DIM; d++)
2042 /* Check if we need periodicity shifts.
2043 * Without PBC or with domain decomposition we don't need them.
2045 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2052 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
2063 /* Loop over charge groups */
2064 for(icg=cg0; (icg < cg1); icg++)
2066 igid = GET_CGINFO_GID(cginfo[icg]);
2067 /* Skip this charge group if all energy groups are excluded! */
2068 if (bExcludeAlleg[igid])
2073 i0 = cgs->index[icg];
2075 if (bMakeQMMMnblist)
2077 /* Skip this charge group if it is not a QM atom while making a
2078 * QM/MM neighbourlist
2080 if (md->bQM[i0]==FALSE)
2082 continue; /* MM particle, go to next particle */
2085 /* Compute the number of charge groups that fall within the control
2088 naaj = calc_naaj(icg,cgsnr);
2095 /* make a normal neighbourlist */
2099 /* Get the j charge-group and dd cell shift ranges */
2100 dd_get_ns_ranges(cr->dd,icg,&jcg0,&jcg1,sh0,sh1);
2105 /* Compute the number of charge groups that fall within the control
2108 naaj = calc_naaj(icg,cgsnr);
2114 /* The i-particle is awlways the test particle,
2115 * so we want all j-particles
2117 max_jcg = cgsnr - 1;
2121 max_jcg = jcg1 - cgsnr;
2126 i_egp_flags = fr->egp_flags + igid*ngid;
2128 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2129 setexcl(i0,cgs->index[icg+1],&top->excls,TRUE,bexcl);
2131 ci2xyz(grid,icg,&cell_x,&cell_y,&cell_z);
2133 /* Changed iicg to icg, DvdS 990115
2134 * (but see consistency check above, DvdS 990330)
2137 fprintf(log,"icg=%5d, naaj=%5d, cell %d %d %d\n",
2138 icg,naaj,cell_x,cell_y,cell_z);
2140 /* Loop over shift vectors in three dimensions */
2141 for (tz=-shp[ZZ]; tz<=shp[ZZ]; tz++)
2143 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2144 /* Calculate range of cells in Z direction that have the shift tz */
2145 zgi = cell_z + tz*Nz;
2148 get_dx(Nz,gridz,rl2,zgi,ZI,&dz0,&dz1,dcz2);
2150 get_dx_dd(Nz,gridz,rl2,zgi,ZI-grid_offset[ZZ],
2151 ncpddc[ZZ],sh0[ZZ],sh1[ZZ],&dz0,&dz1,dcz2);
2157 for (ty=-shp[YY]; ty<=shp[YY]; ty++)
2159 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2160 /* Calculate range of cells in Y direction that have the shift ty */
2163 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2167 ygi = cell_y + ty*Ny;
2170 get_dx(Ny,gridy,rl2,ygi,YI,&dy0,&dy1,dcy2);
2172 get_dx_dd(Ny,gridy,rl2,ygi,YI-grid_offset[YY],
2173 ncpddc[YY],sh0[YY],sh1[YY],&dy0,&dy1,dcy2);
2179 for (tx=-shp[XX]; tx<=shp[XX]; tx++)
2181 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2182 /* Calculate range of cells in X direction that have the shift tx */
2185 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2189 xgi = cell_x + tx*Nx;
2192 get_dx(Nx,gridx,rl2,xgi*Nx,XI,&dx0,&dx1,dcx2);
2194 get_dx_dd(Nx,gridx,rl2,xgi,XI-grid_offset[XX],
2195 ncpddc[XX],sh0[XX],sh1[XX],&dx0,&dx1,dcx2);
2201 /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2202 * from the neigbour list as it will not interact */
2203 if (fr->adress_type != eAdressOff){
2204 if (md->wf[cgs->index[icg]]==0 && egp_explicit(fr, igid)){
2208 /* Get shift vector */
2209 shift=XYZ2IS(tx,ty,tz);
2211 range_check(shift,0,SHIFTS);
2213 for(nn=0; (nn<ngid); nn++)
2220 fprintf(log,"shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2221 shift,dx0,dx1,dy0,dy1,dz0,dz1);
2222 fprintf(log,"cgcm: %8.3f %8.3f %8.3f\n",cgcm[icg][XX],
2223 cgcm[icg][YY],cgcm[icg][ZZ]);
2224 fprintf(log,"xi: %8.3f %8.3f %8.3f\n",XI,YI,ZI);
2226 for (dx=dx0; (dx<=dx1); dx++)
2228 tmp1 = rl2 - dcx2[dx];
2229 for (dy=dy0; (dy<=dy1); dy++)
2231 tmp2 = tmp1 - dcy2[dy];
2234 for (dz=dz0; (dz<=dz1); dz++) {
2235 if (tmp2 > dcz2[dz]) {
2236 /* Find grid-cell cj in which possible neighbours are */
2237 cj = xyz2ci(Ny,Nz,dx,dy,dz);
2239 /* Check out how many cgs (nrj) there in this cell */
2242 /* Find the offset in the cg list */
2245 /* Check if all j's are out of range so we
2246 * can skip the whole cell.
2247 * Should save some time, especially with DD.
2250 (grida[cgj0] >= max_jcg &&
2251 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2257 for (j=0; (j<nrj); j++)
2259 jjcg = grida[cgj0+j];
2261 /* check whether this guy is in range! */
2262 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2265 r2=calc_dx2(XI,YI,ZI,cgcm[jjcg]);
2267 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2268 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2269 /* check energy group exclusions */
2270 if (!(i_egp_flags[jgid] & EGP_EXCL))
2274 if (nsr[jgid] >= MAX_CG)
2276 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2277 nsr[jgid],nl_sr[jgid],
2278 cgs->index,/* cgsatoms, */ bexcl,
2279 shift,fr,FALSE,TRUE,TRUE);
2282 nl_sr[jgid][nsr[jgid]++]=jjcg;
2286 if (nlr_ljc[jgid] >= MAX_CG)
2288 do_longrange(cr,top,fr,ngid,md,icg,jgid,
2290 nl_lr_ljc[jgid],bexcl,shift,x,
2300 nl_lr_ljc[jgid][nlr_ljc[jgid]++]=jjcg;
2304 if (nlr_one[jgid] >= MAX_CG) {
2305 do_longrange(cr,top,fr,ngid,md,icg,jgid,
2307 nl_lr_one[jgid],bexcl,shift,x,
2311 rvdw_lt_rcoul,rcoul_lt_rvdw,FALSE,
2317 nl_lr_one[jgid][nlr_one[jgid]++]=jjcg;
2329 /* CHECK whether there is anything left in the buffers */
2330 for(nn=0; (nn<ngid); nn++)
2334 put_in_list(bHaveVdW,ngid,md,icg,nn,nsr[nn],nl_sr[nn],
2335 cgs->index, /* cgsatoms, */ bexcl,
2336 shift,fr,FALSE,TRUE,TRUE);
2339 if (nlr_ljc[nn] > 0)
2341 do_longrange(cr,top,fr,ngid,md,icg,nn,nlr_ljc[nn],
2342 nl_lr_ljc[nn],bexcl,shift,x,box_size,nrnb,
2343 lambda,dvdlambda,grppener,TRUE,TRUE,FALSE,
2344 put_in_list,bHaveVdW,bDoForces,f);
2347 if (nlr_one[nn] > 0)
2349 do_longrange(cr,top,fr,ngid,md,icg,nn,nlr_one[nn],
2350 nl_lr_one[nn],bexcl,shift,x,box_size,nrnb,
2351 lambda,dvdlambda,grppener,
2352 rvdw_lt_rcoul,rcoul_lt_rvdw,FALSE,
2353 put_in_list,bHaveVdW,bDoForces,f);
2359 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2360 setexcl(cgs->index[icg],cgs->index[icg+1],&top->excls,FALSE,bexcl);
2362 /* Perform any left over force calculations */
2363 for (nn=0; (nn<ngid); nn++)
2367 do_longrange(cr,top,fr,0,md,icg,nn,nlr_ljc[nn],
2368 nl_lr_ljc[nn],bexcl,shift,x,box_size,nrnb,
2369 lambda,dvdlambda,grppener,
2370 TRUE,TRUE,TRUE,put_in_list,bHaveVdW,bDoForces,f);
2373 do_longrange(cr,top,fr,0,md,icg,nn,nlr_one[nn],
2374 nl_lr_one[nn],bexcl,shift,x,box_size,nrnb,
2375 lambda,dvdlambda,grppener,
2376 rvdw_lt_rcoul,rcoul_lt_rvdw,
2377 TRUE,put_in_list,bHaveVdW,bDoForces,f);
2382 /* Close off short range neighbourlists */
2383 close_neighbor_list(fr,FALSE,-1,-1,bMakeQMMMnblist);
2388 void ns_realloc_natoms(gmx_ns_t *ns,int natoms)
2392 if (natoms > ns->nra_alloc)
2394 ns->nra_alloc = over_alloc_dd(natoms);
2395 srenew(ns->bexcl,ns->nra_alloc);
2396 for(i=0; i<ns->nra_alloc; i++)
2403 void init_ns(FILE *fplog,const t_commrec *cr,
2404 gmx_ns_t *ns,t_forcerec *fr,
2405 const gmx_mtop_t *mtop,
2408 int mt,icg,nr_in_cg,maxcg,i,j,jcg,ngid,ncg;
2412 /* Compute largest charge groups size (# atoms) */
2414 for(mt=0; mt<mtop->nmoltype; mt++) {
2415 cgs = &mtop->moltype[mt].cgs;
2416 for (icg=0; (icg < cgs->nr); icg++)
2418 nr_in_cg=max(nr_in_cg,(int)(cgs->index[icg+1]-cgs->index[icg]));
2422 /* Verify whether largest charge group is <= max cg.
2423 * This is determined by the type of the local exclusion type
2424 * Exclusions are stored in bits. (If the type is not large
2425 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2427 maxcg = sizeof(t_excl)*8;
2428 if (nr_in_cg > maxcg)
2430 gmx_fatal(FARGS,"Max #atoms in a charge group: %d > %d\n",
2434 ngid = mtop->groups.grps[egcENER].nr;
2435 snew(ns->bExcludeAlleg,ngid);
2436 for(i=0; i<ngid; i++) {
2437 ns->bExcludeAlleg[i] = TRUE;
2438 for(j=0; j<ngid; j++)
2440 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2442 ns->bExcludeAlleg[i] = FALSE;
2449 ns->grid = init_grid(fplog,fr);
2450 init_nsgrid_lists(fr,ngid,ns);
2455 snew(ns->ns_buf,ngid);
2456 for(i=0; (i<ngid); i++)
2458 snew(ns->ns_buf[i],SHIFTS);
2460 ncg = ncg_mtop(mtop);
2461 snew(ns->simple_aaj,2*ncg);
2462 for(jcg=0; (jcg<ncg); jcg++)
2464 ns->simple_aaj[jcg] = jcg;
2465 ns->simple_aaj[jcg+ncg] = jcg;
2469 /* Create array that determines whether or not atoms have VdW */
2470 snew(ns->bHaveVdW,fr->ntype);
2471 for(i=0; (i<fr->ntype); i++)
2473 for(j=0; (j<fr->ntype); j++)
2475 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2477 ((BHAMA(fr->nbfp,fr->ntype,i,j) != 0) ||
2478 (BHAMB(fr->nbfp,fr->ntype,i,j) != 0) ||
2479 (BHAMC(fr->nbfp,fr->ntype,i,j) != 0)) :
2480 ((C6(fr->nbfp,fr->ntype,i,j) != 0) ||
2481 (C12(fr->nbfp,fr->ntype,i,j) != 0))));
2485 pr_bvec(debug,0,"bHaveVdW",ns->bHaveVdW,fr->ntype,TRUE);
2489 if (!DOMAINDECOMP(cr))
2491 /* This could be reduced with particle decomposition */
2492 ns_realloc_natoms(ns,mtop->natoms);
2495 ns->nblist_initialized=FALSE;
2497 /* nbr list debug dump */
2499 char *ptr=getenv("GMX_DUMP_NL");
2502 ns->dump_nl=strtol(ptr,NULL,10);
2505 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2516 int search_neighbours(FILE *log,t_forcerec *fr,
2517 rvec x[],matrix box,
2518 gmx_localtop_t *top,
2519 gmx_groups_t *groups,
2521 t_nrnb *nrnb,t_mdatoms *md,
2522 real lambda,real *dvdlambda,
2523 gmx_grppairener_t *grppener,
2525 gmx_bool bDoLongRange,
2526 gmx_bool bDoForces,rvec *f)
2528 t_block *cgs=&(top->cgs);
2529 rvec box_size,grid_x0,grid_x1;
2531 real min_size,grid_dens;
2535 gmx_bool *i_egp_flags;
2536 int cg_start,cg_end,start,end;
2539 gmx_domdec_zones_t *dd_zones;
2540 put_in_list_t *put_in_list;
2544 /* Set some local variables */
2546 ngid = groups->grps[egcENER].nr;
2548 for(m=0; (m<DIM); m++)
2550 box_size[m] = box[m][m];
2553 if (fr->ePBC != epbcNONE)
2555 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC,box))
2557 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.");
2561 min_size = min(box_size[XX],min(box_size[YY],box_size[ZZ]));
2562 if (2*fr->rlistlong >= min_size)
2563 gmx_fatal(FARGS,"One of the box diagonal elements has become smaller than twice the cut-off length.");
2567 if (DOMAINDECOMP(cr))
2569 ns_realloc_natoms(ns,cgs->index[cgs->nr]);
2573 /* Reset the neighbourlists */
2574 reset_neighbor_list(fr,FALSE,-1,-1);
2576 if (bGrid && bFillGrid)
2580 if (DOMAINDECOMP(cr))
2582 dd_zones = domdec_zones(cr->dd);
2588 get_nsgrid_boundaries(grid,NULL,box,NULL,NULL,NULL,
2589 cgs->nr,fr->cg_cm,grid_x0,grid_x1,&grid_dens);
2591 grid_first(log,grid,NULL,NULL,fr->ePBC,box,grid_x0,grid_x1,
2592 fr->rlistlong,grid_dens);
2596 /* Don't know why this all is... (DvdS 3/99) */
2602 end = (cgs->nr+1)/2;
2605 if (DOMAINDECOMP(cr))
2608 fill_grid(log,dd_zones,grid,end,-1,end,fr->cg_cm);
2610 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2614 fill_grid(log,NULL,grid,cgs->nr,fr->cg0,fr->hcg,fr->cg_cm);
2615 grid->icg0 = fr->cg0;
2616 grid->icg1 = fr->hcg;
2624 calc_elemnr(log,grid,start,end,cgs->nr);
2626 grid_last(log,grid,start,end,cgs->nr);
2630 check_grid(debug,grid);
2631 print_grid(debug,grid);
2636 /* Set the grid cell index for the test particle only.
2637 * The cell to cg index is not corrected, but that does not matter.
2639 fill_grid(log,NULL,ns->grid,fr->hcg,fr->hcg-1,fr->hcg,fr->cg_cm);
2643 if (!fr->ns.bCGlist)
2645 put_in_list = put_in_list_at;
2649 put_in_list = put_in_list_cg;
2656 nsearch = nsgrid_core(log,cr,fr,box,box_size,ngid,top,
2657 grid,x,ns->bexcl,ns->bExcludeAlleg,
2658 nrnb,md,lambda,dvdlambda,grppener,
2659 put_in_list,ns->bHaveVdW,
2660 bDoLongRange,bDoForces,f,
2663 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2664 * the classical calculation. The charge-charge interaction
2665 * between QM and MM atoms is handled in the QMMM core calculation
2666 * (see QMMM.c). The VDW however, we'd like to compute classically
2667 * and the QM MM atom pairs have just been put in the
2668 * corresponding neighbourlists. in case of QMMM we still need to
2669 * fill a special QMMM neighbourlist that contains all neighbours
2670 * of the QM atoms. If bQMMM is true, this list will now be made:
2672 if (fr->bQMMM && fr->qr->QMMMscheme!=eQMMMschemeoniom)
2674 nsearch += nsgrid_core(log,cr,fr,box,box_size,ngid,top,
2675 grid,x,ns->bexcl,ns->bExcludeAlleg,
2676 nrnb,md,lambda,dvdlambda,grppener,
2677 put_in_list_qmmm,ns->bHaveVdW,
2678 bDoLongRange,bDoForces,f,
2684 nsearch = ns_simple_core(fr,top,md,box,box_size,
2685 ns->bexcl,ns->simple_aaj,
2686 ngid,ns->ns_buf,put_in_list,ns->bHaveVdW);
2694 inc_nrnb(nrnb,eNR_NS,nsearch);
2695 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2700 int natoms_beyond_ns_buffer(t_inputrec *ir,t_forcerec *fr,t_block *cgs,
2701 matrix scale_tot,rvec *x)
2703 int cg0,cg1,cg,a0,a1,a,i,j;
2704 real rint,hbuf2,scale;
2706 gmx_bool bIsotropic;
2711 rint = max(ir->rcoulomb,ir->rvdw);
2712 if (ir->rlist < rint)
2714 gmx_fatal(FARGS,"The neighbor search buffer has negative size: %f nm",
2722 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2724 hbuf2 = sqr(0.5*(ir->rlist - rint));
2725 for(cg=cg0; cg<cg1; cg++)
2727 a0 = cgs->index[cg];
2728 a1 = cgs->index[cg+1];
2729 for(a=a0; a<a1; a++)
2731 if (distance2(cg_cm[cg],x[a]) > hbuf2)
2741 scale = scale_tot[0][0];
2742 for(i=1; i<DIM; i++)
2744 /* With anisotropic scaling, the original spherical ns volumes become
2745 * ellipsoids. To avoid costly transformations we use the minimum
2746 * eigenvalue of the scaling matrix for determining the buffer size.
2747 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2749 scale = min(scale,scale_tot[i][i]);
2750 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2756 if (scale_tot[i][j] != 0)
2762 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2765 for(cg=cg0; cg<cg1; cg++)
2767 svmul(scale,cg_cm[cg],cgsc);
2768 a0 = cgs->index[cg];
2769 a1 = cgs->index[cg+1];
2770 for(a=a0; a<a1; a++)
2772 if (distance2(cgsc,x[a]) > hbuf2)
2781 /* Anistropic scaling */
2782 for(cg=cg0; cg<cg1; cg++)
2784 /* Since scale_tot contains the transpose of the scaling matrix,
2785 * we need to multiply with the transpose.
2787 tmvmul_ur0(scale_tot,cg_cm[cg],cgsc);
2788 a0 = cgs->index[cg];
2789 a1 = cgs->index[cg+1];
2790 for(a=a0; a<a1; a++)
2792 if (distance2(cgsc,x[a]) > hbuf2)